Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PreCompoundModel.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// by V. Lara
27//
28// Modified:
29// 01.04.2008 J.M.Quesada Several changes. Soft cut-off switched off.
30// 01.05.2008 J.M.Quesada Protection against non-physical preeq.
31// transitional regime has been set
32// 03.09.2008 J.M.Quesada for external choice of inverse cross section option
33// 06.09.2008 J.M.Quesada Also external choices have been added for:
34// - superimposed Coulomb barrier (useSICB=true)
35// - "never go back" hipothesis (useNGB=true)
36// - soft cutoff from preeq. to equlibrium (useSCO=true)
37// - CEM transition probabilities (useCEMtr=true)
38// 20.08.2010 V.Ivanchenko Cleanup of the code:
39// - integer Z and A;
40// - emission and transition classes created at
41// initialisation
42// - options are set at initialisation
43// - do not use copy-constructors for G4Fragment
44// 03.01.2012 V.Ivanchenko Added pointer to G4ExcitationHandler to the
45// constructor
46
47#include "G4PreCompoundModel.hh"
49#include "G4SystemOfUnits.hh"
52#include "G4GNASHTransitions.hh"
54#include "G4Proton.hh"
55#include "G4Neutron.hh"
56
57#include "G4NucleiProperties.hh"
58#include "G4NuclearLevelData.hh"
60#include "Randomize.hh"
61#include "G4DynamicParticle.hh"
62#include "G4ParticleTypes.hh"
63#include "G4ParticleTable.hh"
64#include "G4LorentzVector.hh"
65#include "G4Exp.hh"
67
68////////////////////////////////////////////////////////////////////////////////
69
71 : G4VPreCompoundModel(ptr,"PRECO")
72{
73 //G4cout << "### NEW PrecompoundModel " << this << G4endl;
74 if(nullptr == ptr) { SetExcitationHandler(new G4ExcitationHandler()); }
75
77 proton = G4Proton::Proton();
78 neutron = G4Neutron::Neutron();
79 modelID = G4PhysicsModelCatalog::GetModelID("model_PRECO");
80}
81
82////////////////////////////////////////////////////////////////////////////////
83
85{
86 delete theEmission;
87 delete theTransition;
88 delete GetExcitationHandler();
89 theResult.Clear();
90}
91
92////////////////////////////////////////////////////////////////////////////////
93
95{
97}
98
99////////////////////////////////////////////////////////////////////////////////
100
102{
103 if(isInitialised) { return; }
104 isInitialised = true;
105
106 //G4cout << "G4PreCompoundModel::InitialiseModel() started" << G4endl;
107
108 G4DeexPrecoParameters* param = fNuclData->GetParameters();
109
110 fLowLimitExc = param->GetPrecoLowEnergy();
111 fHighLimitExc = param->GetPrecoHighEnergy();
112
113 useSCO = param->UseSoftCutoff();
114
115 minZ = param->GetMinZForPreco();
116 minA = param->GetMinAForPreco();
117
118 theEmission = new G4PreCompoundEmission();
119 if(param->UseHETC()) { theEmission->SetHETCModel(); }
120 theEmission->SetOPTxs(param->GetPrecoModelType());
121
122 if(param->UseGNASH()) { theTransition = new G4GNASHTransitions; }
123 else { theTransition = new G4PreCompoundTransitions(); }
124 theTransition->UseNGB(param->NeverGoBack());
125 theTransition->UseCEMtr(param->UseCEM());
126
127 if(param->PrecoDummy()) { isActive = false; }
128
130}
131
132////////////////////////////////////////////////////////////////////////////////
133
136 G4Nucleus & theNucleus)
137{
138 const G4ParticleDefinition* primary = thePrimary.GetDefinition();
139 if(primary != neutron && primary != proton) {
141 ed << "G4PreCompoundModel is used for ";
142 if(primary) { ed << primary->GetParticleName(); }
143 G4Exception("G4PreCompoundModel::ApplyYourself()","had0033",FatalException,
144 ed,"");
145 return nullptr;
146 }
147
148 G4int Zp = 0;
149 G4int Ap = 1;
150 if(primary == proton) { Zp = 1; }
151
152 G4double timePrimary=thePrimary.GetGlobalTime();
153
154 G4int A = theNucleus.GetA_asInt();
155 G4int Z = theNucleus.GetZ_asInt();
156
157 //G4cout << "### G4PreCompoundModel::ApplyYourself: A= " << A << " Z= " << Z
158 // << " Ap= " << Ap << " Zp= " << Zp << G4endl;
159 // 4-Momentum
160 G4LorentzVector p = thePrimary.Get4Momentum();
162 p += G4LorentzVector(0.0,0.0,0.0,mass);
163 //G4cout << "Primary 4-mom " << p << " mass= " << mass << G4endl;
164
165 // prepare fragment
166 G4Fragment anInitialState(A + Ap, Z + Zp, p);
167 anInitialState.SetNumberOfExcitedParticle(2, 1);
168 anInitialState.SetNumberOfHoles(1,0);
169 anInitialState.SetCreationTime(thePrimary.GetGlobalTime());
170 anInitialState.SetCreatorModelID(modelID);
171
172 // call excitation handler
173 G4ReactionProductVector* result = DeExcite(anInitialState);
174
175 // fill particle change
176 theResult.Clear();
177 theResult.SetStatusChange(stopAndKill);
178 for(auto const & prod : *result) {
179 G4DynamicParticle * aNewDP = new G4DynamicParticle(prod->GetDefinition(),
180 prod->GetTotalEnergy(),
181 prod->GetMomentum());
182 G4HadSecondary aNew = G4HadSecondary(aNewDP);
183 G4double time = std::max(prod->GetFormationTime(), 0.0);
184 aNew.SetTime(timePrimary + time);
185 aNew.SetCreatorModelID(prod->GetCreatorModelID());
186 delete prod;
187 theResult.AddSecondary(aNew);
188 }
189 delete result;
190
191 //return the filled particle change
192 return &theResult;
193}
194
195////////////////////////////////////////////////////////////////////////////////
196
198{
199 if(!isInitialised) { InitialiseModel(); }
200
202 G4double U = aFragment.GetExcitationEnergy();
203 G4int Z = aFragment.GetZ_asInt();
204 G4int A = aFragment.GetA_asInt();
205
206 //G4cout << "### G4PreCompoundModel::DeExcite" << G4endl;
207 //G4cout << aFragment << G4endl;
208
209 // Conditions to skip pre-compound and perform equilibrium emission
210 if (!isActive || (Z < minZ && A < minA) ||
211 U < fLowLimitExc*A || U > A*fHighLimitExc || 0 < aFragment.GetNumberOfLambdas()) {
212 PerformEquilibriumEmission(aFragment, Result);
213 return Result;
214 }
215
216 // main loop
217 G4int count = 0;
218 const G4double ldfact = 12.0/CLHEP::pi2;
219 const G4int countmax = 1000;
220 for (;;) {
221 //G4cout << "### PreCompound loop over fragment" << G4endl;
222 //G4cout << aFragment << G4endl;
223 U = aFragment.GetExcitationEnergy();
224 Z = aFragment.GetZ_asInt();
225 A = aFragment.GetA_asInt();
226 G4int eqExcitonNumber =
227 G4lrint(std::sqrt(ldfact*U*fNuclData->GetLevelDensity(Z, A, U)));
228 //
229 // G4cout<<"Neq="<<EquilibriumExcitonNumber<<G4endl;
230 //
231 // J. M. Quesada (Jan. 08) equilibrium hole number could be used as preeq.
232 // evap. delimiter (IAEA report)
233
234 // Loop for transitions, it is performed while there are
235 // preequilibrium transitions.
236 G4bool isTransition = false;
237
238 // G4cout<<"----------------------------------------"<<G4endl;
239 // G4double NP=aFragment.GetNumberOfParticles();
240 // G4double NH=aFragment.GetNumberOfHoles();
241 // G4double NE=aFragment.GetNumberOfExcitons();
242 // G4cout<<" Ex. Energy="<<aFragment.GetExcitationEnergy()<<G4endl;
243 // G4cout<<"N. excitons="<<NE<<" N. Part="<<NP<<"N. Holes ="<<NH<<G4endl;
244 do {
245 ++count;
246 //G4cout<<"transition number .."<<count
247 // <<" n ="<<aFragment.GetNumberOfExcitons()<<G4endl;
248 // soft cutoff criterium as an "ad-hoc" solution to force
249 // increase in evaporation
250 G4int ne = aFragment.GetNumberOfExcitons();
251 G4bool go_ahead = (ne <= eqExcitonNumber);
252
253 //J. M. Quesada (Apr. 08): soft-cutoff switched off by default
254 if (useSCO && go_ahead) {
255 G4double x = (G4double)(ne - eqExcitonNumber)/(G4double)eqExcitonNumber;
256 if( G4UniformRand() < 1.0 - G4Exp(-x*x/0.32) ) { go_ahead = false; }
257 }
258
259 // JMQ: WARNING: CalculateProbability MUST be called prior to Get!!
260 // (O values would be returned otherwise)
261 G4double transProbability =
262 theTransition->CalculateProbability(aFragment);
263 G4double P1 = theTransition->GetTransitionProb1();
264 G4double P2 = theTransition->GetTransitionProb2();
265 G4double P3 = theTransition->GetTransitionProb3();
266 //G4cout<<"#0 P1="<<P1<<" P2="<<P2<<" P3="<<P3<<G4endl;
267
268 //J.M. Quesada (May 2008) Physical criterium (lamdas) PREVAILS over
269 // approximation (critical exciton number)
270 //V.Ivanchenko (May 2011) added check on number of nucleons
271 // to send a fragment to FermiBreakUp
272 // or check on limits of excitation
273 if(!go_ahead || P1 <= P2+P3 || Z < minZ || A < minA ||
274 U <= fLowLimitExc*A || U > A*fHighLimitExc ||
275 aFragment.GetNumberOfExcitons() <= 0) {
276 //G4cout<<"#4 EquilibriumEmission"<<G4endl;
277 PerformEquilibriumEmission(aFragment,Result);
278 return Result;
279 }
280 G4double emissionProbability =
281 theEmission->GetTotalProbability(aFragment);
282 //G4cout<<"#1 TotalEmissionProbability="<<TotalEmissionProbability
283 // <<" Nex= " <<aFragment.GetNumberOfExcitons()<<G4endl;
284 //J.M.Quesada (May 08) this has already been done in order to decide
285 // what to do (preeq-eq)
286 // Sum of all probabilities
287 G4double TotalProbability = emissionProbability + transProbability;
288
289 // Select subprocess
290 if (TotalProbability*G4UniformRand() > emissionProbability) {
291 //G4cout<<"#2 Transition"<<G4endl;
292 // It will be transition to state with a new number of excitons
293 isTransition = true;
294 // Perform the transition
295 theTransition->PerformTransition(aFragment);
296 } else {
297 //G4cout<<"#3 Emission"<<G4endl;
298 // It will be fragment emission
299 isTransition = false;
300 Result->push_back(theEmission->PerformEmission(aFragment));
301 }
302 // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
303 } while (isTransition); // end of do loop
304
305 // stop if too many iterations
306 if(count >= countmax) {
308 ed << "G4PreCompoundModel loop over " << countmax << " iterations; "
309 << "current G4Fragment: \n" << aFragment;
310 G4Exception("G4PreCompoundModel::DeExcite()","had0034",JustWarning,
311 ed,"");
312 PerformEquilibriumEmission(aFragment, Result);
313 return Result;
314 }
315 } // end of for (;;) loop
316 return Result;
317}
318
319////////////////////////////////////////////////////////////////////////////////
320// Documentation
321////////////////////////////////////////////////////////////////////////////////
322
323void G4PreCompoundModel::ModelDescription(std::ostream& outFile) const
324{
325 outFile
326 << "The GEANT4 precompound model is considered as an extension of the\n"
327 << "hadron kinetic model. It gives a possibility to extend the low energy range\n"
328 << "of the hadron kinetic model for nucleon-nucleus inelastic collision and it \n"
329 << "provides a ”smooth” transition from kinetic stage of reaction described by the\n"
330 << "hadron kinetic model to the equilibrium stage of reaction described by the\n"
331 << "equilibrium deexcitation models.\n"
332 << "The initial information for calculation of pre-compound nuclear stage\n"
333 << "consists of the atomic mass number A, charge Z of residual nucleus, its\n"
334 << "four momentum P0 , excitation energy U and number of excitons n, which equals\n"
335 << "the sum of the number of particles p (from them p_Z are charged) and the number of\n"
336 << "holes h.\n"
337 << "At the preequilibrium stage of reaction, we follow the exciton model approach in ref. [1],\n"
338 << "taking into account the competition among all possible nuclear transitions\n"
339 << "with ∆n = +2, −2, 0 (which are defined by their associated transition probabilities) and\n"
340 << "the emission of neutrons, protons, deuterons, thritium and helium nuclei (also defined by\n"
341 << "their associated emission probabilities according to exciton model)\n"
342 << "\n"
343 << "[1] K.K. Gudima, S.G. Mashnik, V.D. Toneev, Nucl. Phys. A401 329 (1983)\n"
344 << "\n";
345}
346
347////////////////////////////////////////////////////////////////////////////////
348
349void G4PreCompoundModel::DeExciteModelDescription(std::ostream& outFile) const
350{
351 outFile << "description of precompound model as used with DeExcite()" << "\n";
352}
353
354////////////////////////////////////////////////////////////////////////////////
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
@ stopAndKill
CLHEP::HepLorentzVector G4LorentzVector
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4UniformRand()
Definition: Randomize.hh:52
G4double GetPrecoLowEnergy() const
G4double GetPrecoHighEnergy() const
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:312
void SetCreatorModelID(G4int value)
Definition: G4Fragment.hh:433
G4int GetZ_asInt() const
Definition: G4Fragment.hh:289
void SetCreationTime(G4double time)
Definition: G4Fragment.hh:484
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:396
G4int GetNumberOfExcitons() const
Definition: G4Fragment.hh:361
G4int GetNumberOfLambdas() const
Definition: G4Fragment.hh:307
void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP)
Definition: G4Fragment.hh:377
G4int GetA_asInt() const
Definition: G4Fragment.hh:284
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const
G4double GetGlobalTime() const
void SetTime(G4double aT)
void SetCreatorModelID(G4int id)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
G4double GetLevelDensity(G4int Z, G4int A, G4double U)
G4DeexPrecoParameters * GetParameters()
static G4NuclearLevelData * GetInstance()
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 G4int GetModelID(const G4int modelIndex)
G4ReactionProduct * PerformEmission(G4Fragment &aFragment)
G4double GetTotalProbability(const G4Fragment &aFragment)
virtual void ModelDescription(std::ostream &outFile) const final
virtual void InitialiseModel() final
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &thePrimary, G4Nucleus &theNucleus) final
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment) final
virtual void DeExciteModelDescription(std::ostream &outFile) const final
virtual void BuildPhysicsTable(const G4ParticleDefinition &) final
G4PreCompoundModel(G4ExcitationHandler *ptr=nullptr)
static G4Proton * Proton()
Definition: G4Proton.cc:92
G4ExcitationHandler * GetExcitationHandler() const
void SetExcitationHandler(G4ExcitationHandler *ptr)
virtual G4double CalculateProbability(const G4Fragment &aFragment)=0
virtual void PerformTransition(G4Fragment &aFragment)=0
int G4lrint(double ad)
Definition: templates.hh:134