Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4GEMProbability.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$
27 //
28 //---------------------------------------------------------------------
29 //
30 // Geant4 class G4GEMProbability
31 //
32 //
33 // Hadronic Process: Nuclear De-excitations
34 // by V. Lara (Sept 2001)
35 //
36 //
37 // Hadronic Process: Nuclear De-excitations
38 // by V. Lara (Sept 2001)
39 //
40 // J. M. Quesada : several fixes in total GEM width
41 // J. M. Quesada 14/07/2009 bug fixed in total emission width (hbarc)
42 // J. M. Quesada (September 2009) several fixes:
43 // -level density parameter of residual calculated at its right excitation energy.
44 // -InitialLeveldensity calculated according to the right conditions of the
45 // initial excited nucleus.
46 // J. M. Quesada 19/04/2010 fix in emission probability calculation.
47 // V.Ivanchenko 20/04/2010 added usage of G4Pow and use more safe computation
48 // V.Ivanchenko 18/05/2010 trying to speedup the most slow method
49 // by usage of G4Pow, integer Z and A; moved constructor,
50 // destructor and virtual functions to source
51 
52 #include "G4GEMProbability.hh"
53 #include "G4PairingCorrection.hh"
54 #include "G4PhysicalConstants.hh"
55 #include "G4SystemOfUnits.hh"
56 
57 G4GEMProbability:: G4GEMProbability(G4int anA, G4int aZ, G4double aSpin) :
58  theA(anA), theZ(aZ), Spin(aSpin), theCoulombBarrierPtr(0),
59  Normalization(1.0)
60 {
61  theEvapLDPptr = new G4EvaporationLevelDensityParameter;
62  fG4pow = G4Pow::GetInstance();
63  fPlanck= CLHEP::hbar_Planck*fG4pow->logZ(2);
65 }
66 
68 {
69  delete theEvapLDPptr;
70 }
71 
73  G4double MaximalKineticEnergy)
74 {
75  G4double probability = 0.0;
76 
77  if (MaximalKineticEnergy > 0.0 && fragment.GetExcitationEnergy() > 0.0) {
78  G4double CoulombBarrier = GetCoulombBarrier(fragment);
79 
80  probability =
81  CalcProbability(fragment,MaximalKineticEnergy,CoulombBarrier);
82 
83  // Next there is a loop over excited states for this channel
84  // summing probabilities
85  size_t nn = ExcitEnergies.size();
86  if (0 < nn) {
87  G4double SavedSpin = Spin;
88  for (size_t i = 0; i <nn; ++i) {
89  Spin = ExcitSpins[i];
90  // substract excitation energies
91  G4double Tmax = MaximalKineticEnergy - ExcitEnergies[i];
92  if (Tmax > 0.0) {
93  G4double width = CalcProbability(fragment,Tmax,CoulombBarrier);
94  //JMQ April 2010 added condition to prevent reported crash
95  // update probability
96  if (width > 0. && fPlanck < width*ExcitLifetimes[i]) {
97  probability += width;
98  }
99  }
100  }
101  // Restore Spin
102  Spin = SavedSpin;
103  }
104  }
105  Normalization = probability;
106  return probability;
107 }
108 
109 
110 G4double G4GEMProbability::CalcProbability(const G4Fragment & fragment,
111  G4double MaximalKineticEnergy,
112  G4double V)
113 
114 // Calculate integrated probability (width) for evaporation channel
115 {
116  G4int A = fragment.GetA_asInt();
117  G4int Z = fragment.GetZ_asInt();
118 
119  G4int ResidualA = A - theA;
120  G4int ResidualZ = Z - theZ;
121  G4double U = fragment.GetExcitationEnergy();
122 
123  G4double NuclearMass = fragment.ComputeGroundStateMass(theZ, theA);
124 
125  G4double Alpha = CalcAlphaParam(fragment);
126  G4double Beta = CalcBetaParam(fragment);
127 
128  // ***RESIDUAL***
129  //JMQ (September 2009) the following quantities refer to the RESIDUAL:
130 
131  G4double delta0 = fPairCorr->GetPairingCorrection(ResidualA, ResidualZ);
132 
133  G4double a = theEvapLDPptr->
134  LevelDensityParameter(ResidualA,ResidualZ,MaximalKineticEnergy+V-delta0);
135  G4double Ux = (2.5 + 150.0/G4double(ResidualA))*MeV;
136  G4double Ex = Ux + delta0;
137  G4double T = 1.0/(std::sqrt(a/Ux) - 1.5/Ux);
138  //JMQ fixed bug in units
139  G4double E0 = Ex - T*(std::log(T/MeV) - std::log(a*MeV)/4.0
140  - 1.25*std::log(Ux/MeV) + 2.0*std::sqrt(a*Ux));
141  // ***end RESIDUAL ***
142 
143  // ***PARENT***
144  //JMQ (September 2009) the following quantities refer to the PARENT:
145 
146  G4double deltaCN = fPairCorr->GetPairingCorrection(A, Z);
147  G4double aCN = theEvapLDPptr->LevelDensityParameter(A, Z, U-deltaCN);
148  G4double UxCN = (2.5 + 150.0/G4double(A))*MeV;
149  G4double ExCN = UxCN + deltaCN;
150  G4double TCN = 1.0/(std::sqrt(aCN/UxCN) - 1.5/UxCN);
151  // ***end PARENT***
152 
153  G4double Width;
154  G4double InitialLevelDensity;
155  G4double t = MaximalKineticEnergy/T;
156  if ( MaximalKineticEnergy < Ex ) {
157  //JMQ 190709 bug in I1 fixed (T was missing)
158  Width = (I1(t,t)*T + (Beta+V)*I0(t))/std::exp(E0/T);
159  //JMQ 160909 fix: InitialLevelDensity has been taken away
160  //(different conditions for initial CN..)
161  } else {
162 
163  //VI minor speedup
164  G4double expE0T = std::exp(E0/T);
165  const G4double sqrt2 = std::sqrt(2.0);
166 
167  G4double tx = Ex/T;
168  G4double s0 = 2.0*std::sqrt(a*(MaximalKineticEnergy-delta0));
169  G4double sx = 2.0*std::sqrt(a*(Ex-delta0));
170  Width = I1(t,tx)*T/expE0T + I3(s0,sx)*std::exp(s0)/(sqrt2*a);
171  // For charged particles (Beta+V) = 0 beacuse Beta = -V
172  if (theZ == 0) {
173  Width += (Beta+V)*(I0(tx)/expE0T + 2.0*sqrt2*I2(s0,sx)*std::exp(s0));
174  }
175  }
176 
177  //JMQ 14/07/2009 BIG BUG : NuclearMass is in MeV => hbarc instead of hbar_planck must be used
178  // G4double g = (2.0*Spin+1.0)*NuclearMass/(pi2* hbar_Planck*hbar_Planck);
179  G4double gg = (2.0*Spin+1.0)*NuclearMass/(pi2* hbarc*hbarc);
180 
181  //JMQ 190709 fix on Rb and geometrical cross sections according to Furihata's paper
182  // (JAERI-Data/Code 2001-105, p6)
183  // G4double RN = 0.0;
184  G4double Rb = 0.0;
185  if (theA > 4)
186  {
187  G4double Ad = fG4pow->Z13(ResidualA);
188  G4double Aj = fG4pow->Z13(theA);
189  Rb = 1.12*(Aj + Ad) - 0.86*((Aj+Ad)/(Aj*Ad))+2.85;
190  Rb *= fermi;
191  }
192  else if (theA>1)
193  {
194  G4double Ad = fG4pow->Z13(ResidualA);
195  G4double Aj = fG4pow->Z13(theA);
196  Rb=1.5*(Aj+Ad)*fermi;
197  }
198  else
199  {
200  G4double Ad = fG4pow->Z13(ResidualA);
201  Rb = 1.5*Ad*fermi;
202  }
203  // G4double GeometricalXS = pi*RN*RN*std::pow(ResidualA,2./3.);
204  G4double GeometricalXS = pi*Rb*Rb;
205  //end of JMQ fix on Rb by 190709
206 
207  //JMQ 160909 fix: initial level density must be calculated according to the
208  // conditions at the initial compound nucleus
209  // (it has been removed from previous "if" for the residual)
210 
211  if ( U < ExCN )
212  {
213  //JMQ fixed bug in units
214  //VI moved the computation here
215  G4double E0CN = ExCN - TCN*(std::log(TCN/MeV) - std::log(aCN*MeV)/4.0
216  - 1.25*std::log(UxCN/MeV)
217  + 2.0*std::sqrt(aCN*UxCN));
218  InitialLevelDensity = (pi/12.0)*std::exp((U-E0CN)/TCN)/TCN;
219  }
220  else
221  {
222  //VI speedup
223  G4double x = U-deltaCN;
224  G4double x1 = std::sqrt(aCN*x);
225  InitialLevelDensity = (pi/12.0)*std::exp(2*x1)/(x*std::sqrt(x1));
226  }
227 
228  //JMQ 190709 BUG : pi instead of sqrt(pi) must be here according
229  // to Furihata's report:
230  Width *= pi*gg*GeometricalXS*Alpha/(12.0*InitialLevelDensity);
231 
232  return Width;
233 }
234 
235 G4double G4GEMProbability::I3(G4double s0, G4double sx)
236 {
237  G4double s2 = s0*s0;
238  G4double sx2 = sx*sx;
239  G4double S = 1.0/std::sqrt(s0);
240  G4double S2 = S*S;
241  G4double Sx = 1.0/std::sqrt(sx);
242  G4double Sx2 = Sx*Sx;
243 
244  G4double p1 = S *(2.0 + S2 *( 4.0 + S2 *( 13.5 + S2 *( 60.0 + S2 * 325.125 ))));
245  G4double p2 = Sx*Sx2 *(
246  (s2-sx2) + Sx2 *(
247  (1.5*s2+0.5*sx2) + Sx2 *(
248  (3.75*s2+0.25*sx2) + Sx2 *(
249  (12.875*s2+0.625*sx2) + Sx2 *(
250  (59.0625*s2+0.9375*sx2) + Sx2 *(324.8*s2+3.28*sx2))))));
251 
252  p2 *= std::exp(sx-s0);
253 
254  return p1-p2;
255 }