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