Geant4  10.02.p03
G4eCoulombScatteringModel Class Reference

#include <G4eCoulombScatteringModel.hh>

Inheritance diagram for G4eCoulombScatteringModel:
Collaboration diagram for G4eCoulombScatteringModel:

Public Member Functions

 G4eCoulombScatteringModel (G4bool combined=true)
 
virtual ~G4eCoulombScatteringModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle *> *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double)
 
void SetLowEnergyThreshold (G4double val)
 
void SetRecoilThreshold (G4double eth)
 
void SetFixedCut (G4double)
 
G4double GetFixedCut () const
 
- 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 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, 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 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 *> *)
 
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

void DefineMaterial (const G4MaterialCutsCouple *)
 
void SetupParticle (const G4ParticleDefinition *)
 
- 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

G4eCoulombScatteringModeloperator= (const G4eCoulombScatteringModel &right)
 
 G4eCoulombScatteringModel (const G4eCoulombScatteringModel &)
 

Private Attributes

G4IonTabletheIonTable
 
G4ParticleChangeForGamma * fParticleChange
 
G4WentzelOKandVIxSectionwokvi
 
G4NistManagerfNistManager
 
const std::vector< G4double > * pCuts
 
const G4MaterialCutsCouplecurrentCouple
 
const G4MaterialcurrentMaterial
 
G4int currentMaterialIndex
 
G4double cosThetaMin
 
G4double cosThetaMax
 
G4double recoilThreshold
 
G4double elecRatio
 
G4double mass
 
G4double fixedCut
 
const G4ParticleDefinitionparticle
 
const G4ParticleDefinitiontheProton
 
G4bool isCombined
 

Additional Inherited Members

- 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 78 of file G4eCoulombScatteringModel.hh.

Constructor & Destructor Documentation

◆ G4eCoulombScatteringModel() [1/2]

G4eCoulombScatteringModel::G4eCoulombScatteringModel ( G4bool  combined = true)

Definition at line 80 of file G4eCoulombScatteringModel.cc.

81  : G4VEmModel("eCoulombScattering"),
82  cosThetaMin(1.0),
83  cosThetaMax(-1.0),
84  isCombined(combined)
85 {
86  fParticleChange = 0;
90  currentMaterial = 0;
91  fixedCut = -1.0;
92 
93  pCuts = 0;
94 
95  recoilThreshold = 0.*keV; // by default does not work
96 
97  particle = 0;
98  currentCouple = 0;
99  wokvi = new G4WentzelOKandVIxSection(combined);
100 
103  elecRatio = 0.0;
104 }
const G4ParticleDefinition * particle
const std::vector< G4double > * pCuts
const G4MaterialCutsCouple * currentCouple
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
G4IonTable * GetIonTable() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
float proton_mass_c2
Definition: hepunit.py:275
static G4ParticleTable * GetParticleTable()
const G4ParticleDefinition * theProton
static const double keV
Definition: G4SIunits.hh:213
G4WentzelOKandVIxSection * wokvi
G4ParticleChangeForGamma * fParticleChange
Here is the call graph for this function:

◆ ~G4eCoulombScatteringModel()

G4eCoulombScatteringModel::~G4eCoulombScatteringModel ( )
virtual

Definition at line 108 of file G4eCoulombScatteringModel.cc.

109 {
110  delete wokvi;
111 }
G4WentzelOKandVIxSection * wokvi

◆ G4eCoulombScatteringModel() [2/2]

G4eCoulombScatteringModel::G4eCoulombScatteringModel ( const G4eCoulombScatteringModel )
private

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4eCoulombScatteringModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  kinEnergy,
G4double  Z,
G4double  A,
G4double  cut,
G4double  emax 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 191 of file G4eCoulombScatteringModel.cc.

196 {
197  //G4cout << "### G4eCoulombScatteringModel::ComputeCrossSectionPerAtom for "
198  //<< p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
199  G4double cross = 0.0;
200  elecRatio = 0.0;
201  if(p != particle) { SetupParticle(p); }
202 
203  // cross section is set to zero to avoid problems in sample secondary
204  if(kinEnergy <= 0.0) { return cross; }
206  G4double costmin = wokvi->SetupKinematic(kinEnergy, currentMaterial);
207  if(cosThetaMax < costmin) {
208  G4int iz = G4int(Z);
209  G4double cut = cutEnergy;
210  if(fixedCut > 0.0) { cut = fixedCut; }
211  costmin = wokvi->SetupTarget(iz, cut);
212  G4double costmax = cosThetaMax;
213  if(iz == 1 && costmax < 0.0 && particle == theProton) {
214  costmax = 0.0;
215  }
216  if(costmin > costmax) {
217  cross = wokvi->ComputeNuclearCrossSection(costmin, costmax)
218  + wokvi->ComputeElectronCrossSection(costmin, costmax);
219  }
220  /*
221  if(p->GetParticleName() == "mu+")
222  G4cout << "e(MeV)= " << kinEnergy/MeV << " cross(b)= " << cross/barn
223  << " 1-costmin= " << 1-costmin
224  << " 1-costmax= " << 1-costmax
225  << " 1-cosThetaMax= " << 1-cosThetaMax
226  << G4endl;
227  */
228  }
229  return cross;
230 }
const G4ParticleDefinition * particle
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:452
int G4int
Definition: G4Types.hh:78
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
Float_t Z
G4double iz
Definition: TRTMaterials.hh:39
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
void DefineMaterial(const G4MaterialCutsCouple *)
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
const G4ParticleDefinition * theProton
double G4double
Definition: G4Types.hh:76
G4WentzelOKandVIxSection * wokvi
void SetupParticle(const G4ParticleDefinition *)
Here is the call graph for this function:

◆ DefineMaterial()

void G4eCoulombScatteringModel::DefineMaterial ( const G4MaterialCutsCouple cup)
inlineprotected

Definition at line 165 of file G4eCoulombScatteringModel.hh.

166 {
167  if(cup != currentCouple) {
168  currentCouple = cup;
169  currentMaterial = cup->GetMaterial();
171  }
172 }
const G4Material * GetMaterial() const
const G4MaterialCutsCouple * currentCouple
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetFixedCut()

G4double G4eCoulombScatteringModel::GetFixedCut ( ) const
inline

Definition at line 203 of file G4eCoulombScatteringModel.hh.

204 {
205  return fixedCut;
206 }

◆ Initialise()

void G4eCoulombScatteringModel::Initialise ( const G4ParticleDefinition part,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 115 of file G4eCoulombScatteringModel.cc.

117 {
118  SetupParticle(part);
119  currentCouple = 0;
120 
121  if(isCombined) {
122  cosThetaMin = 1.0;
124  if(tet >= pi) { cosThetaMin = -1.0; }
125  else if(tet > 0.0) { cosThetaMin = cos(tet); }
126  }
127 
128  wokvi->Initialise(part, cosThetaMin);
129  /*
130  G4cout << "G4eCoulombScatteringModel: " << particle->GetParticleName()
131  << " 1-cos(ThetaLimit)= " << 1 - cosThetaMin
132  << " cos(thetaMax)= " << cosThetaMax
133  << G4endl;
134  */
135  pCuts = &cuts;
136  //G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
137  /*
138  G4cout << "!!! G4eCoulombScatteringModel::Initialise for "
139  << part->GetParticleName() << " cos(TetMin)= " << cosThetaMin
140  << " cos(TetMax)= " << cosThetaMax <<G4endl;
141  G4cout << "cut= " << (*pCuts)[0] << " cut1= " << (*pCuts)[1] << G4endl;
142  */
143  if(!fParticleChange) {
145  }
146  if(IsMaster() && mass < GeV && part->GetParticleName() != "GenericIon") {
147  InitialiseElementSelectors(part,cuts);
148  }
149 }
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
const std::vector< G4double > * pCuts
const G4MaterialCutsCouple * currentCouple
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:662
static G4double tet[DIM]
static const double pi
Definition: G4SIunits.hh:74
double G4double
Definition: G4Types.hh:76
G4WentzelOKandVIxSection * wokvi
void SetupParticle(const G4ParticleDefinition *)
G4ParticleChangeForGamma * fParticleChange
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
Here is the call graph for this function:

◆ InitialiseLocal()

void G4eCoulombScatteringModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 153 of file G4eCoulombScatteringModel.cc.

155 {
157 }
void SetElementSelectors(std::vector< G4EmElementSelector *> *)
Definition: G4VEmModel.hh:810
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:802
Here is the call graph for this function:

◆ MinPrimaryEnergy()

G4double G4eCoulombScatteringModel::MinPrimaryEnergy ( const G4Material material,
const G4ParticleDefinition part,
G4double   
)
virtual

Reimplemented from G4VEmModel.

Definition at line 162 of file G4eCoulombScatteringModel.cc.

165 {
166  SetupParticle(part);
167 
168  // define cut using cuts for proton
169  G4double cut =
170  std::max(recoilThreshold, (*pCuts)[CurrentCouple()->GetIndex()]);
171 
172  // find out lightest element
173  const G4ElementVector* theElementVector = material->GetElementVector();
174  G4int nelm = material->GetNumberOfElements();
175 
176  // select lightest element
177  G4int Z = 300;
178  for (G4int j=0; j<nelm; ++j) {
179  G4int iz = G4lrint((*theElementVector)[j]->GetZ());
180  if(iz < Z) { Z = iz; }
181  }
183  G4double targetMass = G4NucleiProperties::GetNuclearMass(A, Z);
184  G4double t = std::max(cut, 0.5*(cut + sqrt(2*cut*targetMass)));
185 
186  return t;
187 }
static G4double GetNuclearMass(const G4double A, const G4double Z)
std::vector< G4Element * > G4ElementVector
const std::vector< G4double > * pCuts
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:452
int G4int
Definition: G4Types.hh:78
double A(double temperature)
Float_t Z
G4double iz
Definition: TRTMaterials.hh:39
int G4lrint(double ad)
Definition: templates.hh:163
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
double G4double
Definition: G4Types.hh:76
void SetupParticle(const G4ParticleDefinition *)
G4double GetAtomicMassAmu(const G4String &symb) const
Here is the call graph for this function:

◆ operator=()

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

◆ SampleSecondaries()

void G4eCoulombScatteringModel::SampleSecondaries ( std::vector< G4DynamicParticle *> *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 234 of file G4eCoulombScatteringModel.cc.

240 {
241  G4double kinEnergy = dp->GetKineticEnergy();
243  DefineMaterial(couple);
244  /*
245  G4cout << "G4eCoulombScatteringModel::SampleSecondaries e(MeV)= "
246  << kinEnergy << " " << particle->GetParticleName()
247  << " cut= " << cutEnergy<< G4endl;
248  */
249  // Choose nucleus
250  G4double cut = cutEnergy;
251  if(fixedCut > 0.0) { cut = fixedCut; }
252 
253  wokvi->SetupKinematic(kinEnergy, currentMaterial);
254 
255  const G4Element* currentElement =
256  SelectRandomAtom(couple,particle,kinEnergy,cut,kinEnergy);
257 
258  G4double Z = currentElement->GetZ();
259  G4int iz = G4int(Z);
260 
261  G4double costmin = wokvi->SetupTarget(iz, cut);
262  G4double costmax = cosThetaMax;
263  if(iz == 1 && costmax < 0.0 && particle == theProton) {
264  costmax = 0.0;
265  }
266 
267  if(costmin > costmax) {
268  G4double cross = wokvi->ComputeNuclearCrossSection(costmin, costmax);
269  G4double ecross = wokvi->ComputeElectronCrossSection(costmin, costmax);
270  G4double ratio = ecross/(cross + ecross);
271 
272  G4int ia = SelectIsotopeNumber(currentElement);
273  G4double targetMass = G4NucleiProperties::GetNuclearMass(ia, iz);
274  wokvi->SetTargetMass(targetMass);
275 
276  G4ThreeVector newDirection =
277  wokvi->SampleSingleScattering(costmin, costmax, ratio);
278  G4double cost = newDirection.z();
279  /*
280  G4cout << "SampleSec: e(MeV)= " << kinEnergy/MeV
281  << " 1-costmin= " << 1-costmin
282  << " 1-costmax= " << 1-costmax
283  << " 1-cost= " << 1-cost
284  << " ratio= " << ratio
285  << G4endl;
286  */
287  G4ThreeVector direction = dp->GetMomentumDirection();
288  newDirection.rotateUz(direction);
289 
290  fParticleChange->ProposeMomentumDirection(newDirection);
291 
292  // recoil sampling assuming a small recoil
293  // and first order correction to primary 4-momentum
294  G4double mom2 = wokvi->GetMomentumSquare();
295  G4double trec = mom2*(1.0 - cost)
296  /(targetMass + (mass + kinEnergy)*(1.0 - cost));
297 
298  // the check likely not needed
299  if(trec > kinEnergy) { trec = kinEnergy; }
300  G4double finalT = kinEnergy - trec;
301  G4double edep = 0.0;
302  /*
303  G4cout<<"G4eCoulombScatteringModel: finalT= "<<finalT<<" Trec= "
304  <<trec << " Z= " << iz << " A= " << ia
305  << " tcut(keV)= " << (*pCuts)[currentMaterialIndex]/keV << G4endl;
306  */
307  G4double tcut = recoilThreshold;
308  if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
309 
310  if(trec > tcut) {
312  G4ThreeVector dir = (direction*sqrt(mom2) -
313  newDirection*sqrt(finalT*(2*mass + finalT))).unit();
314  G4DynamicParticle* newdp = new G4DynamicParticle(ion, dir, trec);
315  fvect->push_back(newdp);
316  } else {
317  edep = trec;
318  fParticleChange->ProposeNonIonizingEnergyDeposit(edep);
319  }
320 
321  // finelize primary energy and energy balance
322  // this threshold may be applied only because for low-enegry
323  // e+e- msc model is applied
324  if(finalT < 0.0) {
325  edep += finalT;
326  finalT = 0.0;
327  if(edep < 0.0) { edep = 0.0; }
328  }
329  fParticleChange->SetProposedKineticEnergy(finalT);
330  fParticleChange->ProposeLocalEnergyDeposit(edep);
331  }
332 }
const G4ParticleDefinition * particle
static G4double GetNuclearMass(const G4double A, const G4double Z)
const std::vector< G4double > * pCuts
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:491
TDirectory * dir
int G4int
Definition: G4Types.hh:78
Double_t edep
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
G4double GetKineticEnergy() const
Float_t Z
G4double iz
Definition: TRTMaterials.hh:39
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
void DefineMaterial(const G4MaterialCutsCouple *)
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:585
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
const G4ThreeVector & GetMomentumDirection() const
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
double z() const
const G4ParticleDefinition * theProton
G4ParticleDefinition * GetDefinition() const
double G4double
Definition: G4Types.hh:76
G4WentzelOKandVIxSection * wokvi
void SetupParticle(const G4ParticleDefinition *)
G4double GetZ() const
Definition: G4Element.hh:131
G4ParticleChangeForGamma * fParticleChange
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)
Here is the call graph for this function:

◆ SetFixedCut()

void G4eCoulombScatteringModel::SetFixedCut ( G4double  val)
inline

Definition at line 196 of file G4eCoulombScatteringModel.hh.

197 {
198  fixedCut = val;
199 }

◆ SetLowEnergyThreshold()

void G4eCoulombScatteringModel::SetLowEnergyThreshold ( G4double  val)
inline

◆ SetRecoilThreshold()

void G4eCoulombScatteringModel::SetRecoilThreshold ( G4double  eth)
inline

Definition at line 189 of file G4eCoulombScatteringModel.hh.

190 {
191  recoilThreshold = eth;
192 }

◆ SetupParticle()

void G4eCoulombScatteringModel::SetupParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 177 of file G4eCoulombScatteringModel.hh.

178 {
179  // Initialise mass and charge
180  if(p != particle) {
181  particle = p;
182  mass = particle->GetPDGMass();
183  wokvi->SetupParticle(p);
184  }
185 }
void SetupParticle(const G4ParticleDefinition *)
const G4ParticleDefinition * particle
G4WentzelOKandVIxSection * wokvi
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ cosThetaMax

G4double G4eCoulombScatteringModel::cosThetaMax
private

Definition at line 148 of file G4eCoulombScatteringModel.hh.

◆ cosThetaMin

G4double G4eCoulombScatteringModel::cosThetaMin
private

Definition at line 147 of file G4eCoulombScatteringModel.hh.

◆ currentCouple

const G4MaterialCutsCouple* G4eCoulombScatteringModel::currentCouple
private

Definition at line 143 of file G4eCoulombScatteringModel.hh.

◆ currentMaterial

const G4Material* G4eCoulombScatteringModel::currentMaterial
private

Definition at line 144 of file G4eCoulombScatteringModel.hh.

◆ currentMaterialIndex

G4int G4eCoulombScatteringModel::currentMaterialIndex
private

Definition at line 145 of file G4eCoulombScatteringModel.hh.

◆ elecRatio

G4double G4eCoulombScatteringModel::elecRatio
private

Definition at line 150 of file G4eCoulombScatteringModel.hh.

◆ fixedCut

G4double G4eCoulombScatteringModel::fixedCut
private

Definition at line 153 of file G4eCoulombScatteringModel.hh.

◆ fNistManager

G4NistManager* G4eCoulombScatteringModel::fNistManager
private

Definition at line 139 of file G4eCoulombScatteringModel.hh.

◆ fParticleChange

G4ParticleChangeForGamma* G4eCoulombScatteringModel::fParticleChange
private

Definition at line 137 of file G4eCoulombScatteringModel.hh.

◆ isCombined

G4bool G4eCoulombScatteringModel::isCombined
private

Definition at line 159 of file G4eCoulombScatteringModel.hh.

◆ mass

G4double G4eCoulombScatteringModel::mass
private

Definition at line 151 of file G4eCoulombScatteringModel.hh.

◆ particle

const G4ParticleDefinition* G4eCoulombScatteringModel::particle
private

Definition at line 156 of file G4eCoulombScatteringModel.hh.

◆ pCuts

const std::vector<G4double>* G4eCoulombScatteringModel::pCuts
private

Definition at line 141 of file G4eCoulombScatteringModel.hh.

◆ recoilThreshold

G4double G4eCoulombScatteringModel::recoilThreshold
private

Definition at line 149 of file G4eCoulombScatteringModel.hh.

◆ theIonTable

G4IonTable* G4eCoulombScatteringModel::theIonTable
private

Definition at line 136 of file G4eCoulombScatteringModel.hh.

◆ theProton

const G4ParticleDefinition* G4eCoulombScatteringModel::theProton
private

Definition at line 157 of file G4eCoulombScatteringModel.hh.

◆ wokvi

G4WentzelOKandVIxSection* G4eCoulombScatteringModel::wokvi
private

Definition at line 138 of file G4eCoulombScatteringModel.hh.


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