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