Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QMDReaction Class Reference

#include <G4QMDReaction.hh>

+ Inheritance diagram for G4QMDReaction:

Public Member Functions

 G4QMDReaction ()
 
 ~G4QMDReaction () override
 
std::vector< G4QMDSystem * > GetFinalStates ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) override
 
G4ExcitationHandlerGetExcitationHandler ()
 
void UnUseGEM ()
 
void UseFRAG ()
 
void SetTMAX (G4int i)
 
void SetDT (G4double t)
 
void SetEF (G4double x)
 
void ModelDescription (std::ostream &outFile) const override
 
 G4QMDReaction (const G4QMDReaction &right)=delete
 
const G4QMDReactionoperator= (const G4QMDReaction &right)=delete
 
- 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, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void InitialiseModel ()
 
 G4HadronicInteraction (const G4HadronicInteraction &right)=delete
 
const G4HadronicInteractionoperator= (const G4HadronicInteraction &right)=delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 

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 59 of file G4QMDReaction.hh.

Constructor & Destructor Documentation

◆ G4QMDReaction() [1/2]

G4QMDReaction::G4QMDReaction ( )

Definition at line 54 of file G4QMDReaction.cc.

55: G4HadronicInteraction("QMDModel")
56, system ( NULL )
57, deltaT ( 1 ) // in fsec (c=1)
58, maxTime ( 100 ) // will have maxTime-th time step
59, envelopF ( 1.05 ) // 10% for Peripheral reactions
60, gem ( true )
61, frag ( false )
62, secID( -1 )
63{
65 pipElNucXS = new G4BGGPionElasticXS(G4PionPlus::PionPlus() );
66 pipElNucXS->BuildPhysicsTable(*(G4PionPlus::PionPlus() ) );
67
68 pimElNucXS = new G4BGGPionElasticXS(G4PionMinus::PionMinus() );
69 pimElNucXS->BuildPhysicsTable(*(G4PionMinus::PionMinus() ) );
70
71 pipInelNucXS = new G4BGGPionInelasticXS(G4PionPlus::PionPlus() );
72 pipInelNucXS->BuildPhysicsTable(*(G4PionPlus::PionPlus() ) );
73
74 pimInelNucXS = new G4BGGPionInelasticXS(G4PionMinus::PionMinus() );
75 pimInelNucXS->BuildPhysicsTable(*(G4PionMinus::PionMinus() ) );
76
77 meanField = new G4QMDMeanField();
78 collision = new G4QMDCollision();
79
80 excitationHandler = new G4ExcitationHandler();
81 setEvaporationCh();
82
83 coulomb_collision_gamma_proj = 0.0;
84 coulomb_collision_rx_proj = 0.0;
85 coulomb_collision_rz_proj = 0.0;
86 coulomb_collision_px_proj = 0.0;
87 coulomb_collision_pz_proj = 0.0;
88
89 coulomb_collision_gamma_targ = 0.0;
90 coulomb_collision_rx_targ = 0.0;
91 coulomb_collision_rz_targ = 0.0;
92 coulomb_collision_px_targ = 0.0;
93 coulomb_collision_pz_targ = 0.0;
94
95 secID = G4PhysicsModelCatalog::GetModelID( "model_QMDModel" );
96}
void BuildPhysicsTable(const G4ParticleDefinition &) final
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static G4int GetModelID(const G4int modelIndex)
static G4PionMinus * PionMinus()
static G4PionPlus * PionPlus()
Definition G4PionPlus.cc:93

◆ ~G4QMDReaction()

G4QMDReaction::~G4QMDReaction ( )
override

Definition at line 99 of file G4QMDReaction.cc.

100{
101 delete excitationHandler;
102 delete collision;
103 delete meanField;
104}

◆ G4QMDReaction() [2/2]

G4QMDReaction::G4QMDReaction ( const G4QMDReaction & right)
delete

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 107 of file G4QMDReaction.cc.

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

◆ GetExcitationHandler()

G4ExcitationHandler * G4QMDReaction::GetExcitationHandler ( )
inline

Definition at line 70 of file G4QMDReaction.hh.

70{return excitationHandler;};

◆ GetFinalStates()

std::vector< G4QMDSystem * > G4QMDReaction::GetFinalStates ( )

◆ ModelDescription()

void G4QMDReaction::ModelDescription ( std::ostream & outFile) const
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 842 of file G4QMDReaction.cc.

843{
844 outFile << "Lorentz covarianted Quantum Molecular Dynamics model for nucleus (particle) vs nucleus reactions\n";
845}

◆ operator=()

const G4QMDReaction & G4QMDReaction::operator= ( const G4QMDReaction & right)
delete

◆ SetDT()

void G4QMDReaction::SetDT ( G4double t)
inline

Definition at line 76 of file G4QMDReaction.hh.

76{ deltaT = t; };

◆ SetEF()

void G4QMDReaction::SetEF ( G4double x)
inline

Definition at line 77 of file G4QMDReaction.hh.

77{ envelopF = x; };

◆ SetTMAX()

void G4QMDReaction::SetTMAX ( G4int i)
inline

Definition at line 75 of file G4QMDReaction.hh.

75{ maxTime = i; };

◆ UnUseGEM()

void G4QMDReaction::UnUseGEM ( )
inline

Definition at line 72 of file G4QMDReaction.hh.

72{gem = false; setEvaporationCh();};

◆ UseFRAG()

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: