Geant4 11.2.2
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 50 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: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 volant;
75 delete ablaResult;
76 delete theABLAModel;
77 delete GetExcitationHandler();
78}
G4ExcitationHandler * GetExcitationHandler() const

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 95 of file G4AblaInterface.cc.

97{
98 // This method is adapted from G4PreCompoundModel::ApplyYourself,
99 // and it is used only by Binary Cascade (BIC) when the latter is coupled with Abla
100 // for nuclear de-excitation.
101 // This method allows BIC+ABLA to be used also for proton and neutron projectile
102 // with kinetic energies below 45 MeV, by creating a "compound" nucleus made
103 // by the system "target nucleus + projectile", before calling the DeExcite
104 // method.
105 const G4ParticleDefinition* primary = thePrimary.GetDefinition();
106 if ( primary != G4Neutron::Definition() && primary != G4Proton::Definition() ) {
108 ed << "G4AblaModel is used for ";
109 if ( primary ) ed << primary->GetParticleName();
110 G4Exception( "G4AblaInterface::ApplyYourself()", "had040", FatalException, ed, "" );
111 return nullptr;
112 }
113
114 G4int Zp = 0;
115 G4int Ap = 1;
116 if ( primary == G4Proton::Definition() ) Zp = 1;
117 G4double timePrimary = thePrimary.GetGlobalTime();
118 G4int A = theNucleus.GetA_asInt();
119 G4int Z = theNucleus.GetZ_asInt();
120 G4LorentzVector p = thePrimary.Get4Momentum();
122 p += G4LorentzVector( 0.0, 0.0, 0.0, mass );
123
124 G4Fragment anInitialState(A + Ap, Z + Zp, p);
125 anInitialState.SetNumberOfExcitedParticle(1, Zp);
126 anInitialState.SetNumberOfHoles(1, Zp);
127 anInitialState.SetCreationTime( thePrimary.GetGlobalTime() );
128 anInitialState.SetCreatorModelID( secID );
129
130 G4ReactionProductVector* deExciteResult = DeExcite( anInitialState );
131
132 applyYourselfResult.Clear();
133 applyYourselfResult.SetStatusChange( stopAndKill );
134 for ( auto const & prod : *deExciteResult ) {
135 G4DynamicParticle * aNewDP =
136 new G4DynamicParticle( prod->GetDefinition(), prod->GetTotalEnergy(), prod->GetMomentum() );
137 G4HadSecondary aNew = G4HadSecondary( aNewDP );
138 G4double time = std::max( prod->GetFormationTime(), 0.0 );
139 aNew.SetTime( timePrimary + time );
140 aNew.SetCreatorModelID( prod->GetCreatorModelID() );
141 delete prod;
142 applyYourselfResult.AddSecondary( aNew );
143 }
144 delete deExciteResult;
145 return &applyYourselfResult;
146}
@ 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 SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
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 80 of file G4AblaInterface.cc.

81{
83}

◆ DeExcite()

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

Implements G4VPreCompoundModel.

Definition at line 149 of file G4AblaInterface.cc.

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

Referenced by ApplyYourself().

◆ DeExciteModelDescription()

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

Implements G4VPreCompoundModel.

Definition at line 248 of file G4AblaInterface.cc.

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

◆ InitialiseModel()

void G4AblaInterface::InitialiseModel ( )
finalvirtual

Reimplemented from G4HadronicInteraction.

Definition at line 85 of file G4AblaInterface.cc.

86{
87 if (isInitialised) return;
88 isInitialised = true;
89 theABLAModel->initEvapora();
90 theABLAModel->SetParameters();
92}
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 243 of file G4AblaInterface.cc.

244{
245 outFile << "ABLA++ does not provide an implementation of the ApplyYourself method!\n\n";
246}

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