Geant4 11.4.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 48 of file G4AblaInterface.hh.

Constructor & Destructor Documentation

◆ G4AblaInterface()

G4AblaInterface::G4AblaInterface ( G4ExcitationHandler * ptr = nullptr)

Definition at line 56 of file G4AblaInterface.cc.

57 : G4VPreCompoundModel(ptr, "ABLAXX"),
58 ablaResult(new G4VarNtp),
59 theABLAModel(new G4Abla(ablaResult)),
60 eventNumber(0),
61 secID(-1),
62 isInitialised(false)
63{
65 // G4cout << "### NEW PrecompoundModel " << this << G4endl;
66 if (!ptr) 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 93 of file G4AblaInterface.cc.

95{
96 // This method is adapted from G4PreCompoundModel::ApplyYourself,
97 // and it is used only by Binary Cascade (BIC) when the latter is coupled with
98 // Abla for nuclear de-excitation. This method allows BIC+ABLA to be used also
99 // for proton and neutron projectile with kinetic energies below 45 MeV, by
100 // creating a "compound" nucleus made by the system "target nucleus +
101 // projectile", before calling the DeExcite method.
102 const G4ParticleDefinition* primary = thePrimary.GetDefinition();
103 if (primary != G4Neutron::Definition() && primary != G4Proton::Definition()) {
105 ed << "G4AblaModel is used for ";
106 if (primary) ed << primary->GetParticleName();
107 G4Exception("G4AblaInterface::ApplyYourself()", "had040", FatalException, ed, "");
108 return nullptr;
109 }
110
111 G4int Zp = 0;
112 G4int Ap = 1;
113 if (primary == G4Proton::Definition()) 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 G4DynamicParticle* aNewDP =
133 new G4DynamicParticle(prod->GetDefinition(), prod->GetTotalEnergy(), prod->GetMomentum());
134 G4HadSecondary aNew = G4HadSecondary(aNewDP);
135 G4double time = std::max(prod->GetFormationTime(), 0.0);
136 aNew.SetTime(timePrimary + time);
137 aNew.SetCreatorModelID(prod->GetCreatorModelID());
138 delete prod;
139 applyYourselfResult.AddSecondary(aNew);
140 }
141 delete deExciteResult;
142 return &applyYourselfResult;
143}
@ 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:78
G4int GetZ_asInt() const
Definition G4Nucleus.hh:84
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.

80{
82}

◆ DeExcite()

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

Implements G4VPreCompoundModel.

Definition at line 145 of file G4AblaInterface.cc.

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

Referenced by ApplyYourself().

◆ DeExciteModelDescription()

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

Implements G4VPreCompoundModel.

Definition at line 257 of file G4AblaInterface.cc.

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

◆ InitialiseModel()

void G4AblaInterface::InitialiseModel ( )
finalvirtual

Reimplemented from G4HadronicInteraction.

Definition at line 84 of file G4AblaInterface.cc.

85{
86 if (isInitialised) return;
87 isInitialised = true;
88 theABLAModel->initEvapora();
89 theABLAModel->SetParameters();
91}

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

◆ ModelDescription()

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

Reimplemented from G4HadronicInteraction.

Definition at line 251 of file G4AblaInterface.cc.

252{
253 outFile << "ABLA++ does not provide an implementation of the ApplyYourself "
254 "method!\n\n";
255}

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