Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Attributes | Static Protected Attributes
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 void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
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)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

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
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Static Protected Attributes

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

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

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 100 of file G4eBremParametrizedModel.cc.

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

102  : G4VEmModel(nam),
103  particle(0),
104  isElectron(true),
107  isInitialised(false)
108 {
110 
111  minThreshold = 0.1*keV;
112  lowKinEnergy = 10.*MeV;
113  SetLowEnergyLimit(lowKinEnergy);
114 
116 
118 
120  = densityFactor = densityCorr = fMax = fCoulomb = 0.;
121 
122  InitialiseConstants();
123  if(p) { SetParticle(p); }
124 }
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
electron_Compton_length
Definition: hepunit.py:289
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:585
const G4ParticleDefinition * particle
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
G4eBremParametrizedModel::~G4eBremParametrizedModel ( )
virtual

Definition at line 136 of file G4eBremParametrizedModel.cc.

137 {
138 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 267 of file G4eBremParametrizedModel.cc.

References kinEnergy, G4INCL::Math::min(), and particle.

273 {
274  if(!particle) { SetParticle(p); }
275  if(kineticEnergy < lowKinEnergy) { return 0.0; }
276  G4double cut = std::min(cutEnergy, kineticEnergy);
277  G4double tmax = std::min(maxEnergy, kineticEnergy);
278 
279  if(cut >= tmax) { return 0.0; }
280 
281  SetCurrentElement(Z);
282 
283  G4double cross = ComputeXSectionPerAtom(cut);
284 
285  // allow partial integration
286  if(tmax < kinEnergy) { cross -= ComputeXSectionPerAtom(tmax); }
287 
288  cross *= Z*Z*bremFactor;
289 
290  return cross;
291 }
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
const G4ParticleDefinition * particle
double G4double
Definition: G4Types.hh:76
G4double G4eBremParametrizedModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 201 of file G4eBremParametrizedModel.cc.

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

206 {
207  if(!particle) { SetParticle(p); }
208  if(kineticEnergy < lowKinEnergy) { return 0.0; }
209  G4double cut = std::min(cutEnergy, kineticEnergy);
210  if(cut == 0.0) { return 0.0; }
211 
212  SetupForMaterial(particle, material,kineticEnergy);
213 
214  const G4ElementVector* theElementVector = material->GetElementVector();
215  const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector();
216 
217  G4double dedx = 0.0;
218 
219  // loop for elements in the material
220  for (size_t i=0; i<material->GetNumberOfElements(); i++) {
221 
222  G4VEmModel::SetCurrentElement((*theElementVector)[i]);
223  SetCurrentElement((*theElementVector)[i]->GetZ());
224 
225  dedx += theAtomicNumDensityVector[i]*currentZ*currentZ*ComputeBremLoss(cut);
226  }
227  dedx *= bremFactor;
228 
229  return dedx;
230 }
std::vector< G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:214
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
const G4ParticleDefinition * particle
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
double G4double
Definition: G4Types.hh:76
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:433
void G4eBremParametrizedModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 175 of file G4eBremParametrizedModel.cc.

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

177 {
178  if(p) { SetParticle(p); }
179 
180  lowKinEnergy = LowEnergyLimit();
181 
182  currentZ = 0.;
183 
184  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
185 
186  if(isInitialised) { return; }
188  isInitialised = true;
189 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
G4ParticleChangeForLoss * fParticleChange
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
void G4eBremParametrizedModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 193 of file G4eBremParametrizedModel.cc.

References G4VEmModel::GetElementSelectors(), and G4VEmModel::SetElementSelectors().

195 {
197 }
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:768
G4double G4eBremParametrizedModel::MinEnergyCut ( const G4ParticleDefinition ,
const G4MaterialCutsCouple  
)
virtual

Reimplemented from G4VEmModel.

Definition at line 152 of file G4eBremParametrizedModel.cc.

References minThreshold.

154 {
155  return minThreshold;
156 }
void G4eBremParametrizedModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  vdp,
const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 474 of file G4eBremParametrizedModel.cc.

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

480 {
481  G4double kineticEnergy = dp->GetKineticEnergy();
482  if(kineticEnergy < lowKinEnergy) { return; }
483  G4double cut = std::min(cutEnergy, kineticEnergy);
484  G4double emax = std::min(maxEnergy, kineticEnergy);
485  if(cut >= emax) { return; }
486 
487  SetupForMaterial(particle, couple->GetMaterial(),kineticEnergy);
488 
489  const G4Element* elm =
490  SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
491  SetCurrentElement(elm->GetZ());
492 
493  kinEnergy = kineticEnergy;
494  totalEnergy = kineticEnergy + particleMass;
496 
497  G4double xmin = G4Log(cut*cut + densityCorr);
498  G4double xmax = G4Log(emax*emax + densityCorr);
499  G4double gammaEnergy, f, x;
500 
501  do {
502  x = G4Exp(xmin + G4UniformRand()*(xmax - xmin)) - densityCorr;
503  if(x < 0.0) x = 0.0;
504  gammaEnergy = sqrt(x);
505  f = ComputeDXSectionPerAtom(gammaEnergy);
506 
507  if ( f > fMax ) {
508  G4cout << "### G4eBremParametrizedModel Warning: Majoranta exceeded! "
509  << f << " > " << fMax
510  << " Egamma(MeV)= " << gammaEnergy
511  << " E(mEV)= " << kineticEnergy
512  << G4endl;
513  }
514 
515  } while (f < fMax*G4UniformRand());
516 
517  //
518  // angles of the emitted gamma. ( Z - axis along the parent particle)
519  // use general interface
520  //
521  G4ThreeVector gammaDirection =
522  GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy,
523  G4lrint(currentZ),
524  couple->GetMaterial());
525 
526  // create G4DynamicParticle object for the Gamma
527  G4DynamicParticle* gamma = new G4DynamicParticle(theGamma,gammaDirection,
528  gammaEnergy);
529  vdp->push_back(gamma);
530 
531  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
532  G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
533  - gammaEnergy*gammaDirection).unit();
534 
535  // energy of primary
536  G4double finalE = kineticEnergy - gammaEnergy;
537 
538  // stop tracking and create new secondary instead of primary
539  if(gammaEnergy > SecondaryThreshold()) {
542  G4DynamicParticle* el =
543  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
544  direction, finalE);
545  vdp->push_back(el);
546 
547  // continue tracking
548  } else {
551  }
552 }
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:627
G4double GetKineticEnergy() const
G4ParticleChangeForLoss * fParticleChange
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:578
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
const G4ThreeVector & GetMomentumDirection() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
float electron_mass_c2
Definition: hepunit.py:274
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
int G4lrint(double ad)
Definition: templates.hh:163
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
const G4ParticleDefinition * particle
#define G4endl
Definition: G4ios.hh:61
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:510
const G4Material * GetMaterial() const
void G4eBremParametrizedModel::SetupForMaterial ( const G4ParticleDefinition ,
const G4Material mat,
G4double  kineticEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 160 of file G4eBremParametrizedModel.cc.

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

Referenced by ComputeDEDXPerVolume(), and SampleSecondaries().

163 {
164  densityFactor = mat->GetElectronDensity()*fMigdalConstant;
165 
166  // calculate threshold for density effect
167  kinEnergy = kineticEnergy;
168  totalEnergy = kineticEnergy + particleMass;
170 }
G4double GetElectronDensity() const
Definition: G4Material.hh:215

Field Documentation

G4double G4eBremParametrizedModel::currentZ
protected
G4double G4eBremParametrizedModel::densityCorr
protected
G4double G4eBremParametrizedModel::densityFactor
protected
G4double G4eBremParametrizedModel::facFel
protected

Definition at line 147 of file G4eBremParametrizedModel.hh.

G4double G4eBremParametrizedModel::facFinel
protected

Definition at line 147 of file G4eBremParametrizedModel.hh.

G4double G4eBremParametrizedModel::fCoulomb
protected

Definition at line 148 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

G4double G4eBremParametrizedModel::Fel
protected

Definition at line 146 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

G4double G4eBremParametrizedModel::Finel
protected

Definition at line 146 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

G4double G4eBremParametrizedModel::fMax
protected

Definition at line 148 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), and SampleSecondaries().

G4ParticleChangeForLoss* G4eBremParametrizedModel::fParticleChange
protected

Definition at line 132 of file G4eBremParametrizedModel.hh.

Referenced by Initialise(), and SampleSecondaries().

G4bool G4eBremParametrizedModel::isElectron
protected

Definition at line 150 of file G4eBremParametrizedModel.hh.

G4double G4eBremParametrizedModel::kinEnergy
protected
G4double G4eBremParametrizedModel::lnZ
protected

Definition at line 143 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

G4double G4eBremParametrizedModel::minThreshold
protected

Definition at line 136 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), and MinEnergyCut().

G4NistManager* G4eBremParametrizedModel::nist
protected

Definition at line 129 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

const G4ParticleDefinition* G4eBremParametrizedModel::particle
protected
G4double G4eBremParametrizedModel::particleMass
protected
G4ParticleDefinition* G4eBremParametrizedModel::theGamma
protected

Definition at line 131 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel(), and SampleSecondaries().

G4double G4eBremParametrizedModel::totalEnergy
protected
const G4double G4eBremParametrizedModel::wgi
staticprotected
Initial value:
={ 0.0506, 0.1112, 0.1569, 0.1813,
0.1813, 0.1569, 0.1112, 0.0506 }

Definition at line 134 of file G4eBremParametrizedModel.hh.

const G4double G4eBremParametrizedModel::xgi
staticprotected
Initial value:
={ 0.0199, 0.1017, 0.2372, 0.4083,
0.5917, 0.7628, 0.8983, 0.9801 }

Definition at line 134 of file G4eBremParametrizedModel.hh.

G4double G4eBremParametrizedModel::z13
protected

Definition at line 143 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().

G4double G4eBremParametrizedModel::z23
protected

Definition at line 143 of file G4eBremParametrizedModel.hh.

Referenced by G4eBremParametrizedModel().


The documentation for this class was generated from the following files: