Geant4  10.02.p03
G4hCoulombScatteringModel Class Reference

#include <G4hCoulombScatteringModel.hh>

Inheritance diagram for G4hCoulombScatteringModel:
Collaboration diagram for G4hCoulombScatteringModel:

Public Member Functions

 G4hCoulombScatteringModel (G4bool combined=true)
 
virtual ~G4hCoulombScatteringModel ()
 
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

G4hCoulombScatteringModeloperator= (const G4hCoulombScatteringModel &right)
 
 G4hCoulombScatteringModel (const G4hCoulombScatteringModel &)
 

Private Attributes

G4IonTabletheIonTable
 
G4ParticleChangeForGamma * fParticleChange
 
G4WentzelVIRelXSectionwokvi
 
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 67 of file G4hCoulombScatteringModel.hh.

Constructor & Destructor Documentation

◆ G4hCoulombScatteringModel() [1/2]

G4hCoulombScatteringModel::G4hCoulombScatteringModel ( G4bool  combined = true)

Definition at line 70 of file G4hCoulombScatteringModel.cc.

71  : G4VEmModel("hCoulombScattering"),
72  cosThetaMin(1.0),
73  cosThetaMax(-1.0),
74  isCombined(combined)
75 {
76  fParticleChange = 0;
80  currentMaterial = 0;
81  fixedCut = -1.0;
82 
83  pCuts = 0;
84 
85  recoilThreshold = 0.*keV; // by default does not work
86 
87  particle = 0;
88  currentCouple = 0;
89  wokvi = new G4WentzelVIRelXSection(combined);
90 
93  elecRatio = 0.0;
94 }
const G4ParticleDefinition * theProton
const std::vector< G4double > * pCuts
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
const G4MaterialCutsCouple * currentCouple
const G4ParticleDefinition * particle
G4IonTable * GetIonTable() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
float proton_mass_c2
Definition: hepunit.py:275
G4ParticleChangeForGamma * fParticleChange
static G4ParticleTable * GetParticleTable()
static const double keV
Definition: G4SIunits.hh:213
Here is the call graph for this function:

◆ ~G4hCoulombScatteringModel()

G4hCoulombScatteringModel::~G4hCoulombScatteringModel ( )
virtual

Definition at line 98 of file G4hCoulombScatteringModel.cc.

99 {
100  delete wokvi;
101 }

◆ G4hCoulombScatteringModel() [2/2]

G4hCoulombScatteringModel::G4hCoulombScatteringModel ( const G4hCoulombScatteringModel )
private

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 181 of file G4hCoulombScatteringModel.cc.

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

◆ DefineMaterial()

void G4hCoulombScatteringModel::DefineMaterial ( const G4MaterialCutsCouple cup)
inlineprotected

Definition at line 152 of file G4hCoulombScatteringModel.hh.

153 {
154  if(cup != currentCouple) {
155  currentCouple = cup;
156  currentMaterial = cup->GetMaterial();
158  }
159 }
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 G4hCoulombScatteringModel::GetFixedCut ( ) const
inline

Definition at line 190 of file G4hCoulombScatteringModel.hh.

191 {
192  return fixedCut;
193 }

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 105 of file G4hCoulombScatteringModel.cc.

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

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 143 of file G4hCoulombScatteringModel.cc.

145 {
147 }
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 G4hCoulombScatteringModel::MinPrimaryEnergy ( const G4Material material,
const G4ParticleDefinition part,
G4double   
)
virtual

Reimplemented from G4VEmModel.

Definition at line 152 of file G4hCoulombScatteringModel.cc.

155 {
156  SetupParticle(part);
157 
158  // define cut using cuts for proton
159  G4double cut =
160  std::max(recoilThreshold, (*pCuts)[CurrentCouple()->GetIndex()]);
161 
162  // find out lightest element
163  const G4ElementVector* theElementVector = material->GetElementVector();
164  G4int nelm = material->GetNumberOfElements();
165 
166  // select lightest element
167  G4int Z = 300;
168  for (G4int j=0; j<nelm; ++j) {
169  G4int iz = G4lrint((*theElementVector)[j]->GetZ());
170  if(iz < Z) { Z = iz; }
171  }
173  G4double targetMass = G4NucleiProperties::GetNuclearMass(A, Z);
174  G4double t = std::max(cut, 0.5*(cut + sqrt(2*cut*targetMass)));
175 
176  return t;
177 }
static G4double GetNuclearMass(const G4double A, const G4double Z)
std::vector< G4Element * > G4ElementVector
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:452
const std::vector< G4double > * pCuts
void SetupParticle(const G4ParticleDefinition *)
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
G4double GetAtomicMassAmu(const G4String &symb) const
Here is the call graph for this function:

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 230 of file G4hCoulombScatteringModel.cc.

236 {
237  G4double kinEnergy = dp->GetKineticEnergy();
239  DefineMaterial(couple);
240 
241  // Choose nucleus
242  const G4Element* elm = SelectRandomAtom(couple,particle,
243  kinEnergy,cutEnergy,kinEnergy);
244 
245  G4double Z = elm->GetZ();
246 
247  G4int iz = G4int(Z);
248  G4int ia = SelectIsotopeNumber(elm);
250 
251  wokvi->SetupKinematic(kinEnergy, currentMaterial, cutEnergy, mass2);
252  G4double costmin = wokvi->SetupTarget(iz, cutEnergy);
253  G4double costmax = cosThetaMax;
254  if(iz == 1 && costmax < 0.0 && particle == theProton) {
255  costmax = 0.0;
256  }
257 
258  if(costmin <= costmax) { return; }
259  G4double cross = wokvi->ComputeNuclearCrossSection(costmin, costmax);
260  G4double ecross = wokvi->ComputeElectronCrossSection(costmin, costmax);
261  G4double ratio = ecross/(cross + ecross);
262 
263  G4ThreeVector newDirection =
264  wokvi->SampleSingleScattering(costmin, costmax, ratio);
265 
266  // kinematics in the Lab system
267  G4double ptot = dp->GetTotalMomentum();
268  G4double e1 = dp->GetTotalEnergy();
269 
270  // Lab. system kinematics along projectile direction
271  G4LorentzVector v0 = G4LorentzVector(0, 0, ptot, e1 + mass2);
272  G4double bet = ptot/v0.e();
273  G4double gam = 1.0/sqrt((1.0 - bet)*(1.0 + bet));
274 
275  // CM projectile
276  G4double momCM = gam*(ptot - bet*e1);
277  G4double eCM = gam*(e1 - bet*ptot);
278  // energy & momentum after scattering of incident particle
279  G4double pxCM = momCM*newDirection.x();
280  G4double pyCM = momCM*newDirection.y();
281  G4double pzCM = momCM*newDirection.z();
282 
283  // CM--->Lab
284  G4LorentzVector v1(pxCM , pyCM, gam*(pzCM + bet*eCM), gam*(eCM + bet*pzCM));
285 
287  newDirection = v1.vect().unit();
288  newDirection.rotateUz(dir);
289 
290  fParticleChange->ProposeMomentumDirection(newDirection);
291 
292  // recoil
293  v0 -= v1;
294  G4double trec = v0.e() - mass2;
295  G4double edep = 0.0;
296 
297  G4double tcut = recoilThreshold;
298  if(pCuts) {
299  tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]);
300  //G4cout<<" tcut eV "<<tcut/eV<<endl;
301  }
302 
303  if(trec > tcut) {
305  newDirection = v0.vect().unit();
306  newDirection.rotateUz(dir);
307  G4DynamicParticle* newdp = new G4DynamicParticle(ion, newDirection, trec);
308  fvect->push_back(newdp);
309  } else if(trec > 0.0) {
310  edep = trec;
311  fParticleChange->ProposeNonIonizingEnergyDeposit(edep);
312  }
313 
314  // finelize primary energy and energy balance
315  G4double finalT = v1.e() - mass;
316  if(finalT < 0.0) {
317  edep += finalT;
318  finalT = 0.0;
319  }
320  edep = std::max(edep, 0.0);
321  fParticleChange->SetProposedKineticEnergy(finalT);
322  fParticleChange->ProposeLocalEnergyDeposit(edep);
323 }
const G4ParticleDefinition * theProton
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
static G4double GetNuclearMass(const G4double A, const G4double Z)
void DefineMaterial(const G4MaterialCutsCouple *)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:491
const std::vector< G4double > * pCuts
G4double GetTotalMomentum() const
void SetupParticle(const G4ParticleDefinition *)
TDirectory * dir
int G4int
Definition: G4Types.hh:78
Double_t edep
G4double GetTotalEnergy() const
Hep3Vector vect() const
G4double GetKineticEnergy() const
const G4ParticleDefinition * particle
Float_t Z
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat, G4double cut, G4double tmass)
G4double SetupTarget(G4int Z, G4double cut)
G4double iz
Definition: TRTMaterials.hh:39
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:585
G4ParticleChangeForGamma * fParticleChange
static const G4double e1
double x() const
double y() const
const G4ThreeVector & GetMomentumDirection() const
double z() const
G4ParticleDefinition * GetDefinition() const
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
double G4double
Definition: G4Types.hh:76
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio)
G4double GetZ() const
Definition: G4Element.hh:131
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
CLHEP::HepLorentzVector G4LorentzVector
Here is the call graph for this function:

◆ SetFixedCut()

void G4hCoulombScatteringModel::SetFixedCut ( G4double  val)
inline

Definition at line 183 of file G4hCoulombScatteringModel.hh.

184 {
185  fixedCut = val;
186 }

◆ SetLowEnergyThreshold()

void G4hCoulombScatteringModel::SetLowEnergyThreshold ( G4double  val)
inline

◆ SetRecoilThreshold()

void G4hCoulombScatteringModel::SetRecoilThreshold ( G4double  eth)
inline

Definition at line 176 of file G4hCoulombScatteringModel.hh.

177 {
178  recoilThreshold = eth;
179 }

◆ SetupParticle()

void G4hCoulombScatteringModel::SetupParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 164 of file G4hCoulombScatteringModel.hh.

165 {
166  // Initialise mass and charge
167  if(p != particle) {
168  particle = p;
169  mass = particle->GetPDGMass();
170  wokvi->SetupParticle(p);
171  }
172 }
void SetupParticle(const G4ParticleDefinition *)
const G4ParticleDefinition * particle
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ cosThetaMax

G4double G4hCoulombScatteringModel::cosThetaMax
private

Definition at line 135 of file G4hCoulombScatteringModel.hh.

◆ cosThetaMin

G4double G4hCoulombScatteringModel::cosThetaMin
private

Definition at line 134 of file G4hCoulombScatteringModel.hh.

◆ currentCouple

const G4MaterialCutsCouple* G4hCoulombScatteringModel::currentCouple
private

Definition at line 130 of file G4hCoulombScatteringModel.hh.

◆ currentMaterial

const G4Material* G4hCoulombScatteringModel::currentMaterial
private

Definition at line 131 of file G4hCoulombScatteringModel.hh.

◆ currentMaterialIndex

G4int G4hCoulombScatteringModel::currentMaterialIndex
private

Definition at line 132 of file G4hCoulombScatteringModel.hh.

◆ elecRatio

G4double G4hCoulombScatteringModel::elecRatio
private

Definition at line 137 of file G4hCoulombScatteringModel.hh.

◆ fixedCut

G4double G4hCoulombScatteringModel::fixedCut
private

Definition at line 140 of file G4hCoulombScatteringModel.hh.

◆ fNistManager

G4NistManager* G4hCoulombScatteringModel::fNistManager
private

Definition at line 126 of file G4hCoulombScatteringModel.hh.

◆ fParticleChange

G4ParticleChangeForGamma* G4hCoulombScatteringModel::fParticleChange
private

Definition at line 124 of file G4hCoulombScatteringModel.hh.

◆ isCombined

G4bool G4hCoulombScatteringModel::isCombined
private

Definition at line 146 of file G4hCoulombScatteringModel.hh.

◆ mass

G4double G4hCoulombScatteringModel::mass
private

Definition at line 138 of file G4hCoulombScatteringModel.hh.

◆ particle

const G4ParticleDefinition* G4hCoulombScatteringModel::particle
private

Definition at line 143 of file G4hCoulombScatteringModel.hh.

◆ pCuts

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

Definition at line 128 of file G4hCoulombScatteringModel.hh.

◆ recoilThreshold

G4double G4hCoulombScatteringModel::recoilThreshold
private

Definition at line 136 of file G4hCoulombScatteringModel.hh.

◆ theIonTable

G4IonTable* G4hCoulombScatteringModel::theIonTable
private

Definition at line 123 of file G4hCoulombScatteringModel.hh.

◆ theProton

const G4ParticleDefinition* G4hCoulombScatteringModel::theProton
private

Definition at line 144 of file G4hCoulombScatteringModel.hh.

◆ wokvi

G4WentzelVIRelXSection* G4hCoulombScatteringModel::wokvi
private

Definition at line 125 of file G4hCoulombScatteringModel.hh.


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