Geant4  10.02.p03
G4LivermoreBremsstrahlungModel Class Reference

#include <G4LivermoreBremsstrahlungModel.hh>

Inheritance diagram for G4LivermoreBremsstrahlungModel:
Collaboration diagram for G4LivermoreBremsstrahlungModel:

Public Member Functions

 G4LivermoreBremsstrahlungModel (const G4ParticleDefinition *p=0, const G4String &nam="LowEnBrem")
 
virtual ~G4LivermoreBremsstrahlungModel ()
 
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)
 
G4LivermoreBremsstrahlungModeloperator= (const G4LivermoreBremsstrahlungModel &right)
 
 G4LivermoreBremsstrahlungModel (const G4LivermoreBremsstrahlungModel &)
 

Private Attributes

G4int nwarn
 
size_t idx
 
size_t idy
 
G4bool useBicubicInterpolation
 

Static Private Attributes

static G4Physics2DVectordataSB [101] = {0}
 
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 G4LivermoreBremsstrahlungModel.hh.

Constructor & Destructor Documentation

◆ G4LivermoreBremsstrahlungModel() [1/2]

G4LivermoreBremsstrahlungModel::G4LivermoreBremsstrahlungModel ( const G4ParticleDefinition p = 0,
const G4String nam = "LowEnBrem" 
)

Definition at line 83 of file G4LivermoreBremsstrahlungModel.cc.

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

◆ ~G4LivermoreBremsstrahlungModel()

G4LivermoreBremsstrahlungModel::~G4LivermoreBremsstrahlungModel ( )
virtual

Definition at line 96 of file G4LivermoreBremsstrahlungModel.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] = 0;
103  }
104  }
105  }
106 }
static G4Physics2DVector * dataSB[101]
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
Here is the call graph for this function:

◆ G4LivermoreBremsstrahlungModel() [2/2]

G4LivermoreBremsstrahlungModel::G4LivermoreBremsstrahlungModel ( const G4LivermoreBremsstrahlungModel )
private

Member Function Documentation

◆ ComputeDXSectionPerAtom()

G4double G4LivermoreBremsstrahlungModel::ComputeDXSectionPerAtom ( G4double  gammaEnergy)
protectedvirtual

Reimplemented from G4eBremsstrahlungRelModel.

Definition at line 196 of file G4LivermoreBremsstrahlungModel.cc.

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

Definition at line 139 of file G4LivermoreBremsstrahlungModel.cc.

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

◆ Initialise()

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

Reimplemented from G4eBremsstrahlungRelModel.

Definition at line 110 of file G4LivermoreBremsstrahlungModel.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 = G4int(((*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(!dataSB[Z]) { ReadData(Z, path); }
130  }
131  }
132  }
133 
135 }
static G4Physics2DVector * dataSB[101]
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
int G4int
Definition: G4Types.hh:78
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
Float_t Z
void ReadData(G4int Z, const char *path=0)
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:

◆ InitialiseForElement()

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

Reimplemented from G4VEmModel.

Definition at line 374 of file G4LivermoreBremsstrahlungModel.cc.

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

◆ operator=()

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

◆ ReadData()

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

Definition at line 146 of file G4LivermoreBremsstrahlungModel.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("G4LivermoreBremsstrahlungModel::ReadData()","em0006",
158  FatalException,"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("G4LivermoreBremsstrahlungModel::ReadData()","em0003",
170  FatalException,ed,
171  "G4LEDATA version should be G4EMLOW6.23 or later.");
172  return;
173  }
174  //G4cout << "G4LivermoreBremsstrahlungModel read from <" << ost.str().c_str()
175  // << ">" << G4endl;
177  if(v->Retrieve(fin)) {
179  dataSB[Z] = v;
180  ylimit[Z] = v->Value(0.97, emaxlog, idx, idy);
181  } else {
183  ed << "Bremsstrahlung data file <" << ost.str().c_str()
184  << "> is not retrieved!";
185  G4Exception("G4LivermoreBremsstrahlungModel::ReadData()","em0005",
186  FatalException,ed,
187  "G4LEDATA version should be G4EMLOW6.23 or later.");
188  delete v;
189  }
190  // G4cout << dataSB[Z] << G4endl;
191 }
TString fin
void SetBicubicInterpolation(G4bool)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static G4Physics2DVector * dataSB[101]
static const G4double emaxlog
Float_t Z
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4bool Retrieve(std::ifstream &fIn)
G4double Value(G4double x, G4double y, size_t &lastidx, size_t &lastidy) const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SampleSecondaries()

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

Reimplemented from G4eBremsstrahlungRelModel.

Definition at line 239 of file G4LivermoreBremsstrahlungModel.cc.

245 {
246  G4double kineticEnergy = dp->GetKineticEnergy();
247  G4double cut = std::min(cutEnergy, kineticEnergy);
248  G4double emax = std::min(maxEnergy, kineticEnergy);
249  if(cut >= emax) { return; }
250 
251  SetupForMaterial(particle, couple->GetMaterial(), kineticEnergy);
252 
253  const G4Element* elm =
254  SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
255  SetCurrentElement(elm->GetZ());
256  G4int Z = G4int(currentZ);
257 
258  totalEnergy = kineticEnergy + particleMass;
260  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
261  /*
262  G4cout << "G4LivermoreBremsstrahlungModel::SampleSecondaries E(MeV)= "
263  << kineticEnergy/MeV
264  << " Z= " << Z << " cut(MeV)= " << cut/MeV
265  << " emax(MeV)= " << emax/MeV << " corr= " << densityCorr << G4endl;
266  */
267  G4double xmin = G4Log(cut*cut + densityCorr);
268  G4double xmax = G4Log(emax*emax + densityCorr);
269  G4double y = G4Log(kineticEnergy/MeV);
270 
271  G4double gammaEnergy, v;
272 
273  // majoranta
274  G4double x0 = cut/kineticEnergy;
275  G4double vmax = dataSB[Z]->Value(x0, y, idx, idy)*1.02;
276  // G4double invbeta1 = 0;
277 
278  // majoranta corrected for e-
279  if(isElectron && x0 < 0.97 &&
280  ((kineticEnergy > epeaklimit) || (kineticEnergy < elowlimit))) {
281  G4double ylim = std::min(ylimit[Z],1.1*dataSB[Z]->Value(0.97,y,idx,idy));
282  if(ylim > vmax) { vmax = ylim; }
283  }
284  if(x0 < 0.05) { vmax *= 1.2; }
285 
286  //G4cout<<"y= "<<y<<" xmin= "<<xmin<<" xmax= "<<xmax
287  //<<" vmax= "<<vmax<<G4endl;
288  // G4int ncount = 0;
289  do {
290  //++ncount;
291  G4double x = G4Exp(xmin + G4UniformRand()*(xmax - xmin)) - densityCorr;
292  if(x < 0.0) { x = 0.0; }
293  gammaEnergy = sqrt(x);
294  G4double x1 = gammaEnergy/kineticEnergy;
295  v = dataSB[Z]->Value(x1, y, idx, idy);
296 
297  // correction for positrons
298  if(!isElectron) {
299  G4double e1 = kineticEnergy - cut;
300  G4double invbeta1 = (e1 + particleMass)/sqrt(e1*(e1 + 2*particleMass));
301  G4double e2 = kineticEnergy - gammaEnergy;
302  G4double invbeta2 = (e2 + particleMass)/sqrt(e2*(e2 + 2*particleMass));
303  G4double xxx = twopi*fine_structure_const*currentZ*(invbeta1 - invbeta2);
304 
305  if(xxx < expnumlim) { v = 0.0; }
306  else { v *= G4Exp(xxx); }
307  }
308 
309  if (v > 1.05*vmax && nwarn < 5) {
310  ++nwarn;
312  ed << "### G4LivermoreBremsstrahlungModel Warning: Majoranta exceeded! "
313  << v << " > " << vmax << " by " << v/vmax
314  << " Egamma(MeV)= " << gammaEnergy
315  << " Ee(MeV)= " << kineticEnergy
316  << " Z= " << Z << " " << particle->GetParticleName();
317 
318  if ( 20 == nwarn ) {
319  ed << "\n ### G4LivermoreBremsstrahlungModel Warnings stopped";
320  }
321  G4Exception("G4LivermoreBremsstrahlungModel::SampleScattering","em0044",
322  JustWarning, ed,"");
323 
324  }
325  } while (v < vmax*G4UniformRand());
326 
327  //
328  // angles of the emitted gamma. ( Z - axis along the parent particle)
329  // use general interface
330  //
331 
332  G4ThreeVector gammaDirection =
333  GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy,
334  Z, couple->GetMaterial());
335 
336  // create G4DynamicParticle object for the Gamma
337  G4DynamicParticle* gamma =
338  new G4DynamicParticle(theGamma,gammaDirection,gammaEnergy);
339  vdp->push_back(gamma);
340 
341  G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
342  - gammaEnergy*gammaDirection).unit();
343 
344  /*
345  G4cout << "### G4SBModel: v= "
346  << " Eg(MeV)= " << gammaEnergy
347  << " Ee(MeV)= " << kineticEnergy
348  << " DirE " << direction << " DirG " << gammaDirection
349  << G4endl;
350  */
351  // energy of primary
352  G4double finalE = kineticEnergy - gammaEnergy;
353 
354  // stop tracking and create new secondary instead of primary
355  if(gammaEnergy > SecondaryThreshold()) {
356  fParticleChange->ProposeTrackStatus(fStopAndKill);
357  fParticleChange->SetProposedKineticEnergy(0.0);
358  G4DynamicParticle* el =
359  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
360  direction, finalE);
361  vdp->push_back(el);
362 
363  // continue tracking
364  } else {
365  fParticleChange->SetProposedMomentumDirection(direction);
366  fParticleChange->SetProposedKineticEnergy(finalE);
367  }
368 }
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
static G4Physics2DVector * dataSB[101]
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
#define G4UniformRand()
Definition: Randomize.hh:97
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
Float_t Z
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
const G4ThreeVector & GetMomentumDirection() const
static const G4double epeaklimit
G4ParticleChangeForLoss * fParticleChange
double G4double
Definition: G4Types.hh:76
static const G4double elowlimit
G4double Value(G4double x, G4double y, size_t &lastidx, size_t &lastidy) const
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:

◆ SetBicubicInterpolationFlag()

void G4LivermoreBremsstrahlungModel::SetBicubicInterpolationFlag ( G4bool  val)
inline

Definition at line 106 of file G4LivermoreBremsstrahlungModel.hh.

Member Data Documentation

◆ dataSB

G4Physics2DVector * G4LivermoreBremsstrahlungModel::dataSB = {0}
staticprivate

Definition at line 97 of file G4LivermoreBremsstrahlungModel.hh.

◆ expnumlim

G4double G4LivermoreBremsstrahlungModel::expnumlim = -12.
staticprivate

Definition at line 99 of file G4LivermoreBremsstrahlungModel.hh.

◆ idx

size_t G4LivermoreBremsstrahlungModel::idx
private

Definition at line 101 of file G4LivermoreBremsstrahlungModel.hh.

◆ idy

size_t G4LivermoreBremsstrahlungModel::idy
private

Definition at line 102 of file G4LivermoreBremsstrahlungModel.hh.

◆ nwarn

G4int G4LivermoreBremsstrahlungModel::nwarn
private

Definition at line 100 of file G4LivermoreBremsstrahlungModel.hh.

◆ useBicubicInterpolation

G4bool G4LivermoreBremsstrahlungModel::useBicubicInterpolation
private

Definition at line 103 of file G4LivermoreBremsstrahlungModel.hh.

◆ ylimit

G4double G4LivermoreBremsstrahlungModel::ylimit = {0.0}
staticprivate

Definition at line 98 of file G4LivermoreBremsstrahlungModel.hh.


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