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

#include <G4AblaInterface.hh>

+ Inheritance diagram for G4AblaInterface:

Public Member Functions

 G4AblaInterface (G4ExcitationHandler *ptr=nullptr)
 
virtual ~G4AblaInterface ()
 
virtual G4ReactionProductVectorDeExcite (G4Fragment &aFragment)
 
virtual G4HadFinalStateApplyYourself (G4HadProjectile const &, G4Nucleus &) final
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) final
 
virtual void InitialiseModel () final
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void DeExciteModelDescription (std::ostream &outFile) const
 
- Public Member Functions inherited from G4VPreCompoundModel
 G4VPreCompoundModel (G4ExcitationHandler *ptr=nullptr, const G4String &modelName="PrecompoundModel")
 
virtual ~G4VPreCompoundModel ()
 
void SetExcitationHandler (G4ExcitationHandler *ptr)
 
G4ExcitationHandlerGetExcitationHandler () const
 
 G4VPreCompoundModel (const G4VPreCompoundModel &)=delete
 
const G4VPreCompoundModeloperator= (const G4VPreCompoundModel &right)=delete
 
G4bool operator== (const G4VPreCompoundModel &right) const =delete
 
G4bool operator!= (const G4VPreCompoundModel &right) const =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)
 
 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 49 of file G4AblaInterface.hh.

Constructor & Destructor Documentation

◆ G4AblaInterface()

G4AblaInterface::G4AblaInterface ( G4ExcitationHandler * ptr = nullptr)

Definition at line 55 of file G4AblaInterface.cc.

56 : G4VPreCompoundModel(ptr, "ABLAXX")
57 , ablaResult(new G4VarNtp)
58 , theABLAModel(new G4Abla(ablaResult))
59 , eventNumber(0)
60 , secID(-1)
61 , isInitialised(false)
62{
64 // G4cout << "### NEW PrecompoundModel " << this << G4endl;
65 if (!ptr)
66 SetExcitationHandler(new G4ExcitationHandler);
68 G4cout << G4endl << "G4AblaInterface::InitialiseModel() was right." << G4endl;
69}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
virtual void InitialiseModel() final
const G4String & GetModelName() const
static G4int GetModelID(const G4int modelIndex)
G4VPreCompoundModel(G4ExcitationHandler *ptr=nullptr, const G4String &modelName="PrecompoundModel")
void SetExcitationHandler(G4ExcitationHandler *ptr)

◆ ~G4AblaInterface()

G4AblaInterface::~G4AblaInterface ( )
virtual

Definition at line 71 of file G4AblaInterface.cc.

72{
73 applyYourselfResult.Clear();
74 delete ablaResult;
75 delete theABLAModel;
76 delete GetExcitationHandler();
77}
G4ExcitationHandler * GetExcitationHandler() const

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4AblaInterface::ApplyYourself ( G4HadProjectile const & thePrimary,
G4Nucleus & theNucleus )
finalvirtual

Reimplemented from G4HadronicInteraction.

Definition at line 91 of file G4AblaInterface.cc.

92{
93 // This method is adapted from G4PreCompoundModel::ApplyYourself,
94 // and it is used only by Binary Cascade (BIC) when the latter is coupled with
95 // Abla for nuclear de-excitation. This method allows BIC+ABLA to be used also
96 // for proton and neutron projectile with kinetic energies below 45 MeV, by
97 // creating a "compound" nucleus made by the system "target nucleus +
98 // projectile", before calling the DeExcite method.
99 const G4ParticleDefinition* primary = thePrimary.GetDefinition();
100 if (primary != G4Neutron::Definition() && primary != G4Proton::Definition())
101 {
103 ed << "G4AblaModel is used for ";
104 if (primary)
105 ed << primary->GetParticleName();
106 G4Exception("G4AblaInterface::ApplyYourself()", "had040", FatalException, ed, "");
107 return nullptr;
108 }
109
110 G4int Zp = 0;
111 G4int Ap = 1;
112 if (primary == G4Proton::Definition())
113 Zp = 1;
114 G4double timePrimary = thePrimary.GetGlobalTime();
115 G4int A = theNucleus.GetA_asInt();
116 G4int Z = theNucleus.GetZ_asInt();
117 G4LorentzVector p = thePrimary.Get4Momentum();
119 p += G4LorentzVector(0.0, 0.0, 0.0, mass);
120
121 G4Fragment anInitialState(A + Ap, Z + Zp, p);
122 anInitialState.SetNumberOfExcitedParticle(1, Zp);
123 anInitialState.SetNumberOfHoles(1, Zp);
124 anInitialState.SetCreationTime(thePrimary.GetGlobalTime());
125 anInitialState.SetCreatorModelID(secID);
126
127 G4ReactionProductVector* deExciteResult = DeExcite(anInitialState);
128
129 applyYourselfResult.Clear();
130 applyYourselfResult.SetStatusChange(stopAndKill);
131 for (auto const& prod : *deExciteResult)
132 {
133 G4DynamicParticle* aNewDP =
134 new G4DynamicParticle(prod->GetDefinition(), prod->GetTotalEnergy(), prod->GetMomentum());
135 G4HadSecondary aNew = G4HadSecondary(aNewDP);
136 G4double time = std::max(prod->GetFormationTime(), 0.0);
137 aNew.SetTime(timePrimary + time);
138 aNew.SetCreatorModelID(prod->GetCreatorModelID());
139 delete prod;
140 applyYourselfResult.AddSecondary(aNew);
141 }
142 delete deExciteResult;
143 return &applyYourselfResult;
144}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
@ stopAndKill
CLHEP::HepLorentzVector G4LorentzVector
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double A[17]
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment)
void SetTime(G4double aT)
void SetCreatorModelID(G4int id)
static G4Neutron * Definition()
Definition G4Neutron.cc:50
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetA_asInt() const
Definition G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition G4Nucleus.hh:105
const G4String & GetParticleName() const
static G4Proton * Definition()
Definition G4Proton.cc:45

◆ BuildPhysicsTable()

void G4AblaInterface::BuildPhysicsTable ( const G4ParticleDefinition & )
finalvirtual

Reimplemented from G4HadronicInteraction.

Definition at line 79 of file G4AblaInterface.cc.

79{ InitialiseModel(); }

◆ DeExcite()

G4ReactionProductVector * G4AblaInterface::DeExcite ( G4Fragment & aFragment)
virtual

Implements G4VPreCompoundModel.

Definition at line 146 of file G4AblaInterface.cc.

147{
148 if (!isInitialised)
150
151 ablaResult->clear();
152
153 const G4int ARem = aFragment.GetA_asInt();
154 const G4int ZRem = aFragment.GetZ_asInt();
155 const G4int SRem = -aFragment.GetNumberOfLambdas(); // Strangeness = - (Number of lambdas)
156 const G4double eStarRem = aFragment.GetExcitationEnergy() / MeV;
157 const G4double jRem = aFragment.GetAngularMomentum().mag() / hbar_Planck;
158 const G4LorentzVector& pRem = aFragment.GetMomentum();
159 const G4double pxRem = pRem.x() / MeV;
160 const G4double pyRem = pRem.y() / MeV;
161 const G4double pzRem = pRem.z() / MeV;
162
163 ++eventNumber;
164
165 theABLAModel->DeexcitationAblaxx(ARem, ZRem, eStarRem, jRem, pxRem, pyRem, pzRem, (G4int)eventNumber, SRem);
166
168
169 for (G4int j = 0; j < ablaResult->ntrack; ++j)
170 { // Copy ABLA result to the EventInfo
171 G4ReactionProduct* product = toG4Particle(ablaResult->avv[j],
172 ablaResult->zvv[j],
173 ablaResult->svv[j],
174 ablaResult->enerj[j],
175 ablaResult->pxlab[j],
176 ablaResult->pylab[j],
177 ablaResult->pzlab[j]);
178 if (product)
179 {
180 product->SetCreatorModelID(secID);
181 result->push_back(product);
182 }
183 }
184 return result;
185}
double mag() const
G4double GetExcitationEnergy() const
const G4LorentzVector & GetMomentum() const
G4int GetZ_asInt() const
G4int GetNumberOfLambdas() const
G4ThreeVector GetAngularMomentum() const
G4int GetA_asInt() const
void SetCreatorModelID(const G4int mod)

Referenced by ApplyYourself().

◆ DeExciteModelDescription()

void G4AblaInterface::DeExciteModelDescription ( std::ostream & outFile) const
virtual

Implements G4VPreCompoundModel.

Definition at line 267 of file G4AblaInterface.cc.

268{
269 outFile << "ABLA++ is a statistical model for nuclear de-excitation. It simulates\n"
270 << "the gamma emission and the evaporation of neutrons, light charged\n"
271 << "particles and IMFs, as well as fission where applicable. The code\n"
272 << "included in Geant4 is a C++ translation of the original Fortran\n"
273 << "code ABLA07. Although the model has been recently extended to\n"
274 << "hypernuclei by including the evaporation of lambda particles.\n"
275 << "More details about the physics are available in the Geant4\n"
276 << "Physics Reference Manual and in the reference articles.\n\n"
277 << "References:\n"
278 << "(1) A. Kelic, M. V. Ricciardi, and K. H. Schmidt, in Proceedings of Joint\n"
279 << "ICTP-IAEA Advanced Workshop on Model Codes for Spallation Reactions,\n"
280 << "ICTP Trieste, Italy, 4–8 February 2008, edited by D. Filges, S. "
281 "Leray, Y. Yariv, A. Mengoni, A. Stanculescu, and G. Mank (IAEA "
282 "INDC(NDS)-530, Vienna, 2008), pp. 181–221.\n\n"
283 << "(2) J.L. Rodriguez-Sanchez, J.-C. David et al., Phys. Rev. C 98, 021602R (2018)\n"
284 << "(3) J.L. Rodriguez-Sanchez et al., Phys. Rev. C 105, 014623 (2022)\n"
285 << "(4) J.L. Rodriguez-Sanchez et al., Phys. Rev. Lett. 130, 132501 (2023)\n\n";
286}

◆ InitialiseModel()

void G4AblaInterface::InitialiseModel ( )
finalvirtual

Reimplemented from G4HadronicInteraction.

Definition at line 81 of file G4AblaInterface.cc.

82{
83 if (isInitialised)
84 return;
85 isInitialised = true;
86 theABLAModel->initEvapora();
87 theABLAModel->SetParameters();
89}

Referenced by BuildPhysicsTable(), DeExcite(), and G4AblaInterface().

◆ ModelDescription()

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

Reimplemented from G4HadronicInteraction.

Definition at line 261 of file G4AblaInterface.cc.

262{
263 outFile << "ABLA++ does not provide an implementation of the ApplyYourself "
264 "method!\n\n";
265}

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