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

#include <G4VMscModel.hh>

Inheritance diagram for G4VMscModel:
G4VEmModel G4DummyModel G4GoudsmitSaundersonMscModel G4UrbanMscModel G4WentzelVIModel G4WentzelVIRelModel

Public Member Functions

 G4VMscModel (const G4String &nam)
 
virtual ~G4VMscModel ()
 
virtual G4double ComputeTruePathLengthLimit (const G4Track &track, G4double &stepLimit)
 
virtual G4double ComputeGeomPathLength (G4double truePathLength)
 
virtual G4double ComputeTrueStepLength (G4double geomPathLength)
 
virtual G4ThreeVectorSampleScattering (const G4ThreeVector &, G4double safety)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax)
 
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)
 
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 Initialise (const G4ParticleDefinition *, const G4DataVector &)=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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Member Functions

G4ParticleChangeForMSCGetParticleChangeForMSC (const G4ParticleDefinition *p=0)
 
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

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
 

Detailed Description

Definition at line 68 of file G4VMscModel.hh.

Constructor & Destructor Documentation

G4VMscModel::G4VMscModel ( const G4String nam)

Definition at line 58 of file G4VMscModel.cc.

References DBL_MAX, and G4LossTableManager::Instance().

58  :
59  G4VEmModel(nam),
60  safetyHelper(0),
61  ionisation(0),
62  facrange(0.04),
63  facgeom(2.5),
64  facsafety(0.3),
65  skin(1.0),
66  dtrl(0.05),
67  lambdalimit(mm),
68  geomMin(1.e-6*CLHEP::mm),
69  geomMax(1.e50*CLHEP::mm),
71  samplez(false),
72  latDisplasment(true)
73 {
74  dedx = 2.0*CLHEP::MeV*CLHEP::cm2/CLHEP::g;
75  localrange = DBL_MAX;
76  localtkin = 0.0;
78  currentPart = 0;
79 }
G4double facgeom
Definition: G4VMscModel.hh:177
static G4LossTableManager * Instance()
G4double dtrl
Definition: G4VMscModel.hh:180
G4double facrange
Definition: G4VMscModel.hh:176
G4double lambdalimit
Definition: G4VMscModel.hh:181
G4bool samplez
Definition: G4VMscModel.hh:188
G4double skin
Definition: G4VMscModel.hh:179
G4bool latDisplasment
Definition: G4VMscModel.hh:189
G4double geomMax
Definition: G4VMscModel.hh:183
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
G4double geomMin
Definition: G4VMscModel.hh:182
G4double facsafety
Definition: G4VMscModel.hh:178
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:186
#define DBL_MAX
Definition: templates.hh:83
G4VMscModel::~G4VMscModel ( )
virtual

Definition at line 83 of file G4VMscModel.cc.

84 {}

Member Function Documentation

G4double G4VMscModel::ComputeGeomLimit ( const G4Track track,
G4double presafety,
G4double  limit 
)
inline

Definition at line 256 of file G4VMscModel.hh.

References G4SafetyHelper::CheckNextStep(), geomMax, G4Track::GetMomentumDirection(), G4StepPoint::GetPosition(), G4Step::GetPreStepPoint(), G4Track::GetStep(), G4Track::GetVolume(), and G4SafetyHelper::GetWorldVolume().

Referenced by G4WentzelVIModel::ComputeTruePathLengthLimit(), G4WentzelVIRelModel::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), and G4UrbanMscModel::ComputeTruePathLengthLimit().

259 {
260  G4double res = geomMax;
261  if(track.GetVolume() != safetyHelper->GetWorldVolume()) {
262  res = safetyHelper->CheckNextStep(
263  track.GetStep()->GetPreStepPoint()->GetPosition(),
264  track.GetMomentumDirection(),
265  limit, presafety);
266  }
267  return res;
268 }
G4double geomMax
Definition: G4VMscModel.hh:183
const G4Step * GetStep() const
G4StepPoint * GetPreStepPoint() const
const G4ThreeVector & GetPosition() const
G4VPhysicalVolume * GetWorldVolume()
const G4ThreeVector & GetMomentumDirection() const
G4double CheckNextStep(const G4ThreeVector &position, const G4ThreeVector &direction, const G4double currentMaxStep, G4double &newSafety)
G4VPhysicalVolume * GetVolume() const
double G4double
Definition: G4Types.hh:76
G4double G4VMscModel::ComputeGeomPathLength ( G4double  truePathLength)
virtual

Reimplemented in G4UrbanMscModel, G4GoudsmitSaundersonMscModel, G4WentzelVIModel, and G4WentzelVIRelModel.

Definition at line 152 of file G4VMscModel.cc.

Referenced by ConvertTrueToGeom().

153 {
154  return truePathLength;
155 }
G4double G4VMscModel::ComputeSafety ( const G4ThreeVector position,
G4double  limit 
)
inline

Definition at line 238 of file G4VMscModel.hh.

References G4SafetyHelper::ComputeSafety().

Referenced by G4WentzelVIModel::ComputeTruePathLengthLimit(), G4WentzelVIRelModel::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), and G4UrbanMscModel::ComputeTruePathLengthLimit().

240 {
241  return safetyHelper->ComputeSafety(position);
242 }
G4double ComputeSafety(const G4ThreeVector &pGlobalPoint, G4double maxRadius=DBL_MAX)
G4double G4VMscModel::ComputeTruePathLengthLimit ( const G4Track track,
G4double stepLimit 
)
virtual

Reimplemented in G4UrbanMscModel, G4GoudsmitSaundersonMscModel, G4WentzelVIModel, and G4WentzelVIRelModel.

Definition at line 145 of file G4VMscModel.cc.

References DBL_MAX.

Referenced by G4VMultipleScattering::AlongStepGetPhysicalInteractionLength().

146 {
147  return DBL_MAX;
148 }
#define DBL_MAX
Definition: templates.hh:83
G4double G4VMscModel::ComputeTrueStepLength ( G4double  geomPathLength)
virtual

Reimplemented in G4UrbanMscModel, G4GoudsmitSaundersonMscModel, G4WentzelVIModel, and G4WentzelVIRelModel.

Definition at line 159 of file G4VMscModel.cc.

Referenced by G4VMultipleScattering::AlongStepDoIt().

160 {
161  return geomPathLength;
162 }
G4double G4VMscModel::ConvertTrueToGeom ( G4double tLength,
G4double gLength 
)
inlineprotected

Definition at line 246 of file G4VMscModel.hh.

References ComputeGeomPathLength().

Referenced by G4WentzelVIModel::ComputeTruePathLengthLimit(), G4WentzelVIRelModel::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), and G4UrbanMscModel::ComputeTruePathLengthLimit().

248 {
249  glength = ComputeGeomPathLength(tlength);
250  // should return true length
251  return tlength;
252 }
virtual G4double ComputeGeomPathLength(G4double truePathLength)
Definition: G4VMscModel.cc:152
G4double G4VMscModel::GetDEDX ( const G4ParticleDefinition part,
G4double  kineticEnergy,
const G4MaterialCutsCouple couple 
)
inline

Definition at line 273 of file G4VMscModel.hh.

References G4VEnergyLossProcess::GetDEDX(), G4ParticleDefinition::GetPDGCharge(), and test::x.

Referenced by G4GoudsmitSaundersonMscModel::SampleScattering(), and G4UrbanMscModel::SampleScattering().

275 {
276  G4double x;
277  if(ionisation) { x = ionisation->GetDEDX(kinEnergy, couple); }
278  else {
279  G4double q = part->GetPDGCharge()/CLHEP::eplus;
280  x = dedx*q*q;
281  }
282  return x;
283 }
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
G4double G4VMscModel::GetEnergy ( const G4ParticleDefinition part,
G4double  range,
const G4MaterialCutsCouple couple 
)
inline

Definition at line 308 of file G4VMscModel.hh.

References G4Material::GetDensity(), G4VEnergyLossProcess::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), and G4ParticleDefinition::GetPDGCharge().

Referenced by G4WentzelVIModel::ComputeGeomPathLength(), G4WentzelVIRelModel::ComputeGeomPathLength(), G4GoudsmitSaundersonMscModel::ComputeGeomPathLength(), G4UrbanMscModel::ComputeGeomPathLength(), G4WentzelVIRelModel::ComputeTrueStepLength(), G4WentzelVIModel::ComputeTrueStepLength(), G4GoudsmitSaundersonMscModel::SampleScattering(), and G4UrbanMscModel::SampleScattering().

310 {
311  G4double e;
312  //G4cout << "G4VMscModel::GetEnergy R(mm)= " << range << " " << ionisation
313  // << " Rlocal(mm)= " << localrange << " Elocal(MeV)= " << localtkin
314  // << G4endl;
315  if(ionisation) { e = ionisation->GetKineticEnergy(range, couple); }
316  else {
317  e = localtkin;
318  if(localrange > range) {
319  G4double q = part->GetPDGCharge()/CLHEP::eplus;
320  e -= (localrange - range)*dedx*q*q*couple->GetMaterial()->GetDensity();
321  }
322  }
323  return e;
324 }
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)
G4double GetDensity() const
Definition: G4Material.hh:178
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
const G4Material * GetMaterial() const
G4VEnergyLossProcess * G4VMscModel::GetIonisation ( ) const
inline

Definition at line 328 of file G4VMscModel.hh.

329 {
330  return ionisation;
331 }
G4ParticleChangeForMSC * G4VMscModel::GetParticleChangeForMSC ( const G4ParticleDefinition p = 0)
protected

Definition at line 89 of file G4VMscModel.cc.

References G4LossTableBuilder::BuildTableForModel(), G4PhysicsTable::clearAndDestroy(), G4VEmModel::ForceBuildTableFlag(), G4LossTableBuilder::GetCoupleIndexes(), G4LossTableBuilder::GetDensityFactors(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGMass(), G4TransportationManager::GetSafetyHelper(), G4LossTableManager::GetTableBuilder(), G4TransportationManager::GetTransportationManager(), python.hepunit::GeV, G4VEmModel::HighEnergyActivationLimit(), G4VEmModel::HighEnergyLimit(), G4VEmModel::idxTable, G4SafetyHelper::InitialiseHelper(), G4VEmModel::IsMaster(), G4VEmModel::LowEnergyActivationLimit(), G4VEmModel::LowEnergyLimit(), G4INCL::Math::max(), G4LossTableManager::MaxKinEnergy(), G4INCL::Math::min(), G4LossTableManager::MinKinEnergy(), G4VEmModel::pParticleChange, G4VEmModel::theDensityFactor, G4VEmModel::theDensityIdx, and G4VEmModel::xSectionTable.

Referenced by G4WentzelVIModel::Initialise(), G4WentzelVIRelModel::Initialise(), G4GoudsmitSaundersonMscModel::Initialise(), and G4UrbanMscModel::Initialise().

90 {
91  // recomputed for each new run
92  if(!safetyHelper) {
94  ->GetSafetyHelper();
95  safetyHelper->InitialiseHelper();
96  }
97  G4ParticleChangeForMSC* change = 0;
98  if (pParticleChange) {
99  change = static_cast<G4ParticleChangeForMSC*>(pParticleChange);
100  } else {
101  change = new G4ParticleChangeForMSC();
102  }
103  if(p) {
104 
105  // table is never built for GenericIon
106  if(p->GetParticleName() == "GenericIon") {
107  if(xSectionTable) {
109  delete xSectionTable;
110  xSectionTable = 0;
111  }
112 
113  // table is always built for low mass particles
114  } else if(p->GetPDGMass() < 4.5*GeV || ForceBuildTableFlag()) {
115 
116  idxTable = 0;
117  G4LossTableBuilder* builder = man->GetTableBuilder();
118  if(IsMaster()) {
121  emin = std::max(emin, man->MinKinEnergy());
122  emax = std::min(emax, man->MaxKinEnergy());
123  if(emin < emax) {
124  xSectionTable = builder->BuildTableForModel(xSectionTable, this, p,
125  emin, emax, true);
126  }
127  }
128  theDensityFactor = builder->GetDensityFactors();
129  theDensityIdx = builder->GetCoupleIndexes();
130  }
131  }
132  return change;
133 }
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:613
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:606
size_t idxTable
Definition: G4VEmModel.hh:402
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
G4SafetyHelper * GetSafetyHelper() const
const std::vector< G4double > * GetDensityFactors()
void InitialiseHelper()
G4bool ForceBuildTableFlag() const
Definition: G4VEmModel.hh:648
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
G4PhysicsTable * BuildTableForModel(G4PhysicsTable *table, G4VEmModel *model, const G4ParticleDefinition *, G4double emin, G4double emax, G4bool spline)
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
const G4String & GetParticleName() const
G4LossTableBuilder * GetTableBuilder()
const std::vector< G4int > * theDensityIdx
Definition: G4VEmModel.hh:401
G4double MinKinEnergy() const
const std::vector< G4double > * theDensityFactor
Definition: G4VEmModel.hh:400
static G4TransportationManager * GetTransportationManager()
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:398
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:399
double G4double
Definition: G4Types.hh:76
G4double MaxKinEnergy() const
void clearAndDestroy()
const std::vector< G4int > * GetCoupleIndexes()
G4double G4VMscModel::GetRange ( const G4ParticleDefinition part,
G4double  kineticEnergy,
const G4MaterialCutsCouple couple 
)
inline

Definition at line 288 of file G4VMscModel.hh.

References G4Material::GetDensity(), G4MaterialCutsCouple::GetMaterial(), G4ParticleDefinition::GetPDGCharge(), and G4VEnergyLossProcess::GetRangeForLoss().

Referenced by G4VMultipleScattering::AlongStepDoIt(), G4WentzelVIRelModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), and G4UrbanMscModel::ComputeTruePathLengthLimit().

290 {
291  //G4cout << "G4VMscModel::GetRange E(MeV)= " << kinEnergy << " " << ionisation
292  // << " " << part->GetParticleName()
293  // << G4endl;
294  localtkin = kinEnergy;
295  if(ionisation) {
296  localrange = ionisation->GetRangeForLoss(kinEnergy, couple);
297  } else {
298  G4double q = part->GetPDGCharge()/CLHEP::eplus;
299  localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
300  }
301  //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
302  return localrange;
303 }
G4double GetRangeForLoss(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4double GetDensity() const
Definition: G4Material.hh:178
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
const G4Material * GetMaterial() const
G4double G4VMscModel::GetTransportMeanFreePath ( const G4ParticleDefinition part,
G4double  kinEnergy 
)
inline

Definition at line 345 of file G4VMscModel.hh.

References G4VEmModel::CrossSectionPerVolume(), G4VEmModel::CurrentCouple(), DBL_MAX, G4MaterialCutsCouple::GetIndex(), G4VEmModel::idxTable, test::x, and G4VEmModel::xSectionTable.

Referenced by G4WentzelVIModel::ComputeGeomPathLength(), G4WentzelVIRelModel::ComputeGeomPathLength(), G4GoudsmitSaundersonMscModel::ComputeGeomPathLength(), G4UrbanMscModel::ComputeGeomPathLength(), G4WentzelVIRelModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), G4UrbanMscModel::ComputeTruePathLengthLimit(), G4WentzelVIRelModel::ComputeTrueStepLength(), and G4WentzelVIModel::ComputeTrueStepLength().

347 {
348  G4double x;
349  if(xSectionTable) {
350  G4int idx = CurrentCouple()->GetIndex();
351  x = (*xSectionTable)[(*theDensityIdx)[idx]]->Value(ekin, idxTable)
352  *(*theDensityFactor)[idx]/(ekin*ekin);
353  } else {
354  x = CrossSectionPerVolume(CurrentCouple()->GetMaterial(), part, ekin,
355  0.0, DBL_MAX);
356  }
357  if(0.0 >= x) { x = DBL_MAX; }
358  else { x = 1.0/x; }
359  return x;
360 }
size_t idxTable
Definition: G4VEmModel.hh:402
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:245
int G4int
Definition: G4Types.hh:78
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:426
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:399
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
G4ThreeVector & G4VMscModel::SampleScattering ( const G4ThreeVector ,
G4double  safety 
)
virtual

Reimplemented in G4UrbanMscModel, G4GoudsmitSaundersonMscModel, G4WentzelVIModel, and G4WentzelVIRelModel.

Definition at line 138 of file G4VMscModel.cc.

References fDisplacement.

Referenced by G4VMultipleScattering::AlongStepDoIt().

139 {
140  return fDisplacement;
141 }
G4ThreeVector fDisplacement
Definition: G4VMscModel.hh:185
void G4VMscModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle ,
G4double  tmin,
G4double  tmax 
)
virtual

Implements G4VEmModel.

Reimplemented in G4DummyModel.

Definition at line 166 of file G4VMscModel.cc.

170 {}
void G4VMscModel::SetGeomFactor ( G4double  val)
inline

Definition at line 217 of file G4VMscModel.hh.

References facgeom.

Referenced by G4AdjointhMultipleScattering::InitialiseProcess(), and G4VMultipleScattering::PreparePhysicsTable().

218 {
219  facgeom = val;
220 }
G4double facgeom
Definition: G4VMscModel.hh:177
void G4VMscModel::SetIonisation ( G4VEnergyLossProcess p,
const G4ParticleDefinition part 
)
inline

Definition at line 335 of file G4VMscModel.hh.

Referenced by G4VMultipleScattering::PreparePhysicsTable(), G4VMultipleScattering::SetIonisation(), and G4VMultipleScattering::StartTracking().

337 {
338  ionisation = p;
339  currentPart = part;
340 }
const char * p
Definition: xmltok.h:285
void G4VMscModel::SetLateralDisplasmentFlag ( G4bool  val)
inline

Definition at line 196 of file G4VMscModel.hh.

References latDisplasment.

Referenced by G4AdjointhMultipleScattering::InitialiseProcess(), and G4VMultipleScattering::PreparePhysicsTable().

197 {
198  latDisplasment = val;
199 }
G4bool latDisplasment
Definition: G4VMscModel.hh:189
void G4VMscModel::SetRangeFactor ( G4double  val)
inline

Definition at line 210 of file G4VMscModel.hh.

References facrange.

Referenced by G4AdjointhMultipleScattering::InitialiseProcess(), and G4VMultipleScattering::PreparePhysicsTable().

211 {
212  facrange = val;
213 }
G4double facrange
Definition: G4VMscModel.hh:176
void G4VMscModel::SetSampleZ ( G4bool  val)
inline

Definition at line 231 of file G4VMscModel.hh.

References samplez.

Referenced by G4UrbanMscModel::G4UrbanMscModel().

232 {
233  samplez = val;
234 }
G4bool samplez
Definition: G4VMscModel.hh:188
void G4VMscModel::SetSkin ( G4double  val)
inline

Definition at line 203 of file G4VMscModel.hh.

References skin.

Referenced by G4AdjointhMultipleScattering::InitialiseProcess(), and G4VMultipleScattering::PreparePhysicsTable().

204 {
205  skin = val;
206 }
G4double skin
Definition: G4VMscModel.hh:179
void G4VMscModel::SetStepLimitType ( G4MscStepLimitType  val)
inline

Definition at line 224 of file G4VMscModel.hh.

References steppingAlgorithm.

Referenced by G4AdjointhMultipleScattering::InitialiseProcess(), and G4VMultipleScattering::PreparePhysicsTable().

225 {
226  steppingAlgorithm = val;
227 }
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:186

Field Documentation

G4double G4VMscModel::dtrl
protected
G4double G4VMscModel::facgeom
protected
G4double G4VMscModel::facrange
protected
G4double G4VMscModel::facsafety
protected
G4ThreeVector G4VMscModel::fDisplacement
protected
G4double G4VMscModel::geomMax
protected

Definition at line 183 of file G4VMscModel.hh.

Referenced by ComputeGeomLimit().

G4double G4VMscModel::geomMin
protected

Definition at line 182 of file G4VMscModel.hh.

G4double G4VMscModel::lambdalimit
protected

Definition at line 181 of file G4VMscModel.hh.

G4bool G4VMscModel::latDisplasment
protected
G4bool G4VMscModel::samplez
protected
G4double G4VMscModel::skin
protected
G4MscStepLimitType G4VMscModel::steppingAlgorithm
protected

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