Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes
G4eCoulombScatteringModel Class Reference

#include <G4eCoulombScatteringModel.hh>

Inheritance diagram for G4eCoulombScatteringModel:
G4VEmModel

Public Member Functions

 G4eCoulombScatteringModel (const G4String &nam="eCoulombScattering")
 
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 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 *)
 
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 Member Functions

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

Protected Attributes

G4IonTabletheIonTable
 
G4ParticleChangeForGammafParticleChange
 
G4WentzelOKandVIxSectionwokvi
 
G4NistManagerfNistManager
 
const std::vector< G4double > * pCuts
 
const G4MaterialCutsCouplecurrentCouple
 
const G4MaterialcurrentMaterial
 
G4int currentMaterialIndex
 
G4double cosThetaMin
 
G4double cosThetaMax
 
G4double cosTetMinNuc
 
G4double cosTetMaxNuc
 
G4double recoilThreshold
 
G4double elecRatio
 
G4double mass
 
G4double fixedCut
 
const G4ParticleDefinitionparticle
 
const G4ParticleDefinitiontheProton
 
G4double lowEnergyThreshold
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Detailed Description

Definition at line 78 of file G4eCoulombScatteringModel.hh.

Constructor & Destructor Documentation

G4eCoulombScatteringModel::G4eCoulombScatteringModel ( const G4String nam = "eCoulombScattering")

Definition at line 80 of file G4eCoulombScatteringModel.cc.

References cosTetMaxNuc, cosTetMinNuc, currentCouple, currentMaterial, currentMaterialIndex, elecRatio, fixedCut, fNistManager, fParticleChange, G4ParticleTable::GetIonTable(), G4ParticleTable::GetParticleTable(), G4NistManager::Instance(), python.hepunit::keV, lowEnergyThreshold, mass, particle, pCuts, G4Proton::Proton(), python.hepunit::proton_mass_c2, recoilThreshold, theIonTable, theProton, and wokvi.

81  : G4VEmModel(nam),
82  cosThetaMin(1.0),
83  cosThetaMax(-1.0),
84  isInitialised(false)
85 {
86  fParticleChange = 0;
90  currentMaterial = 0;
91  fixedCut = -1.0;
92 
93  pCuts = 0;
94 
95  lowEnergyThreshold = 1*keV; // particle will be killed for lower energy
96  recoilThreshold = 0.*keV; // by default does not work
97 
98  particle = 0;
99  currentCouple = 0;
101 
103 
104  cosTetMinNuc = 1.0;
105  cosTetMaxNuc = -1.0;
106  elecRatio = 0.0;
108 }
const G4ParticleDefinition * particle
const std::vector< G4double > * pCuts
const G4MaterialCutsCouple * currentCouple
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
G4IonTable * GetIonTable() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
float proton_mass_c2
Definition: hepunit.py:275
static G4ParticleTable * GetParticleTable()
const G4ParticleDefinition * theProton
G4WentzelOKandVIxSection * wokvi
G4ParticleChangeForGamma * fParticleChange
G4eCoulombScatteringModel::~G4eCoulombScatteringModel ( )
virtual

Definition at line 112 of file G4eCoulombScatteringModel.cc.

References wokvi.

113 {
114  delete wokvi;
115 }
G4WentzelOKandVIxSection * wokvi

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 187 of file G4eCoulombScatteringModel.cc.

References G4WentzelOKandVIxSection::ComputeElectronCrossSection(), G4WentzelOKandVIxSection::ComputeNuclearCrossSection(), cosTetMaxNuc, cosTetMinNuc, cosThetaMax, G4VEmModel::CurrentCouple(), currentMaterial, DefineMaterial(), elecRatio, fixedCut, iz, particle, G4WentzelOKandVIxSection::SetupKinematic(), SetupParticle(), G4WentzelOKandVIxSection::SetupTarget(), theProton, and wokvi.

Referenced by SampleSecondaries().

192 {
193  //G4cout << "### G4eCoulombScatteringModel::ComputeCrossSectionPerAtom for "
194  //<< p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
195  G4double cross = 0.0;
196  if(p != particle) { SetupParticle(p); }
197 
198  // cross section is set to zero to avoid problems in sample secondary
199  if(kinEnergy <= 0.0) { return cross; }
202  if(cosThetaMax < cosTetMinNuc) {
203  G4int iz = G4int(Z);
204  G4double cut = cutEnergy;
205  if(fixedCut > 0.0) { cut = fixedCut; }
206  cosTetMinNuc = wokvi->SetupTarget(iz, cut);
208  if(iz == 1 && cosTetMaxNuc < 0.0 && particle == theProton) {
209  cosTetMaxNuc = 0.0;
210  }
213  cross += elecRatio;
214  if(cross > 0.0) { elecRatio /= cross; }
215  }
216  /*
217  if(p->GetParticleName() == "e-")
218  G4cout << "e(MeV)= " << kinEnergy/MeV << " cross(b)= " << cross/barn
219  << " 1-cosTetMinNuc= " << 1-cosTetMinNuc
220  << " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc
221  << G4endl;
222  */
223  return cross;
224 }
const G4ParticleDefinition * particle
int G4int
Definition: G4Types.hh:78
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:426
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 *)
void G4eCoulombScatteringModel::DefineMaterial ( const G4MaterialCutsCouple cup)
inlineprotected
G4double G4eCoulombScatteringModel::GetFixedCut ( ) const
inline

Definition at line 216 of file G4eCoulombScatteringModel.hh.

References fixedCut.

217 {
218  return fixedCut;
219 }
void G4eCoulombScatteringModel::Initialise ( const G4ParticleDefinition part,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 119 of file G4eCoulombScatteringModel.cc.

References cosThetaMin, currentCouple, fParticleChange, G4VEmModel::GetParticleChangeForGamma(), G4WentzelOKandVIxSection::Initialise(), G4VEmModel::InitialiseElementSelectors(), G4VEmModel::IsMaster(), pCuts, G4VEmModel::PolarAngleLimit(), SetupParticle(), and wokvi.

121 {
122  SetupParticle(part);
123  currentCouple = 0;
124  cosThetaMin = cos(PolarAngleLimit());
125  wokvi->Initialise(part, cosThetaMin);
126  /*
127  G4cout << "G4eCoulombScatteringModel: " << particle->GetParticleName()
128  << " 1-cos(ThetaLimit)= " << 1 - cosThetaMin
129  << " cos(thetaMax)= " << cosThetaMax
130  << G4endl;
131  */
132  pCuts = &cuts;
133  // G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
134  /*
135  G4cout << "!!! G4eCoulombScatteringModel::Initialise for "
136  << part->GetParticleName() << " cos(TetMin)= " << cosThetaMin
137  << " cos(TetMax)= " << cosThetaMax <<G4endl;
138  G4cout << "cut= " << pCuts[0] << " cut1= " << pCuts[1] << G4endl;
139  */
140  if(!isInitialised) {
141  isInitialised = true;
143  }
144  if(IsMaster() && mass < GeV && part->GetParticleName() != "GenericIon") {
145  InitialiseElementSelectors(part,cuts);
146  }
147 }
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
const std::vector< G4double > * pCuts
const G4MaterialCutsCouple * currentCouple
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:620
G4WentzelOKandVIxSection * wokvi
void SetupParticle(const G4ParticleDefinition *)
G4ParticleChangeForGamma * fParticleChange
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
void G4eCoulombScatteringModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 151 of file G4eCoulombScatteringModel.cc.

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

153 {
155 }
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:768
G4double G4eCoulombScatteringModel::MinPrimaryEnergy ( const G4Material material,
const G4ParticleDefinition part,
G4double   
)
virtual

Reimplemented from G4VEmModel.

Definition at line 160 of file G4eCoulombScatteringModel.cc.

References G4VEmModel::CurrentCouple(), fNistManager, G4lrint(), G4NistManager::GetAtomicMassAmu(), G4Material::GetElementVector(), G4NucleiProperties::GetNuclearMass(), G4Material::GetNumberOfElements(), iz, lowEnergyThreshold, G4INCL::Math::max(), pCuts, recoilThreshold, and SetupParticle().

163 {
164  SetupParticle(part);
165 
166  // define cut using cuts for proton
167  G4double cut =
168  std::max(recoilThreshold, (*pCuts)[CurrentCouple()->GetIndex()]);
169 
170  // find out lightest element
171  const G4ElementVector* theElementVector = material->GetElementVector();
172  G4int nelm = material->GetNumberOfElements();
173  G4int Z = 300;
174  for (G4int j=0; j<nelm; ++j) {
175  G4int iz = (G4int)(*theElementVector)[j]->GetZ();
176  if(iz < Z) { Z = iz; }
177  }
179  G4double targetMass = G4NucleiProperties::GetNuclearMass(A, Z);
180  G4double t = std::max(cut, 0.5*(cut + sqrt(2*cut*targetMass)));
181 
182  return std::max(lowEnergyThreshold, t);
183 }
static G4double GetNuclearMass(const G4double A, const G4double Z)
std::vector< G4Element * > G4ElementVector
const std::vector< G4double > * pCuts
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:426
G4double iz
Definition: TRTMaterials.hh:39
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double GetAtomicMassAmu(const G4String &symb) const
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
double G4double
Definition: G4Types.hh:76
void SetupParticle(const G4ParticleDefinition *)
void G4eCoulombScatteringModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 228 of file G4eCoulombScatteringModel.cc.

References ComputeCrossSectionPerAtom(), cosTetMinNuc, cosThetaMax, currentMaterialIndex, DefineMaterial(), elecRatio, fixedCut, fParticleChange, G4DynamicParticle::GetDefinition(), G4IonTable::GetIon(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4WentzelOKandVIxSection::GetMomentumSquare(), G4NucleiProperties::GetNuclearMass(), G4Element::GetZ(), iz, lowEnergyThreshold, mass, G4INCL::Math::max(), particle, pCuts, G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VParticleChange::ProposeNonIonizingEnergyDeposit(), recoilThreshold, CLHEP::Hep3Vector::rotateUz(), G4WentzelOKandVIxSection::SampleSingleScattering(), G4VEmModel::SelectIsotopeNumber(), G4VEmModel::SelectRandomAtom(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), G4WentzelOKandVIxSection::SetTargetMass(), SetupParticle(), theIonTable, wokvi, and CLHEP::Hep3Vector::z().

234 {
235  G4double kinEnergy = dp->GetKineticEnergy();
236 
237  // absorb particle below low-energy limit to avoid situation
238  // when a particle has no energy loss
239  if(kinEnergy < lowEnergyThreshold) {
243  return;
244  }
246  DefineMaterial(couple);
247  /*
248  G4cout << "G4eCoulombScatteringModel::SampleSecondaries e(MeV)= "
249  << kinEnergy << " " << particle->GetParticleName()
250  << " cut= " << cutEnergy<< G4endl;
251  */
252  // Choose nucleus
253  G4double cut = cutEnergy;
254  if(fixedCut > 0.0) { cut = fixedCut; }
255 
256  const G4Element* currentElement =
257  SelectRandomAtom(couple,particle,kinEnergy,cut,kinEnergy);
258 
259  G4double Z = currentElement->GetZ();
260 
261  if(ComputeCrossSectionPerAtom(particle,kinEnergy, Z,
262  kinEnergy, cut, kinEnergy) == 0.0)
263  { return; }
264 
265  G4int iz = G4int(Z);
266  G4int ia = SelectIsotopeNumber(currentElement);
267  G4double targetMass = G4NucleiProperties::GetNuclearMass(ia, iz);
268  wokvi->SetTargetMass(targetMass);
269 
270  G4ThreeVector newDirection =
272  G4double cost = newDirection.z();
273 
274  G4ThreeVector direction = dp->GetMomentumDirection();
275  newDirection.rotateUz(direction);
276 
278 
279  // recoil sampling assuming a small recoil
280  // and first order correction to primary 4-momentum
281  G4double mom2 = wokvi->GetMomentumSquare();
282  G4double trec = mom2*(1.0 - cost)
283  /(targetMass + (mass + kinEnergy)*(1.0 - cost));
284 
285  // the check likely not needed
286  if(trec > kinEnergy) { trec = kinEnergy; }
287  G4double finalT = kinEnergy - trec;
288  G4double edep = 0.0;
289  //G4cout<<"G4eCoulombScatteringModel: finalT= "<<finalT<<" Trec= "
290  // <<trec << " Z= " << iz << " A= " << ia<<G4endl;
291 
292  G4double tcut = recoilThreshold;
293  if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
294 
295  if(trec > tcut) {
297  G4ThreeVector dir = (direction*sqrt(mom2) -
298  newDirection*sqrt(finalT*(2*mass + finalT))).unit();
299  G4DynamicParticle* newdp = new G4DynamicParticle(ion, dir, trec);
300  fvect->push_back(newdp);
301  } else {
302  edep = trec;
304  }
305 
306  // finelize primary energy and energy balance
307  // this threshold may be applied only because for low-enegry
308  // e+e- msc model is applied
309  if(finalT <= lowEnergyThreshold) {
310  edep += finalT;
311  finalT = 0.0;
312  }
315 }
const G4ParticleDefinition * particle
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double GetKineticEnergy() const
G4ThreeVector SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)
const std::vector< G4double > * pCuts
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:449
G4double GetZ() const
Definition: G4Element.hh:131
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
double z() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
const G4ThreeVector & GetMomentumDirection() const
G4double iz
Definition: TRTMaterials.hh:39
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
void DefineMaterial(const G4MaterialCutsCouple *)
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:548
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void SetProposedKineticEnergy(G4double proposedKinEnergy)
double G4double
Definition: G4Types.hh:76
G4WentzelOKandVIxSection * wokvi
void SetupParticle(const G4ParticleDefinition *)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
G4ParticleChangeForGamma * fParticleChange
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:510
void G4eCoulombScatteringModel::SetFixedCut ( G4double  val)
inline

Definition at line 209 of file G4eCoulombScatteringModel.hh.

References fixedCut.

210 {
211  fixedCut = val;
212 }
void G4eCoulombScatteringModel::SetLowEnergyThreshold ( G4double  val)
inline

Definition at line 195 of file G4eCoulombScatteringModel.hh.

References lowEnergyThreshold.

Referenced by PhysListEmStandardWVI::ConstructProcess().

196 {
197  lowEnergyThreshold = val;
198 }
void G4eCoulombScatteringModel::SetRecoilThreshold ( G4double  eth)
inline

Definition at line 202 of file G4eCoulombScatteringModel.hh.

References recoilThreshold.

203 {
204  recoilThreshold = eth;
205 }
void G4eCoulombScatteringModel::SetupParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 183 of file G4eCoulombScatteringModel.hh.

References G4ParticleDefinition::GetPDGMass(), mass, particle, G4WentzelOKandVIxSection::SetupParticle(), and wokvi.

Referenced by ComputeCrossSectionPerAtom(), Initialise(), MinPrimaryEnergy(), and SampleSecondaries().

184 {
185  // Initialise mass and charge
186  if(p != particle) {
187  particle = p;
188  mass = particle->GetPDGMass();
189  wokvi->SetupParticle(p);
190  }
191 }
void SetupParticle(const G4ParticleDefinition *)
const G4ParticleDefinition * particle
const char * p
Definition: xmltok.h:285
G4double GetPDGMass() const
G4WentzelOKandVIxSection * wokvi

Field Documentation

G4double G4eCoulombScatteringModel::cosTetMaxNuc
protected
G4double G4eCoulombScatteringModel::cosTetMinNuc
protected
G4double G4eCoulombScatteringModel::cosThetaMax
protected

Definition at line 148 of file G4eCoulombScatteringModel.hh.

Referenced by ComputeCrossSectionPerAtom(), and SampleSecondaries().

G4double G4eCoulombScatteringModel::cosThetaMin
protected

Definition at line 147 of file G4eCoulombScatteringModel.hh.

Referenced by Initialise().

const G4MaterialCutsCouple* G4eCoulombScatteringModel::currentCouple
protected
const G4Material* G4eCoulombScatteringModel::currentMaterial
protected
G4int G4eCoulombScatteringModel::currentMaterialIndex
protected
G4double G4eCoulombScatteringModel::elecRatio
protected
G4double G4eCoulombScatteringModel::fixedCut
protected
G4NistManager* G4eCoulombScatteringModel::fNistManager
protected

Definition at line 139 of file G4eCoulombScatteringModel.hh.

Referenced by G4eCoulombScatteringModel(), and MinPrimaryEnergy().

G4ParticleChangeForGamma* G4eCoulombScatteringModel::fParticleChange
protected
G4double G4eCoulombScatteringModel::lowEnergyThreshold
protected
G4double G4eCoulombScatteringModel::mass
protected
const G4ParticleDefinition* G4eCoulombScatteringModel::particle
protected
const std::vector<G4double>* G4eCoulombScatteringModel::pCuts
protected
G4double G4eCoulombScatteringModel::recoilThreshold
protected
G4IonTable* G4eCoulombScatteringModel::theIonTable
protected

Definition at line 136 of file G4eCoulombScatteringModel.hh.

Referenced by G4eCoulombScatteringModel(), and SampleSecondaries().

const G4ParticleDefinition* G4eCoulombScatteringModel::theProton
protected
G4WentzelOKandVIxSection* G4eCoulombScatteringModel::wokvi
protected

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