G4eBremParametrizedModel Class Reference

#include <G4eBremParametrizedModel.hh>

Inheritance diagram for G4eBremParametrizedModel:

G4VEmModel

Public Member Functions

 G4eBremParametrizedModel (const G4ParticleDefinition *p=0, const G4String &nam="eBremParam")
virtual ~G4eBremParametrizedModel ()
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX)
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy)
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double)

Protected Attributes

G4NistManagernist
const G4ParticleDefinitionparticle
G4ParticleDefinitiontheGamma
G4ParticleChangeForLossfParticleChange
G4double minThreshold
G4double particleMass
G4double kinEnergy
G4double totalEnergy
G4double currentZ
G4double z13
G4double z23
G4double lnZ
G4double densityFactor
G4double densityCorr
G4double Fel
G4double Finel
G4double facFel
G4double facFinel
G4double fMax
G4double fCoulomb
G4bool isElectron

Static Protected Attributes

static const G4double xgi [8]
static const G4double wgi [8]

Detailed Description

Definition at line 61 of file G4eBremParametrizedModel.hh.


Constructor & Destructor Documentation

G4eBremParametrizedModel::G4eBremParametrizedModel ( const G4ParticleDefinition p = 0,
const G4String nam = "eBremParam" 
)

Definition at line 74 of file G4eBremParametrizedModel.cc.

References currentZ, densityCorr, densityFactor, fCoulomb, Fel, Finel, fMax, G4Gamma::Gamma(), G4NistManager::Instance(), kinEnergy, lnZ, minThreshold, nist, particleMass, G4VEmModel::SetAngularDistribution(), G4VEmModel::SetLowEnergyLimit(), theGamma, totalEnergy, z13, and z23.

00076   : G4VEmModel(nam),
00077     particle(0),
00078     isElectron(true),
00079     fMigdalConstant(classic_electr_radius*electron_Compton_length*electron_Compton_length*4.0*pi),
00080     bremFactor(fine_structure_const*classic_electr_radius*classic_electr_radius*16./3.),
00081     isInitialised(false)
00082 {
00083   theGamma = G4Gamma::Gamma();
00084 
00085   minThreshold = 0.1*keV;
00086   lowKinEnergy = 10.*MeV;
00087   SetLowEnergyLimit(lowKinEnergy);  
00088 
00089   nist = G4NistManager::Instance();  
00090 
00091   SetAngularDistribution(new G4ModifiedTsai());
00092 
00093   particleMass = kinEnergy = totalEnergy = currentZ = z13 = z23 = lnZ = Fel = Finel 
00094     = densityFactor = densityCorr = fMax = fCoulomb = 0.;
00095 
00096   InitialiseConstants();
00097   if(p) { SetParticle(p); }
00098 }

G4eBremParametrizedModel::~G4eBremParametrizedModel (  )  [virtual]

Definition at line 110 of file G4eBremParametrizedModel.cc.

00111 {
00112 }


Member Function Documentation

G4double G4eBremParametrizedModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  tkin,
G4double  Z,
G4double  ,
G4double  cutEnergy,
G4double  maxEnergy = DBL_MAX 
) [virtual]

Reimplemented from G4VEmModel.

Definition at line 234 of file G4eBremParametrizedModel.cc.

References kinEnergy, and particle.

00240 {
00241   if(!particle) { SetParticle(p); }
00242   if(kineticEnergy < lowKinEnergy) { return 0.0; }
00243   G4double cut  = std::min(cutEnergy, kineticEnergy);
00244   G4double tmax = std::min(maxEnergy, kineticEnergy);
00245 
00246   if(cut >= tmax) { return 0.0; }
00247 
00248   SetCurrentElement(Z);
00249 
00250   G4double cross = ComputeXSectionPerAtom(cut);
00251 
00252   // allow partial integration
00253   if(tmax < kinEnergy) { cross -= ComputeXSectionPerAtom(tmax); }
00254   
00255   cross *= Z*Z*bremFactor;
00256 
00257   return cross;
00258 }

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

Reimplemented from G4VEmModel.

Definition at line 167 of file G4eBremParametrizedModel.cc.

References currentZ, G4Material::GetAtomicNumDensityVector(), G4Material::GetElementVector(), G4Material::GetNumberOfElements(), particle, G4VEmModel::SetCurrentElement(), and SetupForMaterial().

00172 {
00173   if(!particle) { SetParticle(p); }
00174   if(kineticEnergy < lowKinEnergy) { return 0.0; }
00175   G4double cut = std::min(cutEnergy, kineticEnergy);
00176   if(cut == 0.0) { return 0.0; }
00177 
00178   SetupForMaterial(particle, material,kineticEnergy);
00179 
00180   const G4ElementVector* theElementVector = material->GetElementVector();
00181   const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector();
00182 
00183   G4double dedx = 0.0;
00184 
00185   //  loop for elements in the material
00186   for (size_t i=0; i<material->GetNumberOfElements(); i++) {
00187 
00188     G4VEmModel::SetCurrentElement((*theElementVector)[i]);
00189     SetCurrentElement((*theElementVector)[i]->GetZ());
00190 
00191     dedx += theAtomicNumDensityVector[i]*currentZ*currentZ*ComputeBremLoss(cut);
00192   }
00193   dedx *= bremFactor;
00194 
00195 
00196   return dedx;
00197 }

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

Implements G4VEmModel.

Definition at line 149 of file G4eBremParametrizedModel.cc.

References currentZ, fParticleChange, G4VEmModel::GetParticleChangeForLoss(), G4VEmModel::InitialiseElementSelectors(), and G4VEmModel::LowEnergyLimit().

00151 {
00152   if(p) { SetParticle(p); }
00153 
00154   lowKinEnergy  = LowEnergyLimit();
00155 
00156   currentZ = 0.;
00157 
00158   InitialiseElementSelectors(p, cuts);
00159 
00160   if(isInitialised) { return; }
00161   fParticleChange = GetParticleChangeForLoss();
00162   isInitialised = true;
00163 }

G4double G4eBremParametrizedModel::MinEnergyCut ( const G4ParticleDefinition ,
const G4MaterialCutsCouple  
) [virtual]

Definition at line 126 of file G4eBremParametrizedModel.cc.

References minThreshold.

00128 {
00129   return minThreshold;
00130 }

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

Implements G4VEmModel.

Definition at line 466 of file G4eBremParametrizedModel.cc.

References currentZ, densityCorr, densityFactor, fMax, fParticleChange, fStopAndKill, G4cout, G4endl, G4lrint(), G4UniformRand, G4VEmModel::GetAngularDistribution(), G4DynamicParticle::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), G4DynamicParticle::GetMomentumDirection(), kinEnergy, particle, particleMass, G4VParticleChange::ProposeTrackStatus(), G4VEmAngularDistribution::SampleDirection(), G4VEmModel::SecondaryThreshold(), G4VEmModel::SelectRandomAtom(), G4ParticleChangeForLoss::SetProposedKineticEnergy(), G4ParticleChangeForLoss::SetProposedMomentumDirection(), SetupForMaterial(), theGamma, and totalEnergy.

00472 {
00473   G4double kineticEnergy = dp->GetKineticEnergy();
00474   if(kineticEnergy < lowKinEnergy) { return; }
00475   G4double cut  = std::min(cutEnergy, kineticEnergy);
00476   G4double emax = std::min(maxEnergy, kineticEnergy);
00477   if(cut >= emax) { return; }
00478 
00479   SetupForMaterial(particle, couple->GetMaterial(),kineticEnergy);
00480 
00481   const G4Element* elm = 
00482     SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
00483   SetCurrentElement(elm->GetZ());
00484 
00485   kinEnergy   = kineticEnergy;
00486   totalEnergy = kineticEnergy + particleMass;
00487   densityCorr = densityFactor*totalEnergy*totalEnergy;
00488  
00489   G4double xmin = log(cut*cut + densityCorr);
00490   G4double xmax = log(emax*emax  + densityCorr);
00491   G4double gammaEnergy, f, x; 
00492 
00493   do {
00494     x = exp(xmin + G4UniformRand()*(xmax - xmin)) - densityCorr;
00495     if(x < 0.0) x = 0.0;
00496     gammaEnergy = sqrt(x);
00497     f = ComputeDXSectionPerAtom(gammaEnergy);
00498 
00499     if ( f > fMax ) {
00500       G4cout << "### G4eBremParametrizedModel Warning: Majoranta exceeded! "
00501              << f << " > " << fMax
00502              << " Egamma(MeV)= " << gammaEnergy
00503              << " E(mEV)= " << kineticEnergy
00504              << G4endl;
00505     }
00506 
00507   } while (f < fMax*G4UniformRand());
00508 
00509   //
00510   // angles of the emitted gamma. ( Z - axis along the parent particle)
00511   // use general interface
00512   //
00513   G4ThreeVector gammaDirection = 
00514     GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy,
00515                                               G4lrint(currentZ), 
00516                                               couple->GetMaterial());
00517 
00518   // create G4DynamicParticle object for the Gamma
00519   G4DynamicParticle* gamma = new G4DynamicParticle(theGamma,gammaDirection,
00520                                                    gammaEnergy);
00521   vdp->push_back(gamma);
00522   
00523   G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
00524   G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
00525                              - gammaEnergy*gammaDirection).unit();
00526 
00527   // energy of primary
00528   G4double finalE = kineticEnergy - gammaEnergy;
00529 
00530   // stop tracking and create new secondary instead of primary
00531   if(gammaEnergy > SecondaryThreshold()) {
00532     fParticleChange->ProposeTrackStatus(fStopAndKill);
00533     fParticleChange->SetProposedKineticEnergy(0.0);
00534     G4DynamicParticle* el = 
00535       new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
00536                             direction, finalE);
00537     vdp->push_back(el);
00538 
00539     // continue tracking
00540   } else {
00541     fParticleChange->SetProposedMomentumDirection(direction);
00542     fParticleChange->SetProposedKineticEnergy(finalE);
00543   }
00544 }

void G4eBremParametrizedModel::SetupForMaterial ( const G4ParticleDefinition ,
const G4Material ,
G4double   
) [virtual]

Reimplemented from G4VEmModel.

Definition at line 134 of file G4eBremParametrizedModel.cc.

References densityCorr, densityFactor, G4Material::GetElectronDensity(), kinEnergy, particleMass, and totalEnergy.

Referenced by ComputeDEDXPerVolume(), and SampleSecondaries().

00137 {
00138   densityFactor = mat->GetElectronDensity()*fMigdalConstant;
00139 
00140   // calculate threshold for density effect
00141   kinEnergy   = kineticEnergy;
00142   totalEnergy = kineticEnergy + particleMass;
00143   densityCorr = densityFactor*totalEnergy*totalEnergy;    
00144 }


Field Documentation

G4double G4eBremParametrizedModel::currentZ [protected]

Definition at line 131 of file G4eBremParametrizedModel.hh.

Referenced by ComputeDEDXPerVolume(), G4eBremParametrizedModel(), Initialise(), and SampleSecondaries().

G4double G4eBremParametrizedModel::densityCorr [protected]

Definition at line 134 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), SampleSecondaries(), and SetupForMaterial().

G4double G4eBremParametrizedModel::densityFactor [protected]

Definition at line 133 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), SampleSecondaries(), and SetupForMaterial().

G4double G4eBremParametrizedModel::facFel [protected]

Definition at line 136 of file G4eBremParametrizedModel.hh.

G4double G4eBremParametrizedModel::facFinel [protected]

Definition at line 136 of file G4eBremParametrizedModel.hh.

G4double G4eBremParametrizedModel::fCoulomb [protected]

Definition at line 137 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

G4double G4eBremParametrizedModel::Fel [protected]

Definition at line 135 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

G4double G4eBremParametrizedModel::Finel [protected]

Definition at line 135 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

G4double G4eBremParametrizedModel::fMax [protected]

Definition at line 137 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), and SampleSecondaries().

G4ParticleChangeForLoss* G4eBremParametrizedModel::fParticleChange [protected]

Definition at line 121 of file G4eBremParametrizedModel.hh.

Referenced by Initialise(), and SampleSecondaries().

G4bool G4eBremParametrizedModel::isElectron [protected]

Definition at line 139 of file G4eBremParametrizedModel.hh.

G4double G4eBremParametrizedModel::kinEnergy [protected]

Definition at line 129 of file G4eBremParametrizedModel.hh.

Referenced by ComputeCrossSectionPerAtom(), G4eBremParametrizedModel(), SampleSecondaries(), and SetupForMaterial().

G4double G4eBremParametrizedModel::lnZ [protected]

Definition at line 132 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

G4double G4eBremParametrizedModel::minThreshold [protected]

Definition at line 125 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), and MinEnergyCut().

G4NistManager* G4eBremParametrizedModel::nist [protected]

Definition at line 118 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

const G4ParticleDefinition* G4eBremParametrizedModel::particle [protected]

Definition at line 119 of file G4eBremParametrizedModel.hh.

Referenced by ComputeCrossSectionPerAtom(), ComputeDEDXPerVolume(), and SampleSecondaries().

G4double G4eBremParametrizedModel::particleMass [protected]

Definition at line 128 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), SampleSecondaries(), and SetupForMaterial().

G4ParticleDefinition* G4eBremParametrizedModel::theGamma [protected]

Definition at line 120 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), and SampleSecondaries().

G4double G4eBremParametrizedModel::totalEnergy [protected]

Definition at line 130 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), SampleSecondaries(), and SetupForMaterial().

const G4double G4eBremParametrizedModel::wgi [static, protected]

Initial value:

{ 0.0506, 0.1112, 0.1569, 0.1813,
                                            0.1813, 0.1569, 0.1112, 0.0506 }

Definition at line 123 of file G4eBremParametrizedModel.hh.

const G4double G4eBremParametrizedModel::xgi [static, protected]

Initial value:

{ 0.0199, 0.1017, 0.2372, 0.4083,
                                                  0.5917, 0.7628, 0.8983, 0.9801 }

Definition at line 123 of file G4eBremParametrizedModel.hh.

G4double G4eBremParametrizedModel::z13 [protected]

Definition at line 132 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

G4double G4eBremParametrizedModel::z23 [protected]

Definition at line 132 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().


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