source: trunk/source/processes/hadronic/models/de_excitation/photon_evaporation/src/G4E1Probability001.cc @ 819

Last change on this file since 819 was 819, checked in by garnier, 16 years ago

import all except CVS

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