Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4QMDReaction Class Reference

#include <G4QMDReaction.hh>

Inheritance diagram for G4QMDReaction:
Collaboration diagram for G4QMDReaction:

Public Member Functions

 G4QMDReaction ()
 
 ~G4QMDReaction ()
 
std::vector< G4QMDSystem * > GetFinalStates ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4ExcitationHandlerGetExcitationHandler ()
 
void UnUseGEM ()
 
void UseFRAG ()
 
void SetTMAX (G4int i)
 
void SetDT (G4double t)
 
void SetEF (G4double x)
 
virtual void ModelDescription (std::ostream &outFile) const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
virtual const std::pair
< G4double, G4double
GetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double,
G4double
GetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void InitialiseModel ()
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 60 of file G4QMDReaction.hh.

Constructor & Destructor Documentation

G4QMDReaction::G4QMDReaction ( )

Definition at line 47 of file G4QMDReaction.cc.

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
58  shenXS = new G4IonsShenCrossSection();
59  //genspaXS = new G4GeneralSpaceNNCrossSection();
61  meanField = new G4QMDMeanField();
62  collision = new G4QMDCollision();
63 
64  excitationHandler = new G4ExcitationHandler;
65  excitationHandler->SetDeexChannelsType( fCombined );
66  //fEvaporation - 8 default channels
67  //fCombined - 8 default + 60 GEM
68  //fGEM - 2 default + 66 GEM
69  evaporation = new G4Evaporation;
70  excitationHandler->SetEvaporation( evaporation );
71  setEvaporationCh();
72 
73  coulomb_collision_gamma_proj = 0.0;
74  coulomb_collision_rx_proj = 0.0;
75  coulomb_collision_rz_proj = 0.0;
76  coulomb_collision_px_proj = 0.0;
77  coulomb_collision_pz_proj = 0.0;
78 
79  coulomb_collision_gamma_targ = 0.0;
80  coulomb_collision_rx_targ = 0.0;
81  coulomb_collision_rz_targ = 0.0;
82  coulomb_collision_px_targ = 0.0;
83  coulomb_collision_pz_targ = 0.0;
84 
85 }
G4VCrossSectionDataSet * GetCrossSectionDataSet(const G4String &name, G4bool warning=true)
static const char * Default_Name()
static G4CrossSectionDataSetRegistry * Instance()
G4HadronicInteraction(const G4String &modelName="HadronicModel")
void SetEvaporation(G4VEvaporation *ptr, G4bool isLocal=false)
void SetDeexChannelsType(G4DeexChannelType val)

Here is the call graph for this function:

G4QMDReaction::~G4QMDReaction ( )

Definition at line 89 of file G4QMDReaction.cc.

90 {
91  delete evaporation;
92  delete excitationHandler;
93  delete collision;
94  delete meanField;
95 }

Member Function Documentation

G4HadFinalState * G4QMDReaction::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
)
virtual

Implements G4HadronicInteraction.

Definition at line 99 of file G4QMDReaction.cc.

100 {
101  //G4cout << "G4QMDReaction::ApplyYourself" << G4endl;
102 
104 
105  system = new G4QMDSystem;
106 
107  G4int proj_Z = 0;
108  G4int proj_A = 0;
109  const G4ParticleDefinition* proj_pd = ( const G4ParticleDefinition* ) projectile.GetDefinition();
110  if ( proj_pd->GetParticleType() == "nucleus" )
111  {
112  proj_Z = proj_pd->GetAtomicNumber();
113  proj_A = proj_pd->GetAtomicMass();
114  }
115  else
116  {
117  proj_Z = (int)( proj_pd->GetPDGCharge()/eplus );
118  proj_A = 1;
119  }
120  //G4int targ_Z = int ( target.GetZ() + 0.5 );
121  //G4int targ_A = int ( target.GetN() + 0.5 );
122  //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this)
123  G4int targ_Z = target.GetZ_asInt();
124  G4int targ_A = target.GetA_asInt();
125  const G4ParticleDefinition* targ_pd = G4IonTable::GetIonTable()->GetIon( targ_Z , targ_A , 0.0 );
126 
127 
128  //G4NistManager* nistMan = G4NistManager::Instance();
129 // G4Element* G4NistManager::FindOrBuildElement( targ_Z );
130 
131  const G4DynamicParticle* proj_dp = new G4DynamicParticle ( proj_pd , projectile.Get4Momentum() );
132  //const G4Element* targ_ele = nistMan->FindOrBuildElement( targ_Z );
133  //G4double aTemp = projectile.GetMaterial()->GetTemperature();
134 
135  //090331
136 
137  G4VCrossSectionDataSet* theXS = shenXS;
138 
139  if ( proj_pd->GetParticleType() == "meson" ) theXS = piNucXS;
140 
141  //G4double xs_0 = genspaXS->GetCrossSection ( proj_dp , targ_ele , aTemp );
142  //G4double xs_0 = theXS->GetCrossSection ( proj_dp , targ_ele , aTemp );
143  //110822
144  G4double xs_0 = theXS->GetIsoCrossSection ( proj_dp , targ_Z , targ_A );
145 
146  G4double bmax_0 = std::sqrt( xs_0 / pi );
147  //std::cout << "bmax_0 in fm (fermi) " << bmax_0/fermi << std::endl;
148 
149  //delete proj_dp;
150 
151  G4bool elastic = true;
152 
153  std::vector< G4QMDNucleus* > nucleuses; // Secondary nuceluses
154  G4ThreeVector boostToReac; // ReactionSystem (CM or NN);
155  G4ThreeVector boostBackToLAB; // Reaction System to LAB;
156 
157  G4LorentzVector targ4p( G4ThreeVector( 0.0 ) , targ_pd->GetPDGMass()/GeV );
158  G4ThreeVector boostLABtoCM = targ4p.findBoostToCM( proj_dp->Get4Momentum()/GeV ); // CM of target and proj;
159 
160  G4double p1 = proj_dp->GetMomentum().mag()/GeV/proj_A;
161  G4double m1 = proj_dp->GetDefinition()->GetPDGMass()/GeV/proj_A;
162  G4double e1 = std::sqrt( p1*p1 + m1*m1 );
163  G4double e2 = targ_pd->GetPDGMass()/GeV/targ_A;
164  G4double beta_nn = -p1 / ( e1+e2 );
165 
166  G4ThreeVector boostLABtoNN ( 0. , 0. , beta_nn ); // CM of NN;
167 
168  G4double beta_nncm = ( - boostLABtoCM.beta() + boostLABtoNN.beta() ) / ( 1 - boostLABtoCM.beta() * boostLABtoNN.beta() ) ;
169 
170  //std::cout << targ4p << std::endl;
171  //std::cout << proj_dp->Get4Momentum()<< std::endl;
172  //std::cout << beta_nncm << std::endl;
173  G4ThreeVector boostNNtoCM( 0. , 0. , beta_nncm ); //
174  G4ThreeVector boostCMtoNN( 0. , 0. , -beta_nncm ); //
175 
176  boostToReac = boostLABtoNN;
177  boostBackToLAB = -boostLABtoNN;
178 
179  delete proj_dp;
180 
181  G4int icounter = 0;
182  G4int icounter_max = 1024;
183  while ( elastic ) // Loop checking, 11.03.2015, T. Koi
184  {
185  icounter++;
186  if ( icounter > icounter_max ) {
187  G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
188  break;
189  }
190 
191 // impact parameter
192  //G4double bmax = 1.05*(bmax_0/fermi); // 10% for Peripheral reactions
193  G4double bmax = envelopF*(bmax_0/fermi);
194  G4double b = bmax * std::sqrt ( G4UniformRand() );
195 //071112
196  //G4double b = 0;
197  //G4double b = bmax;
198  //G4double b = bmax/1.05 * 0.7 * G4UniformRand();
199 
200  //G4cout << "G4QMDRESULT bmax_0 = " << bmax_0/fermi << " fm, bmax = " << bmax << " fm , b = " << b << " fm " << G4endl;
201 
202  G4double plab = projectile.GetTotalMomentum()/GeV;
203  G4double elab = ( projectile.GetKineticEnergy() + proj_pd->GetPDGMass() + targ_pd->GetPDGMass() )/GeV;
204 
205  calcOffSetOfCollision( b , proj_pd , targ_pd , plab , elab , bmax , boostCMtoNN );
206 
207 // Projectile
208  G4LorentzVector proj4pLAB = projectile.Get4Momentum()/GeV;
209 
210  G4QMDGroundStateNucleus* proj(NULL);
211  if ( projectile.GetDefinition()->GetParticleType() == "nucleus"
212  || projectile.GetDefinition()->GetParticleName() == "proton"
213  || projectile.GetDefinition()->GetParticleName() == "neutron" )
214  {
215 
216  proj_Z = proj_pd->GetAtomicNumber();
217  proj_A = proj_pd->GetAtomicMass();
218 
219  proj = new G4QMDGroundStateNucleus( proj_Z , proj_A );
220  //proj->ShowParticipants();
221 
222 
223  meanField->SetSystem ( proj );
224  proj->SetTotalPotential( meanField->GetTotalPotential() );
225  proj->CalEnergyAndAngularMomentumInCM();
226 
227  }
228 
229 // Target
230  //G4int iz = int ( target.GetZ() );
231  //G4int ia = int ( target.GetN() );
232  //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this)
233  G4int iz = int ( target.GetZ_asInt() );
234  G4int ia = int ( target.GetA_asInt() );
235 
236  G4QMDGroundStateNucleus* targ = new G4QMDGroundStateNucleus( iz , ia );
237 
238  meanField->SetSystem (targ );
239  targ->SetTotalPotential( meanField->GetTotalPotential() );
240  targ->CalEnergyAndAngularMomentumInCM();
241 
242  //G4LorentzVector targ4p( G4ThreeVector( 0.0 ) , targ->GetNuclearMass()/GeV );
243 // Boost Vector to CM
244  //boostToCM = targ4p.findBoostToCM( proj4pLAB );
245 
246 // Target
247  for ( G4int i = 0 ; i < targ->GetTotalNumberOfParticipant() ; i++ )
248  {
249 
250  G4ThreeVector p0 = targ->GetParticipant( i )->GetMomentum();
251  G4ThreeVector r0 = targ->GetParticipant( i )->GetPosition();
252 
253  G4ThreeVector p ( p0.x() + coulomb_collision_px_targ
254  , p0.y()
255  , p0.z() * coulomb_collision_gamma_targ + coulomb_collision_pz_targ );
256 
257  G4ThreeVector r ( r0.x() + coulomb_collision_rx_targ
258  , r0.y()
259  , r0.z() / coulomb_collision_gamma_targ + coulomb_collision_rz_targ );
260 
261  system->SetParticipant( new G4QMDParticipant( targ->GetParticipant( i )->GetDefinition() , p , r ) );
262  system->GetParticipant( i )->SetTarget();
263 
264  }
265 
266  G4LorentzVector proj4pCM = CLHEP::boostOf ( proj4pLAB , boostToReac );
267  G4LorentzVector targ4pCM = CLHEP::boostOf ( targ4p , boostToReac );
268 
269 // Projectile
270  if ( proj != NULL )
271  {
272 
273 // projectile is nucleus
274 
275  for ( G4int i = 0 ; i < proj->GetTotalNumberOfParticipant() ; i++ )
276  {
277 
278  G4ThreeVector p0 = proj->GetParticipant( i )->GetMomentum();
279  G4ThreeVector r0 = proj->GetParticipant( i )->GetPosition();
280 
281  G4ThreeVector p ( p0.x() + coulomb_collision_px_proj
282  , p0.y()
283  , p0.z() * coulomb_collision_gamma_proj + coulomb_collision_pz_proj );
284 
285  G4ThreeVector r ( r0.x() + coulomb_collision_rx_proj
286  , r0.y()
287  , r0.z() / coulomb_collision_gamma_proj + coulomb_collision_rz_proj );
288 
289  system->SetParticipant( new G4QMDParticipant( proj->GetParticipant( i )->GetDefinition() , p , r ) );
290  system->GetParticipant ( i + targ->GetTotalNumberOfParticipant() )->SetProjectile();
291  }
292 
293  }
294  else
295  {
296 
297 // projectile is particle
298 
299  // avoid multiple set in "elastic" loop
300  if ( system->GetTotalNumberOfParticipant() == targ->GetTotalNumberOfParticipant() )
301  {
302 
303  G4int i = targ->GetTotalNumberOfParticipant();
304 
305  G4ThreeVector p0( 0 );
306  G4ThreeVector r0( 0 );
307 
308  G4ThreeVector p ( p0.x() + coulomb_collision_px_proj
309  , p0.y()
310  , p0.z() * coulomb_collision_gamma_proj + coulomb_collision_pz_proj );
311 
312  G4ThreeVector r ( r0.x() + coulomb_collision_rx_proj
313  , r0.y()
314  , r0.z() / coulomb_collision_gamma_proj + coulomb_collision_rz_proj );
315 
316  system->SetParticipant( new G4QMDParticipant( (G4ParticleDefinition*)projectile.GetDefinition() , p , r ) );
317  // This is not important becase only 1 projectile particle.
318  system->GetParticipant ( i )->SetProjectile();
319  }
320 
321  }
322  //system->ShowParticipants();
323 
324  delete targ;
325  delete proj;
326 
327  meanField->SetSystem ( system );
328  collision->SetMeanField ( meanField );
329 
330 // Time Evolution
331  //std::cout << "Start time evolution " << std::endl;
332  //system->ShowParticipants();
333  for ( G4int i = 0 ; i < maxTime ; i++ )
334  {
335  //G4cout << " do Paropagate " << i << " th time step. " << G4endl;
336  meanField->DoPropagation( deltaT );
337  //system->ShowParticipants();
338  collision->CalKinematicsOfBinaryCollisions( deltaT );
339 
340  if ( i / 10 * 10 == i )
341  {
342  //G4cout << i << " th time step. " << G4endl;
343  //system->ShowParticipants();
344  }
345  //system->ShowParticipants();
346  }
347  //system->ShowParticipants();
348 
349 
350  //std::cout << "Doing Cluster Judgment " << std::endl;
351 
352  nucleuses = meanField->DoClusterJudgment();
353 
354 // Elastic Judgment
355 
356  G4int numberOfSecondary = int ( nucleuses.size() ) + system->GetTotalNumberOfParticipant();
357 
358  G4int sec_a_Z = 0;
359  G4int sec_a_A = 0;
360  const G4ParticleDefinition* sec_a_pd = NULL;
361  G4int sec_b_Z = 0;
362  G4int sec_b_A = 0;
363  const G4ParticleDefinition* sec_b_pd = NULL;
364 
365  if ( numberOfSecondary == 2 )
366  {
367 
368  G4bool elasticLike_system = false;
369  if ( nucleuses.size() == 2 )
370  {
371 
372  sec_a_Z = nucleuses[0]->GetAtomicNumber();
373  sec_a_A = nucleuses[0]->GetMassNumber();
374  sec_b_Z = nucleuses[1]->GetAtomicNumber();
375  sec_b_A = nucleuses[1]->GetMassNumber();
376 
377  if ( ( sec_a_Z == proj_Z && sec_a_A == proj_A && sec_b_Z == targ_Z && sec_b_A == targ_A )
378  || ( sec_a_Z == targ_Z && sec_a_A == targ_A && sec_b_Z == proj_Z && sec_b_A == proj_A ) )
379  {
380  elasticLike_system = true;
381  }
382 
383  }
384  else if ( nucleuses.size() == 1 )
385  {
386 
387  sec_a_Z = nucleuses[0]->GetAtomicNumber();
388  sec_a_A = nucleuses[0]->GetMassNumber();
389  sec_b_pd = system->GetParticipant( 0 )->GetDefinition();
390 
391  if ( ( sec_a_Z == proj_Z && sec_a_A == proj_A && sec_b_pd == targ_pd )
392  || ( sec_a_Z == targ_Z && sec_a_A == targ_A && sec_b_pd == proj_pd ) )
393  {
394  elasticLike_system = true;
395  }
396 
397  }
398  else
399  {
400 
401  sec_a_pd = system->GetParticipant( 0 )->GetDefinition();
402  sec_b_pd = system->GetParticipant( 1 )->GetDefinition();
403 
404  if ( ( sec_a_pd == proj_pd && sec_b_pd == targ_pd )
405  || ( sec_a_pd == targ_pd && sec_b_pd == proj_pd ) )
406  {
407  elasticLike_system = true;
408  }
409 
410  }
411 
412  if ( elasticLike_system == true )
413  {
414 
415  G4bool elasticLike_energy = true;
416 // Cal ExcitationEnergy
417  for ( G4int i = 0 ; i < int ( nucleuses.size() ) ; i++ )
418  {
419 
420  //meanField->SetSystem( nucleuses[i] );
421  meanField->SetNucleus( nucleuses[i] );
422  //nucleuses[i]->SetTotalPotential( meanField->GetTotalPotential() );
423  //nucleuses[i]->CalEnergyAndAngularMomentumInCM();
424 
425  if ( nucleuses[i]->GetExcitationEnergy()*GeV > 1.0*MeV ) elasticLike_energy = false;
426 
427  }
428 
429 // Check Collision
430  G4bool withCollision = true;
431  if ( system->GetNOCollision() == 0 ) withCollision = false;
432 
433 // Final judegement for Inelasitc or Elastic;
434 //
435 // ElasticLike without Collision
436  //if ( elasticLike_energy == true && withCollision == false ) elastic = true; // ielst = 0
437 // ElasticLike with Collision
438  //if ( elasticLike_energy == true && withCollision == true ) elastic = true; // ielst = 1
439 // InelasticLike without Collision
440  //if ( elasticLike_energy == false ) elastic = false; // ielst = 2
441  if ( frag == true )
442  if ( elasticLike_energy == false ) elastic = false;
443 // InelasticLike with Collision
444  if ( elasticLike_energy == false && withCollision == true ) elastic = false; // ielst = 3
445 
446  }
447 
448  }
449  else
450  {
451 
452 // numberOfSecondary != 2
453  elastic = false;
454 
455  }
456 
457 //071115
458  //G4cout << elastic << G4endl;
459  // if elastic is true try again from sampling of impact parameter
460 
461  if ( elastic == true )
462  {
463  // delete this nucleues
464  for ( std::vector< G4QMDNucleus* >::iterator
465  it = nucleuses.begin() ; it != nucleuses.end() ; it++ )
466  {
467  delete *it;
468  }
469  nucleuses.clear();
470  }
471 
472  }
473 
474 
475 // Statical Decay Phase
476 
477  for ( std::vector< G4QMDNucleus* >::iterator it
478  = nucleuses.begin() ; it != nucleuses.end() ; it++ )
479  {
480 
481 /*
482  G4cout << "G4QMDRESULT "
483  << (*it)->GetAtomicNumber()
484  << " "
485  << (*it)->GetMassNumber()
486  << " "
487  << (*it)->Get4Momentum()
488  << " "
489  << (*it)->Get4Momentum().vect()
490  << " "
491  << (*it)->Get4Momentum().restMass()
492  << " "
493  << (*it)->GetNuclearMass()/GeV
494  << G4endl;
495 */
496 
497  meanField->SetNucleus ( *it );
498 
499  if ( (*it)->GetAtomicNumber() == 0 // neutron cluster
500  || (*it)->GetAtomicNumber() == (*it)->GetMassNumber() ) // proton cluster
501  {
502  // push back system
503  for ( G4int i = 0 ; i < (*it)->GetTotalNumberOfParticipant() ; i++ )
504  {
505  G4QMDParticipant* aP = new G4QMDParticipant( ( (*it)->GetParticipant( i ) )->GetDefinition() , ( (*it)->GetParticipant( i ) )->GetMomentum() , ( (*it)->GetParticipant( i ) )->GetPosition() );
506  system->SetParticipant ( aP );
507  }
508  continue;
509  }
510 
511  G4double nucleus_e = std::sqrt ( G4Pow::GetInstance()->powN ( (*it)->GetNuclearMass()/GeV , 2 ) + G4Pow::GetInstance()->powN ( (*it)->Get4Momentum().vect().mag() , 2 ) );
512  G4LorentzVector nucleus_p4CM ( (*it)->Get4Momentum().vect() , nucleus_e );
513 
514 // std::cout << "G4QMDRESULT nucleus deltaQ " << deltaQ << std::endl;
515 
516  G4int ia = (*it)->GetMassNumber();
517  G4int iz = (*it)->GetAtomicNumber();
518 
519  G4LorentzVector lv ( G4ThreeVector( 0.0 ) , (*it)->GetExcitationEnergy()*GeV + G4IonTable::GetIonTable()->GetIonMass( iz , ia ) );
520 
521  G4Fragment* aFragment = new G4Fragment( ia , iz , lv );
522 
524  rv = excitationHandler->BreakItUp( *aFragment );
525  G4bool notBreak = true;
526  for ( G4ReactionProductVector::iterator itt
527  = rv->begin() ; itt != rv->end() ; itt++ )
528  {
529 
530  notBreak = false;
531  // Secondary from this nucleus (*it)
532  const G4ParticleDefinition* pd = (*itt)->GetDefinition();
533 
534  G4LorentzVector p4 ( (*itt)->GetMomentum()/GeV , (*itt)->GetTotalEnergy()/GeV ); //in nucleus(*it) rest system
535  G4LorentzVector p4_CM = CLHEP::boostOf( p4 , -nucleus_p4CM.findBoostToCM() ); // Back to CM
536  G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB ); // Back to LAB
537 
538 
539 //090122
540  //theParticleChange.AddSecondary( dp );
541  if ( !( pd->GetAtomicNumber() == 4 && pd->GetAtomicMass() == 8 ) )
542  {
543  G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV );
545  }
546  else
547  {
548  //Be8 -> Alpha + Alpha + Q
549  G4ThreeVector randomized_direction( G4UniformRand() , G4UniformRand() , G4UniformRand() );
550  randomized_direction = randomized_direction.unit();
551  G4double q_decay = (*itt)->GetMass() - 2*G4Alpha::Alpha()->GetPDGMass();
552  G4double p_decay = std::sqrt ( G4Pow::GetInstance()->powN(G4Alpha::Alpha()->GetPDGMass()+q_decay/2,2) - G4Pow::GetInstance()->powN(G4Alpha::Alpha()->GetPDGMass() , 2 ) );
553  G4LorentzVector p4_a1 ( p_decay*randomized_direction , G4Alpha::Alpha()->GetPDGMass()+q_decay/2 ); //in Be8 rest system
554 
555  G4LorentzVector p4_a1_Be8 = CLHEP::boostOf ( p4_a1/GeV , -p4.findBoostToCM() );
556  G4LorentzVector p4_a1_CM = CLHEP::boostOf ( p4_a1_Be8 , -nucleus_p4CM.findBoostToCM() );
557  G4LorentzVector p4_a1_LAB = CLHEP::boostOf ( p4_a1_CM , boostBackToLAB );
558 
559  G4LorentzVector p4_a2 ( -p_decay*randomized_direction , G4Alpha::Alpha()->GetPDGMass()+q_decay/2 ); //in Be8 rest system
560 
561  G4LorentzVector p4_a2_Be8 = CLHEP::boostOf ( p4_a2/GeV , -p4.findBoostToCM() );
562  G4LorentzVector p4_a2_CM = CLHEP::boostOf ( p4_a2_Be8 , -nucleus_p4CM.findBoostToCM() );
563  G4LorentzVector p4_a2_LAB = CLHEP::boostOf ( p4_a2_CM , boostBackToLAB );
564 
565  G4DynamicParticle* dp1 = new G4DynamicParticle( G4Alpha::Alpha() , p4_a1_LAB*GeV );
566  G4DynamicParticle* dp2 = new G4DynamicParticle( G4Alpha::Alpha() , p4_a2_LAB*GeV );
569  }
570 //090122
571 
572 /*
573  G4cout
574  << "Regist Secondary "
575  << (*itt)->GetDefinition()->GetParticleName()
576  << " "
577  << (*itt)->GetMomentum()/GeV
578  << " "
579  << (*itt)->GetKineticEnergy()/GeV
580  << " "
581  << (*itt)->GetMass()/GeV
582  << " "
583  << (*itt)->GetTotalEnergy()/GeV
584  << " "
585  << (*itt)->GetTotalEnergy()/GeV * (*itt)->GetTotalEnergy()/GeV
586  - (*itt)->GetMomentum()/GeV * (*itt)->GetMomentum()/GeV
587  << " "
588  << nucleus_p4CM.findBoostToCM()
589  << " "
590  << p4
591  << " "
592  << p4_CM
593  << " "
594  << p4_LAB
595  << G4endl;
596 */
597 
598  }
599  if ( notBreak == true )
600  {
601 
602  const G4ParticleDefinition* pd = G4IonTable::GetIonTable()->GetIon( (*it)->GetAtomicNumber() , (*it)->GetMassNumber(), (*it)->GetExcitationEnergy()*GeV );
603  G4LorentzVector p4_CM = nucleus_p4CM;
604  G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB ); // Back to LAB
605  G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV );
607 
608  }
609 
610  for ( G4ReactionProductVector::iterator itt
611  = rv->begin() ; itt != rv->end() ; itt++ )
612  {
613  delete *itt;
614  }
615  delete rv;
616 
617  delete aFragment;
618 
619  }
620 
621 
622 
623  for ( G4int i = 0 ; i < system->GetTotalNumberOfParticipant() ; i++ )
624  {
625 
626  // Secondary particles
627 
628  const G4ParticleDefinition* pd = system->GetParticipant( i )->GetDefinition();
629  G4LorentzVector p4_CM = system->GetParticipant( i )->Get4Momentum();
630  G4LorentzVector p4_LAB = CLHEP::boostOf( p4_CM , boostBackToLAB );
631  G4DynamicParticle* dp = new G4DynamicParticle( pd , p4_LAB*GeV );
633 
634 /*
635  G4cout << "G4QMDRESULT "
636  << "r" << i << " " << system->GetParticipant ( i ) -> GetPosition() << " "
637  << "p" << i << " " << system->GetParticipant ( i ) -> Get4Momentum()
638  << G4endl;
639 */
640 
641  }
642 
643  for ( std::vector< G4QMDNucleus* >::iterator it
644  = nucleuses.begin() ; it != nucleuses.end() ; it++ )
645  {
646  delete *it; // delete nulceuse
647  }
648  nucleuses.clear();
649 
650  system->Clear();
651  delete system;
652 
654 
655  return &theParticleChange;
656 
657 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4int GetNOCollision()
Definition: G4QMDSystem.hh:65
const XML_Char * target
Definition: expat.h:268
G4double powN(G4double x, G4int n) const
Definition: G4Pow.cc:128
void SetParticipant(G4QMDParticipant *particle)
Definition: G4QMDSystem.hh:51
CLHEP::Hep3Vector G4ThreeVector
double x() const
const char * p
Definition: xmltok.h:285
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:503
void Clear()
Definition: G4QMDSystem.cc:68
G4ParticleDefinition * GetDefinition() const
const G4ParticleDefinition * GetDefinition()
void SetNucleus(G4QMDNucleus *aSystem)
int G4int
Definition: G4Types.hh:78
G4ReactionProductVector * BreakItUp(const G4Fragment &theInitialState)
double z() const
G4int GetAtomicNumber() const
void SetStatusChange(G4HadFinalStateStatus aS)
std::vector< G4ReactionProduct * > G4ReactionProductVector
void SetMeanField(G4QMDMeanField *meanfield)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
G4double GetIonMass(G4int Z, G4int A, G4int L=0, G4int lvl=0) const
Definition: G4IonTable.cc:1335
HepLorentzVector boostOf(const HepLorentzVector &vec, const Hep3Vector &betaVector)
static constexpr double eplus
Definition: G4SIunits.hh:199
typedef int(XMLCALL *XML_NotStandaloneHandler)(void *userData)
const G4String & GetParticleType() const
void CalKinematicsOfBinaryCollisions(G4double)
G4int GetTotalNumberOfParticipant()
Definition: G4QMDSystem.hh:60
G4double GetTotalPotential()
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
Hep3Vector findBoostToCM() const
G4double GetPDGMass() const
G4LorentzVector Get4Momentum()
void DoPropagation(G4double)
void SetSystem(G4QMDSystem *aSystem)
double y() const
static constexpr double GeV
Definition: G4SIunits.hh:217
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
static constexpr double pi
Definition: G4SIunits.hh:75
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
double G4double
Definition: G4Types.hh:76
static constexpr double fermi
Definition: G4SIunits.hh:103
G4double GetPDGCharge() const
double mag() const
G4ThreeVector GetMomentum() const
G4double elastic(Particle const *const p1, Particle const *const p2)
G4ExcitationHandler* G4QMDReaction::GetExcitationHandler ( )
inline

Definition at line 70 of file G4QMDReaction.hh.

70 { return excitationHandler; };
std::vector< G4QMDSystem* > G4QMDReaction::GetFinalStates ( )
void G4QMDReaction::ModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 816 of file G4QMDReaction.cc.

817 {
818  outFile << "Lorentz covarianted Quantum Molecular Dynamics model for nucleus (particle) vs nucleus reactions\n";
819 }
void G4QMDReaction::SetDT ( G4double  t)
inline

Definition at line 76 of file G4QMDReaction.hh.

76 { deltaT = t; };
void G4QMDReaction::SetEF ( G4double  x)
inline

Definition at line 77 of file G4QMDReaction.hh.

77 { envelopF = x; };
void G4QMDReaction::SetTMAX ( G4int  i)
inline

Definition at line 75 of file G4QMDReaction.hh.

75 { maxTime = i; };
void G4QMDReaction::UnUseGEM ( )
inline

Definition at line 72 of file G4QMDReaction.hh.

72 { gem = false; setEvaporationCh(); };
void G4QMDReaction::UseFRAG ( )
inline

Definition at line 73 of file G4QMDReaction.hh.

73 { frag = true; };

The documentation for this class was generated from the following files: