G4BraggModel Class Reference

#include <G4BraggModel.hh>

Inheritance diagram for G4BraggModel:

G4VEmModel G4BraggIonGasModel

Public Member Functions

 G4BraggModel (const G4ParticleDefinition *p=0, const G4String &nam="Bragg")
virtual ~G4BraggModel ()
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
virtual G4double GetParticleCharge (const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy)

Protected Member Functions

virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy)
G4double GetChargeSquareRatio () const
void SetChargeSquareRatio (G4double val)

Detailed Description

Definition at line 71 of file G4BraggModel.hh.


Constructor & Destructor Documentation

G4BraggModel::G4BraggModel ( const G4ParticleDefinition p = 0,
const G4String nam = "Bragg" 
)

Definition at line 79 of file G4BraggModel.cc.

References G4Electron::Electron(), G4LossTableManager::EmCorrections(), G4LossTableManager::Instance(), and G4VEmModel::SetHighEnergyLimit().

00080   : G4VEmModel(nam),
00081     particle(0),
00082     currentMaterial(0),
00083     protonMassAMU(1.007276),
00084     iMolecula(-1),
00085     iPSTAR(-1),
00086     isIon(false),
00087     isInitialised(false)
00088 {
00089   fParticleChange = 0;
00090   SetHighEnergyLimit(2.0*MeV);
00091 
00092   lowestKinEnergy  = 1.0*keV;
00093   theZieglerFactor = eV*cm2*1.0e-15;
00094   theElectron = G4Electron::Electron();
00095   expStopPower125 = 0.0;
00096 
00097   corr = G4LossTableManager::Instance()->EmCorrections();
00098   if(p) { SetParticle(p); }
00099   else  { SetParticle(theElectron); }
00100 }

G4BraggModel::~G4BraggModel (  )  [virtual]

Definition at line 104 of file G4BraggModel.cc.

00105 {}


Member Function Documentation

G4double G4BraggModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  Z,
G4double  A,
G4double  cutEnergy,
G4double  maxEnergy 
) [virtual]

Reimplemented from G4VEmModel.

Definition at line 180 of file G4BraggModel.cc.

References ComputeCrossSectionPerElectron().

00186 {
00187   G4double cross = Z*ComputeCrossSectionPerElectron
00188                                          (p,kineticEnergy,cutEnergy,maxEnergy);
00189   return cross;
00190 }

G4double G4BraggModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
) [virtual]

Definition at line 154 of file G4BraggModel.cc.

References MaxSecondaryEnergy().

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

00159 {
00160   G4double cross     = 0.0;
00161   G4double tmax      = MaxSecondaryEnergy(p, kineticEnergy);
00162   G4double maxEnergy = std::min(tmax,maxKinEnergy);
00163   if(cutEnergy < maxEnergy) {
00164 
00165     G4double energy  = kineticEnergy + mass;
00166     G4double energy2 = energy*energy;
00167     G4double beta2   = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
00168     cross = 1.0/cutEnergy - 1.0/maxEnergy - beta2*log(maxEnergy/cutEnergy)/tmax;
00169 
00170     cross *= twopi_mc2_rcl2*chargeSquare/beta2;
00171   }
00172  //   G4cout << "BR: e= " << kineticEnergy << " tmin= " << cutEnergy 
00173  //          << " tmax= " << tmax << " cross= " << cross << G4endl;
00174  
00175   return cross;
00176 }

G4double G4BraggModel::ComputeDEDXPerVolume ( const G4Material ,
const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  cutEnergy 
) [virtual]

Reimplemented from G4VEmModel.

Definition at line 209 of file G4BraggModel.cc.

References G4InuclParticleNames::gam, G4Material::GetElectronDensity(), and MaxSecondaryEnergy().

00213 {
00214   G4double tmax  = MaxSecondaryEnergy(p, kineticEnergy);
00215   G4double tkin  = kineticEnergy/massRate;
00216   G4double dedx  = 0.0;
00217 
00218   if(tkin < lowestKinEnergy) {
00219     dedx = DEDX(material, lowestKinEnergy)*sqrt(tkin/lowestKinEnergy);
00220   } else {
00221     dedx = DEDX(material, tkin); 
00222   }
00223 
00224   if (cutEnergy < tmax) {
00225 
00226     G4double tau   = kineticEnergy/mass;
00227     G4double gam   = tau + 1.0;
00228     G4double bg2   = tau * (tau+2.0);
00229     G4double beta2 = bg2/(gam*gam);
00230     G4double x     = cutEnergy/tmax;
00231 
00232     dedx += (log(x) + (1.0 - x)*beta2) * twopi_mc2_rcl2
00233           * (material->GetElectronDensity())/beta2;
00234   }
00235 
00236   // now compute the total ionization loss
00237 
00238   if (dedx < 0.0) { dedx = 0.0; }
00239 
00240   dedx *= chargeSquare;
00241 
00242   //G4cout << "E(MeV)= " << tkin/MeV << " dedx= " << dedx 
00243   //     << "  " << material->GetName() << G4endl;
00244 
00245   return dedx;
00246 }

G4double G4BraggModel::CrossSectionPerVolume ( const G4Material ,
const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
) [virtual]

Reimplemented from G4VEmModel.

Definition at line 194 of file G4BraggModel.cc.

References ComputeCrossSectionPerElectron(), and G4Material::GetElectronDensity().

00200 {
00201   G4double eDensity = material->GetElectronDensity();
00202   G4double cross = eDensity*ComputeCrossSectionPerElectron
00203                                          (p,kineticEnergy,cutEnergy,maxEnergy);
00204   return cross;
00205 }

G4double G4BraggModel::GetChargeSquareRatio (  )  const [inline, protected]

Definition at line 194 of file G4BraggModel.hh.

00195 {
00196   return chargeSquare;
00197 }

G4double G4BraggModel::GetChargeSquareRatio ( const G4ParticleDefinition ,
const G4Material ,
G4double  kineticEnergy 
) [virtual]

Reimplemented from G4VEmModel.

Definition at line 132 of file G4BraggModel.cc.

References G4EmCorrections::EffectiveChargeCorrection(), G4EmCorrections::EffectiveChargeSquareRatio(), G4VEmModel::GetModelOfFluctuations(), and G4VEmFluctuationModel::SetParticleAndCharge().

00135 {
00136   // this method is called only for ions
00137   G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy);
00138   GetModelOfFluctuations()->SetParticleAndCharge(p, q2);
00139   return q2*corr->EffectiveChargeCorrection(p,mat,kineticEnergy);
00140 }

G4double G4BraggModel::GetParticleCharge ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
) [virtual]

Reimplemented from G4VEmModel.

Reimplemented in G4BraggIonGasModel.

Definition at line 144 of file G4BraggModel.cc.

References G4EmCorrections::GetParticleCharge().

00147 {
00148   // this method is called only for ions, so no check if it is an ion 
00149   return corr->GetParticleCharge(p,mat,kineticEnergy);
00150 }

void G4BraggModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
) [virtual]

Implements G4VEmModel.

Definition at line 109 of file G4BraggModel.cc.

References G4VEmModel::GetParticleChangeForLoss(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetParticleType(), and G4VEmModel::SetDeexcitationFlag().

00111 {
00112   if(p != particle) { SetParticle(p); }
00113 
00114   // always false before the run
00115   SetDeexcitationFlag(false);
00116 
00117   if(!isInitialised) {
00118     isInitialised = true;
00119 
00120     G4String pname = particle->GetParticleName();
00121     if(particle->GetParticleType() == "nucleus" && 
00122        pname != "deuteron" && pname != "triton" &&
00123        pname != "alpha+"   && pname != "helium" &&
00124        pname != "hydrogen") { isIon = true; }
00125 
00126     fParticleChange = GetParticleChangeForLoss();
00127   }
00128 }

G4double G4BraggModel::MaxSecondaryEnergy ( const G4ParticleDefinition ,
G4double  kinEnergy 
) [protected, virtual]

Reimplemented from G4VEmModel.

Definition at line 315 of file G4BraggModel.cc.

Referenced by ComputeCrossSectionPerElectron(), and ComputeDEDXPerVolume().

00317 {
00318   if(pd != particle) { SetParticle(pd); }
00319   G4double tau  = kinEnergy/mass;
00320   G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
00321                   (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
00322   return tmax;
00323 }

virtual void G4BraggModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle ,
G4double  tmin,
G4double  maxEnergy 
) [virtual]

Implements G4VEmModel.

void G4BraggModel::SetChargeSquareRatio ( G4double  val  )  [inline, protected]

Definition at line 199 of file G4BraggModel.hh.

Referenced by G4BraggIonGasModel::ChargeSquareRatio().

00200 {
00201   chargeSquare = val;
00202 }


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