Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4LowEWentzelVIModel Class Reference

#include <G4LowEWentzelVIModel.hh>

Inheritance diagram for G4LowEWentzelVIModel:
Collaboration 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 ComputeGeomPathLength (G4double truePathLength) override
 
virtual G4double ComputeTrueStepLength (G4double geomStepLength) override
 
void SetFixedCut (G4double)
 
G4double GetFixedCut () const
 
G4WentzelOKandVIxSectionGetWVICrossSection ()
 
void SetUseSecondMoment (G4bool)
 
G4bool UseSecondMoment () const
 
G4PhysicsTableGetSecondMomentTable ()
 
G4double SecondMoment (const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double kineticEnergy)
 
void SetSingleScatteringFactor (G4double)
 
- Public Member Functions inherited from G4VMscModel
 G4VMscModel (const G4String &nam)
 
virtual ~G4VMscModel ()
 
virtual 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 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 GetRange (const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetEnergy (const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
 
G4double GetTransportMeanFreePath (const G4ParticleDefinition *part, G4double kinEnergy)
 
- 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 level, const G4ParticleDefinition *, G4double kineticEnergy)
 
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 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)
 
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=nullptr)
 
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)
 

Additional Inherited Members

- Protected Member Functions inherited from G4WentzelVIModel
void DefineMaterial (const G4MaterialCutsCouple *)
 
- 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
 
G4bool singleScatteringMode
 
- 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 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 58 of file G4LowEWentzelVIModel.hh.

Constructor & Destructor Documentation

G4LowEWentzelVIModel::G4LowEWentzelVIModel ( )

Definition at line 58 of file G4LowEWentzelVIModel.cc.

58  :
59  G4WentzelVIModel(false,"LowEnWentzelVI")
60 {
62 }
void SetSingleScatteringFactor(G4double)
G4WentzelVIModel(G4bool comb=true, const G4String &nam="WentzelVIUni")

Here is the call graph for this function:

G4LowEWentzelVIModel::~G4LowEWentzelVIModel ( )
virtual

Definition at line 66 of file G4LowEWentzelVIModel.cc.

67 {}

Member Function Documentation

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

Reimplemented from G4WentzelVIModel.

Definition at line 71 of file G4LowEWentzelVIModel.cc.

74 {
75  G4double tlimit = currentMinimalStep;
76  const G4DynamicParticle* dp = track.GetDynamicParticle();
77  G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
78  G4StepStatus stepStatus = sp->GetStepStatus();
79  singleScatteringMode = false;
80  //G4cout << "G4LowEWentzelVIModel::ComputeTruePathLengthLimit stepStatus= "
81  // << stepStatus << " " << track.GetDefinition()->GetParticleName()
82  // << G4endl;
83 
84  // initialisation for each step, lambda may be computed from scratch
90  /*
91  G4cout << "lambdaeff= " << lambdaeff << " Range= " << currentRange
92  << " tlimit= " << tlimit << " 1-cost= " << 1 - cosTetMaxNuc << G4endl;
93  */
94  // extra check for abnormal situation
95  // this check needed to run MSC with eIoni and eBrem inactivated
96  tlimit = std::min(tlimit, currentRange);
97 
98  // stop here if small range particle
99  if(tlimit < tlimitminfix) {
100  return ConvertTrueToGeom(tlimit, currentMinimalStep);
101  }
102 
103  // pre step
104  G4double presafety = sp->GetSafety();
105  // far from geometry boundary
106  if(currentRange < presafety) {
107  return ConvertTrueToGeom(tlimit, currentMinimalStep);
108  }
109 
110  // compute presafety again if presafety <= 0 and no boundary
111  // i.e. when it is needed for optimization purposes
112  if(stepStatus != fGeomBoundary && presafety < tlimitminfix) {
113  presafety = ComputeSafety(sp->GetPosition(), tlimit);
114  if(currentRange < presafety) {
115  return ConvertTrueToGeom(tlimit, currentMinimalStep);
116  }
117  }
118  /*
119  G4cout << "e(MeV)= " << preKinEnergy/MeV
120  << " " << particle->GetParticleName()
121  << " CurLimit(mm)= " << tlimit/mm <<" safety(mm)= " << presafety/mm
122  << " R(mm)= " <<currentRange/mm
123  << " L0(mm^-1)= " << lambdaeff*mm
124  <<G4endl;
125  */
126  // natural limit for high energy
127  G4double rlimit = std::max(facrange*currentRange, lambdaeff);
128  //G4double rlimit = std::max(facrange*currentRange,
129  // 0.7*(1.0 - cosTetMaxNuc)*lambdaeff);
130 
131  // low-energy e-
132  rlimit = std::max(rlimit, facsafety*presafety);
133 
134  // cut correction
135  //G4double rcut = currentCouple->GetProductionCuts()->GetProductionCut(1);
136  //G4cout << "rcut= " << rcut << " rlimit= " << rlimit << " presafety= "
137  // << presafety << " 1-cosThetaMax= " <<1-cosThetaMax
138  //<< " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc << G4endl;
139  //if(rcut > rlimit) { rlimit = std::min(rlimit, rcut*sqrt(rlimit/rcut)); }
140 
141  tlimit = std::min(tlimit, rlimit);
142  tlimit = std::max(tlimit, tlimitminfix);
143 
144  // step limit in infinite media
145  tlimit = std::min(tlimit, 50*currentMaterial->GetRadlen()/facgeom);
146 
147  //compute geomlimit and force few steps within a volume
149  && stepStatus == fGeomBoundary) {
150 
151  G4double geomlimit = ComputeGeomLimit(track, presafety, currentRange);
152  tlimit = std::min(tlimit, geomlimit/facgeom);
153  }
154  /*
155  G4cout << particle->GetParticleName() << " E(MeV)= " << preKinEnergy
156  << " L0= " << lambdaeff << " R= " << currentRange
157  << " tlimit= " << tlimit
158  << " currentMinimalStep= " << currentMinimalStep << G4endl;
159  */
160  return ConvertTrueToGeom(tlimit, currentMinimalStep);
161 }
G4double facgeom
Definition: G4VMscModel.hh:176
G4WentzelOKandVIxSection * wokvi
G4double GetKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4double facrange
Definition: G4VMscModel.hh:175
G4StepStatus GetStepStatus() const
void DefineMaterial(const G4MaterialCutsCouple *)
G4double ConvertTrueToGeom(G4double &tLength, G4double &gLength)
Definition: G4VMscModel.hh:245
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
const G4MaterialCutsCouple * currentCouple
const G4Material * currentMaterial
const G4Step * GetStep() const
G4double ComputeSafety(const G4ThreeVector &position, G4double limit=DBL_MAX)
Definition: G4VMscModel.hh:237
G4StepStatus
Definition: G4StepStatus.hh:49
G4StepPoint * GetPreStepPoint() const
const G4ThreeVector & GetPosition() const
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:283
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
G4double ComputeGeomLimit(const G4Track &, G4double &presafety, G4double limit)
Definition: G4VMscModel.hh:255
G4double GetTransportMeanFreePath(const G4ParticleDefinition *part, G4double kinEnergy)
Definition: G4VMscModel.hh:340
G4double GetRadlen() const
Definition: G4Material.hh:220
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double facsafety
Definition: G4VMscModel.hh:177
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
const G4ParticleDefinition * particle
G4double GetSafety() const
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:185
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:


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