Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4PhotonEvaporation.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: G4PhotonEvaporation.cc 102726 2017-02-20 13:15:29Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 class file
31 //
32 // CERN, Geneva, Switzerland
33 //
34 // File name: G4PhotonEvaporation
35 //
36 // Author: Vladimir Ivantchenko
37 //
38 // Creation date: 20 December 2011
39 //
40 //Modifications:
41 //
42 //
43 // -------------------------------------------------------------------
44 //
45 
46 #include "G4PhotonEvaporation.hh"
47 
48 #include "Randomize.hh"
49 #include "G4Gamma.hh"
50 #include "G4LorentzVector.hh"
51 #include "G4FragmentVector.hh"
52 #include "G4GammaTransition.hh"
53 #include "G4Pow.hh"
56 
57 G4float G4PhotonEvaporation::GREnergy[] = {0.0f};
58 G4float G4PhotonEvaporation::GRWidth[] = {0.0f};
59 
61  : fLevelManager(nullptr), fTransition(p), fVerbose(0), fPoints(0),
62  vShellNumber(-1), fIndex(0), fMaxLifeTime(DBL_MAX),
63  fICM(false), fRDM(false), fSampleTime(true), isInitialised(false)
64 {
65  //G4cout << "### New G4PhotonEvaporation() " << this << G4endl;
66  fNuclearLevelData = G4NuclearLevelData::GetInstance();
67  LevelDensity = 0.125/CLHEP::MeV;
68  Tolerance = 0.1*CLHEP::keV;
69 
70  if(!fTransition) { fTransition = new G4GammaTransition(); }
71 
72  theA = theZ = fCode = 0;
73  fLevelEnergyMax = fStep = fExcEnergy = fProbability = 0.0;
74 
75  for(G4int i=0; i<MAXDEPOINT; ++i) { fCummProbability[i] = 0.0; }
76  if(0.0f == GREnergy[1]) {
77  G4Pow* g4calc = G4Pow::GetInstance();
78  static const G4float GRWfactor = 0.3f;
79  for (G4int A=1; A<MAXGRDATA; ++A) {
80  GREnergy[A] = (G4float)(40.3*CLHEP::MeV/g4calc->powZ(A,0.2));
81  GRWidth[A] = GRWfactor*GREnergy[A];
82  }
83  }
84 }
85 
87 {
88  delete fTransition;
89 }
90 
92 {
93  if(isInitialised) { return; }
94  isInitialised = true;
95 
96  if(fVerbose > 0) {
97  G4cout << "### G4PhotonEvaporation is initialized " << this << G4endl;
98  }
99  G4DeexPrecoParameters* param = fNuclearLevelData->GetParameters();
100  LevelDensity = param->GetLevelDensity();
101  Tolerance = param->GetMinExcitation();
102  fMaxLifeTime = param->GetMaxLifeTime();
103 
104  fTransition->SetPolarizationFlag(param->CorrelatedGamma());
105  fTransition->SetVerbose(fVerbose);
106 }
107 
108 G4Fragment*
110 {
111  if(fRDM) { fSampleTime = false; }
112  else { fSampleTime = true; }
113 
114  G4Fragment* gamma = GenerateGamma(nucleus);
115  if(fVerbose > 0) {
116  G4cout << "G4PhotonEvaporation::EmittedFragment: RDM= " << fRDM << G4endl;
117  if(gamma) { G4cout << *gamma << G4endl; }
118  G4cout << " Residual: " << *nucleus << G4endl;
119  }
120  return gamma;
121 }
122 
125 {
126  if(fVerbose > 0) {
127  G4cout << "G4PhotonEvaporation::BreakItUp" << G4endl;
128  }
129  G4Fragment* aNucleus = new G4Fragment(nucleus);
130  G4FragmentVector* products = new G4FragmentVector();
131  BreakUpChain(products, aNucleus);
132  products->push_back(aNucleus);
133  return products;
134 }
135 
137  G4Fragment* nucleus)
138 {
139  if(fVerbose > 0) {
140  G4cout << "G4PhotonEvaporation::BreakUpChain RDM= " << fRDM << " "
141  << *nucleus << G4endl;
142  }
143  G4Fragment* gamma = nullptr;
144  if(fRDM) { fSampleTime = false; }
145  else { fSampleTime = true; }
146 
147  do {
148  gamma = GenerateGamma(nucleus);
149  if(gamma) {
150  products->push_back(gamma);
151  if(fVerbose > 0) {
152  G4cout << "G4PhotonEvaporation::BreakUpChain: "
153  << *gamma << G4endl;
154  G4cout << " Residual: " << *nucleus << G4endl;
155  }
156  // for next decays in the chain always sample time
157  fSampleTime = true;
158  }
159  // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
160  } while(gamma);
161  return false;
162 }
163 
164 G4double
166 {
167  if(!isInitialised) { Initialise(); }
168  fProbability = 0.0;
169  fExcEnergy = nucleus->GetExcitationEnergy();
170  G4int Z = nucleus->GetZ_asInt();
171  G4int A = nucleus->GetA_asInt();
172  fCode = 1000*Z + A;
173  if(fVerbose > 1) {
174  G4cout << "G4PhotonEvaporation::GetEmissionProbability: Z="
175  << Z << " A=" << A << " Eexc(MeV)= " << fExcEnergy << G4endl;
176  }
177  // ignore gamma de-excitation for exotic fragments
178  // and for very low excitations
179  if(0 >= Z || 1 >= A || Z == A || Tolerance >= fExcEnergy)
180  { return fProbability; }
181 
182  // ignore gamma de-excitation for highly excited levels
183  if(A >= MAXGRDATA) { A = MAXGRDATA-1; }
184  //G4cout<<" GREnergy= "<< GREnergy[A]<<" GRWidth= "<<GRWidth[A]<<G4endl;
185 
186  static const G4float GREfactor = 5.0f;
187  if(fExcEnergy >= (G4double)(GREfactor*GRWidth[A] + GREnergy[A])) {
188  return fProbability;
189  }
190  // probability computed assuming continium transitions
191  // VI: continium transition are limited only to final states
192  // below Fermi energy (this approach needs further evaluation)
193  G4double emax = std::max(0.0, nucleus->ComputeGroundStateMass(Z, A-1)
195 
196  // max energy level for continues transition
197  emax = std::min(emax, fExcEnergy);
198  const G4double eexcfac = 0.99;
199  if(0.0 == emax || fExcEnergy*eexcfac <= emax) { emax = fExcEnergy*eexcfac; }
200 
201  fStep = emax;
202  static const G4double MaxDeltaEnergy = CLHEP::MeV;
203  fPoints = std::min((G4int)(fStep/MaxDeltaEnergy) + 2, MAXDEPOINT);
204  fStep /= ((G4double)(fPoints - 1));
205  if(fVerbose > 1) {
206  G4cout << "Emax= " << emax << " Npoints= " << fPoints
207  << " Eex= " << fExcEnergy << G4endl;
208  }
209  // integrate probabilities
210  G4double eres = (G4double)GREnergy[A];
211  G4double wres = (G4double)GRWidth[A];
212  G4double eres2= eres*eres;
213  G4double wres2= wres*wres;
214  G4double xsqr = std::sqrt(A*LevelDensity*fExcEnergy);
215 
216  G4double egam = fExcEnergy;
217  G4double gammaE2 = egam*egam;
218  G4double gammaR2 = gammaE2*wres2;
219  G4double egdp2 = gammaE2 - eres2;
220 
221  G4double p0 = G4Exp(-2.0*xsqr)*gammaR2*gammaE2/(egdp2*egdp2 + gammaR2);
222  G4double p1(0.0);
223 
224  for(G4int i=1; i<fPoints; ++i) {
225  egam -= fStep;
226  gammaE2 = egam*egam;
227  gammaR2 = gammaE2*wres2;
228  egdp2 = gammaE2 - eres2;
229  p1 = G4Exp(2.0*(std::sqrt(A*LevelDensity*std::abs(fExcEnergy - egam)) - xsqr))
230  *gammaR2*gammaE2/(egdp2*egdp2 + gammaR2);
231  fProbability += (p1 + p0);
232  fCummProbability[i] = fProbability;
233  //G4cout << "Egamma= " << egam << " Eex= " << fExcEnergy
234  //<< " p0= " << p0 << " p1= " << p1 << " sum= " << fCummProbability[i] <<G4endl;
235  p0 = p1;
236  }
237 
238  static const G4double NormC = 1.25*CLHEP::millibarn
240  fProbability *= fStep*NormC*A;
241  if(fVerbose > 1) { G4cout << "prob= " << fProbability << G4endl; }
242  return fProbability;
243 }
244 
245 G4double
247 {
248  G4double E = energy;
249  InitialiseLevelManager(Z, A);
250  if(fLevelManager) {
251  E = (G4double)fLevelManager->NearestLevelEnergy(energy, fIndex);
252  if(E > fLevelEnergyMax + Tolerance) { E = energy; }
253  }
254  return E;
255 }
256 
258 {
259  InitialiseLevelManager(Z, A);
260  return fLevelEnergyMax;
261 }
262 
263 G4Fragment*
264 G4PhotonEvaporation::GenerateGamma(G4Fragment* nucleus)
265 {
266  if(!isInitialised) { Initialise(); }
267  G4Fragment* result = nullptr;
268  G4double eexc = nucleus->GetExcitationEnergy();
269  if(eexc <= Tolerance) { return result; }
270 
271  InitialiseLevelManager(nucleus->GetZ_asInt(), nucleus->GetA_asInt());
272 
273  G4double time = nucleus->GetCreationTime();
274 
275  G4double efinal = 0.0;
276  G4double ratio = 0.0;
277  vShellNumber = -1;
278  G4int JP1 = 0;
279  G4int JP2 = 0;
280  G4int multiP = 0;
281  G4bool isGamma = true;
282  G4bool isDiscrete = false;
283 
284  const G4NucLevel* level = nullptr;
285  size_t ntrans = 0;
286 
287  if(fLevelManager && eexc <= fLevelEnergyMax + Tolerance) {
288  fIndex = fLevelManager->NearestLevelIndex(eexc, fIndex);
289  if(0 < fIndex) {
290  // for discrete transition
291  level = fLevelManager->GetLevel(fIndex);
292  if(level) {
293  ntrans = level->NumberOfTransitions();
294  JP1 = fLevelManager->SpinTwo(fIndex);
295  if(ntrans > 0) { isDiscrete = true; }
296  else if(fLevelManager->FloatingLevel(fIndex) > 0) {
297  --fIndex;
298  level = fLevelManager->GetLevel(fIndex);
299  ntrans = level->NumberOfTransitions();
300  JP1 = fLevelManager->SpinTwo(fIndex);
301  if(ntrans > 0) { isDiscrete = true; }
302  }
303  }
304  }
305  }
306  if(fVerbose > 1) {
307  G4int prec = G4cout.precision(4);
308  G4cout << "GenerateGamma: Z= " << nucleus->GetZ_asInt()
309  << " A= " << nucleus->GetA_asInt()
310  << " Exc= " << eexc << " Emax= "
311  << fLevelEnergyMax << " idx= " << fIndex
312  << " fCode= " << fCode << " fPoints= " << fPoints
313  << " Ntr= " << ntrans << " discrete: " << isDiscrete
314  << " fProb= " << fProbability << G4endl;
315  G4cout.precision(prec);
316  }
317 
318  // continues part
319  if(!isDiscrete) {
320  // G4cout << "Continues fIndex= " << fIndex << G4endl;
321 
322  // we compare current excitation versus value used for probability
323  // computation and also Z and A used for probability computation
324  if(fCode != 1000*theZ + theA || eexc != fExcEnergy) {
325  GetEmissionProbability(nucleus);
326  }
327  if(fProbability == 0.0) { return result; }
328  G4double y = fCummProbability[fPoints-1]*G4UniformRand();
329  for(G4int i=1; i<fPoints; ++i) {
330  //G4cout << "y= " << y << " cummProb= " << fCummProbability[i] << G4endl;
331  if(y <= fCummProbability[i]) {
332  efinal = fStep*((i - 1) + (y - fCummProbability[i-1])
333  /(fCummProbability[i] - fCummProbability[i-1]));
334  break;
335  }
336  }
337  // final discrete level
338  if(fLevelManager && efinal <= fLevelEnergyMax + Tolerance) {
339  //G4cout << "Efinal= " << efinal << " idx= " << fIndex << G4endl;
340  fIndex = fLevelManager->NearestLevelIndex(efinal, fIndex);
341  efinal = (G4double)fLevelManager->LevelEnergy(fIndex);
342  // protection
343  if(efinal >= eexc && 0 < fIndex) {
344  --fIndex;
345  efinal = (G4double)fLevelManager->LevelEnergy(fIndex);
346  }
347  nucleus->SetFloatingLevelNumber(fLevelManager->FloatingLevel(fIndex));
348  }
349  if(fVerbose > 1) {
350  G4cout << "Continues emission efinal(MeV)= " << efinal << G4endl;
351  }
352  //discrete part
353  } else {
354  if(fVerbose > 1) {
355  G4cout << "Discrete emission from level Index= " << fIndex
356  << " Elevel= " << fLevelManager->LevelEnergy(fIndex)
357  << " RDM= " << fRDM << " ICM= " << fICM << G4endl;
358  }
359  if(0 == fIndex || !level) { return result; }
360 
361  G4double ltime = 0.0;
362  if(fSampleTime) {
363  ltime = (G4double)fLevelManager->LifeTime(fIndex);
364  if(!fRDM && ltime >= fMaxLifeTime) { return result; }
365  }
366 
367  size_t idx = 0;
368  if(1 < ntrans) {
369  idx = level->SampleGammaTransition(G4UniformRand());
370  }
371  if(fVerbose > 1) {
372  G4cout << "Ntrans= " << ntrans << " idx= " << idx
373  << " ICM= " << fICM << G4endl;
374  }
375  G4double prob = (G4double)level->GammaProbability(idx);
376  // prob = 0 means that there is only internal conversion
377  // if((fICM && prob < 1.0) || prob == 0.0) {
378  if(prob < 1.0) {
379  G4double rndm = G4UniformRand();
380  if(rndm > prob) {
381  isGamma = false;
382  rndm = (rndm - prob)/(1.0 - prob);
383  vShellNumber = level->SampleShell(idx, rndm);
384  }
385  }
386  // it is discrete transition with possible gamma correlation
387  isDiscrete = true;
388  ratio = level->MixingRatio(idx);
389  multiP = level->TransitionType(idx);
390  fIndex = level->FinalExcitationIndex(idx);
391  JP2 = fLevelManager->SpinTwo(fIndex);
392 
393  // final energy and time
394  efinal = (G4double)fLevelManager->LevelEnergy(fIndex);
395  if(fSampleTime && ltime > 0.0) {
396  time -= ltime*G4Log(G4UniformRand());
397  }
398  nucleus->SetFloatingLevelNumber(fLevelManager->FloatingLevel(fIndex));
399  }
400  // protection for floating levels
401  if(std::abs(efinal - eexc) <= Tolerance) { return result; }
402 
403  result = fTransition->SampleTransition(nucleus, efinal, ratio, JP1,
404  JP2, multiP, vShellNumber,
405  isDiscrete, isGamma);
406  if(result) { result->SetCreationTime(time); }
407  nucleus->SetCreationTime(time);
408 
409  if(fVerbose > 1) {
410  G4cout << "Final level E= " << efinal << " time= " << time
411  << " idxFinal= " << fIndex << " isDiscrete: " << isDiscrete
412  << " isGamma: " << isGamma << " multiP= " << multiP
413  << " shell= " << vShellNumber << G4endl;
414  }
415  return result;
416 }
417 
419 {}
420 
422 {
423  if(p != fTransition) {
424  delete fTransition;
425  fTransition = p;
426  }
427 }
428 
430 {
431  fICM = val;
432 }
433 
435 {
436  fRDM = val;
437 }
438 
G4double G4ParticleHPJENDLHEData::G4double result
virtual void SetICM(G4bool)
G4float GammaProbability(size_t idx) const
Definition: G4NucLevel.hh:142
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
void SetMaxHalfLife(G4double)
G4FragmentVector * BreakItUp(const G4Fragment &theNucleus)
G4int SpinTwo(size_t i) const
const G4int MAXGRDATA
virtual G4bool BreakUpChain(G4FragmentVector *theResult, G4Fragment *theNucleus) final
const G4int MAXDEPOINT
Definition: G4Pow.hh:56
virtual void RDMForced(G4bool)
static constexpr double pi2
Definition: SystemOfUnits.h:57
static constexpr double millibarn
Definition: SystemOfUnits.h:86
static constexpr double keV
static constexpr double hbarc
G4PhotonEvaporation(G4GammaTransition *ptr=nullptr)
G4float LifeTime(size_t i) const
const char * p
Definition: xmltok.h:285
float G4float
Definition: G4Types.hh:77
virtual G4Fragment * EmittedFragment(G4Fragment *theNucleus) final
virtual void Initialise() final
G4double GetLevelDensity() const
int G4int
Definition: G4Types.hh:78
const G4NucLevel * GetLevel(size_t i) const
G4float MixingRatio(size_t idx) const
Definition: G4NucLevel.hh:134
G4double GetMaxLifeTime() const
void SetPolarizationFlag(G4bool val)
static const double prec
Definition: RanecuEngine.cc:58
virtual G4double GetEmissionProbability(G4Fragment *theNucleus) final
G4float LevelEnergy(size_t i) const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
G4int GetA_asInt() const
Definition: G4Fragment.hh:266
G4double GetCreationTime() const
Definition: G4Fragment.hh:448
bool G4bool
Definition: G4Types.hh:79
static constexpr double neutron_mass_c2
static constexpr double MeV
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:63
G4DeexPrecoParameters * GetParameters()
void SetGammaTransition(G4GammaTransition *)
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:288
size_t SampleGammaTransition(G4double rndm) const
Definition: G4NucLevel.hh:166
G4float NearestLevelEnergy(G4double energy, size_t index=0) const
size_t NearestLevelIndex(G4double energy, size_t index=0) const
static const G4double emax
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
virtual G4Fragment * SampleTransition(G4Fragment *nucleus, G4double newExcEnergy, G4double mpRatio, G4int JP1, G4int JP2, G4int MP, G4int shell, G4bool isDiscrete, G4bool isGamma)
void SetCreationTime(G4double time)
Definition: G4Fragment.hh:453
size_t NumberOfTransitions() const
Definition: G4NucLevel.hh:108
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double energy(const ThreeVector &p, const G4double m)
G4int FloatingLevel(size_t i) const
G4int TransitionType(size_t idx) const
Definition: G4NucLevel.hh:121
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4int GetZ_asInt() const
Definition: G4Fragment.hh:271
void SetFloatingLevelNumber(G4int value)
Definition: G4Fragment.hh:432
G4int SampleShell(size_t idx, G4double rndm) const
Definition: G4NucLevel.hh:176
#define G4endl
Definition: G4ios.hh:61
G4double powZ(G4int Z, G4double y) const
Definition: G4Pow.hh:254
double G4double
Definition: G4Types.hh:76
size_t FinalExcitationIndex(size_t idx) const
Definition: G4NucLevel.hh:113
#define DBL_MAX
Definition: templates.hh:83
virtual G4double GetUpperLevelEnergy(G4int Z, G4int A) final
G4double ComputeGroundStateMass(G4int Z, G4int A) const
Definition: G4Fragment.hh:256
static G4NuclearLevelData * GetInstance()
G4double GetMinExcitation() const
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:283
virtual G4double GetFinalLevelEnergy(G4int Z, G4int A, G4double energy) final
void SetVerbose(G4int val)