00001 // 00002 // ******************************************************************** 00003 // * License and Disclaimer * 00004 // * * 00005 // * The Geant4 software is copyright of the Copyright Holders of * 00006 // * the Geant4 Collaboration. It is provided under the terms and * 00007 // * conditions of the Geant4 Software License, included in the file * 00008 // * LICENSE and available at http://cern.ch/geant4/license . These * 00009 // * include a list of copyright holders. * 00010 // * * 00011 // * Neither the authors of this software system, nor their employing * 00012 // * institutes,nor the agencies providing financial support for this * 00013 // * work make any representation or warranty, express or implied, * 00014 // * regarding this software system or assume any liability for its * 00015 // * use. Please see the license in the file LICENSE and URL above * 00016 // * for the full disclaimer and the limitation of liability. * 00017 // * * 00018 // * This code implementation is the result of the scientific and * 00019 // * technical work of the GEANT4 collaboration. * 00020 // * By using, copying, modifying or distributing the software (or * 00021 // * any work based on the software) you agree to acknowledge its * 00022 // * use in resulting scientific publications, and indicate your * 00023 // * acceptance of all terms of the Geant4 Software license. * 00024 // ******************************************************************** 00025 // 00026 // 00027 00028 #ifndef G4NuclearFermiDensity_h 00029 #define G4NuclearFermiDensity_h 1 00030 00031 #include "globals.hh" 00032 #include "G4ThreeVector.hh" 00033 #include "G4VNuclearDensity.hh" 00034 00035 #include <CLHEP/Units/PhysicalConstants.h> // pi, fermi,.. 00036 #include <cmath> // pow 00037 00038 class G4NuclearFermiDensity : public G4VNuclearDensity 00039 { 00040 public: 00041 G4NuclearFermiDensity(G4int anA, G4int aZ); 00042 ~G4NuclearFermiDensity(); 00043 00044 G4double GetRelativeDensity(const G4ThreeVector & aPosition) const 00045 { 00046 return 1./(1.+std::exp((aPosition.mag()-theR)/a)); 00047 } 00048 00049 G4double GetRadius(const G4double maxRelativeDenisty) const 00050 { 00051 return (maxRelativeDenisty>0 && maxRelativeDenisty <= 1 ) ? 00052 (theR + a*std::log((1-maxRelativeDenisty+std::exp(-1*theR/a))/maxRelativeDenisty)) : DBL_MAX; 00053 } 00054 00055 G4double GetDeriv(const G4ThreeVector & aPosition) const 00056 { 00057 G4double currentR=aPosition.mag(); 00058 if (currentR > 40*theR ) {return 0;} 00059 else return -std::exp((currentR-theR)/a) * sqr(GetDensity(aPosition)) / (a*Getrho0()); 00060 } 00061 00062 private: 00063 G4int theA; 00064 G4double theR; // Nuclear Radius 00065 const G4double a; // Determines the nuclear surface thickness 00066 }; 00067 00068 #endif 00069