Geant4  10.02
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"
41 #include "G4PhysicalConstants.hh"
42 #include "G4SystemOfUnits.hh"
43 #include "G4NistManager.hh"
44 
46 
48 : G4HadronicInteraction("QMDModel")
49 , system ( NULL )
50 , deltaT ( 1 ) // in fsec (c=1)
51 , maxTime ( 100 ) // will have maxTime-th time step
52 , envelopF ( 1.05 ) // 10% for Peripheral reactions
53 , gem ( true )
54 , frag ( false )
55 {
56 
57  //090331
59  //genspaXS = new G4GeneralSpaceNNCrossSection();
61  meanField = new G4QMDMeanField();
62  collision = new G4QMDCollision();
63 
68 
74 
80 
81 }
82 
83 
84 
86 {
87  delete evaporation;
88  delete excitationHandler;
89  delete collision;
90  delete meanField;
91 }
92 
93 
94 
96 {
97  //G4cout << "G4QMDReaction::ApplyYourself" << G4endl;
98 
100 
101  system = new G4QMDSystem;
102 
103  G4int proj_Z = 0;
104  G4int proj_A = 0;
105  const G4ParticleDefinition* proj_pd = ( const G4ParticleDefinition* ) projectile.GetDefinition();
106  if ( proj_pd->GetParticleType() == "nucleus" )
107  {
108  proj_Z = proj_pd->GetAtomicNumber();
109  proj_A = proj_pd->GetAtomicMass();
110  }
111  else
112  {
113  proj_Z = (int)( proj_pd->GetPDGCharge()/eplus );
114  proj_A = 1;
115  }
116  //G4int targ_Z = int ( target.GetZ() + 0.5 );
117  //G4int targ_A = int ( target.GetN() + 0.5 );
118  //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this)
119  G4int targ_Z = target.GetZ_asInt();
120  G4int targ_A = target.GetA_asInt();
121  const G4ParticleDefinition* targ_pd = G4IonTable::GetIonTable()->GetIon( targ_Z , targ_A , 0.0 );
122 
123 
124  //G4NistManager* nistMan = G4NistManager::Instance();
125 // G4Element* G4NistManager::FindOrBuildElement( targ_Z );
126 
127  const G4DynamicParticle* proj_dp = new G4DynamicParticle ( proj_pd , projectile.Get4Momentum() );
128  //const G4Element* targ_ele = nistMan->FindOrBuildElement( targ_Z );
129  //G4double aTemp = projectile.GetMaterial()->GetTemperature();
130 
131  //090331
132 
134 
135  if ( proj_pd->GetParticleType() == "meson" ) theXS = piNucXS;
136 
137  //G4double xs_0 = genspaXS->GetCrossSection ( proj_dp , targ_ele , aTemp );
138  //G4double xs_0 = theXS->GetCrossSection ( proj_dp , targ_ele , aTemp );
139  //110822
140  G4double xs_0 = theXS->GetIsoCrossSection ( proj_dp , targ_Z , targ_A );
141 
142  G4double bmax_0 = std::sqrt( xs_0 / pi );
143  //std::cout << "bmax_0 in fm (fermi) " << bmax_0/fermi << std::endl;
144 
145  //delete proj_dp;
146 
147  G4bool elastic = true;
148 
149  std::vector< G4QMDNucleus* > nucleuses; // Secondary nuceluses
150  G4ThreeVector boostToReac; // ReactionSystem (CM or NN);
151  G4ThreeVector boostBackToLAB; // Reaction System to LAB;
152 
153  G4LorentzVector targ4p( G4ThreeVector( 0.0 ) , targ_pd->GetPDGMass()/GeV );
154  G4ThreeVector boostLABtoCM = targ4p.findBoostToCM( proj_dp->Get4Momentum()/GeV ); // CM of target and proj;
155 
156  G4double p1 = proj_dp->GetMomentum().mag()/GeV/proj_A;
157  G4double m1 = proj_dp->GetDefinition()->GetPDGMass()/GeV/proj_A;
158  G4double e1 = std::sqrt( p1*p1 + m1*m1 );
159  G4double e2 = targ_pd->GetPDGMass()/GeV/targ_A;
160  G4double beta_nn = -p1 / ( e1+e2 );
161 
162  G4ThreeVector boostLABtoNN ( 0. , 0. , beta_nn ); // CM of NN;
163 
164  G4double beta_nncm = ( - boostLABtoCM.beta() + boostLABtoNN.beta() ) / ( 1 - boostLABtoCM.beta() * boostLABtoNN.beta() ) ;
165 
166  //std::cout << targ4p << std::endl;
167  //std::cout << proj_dp->Get4Momentum()<< std::endl;
168  //std::cout << beta_nncm << std::endl;
169  G4ThreeVector boostNNtoCM( 0. , 0. , beta_nncm ); //
170  G4ThreeVector boostCMtoNN( 0. , 0. , -beta_nncm ); //
171 
172  boostToReac = boostLABtoNN;
173  boostBackToLAB = -boostLABtoNN;
174 
175  delete proj_dp;
176 
177  G4int icounter = 0;
178  G4int icounter_max = 1024;
179  while ( elastic ) // Loop checking, 11.03.2015, T. Koi
180  {
181  icounter++;
182  if ( icounter > icounter_max ) {
183  G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
184  break;
185  }
186 
187 // impact parameter
188  //G4double bmax = 1.05*(bmax_0/fermi); // 10% for Peripheral reactions
189  G4double bmax = envelopF*(bmax_0/fermi);
190  G4double b = bmax * std::sqrt ( G4UniformRand() );
191 //071112
192  //G4double b = 0;
193  //G4double b = bmax;
194  //G4double b = bmax/1.05 * 0.7 * G4UniformRand();
195 
196  //G4cout << "G4QMDRESULT bmax_0 = " << bmax_0/fermi << " fm, bmax = " << bmax << " fm , b = " << b << " fm " << G4endl;
197 
198  G4double plab = projectile.GetTotalMomentum()/GeV;
199  G4double elab = ( projectile.GetKineticEnergy() + proj_pd->GetPDGMass() + targ_pd->GetPDGMass() )/GeV;
200 
201  calcOffSetOfCollision( b , proj_pd , targ_pd , plab , elab , bmax , boostCMtoNN );
202 
203 // Projectile
204  G4LorentzVector proj4pLAB = projectile.Get4Momentum()/GeV;
205 
206  G4QMDGroundStateNucleus* proj(NULL);
207  if ( projectile.GetDefinition()->GetParticleType() == "nucleus"
208  || projectile.GetDefinition()->GetParticleName() == "proton"
209  || projectile.GetDefinition()->GetParticleName() == "neutron" )
210  {
211 
212  proj_Z = proj_pd->GetAtomicNumber();
213  proj_A = proj_pd->GetAtomicMass();
214 
215  proj = new G4QMDGroundStateNucleus( proj_Z , proj_A );
216  //proj->ShowParticipants();
217 
218 
219  meanField->SetSystem ( proj );
222 
223  }
224 
225 // Target
226  //G4int iz = int ( target.GetZ() );
227  //G4int ia = int ( target.GetN() );
228  //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this)
229  G4int iz = int ( target.GetZ_asInt() );
230  G4int ia = int ( target.GetA_asInt() );
231 
232  G4QMDGroundStateNucleus* targ = new G4QMDGroundStateNucleus( iz , ia );
233 
234  meanField->SetSystem (targ );
235  targ->SetTotalPotential( meanField->GetTotalPotential() );
236  targ->CalEnergyAndAngularMomentumInCM();
237 
238  //G4LorentzVector targ4p( G4ThreeVector( 0.0 ) , targ->GetNuclearMass()/GeV );
239 // Boost Vector to CM
240  //boostToCM = targ4p.findBoostToCM( proj4pLAB );
241 
242 // Target
243  for ( G4int i = 0 ; i < targ->GetTotalNumberOfParticipant() ; i++ )
244  {
245 
246  G4ThreeVector p0 = targ->GetParticipant( i )->GetMomentum();
247  G4ThreeVector r0 = targ->GetParticipant( i )->GetPosition();
248 
250  , p0.y()
252 
254  , r0.y()
256 
257  system->SetParticipant( new G4QMDParticipant( targ->GetParticipant( i )->GetDefinition() , p , r ) );
259 
260  }
261 
262  G4LorentzVector proj4pCM = CLHEP::boostOf ( proj4pLAB , boostToReac );
263  G4LorentzVector targ4pCM = CLHEP::boostOf ( targ4p , boostToReac );
264 
265 // Projectile
266  if ( proj != NULL )
267  {
268 
269 // projectile is nucleus
270 
271  for ( G4int i = 0 ; i < proj->GetTotalNumberOfParticipant() ; i++ )
272  {
273 
274  G4ThreeVector p0 = proj->GetParticipant( i )->GetMomentum();
275  G4ThreeVector r0 = proj->GetParticipant( i )->GetPosition();
276 
278  , p0.y()
280 
282  , r0.y()
284 
285  system->SetParticipant( new G4QMDParticipant( proj->GetParticipant( i )->GetDefinition() , p , r ) );
286  system->GetParticipant ( i + targ->GetTotalNumberOfParticipant() )->SetProjectile();
287  }
288 
289  }
290  else
291  {
292 
293 // projectile is particle
294 
295  // avoid multiple set in "elastic" loop
296  if ( system->GetTotalNumberOfParticipant() == targ->GetTotalNumberOfParticipant() )
297  {
298 
299  G4int i = targ->GetTotalNumberOfParticipant();
300 
301  G4ThreeVector p0( 0 );
302  G4ThreeVector r0( 0 );
303 
305  , p0.y()
307 
309  , r0.y()
311 
312  system->SetParticipant( new G4QMDParticipant( (G4ParticleDefinition*)projectile.GetDefinition() , p , r ) );
313  // This is not important becase only 1 projectile particle.
315  }
316 
317  }
318  //system->ShowParticipants();
319 
320  delete targ;
321  delete proj;
322 
325 
326 // Time Evolution
327  //std::cout << "Start time evolution " << std::endl;
328  //system->ShowParticipants();
329  for ( G4int i = 0 ; i < maxTime ; i++ )
330  {
331  //G4cout << " do Paropagate " << i << " th time step. " << G4endl;
333  //system->ShowParticipants();
335 
336  if ( i / 10 * 10 == i )
337  {
338  //G4cout << i << " th time step. " << G4endl;
339  //system->ShowParticipants();
340  }
341  //system->ShowParticipants();
342  }
343  //system->ShowParticipants();
344 
345 
346  //std::cout << "Doing Cluster Judgment " << std::endl;
347 
348  nucleuses = meanField->DoClusterJudgment();
349 
350 // Elastic Judgment
351 
352  G4int numberOfSecondary = int ( nucleuses.size() ) + system->GetTotalNumberOfParticipant();
353 
354  G4int sec_a_Z = 0;
355  G4int sec_a_A = 0;
356  const G4ParticleDefinition* sec_a_pd = NULL;
357  G4int sec_b_Z = 0;
358  G4int sec_b_A = 0;
359  const G4ParticleDefinition* sec_b_pd = NULL;
360 
361  if ( numberOfSecondary == 2 )
362  {
363 
364  G4bool elasticLike_system = false;
365  if ( nucleuses.size() == 2 )
366  {
367 
368  sec_a_Z = nucleuses[0]->GetAtomicNumber();
369  sec_a_A = nucleuses[0]->GetMassNumber();
370  sec_b_Z = nucleuses[1]->GetAtomicNumber();
371  sec_b_A = nucleuses[1]->GetMassNumber();
372 
373  if ( ( sec_a_Z == proj_Z && sec_a_A == proj_A && sec_b_Z == targ_Z && sec_b_A == targ_A )
374  || ( sec_a_Z == targ_Z && sec_a_A == targ_A && sec_b_Z == proj_Z && sec_b_A == proj_A ) )
375  {
376  elasticLike_system = true;
377  }
378 
379  }
380  else if ( nucleuses.size() == 1 )
381  {
382 
383  sec_a_Z = nucleuses[0]->GetAtomicNumber();
384  sec_a_A = nucleuses[0]->GetMassNumber();
385  sec_b_pd = system->GetParticipant( 0 )->GetDefinition();
386 
387  if ( ( sec_a_Z == proj_Z && sec_a_A == proj_A && sec_b_pd == targ_pd )
388  || ( sec_a_Z == targ_Z && sec_a_A == targ_A && sec_b_pd == proj_pd ) )
389  {
390  elasticLike_system = true;
391  }
392 
393  }
394  else
395  {
396 
397  sec_a_pd = system->GetParticipant( 0 )->GetDefinition();
398  sec_b_pd = system->GetParticipant( 1 )->GetDefinition();
399 
400  if ( ( sec_a_pd == proj_pd && sec_b_pd == targ_pd )
401  || ( sec_a_pd == targ_pd && sec_b_pd == proj_pd ) )
402  {
403  elasticLike_system = true;
404  }
405 
406  }
407 
408  if ( elasticLike_system == true )
409  {
410 
411  G4bool elasticLike_energy = true;
412 // Cal ExcitationEnergy
413  for ( G4int i = 0 ; i < int ( nucleuses.size() ) ; i++ )
414  {
415 
416  //meanField->SetSystem( nucleuses[i] );
417  meanField->SetNucleus( nucleuses[i] );
418  //nucleuses[i]->SetTotalPotential( meanField->GetTotalPotential() );
419  //nucleuses[i]->CalEnergyAndAngularMomentumInCM();
420 
421  if ( nucleuses[i]->GetExcitationEnergy()*GeV > 1.0*MeV ) elasticLike_energy = false;
422 
423  }
424 
425 // Check Collision
426  G4bool withCollision = true;
427  if ( system->GetNOCollision() == 0 ) withCollision = false;
428 
429 // Final judegement for Inelasitc or Elastic;
430 //
431 // ElasticLike without Collision
432  //if ( elasticLike_energy == true && withCollision == false ) elastic = true; // ielst = 0
433 // ElasticLike with Collision
434  //if ( elasticLike_energy == true && withCollision == true ) elastic = true; // ielst = 1
435 // InelasticLike without Collision
436  //if ( elasticLike_energy == false ) elastic = false; // ielst = 2
437  if ( frag == true )
438  if ( elasticLike_energy == false ) elastic = false;
439 // InelasticLike with Collision
440  if ( elasticLike_energy == false && withCollision == true ) elastic = false; // ielst = 3
441 
442  }
443 
444  }
445  else
446  {
447 
448 // numberOfSecondary != 2
449  elastic = false;
450 
451  }
452 
453 //071115
454  //G4cout << elastic << G4endl;
455  // if elastic is true try again from sampling of impact parameter
456 
457  if ( elastic == true )
458  {
459  // delete this nucleues
460  for ( std::vector< G4QMDNucleus* >::iterator
461  it = nucleuses.begin() ; it != nucleuses.end() ; it++ )
462  {
463  delete *it;
464  }
465  nucleuses.clear();
466  }
467 
468  }
469 
470 
471 // Statical Decay Phase
472 
473  for ( std::vector< G4QMDNucleus* >::iterator it
474  = nucleuses.begin() ; it != nucleuses.end() ; it++ )
475  {
476 
477 /*
478  G4cout << "G4QMDRESULT "
479  << (*it)->GetAtomicNumber()
480  << " "
481  << (*it)->GetMassNumber()
482  << " "
483  << (*it)->Get4Momentum()
484  << " "
485  << (*it)->Get4Momentum().vect()
486  << " "
487  << (*it)->Get4Momentum().restMass()
488  << " "
489  << (*it)->GetNuclearMass()/GeV
490  << G4endl;
491 */
492 
493  meanField->SetNucleus ( *it );
494 
495  if ( (*it)->GetAtomicNumber() == 0 // neutron cluster
496  || (*it)->GetAtomicNumber() == (*it)->GetMassNumber() ) // proton cluster
497  {
498  // push back system
499  for ( G4int i = 0 ; i < (*it)->GetTotalNumberOfParticipant() ; i++ )
500  {
501  G4QMDParticipant* aP = new G4QMDParticipant( ( (*it)->GetParticipant( i ) )->GetDefinition() , ( (*it)->GetParticipant( i ) )->GetMomentum() , ( (*it)->GetParticipant( i ) )->GetPosition() );
502  system->SetParticipant ( aP );
503  }
504  continue;
505  }
506 
507  G4double nucleus_e = std::sqrt ( G4Pow::GetInstance()->powN ( (*it)->GetNuclearMass()/GeV , 2 ) + G4Pow::GetInstance()->powN ( (*it)->Get4Momentum().vect().mag() , 2 ) );
508  G4LorentzVector nucleus_p4CM ( (*it)->Get4Momentum().vect() , nucleus_e );
509 
510 // std::cout << "G4QMDRESULT nucleus deltaQ " << deltaQ << std::endl;
511 
512  G4int ia = (*it)->GetMassNumber();
513  G4int iz = (*it)->GetAtomicNumber();
514 
515  G4LorentzVector lv ( G4ThreeVector( 0.0 ) , (*it)->GetExcitationEnergy()*GeV + G4IonTable::GetIonTable()->GetIonMass( iz , ia ) );
516 
517  G4Fragment* aFragment = new G4Fragment( ia , iz , lv );
518 
520  rv = excitationHandler->BreakItUp( *aFragment );
521  G4bool notBreak = true;
522  for ( G4ReactionProductVector::iterator itt
523  = rv->begin() ; itt != rv->end() ; itt++ )
524  {
525 
526  notBreak = false;
527  // Secondary from this nucleus (*it)
528  const G4ParticleDefinition* pd = (*itt)->GetDefinition();
529 
530  G4LorentzVector p4 ( (*itt)->GetMomentum()/GeV , (*itt)->GetTotalEnergy()/GeV ); //in nucleus(*it) rest system
531  G4LorentzVector p4_CM = CLHEP::boostOf( p4 , -nucleus_p4CM.findBoostToCM() ); // Back to CM
532  G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB ); // Back to LAB
533 
534 
535 //090122
536  //theParticleChange.AddSecondary( dp );
537  if ( !( pd->GetAtomicNumber() == 4 && pd->GetAtomicMass() == 8 ) )
538  {
539  G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV );
541  }
542  else
543  {
544  //Be8 -> Alpha + Alpha + Q
545  G4ThreeVector randomized_direction( G4UniformRand() , G4UniformRand() , G4UniformRand() );
546  randomized_direction = randomized_direction.unit();
547  G4double q_decay = (*itt)->GetMass() - 2*G4Alpha::Alpha()->GetPDGMass();
548  G4double p_decay = std::sqrt ( G4Pow::GetInstance()->powN(G4Alpha::Alpha()->GetPDGMass()+q_decay/2,2) - G4Pow::GetInstance()->powN(G4Alpha::Alpha()->GetPDGMass() , 2 ) );
549  G4LorentzVector p4_a1 ( p_decay*randomized_direction , G4Alpha::Alpha()->GetPDGMass()+q_decay/2 ); //in Be8 rest system
550 
551  G4LorentzVector p4_a1_Be8 = CLHEP::boostOf ( p4_a1/GeV , -p4.findBoostToCM() );
552  G4LorentzVector p4_a1_CM = CLHEP::boostOf ( p4_a1_Be8 , -nucleus_p4CM.findBoostToCM() );
553  G4LorentzVector p4_a1_LAB = CLHEP::boostOf ( p4_a1_CM , boostBackToLAB );
554 
555  G4LorentzVector p4_a2 ( -p_decay*randomized_direction , G4Alpha::Alpha()->GetPDGMass()+q_decay/2 ); //in Be8 rest system
556 
557  G4LorentzVector p4_a2_Be8 = CLHEP::boostOf ( p4_a2/GeV , -p4.findBoostToCM() );
558  G4LorentzVector p4_a2_CM = CLHEP::boostOf ( p4_a2_Be8 , -nucleus_p4CM.findBoostToCM() );
559  G4LorentzVector p4_a2_LAB = CLHEP::boostOf ( p4_a2_CM , boostBackToLAB );
560 
561  G4DynamicParticle* dp1 = new G4DynamicParticle( G4Alpha::Alpha() , p4_a1_LAB*GeV );
562  G4DynamicParticle* dp2 = new G4DynamicParticle( G4Alpha::Alpha() , p4_a2_LAB*GeV );
565  }
566 //090122
567 
568 /*
569  G4cout
570  << "Regist Secondary "
571  << (*itt)->GetDefinition()->GetParticleName()
572  << " "
573  << (*itt)->GetMomentum()/GeV
574  << " "
575  << (*itt)->GetKineticEnergy()/GeV
576  << " "
577  << (*itt)->GetMass()/GeV
578  << " "
579  << (*itt)->GetTotalEnergy()/GeV
580  << " "
581  << (*itt)->GetTotalEnergy()/GeV * (*itt)->GetTotalEnergy()/GeV
582  - (*itt)->GetMomentum()/GeV * (*itt)->GetMomentum()/GeV
583  << " "
584  << nucleus_p4CM.findBoostToCM()
585  << " "
586  << p4
587  << " "
588  << p4_CM
589  << " "
590  << p4_LAB
591  << G4endl;
592 */
593 
594  }
595  if ( notBreak == true )
596  {
597 
598  const G4ParticleDefinition* pd = G4IonTable::GetIonTable()->GetIon( (*it)->GetAtomicNumber() , (*it)->GetMassNumber(), (*it)->GetExcitationEnergy()*GeV );
599  G4LorentzVector p4_CM = nucleus_p4CM;
600  G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB ); // Back to LAB
601  G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV );
603 
604  }
605 
606  for ( G4ReactionProductVector::iterator itt
607  = rv->begin() ; itt != rv->end() ; itt++ )
608  {
609  delete *itt;
610  }
611  delete rv;
612 
613  delete aFragment;
614 
615  }
616 
617 
618 
619  for ( G4int i = 0 ; i < system->GetTotalNumberOfParticipant() ; i++ )
620  {
621 
622  // Secondary particles
623 
626  G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB );
627  G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV );
629 
630 /*
631  G4cout << "G4QMDRESULT "
632  << "r" << i << " " << system->GetParticipant ( i ) -> GetPosition() << " "
633  << "p" << i << " " << system->GetParticipant ( i ) -> Get4Momentum()
634  << G4endl;
635 */
636 
637  }
638 
639  for ( std::vector< G4QMDNucleus* >::iterator it
640  = nucleuses.begin() ; it != nucleuses.end() ; it++ )
641  {
642  delete *it; // delete nulceuse
643  }
644  nucleuses.clear();
645 
646  system->Clear();
647  delete system;
648 
650 
651  return &theParticleChange;
652 
653 }
654 
655 
656 
658 const G4ParticleDefinition* pd_proj ,
659 const G4ParticleDefinition* pd_targ ,
660 G4double ptot , G4double etot , G4double bmax , G4ThreeVector boostToCM )
661 {
662 
663  G4double mass_proj = pd_proj->GetPDGMass()/GeV;
664  G4double mass_targ = pd_targ->GetPDGMass()/GeV;
665 
666  G4double stot = std::sqrt ( etot*etot - ptot*ptot );
667 
668  G4double pstt = std::sqrt ( ( stot*stot - ( mass_proj + mass_targ ) * ( mass_proj + mass_targ )
669  ) * ( stot*stot - ( mass_proj - mass_targ ) * ( mass_proj - mass_targ ) ) )
670  / ( 2.0 * stot );
671 
672  G4double pzcc = pstt;
673  G4double eccm = stot - ( mass_proj + mass_targ );
674 
675  G4int zp = 1;
676  G4int ap = 1;
677  if ( pd_proj->GetParticleType() == "nucleus" )
678  {
679  zp = pd_proj->GetAtomicNumber();
680  ap = pd_proj->GetAtomicMass();
681  }
682  else
683  {
684  // proton, neutron, mesons
685  zp = int ( pd_proj->GetPDGCharge()/eplus + 0.5 );
686  // ap = 1;
687  }
688 
689 
690  G4int zt = pd_targ->GetAtomicNumber();
691  G4int at = pd_targ->GetAtomicMass();
692 
693 
694  //G4double rmax0 = 8.0; // T.K dicide parameter value // for low energy
695  G4double rmax0 = bmax + 4.0;
696  G4double rmax = std::sqrt( rmax0*rmax0 + b*b );
697 
698  G4double ccoul = 0.001439767;
699  G4double pcca = 1.0 - double ( zp * zt ) * ccoul / eccm / rmax - ( b / rmax )*( b / rmax );
700 
701  G4double pccf = std::sqrt( pcca );
702 
703  //Fix for neutral particles
704  G4double aas1 = 0.0;
705  G4double bbs = 0.0;
706 
707  if ( zp != 0 )
708  {
709  G4double aas = 2.0 * eccm * b / double ( zp * zt ) / ccoul;
710  bbs = 1.0 / std::sqrt ( 1.0 + aas*aas );
711  aas1 = ( 1.0 + aas * b / rmax ) * bbs;
712  }
713 
714  G4double cost = 0.0;
715  G4double sint = 0.0;
716  G4double thet1 = 0.0;
717  G4double thet2 = 0.0;
718  if ( 1.0 - aas1*aas1 <= 0 || 1.0 - bbs*bbs <= 0.0 )
719  {
720  cost = 1.0;
721  sint = 0.0;
722  }
723  else
724  {
725  G4double aat1 = aas1 / std::sqrt ( 1.0 - aas1*aas1 );
726  G4double aat2 = bbs / std::sqrt ( 1.0 - bbs*bbs );
727 
728  thet1 = std::atan ( aat1 );
729  thet2 = std::atan ( aat2 );
730 
731 // TK enter to else block
732  G4double theta = thet1 - thet2;
733  cost = std::cos( theta );
734  sint = std::sin( theta );
735  }
736 
737  G4double rzpr = -rmax * cost * ( mass_targ ) / ( mass_proj + mass_targ );
738  G4double rzta = rmax * cost * ( mass_proj ) / ( mass_proj + mass_targ );
739 
740  G4double rxpr = rmax / 2.0 * sint;
741 
742  G4double rxta = -rxpr;
743 
744 
745  G4double pzpc = pzcc * ( cost * pccf + sint * b / rmax );
746  G4double pxpr = pzcc * ( -sint * pccf + cost * b / rmax );
747 
748  G4double pztc = - pzpc;
749  G4double pxta = - pxpr;
750 
751  G4double epc = std::sqrt ( pzpc*pzpc + pxpr*pxpr + mass_proj*mass_proj );
752  G4double etc = std::sqrt ( pztc*pztc + pxta*pxta + mass_targ*mass_targ );
753 
754  G4double pzpr = pzpc;
755  G4double pzta = pztc;
756  G4double epr = epc;
757  G4double eta = etc;
758 
759 // CM -> NN
760  G4double gammacm = boostToCM.gamma();
761  //G4double betacm = -boostToCM.beta();
762  G4double betacm = boostToCM.z();
763  pzpr = pzpc + betacm * gammacm * ( gammacm / ( 1. + gammacm ) * pzpc * betacm + epc );
764  pzta = pztc + betacm * gammacm * ( gammacm / ( 1. + gammacm ) * pztc * betacm + etc );
765  epr = gammacm * ( epc + betacm * pzpc );
766  eta = gammacm * ( etc + betacm * pztc );
767 
768  //G4double betpr = pzpr / epr;
769  //G4double betta = pzta / eta;
770 
771  G4double gammpr = epr / ( mass_proj );
772  G4double gammta = eta / ( mass_targ );
773 
774  pzta = pzta / double ( at );
775  pxta = pxta / double ( at );
776 
777  pzpr = pzpr / double ( ap );
778  pxpr = pxpr / double ( ap );
779 
780  G4double zeroz = 0.0;
781 
782  rzpr = rzpr -zeroz;
783  rzta = rzta -zeroz;
784 
785  // Set results
791 
797 
798 }
799 
800 
801 
803 {
804 
805  if ( gem == true )
807  else
809 
810 }
811 
812 void G4QMDReaction::ModelDescription(std::ostream& outFile) const
813 {
814  outFile << "Lorentz covarianted Quantum Molecular Dynamics model for nucleus (particle) vs nucleus reactions\n";
815 }
G4VCrossSectionDataSet * GetCrossSectionDataSet(const G4String &name, G4bool warning=true)
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4double coulomb_collision_rx_proj
G4int GetNOCollision()
Definition: G4QMDSystem.hh:65
G4ThreeVector GetPosition()
void SetDefaultChannel()
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
static const double MeV
Definition: G4SIunits.hh:211
G4double coulomb_collision_px_proj
G4double coulomb_collision_rx_targ
G4double coulomb_collision_pz_targ
G4Evaporation * evaporation
G4double powN(G4double x, G4int n) const
Definition: G4Pow.cc:128
void SetParticipant(G4QMDParticipant *particle)
Definition: G4QMDSystem.hh:51
CLHEP::Hep3Vector G4ThreeVector
void setEvaporationCh()
G4PiNuclearCrossSection * piNucXS
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:491
G4QMDSystem * system
static const G4double e2
void Clear()
Definition: G4QMDSystem.cc:68
G4ParticleDefinition * GetDefinition() const
const G4ParticleDefinition * GetDefinition()
void SetNucleus(G4QMDNucleus *aSystem)
int G4int
Definition: G4Types.hh:78
static const char * Default_Name()
G4ReactionProductVector * BreakItUp(const G4Fragment &theInitialState)
G4double deltaT
const G4String & GetParticleName() const
G4double envelopF
G4int GetAtomicNumber() const
void SetStatusChange(G4HadFinalStateStatus aS)
G4ThreeVector GetMomentum()
std::vector< G4ReactionProduct * > G4ReactionProductVector
void SetGEMChannel()
void SetMeanField(G4QMDMeanField *meanfield)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const G4ParticleDefinition * GetDefinition() const
bool G4bool
Definition: G4Types.hh:79
G4double GetIonMass(G4int Z, G4int A, G4int L=0, G4int lvl=0) const
Definition: G4IonTable.cc:1277
G4double iz
Definition: TRTMaterials.hh:39
G4double GetKineticEnergy() const
G4double coulomb_collision_px_targ
static G4CrossSectionDataSetRegistry * Instance()
static const double GeV
Definition: G4SIunits.hh:214
const G4String & GetParticleType() const
void CalKinematicsOfBinaryCollisions(G4double)
G4int GetTotalNumberOfParticipant()
Definition: G4QMDSystem.hh:60
G4double GetTotalPotential()
const G4LorentzVector & Get4Momentum() const
G4int GetAtomicMass() const
G4LorentzVector Get4Momentum() const
std::vector< G4QMDNucleus * > DoClusterJudgment()
G4QMDParticipant * GetParticipant(G4int i)
Definition: G4QMDSystem.hh:62
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:78
static const G4double e1
G4QMDMeanField * meanField
G4double coulomb_collision_rz_targ
void calcOffSetOfCollision(G4double, const G4ParticleDefinition *, const G4ParticleDefinition *, G4double, G4double, G4double, G4ThreeVector)
G4IonsShenCrossSection * shenXS
G4double GetPDGMass() const
static const double pi
Definition: G4SIunits.hh:74
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4double coulomb_collision_gamma_targ
G4LorentzVector Get4Momentum()
void DoPropagation(G4double)
void SetSystem(G4QMDSystem *aSystem)
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
void SetTotalPotential(G4double x)
Definition: G4QMDNucleus.hh:62
G4double coulomb_collision_rz_proj
G4ExcitationHandler * excitationHandler
void SetEvaporation(G4VEvaporation *ptr)
G4QMDCollision * collision
#define G4endl
Definition: G4ios.hh:61
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
G4double coulomb_collision_gamma_proj
double G4double
Definition: G4Types.hh:76
void CalEnergyAndAngularMomentumInCM()
static const double eplus
Definition: G4SIunits.hh:196
G4double GetPDGCharge() const
static const double fermi
Definition: G4SIunits.hh:102
const G4double r0
G4ThreeVector GetMomentum() const
G4double coulomb_collision_pz_proj
G4double GetTotalMomentum() const
G4double elastic(Particle const *const p1, Particle const *const p2)
virtual void ModelDescription(std::ostream &outFile) const
CLHEP::HepLorentzVector G4LorentzVector