Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QMDReaction.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// 080505 Fixed and changed sampling method of impact parameter by T. Koi
27// 080602 Fix memory leaks by T. Koi
28// 080612 Delete unnecessary dependency and unused functions
29// Change criterion of reaction by T. Koi
30// 081107 Add UnUseGEM (then use the default channel of G4Evaporation)
31// UseFrag (chage criterion of a inelastic reaction)
32// Fix bug in nucleon projectiles by T. Koi
33// 090122 Be8 -> Alpha + Alpha
34// 090331 Change member shenXS and genspaXS object to pointer
35// 091119 Fix for incidence of neutral particles
36//
37#include "G4QMDReaction.hh"
38#include "G4QMDNucleus.hh"
40#include "G4Pow.hh"
42#include "G4SystemOfUnits.hh"
43#include "G4NistManager.hh"
44
46#include "G4BGGPionElasticXS.hh"
52
53
55: G4HadronicInteraction("QMDModel")
56, system ( NULL )
57, deltaT ( 1 ) // in fsec (c=1)
58, maxTime ( 100 ) // will have maxTime-th time step
59, envelopF ( 1.05 ) // 10% for Peripheral reactions
60, gem ( true )
61, frag ( false )
62, secID( -1 )
63{
65 pipElNucXS = new G4BGGPionElasticXS(G4PionPlus::PionPlus() );
66 pipElNucXS->BuildPhysicsTable(*(G4PionPlus::PionPlus() ) );
67
68 pimElNucXS = new G4BGGPionElasticXS(G4PionMinus::PionMinus() );
69 pimElNucXS->BuildPhysicsTable(*(G4PionMinus::PionMinus() ) );
70
71 pipInelNucXS = new G4BGGPionInelasticXS(G4PionPlus::PionPlus() );
72 pipInelNucXS->BuildPhysicsTable(*(G4PionPlus::PionPlus() ) );
73
74 pimInelNucXS = new G4BGGPionInelasticXS(G4PionMinus::PionMinus() );
75 pimInelNucXS->BuildPhysicsTable(*(G4PionMinus::PionMinus() ) );
76
77 meanField = new G4QMDMeanField();
78 collision = new G4QMDCollision();
79
80 excitationHandler = new G4ExcitationHandler();
81 setEvaporationCh();
82
83 coulomb_collision_gamma_proj = 0.0;
84 coulomb_collision_rx_proj = 0.0;
85 coulomb_collision_rz_proj = 0.0;
86 coulomb_collision_px_proj = 0.0;
87 coulomb_collision_pz_proj = 0.0;
88
89 coulomb_collision_gamma_targ = 0.0;
90 coulomb_collision_rx_targ = 0.0;
91 coulomb_collision_rz_targ = 0.0;
92 coulomb_collision_px_targ = 0.0;
93 coulomb_collision_pz_targ = 0.0;
94
95 secID = G4PhysicsModelCatalog::GetModelID( "model_QMDModel" );
96}
97
98
100{
101 delete excitationHandler;
102 delete collision;
103 delete meanField;
104}
105
106
108{
109 //G4cout << "G4QMDReaction::ApplyYourself" << G4endl;
110
111 theParticleChange.Clear();
112
113 system = new G4QMDSystem;
114
115 G4int proj_Z = 0;
116 G4int proj_A = 0;
117 const G4ParticleDefinition* proj_pd = ( const G4ParticleDefinition* ) projectile.GetDefinition();
118 if ( proj_pd->GetParticleType() == "nucleus" )
119 {
120 proj_Z = proj_pd->GetAtomicNumber();
121 proj_A = proj_pd->GetAtomicMass();
122 }
123 else
124 {
125 proj_Z = (int)( proj_pd->GetPDGCharge()/eplus );
126 proj_A = 1;
127 }
128 //G4int targ_Z = int ( target.GetZ() + 0.5 );
129 //G4int targ_A = int ( target.GetN() + 0.5 );
130 //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this)
131 G4int targ_Z = target.GetZ_asInt();
132 G4int targ_A = target.GetA_asInt();
133 const G4ParticleDefinition* targ_pd = G4IonTable::GetIonTable()->GetIon( targ_Z , targ_A , 0.0 );
134
135
136 //G4NistManager* nistMan = G4NistManager::Instance();
137// G4Element* G4NistManager::FindOrBuildElement( targ_Z );
138
139 const G4DynamicParticle* proj_dp = new G4DynamicParticle ( proj_pd , projectile.Get4Momentum() );
140 //const G4Element* targ_ele = nistMan->FindOrBuildElement( targ_Z );
141 //G4double aTemp = projectile.GetMaterial()->GetTemperature();
142
143 // Glauber-Gribov nucleus-nucleus cross section does not have GetIsoCrossSection,
144 // therefore call GetElementCrossSection instead.
145 //G4double xs_0 = theXS->GetIsoCrossSection ( proj_dp , targ_Z , targ_A );
146 G4double xs_0 = theXS->GetElementCrossSection( proj_dp , targ_Z , projectile.GetMaterial() );
147
148 // When the projectile is a pion
149 if (proj_pd == G4PionPlus::PionPlus() ) {
150 xs_0 = pipElNucXS->GetElementCrossSection(proj_dp, targ_Z, projectile.GetMaterial() ) +
151 pipInelNucXS->GetElementCrossSection(proj_dp, targ_Z, projectile.GetMaterial() );
152 } else if (proj_pd == G4PionMinus::PionMinus() ) {
153 xs_0 = pimElNucXS->GetElementCrossSection(proj_dp, targ_Z, projectile.GetMaterial() ) +
154 pimInelNucXS->GetElementCrossSection(proj_dp, targ_Z, projectile.GetMaterial() );
155 }
156
157 //G4double xs_0 = genspaXS->GetCrossSection ( proj_dp , targ_ele , aTemp );
158 //G4double xs_0 = theXS->GetCrossSection ( proj_dp , targ_ele , aTemp );
159 //110822
160
161 G4double bmax_0 = std::sqrt( xs_0 / pi );
162 //std::cout << "bmax_0 in fm (fermi) " << bmax_0/fermi << std::endl;
163
164 //delete proj_dp;
165
166 G4bool elastic = true;
167
168 std::vector< G4QMDNucleus* > nucleuses; // Secondary nuceluses
169 G4ThreeVector boostToReac; // ReactionSystem (CM or NN);
170 G4ThreeVector boostBackToLAB; // Reaction System to LAB;
171
172 G4LorentzVector targ4p( G4ThreeVector( 0.0 ) , targ_pd->GetPDGMass()/GeV );
173 G4ThreeVector boostLABtoCM = targ4p.findBoostToCM( proj_dp->Get4Momentum()/GeV ); // CM of target and proj;
174
175 G4double p1 = proj_dp->GetMomentum().mag()/GeV/proj_A;
176 G4double m1 = proj_dp->GetDefinition()->GetPDGMass()/GeV/proj_A;
177 G4double e1 = std::sqrt( p1*p1 + m1*m1 );
178 G4double e2 = targ_pd->GetPDGMass()/GeV/targ_A;
179 G4double beta_nn = -p1 / ( e1+e2 );
180
181 G4ThreeVector boostLABtoNN ( 0. , 0. , beta_nn ); // CM of NN;
182
183 G4double beta_nncm = ( - boostLABtoCM.beta() + boostLABtoNN.beta() ) / ( 1 - boostLABtoCM.beta() * boostLABtoNN.beta() ) ;
184
185 //std::cout << targ4p << std::endl;
186 //std::cout << proj_dp->Get4Momentum()<< std::endl;
187 //std::cout << beta_nncm << std::endl;
188 G4ThreeVector boostNNtoCM( 0. , 0. , beta_nncm ); //
189 G4ThreeVector boostCMtoNN( 0. , 0. , -beta_nncm ); //
190
191 boostToReac = boostLABtoNN;
192 boostBackToLAB = -boostLABtoNN;
193
194 delete proj_dp;
195
196 G4int icounter = 0;
197 G4int icounter_max = 1024;
198 while ( elastic ) // Loop checking, 11.03.2015, T. Koi
199 {
200 icounter++;
201 if ( icounter > icounter_max ) {
202 G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
203 break;
204 }
205
206// impact parameter
207 //G4double bmax = 1.05*(bmax_0/fermi); // 10% for Peripheral reactions
208 G4double bmax = envelopF*(bmax_0/fermi);
209 G4double b = bmax * std::sqrt ( G4UniformRand() );
210//071112
211 //G4double b = 0;
212 //G4double b = bmax;
213 //G4double b = bmax/1.05 * 0.7 * G4UniformRand();
214
215 //G4cout << "G4QMDRESULT bmax_0 = " << bmax_0/fermi << " fm, bmax = " << bmax << " fm , b = " << b << " fm " << G4endl;
216
217 G4double plab = projectile.GetTotalMomentum()/GeV;
218 G4double elab = ( projectile.GetKineticEnergy() + proj_pd->GetPDGMass() + targ_pd->GetPDGMass() )/GeV;
219
220 calcOffSetOfCollision( b , proj_pd , targ_pd , plab , elab , bmax , boostCMtoNN );
221
222// Projectile
223 G4LorentzVector proj4pLAB = projectile.Get4Momentum()/GeV;
224
225 G4QMDGroundStateNucleus* proj(NULL);
226 if ( projectile.GetDefinition()->GetParticleType() == "nucleus"
227 || projectile.GetDefinition()->GetParticleName() == "proton"
228 || projectile.GetDefinition()->GetParticleName() == "neutron" )
229 {
230
231 proj_Z = proj_pd->GetAtomicNumber();
232 proj_A = proj_pd->GetAtomicMass();
233
234 proj = new G4QMDGroundStateNucleus( proj_Z , proj_A );
235 //proj->ShowParticipants();
236
237
238 meanField->SetSystem ( proj );
239 proj->SetTotalPotential( meanField->GetTotalPotential() );
241
242 }
243
244// Target
245 //G4int iz = int ( target.GetZ() );
246 //G4int ia = int ( target.GetN() );
247 //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this)
248 G4int iz = int ( target.GetZ_asInt() );
249 G4int ia = int ( target.GetA_asInt() );
250
252
253 meanField->SetSystem (targ );
254 targ->SetTotalPotential( meanField->GetTotalPotential() );
256
257 //G4LorentzVector targ4p( G4ThreeVector( 0.0 ) , targ->GetNuclearMass()/GeV );
258// Boost Vector to CM
259 //boostToCM = targ4p.findBoostToCM( proj4pLAB );
260
261// Target
262 for ( G4int i = 0 ; i < targ->GetTotalNumberOfParticipant() ; i++ )
263 {
264
265 G4ThreeVector p0 = targ->GetParticipant( i )->GetMomentum();
266 G4ThreeVector r0 = targ->GetParticipant( i )->GetPosition();
267
268 G4ThreeVector p ( p0.x() + coulomb_collision_px_targ
269 , p0.y()
270 , p0.z() * coulomb_collision_gamma_targ + coulomb_collision_pz_targ );
271
272 G4ThreeVector r ( r0.x() + coulomb_collision_rx_targ
273 , r0.y()
274 , r0.z() / coulomb_collision_gamma_targ + coulomb_collision_rz_targ );
275
276 system->SetParticipant( new G4QMDParticipant( targ->GetParticipant( i )->GetDefinition() , p , r ) );
277 system->GetParticipant( i )->SetTarget();
278
279 }
280
281 G4LorentzVector proj4pCM = CLHEP::boostOf ( proj4pLAB , boostToReac );
282 G4LorentzVector targ4pCM = CLHEP::boostOf ( targ4p , boostToReac );
283
284// Projectile
285 //G4cout << "proj : " << proj << G4endl;
286 //if ( proj != NULL )
287 if ( proj_A != 1 )
288 {
289
290// projectile is nucleus
291
292 for ( G4int i = 0 ; i < proj->GetTotalNumberOfParticipant() ; i++ )
293 {
294
295 G4ThreeVector p0 = proj->GetParticipant( i )->GetMomentum();
296 G4ThreeVector r0 = proj->GetParticipant( i )->GetPosition();
297
298 G4ThreeVector p ( p0.x() + coulomb_collision_px_proj
299 , p0.y()
300 , p0.z() * coulomb_collision_gamma_proj + coulomb_collision_pz_proj );
301
302 G4ThreeVector r ( r0.x() + coulomb_collision_rx_proj
303 , r0.y()
304 , r0.z() / coulomb_collision_gamma_proj + coulomb_collision_rz_proj );
305
306 system->SetParticipant( new G4QMDParticipant( proj->GetParticipant( i )->GetDefinition() , p , r ) );
307 system->GetParticipant ( i + targ->GetTotalNumberOfParticipant() )->SetProjectile();
308 }
309
310 }
311 else
312 {
313
314// projectile is particle
315
316 // avoid multiple set in "elastic" loop
317 //G4cout << "system Total Participants : " << system->GetTotalNumberOfParticipant() << ", target : " << targ->GetTotalNumberOfParticipant() << G4endl;
318 if ( system->GetTotalNumberOfParticipant() == targ->GetTotalNumberOfParticipant() )
319 {
320
322
323 G4ThreeVector p0( 0 );
324 G4ThreeVector r0( 0 );
325
326 G4ThreeVector p ( p0.x() + coulomb_collision_px_proj
327 , p0.y()
328 , p0.z() * coulomb_collision_gamma_proj + coulomb_collision_pz_proj );
329
330 G4ThreeVector r ( r0.x() + coulomb_collision_rx_proj
331 , r0.y()
332 , r0.z() / coulomb_collision_gamma_proj + coulomb_collision_rz_proj );
333
334 system->SetParticipant( new G4QMDParticipant( (G4ParticleDefinition*)projectile.GetDefinition() , p , r ) );
335 // This is not important becase only 1 projectile particle.
336 system->GetParticipant ( i )->SetProjectile();
337 }
338
339 }
340 //system->ShowParticipants();
341
342 delete targ;
343 delete proj;
344
345 meanField->SetSystem ( system );
346 collision->SetMeanField ( meanField );
347
348// Time Evolution
349 //std::cout << "Start time evolution " << std::endl;
350 //system->ShowParticipants();
351 for ( G4int i = 0 ; i < maxTime ; i++ )
352 {
353 //G4cout << " do Paropagate " << i << " th time step. " << G4endl;
354 meanField->DoPropagation( deltaT );
355 //system->ShowParticipants();
356 collision->CalKinematicsOfBinaryCollisions( deltaT );
357
358 if ( i / 10 * 10 == i )
359 {
360 //G4cout << i << " th time step. " << G4endl;
361 //system->ShowParticipants();
362 }
363 //system->ShowParticipants();
364 }
365 //system->ShowParticipants();
366
367
368 //std::cout << "Doing Cluster Judgment " << std::endl;
369
370 nucleuses = meanField->DoClusterJudgment();
371
372// Elastic Judgment
373
374 G4int numberOfSecondary = int ( nucleuses.size() ) + system->GetTotalNumberOfParticipant();
375
376 G4int sec_a_Z = 0;
377 G4int sec_a_A = 0;
378 const G4ParticleDefinition* sec_a_pd = NULL;
379 G4int sec_b_Z = 0;
380 G4int sec_b_A = 0;
381 const G4ParticleDefinition* sec_b_pd = NULL;
382
383 if ( numberOfSecondary == 2 )
384 {
385
386 G4bool elasticLike_system = false;
387 if ( nucleuses.size() == 2 )
388 {
389
390 sec_a_Z = nucleuses[0]->GetAtomicNumber();
391 sec_a_A = nucleuses[0]->GetMassNumber();
392 sec_b_Z = nucleuses[1]->GetAtomicNumber();
393 sec_b_A = nucleuses[1]->GetMassNumber();
394
395 if ( ( sec_a_Z == proj_Z && sec_a_A == proj_A && sec_b_Z == targ_Z && sec_b_A == targ_A )
396 || ( sec_a_Z == targ_Z && sec_a_A == targ_A && sec_b_Z == proj_Z && sec_b_A == proj_A ) )
397 {
398 elasticLike_system = true;
399 }
400
401 }
402 else if ( nucleuses.size() == 1 )
403 {
404
405 sec_a_Z = nucleuses[0]->GetAtomicNumber();
406 sec_a_A = nucleuses[0]->GetMassNumber();
407 sec_b_pd = system->GetParticipant( 0 )->GetDefinition();
408
409 if ( ( sec_a_Z == proj_Z && sec_a_A == proj_A && sec_b_pd == targ_pd )
410 || ( sec_a_Z == targ_Z && sec_a_A == targ_A && sec_b_pd == proj_pd ) )
411 {
412 elasticLike_system = true;
413 }
414
415 }
416 else
417 {
418
419 sec_a_pd = system->GetParticipant( 0 )->GetDefinition();
420 sec_b_pd = system->GetParticipant( 1 )->GetDefinition();
421
422 if ( ( sec_a_pd == proj_pd && sec_b_pd == targ_pd )
423 || ( sec_a_pd == targ_pd && sec_b_pd == proj_pd ) )
424 {
425 elasticLike_system = true;
426 }
427
428 }
429
430 if ( elasticLike_system == true )
431 {
432
433 G4bool elasticLike_energy = true;
434// Cal ExcitationEnergy
435 for ( G4int i = 0 ; i < int ( nucleuses.size() ) ; i++ )
436 {
437
438 //meanField->SetSystem( nucleuses[i] );
439 meanField->SetNucleus( nucleuses[i] );
440 //nucleuses[i]->SetTotalPotential( meanField->GetTotalPotential() );
441 //nucleuses[i]->CalEnergyAndAngularMomentumInCM();
442
443 if ( nucleuses[i]->GetExcitationEnergy()*GeV > 1.0*MeV ) elasticLike_energy = false;
444
445 }
446
447// Check Collision
448 G4bool withCollision = true;
449 if ( system->GetNOCollision() == 0 ) withCollision = false;
450
451// Final judegement for Inelasitc or Elastic;
452//
453// ElasticLike without Collision
454 //if ( elasticLike_energy == true && withCollision == false ) elastic = true; // ielst = 0
455// ElasticLike with Collision
456 //if ( elasticLike_energy == true && withCollision == true ) elastic = true; // ielst = 1
457// InelasticLike without Collision
458 //if ( elasticLike_energy == false ) elastic = false; // ielst = 2
459 if ( frag == true )
460 if ( elasticLike_energy == false ) elastic = false;
461// InelasticLike with Collision
462 if ( elasticLike_energy == false && withCollision == true ) elastic = false; // ielst = 3
463
464 }
465
466 }
467 else
468 {
469
470// numberOfSecondary != 2
471 elastic = false;
472
473 }
474
475//071115
476 //G4cout << elastic << G4endl;
477 // if elastic is true try again from sampling of impact parameter
478
479 if ( elastic == true )
480 {
481 // delete this nucleues
482 for ( std::vector< G4QMDNucleus* >::iterator
483 it = nucleuses.begin() ; it != nucleuses.end() ; it++ )
484 {
485 delete *it;
486 }
487 nucleuses.clear();
488 }
489
490 }
491
492// random rotation
493 G4double randAzimuthAngle = CLHEP::twopi*G4UniformRand();
494 G4ThreeVector zAxis(0.,0.,1.);
495
496// Statical Decay Phase
497
498 for ( std::vector< G4QMDNucleus* >::iterator it
499 = nucleuses.begin() ; it != nucleuses.end() ; it++ )
500 {
501
502/*
503 G4cout << "G4QMDRESULT "
504 << (*it)->GetAtomicNumber()
505 << " "
506 << (*it)->GetMassNumber()
507 << " "
508 << (*it)->Get4Momentum()
509 << " "
510 << (*it)->Get4Momentum().vect()
511 << " "
512 << (*it)->Get4Momentum().restMass()
513 << " "
514 << (*it)->GetNuclearMass()/GeV
515 << G4endl;
516*/
517
518 meanField->SetNucleus ( *it );
519
520 if ( (*it)->GetAtomicNumber() == 0 // neutron cluster
521 || (*it)->GetAtomicNumber() == (*it)->GetMassNumber() ) // proton cluster
522 {
523 // push back system
524 for ( G4int i = 0 ; i < (*it)->GetTotalNumberOfParticipant() ; i++ )
525 {
526 G4QMDParticipant* aP = new G4QMDParticipant( ( (*it)->GetParticipant( i ) )->GetDefinition() , ( (*it)->GetParticipant( i ) )->GetMomentum() , ( (*it)->GetParticipant( i ) )->GetPosition() );
527 system->SetParticipant ( aP );
528 }
529 continue;
530 }
531
532 G4double nucleus_e = std::sqrt ( G4Pow::GetInstance()->powN ( (*it)->GetNuclearMass()/GeV , 2 ) + G4Pow::GetInstance()->powN ( (*it)->Get4Momentum().vect().mag() , 2 ) );
533 G4LorentzVector nucleus_p4CM ( (*it)->Get4Momentum().vect() , nucleus_e );
534
535// std::cout << "G4QMDRESULT nucleus deltaQ " << deltaQ << std::endl;
536
537 G4int ia = (*it)->GetMassNumber();
538 G4int iz = (*it)->GetAtomicNumber();
539
540 G4LorentzVector lv ( G4ThreeVector( 0.0 ) , (*it)->GetExcitationEnergy()*GeV + G4IonTable::GetIonTable()->GetIonMass( iz , ia ) );
541
542 G4Fragment* aFragment = new G4Fragment( ia , iz , lv );
543
545 rv = excitationHandler->BreakItUp( *aFragment );
546 G4bool notBreak = true;
547 for ( G4ReactionProductVector::iterator itt
548 = rv->begin() ; itt != rv->end() ; itt++ )
549 {
550
551 notBreak = false;
552 // Secondary from this nucleus (*it)
553 const G4ParticleDefinition* pd = (*itt)->GetDefinition();
554
555 G4LorentzVector p4 ( (*itt)->GetMomentum()/GeV , (*itt)->GetTotalEnergy()/GeV ); //in nucleus(*it) rest system
556 G4LorentzVector p4_CM = CLHEP::boostOf( p4 , -nucleus_p4CM.findBoostToCM() ); // Back to CM
557 G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB ); // Back to LAB
558 p4_LAB.rotate(randAzimuthAngle, zAxis); // random rotation
559
560
561//090122
562 //theParticleChange.AddSecondary( dp );
563 if ( !( pd->GetAtomicNumber() == 4 && pd->GetAtomicMass() == 8 ) )
564 {
565 //G4cout << "pd out of notBreak loop : " << pd->GetParticleName() << G4endl;
566 G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV );
567 theParticleChange.AddSecondary( dp );
568 }
569 else
570 {
571 //Be8 -> Alpha + Alpha + Q
572 G4ThreeVector randomized_direction( G4UniformRand() , G4UniformRand() , G4UniformRand() );
573 randomized_direction = randomized_direction.unit();
574 G4double q_decay = (*itt)->GetMass() - 2*G4Alpha::Alpha()->GetPDGMass();
575 G4double p_decay = std::sqrt ( G4Pow::GetInstance()->powN(G4Alpha::Alpha()->GetPDGMass()+q_decay/2,2) - G4Pow::GetInstance()->powN(G4Alpha::Alpha()->GetPDGMass() , 2 ) );
576 G4LorentzVector p4_a1 ( p_decay*randomized_direction , G4Alpha::Alpha()->GetPDGMass()+q_decay/2 ); //in Be8 rest system
577
578 G4LorentzVector p4_a1_Be8 = CLHEP::boostOf ( p4_a1/GeV , -p4.findBoostToCM() );
579 G4LorentzVector p4_a1_CM = CLHEP::boostOf ( p4_a1_Be8 , -nucleus_p4CM.findBoostToCM() );
580 G4LorentzVector p4_a1_LAB = CLHEP::boostOf ( p4_a1_CM , boostBackToLAB );
581 p4_a1_LAB.rotate(randAzimuthAngle, zAxis); // random rotation
582
583 G4LorentzVector p4_a2 ( -p_decay*randomized_direction , G4Alpha::Alpha()->GetPDGMass()+q_decay/2 ); //in Be8 rest system
584
585 G4LorentzVector p4_a2_Be8 = CLHEP::boostOf ( p4_a2/GeV , -p4.findBoostToCM() );
586 G4LorentzVector p4_a2_CM = CLHEP::boostOf ( p4_a2_Be8 , -nucleus_p4CM.findBoostToCM() );
587 G4LorentzVector p4_a2_LAB = CLHEP::boostOf ( p4_a2_CM , boostBackToLAB );
588 p4_a2_LAB.rotate(randAzimuthAngle, zAxis); // random rotation
589
590 G4DynamicParticle* dp1 = new G4DynamicParticle( G4Alpha::Alpha() , p4_a1_LAB*GeV );
591 G4DynamicParticle* dp2 = new G4DynamicParticle( G4Alpha::Alpha() , p4_a2_LAB*GeV );
592 theParticleChange.AddSecondary( dp1 );
593 theParticleChange.AddSecondary( dp2 );
594 }
595//090122
596
597/*
598 G4cout
599 << "Regist Secondary "
600 << (*itt)->GetDefinition()->GetParticleName()
601 << " "
602 << (*itt)->GetMomentum()/GeV
603 << " "
604 << (*itt)->GetKineticEnergy()/GeV
605 << " "
606 << (*itt)->GetMass()/GeV
607 << " "
608 << (*itt)->GetTotalEnergy()/GeV
609 << " "
610 << (*itt)->GetTotalEnergy()/GeV * (*itt)->GetTotalEnergy()/GeV
611 - (*itt)->GetMomentum()/GeV * (*itt)->GetMomentum()/GeV
612 << " "
613 << nucleus_p4CM.findBoostToCM()
614 << " "
615 << p4
616 << " "
617 << p4_CM
618 << " "
619 << p4_LAB
620 << G4endl;
621*/
622
623 }
624 if ( notBreak == true )
625 {
626
627 const G4ParticleDefinition* pd = G4IonTable::GetIonTable()->GetIon( (*it)->GetAtomicNumber() , (*it)->GetMassNumber(), (*it)->GetExcitationEnergy()*GeV );
628 //G4cout << "pd in notBreak loop : " << pd->GetParticleName() << G4endl;
629 G4LorentzVector p4_CM = nucleus_p4CM;
630 G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB ); // Back to LAB
631 p4_LAB.rotate(randAzimuthAngle, zAxis); // random rotation
632 G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV );
633 theParticleChange.AddSecondary( dp );
634
635 }
636
637 for ( G4ReactionProductVector::iterator itt
638 = rv->begin() ; itt != rv->end() ; itt++ )
639 {
640 delete *itt;
641 }
642 delete rv;
643
644 delete aFragment;
645
646 }
647
648
649
650 for ( G4int i = 0 ; i < system->GetTotalNumberOfParticipant() ; i++ )
651 {
652 // Secondary particles
653
654 const G4ParticleDefinition* pd = system->GetParticipant( i )->GetDefinition();
655 G4LorentzVector p4_CM = system->GetParticipant( i )->Get4Momentum();
656 G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB );
657 p4_LAB.rotate(randAzimuthAngle, zAxis); // random rotation
658 G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV );
659 theParticleChange.AddSecondary( dp );
660 //G4cout << "In the last theParticleChange loop : " << pd->GetParticleName() << G4endl;
661
662/*
663 G4cout << "G4QMDRESULT "
664 << "r" << i << " " << system->GetParticipant ( i ) -> GetPosition() << " "
665 << "p" << i << " " << system->GetParticipant ( i ) -> Get4Momentum()
666 << G4endl;
667*/
668
669 }
670
671 for ( std::vector< G4QMDNucleus* >::iterator it
672 = nucleuses.begin() ; it != nucleuses.end() ; it++ )
673 {
674 delete *it; // delete nulceuse
675 }
676 nucleuses.clear();
677
678 system->Clear();
679 delete system;
680
681 theParticleChange.SetStatusChange( stopAndKill );
682
683 for (G4int i = 0; i < G4int(theParticleChange.GetNumberOfSecondaries() ); i++)
684 {
685 //G4cout << "Particle : " << theParticleChange.GetSecondary(i)->GetParticle()->GetParticleDefinition()->GetParticleName() << G4endl;
686 //G4cout << "KEnergy : " << theParticleChange.GetSecondary(i)->GetParticle()->GetKineticEnergy() << G4endl;
687 //G4cout << "modelID : " << theParticleChange.GetSecondary(i)->GetCreatorModelID() << G4endl;
688 theParticleChange.GetSecondary(i)->SetCreatorModelID(secID);
689 }
690
691 return &theParticleChange;
692
693}
694
695
696
697void G4QMDReaction::calcOffSetOfCollision( G4double b ,
698const G4ParticleDefinition* pd_proj ,
699const G4ParticleDefinition* pd_targ ,
700G4double ptot , G4double etot , G4double bmax , G4ThreeVector boostToCM )
701{
702
703 G4double mass_proj = pd_proj->GetPDGMass()/GeV;
704 G4double mass_targ = pd_targ->GetPDGMass()/GeV;
705
706 G4double stot = std::sqrt ( etot*etot - ptot*ptot );
707
708 G4double pstt = std::sqrt ( ( stot*stot - ( mass_proj + mass_targ ) * ( mass_proj + mass_targ )
709 ) * ( stot*stot - ( mass_proj - mass_targ ) * ( mass_proj - mass_targ ) ) )
710 / ( 2.0 * stot );
711
712 G4double pzcc = pstt;
713 G4double eccm = stot - ( mass_proj + mass_targ );
714
715 G4int zp = 1;
716 G4int ap = 1;
717 if ( pd_proj->GetParticleType() == "nucleus" )
718 {
719 zp = pd_proj->GetAtomicNumber();
720 ap = pd_proj->GetAtomicMass();
721 }
722 else
723 {
724 // proton, neutron, mesons
725 zp = int ( pd_proj->GetPDGCharge()/eplus + 0.5 );
726 // ap = 1;
727 }
728
729
730 G4int zt = pd_targ->GetAtomicNumber();
731 G4int at = pd_targ->GetAtomicMass();
732
733
734 // Check the ramx0 value
735 //G4double rmax0 = 8.0; // T.K dicide parameter value // for low energy
736 G4double rmax0 = bmax + 4.0;
737 G4double rmax = std::sqrt( rmax0*rmax0 + b*b );
738
739 G4double ccoul = 0.001439767;
740 G4double pcca = 1.0 - double ( zp * zt ) * ccoul / eccm / rmax - ( b / rmax )*( b / rmax );
741
742 G4double pccf = std::sqrt( pcca );
743
744 //Fix for neutral particles
745 G4double aas1 = 0.0;
746 G4double bbs = 0.0;
747
748 if ( zp != 0 )
749 {
750 G4double aas = 2.0 * eccm * b / double ( zp * zt ) / ccoul;
751 bbs = 1.0 / std::sqrt ( 1.0 + aas*aas );
752 aas1 = ( 1.0 + aas * b / rmax ) * bbs;
753 }
754
755 G4double cost = 0.0;
756 G4double sint = 0.0;
757 G4double thet1 = 0.0;
758 G4double thet2 = 0.0;
759 if ( 1.0 - aas1*aas1 <= 0 || 1.0 - bbs*bbs <= 0.0 )
760 {
761 cost = 1.0;
762 sint = 0.0;
763 }
764 else
765 {
766 G4double aat1 = aas1 / std::sqrt ( 1.0 - aas1*aas1 );
767 G4double aat2 = bbs / std::sqrt ( 1.0 - bbs*bbs );
768
769 thet1 = std::atan ( aat1 );
770 thet2 = std::atan ( aat2 );
771
772// TK enter to else block
773 G4double theta = thet1 - thet2;
774 cost = std::cos( theta );
775 sint = std::sin( theta );
776 }
777
778 G4double rzpr = -rmax * cost * ( mass_targ ) / ( mass_proj + mass_targ );
779 G4double rzta = rmax * cost * ( mass_proj ) / ( mass_proj + mass_targ );
780
781 G4double rxpr = rmax / 2.0 * sint;
782
783 G4double rxta = -rxpr;
784
785
786 G4double pzpc = pzcc * ( cost * pccf + sint * b / rmax );
787 G4double pxpr = pzcc * ( -sint * pccf + cost * b / rmax );
788
789 G4double pztc = - pzpc;
790 G4double pxta = - pxpr;
791
792 G4double epc = std::sqrt ( pzpc*pzpc + pxpr*pxpr + mass_proj*mass_proj );
793 G4double etc = std::sqrt ( pztc*pztc + pxta*pxta + mass_targ*mass_targ );
794
795 G4double pzpr = pzpc;
796 G4double pzta = pztc;
797 G4double epr = epc;
798 G4double eta = etc;
799
800// CM -> NN
801 G4double gammacm = boostToCM.gamma();
802 //G4double betacm = -boostToCM.beta();
803 G4double betacm = boostToCM.z();
804 pzpr = pzpc + betacm * gammacm * ( gammacm / ( 1. + gammacm ) * pzpc * betacm + epc );
805 pzta = pztc + betacm * gammacm * ( gammacm / ( 1. + gammacm ) * pztc * betacm + etc );
806 epr = gammacm * ( epc + betacm * pzpc );
807 eta = gammacm * ( etc + betacm * pztc );
808
809 //G4double betpr = pzpr / epr;
810 //G4double betta = pzta / eta;
811
812 G4double gammpr = epr / ( mass_proj );
813 G4double gammta = eta / ( mass_targ );
814
815 pzta = pzta / double ( at );
816 pxta = pxta / double ( at );
817
818 pzpr = pzpr / double ( ap );
819 pxpr = pxpr / double ( ap );
820
821 G4double zeroz = 0.0;
822
823 rzpr = rzpr -zeroz;
824 rzta = rzta -zeroz;
825
826 // Set results
827 coulomb_collision_gamma_proj = gammpr;
828 coulomb_collision_rx_proj = rxpr;
829 coulomb_collision_rz_proj = rzpr;
830 coulomb_collision_px_proj = pxpr;
831 coulomb_collision_pz_proj = pzpr;
832
833 coulomb_collision_gamma_targ = gammta;
834 coulomb_collision_rx_targ = rxta;
835 coulomb_collision_rz_targ = rzta;
836 coulomb_collision_px_targ = pxta;
837 coulomb_collision_pz_targ = pzta;
838
839}
840
841void G4QMDReaction::setEvaporationCh()
842{
843 //fEvaporation - 8 default channels
844 //fCombined - 8 default + 60 GEM
845 //fGEM - 2 default + 66 GEM
846 G4DeexChannelType ctype = gem ? fGEM : fCombined;
847 excitationHandler->SetDeexChannelsType(ctype);
848}
849
850void G4QMDReaction::ModelDescription(std::ostream& outFile) const
851{
852 outFile << "Lorentz covarianted Quantum Molecular Dynamics model for nucleus (particle) vs nucleus reactions\n";
853}
@ stopAndKill
CLHEP::HepLorentzVector G4LorentzVector
std::vector< G4ReactionProduct * > G4ReactionProductVector
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
double beta() const
double z() const
Hep3Vector unit() const
double x() const
double y() const
double mag() const
double gamma() const
HepLorentzVector & rotate(double, const Hep3Vector &)
Hep3Vector findBoostToCM() const
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
G4ParticleDefinition * GetDefinition() const
G4LorentzVector Get4Momentum() const
G4ThreeVector GetMomentum() const
const G4Material * GetMaterial() const
G4double GetTotalMomentum() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4HadronicInteraction(const G4String &modelName="HadronicModel")
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
static G4IonTable * GetIonTable()
G4int GetA_asInt() const
Definition G4Nucleus.hh:78
G4int GetZ_asInt() const
Definition G4Nucleus.hh:84
G4int GetAtomicNumber() const
const G4String & GetParticleType() const
G4int GetAtomicMass() const
const G4String & GetParticleName() const
static G4int GetModelID(const G4int modelIndex)
static G4PionMinus * PionMinus()
static G4PionPlus * PionPlus()
Definition G4PionPlus.cc:93
static G4Pow * GetInstance()
Definition G4Pow.cc:41
void SetTotalPotential(G4double x)
void CalEnergyAndAngularMomentumInCM()
G4ThreeVector GetPosition()
const G4ParticleDefinition * GetDefinition()
G4ThreeVector GetMomentum()
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) override
void ModelDescription(std::ostream &outFile) const override
~G4QMDReaction() override
G4QMDParticipant * GetParticipant(G4int i)
G4int GetTotalNumberOfParticipant()
void SetSystem(G4QMDSystem *, G4ThreeVector, G4ThreeVector)
HepLorentzVector boostOf(const HepLorentzVector &vec, const Hep3Vector &betaVector)