Geant4 11.1.1
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 &)
 
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 ()
 
virtual G4ReactionProductVectorDeExcite (G4Fragment &aFragment)=0
 
virtual void DeExciteModelDescription (std::ostream &outFile) const =0
 
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 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 49 of file G4AblaInterface.hh.

Constructor & Destructor Documentation

◆ G4AblaInterface()

G4AblaInterface::G4AblaInterface ( G4ExcitationHandler ptr = nullptr)

Definition at line 55 of file G4AblaInterface.cc.

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

◆ ~G4AblaInterface()

G4AblaInterface::~G4AblaInterface ( )
virtual

Definition at line 71 of file G4AblaInterface.cc.

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

Member Function Documentation

◆ ApplyYourself()

virtual G4HadFinalState * G4AblaInterface::ApplyYourself ( G4HadProjectile const &  ,
G4Nucleus  
)
inlinevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 57 of file G4AblaInterface.hh.

58 {
59return nullptr;
60 }

◆ 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 93 of file G4AblaInterface.cc.

93 {
94 if (!isInitialised) InitialiseModel();
95
96 volant->clear();
97 ablaResult->clear();
98
99 const G4int ARem = aFragment.GetA_asInt();
100 const G4int ZRem = aFragment.GetZ_asInt();
101 const G4int SRem = -aFragment.GetNumberOfLambdas(); // Strangeness = - (Number of lambdas)
102 const G4double eStarRem = aFragment.GetExcitationEnergy() / MeV;
103 const G4double jRem = aFragment.GetAngularMomentum().mag() / hbar_Planck;
104 const G4LorentzVector& pRem = aFragment.GetMomentum();
105 const G4double pxRem = pRem.x() / MeV;
106 const G4double pyRem = pRem.y() / MeV;
107 const G4double pzRem = pRem.z() / MeV;
108
109 ++eventNumber;
110
111 theABLAModel->DeexcitationAblaxx(ARem, ZRem, eStarRem, jRem, pxRem, pyRem,
112 pzRem, (G4int)eventNumber, SRem);
113
115
116 for(G4int j = 0; j < ablaResult->ntrack; ++j)
117 { // Copy ABLA result to the EventInfo
118 G4ReactionProduct* product =
119 toG4Particle(ablaResult->avv[j], ablaResult->zvv[j], ablaResult->svv[j],
120 ablaResult->enerj[j], ablaResult->pxlab[j],
121 ablaResult->pylab[j], ablaResult->pzlab[j]);
122 if(product)
123 {
124 product->SetCreatorModelID(secID);
125 result->push_back(product);
126 }
127 }
128 return result;
129}
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
double mag() const
void DeexcitationAblaxx(G4int nucleusA, G4int nucleusZ, G4double excitationEnergy, G4double angularMomentum, G4double momX, G4double momY, G4double momZ, G4int eventnumber)
Definition: G4Abla.cc:96
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:312
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:322
G4int GetZ_asInt() const
Definition: G4Fragment.hh:289
G4int GetNumberOfLambdas() const
Definition: G4Fragment.hh:307
G4ThreeVector GetAngularMomentum() const
Definition: G4Fragment.cc:330
G4int GetA_asInt() const
Definition: G4Fragment.hh:284
void SetCreatorModelID(const G4int mod)
void clear()
G4double enerj[VARNTPSIZE]
G4double pylab[VARNTPSIZE]
G4int svv[VARNTPSIZE]
G4int avv[VARNTPSIZE]
G4double pzlab[VARNTPSIZE]
G4int zvv[VARNTPSIZE]
G4double pxlab[VARNTPSIZE]
void clear()

◆ DeExciteModelDescription()

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

Implements G4VPreCompoundModel.

Definition at line 192 of file G4AblaInterface.cc.

193{
194 outFile
195 << "ABLA++ is a statistical model for nuclear de-excitation. It simulates\n"
196 << "the gamma emission and the evaporation of neutrons, light charged\n"
197 << "particles and IMFs, as well as fission where applicable. The code\n"
198 << "included in Geant4 is a C++ translation of the original Fortran\n"
199 << "code ABLA07. Although the model has been recently extended to\n"
200 << "hypernuclei by including the evaporation of lambda particles.\n"
201 << "More details about the physics are available in the Geant4\n"
202 << "Physics Reference Manual and in the reference articles.\n\n"
203 << "References:\n"
204 << "(1) A. Kelic, M. V. Ricciardi, and K. H. Schmidt, in Proceedings of "
205 "Joint\n"
206 << "ICTP-IAEA Advanced Workshop on Model Codes for Spallation Reactions,\n"
207 << "ICTP Trieste, Italy, 4–8 February 2008, edited by D. Filges, S. Leray, "
208 "Y. Yariv,\n"
209 << "A. Mengoni, A. Stanculescu, and G. Mank (IAEA INDC(NDS)-530, Vienna, "
210 "2008), pp. 181–221.\n\n"
211 << "(2) J.L. Rodriguez-Sanchez, J.-C. David et al., Phys. Rev. C 98, "
212 "021602 (2018)\n\n";
213}

◆ 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}
void initEvapora()
Definition: G4Abla.cc:2133
void SetParameters()
Definition: G4Abla.cc:2320

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

◆ ModelDescription()

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

Reimplemented from G4HadronicInteraction.

Definition at line 187 of file G4AblaInterface.cc.

188{
189 outFile << "ABLA++ does not provide an implementation of the ApplyYourself method!\n\n";
190}

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