Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4SeltzerBergerModel Class Reference

#include <G4SeltzerBergerModel.hh>

+ Inheritance diagram for G4SeltzerBergerModel:

Public Member Functions

 G4SeltzerBergerModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="eBremSB")
 
 ~G4SeltzerBergerModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
 
void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double) override
 
void SetBicubicInterpolationFlag (G4bool val)
 
G4SeltzerBergerModeloperator= (const G4SeltzerBergerModel &right)=delete
 
 G4SeltzerBergerModel (const G4SeltzerBergerModel &)=delete
 
- Public Member Functions inherited from G4eBremsstrahlungRelModel
 G4eBremsstrahlungRelModel (const G4ParticleDefinition *p=0, const G4String &nam="eBremLPM")
 
virtual ~G4eBremsstrahlungRelModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double ekin, G4double cutEnergy) override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double ekin, G4double zet, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double) override
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cutEnergy) override
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, 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 G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, 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 *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
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) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

G4double ComputeDXSectionPerAtom (G4double gammaEnergy) override
 
- Protected Member Functions inherited from G4eBremsstrahlungRelModel
virtual G4double ComputeDXSectionPerAtom (G4double gammaEnergy)
 
void SetParticle (const G4ParticleDefinition *p)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Additional Inherited Members

- Protected Attributes inherited from G4eBremsstrahlungRelModel
G4bool fIsElectron
 
G4bool fIsScatOffElectron
 
G4bool fIsLPMActive
 
G4int fCurrentIZ
 
G4double fPrimaryParticleMass
 
G4double fPrimaryKinEnergy
 
G4double fPrimaryTotalEnergy
 
G4double fDensityFactor
 
G4double fDensityCorr
 
G4double fLowestKinEnergy
 
G4double fNucTerm
 
G4double fSumTerm
 
const G4ParticleDefinitionfPrimaryParticle
 
G4ParticleDefinitionfGammaParticle
 
G4ParticleChangeForLossfParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 
- Static Protected Attributes inherited from G4eBremsstrahlungRelModel
static const G4double gBremFactor
 
static const G4double gMigdalConstant
 

Detailed Description

Definition at line 69 of file G4SeltzerBergerModel.hh.

Constructor & Destructor Documentation

◆ G4SeltzerBergerModel() [1/2]

G4SeltzerBergerModel::G4SeltzerBergerModel ( const G4ParticleDefinition p = nullptr,
const G4String nam = "eBremSB" 
)
explicit

Definition at line 88 of file G4SeltzerBergerModel.cc.

90: G4eBremsstrahlungRelModel(p,nam), fIsUseBicubicInterpolation(false),
91 fIsUseSamplingTables(true), fNumWarnings(0), fIndx(0), fIndy(0)
92{
93 fLowestKinEnergy = 1.0*keV;
95 SetLPMFlag(false);
97}
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:806
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618

◆ ~G4SeltzerBergerModel()

G4SeltzerBergerModel::~G4SeltzerBergerModel ( )
override

Definition at line 99 of file G4SeltzerBergerModel.cc.

100{
101 // delete SB-DCS data per Z
102 if (IsMaster()) {
103 for (size_t iz = 0; iz < gMaxZet; ++iz) {
104 if (gSBDCSData[iz]) {
105 delete gSBDCSData[iz];
106 gSBDCSData[iz] = nullptr;
107 }
108 }
109 if (gSBSamplingTable) {
110 delete gSBSamplingTable;
111 gSBSamplingTable = nullptr;
112 }
113 }
114}
G4bool IsMaster() const
Definition: G4VEmModel.hh:736

◆ G4SeltzerBergerModel() [2/2]

G4SeltzerBergerModel::G4SeltzerBergerModel ( const G4SeltzerBergerModel )
delete

Member Function Documentation

◆ ComputeDXSectionPerAtom()

G4double G4SeltzerBergerModel::ComputeDXSectionPerAtom ( G4double  gammaEnergy)
overrideprotectedvirtual

Reimplemented from G4eBremsstrahlungRelModel.

Definition at line 218 of file G4SeltzerBergerModel.cc.

219{
220 G4double dxsec = 0.0;
221 if (gammaEnergy < 0.0 || fPrimaryKinEnergy <= 0.0) {
222 return dxsec;
223 }
224 // reduced photon energy
225 const G4double x = gammaEnergy/fPrimaryKinEnergy;
226 // l-kinetic energy of the e-/e+
227 const G4double y = G4Log(fPrimaryKinEnergy/CLHEP::MeV);
228 // make sure that the Z-related SB-DCS are loaded
229 // NOTE: fCurrentIZ should have been set before.
230 fCurrentIZ = std::max(std::min(fCurrentIZ, gMaxZet-1), 1);
231 if (nullptr == gSBDCSData[fCurrentIZ]) {
232 ReadData(fCurrentIZ);
233 }
234 // NOTE: SetupForMaterial should have been called before!
235 const G4double pt2 = fPrimaryKinEnergy*(fPrimaryKinEnergy+2.*kMC2);
237 G4double val = gSBDCSData[fCurrentIZ]->Value(x,y,fIndx,fIndy);
238 dxsec = val*invb2*CLHEP::millibarn/gBremFactor;
239 // e+ correction
240 if (!fIsElectron) {
241 const G4double invbeta1 = std::sqrt(invb2);
242 const G4double e2 = fPrimaryKinEnergy-gammaEnergy;
243 if (e2 > 0.0) {
244 const G4double invbeta2 = (e2+kMC2)/std::sqrt(e2*(e2+2.0*kMC2));
245 const G4double dum0 = kAlpha*fCurrentIZ*(invbeta1-invbeta2);
246 if (dum0 < gExpNumLimit) {
247 dxsec = 0.0;
248 } else {
249 dxsec *= G4Exp(dum0);
250 }
251 } else {
252 dxsec = 0.0;
253 }
254 }
255 return dxsec;
256}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
double G4double
Definition: G4Types.hh:83
G4double Value(G4double x, G4double y, std::size_t &lastidx, std::size_t &lastidy) const

◆ Initialise()

void G4SeltzerBergerModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
overridevirtual

Reimplemented from G4eBremsstrahlungRelModel.

Definition at line 116 of file G4SeltzerBergerModel.cc.

118{
119 if (p) {
120 SetParticle(p);
121 }
122 fIsUseSamplingTables = G4EmParameters::Instance()->EnableSamplingTable();
123 // Access to elements
124 if (IsMaster()) {
125
126 auto theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
127 size_t numOfCouples = theCoupleTable->GetTableSize();
128 for(size_t j=0; j<numOfCouples; ++j) {
129 auto mat = theCoupleTable->GetMaterialCutsCouple(j)->GetMaterial();
130 auto elmVec = mat->GetElementVector();
131 for (auto & elm : *elmVec) {
132 G4int Z = std::max(1,std::min(elm->GetZasInt(), gMaxZet-1));
133 // load SB-DCS data for this atomic number if it has not been loaded yet
134 if (gSBDCSData[Z] == nullptr) ReadData(Z);
135 }
136 }
137 // elem.selectr. only for master: base class init-local will set for workers
140 }
141 // init sampling tables if it was requested
142 if (fIsUseSamplingTables) {
143 if (!gSBSamplingTable) {
144 gSBSamplingTable = new G4SBBremTable();
145 }
146 gSBSamplingTable->Initialize(std::max(fLowestKinEnergy,LowEnergyLimit()),
148 }
149 }
150 //
152 if (GetTripletModel()) {
153 GetTripletModel()->Initialise(p, cuts);
154 fIsScatOffElectron = true;
155 }
156}
int G4int
Definition: G4Types.hh:85
static G4EmParameters * Instance()
G4bool EnableSamplingTable() const
static G4ProductionCutsTable * GetProductionCutsTable()
void Initialize(const G4double lowe, const G4double highe)
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4VEmModel * GetTripletModel()
Definition: G4VEmModel.hh:628
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118
void SetParticle(const G4ParticleDefinition *p)
G4ParticleChangeForLoss * fParticleChange

Referenced by G4ePolarizedBremsstrahlungModel::Initialise().

◆ operator=()

G4SeltzerBergerModel & G4SeltzerBergerModel::operator= ( const G4SeltzerBergerModel right)
delete

◆ SampleSecondaries()

void G4SeltzerBergerModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  vdp,
const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  cutEnergy,
G4double  maxEnergy 
)
overridevirtual

Reimplemented from G4eBremsstrahlungRelModel.

Definition at line 259 of file G4SeltzerBergerModel.cc.

264{
265 const G4double kinEnergy = dp->GetKineticEnergy();
266 const G4double logKinEnergy = dp->GetLogKineticEnergy();
267 const G4double tmin = std::min(cutEnergy, kinEnergy);
268 const G4double tmax = std::min(maxEnergy, kinEnergy);
269 if (tmin >= tmax) {
270 return;
271 }
272 // set local variables and select target element
273 SetupForMaterial(fPrimaryParticle, couple->GetMaterial(), kinEnergy);
274 const G4Element* elm = SelectTargetAtom(couple, fPrimaryParticle, kinEnergy,
275 logKinEnergy, tmin, tmax);
276 fCurrentIZ = std::max(std::min(elm->GetZasInt(),gMaxZet-1), 1);
277 //
278 const G4double totMomentum = std::sqrt(kinEnergy*(fPrimaryTotalEnergy+kMC2));
279 /*
280 G4cout << "G4SeltzerBergerModel::SampleSecondaries E(MeV)= "
281 << kinEnergy/MeV
282 << " Z= " << fCurrentIZ << " cut(MeV)= " << tmin/MeV
283 << " emax(MeV)= " << tmax/MeV << " corr= " << fDensityCorr << G4endl;
284 */
285 // sample emitted photon energy either by rejection or from samplign tables
286 const G4double gammaEnergy = !fIsUseSamplingTables
287 ? SampleEnergyTransfer(kinEnergy, logKinEnergy, tmin, tmax)
288 : gSBSamplingTable->SampleEnergyTransfer(kinEnergy, logKinEnergy, tmin,
289 fDensityCorr, fCurrentIZ, couple->GetIndex(), fIsElectron);
290 // should never happen under normal conditions but protect it
291 if (gammaEnergy <= 0.) {
292 return;
293 }
294 //
295 // angles of the emitted gamma. ( Z - axis along the parent particle) use
296 // general interface
298 fPrimaryTotalEnergy-gammaEnergy, fCurrentIZ, couple->GetMaterial());
299 // create G4DynamicParticle object for the emitted Gamma
301 gammaEnergy);
302 vdp->push_back(gamma);
303 //
304 // compute post-interaction kinematics of the primary e-/e+
305 G4ThreeVector dir =
306 (totMomentum*dp->GetMomentumDirection()-gammaEnergy*gamDir).unit();
307 const G4double finalE = kinEnergy - gammaEnergy;
308 /*
309 G4cout << "### G4SBModel: v= "
310 << " Eg(MeV)= " << gammaEnergy
311 << " Ee(MeV)= " << kineticEnergy
312 << " DirE " << direction << " DirG " << gammaDirection
313 << G4endl;
314 */
315 // if secondary gamma energy is higher than threshold(very high by default)
316 // then stop tracking the primary particle and create new secondary e-/e+
317 // instead of the primary
318 if (gammaEnergy > SecondaryThreshold()) {
322 const_cast<G4ParticleDefinition*>(fPrimaryParticle), dir, finalE);
323 vdp->push_back(el);
324 } else { // continue tracking the primary e-/e+ otherwise
327 }
328}
@ fStopAndKill
const G4ThreeVector & GetMomentumDirection() const
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
G4int GetZasInt() const
Definition: G4Element.hh:131
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:611
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:587
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:680
void ProposeTrackStatus(G4TrackStatus status)
const G4ParticleDefinition * fPrimaryParticle
G4ParticleDefinition * fGammaParticle

Referenced by G4ePolarizedBremsstrahlungModel::SampleSecondaries().

◆ SetBicubicInterpolationFlag()

void G4SeltzerBergerModel::SetBicubicInterpolationFlag ( G4bool  val)
inline

Definition at line 90 of file G4SeltzerBergerModel.hh.

91 { fIsUseBicubicInterpolation = val; };

◆ SetupForMaterial()

void G4SeltzerBergerModel::SetupForMaterial ( const G4ParticleDefinition ,
const G4Material mat,
G4double  kineticEnergy 
)
overridevirtual

Reimplemented from G4eBremsstrahlungRelModel.

Definition at line 409 of file G4SeltzerBergerModel.cc.

412{
414 // calculate threshold for density effect: gamma*k_p = sqrt(fDensityCorr)
415 fPrimaryKinEnergy = kineticEnergy;
416 fPrimaryTotalEnergy = kineticEnergy+CLHEP::electron_mass_c2;
419}
G4double GetElectronDensity() const
Definition: G4Material.hh:215
G4bool LPMFlag() const
Definition: G4VEmModel.hh:687

Referenced by SampleSecondaries().


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