Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4WentzelVIModel Class Reference

#include <G4WentzelVIModel.hh>

+ Inheritance diagram for G4WentzelVIModel:

Public Member Functions

 G4WentzelVIModel (G4bool comb=true, const G4String &nam="WentzelVIUni")
 
virtual ~G4WentzelVIModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
virtual void StartTracking (G4Track *) override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=DBL_MAX, G4double emax=DBL_MAX) override
 
virtual G4ThreeVectorSampleScattering (const G4ThreeVector &, G4double safety) override
 
virtual G4double ComputeTruePathLengthLimit (const G4Track &track, G4double &currentMinimalStep) override
 
virtual G4double ComputeGeomPathLength (G4double truePathLength) override
 
virtual G4double ComputeTrueStepLength (G4double geomStepLength) override
 
void SetFixedCut (G4double)
 
G4double GetFixedCut () const
 
void SetWVICrossSection (G4WentzelOKandVIxSection *)
 
G4WentzelOKandVIxSectionGetWVICrossSection ()
 
void SetUseSecondMoment (G4bool)
 
G4bool UseSecondMoment () const
 
G4PhysicsTableGetSecondMomentTable ()
 
G4double SecondMoment (const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double kineticEnergy)
 
void SetSingleScatteringFactor (G4double)
 
void DefineMaterial (const G4MaterialCutsCouple *)
 
- Public Member Functions inherited from G4VMscModel
 G4VMscModel (const G4String &nam)
 
 ~G4VMscModel () override
 
virtual G4double ComputeTruePathLengthLimit (const G4Track &track, G4double &stepLimit)=0
 
virtual G4double ComputeGeomPathLength (G4double truePathLength)=0
 
virtual G4double ComputeTrueStepLength (G4double geomPathLength)=0
 
virtual G4ThreeVectorSampleScattering (const G4ThreeVector &, G4double safety)=0
 
void InitialiseParameters (const G4ParticleDefinition *)
 
void DumpParameters (std::ostream &out) const
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax) override
 
void SetStepLimitType (G4MscStepLimitType)
 
void SetLateralDisplasmentFlag (G4bool val)
 
void SetRangeFactor (G4double)
 
void SetGeomFactor (G4double)
 
void SetSkin (G4double)
 
void SetLambdaLimit (G4double)
 
void SetSafetyFactor (G4double)
 
void SetSampleZ (G4bool)
 
G4VEnergyLossProcessGetIonisation () const
 
void SetIonisation (G4VEnergyLossProcess *, const G4ParticleDefinition *part)
 
G4double ComputeSafety (const G4ThreeVector &position, G4double limit=DBL_MAX)
 
G4double ComputeGeomLimit (const G4Track &, G4double &presafety, G4double limit)
 
G4double GetDEDX (const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetDEDX (const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple, G4double logKineticEnergy)
 
G4double GetRange (const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetRange (const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple, G4double logKineticEnergy)
 
G4double GetEnergy (const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
 
G4double GetTransportMeanFreePath (const G4ParticleDefinition *part, G4double kinEnergy)
 
G4double GetTransportMeanFreePath (const G4ParticleDefinition *part, G4double kinEnergy, G4double logKinEnergy)
 
G4VMscModeloperator= (const G4VMscModel &right)=delete
 
 G4VMscModel (const G4VMscModel &)=delete
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
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 level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
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 MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
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 * > *)
 
virtual 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)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, 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 *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
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) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () 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)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

G4double ComputeTransportXSectionPerVolume (G4double cosTheta)
 
void SetupParticle (const G4ParticleDefinition *)
 
- Protected Member Functions inherited from G4VMscModel
G4ParticleChangeForMSCGetParticleChangeForMSC (const G4ParticleDefinition *p=nullptr)
 
G4double ConvertTrueToGeom (G4double &tLength, G4double &gLength)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

G4WentzelOKandVIxSectionwokvi
 
G4double tlimitminfix
 
G4double ssFactor
 
G4double invssFactor
 
G4double preKinEnergy
 
G4double tPathLength
 
G4double zPathLength
 
G4double lambdaeff
 
G4double currentRange
 
G4double cosTetMaxNuc
 
G4int currentMaterialIndex
 
const G4MaterialCutsCouplecurrentCouple
 
const G4MaterialcurrentMaterial
 
const G4ParticleDefinitionparticle
 
G4ParticleChangeForMSCfParticleChange
 
const G4DataVectorcurrentCuts
 
G4double invsqrt12
 
G4double fixedCut
 
G4double effKinEnergy
 
G4double cosThetaMin
 
G4double cosThetaMax
 
G4PhysicsTablefSecondMoments
 
size_t idx2
 
G4int minNCollisions
 
G4int nelments
 
std::vector< G4doublexsecn
 
std::vector< G4doubleprob
 
G4double xtsec
 
G4double numlimit
 
G4double lowEnergyLimit
 
G4bool singleScatteringMode
 
G4bool isCombined
 
G4bool useSecondMoment
 
- Protected Attributes inherited from G4VMscModel
G4double facrange
 
G4double facgeom
 
G4double facsafety
 
G4double skin
 
G4double dtrl
 
G4double lambdalimit
 
G4double geomMin
 
G4double geomMax
 
G4ThreeVector fDisplacement
 
G4MscStepLimitType steppingAlgorithm
 
G4bool samplez
 
G4bool latDisplasment
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Detailed Description

Definition at line 65 of file G4WentzelVIModel.hh.

Constructor & Destructor Documentation

◆ G4WentzelVIModel()

G4WentzelVIModel::G4WentzelVIModel ( G4bool  comb = true,
const G4String nam = "WentzelVIUni" 
)
explicit

Definition at line 73 of file G4WentzelVIModel.cc.

74 : G4VMscModel(nam),
75 ssFactor(1.05),
76 invssFactor(1.0),
77 currentCouple(nullptr),
78 cosThetaMin(1.0),
79 cosThetaMax(-1.0),
80 fSecondMoments(nullptr),
81 idx2(0),
82 numlimit(0.1),
84 isCombined(comb),
85 useSecondMoment(false)
86{
88 invsqrt12 = 1./sqrt(12.);
89 tlimitminfix = 1.e-6*mm;
90 lowEnergyLimit = 1.0*eV;
91 particle = nullptr;
92 nelments = 5;
93 xsecn.resize(nelments);
94 prob.resize(nelments);
96 fixedCut = -1.0;
97
98 minNCollisions = 10;
99
101 = currentRange = xtsec = cosTetMaxNuc = 0.0;
103
104 fParticleChange = nullptr;
105 currentCuts = nullptr;
106 currentMaterial = nullptr;
107}
const G4MaterialCutsCouple * currentCouple
std::vector< G4double > prob
G4PhysicsTable * fSecondMoments
std::vector< G4double > xsecn
G4ParticleChangeForMSC * fParticleChange
const G4ParticleDefinition * particle
const G4DataVector * currentCuts
const G4Material * currentMaterial
void SetSingleScatteringFactor(G4double)
G4WentzelOKandVIxSection * wokvi

◆ ~G4WentzelVIModel()

G4WentzelVIModel::~G4WentzelVIModel ( )
virtual

Definition at line 111 of file G4WentzelVIModel.cc.

112{
113 delete wokvi;
114 if(fSecondMoments && IsMaster()) {
115 delete fSecondMoments;
116 fSecondMoments = nullptr;
117 }
118}
G4bool IsMaster() const
Definition: G4VEmModel.hh:736

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4WentzelVIModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  KineticEnergy,
G4double  AtomicNumber,
G4double  AtomicWeight = 0.,
G4double  cut = DBL_MAX,
G4double  emax = DBL_MAX 
)
overridevirtual

Reimplemented from G4VEmModel.

Reimplemented in G4WentzelVIRelModel.

Definition at line 223 of file G4WentzelVIModel.cc.

228{
229 G4double cross = 0.0;
230 SetupParticle(p);
231 if(kinEnergy < lowEnergyLimit) { return cross; }
232 if(!CurrentCouple()) {
233 G4Exception("G4WentzelVIModel::ComputeCrossSectionPerAtom", "em0011",
234 FatalException, " G4MaterialCutsCouple is not defined");
235 return 0.0;
236 }
239 if(cosTetMaxNuc < 1.0) {
240 G4double cut = (0.0 < fixedCut) ? fixedCut : cutEnergy;
241 G4double cost = wokvi->SetupTarget(G4lrint(Z), cut);
243 /*
244 if(p->GetParticleName() == "e-")
245 G4cout << "G4WentzelVIModel::CS: Z= " << G4int(Z) << " e(MeV)= "<<kinEnergy
246 << " 1-cosN= " << 1 - cosTetMaxNuc << " cross(bn)= " << cross/barn
247 << " " << particle->GetParticleName() << G4endl;
248 */
249 }
250 return cross;
251}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:480
G4double SetupTarget(G4int Z, G4double cut)
G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax)
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
void DefineMaterial(const G4MaterialCutsCouple *)
void SetupParticle(const G4ParticleDefinition *)
int G4lrint(double ad)
Definition: templates.hh:134

◆ ComputeGeomPathLength()

G4double G4WentzelVIModel::ComputeGeomPathLength ( G4double  truePathLength)
overridevirtual

Implements G4VMscModel.

Definition at line 365 of file G4WentzelVIModel.cc.

366{
367 zPathLength = tPathLength = truelength;
368
369 // small step use only single scattering
370 cosThetaMin = 1.0;
372 //G4cout << "xtsec= " << xtsec << " Nav= "
373 // << zPathLength*xtsec << G4endl;
377
378 } else {
379 //G4cout << "ComputeGeomPathLength: tLength= " << tPathLength
380 // << " Leff= " << lambdaeff << G4endl;
381 // small step
384 zPathLength *= (1.0 - 0.5*tau + tau*tau/6.0);
385
386 // medium step
387 } else {
388 G4double e1 = 0.0;
391 }
392 effKinEnergy = 0.5*(e1 + preKinEnergy);
395 //G4cout << " tLength= "<< tPathLength<< " Leff= " << lambdaeff << G4endl;
399 }
400 }
401 }
402 //G4cout << "Comp.geom: zLength= "<<zPathLength<<" tLength= "
403 // << tPathLength<< " Leff= " << lambdaeff << G4endl;
404 return zPathLength;
405}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
int G4int
Definition: G4Types.hh:85
G4double GetTransportMeanFreePath(const G4ParticleDefinition *part, G4double kinEnergy)
Definition: G4VMscModel.hh:405
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:368
G4double ComputeTransportXSectionPerVolume(G4double cosTheta)
#define DBL_MAX
Definition: templates.hh:62

◆ ComputeTransportXSectionPerVolume()

G4double G4WentzelVIModel::ComputeTransportXSectionPerVolume ( G4double  cosTheta)
protected

Definition at line 716 of file G4WentzelVIModel.cc.

717{
718 // prepare recomputation of x-sections
719 const G4ElementVector* theElementVector = currentMaterial->GetElementVector();
720 const G4double* theAtomNumDensityVector =
723 if(nelm > nelments) {
724 nelments = nelm;
725 xsecn.resize(nelm);
726 prob.resize(nelm);
727 }
728
729 // check consistency
730 xtsec = 0.0;
731 if(cosTetMaxNuc >= cosTheta) { return 0.0; }
732
733 G4double cut = (*currentCuts)[currentMaterialIndex];
734 if(fixedCut > 0.0) { cut = fixedCut; }
735
736 // loop over elements
737 G4double xs = 0.0;
738 for (G4int i=0; i<nelm; ++i) {
739 G4double costm =
740 wokvi->SetupTarget((*theElementVector)[i]->GetZasInt(), cut);
741 G4double density = theAtomNumDensityVector[i];
742
743 G4double esec = 0.0;
744 if(costm < cosTheta) {
745
746 // recompute the transport x-section
747 if(1.0 > cosTheta) {
748 xs += density*wokvi->ComputeTransportCrossSectionPerAtom(cosTheta);
749 }
750 // recompute the total x-section
751 G4double nucsec = wokvi->ComputeNuclearCrossSection(cosTheta, costm);
752 esec = wokvi->ComputeElectronCrossSection(cosTheta, costm);
753 nucsec += esec;
754 if(nucsec > 0.0) { esec /= nucsec; }
755 xtsec += nucsec*density;
756 }
757 xsecn[i] = xtsec;
758 prob[i] = esec;
759 //G4cout << i << " xs= " << xs << " xtsec= " << xtsec
760 // << " 1-cosTheta= " << 1-cosTheta
761 // << " 1-cosTetMaxNuc2= " <<1-cosTetMaxNuc2<< G4endl;
762 }
763
764 //G4cout << "ComputeXS result: xsec(1/mm)= " << xs
765 // << " txsec(1/mm)= " << xtsec <<G4endl;
766 return xs;
767}
std::vector< G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)

Referenced by ComputeGeomPathLength(), and ComputeTrueStepLength().

◆ ComputeTruePathLengthLimit()

G4double G4WentzelVIModel::ComputeTruePathLengthLimit ( const G4Track track,
G4double currentMinimalStep 
)
overridevirtual

Implements G4VMscModel.

Reimplemented in G4LowEWentzelVIModel.

Definition at line 267 of file G4WentzelVIModel.cc.

270{
271 G4double tlimit = currentMinimalStep;
272 const G4DynamicParticle* dp = track.GetDynamicParticle();
274 G4StepStatus stepStatus = sp->GetStepStatus();
275 singleScatteringMode = false;
276
277 //G4cout << "G4WentzelVIModel::ComputeTruePathLengthLimit stepStatus= "
278 // << stepStatus << " " << track.GetDefinition()->GetParticleName()
279 // << G4endl;
280
281 // initialisation for each step, lambda may be computed from scratch
285 const G4double logPreKinEnergy = dp->GetLogKineticEnergy();
289
290 //G4cout << "lambdaeff= " << lambdaeff << " Range= " << currentRange
291 // << " tlimit= " << tlimit << " 1-cost= " << 1 - cosTetMaxNuc << G4endl;
292
293 // extra check for abnormal situation
294 // this check needed to run MSC with eIoni and eBrem inactivated
295 if(tlimit > currentRange) { tlimit = currentRange; }
296
297 // stop here if small range particle
298 if(tlimit < tlimitminfix) {
299 return ConvertTrueToGeom(tlimit, currentMinimalStep);
300 }
301
302 // pre step
303 G4double presafety = sp->GetSafety();
304 // far from geometry boundary
305 if(currentRange < presafety) {
306 return ConvertTrueToGeom(tlimit, currentMinimalStep);
307 }
308
309 // compute presafety again if presafety <= 0 and no boundary
310 // i.e. when it is needed for optimization purposes
311 if(stepStatus != fGeomBoundary && presafety < tlimitminfix) {
312 presafety = ComputeSafety(sp->GetPosition(), tlimit);
313 if(currentRange < presafety) {
314 return ConvertTrueToGeom(tlimit, currentMinimalStep);
315 }
316 }
317 /*
318 G4cout << "e(MeV)= " << preKinEnergy/MeV
319 << " " << particle->GetParticleName()
320 << " CurLimit(mm)= " << tlimit/mm <<" safety(mm)= " << presafety/mm
321 << " R(mm)= " <<currentRange/mm
322 << " L0(mm^-1)= " << lambdaeff*mm
323 << G4endl;
324 */
325 // natural limit for high energy
326 G4double rlimit = std::max(facrange*currentRange,
328
329 // low-energy e-
331 rlimit = std::min(rlimit, facsafety*presafety);
332 }
333
334 // cut correction
336 //G4cout << "rcut= " << rcut << " rlimit= " << rlimit << " presafety= "
337 // << presafety << " 1-cosThetaMax= " <<1-cosThetaMax
338 //<< " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc << G4endl;
339 if(rcut > rlimit) { rlimit = std::min(rlimit, rcut*sqrt(rlimit/rcut)); }
340
341 tlimit = std::min(tlimit, rlimit);
342 tlimit = std::max(tlimit, tlimitminfix);
343
344 // step limit in infinite media
345 tlimit = std::min(tlimit, 50*currentMaterial->GetRadlen()/facgeom);
346
347 //compute geomlimit and force few steps within a volume
349 && stepStatus == fGeomBoundary) {
350
351 G4double geomlimit = ComputeGeomLimit(track, presafety, currentRange);
352 tlimit = std::min(tlimit, geomlimit/facgeom);
353 }
354 /*
355 G4cout << particle->GetParticleName() << " e= " << preKinEnergy
356 << " L0= " << lambdaeff << " R= " << currentRange
357 << " tlimit= " << tlimit
358 << " currentMinimalStep= " << currentMinimalStep << G4endl;
359 */
360 return ConvertTrueToGeom(tlimit, currentMinimalStep);
361}
@ fUseDistanceToBoundary
G4StepStatus
Definition: G4StepStatus.hh:40
@ fGeomBoundary
Definition: G4StepStatus.hh:43
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
G4ProductionCuts * GetProductionCuts() const
G4double GetRadlen() const
Definition: G4Material.hh:218
G4double GetProductionCut(G4int index) const
G4StepPoint * GetPreStepPoint() const
const G4DynamicParticle * GetDynamicParticle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
const G4Step * GetStep() const
G4double facrange
Definition: G4VMscModel.hh:193
G4double ComputeGeomLimit(const G4Track &, G4double &presafety, G4double limit)
Definition: G4VMscModel.hh:287
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:330
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:203
G4double ConvertTrueToGeom(G4double &tLength, G4double &gLength)
Definition: G4VMscModel.hh:277
G4double ComputeSafety(const G4ThreeVector &position, G4double limit=DBL_MAX)
Definition: G4VMscModel.hh:269
G4double facsafety
Definition: G4VMscModel.hh:195
G4double facgeom
Definition: G4VMscModel.hh:194

◆ ComputeTrueStepLength()

G4double G4WentzelVIModel::ComputeTrueStepLength ( G4double  geomStepLength)
overridevirtual

Implements G4VMscModel.

Definition at line 409 of file G4WentzelVIModel.cc.

410{
411 // initialisation of single scattering x-section
412 /*
413 G4cout << "ComputeTrueStepLength: Step= " << geomStepLength
414 << " geomL= " << zPathLength
415 << " Lambda= " << lambdaeff
416 << " 1-cosThetaMaxNuc= " << 1 - cosTetMaxNuc << G4endl;
417 */
419 zPathLength = tPathLength = geomStepLength;
420
421 } else {
422
423 // step defined by transportation
424 // change both geom and true step lengths
425 if(geomStepLength < zPathLength) {
426
427 // single scattering
428 if(G4int(geomStepLength*xtsec) < minNCollisions) {
429 zPathLength = tPathLength = geomStepLength;
432
433 // multiple scattering
434 } else {
435 // small step
436 if(geomStepLength < numlimit*lambdaeff) {
437 G4double tau = geomStepLength/lambdaeff;
438 tPathLength = geomStepLength*(1.0 + 0.5*tau + tau*tau/3.0);
439
440 // energy correction for a big step
441 } else {
442 tPathLength *= geomStepLength/zPathLength;
443 G4double e1 = 0.0;
446 }
447 effKinEnergy = 0.5*(e1 + preKinEnergy);
450 G4double tau = geomStepLength/lambdaeff;
451
452 if(tau < 0.999999) { tPathLength = -lambdaeff*G4Log(1.0 - tau); }
453 else { tPathLength = currentRange; }
454 }
455 zPathLength = geomStepLength;
456 }
457 }
458 }
459 // check of step length
460 // define threshold angle between single and multiple scattering
463 xtsec = 0.0;
464
465 // recompute transport cross section - do not change energy
466 // anymore - cannot be applied for big steps
468 // new computation
470 //G4cout << "%%%% cross= " << cross << " xtsec= " << xtsec
471 // << " 1-cosTMin= " << 1.0 - cosThetaMin << G4endl;
472 if(cross <= 0.0) {
476 cosThetaMin = 1.0;
477 } else if(xtsec > 0.0) {
478
479 lambdaeff = 1./cross;
480 G4double tau = zPathLength*cross;
481 if(tau < numlimit) {
482 tPathLength = zPathLength*(1.0 + 0.5*tau + tau*tau/3.0);
483 } else if(tau < 0.999999) {
484 tPathLength = -lambdaeff*G4Log(1.0 - tau);
485 } else {
487 }
488 }
489 }
490 }
492 /*
493 G4cout <<"Comp.true: zLength= "<<zPathLength<<" tLength= "<<tPathLength
494 <<" Leff(mm)= "<<lambdaeff/mm<<" sig0(1/mm)= " << xtsec <<G4endl;
495 G4cout << particle->GetParticleName() << " 1-cosThetaMin= " << 1-cosThetaMin
496 << " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc
497 << " e(MeV)= " << preKinEnergy/MeV << " "
498 << " SSmode= " << singleScatteringMode << G4endl;
499 */
500 return tPathLength;
501}
G4double G4Log(G4double x)
Definition: G4Log.hh:226

◆ DefineMaterial()

void G4WentzelVIModel::DefineMaterial ( const G4MaterialCutsCouple cup)

Definition at line 211 of file G4WentzelVIModel.cc.

212{
213 if(cup != currentCouple) {
214 currentCouple = cup;
215 SetCurrentCouple(cup);
218 }
219}
const G4Material * GetMaterial() const
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:465

Referenced by ComputeCrossSectionPerAtom(), G4LowEWentzelVIModel::ComputeTruePathLengthLimit(), ComputeTruePathLengthLimit(), Initialise(), and SecondMoment().

◆ GetFixedCut()

G4double G4WentzelVIModel::GetFixedCut ( ) const
inline

Definition at line 217 of file G4WentzelVIModel.hh.

218{
219 return fixedCut;
220}

◆ GetSecondMomentTable()

G4PhysicsTable * G4WentzelVIModel::GetSecondMomentTable ( )
inline

Definition at line 255 of file G4WentzelVIModel.hh.

256{
257 return fSecondMoments;
258}

Referenced by InitialiseLocal().

◆ GetWVICrossSection()

G4WentzelOKandVIxSection * G4WentzelVIModel::GetWVICrossSection ( )
inline

Definition at line 234 of file G4WentzelVIModel.hh.

235{
236 return wokvi;
237}

Referenced by G4WentzelVIRelModel::DefineMaterial().

◆ Initialise()

void G4WentzelVIModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
overridevirtual

Implements G4VEmModel.

Reimplemented in G4WentzelVIRelModel.

Definition at line 122 of file G4WentzelVIModel.cc.

124{
125 // reset parameters
126 SetupParticle(p);
128 currentRange = 0.0;
129
130 if(isCombined) {
132 if(tet <= 0.0) { cosThetaMax = 1.0; }
133 else if(tet < CLHEP::pi) { cosThetaMax = cos(tet); }
134 }
135 //G4cout << "G4WentzelVIModel::Initialise " << p->GetParticleName()
136 // << " " << this << " " << wokvi << G4endl;
137
139 /*
140 G4cout << "G4WentzelVIModel: " << particle->GetParticleName()
141 << " 1-cos(ThetaLimit)= " << 1 - cosThetaMax
142 << " SingScatFactor= " << ssFactor
143 << G4endl;
144 */
145 currentCuts = &cuts;
146
147 // set values of some data members
149
150 // build second moment table only if transport table is build
152 if(useSecondMoment && IsMaster() && table) {
153
154 //G4cout << "### G4WentzelVIModel::Initialise: build 2nd moment table "
155 // << table << G4endl;
158 // Access to materials
159 const G4ProductionCutsTable* theCoupleTable =
161 size_t numOfCouples = theCoupleTable->GetTableSize();
162
163 G4bool splineFlag = true;
164 G4PhysicsVector* aVector = nullptr;
165 G4PhysicsVector* bVector = nullptr;
168 if(emin < emax) {
170 *G4lrint(std::log10(emax/emin));
171 if(n < 3) { n = 3; }
172
173 for(size_t i=0; i<numOfCouples; ++i) {
174
175 //G4cout<< "i= " << i << " Flag= " << fSecondMoments->GetFlag(i)
176 // << G4endl;
177 if(fSecondMoments->GetFlag(i)) {
178 DefineMaterial(theCoupleTable->GetMaterialCutsCouple(i));
179
180 delete (*fSecondMoments)[i];
181 if(!aVector) {
182 aVector = new G4PhysicsLogVector(emin, emax, n);
183 bVector = aVector;
184 } else {
185 bVector = new G4PhysicsVector(*aVector);
186 }
187 for(size_t j=0; j<n; ++j) {
188 G4double e = bVector->Energy(j);
189 bVector->PutValue(j, ComputeSecondMoment(p, e)*e*e);
190 }
191 if(splineFlag) { bVector->FillSecondDerivatives(); }
192 (*fSecondMoments)[i] = bVector;
193 }
194 }
195 }
196 //G4cout << *fSecondMoments << G4endl;
197 }
198}
bool G4bool
Definition: G4Types.hh:86
static G4EmParameters * Instance()
G4int NumberOfBinsPerDecade() const
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
G4bool GetFlag(std::size_t i) const
G4double Energy(std::size_t index) const
void PutValue(std::size_t index, G4double theValue)
void FillSecondDerivatives()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:673
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:659
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:860
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:666
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=nullptr)
Definition: G4VMscModel.cc:91
void InitialiseParameters(const G4ParticleDefinition *)
Definition: G4VMscModel.cc:139
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)

Referenced by G4WentzelVIRelModel::Initialise().

◆ InitialiseLocal()

void G4WentzelVIModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 202 of file G4WentzelVIModel.cc.

204{
205 fSecondMoments = static_cast<G4WentzelVIModel*>(masterModel)
207}
G4PhysicsTable * GetSecondMomentTable()

◆ SampleScattering()

G4ThreeVector & G4WentzelVIModel::SampleScattering ( const G4ThreeVector oldDirection,
G4double  safety 
)
overridevirtual

Implements G4VMscModel.

Definition at line 506 of file G4WentzelVIModel.cc.

508{
509 fDisplacement.set(0.0,0.0,0.0);
510 //G4cout << "!##! G4WentzelVIModel::SampleScattering for "
511 // << particle->GetParticleName() << G4endl;
512
513 // ignore scattering for zero step length and energy below the limit
515 { return fDisplacement; }
516
517 G4double invlambda = 0.0;
518 if(lambdaeff < DBL_MAX) { invlambda = 0.5/lambdaeff; }
519
520 // use average kinetic energy over the step
521 G4double cut = (*currentCuts)[currentMaterialIndex];
522 if(fixedCut > 0.0) { cut = fixedCut; }
523 /*
524 G4cout <<"SampleScat: E0(MeV)= "<< preKinEnergy/MeV
525 << " Leff= " << lambdaeff <<" sig0(1/mm)= " << xtsec
526 << " xmsc= " << tPathLength*invlambda
527 << " safety= " << safety << G4endl;
528 */
529 // step limit due msc
530 G4int nMscSteps = 1;
532 G4double z0 = x0*invlambda;
533 //G4double zzz = 0.0;
534 G4double prob2 = 0.0;
535
536 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
537
538 // large scattering angle case - two step approach
540 static const G4double zzmin = 0.05;
541 if(useSecondMoment) {
542 G4double z1 = invlambda*invlambda;
544 prob2 = (z2 - z1)/(1.5*z1 - z2);
545 }
546 // if(z0 > zzmin && safety > tlimitminfix) {
547 if(z0 > zzmin) {
548 x0 *= 0.5;
549 z0 *= 0.5;
550 nMscSteps = 2;
551 }
552 //if(z0 > zzmin) { zzz = G4Exp(-1.0/z0); }
553 G4double zzz = 0.0;
554 if(z0 > zzmin) {
555 zzz = G4Exp(-1.0/z0);
556 z0 += zzz;
557 prob2 *= (1 + zzz);
558 }
559 prob2 /= (1 + prob2);
560 }
561
562 // step limit due to single scattering
563 G4double x1 = 2*tPathLength;
564 if(0.0 < xtsec) { x1 = -G4Log(rndmEngine->flat())/xtsec; }
565
566 // no scattering case
568 { return fDisplacement; }
569
570 const G4ElementVector* theElementVector =
573
574 // geometry
575 G4double sint, cost, phi;
576 G4ThreeVector temp(0.0,0.0,1.0);
577
578 // current position and direction relative to the end point
579 // because of magnetic field geometry is computed relatively to the
580 // end point of the step
581 G4ThreeVector dir(0.0,0.0,1.0);
582 fDisplacement.set(0.0,0.0,-zPathLength);
583
585
586 // start a loop
587 G4double x2 = x0;
588 G4double step, z;
589 G4bool singleScat;
590 /*
591 G4cout << "Start of the loop x1(mm)= " << x1 << " x2(mm)= " << x2
592 << " 1-cost1= " << 1 - cosThetaMin << " SSmode= " << singleScatteringMode
593 << " xtsec= " << xtsec << " Nst= " << nMscSteps << G4endl;
594 */
595 do {
596
597 //G4cout << "# x1(mm)= "<< x1<< " x2(mm)= "<< x2 << G4endl;
598 // single scattering case
599 if(singleScatteringMode && x1 > x2) {
600 fDisplacement += x2*mscfac*dir;
601 break;
602 }
603
604 // what is next single of multiple?
605 if(x1 <= x2) {
606 step = x1;
607 singleScat = true;
608 } else {
609 step = x2;
610 singleScat = false;
611 }
612
613 //G4cout << "# step(mm)= "<< step<< " singlScat= "<< singleScat << G4endl;
614
615 // new position
616 fDisplacement += step*mscfac*dir;
617
618 if(singleScat) {
619
620 // select element
621 G4int i = 0;
622 if(nelm > 1) {
623 G4double qsec = rndmEngine->flat()*xtsec;
624 for (; i<nelm; ++i) { if(xsecn[i] >= qsec) { break; } }
625 }
626 G4double cosTetM =
627 wokvi->SetupTarget((*theElementVector)[i]->GetZasInt(), cut);
628 //G4cout << "!!! " << cosThetaMin << " " << cosTetM << " "
629 // << prob[i] << G4endl;
630 temp = wokvi->SampleSingleScattering(cosThetaMin, cosTetM, prob[i]);
631
632 // direction is changed
633 temp.rotateUz(dir);
634 dir = temp;
635 //G4cout << dir << G4endl;
636
637 // new proposed step length
638 x2 -= step;
639 x1 = -G4Log(rndmEngine->flat())/xtsec;
640
641 // multiple scattering
642 } else {
643 --nMscSteps;
644 x1 -= step;
645 x2 = x0;
646
647 // sample z in interval 0 - 1
648 G4bool isFirst = true;
649 if(prob2 > 0.0 && rndmEngine->flat() < prob2) { isFirst = false; }
650 do {
651 //z = -z0*G4Log(1.0 - (1.0 - zzz)*rndmEngine->flat());
652 if(isFirst) { z = -G4Log(rndmEngine->flat()); }
653 else { z = G4RandGamma::shoot(rndmEngine, 2.0, 2.0); }
654 z *= z0;
655 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
656 } while(z > 1.0);
657
658 cost = 1.0 - 2.0*z/*factCM*/;
659 if(cost > 1.0) { cost = 1.0; }
660 else if(cost < -1.0) { cost =-1.0; }
661 sint = sqrt((1.0 - cost)*(1.0 + cost));
662 phi = twopi*rndmEngine->flat();
663 G4double vx1 = sint*cos(phi);
664 G4double vy1 = sint*sin(phi);
665
666 // lateral displacement
667 if (latDisplasment) {
668 G4double rms = invsqrt12*sqrt(2*z0);
669 G4double r = x0*mscfac;
670 G4double dx = r*(0.5*vx1 + rms*G4RandGauss::shoot(rndmEngine,0.0,1.0));
671 G4double dy = r*(0.5*vy1 + rms*G4RandGauss::shoot(rndmEngine,0.0,1.0));
672 G4double d = r*r - dx*dx - dy*dy;
673
674 // change position
675 if(d >= 0.0) {
676 temp.set(dx,dy,sqrt(d) - r);
677 temp.rotateUz(dir);
678 fDisplacement += temp;
679 }
680 }
681 // change direction
682 temp.set(vx1,vy1,cost);
683 temp.rotateUz(dir);
684 dir = temp;
685 }
686 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
687 } while (0 < nMscSteps);
688
689 dir.rotateUz(oldDirection);
690
691 //G4cout<<"G4WentzelVIModel sampling is done 1-cost= "<< 1.-dir.z()<<G4endl;
692 // end of sampling -------------------------------
693
695
696 // lateral displacement
697 fDisplacement.rotateUz(oldDirection);
698
699 /*
700 G4cout << " r(mm)= " << fDisplacement.mag()
701 << " safety= " << safety
702 << " trueStep(mm)= " << tPathLength
703 << " geomStep(mm)= " << zPathLength
704 << " x= " << fDisplacement.x()
705 << " y= " << fDisplacement.y()
706 << " z= " << fDisplacement.z()
707 << G4endl;
708 */
709
710 //G4cout<< "G4WentzelVIModel::SampleScattering end NewDir= " << dir<< G4endl;
711 return fDisplacement;
712}
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
virtual double flat()=0
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4bool latDisplasment
Definition: G4VMscModel.hh:206
G4ThreeVector fDisplacement
Definition: G4VMscModel.hh:202
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio)
G4double SecondMoment(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double kineticEnergy)

◆ SecondMoment()

G4double G4WentzelVIModel::SecondMoment ( const G4ParticleDefinition part,
const G4MaterialCutsCouple couple,
G4double  kineticEnergy 
)
inline

Definition at line 263 of file G4WentzelVIModel.hh.

266{
267 G4double x = 0.0;
268 if(useSecondMoment) {
269 DefineMaterial(couple);
270 x = (fSecondMoments) ?
271 (*fSecondMoments)[(*theDensityIdx)[currentMaterialIndex]]->Value(ekin, idx2)
272 *(*theDensityFactor)[currentMaterialIndex]/(ekin*ekin)
273 : ComputeSecondMoment(part, ekin);
274 }
275 return x;
276}

Referenced by SampleScattering().

◆ SetFixedCut()

void G4WentzelVIModel::SetFixedCut ( G4double  val)
inline

Definition at line 210 of file G4WentzelVIModel.hh.

211{
212 fixedCut = val;
213}

◆ SetSingleScatteringFactor()

void G4WentzelVIModel::SetSingleScatteringFactor ( G4double  val)

Definition at line 801 of file G4WentzelVIModel.cc.

802{
803 if(val > 0.05) {
804 ssFactor = val;
805 invssFactor = 1.0/(val - 0.05);
806 }
807}

Referenced by G4LowEWentzelVIModel::G4LowEWentzelVIModel(), and G4WentzelVIModel().

◆ SetupParticle()

void G4WentzelVIModel::SetupParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 199 of file G4WentzelVIModel.hh.

200{
201 // Initialise mass and charge
202 if(p != particle) {
203 particle = p;
205 }
206}
void SetupParticle(const G4ParticleDefinition *)

Referenced by ComputeCrossSectionPerAtom(), G4WentzelVIRelModel::ComputeCrossSectionPerAtom(), Initialise(), and StartTracking().

◆ SetUseSecondMoment()

void G4WentzelVIModel::SetUseSecondMoment ( G4bool  val)
inline

Definition at line 241 of file G4WentzelVIModel.hh.

242{
243 useSecondMoment = val;
244}

◆ SetWVICrossSection()

void G4WentzelVIModel::SetWVICrossSection ( G4WentzelOKandVIxSection ptr)
inline

Definition at line 224 of file G4WentzelVIModel.hh.

225{
226 if(ptr != wokvi) {
227 delete wokvi;
228 wokvi = ptr;
229 }
230}

Referenced by G4WentzelVIRelModel::G4WentzelVIRelModel().

◆ StartTracking()

void G4WentzelVIModel::StartTracking ( G4Track track)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 255 of file G4WentzelVIModel.cc.

256{
257 /*
258 G4cout << "G4WentzelVIModel::StartTracking " << track << " " << this << " "
259 << track->GetParticleDefinition()->GetParticleName()
260 << " workvi: " << wokvi << G4endl;
261 */
263}
const G4ParticleDefinition * GetParticleDefinition() const

◆ UseSecondMoment()

G4bool G4WentzelVIModel::UseSecondMoment ( ) const
inline

Definition at line 248 of file G4WentzelVIModel.hh.

249{
250 return useSecondMoment;
251}

Member Data Documentation

◆ cosTetMaxNuc

◆ cosThetaMax

G4double G4WentzelVIModel::cosThetaMax
protected

Definition at line 173 of file G4WentzelVIModel.hh.

Referenced by ComputeTruePathLengthLimit(), and Initialise().

◆ cosThetaMin

G4double G4WentzelVIModel::cosThetaMin
protected

◆ currentCouple

◆ currentCuts

const G4DataVector* G4WentzelVIModel::currentCuts
protected

Definition at line 163 of file G4WentzelVIModel.hh.

Referenced by G4WentzelVIModel(), and Initialise().

◆ currentMaterial

◆ currentMaterialIndex

G4int G4WentzelVIModel::currentMaterialIndex
protected

◆ currentRange

◆ effKinEnergy

G4double G4WentzelVIModel::effKinEnergy
protected

◆ fixedCut

◆ fParticleChange

G4ParticleChangeForMSC* G4WentzelVIModel::fParticleChange
protected

Definition at line 162 of file G4WentzelVIModel.hh.

Referenced by G4WentzelVIModel(), Initialise(), and SampleScattering().

◆ fSecondMoments

G4PhysicsTable* G4WentzelVIModel::fSecondMoments
protected

◆ idx2

size_t G4WentzelVIModel::idx2
protected

Definition at line 176 of file G4WentzelVIModel.hh.

Referenced by SecondMoment().

◆ invsqrt12

G4double G4WentzelVIModel::invsqrt12
protected

Definition at line 165 of file G4WentzelVIModel.hh.

Referenced by G4WentzelVIModel(), and SampleScattering().

◆ invssFactor

G4double G4WentzelVIModel::invssFactor
protected

Definition at line 146 of file G4WentzelVIModel.hh.

Referenced by ComputeTruePathLengthLimit(), and SetSingleScatteringFactor().

◆ isCombined

G4bool G4WentzelVIModel::isCombined
protected

Definition at line 192 of file G4WentzelVIModel.hh.

Referenced by G4WentzelVIModel(), and Initialise().

◆ lambdaeff

◆ lowEnergyLimit

G4double G4WentzelVIModel::lowEnergyLimit
protected

◆ minNCollisions

G4int G4WentzelVIModel::minNCollisions
protected

◆ nelments

G4int G4WentzelVIModel::nelments
protected

Definition at line 180 of file G4WentzelVIModel.hh.

Referenced by ComputeTransportXSectionPerVolume(), and G4WentzelVIModel().

◆ numlimit

G4double G4WentzelVIModel::numlimit
protected

Definition at line 185 of file G4WentzelVIModel.hh.

Referenced by ComputeGeomPathLength(), and ComputeTrueStepLength().

◆ particle

◆ preKinEnergy

◆ prob

std::vector<G4double> G4WentzelVIModel::prob
protected

◆ singleScatteringMode

◆ ssFactor

G4double G4WentzelVIModel::ssFactor
protected

Definition at line 145 of file G4WentzelVIModel.hh.

Referenced by ComputeTrueStepLength(), and SetSingleScatteringFactor().

◆ tlimitminfix

G4double G4WentzelVIModel::tlimitminfix
protected

◆ tPathLength

G4double G4WentzelVIModel::tPathLength
protected

◆ useSecondMoment

G4bool G4WentzelVIModel::useSecondMoment
protected

◆ wokvi

◆ xsecn

std::vector<G4double> G4WentzelVIModel::xsecn
protected

◆ xtsec

◆ zPathLength

G4double G4WentzelVIModel::zPathLength
protected

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