// // ******************************************************************** // * License and Disclaimer * // * * // * The Geant4 software is copyright of the Copyright Holders of * // * the Geant4 Collaboration. It is provided under the terms and * // * conditions of the Geant4 Software License, included in the file * // * LICENSE and available at http://cern.ch/geant4/license . These * // * include a list of copyright holders. * // * * // * Neither the authors of this software system, nor their employing * // * institutes,nor the agencies providing financial support for this * // * work make any representation or warranty, express or implied, * // * regarding this software system or assume any liability for its * // * use. Please see the license in the file LICENSE and URL above * // * for the full disclaimer and the limitation of liability. * // * * // * This code implementation is the result of the scientific and * // * technical work of the GEANT4 collaboration. * // * By using, copying, modifying or distributing the software (or * // * any work based on the software) you agree to acknowledge its * // * use in resulting scientific publications, and indicate your * // * acceptance of all terms of the Geant4 Software license. * // ******************************************************************** // // $Id: G4E1SingleProbability1.cc,v 1.5 2010/11/17 16:50:53 vnivanch Exp $ // GEANT4 tag $Name: geant4-09-04-ref-00 $ // // Class G4E1SingleProbability1.cc // #include "G4E1SingleProbability1.hh" #include "G4ConstantLevelDensityParameter.hh" #include "Randomize.hh" #include "G4Pow.hh" // Constructors and operators // G4E1SingleProbability1::G4E1SingleProbability1() {} G4E1SingleProbability1::~G4E1SingleProbability1() {} // Calculate the emission probability // G4double G4E1SingleProbability1::EmissionProbDensity(const G4Fragment& frag, G4double exciteE) { // Calculate the probability density here // From nuclear fragment properties and the excitation energy, calculate // the probability density for photon evaporation from U to U - exciteE // (U = nucleus excitation energy, exciteE = total evaporated photon // energy). // fragment = nuclear fragment BEFORE de-excitation G4double theProb = 0.0; G4int Afrag = frag.GetA_asInt(); G4int Zfrag = frag.GetZ_asInt(); G4double Uexcite = frag.GetExcitationEnergy(); if( (Uexcite-exciteE) < 0.0 || exciteE < 0 || Uexcite <= 0) return theProb; // Need a level density parameter. // For now, just use the constant approximation (not reliable near magic // nuclei). G4ConstantLevelDensityParameter a; G4double aLevelDensityParam = a.LevelDensityParameter(Afrag,Zfrag,Uexcite); G4double levelDensBef = std::exp(2.0*std::sqrt(aLevelDensityParam*Uexcite)); G4double levelDensAft = std::exp(2.0*std::sqrt(aLevelDensityParam*(Uexcite-exciteE))); // Now form the probability density // Define constants for the photoabsorption cross-section (the reverse // process of our de-excitation) G4double sigma0 = 2.5 * Afrag * millibarn; // millibarns G4double Egdp = (40.3 / G4Pow::GetInstance()->powZ(Afrag,0.2) )*MeV; G4double GammaR = 0.30 * Egdp; const G4double normC = 1.0 / ((pi * hbarc)*(pi * hbarc)); // CD //cout<<" PROB TESTS "<