source: trunk/source/processes/hadronic/models/de_excitation/evaporation/src/G4EvaporationProbability.cc @ 1007

Last change on this file since 1007 was 1007, checked in by garnier, 15 years ago

update to geant4.9.2

File size: 8.6 KB
Line 
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//J.M. Quesada (August2008). Based on:
27//
28// Hadronic Process: Nuclear De-excitations
29// by V. Lara (Oct 1998)
30//
31// Modif (03 September 2008) by J. M. Quesada for external choice of inverse
32// cross section option
33// JMQ (06 September 2008) Also external choices have been added for
34// superimposed Coulomb barrier (if useSICB is set true, by default is false)
35
36#include <iostream>
37using namespace std;
38
39#include "G4EvaporationProbability.hh"
40#include "G4PairingCorrection.hh"
41
42
43
44G4EvaporationProbability::G4EvaporationProbability(const G4EvaporationProbability &) : G4VEmissionProbability()
45{
46    throw G4HadronicException(__FILE__, __LINE__, "G4EvaporationProbability::copy_constructor meant to not be accessable");
47}
48
49
50
51
52const G4EvaporationProbability & G4EvaporationProbability::
53operator=(const G4EvaporationProbability &)
54{
55    throw G4HadronicException(__FILE__, __LINE__, "G4EvaporationProbability::operator= meant to not be accessable");
56    return *this;
57}
58
59
60G4bool G4EvaporationProbability::operator==(const G4EvaporationProbability &) const
61{
62    return false;
63}
64
65G4bool G4EvaporationProbability::operator!=(const G4EvaporationProbability &) const
66{
67    return true;
68}
69 
70G4double G4EvaporationProbability::EmissionProbability(const G4Fragment & fragment, const G4double anEnergy)
71{
72    G4double EmissionProbability = 0.0;
73    G4double MaximalKineticEnergy = anEnergy;
74
75    if (MaximalKineticEnergy > 0.0 && fragment.GetExcitationEnergy() > 0.0) {
76        EmissionProbability = CalculateProbability(fragment, MaximalKineticEnergy);
77
78    }
79    return EmissionProbability;
80}
81
82////////////////////////////////////
83
84// Computes the integrated probability of evaporation channel
85G4double G4EvaporationProbability::CalculateProbability(const G4Fragment & fragment, const G4double MaximalKineticEnergy)
86{
87    G4double ResidualA = fragment.GetA() - theA;
88    G4double ResidualZ = fragment.GetZ() - theZ;
89    G4double U = fragment.GetExcitationEnergy();
90   
91 if (OPTxs==0) {
92
93       
94    G4double NuclearMass = G4ParticleTable::GetParticleTable()->GetIonTable()->GetNucleusMass(theZ,theA);
95
96
97    G4double delta0 = G4PairingCorrection::GetInstance()->GetPairingCorrection(static_cast<G4int>(fragment.GetA()),
98                                                                               static_cast<G4int>(fragment.GetZ()));
99
100    G4double SystemEntropy = 2.0*std::sqrt(theEvapLDPptr->LevelDensityParameter(static_cast<G4int>(fragment.GetA()),
101                                                                           static_cast<G4int>(fragment.GetZ()),U)*
102                                      (U-delta0));
103                                                                 
104
105    G4double RN = 1.5*fermi;
106
107    G4double Alpha = CalcAlphaParam(fragment);
108    G4double Beta = CalcBetaParam(fragment);
109       
110    G4double Rmax = MaximalKineticEnergy;
111    G4double a = theEvapLDPptr->LevelDensityParameter(static_cast<G4int>(ResidualA),
112                                                      static_cast<G4int>(ResidualZ),
113                                                      Rmax);
114    G4double GlobalFactor = static_cast<G4double>(Gamma) * (Alpha/(a*a)) *
115        (NuclearMass*RN*RN*std::pow(ResidualA,2./3.))/
116        (2.*pi* hbar_Planck*hbar_Planck);
117    G4double Term1 = (2.0*Beta*a-3.0)/2.0 + Rmax*a;
118    G4double Term2 = (2.0*Beta*a-3.0)*std::sqrt(Rmax*a) + 2.0*a*Rmax;
119       
120    G4double ExpTerm1 = 0.0;
121    if (SystemEntropy <= 600.0) ExpTerm1 = std::exp(-SystemEntropy);
122       
123    G4double ExpTerm2 = 2.*std::sqrt(a*Rmax) - SystemEntropy;
124    if (ExpTerm2 > 700.0) ExpTerm2 = 700.0;
125    ExpTerm2 = std::exp(ExpTerm2);
126       
127    G4double Width = GlobalFactor*(Term1*ExpTerm1 + Term2*ExpTerm2);
128       
129    return Width;
130             
131 } else if (OPTxs==1 || OPTxs==2 ||OPTxs==3 || OPTxs==4) {
132
133   G4double limit;
134   if (useSICB) limit=theCoulombBarrierptr->GetCoulombBarrier(G4lrint(ResidualA),G4lrint(ResidualZ),U);
135   else limit=0.;
136
137   if (MaximalKineticEnergy <= limit)  return 0.0;
138
139
140// if Coulomb barrier cutoff is superimposed for all cross sections the limit is the Coulomb Barrier
141        G4double LowerLimit= limit;
142
143//  MaximalKineticEnergy: asimptotic value (already accounted for in G4EvaporationChannel)     
144
145  G4double UpperLimit = MaximalKineticEnergy;
146
147
148  G4double Width = IntegrateEmissionProbability(fragment,LowerLimit,UpperLimit);
149
150  return Width;
151 } else{
152   std::ostringstream errOs;
153   errOs << "Bad option for cross sections at evaporation"  <<G4endl;
154   throw G4HadronicException(__FILE__, __LINE__, errOs.str());
155 }
156 
157}
158
159/////////////////////////////////////////////////////////////////////
160
161G4double G4EvaporationProbability::
162IntegrateEmissionProbability(const G4Fragment & fragment, const G4double & Low, const G4double & Up )
163{
164
165  static const G4int N = 10;
166  // 10-Points Gauss-Legendre abcisas and weights
167  static const G4double w[N] = {
168    0.0666713443086881,
169    0.149451349150581,
170    0.219086362515982,
171    0.269266719309996,
172    0.295524224714753,
173    0.295524224714753,
174    0.269266719309996,
175    0.219086362515982,
176    0.149451349150581,
177    0.0666713443086881
178  };
179  static const G4double x[N] = {
180    -0.973906528517172,
181    -0.865063366688985,
182    -0.679409568299024,
183    -0.433395394129247,
184    -0.148874338981631,
185    0.148874338981631,
186    0.433395394129247,
187    0.679409568299024,
188    0.865063366688985,
189    0.973906528517172
190  };
191
192  G4double Total = 0.0;
193
194
195  for (G4int i = 0; i < N; i++) 
196    {
197
198      G4double KineticE = ((Up-Low)*x[i]+(Up+Low))/2.0;
199
200      Total += w[i]*ProbabilityDistributionFunction(fragment, KineticE);
201
202    }
203  Total *= (Up-Low)/2.0;
204  return Total;
205}
206
207
208/////////////////////////////////////////////////////////
209//New method (OPT=1,2,3,4)
210
211G4double G4EvaporationProbability::ProbabilityDistributionFunction( const G4Fragment & fragment, const G4double K)
212{ 
213
214 
215
216
217  G4double ResidualA = fragment.GetA() - theA;
218  G4double ResidualZ = fragment.GetZ() - theZ; 
219  G4double U = fragment.GetExcitationEnergy();
220
221  //        if(K <= theCoulombBarrierptr->GetCoulombBarrier(G4lrint(ResidualA),G4lrint(ResidualZ),U)) return 0.0;   
222
223  G4double delta1 = G4PairingCorrection::GetInstance()->GetPairingCorrection(static_cast<G4int>(ResidualA),static_cast<G4int>(ResidualZ));
224
225 
226  G4double delta0 = G4PairingCorrection::GetInstance()->GetPairingCorrection(static_cast<G4int>(fragment.GetA()),static_cast<G4int>(fragment.GetZ()));
227
228 
229  G4double ParticleMass = G4ParticleTable::GetParticleTable()->GetIonTable()->GetNucleusMass(theZ,theA);
230
231  G4double theSeparationEnergy= G4NucleiProperties::GetMassExcess(static_cast<G4int>(theA),static_cast<G4int>(theZ)) +
232    G4NucleiProperties::GetMassExcess(static_cast<G4int>(ResidualA),static_cast<G4int>(ResidualZ)) -
233    G4NucleiProperties::GetMassExcess(static_cast<G4int>(fragment.GetA()),static_cast<G4int>(fragment.GetZ()));
234
235  G4double a0 = theEvapLDPptr->LevelDensityParameter(static_cast<G4int>(fragment.GetA()),static_cast<G4int>(fragment.GetZ()),U - delta0);
236
237  G4double a1 = theEvapLDPptr->LevelDensityParameter(static_cast<G4int>(ResidualA),static_cast<G4int>(ResidualZ),U - theSeparationEnergy - delta1);
238 
239 
240  G4double E0=U-delta0;
241
242  G4double E1=U-theSeparationEnergy-delta1-K;
243
244  if (E1<0.) return 0.;
245
246
247//Without 1/hbar_Panck remains as a width
248
249  G4double  Prob=Gamma*ParticleMass/((pi*hbar_Planck)*(pi*hbar_Planck)*std::exp(2*std::sqrt(a0*E0)))*K*CrossSection(fragment,K)*std::exp(2*std::sqrt(a1*E1))*millibarn;
250
251return Prob;
252
253}
254
255
Note: See TracBrowser for help on using the repository browser.