Geant4  10.02
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 94676 2015-12-02 09:51:20Z gunter $
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"
54 #include "G4Pow.hh"
55 #include <CLHEP/Units/SystemOfUnits.h>
56 #include <CLHEP/Units/PhysicalConstants.h>
57 
60 const G4float GREfactor = 5.0f;
61 const G4float GRWfactor = 0.3f;
65 const G4double NormC = 2.5*CLHEP::millibarn/(CLHEP::pi2*CLHEP::hbarc*CLHEP::hbarc);
67 
69  : fLevelManager(nullptr), fTransition(p), fVerbose(0), fPoints(0),
70  vShellNumber(-1), fIndex(0), fTimeLimit(DBL_MAX), fMaxLifeTime(DBL_MAX),
71  fICM(false), fRDM(false), fSampleTime(true)
72 {
74  if(!fTransition) {
75  char* en = getenv("G4UseNuclearPolarization");
76  if(en) { fTransition = new G4PolarizedGammaTransition(); }
77  else { fTransition = new G4GammaTransition(); }
78  }
79 
80  char* env = getenv("G4AddTimeLimitToPhotonEvaporation");
81  if(env) { fTimeLimit = 1.e-16*CLHEP::second; }
82 
83  theA = theZ = fCode = 0;
85 
86  for(G4int i=0; i<MAXDEPOINT; ++i) { fCummProbability[i] = 0.0; }
87  if(0.0f == GREnergy[1]) {
88  G4Pow* g4pow = G4Pow::GetInstance();
89  for (G4int A=1; A<MAXGRDATA; ++A) {
90  GREnergy[A] = (G4float)(40.3*CLHEP::MeV/g4pow->powZ(A,0.2));
92  }
93  }
94 }
95 
97 {
98  delete fTransition;
99 }
100 
101 G4Fragment*
103 {
104  if(fRDM) { fSampleTime = false; }
105  else { fSampleTime = true; }
106 
107  G4Fragment* gamma = GenerateGamma(nucleus);
108  if(fVerbose > 0) {
109  G4cout << "G4PhotonEvaporation::EmittedFragment: RDM= " << fRDM << G4endl;
110  if(gamma) { G4cout << *gamma << G4endl; }
111  G4cout << " Residual: " << *nucleus << G4endl;
112  }
113  return gamma;
114 }
115 
118 {
119  //G4cout << "G4PhotonEvaporation::BreakUpFragment" << G4endl;
120  G4FragmentVector* products = new G4FragmentVector();
121  BreakUpChain(products, nucleus);
122  return products;
123 }
124 
127 {
128  //G4cout << "G4PhotonEvaporation::BreakUp" << G4endl;
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 
138 {
139  //G4cout << "G4PhotonEvaporation::BreakItUp" << G4endl;
140  G4Fragment* aNucleus = new G4Fragment(nucleus);
141  G4FragmentVector* products = new G4FragmentVector();
142  BreakUpChain(products, aNucleus);
143  products->push_back(aNucleus);
144  return products;
145 }
146 
148  G4Fragment* nucleus)
149 {
150  if(fVerbose > 0) {
151  G4cout << "G4PhotonEvaporation::BreakUpChain RDM= " << fRDM << " "
152  << *nucleus << G4endl;
153  }
154  G4Fragment* gamma = 0;
155  if(fRDM) { fSampleTime = false; }
156  else { fSampleTime = true; }
157 
158  do {
159  gamma = GenerateGamma(nucleus);
160  if(gamma) {
161  products->push_back(gamma);
162  if(fVerbose > 0) {
163  G4cout << "G4PhotonEvaporation::BreakUpChain: "
164  << *gamma << G4endl;
165  G4cout << " Residual: " << *nucleus << G4endl;
166  }
167  // for next decays in the chain always sample time
168  fSampleTime = true;
169  }
170  // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
171  } while(gamma);
172  return false;
173 }
174 
175 G4double
177 {
178  fProbability = 0.0;
179  fExcEnergy = nucleus->GetExcitationEnergy();
180  G4int Z = nucleus->GetZ_asInt();
181  G4int A = nucleus->GetA_asInt();
182  fCode = 1000*Z + A;
183  if(fVerbose > 1) {
184  G4cout << "G4PhotonEvaporation::GetEmissionProbability: Z="
185  << Z << " A=" << A << " Eexc(MeV)= " << fExcEnergy << G4endl;
186  }
187  // ignore gamma de-excitation for exotic fragments
188  // and for very low excitations
189  if(0 >= Z || 1 >= A || Z == A || Tolerance >= fExcEnergy)
190  { return fProbability; }
191 
192  // ignore gamma de-excitation for highly excited levels
193  if(A >= MAXGRDATA) { A = MAXGRDATA-1; }
194  //G4cout << " GREnergy= " << GREnergy[A] << " GRWidth= " << GRWidth[A] << G4endl;
195  if(fExcEnergy >= (G4double)(GREfactor*GRWidth[A] + GREnergy[A])) {
196  return fProbability;
197  }
198 
199  // probability computed assuming continium transitions
200  // VI: continium transition are limited only to final states
201  // below Fermi energy (this approach needs further evaluation)
202  fFermiEnergy = std::max(0.0, nucleus->ComputeGroundStateMass(Z, A-1)
203  + CLHEP::neutron_mass_c2 - nucleus->GetGroundStateMass());
204 
205  // max energy level for continues transition
207  const G4double eexcfac = 0.99;
208  if(0.0 == emax || fExcEnergy*eexcfac <= emax) { emax = fExcEnergy*eexcfac; }
209 
210  fStep = emax;
212  fStep /= ((G4double)(fPoints - 1));
213  if(fVerbose > 1) {
214  G4cout << "Emax= " << emax << " Npoints= " << fPoints
215  <<" Efermi= " << fFermiEnergy << " Eex= " << fExcEnergy << G4endl;
216  }
217  // integrate probabilities
218  G4double eres = (G4double)GREnergy[A];
219  G4double wres = (G4double)GRWidth[A];
220  G4double eres2= eres*eres;
221  G4double wres2= wres*wres;
222  G4double xsqr = std::sqrt(A*LevelDensity*fExcEnergy);
223 
224  G4double egam = fExcEnergy - emax;
225  G4double gammaE2 = egam*egam;
226  G4double gammaR2 = gammaE2*wres*wres;
227  G4double egdp2 = gammaE2 - eres*eres;
228 
229  G4double p0 = G4Exp(2*(std::sqrt(A*LevelDensity*emax) - xsqr))
230  *gammaR2*gammaE2/(egdp2*egdp2 + gammaR2);
231  G4double p1(0.0);
232 
233  for(G4int i=1; i<fPoints; ++i) {
234  egam += fStep;
235  gammaE2 = egam*egam;
236  gammaR2 = gammaE2*wres2;
237  egdp2 = gammaE2 - eres2;
238  //G4cout << "Egamma= " << egam << " Eex= " << fExcEnergy
239  //<< " p0= " << p0 << " p1= " << p1 << G4endl;
240  p1 = G4Exp(2*(std::sqrt(A*LevelDensity*std::abs(fExcEnergy - egam)) - xsqr))
241  *gammaR2*gammaE2/(egdp2*egdp2 + gammaR2);
242  fProbability += (p1 + p0);
244  p0 = p1;
245  }
246  fProbability *= 0.5*fStep*NormC*A;
247  if(fVerbose > 1) { G4cout << "prob= " << fProbability << G4endl; }
248  return fProbability;
249 }
250 
251 G4double
253 {
254  G4double E = energy;
256  if(fLevelManager) {
258  if(E > fLevelEnergyMax + Tolerance) { E = energy; }
259  }
260  return E;
261 }
262 
264 {
266  return fLevelEnergyMax;
267 }
268 
269 G4Fragment*
271 {
272  G4Fragment* result = 0;
273  G4double eexc = nucleus->GetExcitationEnergy();
274  if(eexc < MinGammaEnergy) { return result; }
275 
276  InitialiseLevelManager(nucleus->GetZ_asInt(), nucleus->GetA_asInt());
277 
278  G4double time = nucleus->GetCreationTime();
279 
280  G4double efinal = 0.0;
281  vShellNumber = -1;
282  size_t shell = 0;
283  G4int deltaS = 1;
284  G4bool isGamma = true;
285  G4bool isLongLived = false;
286  G4bool isX = false;
287  G4bool icm = fICM;
288 
289  if(fVerbose > 1) {
290  G4cout << "GenerateGamma: Exc= " << eexc << " Emax= "
291  << fLevelEnergyMax << G4endl;
292  }
293  // continues part
294  if(!fLevelManager || eexc > fLevelEnergyMax + Tolerance) {
295  //G4cout << "Continues fPoints= " << fPoints << " " << fLevelManager << G4endl;
296 
297  // we compare current excitation versus value used for probability
298  // computation and also Z and A used for probability computation
299  if(fCode != 1000*theZ + theA && eexc != fExcEnergy) {
300  GetEmissionProbability(nucleus);
301  }
302  if(fProbability == 0.0) { return result; }
304  for(G4int i=1; i<fPoints; ++i) {
305  //G4cout << "y= " << y << " cummProb= " << fCummProbability[i] << G4endl;
306  if(y <= fCummProbability[i]) {
307  efinal = fStep*((i - 1) + (y - fCummProbability[i-1])
308  /(fCummProbability[i] - fCummProbability[i-1]));
309  break;
310  }
311  }
312  // final discrete level
313  if(fLevelManager && efinal <= fLevelEnergyMax + Tolerance) {
314  //G4cout << "Efinal= " << efinal << " idx= " << fIndex << G4endl;
317  }
318  //discrete part
319  } else {
321  if(fVerbose > 1) {
322  G4cout << "Discrete emission from level Index= " << fIndex
323  << " Elevel= " << fLevelManager->LevelEnergy(fIndex)
324  << " ICM= " << fICM << G4endl;
325  }
326  if(0 == fIndex) { return result; }
327  G4double ltime = 0.0;
328  if(fICM) { ltime = (G4double)fLevelManager->LifeTime(fIndex); }
329  else { ltime = (G4double)fLevelManager->LifeTimeGamma(fIndex); }
330 
331  if(ltime >= fMaxLifeTime) { return result; }
332  if(ltime > fTimeLimit) {
333  icm = true;
334  isLongLived = true;
335  }
336  const G4NucLevel* level = fLevelManager->GetLevel(fIndex);
337  size_t ntrans = level->NumberOfTransitions();
338  size_t idx = 0;
339  isX = level->IsXLevel();
340  //G4cout << "Ntrans= " << ntrans << " idx= " << idx << " isX= " << isX
341  // << " icm= " << icm << G4endl;
342  if(!isX) {
343  if(1 < ntrans) {
344  G4double rndm = G4UniformRand();
345  if(icm) { idx = level->SampleGammaETransition(rndm); }
346  else { idx = level->SampleGammaTransition(rndm); }
347  //G4cout << "Sampled idx= " << idx << " rndm= " << rndm << G4endl;
348  }
349  if(icm) {
350  G4double rndm = G4UniformRand();
351  G4double prob = level->GammaProbability(idx);
352  if(rndm > prob) {
353  rndm = (rndm - prob)/(1.0 - prob);
354  shell = level->SampleShell(idx, rndm);
355  vShellNumber = shell;
356  isGamma = false;
357  }
358  }
359  }
360  efinal = (G4double)level->FinalExcitationEnergy(idx);
361  fIndex = idx;
362 
363  if(fSampleTime && ltime > 0.0) {
364  time -= ltime*G4Log(G4UniformRand());
365  }
366  }
367  // no gamma emission for X-level
368  if(isX) {
369  G4ThreeVector v = nucleus->GetMomentum().vect().unit();
370  G4double mass = nucleus->GetGroundStateMass() + efinal;
371  G4double e = std::max(mass,nucleus->GetMomentum().e());
372  G4double mom = std::sqrt((e - mass)*(e + mass));
373  v *= mom;
374  nucleus->SetMomentum(G4LorentzVector(v.x(),v.y(),v.z(),efinal+mass));
375  // normal case
376  } else {
377  result = fTransition->SampleTransition(nucleus, efinal,
378  deltaS, shell,
379  isGamma, isLongLived);
380  if(result) { result->SetCreationTime(time); }
381  }
382  nucleus->SetCreationTime(time);
383 
384  if(fVerbose > 1) {
385  G4cout << "Final level E= " << efinal << " time= " << time
386  << " idx= " << fIndex << " isX " << isX
387  << " isGamma: " << isGamma << " isLongLived: " << isLongLived
388  << " deltaS= " << deltaS << " shell= " << shell << G4endl;
389  }
390  return result;
391 }
392 
394 {
395  static const G4double tfact = G4Pow::GetInstance()->logZ(2);
396  fMaxLifeTime = val/tfact;
397 }
398 
400 {
401  if(p != fTransition) {
402  delete fTransition;
403  fTransition = p;
404  }
405 }
406 
408 {
409  fICM = val;
410 }
411 
413 {
414  fRDM = val;
415 }
416 
virtual void SetICM(G4bool)
virtual G4bool BreakUpChain(G4FragmentVector *theResult, G4Fragment *theNucleus)
G4float GammaProbability(size_t idx) const
Definition: G4NucLevel.hh:115
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
void SetMaxHalfLife(G4double)
virtual G4FragmentVector * BreakItUp(const G4Fragment &theNucleus)
static const double MeV
Definition: G4SIunits.hh:211
G4GammaTransition * fTransition
const G4int MAXGRDATA
static G4float GREnergy[MAXGRDATA]
virtual G4double GetUpperLevelEnergy(G4int Z, G4int A)
CLHEP::Hep3Vector G4ThreeVector
const G4int MAXDEPOINT
const G4double Tolerance
Definition: G4Pow.hh:56
virtual void RDMForced(G4bool)
G4float LifeTime(size_t i) const
float G4float
Definition: G4Types.hh:77
virtual G4double GetFinalLevelEnergy(G4int Z, G4int A, G4double energy)
static G4float GRWidth[MAXGRDATA]
G4double fCummProbability[MAXDEPOINT]
G4PhotonEvaporation(G4GammaTransition *ptr=0)
const G4LevelManager * fLevelManager
static const double pi2
Definition: G4SIunits.hh:77
int G4int
Definition: G4Types.hh:78
const G4NucLevel * GetLevel(size_t i) const
virtual G4FragmentVector * BreakUp(const G4Fragment &theNucleus)
G4float FinalExcitationEnergy(size_t idx) const
Definition: G4NucLevel.hh:109
G4double logZ(G4int Z) const
Definition: G4Pow.hh:166
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:251
G4double GetCreationTime() const
Definition: G4Fragment.hh:405
G4bool IsXLevel() const
Definition: G4NucLevel.hh:104
bool G4bool
Definition: G4Types.hh:79
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:284
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:289
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:63
static const double second
Definition: G4SIunits.hh:156
virtual G4double GetEmissionProbability(G4Fragment *theNucleus)
void SetGammaTransition(G4GammaTransition *)
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:273
size_t SampleGammaTransition(G4double rndm) const
Definition: G4NucLevel.hh:127
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
const G4double MinGammaEnergy
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
const G4float GRWfactor
void SetCreationTime(G4double time)
Definition: G4Fragment.hh:410
size_t NumberOfTransitions() const
Definition: G4NucLevel.hh:99
virtual G4Fragment * SampleTransition(G4Fragment *nucleus, G4double newExcEnergy, G4int deltaS, size_t shell, G4bool isGamma, G4bool isLongLived)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double energy(const ThreeVector &p, const G4double m)
G4Fragment * GenerateGamma(G4Fragment *nucleus)
const G4double MaxDeltaEnergy
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4NuclearLevelData * fNuclearLevelData
size_t SampleShell(size_t idx, G4double rndm) const
Definition: G4NucLevel.hh:147
static const double millibarn
Definition: G4SIunits.hh:105
G4int GetZ_asInt() const
Definition: G4Fragment.hh:256
const G4float GREfactor
#define G4endl
Definition: G4ios.hh:61
const G4double LevelDensity
virtual G4Fragment * EmittedFragment(G4Fragment *theNucleus)
size_t SampleGammaETransition(G4double rndm) const
Definition: G4NucLevel.hh:137
static const double keV
Definition: G4SIunits.hh:213
virtual G4FragmentVector * BreakUpFragment(G4Fragment *theNucleus)
G4double powZ(G4int Z, G4double y) const
Definition: G4Pow.hh:254
double G4double
Definition: G4Types.hh:76
G4float LifeTimeGamma(size_t i) const
void InitialiseLevelManager(G4int Z, G4int A)
#define DBL_MAX
Definition: templates.hh:83
G4double ComputeGroundStateMass(G4int Z, G4int A) const
Definition: G4Fragment.hh:241
static G4NuclearLevelData * GetInstance()
const G4double NormC
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:268
CLHEP::HepLorentzVector G4LorentzVector