Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4PreCompoundModel.cc
이 파일의 문서화 페이지로 가기
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 108685 2018-02-27 07:58:38Z 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
43 // initialisation
44 // - options are set at initialisation
45 // - do not use copy-constructors for G4Fragment
46 // 03.01.2012 V.Ivanchenko Added pointer to G4ExcitationHandler to the
47 // constructor
48 
49 #include "G4PreCompoundModel.hh"
50 #include "G4PhysicalConstants.hh"
51 #include "G4SystemOfUnits.hh"
52 #include "G4PreCompoundEmission.hh"
54 #include "G4GNASHTransitions.hh"
55 #include "G4ParticleDefinition.hh"
56 #include "G4Proton.hh"
57 #include "G4Neutron.hh"
58 
59 #include "G4NucleiProperties.hh"
60 #include "G4NuclearLevelData.hh"
61 #include "G4DeexPrecoParameters.hh"
62 #include "Randomize.hh"
63 #include "G4DynamicParticle.hh"
64 #include "G4ParticleTypes.hh"
65 #include "G4ParticleTable.hh"
66 #include "G4LorentzVector.hh"
67 #include "G4Exp.hh"
68 
70 
72  : G4VPreCompoundModel(ptr,"PRECO"),theEmission(nullptr),theTransition(nullptr),
73  useSCO(false),isInitialised(false),isActive(true),minZ(3),minA(5)
74 {
75  //G4cout << "### NEW PrecompoundModel " << this << G4endl;
76  if(!ptr) { SetExcitationHandler(new G4ExcitationHandler()); }
77 
81 }
82 
84 
86 {
87  delete theEmission;
88  delete theTransition;
89  delete GetExcitationHandler();
90 }
91 
93 
95 {
97 }
98 
100 
102 {
103  if(isInitialised) { return; }
104  isInitialised = true;
105 
106  //G4cout << "G4PreCompoundModel::InitialiseModel() started" << G4endl;
107 
108  G4DeexPrecoParameters* param =
110 
111  // 12/pi2 factor is used in real computation
112  fLevelDensity = param->GetLevelDensity()*12.0/CLHEP::pi2;
113  fLimitEnergy = param->GetPrecoLowEnergy();
114 
115  useSCO = param->UseSoftCutoff();
116 
117  minZ = param->GetMinZForPreco();
118  minA = param->GetMinAForPreco();
119 
121  if(param->UseHETC()) { theEmission->SetHETCModel(); }
122  //else { theEmission->SetDefaultModel(); }
124 
125  if(param->UseGNASH()) { theTransition = new G4GNASHTransitions; }
126  else { theTransition = new G4PreCompoundTransitions(); }
127  theTransition->UseNGB(param->NeverGoBack());
128  theTransition->UseCEMtr(param->UseCEM());
129 
130  if(param->PrecoDummy()) { isActive = false; }
131 
133 }
134 
136 
139  G4Nucleus & theNucleus)
140 {
141  const G4ParticleDefinition* primary = thePrimary.GetDefinition();
142  if(primary != neutron && primary != proton) {
144  ed << "G4PreCompoundModel is used for ";
145  if(primary) { ed << primary->GetParticleName(); }
146  G4Exception("G4PreCompoundModel::ApplyYourself()","had0033",FatalException,
147  ed,"");
148  return 0;
149  }
150 
151  G4int Zp = 0;
152  G4int Ap = 1;
153  if(primary == proton) { Zp = 1; }
154 
155  G4double timePrimary=thePrimary.GetGlobalTime();
156 
157  G4int A = theNucleus.GetA_asInt();
158  G4int Z = theNucleus.GetZ_asInt();
159 
160  //G4cout << "### G4PreCompoundModel::ApplyYourself: A= " << A << " Z= " << Z
161  // << " Ap= " << Ap << " Zp= " << Zp << G4endl;
162  // 4-Momentum
163  G4LorentzVector p = thePrimary.Get4Momentum();
165  p += G4LorentzVector(0.0,0.0,0.0,mass);
166  //G4cout << "Primary 4-mom " << p << " mass= " << mass << G4endl;
167 
168  // prepare fragment
169  G4Fragment anInitialState(A + Ap, Z + Zp, p);
170  anInitialState.SetNumberOfExcitedParticle(2, 1);
171  anInitialState.SetNumberOfHoles(1,0);
172  anInitialState.SetCreationTime(thePrimary.GetGlobalTime());
173 
174  // call excitation handler
175  G4ReactionProductVector * result = DeExcite(anInitialState);
176 
177  // fill particle change
178  theResult.Clear();
180  for(G4ReactionProductVector::iterator i= result->begin();
181  i != result->end(); ++i)
182  {
183  G4DynamicParticle * aNewDP =
184  new G4DynamicParticle((*i)->GetDefinition(),
185  (*i)->GetTotalEnergy(),
186  (*i)->GetMomentum());
187  G4HadSecondary aNew = G4HadSecondary(aNewDP);
188  G4double time=(*i)->GetFormationTime();
189  if(time < 0.0) { time = 0.0; }
190  aNew.SetTime(timePrimary + time);
191  aNew.SetCreatorModelType((*i)->GetCreatorModel());
192  delete (*i);
193  theResult.AddSecondary(aNew);
194  }
195  delete result;
196 
197  //return the filled particle change
198  return &theResult;
199 }
200 
202 
204 {
205  if(!isInitialised) { InitialiseModel(); }
206 
208  G4double Eex = aFragment.GetExcitationEnergy();
209  G4int Z = aFragment.GetZ_asInt();
210  G4int A = aFragment.GetA_asInt();
211 
212  //G4cout << "### G4PreCompoundModel::DeExcite" << G4endl;
213  //G4cout << aFragment << G4endl;
214 
215  // Perform Equilibrium Emission
216  if (!isActive || (Z < minZ && A < minA) || Eex < fLimitEnergy*A) {
217  PerformEquilibriumEmission(aFragment, Result);
218  return Result;
219  }
220 
221  // main loop
222  G4int count = 0;
223  static const G4int countmax = 1000;
224  for (;;) {
225  //G4cout << "### PreCompound loop over fragment" << G4endl;
226  //G4cout << aFragment << G4endl;
227  G4int EquilibriumExcitonNumber =
228  G4lrint(std::sqrt(aFragment.GetExcitationEnergy()
229  *aFragment.GetA_asInt()*fLevelDensity));
230  //
231  // G4cout<<"Neq="<<EquilibriumExcitonNumber<<G4endl;
232  //
233  // J. M. Quesada (Jan. 08) equilibrium hole number could be used as preeq.
234  // evap. delimiter (IAEA report)
235 
236  // Loop for transitions, it is performed while there are
237  // preequilibrium transitions.
238  G4bool ThereIsTransition = false;
239 
240  // G4cout<<"----------------------------------------"<<G4endl;
241  // G4double NP=aFragment.GetNumberOfParticles();
242  // G4double NH=aFragment.GetNumberOfHoles();
243  // G4double NE=aFragment.GetNumberOfExcitons();
244  // G4cout<<" Ex. Energy="<<aFragment.GetExcitationEnergy()<<G4endl;
245  // G4cout<<"N. excitons="<<NE<<" N. Part="<<NP<<"N. Holes ="<<NH<<G4endl;
246  do {
247  ++count;
248  //G4cout<<"transition number .."<<count
249  // <<" n ="<<aFragment.GetNumberOfExcitons()<<G4endl;
250  G4bool go_ahead = false;
251  // soft cutoff criterium as an "ad-hoc" solution to force
252  // increase in evaporation
253  G4int test = aFragment.GetNumberOfExcitons();
254  if (test <= EquilibriumExcitonNumber) { go_ahead=true; }
255 
256  //J. M. Quesada (Apr. 08): soft-cutoff switched off by default
257  if (useSCO && go_ahead)
258  {
259  G4double x = G4double(test)/G4double(EquilibriumExcitonNumber) - 1;
260  if( G4UniformRand() < 1.0 - G4Exp(-x*x/0.32) ) { go_ahead = false; }
261  }
262 
263  // JMQ: WARNING: CalculateProbability MUST be called prior to Get!!
264  // (O values would be returned otherwise)
265  G4double TotalTransitionProbability =
270  //G4cout<<"#0 P1="<<P1<<" P2="<<P2<<" P3="<<P3<<G4endl;
271 
272  //J.M. Quesada (May 2008) Physical criterium (lamdas) PREVAILS over
273  // approximation (critical exciton number)
274  //V.Ivanchenko (May 2011) added check on number of nucleons
275  // to send a fragment to FermiBreakUp
276  if(!go_ahead || P1 <= P2+P3 ||
277  aFragment.GetZ_asInt() < minZ || aFragment.GetA_asInt() < minA ||
278  (aFragment.GetExcitationEnergy() <= fLimitEnergy*aFragment.GetA_asInt()))
279  {
280  //G4cout<<"#4 EquilibriumEmission"<<G4endl;
281  PerformEquilibriumEmission(aFragment,Result);
282  return Result;
283  }
284  else
285  {
286  //
287  // Check if number of excitons is greater than 0
288  // else perform equilibrium emission
289  if (aFragment.GetNumberOfExcitons() <= 0)
290  {
291  PerformEquilibriumEmission(aFragment,Result);
292  return Result;
293  }
294 
295  G4double TotalEmissionProbability =
296  theEmission->GetTotalProbability(aFragment);
297  //
298  //G4cout<<"#1 TotalEmissionProbability="<<TotalEmissionProbability
299  // <<" Nex= " <<aFragment.GetNumberOfExcitons()<<G4endl;
300  //J.M.Quesada (May 08) this has already been done in order to decide
301  // what to do (preeq-eq)
302  // Sum of all probabilities
303  G4double TotalProbability = TotalEmissionProbability
304  + TotalTransitionProbability;
305 
306  // Select subprocess
307  if (TotalProbability*G4UniformRand() > TotalEmissionProbability)
308  {
309  //G4cout<<"#2 Transition"<<G4endl;
310  // It will be transition to state with a new number of excitons
311  ThereIsTransition = true;
312  // Perform the transition
313  theTransition->PerformTransition(aFragment);
314  }
315  else
316  {
317  //G4cout<<"#3 Emission"<<G4endl;
318  // It will be fragment emission
319  ThereIsTransition = false;
320  Result->push_back(theEmission->PerformEmission(aFragment));
321  }
322  }
323  // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
324  } while (ThereIsTransition); // end of do loop
325 
326  // stop if too many iterations
327  if(count >= countmax) {
329  ed << "G4PreCompoundModel loop over " << countmax << " iterations; "
330  << "current G4Fragment: \n" << aFragment;
331  G4Exception("G4PreCompoundModel::DeExcite()","had0034",JustWarning,
332  ed,"");
333  PerformEquilibriumEmission(aFragment, Result);
334  return Result;
335  }
336  } // end of for (;;) loop
337  return Result;
338 }
339 
341 // Initialisation
343 
345 {
346  PrintWarning("UseHETCEmission");
347 }
348 
350 {
351  PrintWarning("UseDefaultEmission");
352 }
353 
355 {
356  PrintWarning("UseGNASHTransition");
357 }
358 
360 {
361  PrintWarning("UseDefaultTransition");
362 }
363 
365 {
366  PrintWarning("UseOPTxs");
367 }
368 
370 {
371  PrintWarning("UseSICB");
372 }
373 
375 {
376  PrintWarning("UseNGB");
377 }
378 
380 {
381  PrintWarning("UseSCO");
382 }
383 
385 {
386  PrintWarning("UseCEMtr");
387 }
388 
390 {
392  ed << "Obsolete method of the preCompound model is called: "
393  << mname << "() \n Instead a corresponding method of "
394  << "G4DeexPrecoParameters class should be used";
395 
396  G4Exception("G4PreCompoundModel::ReadData()","had0803",JustWarning,ed);
397 }
398 
400 // Documentation
402 
403 void G4PreCompoundModel::ModelDescription(std::ostream& outFile) const
404 {
405  outFile
406  << "The GEANT4 precompound model is considered as an extension of the\n"
407  << "hadron kinetic model. It gives a possibility to extend the low energy range\n"
408  << "of the hadron kinetic model for nucleon-nucleus inelastic collision and it \n"
409  << "provides a ”smooth” transition from kinetic stage of reaction described by the\n"
410  << "hadron kinetic model to the equilibrium stage of reaction described by the\n"
411  << "equilibrium deexcitation models.\n"
412  << "The initial information for calculation of pre-compound nuclear stage\n"
413  << "consists of the atomic mass number A, charge Z of residual nucleus, its\n"
414  << "four momentum P0 , excitation energy U and number of excitons n, which equals\n"
415  << "the sum of the number of particles p (from them p_Z are charged) and the number of\n"
416  << "holes h.\n"
417  << "At the preequilibrium stage of reaction, we follow the exciton model approach in ref. [1],\n"
418  << "taking into account the competition among all possible nuclear transitions\n"
419  << "with ∆n = +2, −2, 0 (which are defined by their associated transition probabilities) and\n"
420  << "the emission of neutrons, protons, deutrons, thritium and helium nuclei (also defined by\n"
421  << "their associated emission probabilities according to exciton model)\n"
422  << "\n"
423  << "[1] K.K. Gudima, S.G. Mashnik, V.D. Toneev, Nucl. Phys. A401 329 (1983)\n"
424  << "\n";
425 }
426 
427 void G4PreCompoundModel::DeExciteModelDescription(std::ostream& outFile) const
428 {
429  outFile << "description of precompound model as used with DeExcite()" << "\n";
430 }
431 
Float_t x
Definition: compare.C:6
void SetExcitationHandler(G4ExcitationHandler *ptr)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4HadFinalState theResult
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:45
G4ReactionProduct * PerformEmission(G4Fragment &aFragment)
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:368
void SetCreationTime(G4double time)
Definition: G4Fragment.hh:446
G4DeexPrecoParameters * GetParameters()
G4ExcitationHandler * GetExcitationHandler() const
G4PreCompoundEmission * theEmission
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &thePrimary, G4Nucleus &theNucleus) final
const char * p
Definition: xmltok.h:285
static G4NuclearLevelData * GetInstance()
const G4String & GetParticleName() const
static const G4double * P2[nN]
G4int GetA_asInt() const
Definition: G4Fragment.hh:259
G4double GetPrecoLowEnergy() const
static constexpr double pi2
Definition: SystemOfUnits.h:57
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
static G4Proton * Proton()
Definition: G4Proton.cc:93
virtual G4double CalculateProbability(const G4Fragment &aFragment)=0
void PerformEquilibriumEmission(const G4Fragment &aFragment, G4ReactionProductVector *theResult) const
virtual void ModelDescription(std::ostream &outFile) const final
Float_t Z
virtual void BuildPhysicsTable(const G4ParticleDefinition &) final
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4double GetGlobalTime() const
virtual void DeExciteModelDescription(std::ostream &outFile) const final
#define G4UniformRand()
Definition: Randomize.hh:53
double A(double temperature)
static G4double GetNuclearMass(const G4double A, const G4double Z)
void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP)
Definition: G4Fragment.hh:349
G4double GetTotalProbability(const G4Fragment &aFragment)
virtual void PerformTransition(G4Fragment &aFragment)=0
G4PreCompoundModel(G4ExcitationHandler *ptr=nullptr)
std::vector< G4ReactionProduct * > G4ReactionProductVector
G4int GetZ_asInt() const
Definition: G4Fragment.hh:264
const G4LorentzVector & Get4Momentum() const
G4double G4ParticleHPJENDLHEData::G4double result
void PrintWarning(const G4String &mname)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
const G4ParticleDefinition * GetDefinition() const
G4double GetLevelDensity() const
int G4lrint(double ad)
Definition: templates.hh:151
int G4int
Definition: G4Types.hh:78
const G4ParticleDefinition * proton
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
virtual void InitialiseModel() final
const G4ParticleDefinition * neutron
static const G4double * P1[nN]
CLHEP::HepLorentzVector G4LorentzVector
G4int GetNumberOfExcitons() const
Definition: G4Fragment.hh:333
void SetOPTxs(G4int opt)
G4VPreCompoundTransitions * theTransition
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:276
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
void SetStatusChange(G4HadFinalStateStatus aS)
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment) final