G4IonisParamElm Class Reference

#include <G4IonisParamElm.hh>


Public Member Functions

 G4IonisParamElm (G4double Z)
virtual ~G4IonisParamElm ()
G4double GetZ () const
G4double GetZ3 () const
G4double GetZZ3 () const
G4double GetlogZ3 () const
G4double GetTau0 () const
G4double GetTaul () const
G4double GetAlow () const
G4double GetBlow () const
G4double GetClow () const
G4double GetMeanExcitationEnergy () const
G4double GetFermiVelocity () const
G4double GetLFactor () const
G4doubleGetShellCorrectionVector () const
 G4IonisParamElm (G4IonisParamElm &)
G4IonisParamElmoperator= (const G4IonisParamElm &)
G4int operator== (const G4IonisParamElm &) const
G4int operator!= (const G4IonisParamElm &) const
 G4IonisParamElm (__void__ &)


Detailed Description

Definition at line 51 of file G4IonisParamElm.hh.


Constructor & Destructor Documentation

G4IonisParamElm::G4IonisParamElm ( G4double  Z  ) 

Definition at line 52 of file G4IonisParamElm.cc.

References FatalException, G4Exception(), G4Pow::GetInstance(), G4NistManager::GetMeanIonisationEnergy(), G4NistManager::Instance(), G4Pow::logZ(), and G4Pow::Z13().

00053 {
00054   G4int Z = G4int(AtomNumber + 0.5);
00055   if (Z < 1) {
00056     G4Exception("G4IonisParamElm::G4IonisParamElm()",  "mat501", FatalException,
00057                 "It is not allowed to create an Element with Z<1");
00058   }
00059   G4Pow* g4pow = G4Pow::GetInstance();
00060 
00061   // some basic functions of the atomic number
00062   fZ     = Z;
00063   fZ3    = g4pow->Z13(Z);
00064   fZZ3   = fZ3*g4pow->Z13(Z+1);
00065   flogZ3 = g4pow->logZ(Z)/3.;
00066    
00067   fMeanExcitationEnergy = 
00068     G4NistManager::Instance()->GetMeanIonisationEnergy(Z);
00069 
00070   // compute parameters for ion transport
00071   // The aproximation from:
00072   // J.F.Ziegler, J.P. Biersack, U. Littmark
00073   // The Stopping and Range of Ions in Matter,
00074   // Vol.1, Pergamon Press, 1985
00075   // Fast ions or hadrons
00076 
00077   G4int iz = Z - 1;
00078   if(91 < iz) { iz = 91; }
00079 
00080   static G4double vFermi[92] = {
00081     1.0309,  0.15976, 0.59782, 1.0781,  1.0486,  1.0,     1.058,   0.93942, 0.74562, 0.3424,
00082     0.45259, 0.71074, 0.90519, 0.97411, 0.97184, 0.89852, 0.70827, 0.39816, 0.36552, 0.62712,
00083     0.81707, 0.9943,  1.1423,  1.2381,  1.1222,  0.92705, 1.0047,  1.2,     1.0661,  0.97411,
00084     0.84912, 0.95,    1.0903,  1.0429,  0.49715, 0.37755, 0.35211, 0.57801, 0.77773, 1.0207,
00085     1.029,   1.2542,  1.122,   1.1241,  1.0882,  1.2709,  1.2542,  0.90094, 0.74093, 0.86054,
00086     0.93155, 1.0047,  0.55379, 0.43289, 0.32636, 0.5131,  0.695,   0.72591, 0.71202, 0.67413,
00087     0.71418, 0.71453, 0.5911,  0.70263, 0.68049, 0.68203, 0.68121, 0.68532, 0.68715, 0.61884,
00088     0.71801, 0.83048, 1.1222,  1.2381,  1.045,   1.0733,  1.0953,  1.2381,  1.2879,  0.78654,
00089     0.66401, 0.84912, 0.88433, 0.80746, 0.43357, 0.41923, 0.43638, 0.51464, 0.73087, 0.81065,
00090     1.9578,  1.0257} ;
00091 
00092   static G4double lFactor[92] = {
00093     1.0,  1.0,  1.1,  1.06, 1.01, 1.03, 1.04, 0.99, 0.95, 0.9,
00094     0.82, 0.81, 0.83, 0.88, 1.0,  0.95, 0.97, 0.99, 0.98, 0.97,
00095     0.98, 0.97, 0.96, 0.93, 0.91, 0.9,  0.88, 0.9,  0.9,  0.9,
00096     0.9,  0.85, 0.9,  0.9,  0.91, 0.92, 0.9,  0.9,  0.9,  0.9,
00097     0.9,  0.88, 0.9,  0.88, 0.88, 0.9,  0.9,  0.88, 0.9,  0.9,
00098     0.9,  0.9,  0.96, 1.2,  0.9,  0.88, 0.88, 0.85, 0.9,  0.9,
00099     0.92, 0.95, 0.99, 1.03, 1.05, 1.07, 1.08, 1.1,  1.08, 1.08,
00100     1.08, 1.08, 1.09, 1.09, 1.1,  1.11, 1.12, 1.13, 1.14, 1.15,
00101     1.17, 1.2,  1.18, 1.17, 1.17, 1.16, 1.16, 1.16, 1.16, 1.16,
00102     1.16, 1.16} ;
00103 
00104   fVFermi  = vFermi[iz];
00105   fLFactor = lFactor[iz];
00106 
00107   // obsolete parameters for ionisation
00108   fTau0 = 0.1*fZ3*MeV/proton_mass_c2;
00109   fTaul = 2.*MeV/proton_mass_c2;
00110 
00111   // compute the Bethe-Bloch formula for energy = fTaul*particle mass
00112   G4double rate = fMeanExcitationEnergy/electron_mass_c2 ;
00113   G4double w = fTaul*(fTaul+2.) ;
00114   fBetheBlochLow = (fTaul+1.)*(fTaul+1.)*std::log(2.*w/rate)/w - 1. ;
00115   fBetheBlochLow = 2.*fZ*twopi_mc2_rcl2*fBetheBlochLow ; 
00116   
00117   fClow = std::sqrt(fTaul)*fBetheBlochLow;
00118   fAlow = 6.458040 * fClow/fTau0;
00119   G4double Taum = 0.035*fZ3*MeV/proton_mass_c2;
00120   fBlow =-3.229020*fClow/(fTau0*std::sqrt(Taum));
00121 
00122   // Shell correction parameterization
00123   fShellCorrectionVector = new G4double[3]; //[3]
00124   rate = 0.001*fMeanExcitationEnergy/eV;
00125   G4double rate2 = rate*rate;
00126     /*    
00127     fShellCorrectionVector[0] = ( 1.10289e5 + 5.14781e8*rate)*rate2 ;
00128     fShellCorrectionVector[1] = ( 7.93805e3 - 2.22565e7*rate)*rate2 ;
00129     fShellCorrectionVector[2] = (-9.92256e1 + 2.10823e5*rate)*rate2 ;
00130     */
00131   fShellCorrectionVector[0] = ( 0.422377   + 3.858019*rate)*rate2 ;
00132   fShellCorrectionVector[1] = ( 0.0304043  - 0.1667989*rate)*rate2 ;
00133   fShellCorrectionVector[2] = (-0.00038106 + 0.00157955*rate)*rate2 ;
00134 }

G4IonisParamElm::~G4IonisParamElm (  )  [virtual]

Definition at line 150 of file G4IonisParamElm.cc.

00151 {
00152   if (fShellCorrectionVector) { delete [] fShellCorrectionVector; }
00153 }

G4IonisParamElm::G4IonisParamElm ( G4IonisParamElm  ) 

Definition at line 157 of file G4IonisParamElm.cc.

00158 {
00159   fShellCorrectionVector = 0;
00160   *this = right;
00161 }

G4IonisParamElm::G4IonisParamElm ( __void__ &   ) 

Definition at line 141 of file G4IonisParamElm.cc.

00142   : fShellCorrectionVector(0)
00143 {
00144   fZ=fZ3=fZZ3=flogZ3=fTau0=fTaul=fBetheBlochLow=fAlow=fBlow=fClow
00145     =fMeanExcitationEnergy=fVFermi=fLFactor=0.0;
00146 }


Member Function Documentation

G4double G4IonisParamElm::GetAlow (  )  const [inline]

Definition at line 70 of file G4IonisParamElm.hh.

00070 {return fAlow;}

G4double G4IonisParamElm::GetBlow (  )  const [inline]

Definition at line 71 of file G4IonisParamElm.hh.

00071 {return fBlow;}

G4double G4IonisParamElm::GetClow (  )  const [inline]

Definition at line 72 of file G4IonisParamElm.hh.

00072 {return fClow;}

G4double G4IonisParamElm::GetFermiVelocity (  )  const [inline]

Definition at line 78 of file G4IonisParamElm.hh.

00078 {return fVFermi;};

G4double G4IonisParamElm::GetLFactor (  )  const [inline]

Definition at line 79 of file G4IonisParamElm.hh.

00079 {return fLFactor;};

G4double G4IonisParamElm::GetlogZ3 (  )  const [inline]

Definition at line 63 of file G4IonisParamElm.hh.

Referenced by G4PairProductionRelModel::SampleSecondaries(), G4LivermorePolarizedGammaConversionModel::SampleSecondaries(), G4LivermoreNuclearGammaConversionModel::SampleSecondaries(), G4LivermoreGammaConversionModelRC::SampleSecondaries(), G4LivermoreGammaConversionModel::SampleSecondaries(), G4eBremsstrahlungModel::SampleSecondaries(), and G4BetheHeitlerModel::SampleSecondaries().

00063 {return flogZ3;}   // std::log(Z)/3

G4double G4IonisParamElm::GetMeanExcitationEnergy (  )  const [inline]

Definition at line 75 of file G4IonisParamElm.hh.

00075 {return fMeanExcitationEnergy;}

G4double* G4IonisParamElm::GetShellCorrectionVector (  )  const [inline]

Definition at line 81 of file G4IonisParamElm.hh.

00081 {return fShellCorrectionVector;}

G4double G4IonisParamElm::GetTau0 (  )  const [inline]

Definition at line 65 of file G4IonisParamElm.hh.

00065 {return fTau0;};

G4double G4IonisParamElm::GetTaul (  )  const [inline]

Definition at line 68 of file G4IonisParamElm.hh.

00068 {return fTaul;}        // 2*MeV/proton mass

G4double G4IonisParamElm::GetZ (  )  const [inline]

Definition at line 60 of file G4IonisParamElm.hh.

00060 {return fZ;}       // atomic number Z

G4double G4IonisParamElm::GetZ3 (  )  const [inline]

Definition at line 61 of file G4IonisParamElm.hh.

Referenced by G4PairProductionRelModel::SampleSecondaries(), G4LivermorePolarizedGammaConversionModel::SampleSecondaries(), G4LivermoreNuclearGammaConversionModel::SampleSecondaries(), G4LivermoreGammaConversionModelRC::SampleSecondaries(), G4LivermoreGammaConversionModel::SampleSecondaries(), G4eBremsstrahlungModel::SampleSecondaries(), and G4BetheHeitlerModel::SampleSecondaries().

00061 {return fZ3;}      // std::pow (Z,1/3)

G4double G4IonisParamElm::GetZZ3 (  )  const [inline]

Definition at line 62 of file G4IonisParamElm.hh.

Referenced by G4eBremsstrahlungModel::SampleSecondaries().

00062 {return fZZ3;}     // std::pow (Z(Z+1),1/3)

G4int G4IonisParamElm::operator!= ( const G4IonisParamElm  )  const

Definition at line 200 of file G4IonisParamElm.cc.

00201 {
00202   return (this != (G4IonisParamElm *) &right);
00203 }

G4IonisParamElm & G4IonisParamElm::operator= ( const G4IonisParamElm  ) 

Definition at line 165 of file G4IonisParamElm.cc.

References fAlow, fBetheBlochLow, fBlow, fClow, fLFactor, flogZ3, fMeanExcitationEnergy, fShellCorrectionVector, fTau0, fTaul, fVFermi, fZ, fZ3, and fZZ3.

00166 {
00167   if (this != &right)
00168     {
00169       fZ                     = right.fZ;
00170       fZ3                    = right.fZ3;
00171       fZZ3                   = right.fZZ3;
00172       flogZ3                 = right.flogZ3;
00173       fTau0                  = right.fTau0;
00174       fTaul                  = right.fTaul;
00175       fBetheBlochLow         = right.fBetheBlochLow;
00176       fAlow                  = right.fAlow;
00177       fBlow                  = right.fBlow;
00178       fClow                  = right.fClow;
00179       fMeanExcitationEnergy  = right.fMeanExcitationEnergy;
00180       fVFermi                = right.fVFermi;
00181       fLFactor               = right.fLFactor;
00182       if (fShellCorrectionVector) { delete [] fShellCorrectionVector; } 
00183       fShellCorrectionVector = new G4double[3];            
00184       fShellCorrectionVector[0] = right.fShellCorrectionVector[0];
00185       fShellCorrectionVector[1] = right.fShellCorrectionVector[1];
00186       fShellCorrectionVector[2] = right.fShellCorrectionVector[2];      
00187     } 
00188   return *this;
00189 }

G4int G4IonisParamElm::operator== ( const G4IonisParamElm  )  const

Definition at line 193 of file G4IonisParamElm.cc.

00194 {
00195   return (this == (G4IonisParamElm *) &right);
00196 }


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