Geant4 10.7.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 ()
 
virtual ~G4AblaInterface ()
 
virtual G4ReactionProductVectorDeExcite (G4Fragment &aFragment)
 
virtual G4HadFinalStateApplyYourself (G4HadProjectile const &, G4Nucleus &)
 
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 50 of file G4AblaInterface.hh.

Constructor & Destructor Documentation

◆ G4AblaInterface()

G4AblaInterface::G4AblaInterface ( )

Definition at line 51 of file G4AblaInterface.cc.

51 :
52 G4VPreCompoundModel(NULL, "ABLA"),
53 ablaResult(new G4VarNtp),
54 volant(new G4Volant),
55 theABLAModel(new G4Abla(volant, ablaResult)),
56 eventNumber(0)
57{
58 theABLAModel->initEvapora();
59 theABLAModel->SetParameters();
60}
Definition: G4Abla.hh:54
void initEvapora()
Definition: G4Abla.cc:2133
void SetParameters()
Definition: G4Abla.cc:2320

◆ ~G4AblaInterface()

G4AblaInterface::~G4AblaInterface ( )
virtual

Definition at line 62 of file G4AblaInterface.cc.

62 {
63 delete volant;
64 delete ablaResult;
65 delete theABLAModel;
66}

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 57 of file G4AblaInterface.hh.

57 {
58 return NULL;
59 }

◆ DeExcite()

G4ReactionProductVector * G4AblaInterface::DeExcite ( G4Fragment aFragment)
virtual

Implements G4VPreCompoundModel.

Definition at line 68 of file G4AblaInterface.cc.

68 {
69 volant->clear();
70 ablaResult->clear();
71
72 const G4int ARem = aFragment.GetA_asInt();
73 const G4int ZRem = aFragment.GetZ_asInt();
74 const G4double eStarRem = aFragment.GetExcitationEnergy() / MeV;
75 const G4double jRem = aFragment.GetAngularMomentum().mag() / hbar_Planck;
76 const G4LorentzVector &pRem = aFragment.GetMomentum();
77 const G4double pxRem = pRem.x() / MeV;
78 const G4double pyRem = pRem.y() / MeV;
79 const G4double pzRem = pRem.z() / MeV;
80
81 eventNumber++;
82
83 theABLAModel->DeexcitationAblaxx(ARem, ZRem,
84 eStarRem,
85 jRem,
86 pxRem,
87 pyRem,
88 pzRem,
89 eventNumber);
90
92
93 for(int j = 0; j < ablaResult->ntrack; ++j) { // Copy ABLA result to the EventInfo
94
95 G4ReactionProduct *product = toG4Particle(ablaResult->avv[j],
96 ablaResult->zvv[j],
97 ablaResult->svv[j],
98 ablaResult->enerj[j],
99 ablaResult->pxlab[j],
100 ablaResult->pylab[j],
101 ablaResult->pzlab[j]);
102 if(product)
103 result->push_back(product);
104 }
105 return result;
106}
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:275
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:299
G4int GetZ_asInt() const
Definition: G4Fragment.hh:263
G4ThreeVector GetAngularMomentum() const
Definition: G4Fragment.cc:254
G4int GetA_asInt() const
Definition: G4Fragment.hh:258
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 150 of file G4AblaInterface.cc.

150 {
151 outFile
152 << "ABLA++ is a statistical model for nuclear de-excitation. It simulates\n"
153 << "the gamma emission and the evaporation of neutrons, light charged particles\n"
154 << "and IMFs, as well as fission where applicable. The code included in Geant4\n"
155 << "is a C++ translation of the original Fortran code ABLA07. Although the model\n"
156 << "has been recently extended to hypernuclei by including the evaporation of lambda\n"
157 << "particles. More details about the physics are available in the\n"
158 << "Geant4 Physics Reference Manual and in the reference articles.\n\n"
159 << "References:\n"
160 << "(1) A. Kelic, M. V. Ricciardi, and K. H. Schmidt, in Proceedings of Joint\n"
161 << "ICTP-IAEA Advanced Workshop on Model Codes for Spallation Reactions,\n"
162 << "ICTP Trieste, Italy, 4–8 February 2008, edited by D. Filges, S. Leray, Y. Yariv,\n"
163 << "A. Mengoni, A. Stanculescu, and G. Mank (IAEA INDC(NDS)-530, Vienna, 2008), pp. 181–221.\n\n"
164 << "(2) J.L. Rodriguez-Sanchez, J.-C. David et al., Phys. Rev. C 98, 021602 (2018)\n\n";
165}

◆ ModelDescription()

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

Reimplemented from G4HadronicInteraction.

Definition at line 146 of file G4AblaInterface.cc.

146 {
147 outFile << "ABLA++ does not provide an implementation of the ApplyYourself method!\n\n";
148}

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