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