Geant4  10.01.p01
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 // $Id: G4PreCompoundModel.cc 80062 2014-03-31 13:41:30Z gcosmo $
27 //
28 // by V. Lara
29 //
30 // Modified:
31 // 01.04.2008 J.M.Quesada Several changes. Soft cut-off switched off.
32 // 01.05.2008 J.M.Quesada Protection against non-physical preeq.
33 // transitional regime has been set
34 // 03.09.2008 J.M.Quesada for external choice of inverse cross section option
35 // 06.09.2008 J.M.Quesada Also external choices have been added for:
36 // - superimposed Coulomb barrier (useSICB=true)
37 // - "never go back" hipothesis (useNGB=true)
38 // - soft cutoff from preeq. to equlibrium (useSCO=true)
39 // - CEM transition probabilities (useCEMtr=true)
40 // 20.08.2010 V.Ivanchenko Cleanup of the code:
41 // - integer Z and A;
42 // - emission and transition classes created at initialisation
43 // - options are set at initialisation
44 // - do not use copy-constructors for G4Fragment
45 // 03.01.2012 V.Ivanchenko Added pointer to G4ExcitationHandler to the
46 // constructor
47 
48 #include "G4PreCompoundModel.hh"
49 #include "G4PhysicalConstants.hh"
50 #include "G4SystemOfUnits.hh"
51 #include "G4PreCompoundEmission.hh"
53 #include "G4GNASHTransitions.hh"
54 #include "G4ParticleDefinition.hh"
55 #include "G4Proton.hh"
56 #include "G4Neutron.hh"
57 
58 #include "G4NucleiProperties.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"
66 
68  : G4VPreCompoundModel(ptr,"PRECO"), useHETCEmission(false),
69  useGNASHTransition(false), OPTxs(3), useSICB(false),
70  useNGB(false), useSCO(false), useCEMtr(true), maxZ(3), maxA(5)
71  //maxZ(9), maxA(17)
72 {
73  if(!ptr) { SetExcitationHandler(new G4ExcitationHandler()); }
76 
79  else { theEmission->SetDefaultModel(); }
82 
84  else { theTransition = new G4PreCompoundTransitions(); }
87 
90 }
91 
93 {
94  delete theEmission;
95  delete theTransition;
96  delete GetExcitationHandler();
97 }
98 
99 void G4PreCompoundModel::ModelDescription(std::ostream& outFile) const
100 {
101  outFile << "The GEANT4 precompound model is considered as an extension of the\n"
102  << "hadron kinetic model. It gives a possibility to extend the low energy range\n"
103  << "of the hadron kinetic model for nucleon-nucleus inelastic collision and it \n"
104  << "provides a ”smooth” transition from kinetic stage of reaction described by the\n"
105  << "hadron kinetic model to the equilibrium stage of reaction described by the\n"
106  << "equilibrium deexcitation models.\n"
107  << "The initial information for calculation of pre-compound nuclear stage\n"
108  << "consists of the atomic mass number A, charge Z of residual nucleus, its\n"
109  << "four momentum P0 , excitation energy U and number of excitons n, which equals\n"
110  << "the sum of the number of particles p (from them p_Z are charged) and the number of\n"
111  << "holes h.\n"
112  << "At the preequilibrium stage of reaction, we follow the exciton model approach in ref. [1],\n"
113  << "taking into account the competition among all possible nuclear transitions\n"
114  << "with ∆n = +2, −2, 0 (which are defined by their associated transition probabilities) and\n"
115  << "the emission of neutrons, protons, deutrons, thritium and helium nuclei (also defined by\n"
116  << "their associated emission probabilities according to exciton model)\n"
117  << "\n"
118  << "[1] K.K. Gudima, S.G. Mashnik, V.D. Toneev, Nucl. Phys. A401 329 (1983)\n"
119  << "\n";
120 }
121 void G4PreCompoundModel::DeExciteModelDescription(std::ostream& outFile) const
122 {
123  outFile << "description of precompound model as used with DeExcite()"
124  << "\n";
125 }
127 
129  G4Nucleus & theNucleus)
130 {
131  const G4ParticleDefinition* primary = thePrimary.GetDefinition();
132  if(primary != neutron && primary != proton) {
133  std::ostringstream errOs;
134  errOs << "BAD primary type in G4PreCompoundModel: "
135  << primary->GetParticleName() <<G4endl;
136  throw G4HadronicException(__FILE__, __LINE__, errOs.str());
137  }
138 
139  G4int Zp = 0;
140  G4int Ap = 1;
141  if(primary == proton) { Zp = 1; }
142 
143  G4int A = theNucleus.GetA_asInt();
144  G4int Z = theNucleus.GetZ_asInt();
145 
146  //G4cout << "### G4PreCompoundModel::ApplyYourself: A= " << A << " Z= " << Z
147  // << " Ap= " << Ap << " Zp= " << Zp << G4endl;
148  // 4-Momentum
149  G4LorentzVector p = thePrimary.Get4Momentum();
151  p += G4LorentzVector(0.0,0.0,0.0,mass);
152  //G4cout << "Primary 4-mom " << p << " mass= " << mass << G4endl;
153 
154  // prepare fragment
155  G4Fragment anInitialState(A + Ap, Z + Zp, p);
156 
157  // projectile and target nucleons
158  // Add nucleon on which interaction happens
159  //++Ap;
160  //if(A*G4UniformRand() <= G4double(Z)) { Zp += 1; }
161  anInitialState.SetNumberOfExcitedParticle(2, 1);
162  anInitialState.SetNumberOfHoles(1,0);
163  // anInitialState.SetNumberOfExcitedParticle(Ap, Zp);
164  // anInitialState.SetNumberOfHoles(Ap,Zp);
165 
166  anInitialState.SetCreationTime(thePrimary.GetGlobalTime());
167 
168  // call excitation handler
169  G4ReactionProductVector * result = DeExcite(anInitialState);
170 
171  // fill particle change
172  theResult.Clear();
174  for(G4ReactionProductVector::iterator i= result->begin(); i != result->end(); ++i)
175  {
176  G4DynamicParticle * aNew =
177  new G4DynamicParticle((*i)->GetDefinition(),
178  (*i)->GetTotalEnergy(),
179  (*i)->GetMomentum());
180  delete (*i);
181  theResult.AddSecondary(aNew);
182  }
183  delete result;
184 
185  //return the filled particle change
186  return &theResult;
187 }
188 
190 
192 {
194  G4double Eex = aFragment.GetExcitationEnergy();
195  G4int Z = aFragment.GetZ_asInt();
196  G4int A = aFragment.GetA_asInt();
197 
198  //G4cout << "### G4PreCompoundModel::DeExcite" << G4endl;
199  //G4cout << aFragment << G4endl;
200 
201  // Perform Equilibrium Emission
202  if ((Z < maxZ && A < maxA) || Eex < MeV /*|| Eex > 3.*MeV*A*/) {
203  PerformEquilibriumEmission(aFragment, Result);
204  return Result;
205  }
206 
207  // main loop
208  for (;;) {
209 
210  //fragment++;
211  //G4cout<<"-------------------------------------------------------------------"<<G4endl;
212  //G4cout<<"Fragment number .. "<<fragment<<G4endl;
213 
214  // Initialize fragment according with the nucleus parameters
215  //G4cout << "### Loop over fragment" << G4endl;
216  //G4cout << aFragment << G4endl;
217 
218  theEmission->Initialize(aFragment);
219 
220  G4double gg = (6.0/pi2)*aFragment.GetA_asInt()*fLevelDensity;
221 
222  G4int EquilibriumExcitonNumber =
223  G4lrint(std::sqrt(2*gg*aFragment.GetExcitationEnergy()));
224  //
225  // G4cout<<"Neq="<<EquilibriumExcitonNumber<<G4endl;
226  //
227  // J. M. Quesada (Jan. 08) equilibrium hole number could be used as preeq.
228  // evap. delimiter (IAEA report)
229 
230  // Loop for transitions, it is performed while there are
231  // preequilibrium transitions.
232  G4bool ThereIsTransition = false;
233 
234  // G4cout<<"----------------------------------------"<<G4endl;
235  // G4double NP=aFragment.GetNumberOfParticles();
236  // G4double NH=aFragment.GetNumberOfHoles();
237  // G4double NE=aFragment.GetNumberOfExcitons();
238  // G4cout<<" Ex. Energy="<<aFragment.GetExcitationEnergy()<<G4endl;
239  // G4cout<<"N. excitons="<<NE<<" N. Part="<<NP<<"N. Holes ="<<NH<<G4endl;
240  //G4int transition=0;
241  do {
242  //transition++;
243  //G4cout<<"transition number .."<<transition<<G4endl;
244  //G4cout<<" n ="<<aFragment.GetNumberOfExcitons()<<G4endl;
245  G4bool go_ahead = false;
246  // soft cutoff criterium as an "ad-hoc" solution to force
247  // increase in evaporation
248  G4int test = aFragment.GetNumberOfExcitons();
249  if (test <= EquilibriumExcitonNumber) { go_ahead=true; }
250 
251  //J. M. Quesada (Apr. 08): soft-cutoff switched off by default
252  if (useSCO && go_ahead)
253  {
254  G4double x = G4double(test)/G4double(EquilibriumExcitonNumber) - 1;
255  if( G4UniformRand() < 1.0 - G4Exp(-x*x/0.32) ) { go_ahead = false; }
256  }
257 
258  // JMQ: WARNING: CalculateProbability MUST be called prior to Get methods !!
259  // (O values would be returned otherwise)
260  G4double TotalTransitionProbability =
265  //G4cout<<"#0 P1="<<P1<<" P2="<<P2<<" P3="<<P3<<G4endl;
266 
267  //J.M. Quesada (May 2008) Physical criterium (lamdas) PREVAILS over
268  // approximation (critical exciton number)
269  //V.Ivanchenko (May 2011) added check on number of nucleons
270  // to send a fragment to FermiBreakUp
271  if(!go_ahead || P1 <= P2+P3 ||
272  (aFragment.GetZ_asInt() < maxZ && aFragment.GetA_asInt() < maxA) )
273  {
274  //G4cout<<"#4 EquilibriumEmission"<<G4endl;
275  PerformEquilibriumEmission(aFragment,Result);
276  return Result;
277  }
278  else
279  {
280  G4double TotalEmissionProbability =
281  theEmission->GetTotalProbability(aFragment);
282  //
283  //G4cout<<"#1 TotalEmissionProbability="<<TotalEmissionProbability
284  // <<" Nex= " <<aFragment.GetNumberOfExcitons()<<G4endl;
285  //
286  // Check if number of excitons is greater than 0
287  // else perform equilibrium emission
288  if (aFragment.GetNumberOfExcitons() <= 0)
289  {
290  PerformEquilibriumEmission(aFragment,Result);
291  return Result;
292  }
293 
294  //J.M.Quesada (May 08) this has already been done in order to decide
295  // what to do (preeq-eq)
296  // Sum of all probabilities
297  G4double TotalProbability = TotalEmissionProbability
298  + TotalTransitionProbability;
299 
300  // Select subprocess
301  if (TotalProbability*G4UniformRand() > TotalEmissionProbability)
302  {
303  //G4cout<<"#2 Transition"<<G4endl;
304  // It will be transition to state with a new number of excitons
305  ThereIsTransition = true;
306  // Perform the transition
307  theTransition->PerformTransition(aFragment);
308  }
309  else
310  {
311  //G4cout<<"#3 Emission"<<G4endl;
312  // It will be fragment emission
313  ThereIsTransition = false;
314  Result->push_back(theEmission->PerformEmission(aFragment));
315  }
316  }
317  } while (ThereIsTransition); // end of do loop
318  } // end of for (;;) loop
319  return Result;
320 }
321 
323 // Initialisation
325 
327 {
328  useHETCEmission = true;
330 }
331 
333 {
334  useHETCEmission = false;
336 }
337 
339  useGNASHTransition = true;
340  delete theTransition;
344 }
345 
347  useGNASHTransition = false;
348  delete theTransition;
352 }
353 
355 {
356  OPTxs = opt;
358 }
359 
361 {
362  useSICB = true;
364 }
365 
367 {
368  useNGB = true;
369 }
370 
372 {
373  useSCO = true;
374 }
375 
377 {
378  useCEMtr = true;
379 }
380 
virtual void PerformTransition(G4Fragment &aFragment)=0
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
void PerformEquilibriumEmission(const G4Fragment &aFragment, G4ReactionProductVector *theResult) const
static const double MeV
Definition: G4SIunits.hh:193
static G4double GetNuclearMass(const G4double A, const G4double Z)
static const G4double * P1[nN]
G4PreCompoundEmission * theEmission
G4PreCompoundModel(G4ExcitationHandler *ptr=0)
const G4ParticleDefinition * neutron
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:360
int G4int
Definition: G4Types.hh:78
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &thePrimary, G4Nucleus &theNucleus)
const G4String & GetParticleName() const
G4HadFinalState theResult
void SetStatusChange(G4HadFinalStateStatus aS)
void SetExcitationHandler(G4ExcitationHandler *ptr)
std::vector< G4ReactionProduct * > G4ReactionProductVector
void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP)
Definition: G4Fragment.hh:341
G4ExcitationHandler * GetExcitationHandler() const
#define G4UniformRand()
Definition: Randomize.hh:95
G4int GetA_asInt() const
Definition: G4Fragment.hh:243
const G4ParticleDefinition * GetDefinition() const
G4VPreCompoundTransitions * theTransition
bool G4bool
Definition: G4Types.hh:79
G4double GetGlobalTime() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static const G4double A[nN]
const G4LorentzVector & Get4Momentum() const
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment)
virtual G4double CalculateProbability(const G4Fragment &aFragment)=0
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4ReactionProduct * PerformEmission(G4Fragment &aFragment)
void SetCreationTime(G4double time)
Definition: G4Fragment.hh:418
G4int GetNumberOfExcitons() const
Definition: G4Fragment.hh:325
int G4lrint(double ad)
Definition: templates.hh:163
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
virtual void ModelDescription(std::ostream &outFile) const
G4int GetZ_asInt() const
Definition: G4Fragment.hh:248
static const G4double * P2[nN]
#define G4endl
Definition: G4ios.hh:61
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
const G4ParticleDefinition * proton
double G4double
Definition: G4Types.hh:76
virtual void DeExciteModelDescription(std::ostream &outFile) const
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:260
void SetOPTxs(G4int opt)
void Initialize(const G4Fragment &aFragment)
G4double GetTotalProbability(const G4Fragment &aFragment)
CLHEP::HepLorentzVector G4LorentzVector