Geant4  10.02.p03
G4SeltzerBergerModel Class Reference

#include <G4SeltzerBergerModel.hh>

Inheritance diagram for G4SeltzerBergerModel:
Collaboration diagram for G4SeltzerBergerModel:

Public Member Functions

 G4SeltzerBergerModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="eBremSB")
 
virtual ~G4SeltzerBergerModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle *> *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy)
 
void SetBicubicInterpolationFlag (G4bool)
 
- Public Member Functions inherited from G4eBremsstrahlungRelModel
 G4eBremsstrahlungRelModel (const G4ParticleDefinition *p=0, const G4String &nam="eBremLPM")
 
virtual ~G4eBremsstrahlungRelModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
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 SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut)
 
void SetLPMconstant (G4double val)
 
G4double LPMconstant () const
 
void SetLowestKinEnergy (G4double)
 
G4double LowestKinEnergy () const
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int, const G4ParticleDefinition *, G4double)
 
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 MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
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 *> *)
 
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
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 

Protected Member Functions

virtual G4double ComputeDXSectionPerAtom (G4double gammaEnergy)
 
virtual G4String DirectoryPath () const
 
- Protected Member Functions inherited from G4eBremsstrahlungRelModel
void SetCurrentElement (const G4double)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Private Member Functions

void ReadData (G4int Z, const char *path=0)
 
G4SeltzerBergerModeloperator= (const G4SeltzerBergerModel &right)
 
 G4SeltzerBergerModel (const G4SeltzerBergerModel &)
 

Private Attributes

G4int nwarn
 
size_t idx
 
size_t idy
 
G4bool useBicubicInterpolation
 

Static Private Attributes

static G4Physics2DVectordataSB [101] = {nullptr}
 
static G4double ylimit [101] = {0.0}
 
static G4double expnumlim = -12.
 

Additional Inherited Members

- Protected Attributes inherited from G4eBremsstrahlungRelModel
G4NistManagernist
 
const G4ParticleDefinitionparticle
 
G4ParticleDefinitiontheGamma
 
G4ParticleChangeForLoss * fParticleChange
 
G4double bremFactor
 
G4double particleMass
 
G4double kinEnergy
 
G4double totalEnergy
 
G4double currentZ
 
G4double densityFactor
 
G4double densityCorr
 
G4bool isElectron
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChange * pParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 61 of file G4SeltzerBergerModel.hh.

Constructor & Destructor Documentation

◆ G4SeltzerBergerModel() [1/2]

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

Definition at line 83 of file G4SeltzerBergerModel.cc.

86 {
89  SetLPMFlag(false);
90  nwarn = 0;
91  idx = idy = 0;
92 }
G4eBremsstrahlungRelModel(const G4ParticleDefinition *p=0, const G4String &nam="eBremLPM")
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:774
static const double keV
Definition: G4SIunits.hh:213
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
Here is the call graph for this function:

◆ ~G4SeltzerBergerModel()

G4SeltzerBergerModel::~G4SeltzerBergerModel ( )
virtual

Definition at line 96 of file G4SeltzerBergerModel.cc.

97 {
98  if(IsMaster()) {
99  for(size_t i=0; i<101; ++i) {
100  if(dataSB[i]) {
101  delete dataSB[i];
102  dataSB[i] = nullptr;
103  }
104  }
105  }
106 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
static G4Physics2DVector * dataSB[101]
Here is the call graph for this function:

◆ G4SeltzerBergerModel() [2/2]

G4SeltzerBergerModel::G4SeltzerBergerModel ( const G4SeltzerBergerModel )
private

Member Function Documentation

◆ ComputeDXSectionPerAtom()

G4double G4SeltzerBergerModel::ComputeDXSectionPerAtom ( G4double  gammaEnergy)
protectedvirtual

Reimplemented from G4eBremsstrahlungRelModel.

Definition at line 193 of file G4SeltzerBergerModel.cc.

194 {
195 
196  if(gammaEnergy < 0.0 || kinEnergy <= 0.0) { return 0.0; }
197  G4double x = gammaEnergy/kinEnergy;
199  G4int Z = G4lrint(currentZ);
200 
201  //G4cout << "G4SeltzerBergerModel::ComputeDXSectionPerAtom Z= " << Z
202  // << " x= " << x << " y= " << y << " " << dataSB[Z] << G4endl;
203  if(nullptr == dataSB[Z]) { InitialiseForElement(0, Z); }
204  /*
205  G4ExceptionDescription ed;
206  ed << "Bremsstrahlung data for Z= " << Z
207  << " are not initialized!";
208  G4Exception("G4SeltzerBergerModel::ComputeDXSectionPerAtom()","em0005",
209  FatalException, ed,
210  "G4LEDATA version should be G4EMLOW6.23 or later.");
211  }
212  */
213  G4double invb2 =
215  G4double cross = dataSB[Z]->Value(x,y,idx,idy)*invb2*millibarn/bremFactor;
216 
217  if(!isElectron) {
218  G4double invbeta1 = sqrt(invb2);
219  G4double e2 = kinEnergy - gammaEnergy;
220  if(e2 > 0.0) {
221  G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
222  G4double xxx = alpha*currentZ*(invbeta1 - invbeta2);
223  if(xxx < expnumlim) { cross = 0.0; }
224  else { cross *= G4Exp(xxx); }
225  } else {
226  cross = 0.0;
227  }
228  }
229 
230  return cross;
231 }
static const double MeV
Definition: G4SIunits.hh:211
static const G4double e2
int G4int
Definition: G4Types.hh:78
Double_t y
Float_t Z
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
int G4lrint(double ad)
Definition: templates.hh:163
static G4Physics2DVector * dataSB[101]
static const double millibarn
Definition: G4SIunits.hh:105
double G4double
Definition: G4Types.hh:76
G4double Value(G4double x, G4double y, size_t &lastidx, size_t &lastidy) const
static const G4double alpha
Here is the call graph for this function:

◆ DirectoryPath()

G4String G4SeltzerBergerModel::DirectoryPath ( ) const
protectedvirtual

Definition at line 139 of file G4SeltzerBergerModel.cc.

140 {
141  return "/brem_SB/br";
142 }
Here is the caller graph for this function:

◆ Initialise()

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

Reimplemented from G4eBremsstrahlungRelModel.

Reimplemented in G4ePolarizedBremsstrahlungModel.

Definition at line 110 of file G4SeltzerBergerModel.cc.

112 {
113  // Access to elements
114  if(IsMaster()) {
115 
116  // check environment variable
117  // Build the complete string identifying the file with the data set
118  char* path = getenv("G4LEDATA");
119 
120  const G4ElementTable* theElmTable = G4Element::GetElementTable();
121  size_t numOfElm = G4Element::GetNumberOfElements();
122  if(numOfElm > 0) {
123  for(size_t i=0; i<numOfElm; ++i) {
124  G4int Z = G4lrint(((*theElmTable)[i])->GetZ());
125  if(Z < 1) { Z = 1; }
126  else if(Z > 100) { Z = 100; }
127  //G4cout << "Z= " << Z << G4endl;
128  // Initialisation
129  if(nullptr == dataSB[Z]) { ReadData(Z, path); }
130  }
131  }
132  }
133 
135 }
void ReadData(G4int Z, const char *path=0)
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
int G4int
Definition: G4Types.hh:78
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
Float_t Z
int G4lrint(double ad)
Definition: templates.hh:163
static G4Physics2DVector * dataSB[101]
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InitialiseForElement()

void G4SeltzerBergerModel::InitialiseForElement ( const G4ParticleDefinition ,
G4int  Z 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 375 of file G4SeltzerBergerModel.cc.

377 {
378  G4AutoLock l(&SeltzerBergerModelMutex);
379  // G4cout << "G4SeltzerBergerModel::InitialiseForElement Z= " << Z << G4endl;
380  if(nullptr == dataSB[Z]) { ReadData(Z); }
381 }
void ReadData(G4int Z, const char *path=0)
Float_t Z
static G4Physics2DVector * dataSB[101]
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

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

◆ ReadData()

void G4SeltzerBergerModel::ReadData ( G4int  Z,
const char *  path = 0 
)
private

Definition at line 146 of file G4SeltzerBergerModel.cc.

147 {
148  // G4cout << "ReadData Z= " << Z << G4endl;
149  // G4cout << "Status for Z= " << dataSB[Z] << G4endl;
150  //if(path) { G4cout << path << G4endl; }
151  if(dataSB[Z]) { return; }
152  const char* datadir = path;
153 
154  if(!datadir) {
155  datadir = getenv("G4LEDATA");
156  if(!datadir) {
157  G4Exception("G4SeltzerBergerModel::ReadData()","em0006",FatalException,
158  "Environment variable G4LEDATA not defined");
159  return;
160  }
161  }
162  std::ostringstream ost;
163  ost << datadir << DirectoryPath() << Z;
164  std::ifstream fin(ost.str().c_str());
165  if( !fin.is_open()) {
167  ed << "Bremsstrahlung data file <" << ost.str().c_str()
168  << "> is not opened!";
169  G4Exception("G4SeltzerBergerModel::ReadData()","em0003",FatalException,
170  ed,"G4LEDATA version should be G4EMLOW6.23 or later.");
171  return;
172  }
173  //G4cout << "G4SeltzerBergerModel read from <" << ost.str().c_str()
174  // << ">" << G4endl;
176  if(v->Retrieve(fin)) {
178  dataSB[Z] = v;
179  ylimit[Z] = v->Value(0.97, emaxlog, idx, idy);
180  } else {
182  ed << "Bremsstrahlung data file <" << ost.str().c_str()
183  << "> is not retrieved!";
184  G4Exception("G4SeltzerBergerModel::ReadData()","em0005",FatalException,
185  ed,"G4LEDATA version should be G4EMLOW6.23 or later.");
186  delete v;
187  }
188  // G4cout << dataSB[Z] << G4endl;
189 }
TString fin
void SetBicubicInterpolation(G4bool)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
virtual G4String DirectoryPath() const
Float_t Z
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4double ylimit[101]
G4bool Retrieve(std::ifstream &fIn)
static G4Physics2DVector * dataSB[101]
G4double Value(G4double x, G4double y, size_t &lastidx, size_t &lastidy) const
static const G4double emaxlog
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SampleSecondaries()

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

Reimplemented from G4eBremsstrahlungRelModel.

Reimplemented in G4ePolarizedBremsstrahlungModel.

Definition at line 236 of file G4SeltzerBergerModel.cc.

241 {
242  G4double kineticEnergy = dp->GetKineticEnergy();
243  G4double cut = std::min(cutEnergy, kineticEnergy);
244  G4double emax = std::min(maxEnergy, kineticEnergy);
245  if(cut >= emax) { return; }
246 
247  SetupForMaterial(particle, couple->GetMaterial(), kineticEnergy);
248 
249  const G4Element* elm =
250  SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
251  SetCurrentElement(elm->GetZ());
252  G4int Z = G4int(currentZ);
253 
254  totalEnergy = kineticEnergy + particleMass;
256  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
257  /*
258  G4cout << "G4SeltzerBergerModel::SampleSecondaries E(MeV)= "
259  << kineticEnergy/MeV
260  << " Z= " << Z << " cut(MeV)= " << cut/MeV
261  << " emax(MeV)= " << emax/MeV << " corr= " << densityCorr << G4endl;
262  */
263  G4double xmin = G4Log(cut*cut + densityCorr);
264  G4double xmax = G4Log(emax*emax + densityCorr);
265  G4double y = G4Log(kineticEnergy/MeV);
266 
267  G4double gammaEnergy, v;
268 
269  // majoranta
270  G4double x0 = cut/kineticEnergy;
271  G4double vmax = dataSB[Z]->Value(x0, y, idx, idy)*1.02;
272  // G4double invbeta1 = 0;
273 
274  // majoranta corrected for e-
275  if(isElectron && x0 < 0.97 &&
276  ((kineticEnergy > epeaklimit) || (kineticEnergy < elowlimit))) {
277  G4double ylim = std::min(ylimit[Z],1.1*dataSB[Z]->Value(0.97,y,idx,idy));
278  if(ylim > vmax) { vmax = ylim; }
279  }
280  if(x0 < 0.05) { vmax *= 1.2; }
281 
282  //G4cout<<"y= "<<y<<" xmin= "<<xmin<<" xmax= "<<xmax
283  //<<" vmax= "<<vmax<<G4endl;
284  static const G4int ncountmax = 100;
285  CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
286  G4double rndm[2];
287 
288  for(G4int nn=0; nn<ncountmax; ++nn) {
289  rndmEngine->flatArray(2, rndm);
290  G4double x = G4Exp(xmin + rndm[0]*(xmax - xmin)) - densityCorr;
291  if(x < 0.0) { x = 0.0; }
292  gammaEnergy = sqrt(x);
293  G4double x1 = gammaEnergy/kineticEnergy;
294  v = dataSB[Z]->Value(x1, y, idx, idy);
295 
296  // correction for positrons
297  if(!isElectron) {
298  G4double e1 = kineticEnergy - cut;
299  G4double invbeta1 = (e1 + particleMass)/sqrt(e1*(e1 + 2*particleMass));
300  G4double e2 = kineticEnergy - gammaEnergy;
301  G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
302  G4double xxx = twopi*fine_structure_const*currentZ*(invbeta1 - invbeta2);
303 
304  if(xxx < expnumlim) { v = 0.0; }
305  else { v *= G4Exp(xxx); }
306  }
307 
308  if (v > 1.05*vmax && nwarn < 5) {
309  ++nwarn;
311  ed << "### G4SeltzerBergerModel Warning: Majoranta exceeded! "
312  << v << " > " << vmax << " by " << v/vmax
313  << " Niter= " << nn
314  << " Egamma(MeV)= " << gammaEnergy
315  << " Ee(MeV)= " << kineticEnergy
316  << " Z= " << Z << " " << particle->GetParticleName();
317 
318  if ( 20 == nwarn ) {
319  ed << "\n ### G4SeltzerBergerModel Warnings stopped";
320  }
321  G4Exception("G4SeltzerBergerModel::SampleScattering","em0044",
322  JustWarning, ed,"");
323 
324  }
325  if(v >= vmax*rndm[1]) { break; }
326  }
327 
328  //
329  // angles of the emitted gamma. ( Z - axis along the parent particle)
330  // use general interface
331  //
332 
333  G4ThreeVector gammaDirection =
334  GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy,
335  Z, couple->GetMaterial());
336 
337  // create G4DynamicParticle object for the Gamma
338  G4DynamicParticle* gamma =
339  new G4DynamicParticle(theGamma,gammaDirection,gammaEnergy);
340  vdp->push_back(gamma);
341 
342  G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
343  - gammaEnergy*gammaDirection).unit();
344 
345  /*
346  G4cout << "### G4SBModel: v= "
347  << " Eg(MeV)= " << gammaEnergy
348  << " Ee(MeV)= " << kineticEnergy
349  << " DirE " << direction << " DirG " << gammaDirection
350  << G4endl;
351  */
352  // energy of primary
353  G4double finalE = kineticEnergy - gammaEnergy;
354 
355  // stop tracking and create new secondary instead of primary
356  if(gammaEnergy > SecondaryThreshold()) {
357  fParticleChange->ProposeTrackStatus(fStopAndKill);
358  fParticleChange->SetProposedKineticEnergy(0.0);
359  G4DynamicParticle* el =
360  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
361  direction, finalE);
362  vdp->push_back(el);
363 
364  // continue tracking
365  } else {
366  fParticleChange->SetProposedMomentumDirection(direction);
367  fParticleChange->SetProposedKineticEnergy(finalE);
368  }
369 }
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:669
static const double MeV
Definition: G4SIunits.hh:211
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const G4Material * GetMaterial() const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:617
static const G4double e2
int G4int
Definition: G4Types.hh:78
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
const G4ParticleDefinition * particle
Double_t y
int fine_structure_const
Definition: hepunit.py:287
G4double GetKineticEnergy() const
const G4String & GetParticleName() const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
Float_t Z
static const G4double elowlimit
static const double twopi
Definition: G4SIunits.hh:75
Double_t x1[nxs]
float electron_mass_c2
Definition: hepunit.py:274
static const G4double e1
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const G4double emax
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static G4double ylimit[101]
const G4ThreeVector & GetMomentumDirection() const
static G4Physics2DVector * dataSB[101]
G4ParticleChangeForLoss * fParticleChange
static const G4double epeaklimit
double G4double
Definition: G4Types.hh:76
G4double Value(G4double x, G4double y, size_t &lastidx, size_t &lastidy) const
virtual void flatArray(const int size, double *vect)=0
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:369
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetBicubicInterpolationFlag()

void G4SeltzerBergerModel::SetBicubicInterpolationFlag ( G4bool  val)
inline

Definition at line 106 of file G4SeltzerBergerModel.hh.

107 {
109 }

Member Data Documentation

◆ dataSB

G4Physics2DVector * G4SeltzerBergerModel::dataSB = {nullptr}
staticprivate

Definition at line 97 of file G4SeltzerBergerModel.hh.

◆ expnumlim

G4double G4SeltzerBergerModel::expnumlim = -12.
staticprivate

Definition at line 99 of file G4SeltzerBergerModel.hh.

◆ idx

size_t G4SeltzerBergerModel::idx
private

Definition at line 101 of file G4SeltzerBergerModel.hh.

◆ idy

size_t G4SeltzerBergerModel::idy
private

Definition at line 102 of file G4SeltzerBergerModel.hh.

◆ nwarn

G4int G4SeltzerBergerModel::nwarn
private

Definition at line 100 of file G4SeltzerBergerModel.hh.

◆ useBicubicInterpolation

G4bool G4SeltzerBergerModel::useBicubicInterpolation
private

Definition at line 103 of file G4SeltzerBergerModel.hh.

◆ ylimit

G4double G4SeltzerBergerModel::ylimit = {0.0}
staticprivate

Definition at line 98 of file G4SeltzerBergerModel.hh.


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