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