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

#include <G4PreCompoundModel.hh>

Inheritance diagram for G4PreCompoundModel:

Public Member Functions

 G4PreCompoundModel (G4ExcitationHandler *ptr=nullptr)
 ~G4PreCompoundModel () override
G4HadFinalStateApplyYourself (const G4HadProjectile &thePrimary, G4Nucleus &theNucleus) override
G4ReactionProductVectorDeExcite (G4Fragment &aFragment) override
void BuildPhysicsTable (const G4ParticleDefinition &) override
void InitialiseModel () override
void ModelDescription (std::ostream &outFile) const override
void DeExciteModelDescription (std::ostream &outFile) const override
 G4PreCompoundModel (const G4PreCompoundModel &)=delete
const G4PreCompoundModeloperator= (const G4PreCompoundModel &right)=delete
G4bool operator== (const G4PreCompoundModel &right) const =delete
G4bool operator!= (const G4PreCompoundModel &right) const =delete
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 64 of file G4PreCompoundModel.hh.

Constructor & Destructor Documentation

◆ G4PreCompoundModel() [1/2]

G4PreCompoundModel::G4PreCompoundModel ( G4ExcitationHandler * ptr = nullptr)
explicit

Definition at line 74 of file G4PreCompoundModel.cc.

75 : G4VPreCompoundModel(ptr,"PRECO")
76{
77 //G4cout << "### NEW PrecompoundModel " << this << G4endl;
78 if(nullptr == ptr) { SetExcitationHandler(new G4ExcitationHandler()); }
79
81 proton = G4Proton::Proton();
82 neutron = G4Neutron::Neutron();
83 modelID = G4PhysicsModelCatalog::GetModelID("model_PRECO");
84}
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
static G4NuclearLevelData * GetInstance()
static G4int GetModelID(const G4int modelIndex)
static G4Proton * Proton()
Definition G4Proton.cc:90
G4VPreCompoundModel(G4ExcitationHandler *ptr=nullptr, const G4String &modelName="PrecompoundModel")
void SetExcitationHandler(G4ExcitationHandler *ptr)

Referenced by G4PreCompoundModel(), operator!=(), operator=(), and operator==().

◆ ~G4PreCompoundModel()

G4PreCompoundModel::~G4PreCompoundModel ( )
override

Definition at line 88 of file G4PreCompoundModel.cc.

89{
90 delete theEmission;
91 delete theTransition;
92 delete GetExcitationHandler();
93 delete theMultiFrag;
94 theResult.Clear();
95}
G4ExcitationHandler * GetExcitationHandler() const

◆ G4PreCompoundModel() [2/2]

G4PreCompoundModel::G4PreCompoundModel ( const G4PreCompoundModel & )
delete

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4PreCompoundModel::ApplyYourself ( const G4HadProjectile & thePrimary,
G4Nucleus & theNucleus )
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 149 of file G4PreCompoundModel.cc.

151{
152 const G4ParticleDefinition* primary = thePrimary.GetDefinition();
153 if(primary != neutron && primary != proton) {
155 ed << "G4PreCompoundModel is used for ";
156 if(primary) { ed << primary->GetParticleName(); }
157 G4Exception("G4PreCompoundModel::ApplyYourself()","had0033",FatalException,
158 ed,"");
159 return nullptr;
160 }
161
162 G4int Zp = 0;
163 G4int Ap = 1;
164 if(primary == proton) { Zp = 1; }
165
166 G4double timePrimary=thePrimary.GetGlobalTime();
167 G4int A = theNucleus.GetA_asInt();
168 G4int Z = theNucleus.GetZ_asInt();
169
170 //G4cout << "### G4PreCompoundModel::ApplyYourself: A= " << A << " Z= " << Z
171 // << " Ap= " << Ap << " Zp= " << Zp << G4endl;
172
173 // 4-Momentum
174 G4LorentzVector p = thePrimary.Get4Momentum();
176 p += G4LorentzVector(0.0,0.0,0.0,mass);
177 //G4cout << "Primary 4-mom " << p << " mass= " << mass << G4endl;
178
179 // prepare fragment
180 G4Fragment anInitialState(A + Ap, Z + Zp, p);
181 anInitialState.SetNumberOfExcitedParticle(2, 1);
182 anInitialState.SetNumberOfHoles(1,0);
183 anInitialState.SetCreationTime(thePrimary.GetGlobalTime());
184 anInitialState.SetCreatorModelID(modelID);
185
186 // call excitation handler
187 G4ReactionProductVector* result = DeExcite(anInitialState);
188
189 // fill particle change
190 theResult.Clear();
191 theResult.SetStatusChange(stopAndKill);
192 for(auto const & prod : *result) {
193 G4DynamicParticle * aNewDP = new G4DynamicParticle(prod->GetDefinition(),
194 prod->GetTotalEnergy(),
195 prod->GetMomentum());
196 G4HadSecondary aNew = G4HadSecondary(aNewDP);
197 G4double time = std::max(prod->GetFormationTime(), 0.0);
198 aNew.SetTime(timePrimary + time);
199 aNew.SetCreatorModelID(prod->GetCreatorModelID());
200 delete prod;
201 theResult.AddSecondary(aNew);
202 }
203 delete result;
204 return &theResult;
205}
@ 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]
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const
G4double GetGlobalTime() const
void SetTime(G4double aT)
void SetCreatorModelID(G4int id)
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
G4ReactionProductVector * DeExcite(G4Fragment &aFragment) override

◆ BuildPhysicsTable()

void G4PreCompoundModel::BuildPhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 99 of file G4PreCompoundModel.cc.

100{
102}
void InitialiseModel() override

◆ DeExcite()

G4ReactionProductVector * G4PreCompoundModel::DeExcite ( G4Fragment & aFragment)
overridevirtual

Implements G4VPreCompoundModel.

Definition at line 209 of file G4PreCompoundModel.cc.

210{
211 if (!isInitialised) { InitialiseModel(); }
212 if (usePrecoInterface) { return fInterface->DeExcite(aFragment); }
213
215
216 G4int Z = aFragment.GetZ_asInt();
217 G4int A = aFragment.GetA_asInt();
218 G4double U = aFragment.GetExcitationEnergy();
219
220 if (1 < fVerbose) {
221 G4cout << "### G4PreCompoundModel::DeExcite Z=" << Z << " A=" << A
222 << " U(MeV)=" << U << " Umin=" << fLowLimitExc*A << G4endl;
223 G4cout << aFragment << G4endl;
224 }
225
226 // Conditions to skip pre-compound and perform equilibrium emission
227 if (!isActive || (Z < minZ && A < minA) || U < fLowLimitExc*A ||
228 U > fHighLimitExc*A || 0 < aFragment.GetNumberOfLambdas() ||
229 0 == Z || Z >= A) {
230 PerformEquilibriumEmission(aFragment, result);
231 return result;
232 }
233
234 // apply multi-fragmentation above threshold
235 G4FragmentVector* theTempResult;
236 if (U > A*fMinEForMultiFrag) {
237 theTempResult = theMultiFrag->BreakItUp(aFragment);
238 }
239 else {
240 theTempResult = new G4FragmentVector();
241 theTempResult->push_back(new G4Fragment(aFragment));
242 }
243
244 // apply pre-compound on the list of fragments
245 for (auto & ptr : *theTempResult) {
246 DoIt(result, *ptr);
247 delete ptr;
248 }
249 delete theTempResult;
250 return result;
251}
std::vector< G4Fragment * > G4FragmentVector
Definition G4Fragment.hh:65
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double GetExcitationEnergy() const
G4int GetZ_asInt() const
G4int GetNumberOfLambdas() const
G4int GetA_asInt() const

Referenced by ApplyYourself().

◆ DeExciteModelDescription()

void G4PreCompoundModel::DeExciteModelDescription ( std::ostream & outFile) const
overridevirtual

Implements G4VPreCompoundModel.

Definition at line 384 of file G4PreCompoundModel.cc.

385{
386 outFile << "description of precompound model as used with DeExcite()" << "\n";
387}

◆ InitialiseModel()

void G4PreCompoundModel::InitialiseModel ( )
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 106 of file G4PreCompoundModel.cc.

107{
108 if (isInitialised) { return; }
109 isInitialised = true;
110
111 G4DeexPrecoParameters* param = fNuclData->GetParameters();
112 fVerbose = param->GetVerbose();
113
114 if (param->PrecoDummy() || eDeexcitation == param->GetPreCompoundType()) {
115 isActive = false;
116 }
117 else if (ePrecoInterface == param->GetPreCompoundType()) {
118 usePrecoInterface = true;
119 fInterface = new G4PreCompoundInterface();
120 }
121 else {
122 theEmission = new G4PreCompoundEmission();
123 if (param->UseHETC()) { theEmission->SetHETCModel(); }
124 theEmission->SetOPTxs(param->GetPrecoModelType());
125
126 if (param->UseGNASH()) { theTransition = new G4GNASHTransitions; }
127 else { theTransition = new G4PreCompoundTransitions(); }
128 theTransition->UseNGB(param->NeverGoBack());
129 theTransition->UseCEMtr(param->UseCEM());
130 }
131 if (isActive) {
132 theMultiFrag = new G4StatMF();
133 fLowLimitExc = param->GetPrecoLowEnergy();
134 fHighLimitExc = param->GetPrecoHighEnergy();
135 fMinEForMultiFrag = param->GetMinExPerNucleounForMF();
136
137 useSCO = param->UseSoftCutoff();
138
139 minZ = param->GetMinZForPreco();
140 minA = param->GetMinAForPreco();
141 }
142
144}
@ ePrecoInterface
G4double GetPrecoHighEnergy() const
G4double GetMinExPerNucleounForMF() const
G4PreCompoundType GetPreCompoundType() const

Referenced by BuildPhysicsTable(), and DeExcite().

◆ ModelDescription()

void G4PreCompoundModel::ModelDescription ( std::ostream & outFile) const
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 358 of file G4PreCompoundModel.cc.

359{
360 outFile
361 << "The GEANT4 precompound model is considered as an extension of the\n"
362 << "hadron kinetic model. It gives a possibility to extend the low energy range\n"
363 << "of the hadron kinetic model for nucleon-nucleus inelastic collision and it \n"
364 << "provides a ”smooth” transition from kinetic stage of reaction described by the\n"
365 << "hadron kinetic model to the equilibrium stage of reaction described by the\n"
366 << "equilibrium deexcitation models.\n"
367 << "The initial information for calculation of pre-compound nuclear stage\n"
368 << "consists of the atomic mass number A, charge Z of residual nucleus, its\n"
369 << "four momentum P0, excitation energy U and number of excitons n, which equals\n"
370 << "the sum of the number of particles p (from them p_Z are charged) and the number of\n"
371 << "holes h.\n"
372 << "At the preequilibrium stage of reaction, we follow the exciton model approach in ref. [1],\n"
373 << "taking into account the competition among all possible nuclear transitions\n"
374 << "with ∆n = +2, −2, 0 (which are defined by their associated transition probabilities) and\n"
375 << "the emission of neutrons, protons, deuterons, thritium and helium nuclei (also defined by\n"
376 << "their associated emission probabilities according to exciton model)\n"
377 << "\n"
378 << "[1] K.K. Gudima, S.G. Mashnik, V.D. Toneev, Nucl. Phys. A401 329 (1983)\n"
379 << "\n";
380}

◆ operator!=()

G4bool G4PreCompoundModel::operator!= ( const G4PreCompoundModel & right) const
delete

◆ operator=()

const G4PreCompoundModel & G4PreCompoundModel::operator= ( const G4PreCompoundModel & right)
delete

◆ operator==()

G4bool G4PreCompoundModel::operator== ( const G4PreCompoundModel & right) const
delete

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