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

#include <G4WilsonAbrasionModel.hh>

+ Inheritance diagram for G4WilsonAbrasionModel:

Public Member Functions

 G4WilsonAbrasionModel (G4bool useAblation1=false)
 
 G4WilsonAbrasionModel (G4ExcitationHandler *)
 
 ~G4WilsonAbrasionModel ()
 
 G4WilsonAbrasionModel (const G4WilsonAbrasionModel &right)
 
const G4WilsonAbrasionModeloperator= (G4WilsonAbrasionModel &right)
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &, G4Nucleus &)
 
void SetVerboseLevel (G4int)
 
void SetUseAblation (G4bool)
 
G4bool GetUseAblation ()
 
void SetConserveMomentum (G4bool)
 
G4bool GetConserveMomentum ()
 
void SetExcitationHandler (G4ExcitationHandler *)
 
G4ExcitationHandlerGetExcitationHandler ()
 
virtual void ModelDescription (std::ostream &) const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
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 ModelDescription (std::ostream &outFile) const
 
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 77 of file G4WilsonAbrasionModel.hh.

Constructor & Destructor Documentation

◆ G4WilsonAbrasionModel() [1/3]

G4WilsonAbrasionModel::G4WilsonAbrasionModel ( G4bool  useAblation1 = false)

Definition at line 115 of file G4WilsonAbrasionModel.cc.

116 :G4HadronicInteraction("G4WilsonAbrasion")
117{
118 // Send message to stdout to advise that the G4Abrasion model is being used.
119 PrintWelcomeMessage();
120
121 // Set the default verbose level to 0 - no output.
122 verboseLevel = 0;
123 useAblation = useAblation1;
124 theAblation = nullptr;
125
126 // No de-excitation handler has been supplied - define the default handler.
127
128 theExcitationHandler = new G4ExcitationHandler();
129 if (useAblation)
130 {
131 theAblation = new G4WilsonAblationModel;
132 theAblation->SetVerboseLevel(verboseLevel);
133 theExcitationHandler->SetEvaporation(theAblation, true);
134 }
135
136 // Set the minimum and maximum range for the model (despite nomanclature,
137 // this is in energy per nucleon number).
138
139 SetMinEnergy(70.0*MeV);
140 SetMaxEnergy(10.1*GeV);
141 isBlocked = false;
142
143 // npK, when mutiplied by the nuclear Fermi momentum, determines the range of
144 // momentum over which the secondary nucleon momentum is sampled.
145
146 r0sq = 0.0;
147 npK = 5.0;
148 B = 10.0 * MeV;
149 third = 1.0 / 3.0;
150 fradius = 0.99;
151 conserveEnergy = false;
152 conserveMomentum = true;
153}
void SetEvaporation(G4VEvaporation *ptr, G4bool isLocal=false)
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)

◆ G4WilsonAbrasionModel() [2/3]

G4WilsonAbrasionModel::G4WilsonAbrasionModel ( G4ExcitationHandler aExcitationHandler)

Definition at line 167 of file G4WilsonAbrasionModel.cc.

168{
169// Send message to stdout to advise that the G4Abrasion model is being used.
170
171 PrintWelcomeMessage();
172
173// Set the default verbose level to 0 - no output.
174
175 verboseLevel = 0;
176
177 theAblation = nullptr; //A.R. 26-Jul-2012 Coverity fix.
178 useAblation = false; //A.R. 14-Aug-2012 Coverity fix.
179
180//
181// The user is able to provide the excitation handler as well as an argument
182// which is provided in this instantiation is used to determine
183// whether the spectators of the interaction are free following the abrasion.
184//
185 theExcitationHandler = aExcitationHandler;
186//
187//
188// Set the minimum and maximum range for the model (despite nomanclature, this
189// is in energy per nucleon number).
190//
191 SetMinEnergy(70.0*MeV);
192 SetMaxEnergy(10.1*GeV);
193 isBlocked = false;
194//
195//
196// npK, when mutiplied by the nuclear Fermi momentum, determines the range of
197// momentum over which the secondary nucleon momentum is sampled.
198//
199 r0sq = 0.0; //A.R. 14-Aug-2012 Coverity fix.
200 npK = 5.0;
201 B = 10.0 * MeV;
202 third = 1.0 / 3.0;
203 fradius = 0.99;
204 conserveEnergy = false;
205 conserveMomentum = true;
206}

◆ ~G4WilsonAbrasionModel()

G4WilsonAbrasionModel::~G4WilsonAbrasionModel ( )

Definition at line 209 of file G4WilsonAbrasionModel.cc.

210{
211 delete theExcitationHandler;
212}

◆ G4WilsonAbrasionModel() [3/3]

G4WilsonAbrasionModel::G4WilsonAbrasionModel ( const G4WilsonAbrasionModel right)

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4WilsonAbrasionModel::ApplyYourself ( const G4HadProjectile theTrack,
G4Nucleus theTarget 
)
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 215 of file G4WilsonAbrasionModel.cc.

217{
218//
219//
220// The secondaries will be returned in G4HadFinalState &theParticleChange -
221// initialise this. The original track will always be discontinued and
222// secondaries followed.
223//
226//
227//
228// Get relevant information about the projectile and target (A, Z, energy/nuc,
229// momentum, etc).
230//
231 const G4ParticleDefinition *definitionP = theTrack.GetDefinition();
232 const G4double AP = definitionP->GetBaryonNumber();
233 const G4double ZP = definitionP->GetPDGCharge();
234 G4LorentzVector pP = theTrack.Get4Momentum();
235 G4double E = theTrack.GetKineticEnergy()/AP;
236 G4double AT = theTarget.GetA_asInt();
237 G4double ZT = theTarget.GetZ_asInt();
238 G4double TotalEPre = theTrack.GetTotalEnergy() +
239 theTarget.AtomicMass(AT, ZT) + theTarget.GetEnergyDeposit();
240 G4double TotalEPost = 0.0;
241//
242//
243// Determine the radii of the projectile and target nuclei.
244//
246 G4double rP = aR.GetWilsonRadius(AP);
247 G4double rT = aR.GetWilsonRadius(AT);
248 G4double rPsq = rP * rP;
249 G4double rTsq = rT * rT;
250 if (verboseLevel >= 2)
251 {
252 G4cout <<"########################################"
253 <<"########################################"
254 <<G4endl;
255 G4cout.precision(6);
256 G4cout <<"IN G4WilsonAbrasionModel" <<G4endl;
257 G4cout <<"Initial projectile A=" <<AP
258 <<", Z=" <<ZP
259 <<", radius = " <<rP/fermi <<" fm"
260 <<G4endl;
261 G4cout <<"Initial target A=" <<AT
262 <<", Z=" <<ZT
263 <<", radius = " <<rT/fermi <<" fm"
264 <<G4endl;
265 G4cout <<"Projectile momentum and Energy/nuc = " <<pP <<" ," <<E <<G4endl;
266 }
267//
268//
269// The following variables are used to determine the impact parameter in the
270// near-field (i.e. taking into consideration the electrostatic repulsion).
271//
272 G4double rm = ZP * ZT * elm_coupling / (E * AP);
273 G4double r = 0.0;
274 G4double rsq = 0.0;
275//
276//
277// Initialise some of the variables which wll be used to calculate the chord-
278// length for nucleons in the projectile and target, and hence calculate the
279// number of abraded nucleons and the excitation energy.
280//
281 G4NuclearAbrasionGeometry *theAbrasionGeometry = nullptr;
282 G4double CT = 0.0;
283 G4double F = 0.0;
284 G4int Dabr = 0;
285//
286//
287// The following loop is performed until the number of nucleons which are
288// abraded by the process is >1, i.e. an interaction MUST occur.
289//
290 G4bool skipInteraction = false; // It will be set true if the two nuclei fail to collide
291 const G4int maxNumberOfLoops = 1000;
292 G4int loopCounter = -1;
293 while (Dabr == 0 && ++loopCounter < maxNumberOfLoops) /* Loop checking, 07.08.2015, A.Ribon */
294 {
295//
296//
297// Sample the impact parameter. For the moment, this class takes account of
298// electrostatic effects on the impact parameter, but (like HZETRN AND NUCFRG2)
299// does not make any correction for the effects of nuclear-nuclear repulsion.
300//
301 G4double rPT = rP + rT;
302 G4double rPTsq = rPT * rPT;
303//
304//
305// This is a "catch" to make sure we don't go into an infinite loop because the
306// energy is too low to overcome nuclear repulsion. PRT 20091023. If the
307// value of rm < fradius * rPT then we're unlikely to sample a small enough
308// impact parameter (energy of incident particle is too low).
309//
310 if (rm >= fradius * rPT) {
311 skipInteraction = true;
312 }
313//
314//
315// Now sample impact parameter until the criterion is met that projectile
316// and target overlap, but repulsion is taken into consideration.
317//
318 G4int evtcnt = 0;
319 r = 1.1 * rPT;
320 while (r > rPT && ++evtcnt < 1000) /* Loop checking, 07.08.2015, A.Ribon */
321 {
322 G4double bsq = rPTsq * G4UniformRand();
323 r = (rm + std::sqrt(rm*rm + 4.0*bsq)) / 2.0;
324 }
325//
326//
327// We've tried to sample this 1000 times, but failed.
328//
329 if (evtcnt >= 1000) {
330 skipInteraction = true;
331 }
332
333 rsq = r * r;
334//
335//
336// Now determine the chord-length through the target nucleus.
337//
338 if (rT > rP)
339 {
340 G4double x = (rPsq + rsq - rTsq) / 2.0 / r;
341 if (x > 0.0) CT = 2.0 * std::sqrt(rTsq - x*x);
342 else CT = 2.0 * std::sqrt(rTsq - rsq);
343 }
344 else
345 {
346 G4double x = (rTsq + rsq - rPsq) / 2.0 / r;
347 if (x > 0.0) CT = 2.0 * std::sqrt(rTsq - x*x);
348 else CT = 2.0 * rT;
349 }
350//
351//
352// Determine the number of abraded nucleons. Note that the mean number of
353// abraded nucleons is used to sample the Poisson distribution. The Poisson
354// distribution is sampled only ten times with the current impact parameter,
355// and if it fails after this to find a case for which the number of abraded
356// nucleons >1, the impact parameter is re-sampled.
357//
358 delete theAbrasionGeometry;
359 theAbrasionGeometry = new G4NuclearAbrasionGeometry(AP,AT,r);
360 F = theAbrasionGeometry->F();
361 G4double lambda = 16.6*fermi / G4Pow::GetInstance()->powA(E/MeV,0.26);
362 G4double Mabr = F * AP * (1.0 - G4Exp(-CT/lambda));
363 G4long n = 0;
364 for (G4int i = 0; i<10; ++i)
365 {
366 n = G4Poisson(Mabr);
367 if (n > 0)
368 {
369 if (n>AP) Dabr = (G4int) AP;
370 else Dabr = (G4int) n;
371 break;
372 }
373 }
374 } // End of while loop
375
376 if ( loopCounter >= maxNumberOfLoops || skipInteraction ) {
377 // Assume nuclei do not collide and return unchanged primary.
381 if (verboseLevel >= 2) {
382 G4cout <<"Particle energy too low to overcome repulsion." <<G4endl;
383 G4cout <<"Event rejected and original track maintained" <<G4endl;
384 G4cout <<"########################################"
385 <<"########################################"
386 <<G4endl;
387 }
388 delete theAbrasionGeometry;
389 return &theParticleChange;
390 }
391
392 if (verboseLevel >= 2)
393 {
394 G4cout <<G4endl;
395 G4cout <<"Impact parameter = " <<r/fermi <<" fm" <<G4endl;
396 G4cout <<"# Abraded nucleons = " <<Dabr <<G4endl;
397 }
398//
399//
400// The number of abraded nucleons must be no greater than the number of
401// nucleons in either the projectile or the target. If AP - Dabr < 2 or
402// AT - Dabr < 2 then either we have only a nucleon left behind in the
403// projectile/target or we've tried to abrade too many nucleons - and Dabr
404// should be limited.
405//
406 if (AP - (G4double) Dabr < 2.0) Dabr = (G4int) AP;
407 if (AT - (G4double) Dabr < 2.0) Dabr = (G4int) AT;
408//
409//
410// Determine the abraded secondary nucleons from the projectile. *fragmentP
411// is a pointer to the prefragment from the projectile and nSecP is the number
412// of nucleons in theParticleChange which have been abraded. The total energy
413// from these is determined.
414//
415 G4ThreeVector boost = pP.findBoostToCM();
416 G4Fragment *fragmentP = GetAbradedNucleons (Dabr, AP, ZP, rP);
418 G4int i = 0;
419 for (i=0; i<nSecP; ++i)
420 {
421 TotalEPost += theParticleChange.GetSecondary(i)->
422 GetParticle()->GetTotalEnergy();
423 }
424//
425//
426// Determine the number of spectators in the interaction region for the
427// projectile.
428//
429 G4int DspcP = (G4int) (AP*F) - Dabr;
430 if (DspcP <= 0) DspcP = 0;
431 else if (DspcP > AP-Dabr) DspcP = ((G4int) AP) - Dabr;
432//
433//
434// Determine excitation energy associated with excess surface area of the
435// projectile (EsP) and the excitation due to scattering of nucleons which are
436// retained within the projectile (ExP). Add the total energy from the excited
437// nucleus to the total energy of the secondaries.
438//
439 G4bool excitationAbsorbedByProjectile = false;
440 if (fragmentP != nullptr)
441 {
442 G4double EsP = theAbrasionGeometry->GetExcitationEnergyOfProjectile();
443 G4double ExP = 0.0;
444 if (Dabr < AT)
445 excitationAbsorbedByProjectile = G4UniformRand() < 0.5;
446 if (excitationAbsorbedByProjectile)
447 ExP = GetNucleonInducedExcitation(rP, rT, r);
448 G4double xP = EsP + ExP;
449 if (xP > B*(AP-Dabr)) xP = B*(AP-Dabr);
450 G4LorentzVector lorentzVector = fragmentP->GetMomentum();
451 lorentzVector.setE(lorentzVector.e()+xP);
452 fragmentP->SetMomentum(lorentzVector);
453 TotalEPost += lorentzVector.e();
454 }
455 G4double EMassP = TotalEPost;
456//
457//
458// Determine the abraded secondary nucleons from the target. Note that it's
459// assumed that the same number of nucleons are abraded from the target as for
460// the projectile, and obviously no boost is applied to the products. *fragmentT
461// is a pointer to the prefragment from the target and nSec is the total number
462// of nucleons in theParticleChange which have been abraded. The total energy
463// from these is determined.
464//
465 G4Fragment *fragmentT = GetAbradedNucleons (Dabr, AT, ZT, rT);
467 for (i=nSecP; i<nSec; ++i)
468 {
469 TotalEPost += theParticleChange.GetSecondary(i)->
470 GetParticle()->GetTotalEnergy();
471 }
472//
473//
474// Determine the number of spectators in the interaction region for the
475// target.
476//
477 G4int DspcT = (G4int) (AT*F) - Dabr;
478 if (DspcT <= 0) DspcT = 0;
479 else if (DspcT > AP-Dabr) DspcT = ((G4int) AT) - Dabr;
480//
481//
482// Determine excitation energy associated with excess surface area of the
483// target (EsT) and the excitation due to scattering of nucleons which are
484// retained within the target (ExT). Add the total energy from the excited
485// nucleus to the total energy of the secondaries.
486//
487 if (fragmentT != nullptr)
488 {
489 G4double EsT = theAbrasionGeometry->GetExcitationEnergyOfTarget();
490 G4double ExT = 0.0;
491 if (!excitationAbsorbedByProjectile)
492 ExT = GetNucleonInducedExcitation(rT, rP, r);
493 G4double xT = EsT + ExT;
494 if (xT > B*(AT-Dabr)) xT = B*(AT-Dabr);
495 G4LorentzVector lorentzVector = fragmentT->GetMomentum();
496 lorentzVector.setE(lorentzVector.e()+xT);
497 fragmentT->SetMomentum(lorentzVector);
498 TotalEPost += lorentzVector.e();
499 }
500//
501//
502// Now determine the difference between the pre and post interaction
503// energy - this will be used to determine the Lorentz boost if conservation
504// of energy is to be imposed/attempted.
505//
506 G4double deltaE = TotalEPre - TotalEPost;
507 if (deltaE > 0.0 && conserveEnergy)
508 {
509 G4double beta = std::sqrt(1.0 - EMassP*EMassP/G4Pow::GetInstance()->powN(deltaE+EMassP,2));
510 boost = boost / boost.mag() * beta;
511 }
512//
513//
514// Now boost the secondaries from the projectile.
515//
516 G4ThreeVector pBalance = pP.vect();
517 for (i=0; i<nSecP; ++i)
518 {
520 GetParticle();
521 G4LorentzVector lorentzVector = dynamicP->Get4Momentum();
522 lorentzVector.boost(-boost);
523 dynamicP->Set4Momentum(lorentzVector);
524 pBalance -= lorentzVector.vect();
525 }
526//
527//
528// Set the boost for the projectile prefragment. This is now based on the
529// conservation of momentum. However, if the user selected momentum of the
530// prefragment is not to be conserved this simply boosted to the velocity of the
531// original projectile times the ratio of the unexcited to the excited mass
532// of the prefragment (the excitation increases the effective mass of the
533// prefragment, and therefore modifying the boost is an attempt to prevent
534// the momentum of the prefragment being excessive).
535//
536 if (fragmentP != nullptr)
537 {
538 G4LorentzVector lorentzVector = fragmentP->GetMomentum();
539 G4double fragmentM = lorentzVector.m();
540 if (conserveMomentum)
541 fragmentP->SetMomentum
542 (G4LorentzVector(pBalance,std::sqrt(pBalance.mag2()+fragmentM*fragmentM+1.0*eV*eV)));
543 else
544 {
545 G4double fragmentGroundStateM = fragmentP->GetGroundStateMass();
546 fragmentP->SetMomentum(lorentzVector.boost(-boost * fragmentGroundStateM/fragmentM));
547 }
548 }
549//
550//
551// Output information to user if verbose information requested.
552//
553 if (verboseLevel >= 2)
554 {
555 G4cout <<G4endl;
556 G4cout <<"-----------------------------------" <<G4endl;
557 G4cout <<"Secondary nucleons from projectile:" <<G4endl;
558 G4cout <<"-----------------------------------" <<G4endl;
559 G4cout.precision(7);
560 for (i=0; i<nSecP; ++i)
561 {
562 G4cout <<"Particle # " <<i <<G4endl;
565 G4cout <<"New nucleon (P) " <<dyn->GetDefinition()->GetParticleName()
566 <<" : " <<dyn->Get4Momentum()
567 <<G4endl;
568 }
569 G4cout <<"---------------------------" <<G4endl;
570 G4cout <<"The projectile prefragment:" <<G4endl;
571 G4cout <<"---------------------------" <<G4endl;
572 if (fragmentP != nullptr)
573 G4cout <<*fragmentP <<G4endl;
574 else
575 G4cout <<"(No residual prefragment)" <<G4endl;
576 G4cout <<G4endl;
577 G4cout <<"-------------------------------" <<G4endl;
578 G4cout <<"Secondary nucleons from target:" <<G4endl;
579 G4cout <<"-------------------------------" <<G4endl;
580 G4cout.precision(7);
581 for (i=nSecP; i<nSec; ++i)
582 {
583 G4cout <<"Particle # " <<i <<G4endl;
586 G4cout <<"New nucleon (T) " <<dyn->GetDefinition()->GetParticleName()
587 <<" : " <<dyn->Get4Momentum()
588 <<G4endl;
589 }
590 G4cout <<"-----------------------" <<G4endl;
591 G4cout <<"The target prefragment:" <<G4endl;
592 G4cout <<"-----------------------" <<G4endl;
593 if (fragmentT != nullptr)
594 G4cout <<*fragmentT <<G4endl;
595 else
596 G4cout <<"(No residual prefragment)" <<G4endl;
597 }
598//
599//
600// Now we can decay the nuclear fragments if present. The secondaries are
601// collected and boosted as well. This is performed first for the projectile...
602//
603 if (fragmentP !=nullptr)
604 {
605 G4ReactionProductVector *products = nullptr;
606 // if (fragmentP->GetZ_asInt() != fragmentP->GetA_asInt())
607 products = theExcitationHandler->BreakItUp(*fragmentP);
608 // else
609 // products = theExcitationHandlerx->BreakItUp(*fragmentP);
610 delete fragmentP;
611 fragmentP = nullptr;
612
613 G4ReactionProductVector::iterator iter;
614 for (iter = products->begin(); iter != products->end(); ++iter)
615 {
616 G4DynamicParticle *secondary =
617 new G4DynamicParticle((*iter)->GetDefinition(),
618 (*iter)->GetTotalEnergy(), (*iter)->GetMomentum());
619 theParticleChange.AddSecondary (secondary); // Added MHM 20050118
620 G4String particleName = (*iter)->GetDefinition()->GetParticleName();
621 delete (*iter); // get rid of leftover particle def! // Added MHM 20050118
622 if (verboseLevel >= 2 && particleName.find("[",0) < particleName.size())
623 {
624 G4cout <<"------------------------" <<G4endl;
625 G4cout <<"The projectile fragment:" <<G4endl;
626 G4cout <<"------------------------" <<G4endl;
627 G4cout <<" fragmentP = " <<particleName
628 <<" Energy = " <<secondary->GetKineticEnergy()
629 <<G4endl;
630 }
631 }
632 delete products; // Added MHM 20050118
633 }
634//
635//
636// Now decay the target nucleus - no boost is applied since in this
637// approximation it is assumed that there is negligible momentum transfer from
638// the projectile.
639//
640 if (fragmentT != nullptr)
641 {
642 G4ReactionProductVector *products = nullptr;
643 // if (fragmentT->GetZ_asInt() != fragmentT->GetA_asInt())
644 products = theExcitationHandler->BreakItUp(*fragmentT);
645 // else
646 // products = theExcitationHandlerx->BreakItUp(*fragmentT);
647 // delete fragmentT;
648 fragmentT = nullptr;
649
650 G4ReactionProductVector::iterator iter;
651 for (iter = products->begin(); iter != products->end(); ++iter)
652 {
653 G4DynamicParticle *secondary =
654 new G4DynamicParticle((*iter)->GetDefinition(),
655 (*iter)->GetTotalEnergy(), (*iter)->GetMomentum());
657 G4String particleName = (*iter)->GetDefinition()->GetParticleName();
658 delete (*iter); // get rid of leftover particle def! // Added MHM 20050118
659 if (verboseLevel >= 2 && particleName.find("[",0) < particleName.size())
660 {
661 G4cout <<"--------------------" <<G4endl;
662 G4cout <<"The target fragment:" <<G4endl;
663 G4cout <<"--------------------" <<G4endl;
664 G4cout <<" fragmentT = " <<particleName
665 <<" Energy = " <<secondary->GetKineticEnergy()
666 <<G4endl;
667 }
668 }
669 delete products; // Added MHM 20050118
670 }
671
672 if (verboseLevel >= 2)
673 G4cout <<"########################################"
674 <<"########################################"
675 <<G4endl;
676
677 delete theAbrasionGeometry;
678 return &theParticleChange;
679}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
@ isAlive
@ stopAndKill
CLHEP::HepLorentzVector G4LorentzVector
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:50
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:83
long G4long
Definition: G4Types.hh:87
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
double mag2() const
double mag() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
Hep3Vector findBoostToCM() const
void DumpInfo(G4int mode=0) const
G4ParticleDefinition * GetDefinition() const
G4LorentzVector Get4Momentum() const
G4double GetKineticEnergy() const
void Set4Momentum(const G4LorentzVector &momentum)
G4ReactionProductVector * BreakItUp(const G4Fragment &theInitialState)
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:280
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:299
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:304
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
std::size_t GetNumberOfSecondaries() const
void SetEnergyChange(G4double anEnergy)
G4HadSecondary * GetSecondary(size_t i)
void SetMomentumChange(const G4ThreeVector &aV)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetTotalEnergy() const
G4DynamicParticle * GetParticle()
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4double GetEnergyDeposit()
Definition: G4Nucleus.hh:184
G4double AtomicMass(const G4double A, const G4double Z) const
Definition: G4Nucleus.cc:254
G4double GetPDGCharge() const
const G4String & GetParticleName() const
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:230
G4double GetWilsonRadius(G4double A)

◆ GetConserveMomentum()

G4bool G4WilsonAbrasionModel::GetConserveMomentum ( )
inline

Definition at line 138 of file G4WilsonAbrasionModel.hh.

139 {return conserveMomentum;}

◆ GetExcitationHandler()

G4ExcitationHandler * G4WilsonAbrasionModel::GetExcitationHandler ( )
inline

Definition at line 123 of file G4WilsonAbrasionModel.hh.

124 {return theExcitationHandler;}

◆ GetUseAblation()

G4bool G4WilsonAbrasionModel::GetUseAblation ( )
inline

Definition at line 126 of file G4WilsonAbrasionModel.hh.

127 {return useAblation;}

◆ ModelDescription()

void G4WilsonAbrasionModel::ModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 155 of file G4WilsonAbrasionModel.cc.

156{
157 outFile << "G4WilsonAbrasionModel is a macroscopic treatment of\n"
158 << "nucleus-nucleus collisions using simple geometric arguments.\n"
159 << "The smaller projectile nucleus gouges out a part of the larger\n"
160 << "target nucleus, leaving a residual nucleus and a fireball\n"
161 << "region where the projectile and target intersect. The fireball"
162 << "is then treated as a highly excited nuclear fragment. This\n"
163 << "model is based on the NUCFRG2 model and is valid for all\n"
164 << "projectile energies between 70 MeV/n and 10.1 GeV/n. \n";
165}

◆ operator=()

const G4WilsonAbrasionModel & G4WilsonAbrasionModel::operator= ( G4WilsonAbrasionModel right)

◆ SetConserveMomentum()

void G4WilsonAbrasionModel::SetConserveMomentum ( G4bool  conserveMomentum1)
inline

Definition at line 135 of file G4WilsonAbrasionModel.hh.

136 {conserveMomentum = conserveMomentum1;}

◆ SetExcitationHandler()

void G4WilsonAbrasionModel::SetExcitationHandler ( G4ExcitationHandler aExcitationHandler)
inline

Definition at line 120 of file G4WilsonAbrasionModel.hh.

121 {theExcitationHandler = aExcitationHandler;}

◆ SetUseAblation()

void G4WilsonAbrasionModel::SetUseAblation ( G4bool  useAblation1)

Definition at line 848 of file G4WilsonAbrasionModel.cc.

849{
850 if (useAblation != useAblation1)
851 {
852 useAblation = useAblation1;
853 if (useAblation)
854 {
855 theAblation = new G4WilsonAblationModel;
856 theAblation->SetVerboseLevel(verboseLevel);
857 theExcitationHandler->SetEvaporation(theAblation);
858 }
859 else
860 {
861 delete theExcitationHandler;
862 theAblation = nullptr;
863 theExcitationHandler = new G4ExcitationHandler();
864 }
865 }
866 return;
867}

◆ SetVerboseLevel()

void G4WilsonAbrasionModel::SetVerboseLevel ( G4int  verboseLevel1)
inline

Definition at line 141 of file G4WilsonAbrasionModel.hh.

142{
143 verboseLevel = verboseLevel1;
144 if (useAblation) theAblation->SetVerboseLevel(verboseLevel);
145}

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