Geant4 11.2.2
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)=delete
 
const G4WilsonAbrasionModeloperator= (G4WilsonAbrasionModel &right)=delete
 
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 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 77 of file G4WilsonAbrasionModel.hh.

Constructor & Destructor Documentation

◆ G4WilsonAbrasionModel() [1/3]

G4WilsonAbrasionModel::G4WilsonAbrasionModel ( G4bool useAblation1 = false)

Definition at line 117 of file G4WilsonAbrasionModel.cc.

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

◆ G4WilsonAbrasionModel() [2/3]

G4WilsonAbrasionModel::G4WilsonAbrasionModel ( G4ExcitationHandler * aExcitationHandler)

Definition at line 172 of file G4WilsonAbrasionModel.cc.

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

◆ ~G4WilsonAbrasionModel()

G4WilsonAbrasionModel::~G4WilsonAbrasionModel ( )

Definition at line 218 of file G4WilsonAbrasionModel.cc.

219{
220 delete theExcitationHandler;
221}

◆ G4WilsonAbrasionModel() [3/3]

G4WilsonAbrasionModel::G4WilsonAbrasionModel ( const G4WilsonAbrasionModel & right)
delete

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 224 of file G4WilsonAbrasionModel.cc.

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

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

◆ operator=()

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

◆ 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 857 of file G4WilsonAbrasionModel.cc.

858{
859 if (useAblation != useAblation1)
860 {
861 useAblation = useAblation1;
862 if (useAblation)
863 {
864 theAblation = new G4WilsonAblationModel;
865 theAblation->SetVerboseLevel(verboseLevel);
866 theExcitationHandler->SetEvaporation(theAblation);
867 }
868 else
869 {
870 delete theExcitationHandler;
871 theAblation = nullptr;
872 theExcitationHandler = new G4ExcitationHandler();
873 }
874 }
875 return;
876}

◆ 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: