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

#include <G4LowEWentzelVIModel.hh>

+ Inheritance diagram for G4LowEWentzelVIModel:

Public Member Functions

 G4LowEWentzelVIModel ()
 
virtual ~G4LowEWentzelVIModel ()
 
virtual G4double ComputeTruePathLengthLimit (const G4Track &track, G4double &currentMinimalStep)
 
- Public Member Functions inherited from G4WentzelVIModel
 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
 

Additional Inherited Members

- Protected Member Functions inherited from G4WentzelVIModel
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 inherited from G4WentzelVIModel
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 57 of file G4LowEWentzelVIModel.hh.

Constructor & Destructor Documentation

◆ G4LowEWentzelVIModel()

G4LowEWentzelVIModel::G4LowEWentzelVIModel ( )

Definition at line 57 of file G4LowEWentzelVIModel.cc.

57 :
58 G4WentzelVIModel(false,"LowEnWentzelVI")
59{
61}
void SetSingleScatteringFactor(G4double)

◆ ~G4LowEWentzelVIModel()

G4LowEWentzelVIModel::~G4LowEWentzelVIModel ( )
virtual

Definition at line 65 of file G4LowEWentzelVIModel.cc.

66{}

Member Function Documentation

◆ ComputeTruePathLengthLimit()

G4double G4LowEWentzelVIModel::ComputeTruePathLengthLimit ( const G4Track track,
G4double currentMinimalStep 
)
virtual

Reimplemented from G4WentzelVIModel.

Definition at line 70 of file G4LowEWentzelVIModel.cc.

73{
74 G4double tlimit = currentMinimalStep;
75 const G4DynamicParticle* dp = track.GetDynamicParticle();
77 G4StepStatus stepStatus = sp->GetStepStatus();
79 //G4cout << "G4LowEWentzelVIModel::ComputeTruePathLengthLimit stepStatus= "
80 // << stepStatus << " " << track.GetDefinition()->GetParticleName()
81 // << G4endl;
82
83 // initialisation for each step, lambda may be computed from scratch
89 /*
90 G4cout << "lambdaeff= " << lambdaeff << " Range= " << currentRange
91 << " tlimit= " << tlimit << " 1-cost= " << 1 - cosTetMaxNuc << G4endl;
92 */
93 // extra check for abnormal situation
94 // this check needed to run MSC with eIoni and eBrem inactivated
95 tlimit = std::min(tlimit, currentRange);
96
97 // stop here if small range particle
98 if(tlimit < tlimitminfix) {
99 return ConvertTrueToGeom(tlimit, currentMinimalStep);
100 }
101
102 // pre step
103 G4double presafety = sp->GetSafety();
104 // far from geometry boundary
105 if(currentRange < presafety) {
106 return ConvertTrueToGeom(tlimit, currentMinimalStep);
107 }
108
109 // compute presafety again if presafety <= 0 and no boundary
110 // i.e. when it is needed for optimization purposes
111 if(stepStatus != fGeomBoundary && presafety < tlimitminfix) {
112 presafety = ComputeSafety(sp->GetPosition(), tlimit);
113 if(currentRange < presafety) {
114 return ConvertTrueToGeom(tlimit, currentMinimalStep);
115 }
116 }
117 /*
118 G4cout << "e(MeV)= " << preKinEnergy/MeV
119 << " " << particle->GetParticleName()
120 << " CurLimit(mm)= " << tlimit/mm <<" safety(mm)= " << presafety/mm
121 << " R(mm)= " <<currentRange/mm
122 << " L0(mm^-1)= " << lambdaeff*mm
123 <<G4endl;
124 */
125 // natural limit for high energy
126 G4double rlimit = std::max(facrange*currentRange, lambdaeff);
127 //G4double rlimit = std::max(facrange*currentRange,
128 // 0.7*(1.0 - cosTetMaxNuc)*lambdaeff);
129
130 // low-energy e-
131 rlimit = std::max(rlimit, facsafety*presafety);
132
133 // cut correction
134 //G4double rcut = currentCouple->GetProductionCuts()->GetProductionCut(1);
135 //G4cout << "rcut= " << rcut << " rlimit= " << rlimit << " presafety= "
136 // << presafety << " 1-cosThetaMax= " <<1-cosThetaMax
137 //<< " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc << G4endl;
138 //if(rcut > rlimit) { rlimit = std::min(rlimit, rcut*sqrt(rlimit/rcut)); }
139
140 tlimit = std::min(tlimit, rlimit);
141 tlimit = std::max(tlimit, tlimitminfix);
142
143 // step limit in infinite media
144 tlimit = std::min(tlimit, 50*currentMaterial->GetRadlen()/facgeom);
145
146 //compute geomlimit and force few steps within a volume
148 && stepStatus == fGeomBoundary) {
149
150 G4double geomlimit = ComputeGeomLimit(track, presafety, currentRange);
151 tlimit = std::min(tlimit, geomlimit/facgeom);
152 }
153 /*
154 G4cout << particle->GetParticleName() << " E(MeV)= " << preKinEnergy
155 << " L0= " << lambdaeff << " R= " << currentRange
156 << " tlimit= " << tlimit
157 << " currentMinimalStep= " << currentMinimalStep << G4endl;
158 */
159 return ConvertTrueToGeom(tlimit, currentMinimalStep);
160}
@ fUseDistanceToBoundary
G4StepStatus
Definition: G4StepStatus.hh:40
@ fGeomBoundary
Definition: G4StepStatus.hh:43
double G4double
Definition: G4Types.hh:83
G4double GetKineticEnergy() const
G4double GetRadlen() const
Definition: G4Material.hh:218
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 GetTransportMeanFreePath(const G4ParticleDefinition *part, G4double kinEnergy)
Definition: G4VMscModel.hh:405
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
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
const G4MaterialCutsCouple * currentCouple
void DefineMaterial(const G4MaterialCutsCouple *)
const G4ParticleDefinition * particle
const G4Material * currentMaterial
G4WentzelOKandVIxSection * wokvi

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