source: trunk/source/processes/hadronic/models/de_excitation/photon_evaporation/src/G4E1Probability.cc @ 1340

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 6.9 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//
27// $Id: G4E1Probability.cc,v 1.9 2010/05/25 10:47:24 vnivanch Exp $
28// GEANT4 tag $Name: geant4-09-03-ref-09 $
29//
30//---------------------------------------------------------------------
31//
32// Geant4 class G4E1Probability
33//
34// by V. Lara (May 2003)
35//
36// Modifications:
37// 18.05.2010 V.Ivanchenko trying to speedup the most slow method
38//            by usage of G4Pow, integer A and introduction of const members
39//
40//
41
42#include "G4E1Probability.hh"
43//#include "G4ConstantLevelDensityParameter.hh"
44#include "Randomize.hh"
45#include "G4Pow.hh"
46
47// Constructors and operators
48//
49
50G4E1Probability::G4E1Probability():G4VEmissionProbability()
51{
52  G4double x = CLHEP::pi*CLHEP::hbarc;
53  normC = 1.0 / (x*x);
54  theLevelDensityParameter = 0.125/MeV;
55  fG4pow = G4Pow::GetInstance(); 
56}
57
58G4E1Probability::~G4E1Probability()
59{}
60
61// Calculate the emission probability
62//
63
64G4double G4E1Probability::EmissionProbDensity(const G4Fragment& frag, 
65                                              const G4double gammaE)
66{
67
68  // Calculate the probability density here
69
70  // From nuclear fragment properties and the excitation energy, calculate
71  // the probability density for photon evaporation from U to U - gammaE
72  // (U = nucleus excitation energy, gammaE = total evaporated photon
73  // energy). Fragment = nuclear fragment BEFORE de-excitation
74
75  G4double theProb = 0.0;
76
77  G4int Afrag = frag.GetA_asInt();
78  G4double Uexcite = frag.GetExcitationEnergy();
79
80  if( (Uexcite-gammaE) < 0.0 || gammaE < 0) { return theProb; }
81
82  // Need a level density parameter.
83  // For now, just use the constant approximation (not reliable near magic
84  // nuclei).
85
86  G4double aLevelDensityParam = Afrag*theLevelDensityParameter;
87
88  //  G4double levelDensBef = std::exp(2*std::sqrt(aLevelDensityParam*Uexcite));
89  //  G4double levelDensAft = std::exp(2*std::sqrt(aLevelDensityParam*(Uexcite-gammaE)));
90  // VI reduce number of calls to exp
91  G4double levelDens = 1.0;
92  if( aLevelDensityParam > 0.0 ) {
93    levelDens = std::exp(2*(std::sqrt(aLevelDensityParam*(Uexcite-gammaE)) 
94      - std::sqrt(aLevelDensityParam*Uexcite))); 
95  }
96  // Now form the probability density
97
98  // Define constants for the photoabsorption cross-section (the reverse
99  // process of our de-excitation)
100
101  G4double sigma0 = 2.5 * Afrag * millibarn;  // millibarns
102
103  G4double Egdp   = (40.3 / fG4pow->powZ(Afrag,0.2) )*MeV;
104  G4double GammaR = 0.30 * Egdp;
105 
106  // CD
107  //cout<<"  PROB TESTS "<<G4endl;
108  //cout<<" hbarc = "<<hbarc<<G4endl;
109  //cout<<" pi = "<<pi<<G4endl;
110  //cout<<" Uexcite, gammaE = "<<Uexcite<<"  "<<gammaE<<G4endl;
111  //cout<<" Uexcite, gammaE = "<<Uexcite*MeV<<"  "<<gammaE*MeV<<G4endl;
112  //cout<<" lev density param = "<<aLevelDensityParam<<G4endl;
113  //cout<<" level densities = "<<levelDensBef<<"  "<<levelDensAft<<G4endl;
114  //cout<<" sigma0 = "<<sigma0<<G4endl;
115  //cout<<" Egdp, GammaR = "<<Egdp<<"  "<<GammaR<<G4endl;
116  //cout<<" normC = "<<normC<<G4endl;
117
118  // VI implementation 18.05.2010
119  G4double gammaE2 = gammaE*gammaE;
120  G4double gammaR2 = gammaE2*GammaR*GammaR;
121  G4double egdp2   = gammaE2 - Egdp*Egdp;
122  G4double sigmaAbs = sigma0*gammaR2/(egdp2*egdp2 + gammaR2); 
123  theProb = normC * sigmaAbs * gammaE2 * levelDens;
124
125  // old implementation
126  //  G4double numerator = sigma0 * gammaE*gammaE * GammaR*GammaR;
127  // G4double denominator = (gammaE*gammaE - Egdp*Egdp)*
128  //         (gammaE*gammaE - Egdp*Egdp) + GammaR*GammaR*gammaE*gammaE;
129
130  //G4double sigmaAbs = numerator/denominator;
131  //theProb = normC * sigmaAbs * gammaE2 * levelDensAft/levelDensBef;
132
133  // CD
134  //cout<<" sigmaAbs = "<<sigmaAbs<<G4endl;
135  //cout<<" Probability = "<<theProb<<G4endl;
136
137  return theProb;
138
139}
140
141G4double G4E1Probability::EmissionProbability(const G4Fragment& frag, 
142                                                 const G4double gammaE)
143{
144
145  // From nuclear fragment properties and the excitation energy, calculate
146  // the probability for photon evaporation down to last ground level.
147  // fragment = nuclear fragment BEFORE de-excitation
148
149  G4double theProb = 0.0;
150
151  G4double Uafter = 0.0;
152  const G4double Uexcite = frag.GetExcitationEnergy();
153
154  G4double normC = 3.0;
155
156  const G4double upperLim = Uexcite;
157  const G4double lowerLim = Uafter;
158  const G4int numIters = 100;
159
160  // Fall-back is a uniform random number
161
162  //G4double uniformNum = G4UniformRand();
163  //theProb = uniformNum;
164
165  // Need to integrate EmissionProbDensity from lowerLim to upperLim
166  // and multiply by normC
167
168  G4double integ = normC *
169           EmissionIntegration(frag,gammaE,lowerLim,upperLim,numIters);
170  if(integ > 0.0) theProb = integ/(upperLim-lowerLim);
171
172  return theProb;
173
174}
175
176G4double G4E1Probability::EmissionIntegration(const G4Fragment& frag, 
177                             const G4double ,
178                             const G4double lowLim, const G4double upLim,
179                             const G4int numIters)
180
181{
182
183  // Simple Gaussian quadrature integration
184
185  G4double x;
186  G4double root3 = 1.0/std::sqrt(3.0);
187
188  G4double Step = (upLim-lowLim)/(2.0*numIters);
189  G4double Delta = Step*root3;
190
191  G4double mean = 0.0;
192
193  G4double theInt = 0.0;
194
195  for(G4int i = 0; i < numIters; i++) {
196
197    x = (2*i + 1)*Step;
198    G4double E1ProbDensityA = EmissionProbDensity(frag,x+Delta);
199    G4double E1ProbDensityB = EmissionProbDensity(frag,x-Delta);
200
201    mean += E1ProbDensityA + E1ProbDensityB;
202
203  }
204
205  if(mean*Step > 0.0) theInt = mean*Step;
206
207  return theInt;
208
209}
210
211
Note: See TracBrowser for help on using the repository browser.