G4E1Probability Class Reference

#include <G4E1Probability.hh>

Inheritance diagram for G4E1Probability:

G4VEmissionProbability

Public Member Functions

 G4E1Probability ()
virtual ~G4E1Probability ()
G4double EmissionProbability (const G4Fragment &frag, G4double excite)
G4double EmissionProbDensity (const G4Fragment &frag, G4double ePhoton)

Detailed Description

Definition at line 49 of file G4E1Probability.hh.


Constructor & Destructor Documentation

G4E1Probability::G4E1Probability (  ) 

Definition at line 50 of file G4E1Probability.cc.

References G4Pow::GetInstance(), and G4INCL::Math::pi.

00050                                 :G4VEmissionProbability()
00051 {
00052   G4double x = CLHEP::pi*CLHEP::hbarc;
00053   normC = 1.0 / (x*x);
00054   theLevelDensityParameter = 0.125/MeV;
00055   fG4pow = G4Pow::GetInstance(); 
00056 }

G4E1Probability::~G4E1Probability (  )  [virtual]

Definition at line 58 of file G4E1Probability.cc.

00059 {}


Member Function Documentation

G4double G4E1Probability::EmissionProbability ( const G4Fragment frag,
G4double  excite 
) [virtual]

Implements G4VEmissionProbability.

Definition at line 139 of file G4E1Probability.cc.

00141 {
00142   // From nuclear fragment properties and the excitation energy, calculate
00143   // the probability for photon evaporation down to last ground level.
00144   // fragment = nuclear fragment BEFORE de-excitation
00145 
00146   G4double upperLim = gammaE;
00147   G4double lowerLim = 0.0; 
00148 
00149   //G4cout << "G4E1Probability::EmissionProbability:  Emin= " << lowerLim
00150   //     << " Emax= " << upperLim << G4endl;
00151   if( upperLim - lowerLim <= CLHEP::keV ) { return 0.0; } 
00152 
00153   // Need to integrate EmissionProbDensity from lowerLim to upperLim 
00154   // and multiply by factor 3 (?!)
00155 
00156   G4double integ = 3.0 * EmissionIntegration(frag,lowerLim,upperLim);
00157 
00158   return integ;
00159 
00160 }

G4double G4E1Probability::EmissionProbDensity ( const G4Fragment frag,
G4double  ePhoton 
)

Definition at line 64 of file G4E1Probability.cc.

References G4Fragment::GetA_asInt(), G4Fragment::GetExcitationEnergy(), and G4Pow::powZ().

00066 {
00067 
00068   // Calculate the probability density here
00069 
00070   // From nuclear fragment properties and the excitation energy, calculate
00071   // the probability density for photon evaporation from U to U - gammaE
00072   // (U = nucleus excitation energy, gammaE = total evaporated photon
00073   // energy). Fragment = nuclear fragment BEFORE de-excitation
00074 
00075   G4double theProb = 0.0;
00076 
00077   G4int Afrag = frag.GetA_asInt();
00078   G4double Uexcite = frag.GetExcitationEnergy();
00079   G4double U = std::max(0.0,Uexcite-gammaE);
00080 
00081   if(gammaE < 0.0) { return theProb; }
00082 
00083   // Need a level density parameter.
00084   // For now, just use the constant approximation (not reliable near magic
00085   // nuclei).
00086 
00087   G4double aLevelDensityParam = Afrag*theLevelDensityParameter;
00088 
00089   //  G4double levelDensBef = std::exp(2*std::sqrt(aLevelDensityParam*Uexcite));
00090   //  G4double levelDensAft = std::exp(2*std::sqrt(aLevelDensityParam*(Uexcite-gammaE)));
00091   // VI reduce number of calls to exp 
00092   G4double levelDens = 
00093     std::exp(2*(std::sqrt(aLevelDensityParam*U)-std::sqrt(aLevelDensityParam*Uexcite)));
00094   // Now form the probability density
00095 
00096   // Define constants for the photoabsorption cross-section (the reverse
00097   // process of our de-excitation)
00098 
00099   G4double sigma0 = 2.5 * Afrag * millibarn;  // millibarns
00100 
00101   G4double Egdp   = (40.3 / fG4pow->powZ(Afrag,0.2) )*MeV;
00102   G4double GammaR = 0.30 * Egdp;
00103  
00104   // CD
00105   //cout<<"  PROB TESTS "<<G4endl;
00106   //cout<<" hbarc = "<<hbarc<<G4endl;
00107   //cout<<" pi = "<<pi<<G4endl;
00108   //cout<<" Uexcite, gammaE = "<<Uexcite<<"  "<<gammaE<<G4endl;
00109   //cout<<" Uexcite, gammaE = "<<Uexcite*MeV<<"  "<<gammaE*MeV<<G4endl;
00110   //cout<<" lev density param = "<<aLevelDensityParam<<G4endl;
00111   //cout<<" level densities = "<<levelDensBef<<"  "<<levelDensAft<<G4endl;
00112   //cout<<" sigma0 = "<<sigma0<<G4endl;
00113   //cout<<" Egdp, GammaR = "<<Egdp<<"  "<<GammaR<<G4endl;
00114   //cout<<" normC = "<<normC<<G4endl;
00115 
00116   // VI implementation 18.05.2010
00117   G4double gammaE2 = gammaE*gammaE;
00118   G4double gammaR2 = gammaE2*GammaR*GammaR;
00119   G4double egdp2   = gammaE2 - Egdp*Egdp;
00120   G4double sigmaAbs = sigma0*gammaR2/(egdp2*egdp2 + gammaR2); 
00121   theProb = normC * sigmaAbs * gammaE2 * levelDens;
00122 
00123   // old implementation
00124   //  G4double numerator = sigma0 * gammaE*gammaE * GammaR*GammaR;
00125   // G4double denominator = (gammaE*gammaE - Egdp*Egdp)*
00126   //         (gammaE*gammaE - Egdp*Egdp) + GammaR*GammaR*gammaE*gammaE;
00127 
00128   //G4double sigmaAbs = numerator/denominator; 
00129   //theProb = normC * sigmaAbs * gammaE2 * levelDensAft/levelDensBef;
00130 
00131   // CD
00132   //cout<<" sigmaAbs = "<<sigmaAbs<<G4endl;
00133   //cout<<" Probability = "<<theProb<<G4endl;
00134 
00135   return theProb;
00136 
00137 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:51:49 2013 for Geant4 by  doxygen 1.4.7