76 currentCouple(nullptr),
80 singleScatteringMode(false)
82 invsqrt12 = 1./sqrt(12.);
83 tlimitminfix = 1.e-6*
mm;
84 lowEnergyLimit = 1.0*
eV;
87 xsecn.resize(nelments);
88 prob.resize(nelments);
94 preKinEnergy = tPathLength = zPathLength = lambdaeff = currentRange =
96 currentMaterialIndex = 0;
97 cosThetaMax = cosTetMaxNuc = 1.0;
99 fParticleChange =
nullptr;
100 currentCuts =
nullptr;
101 currentMaterial =
nullptr;
122 if(tet >=
pi) { cosThetaMax = -1.0; }
123 else if(tet > 0.0) { cosThetaMax = cos(tet); }
147 if(p != particle) { SetupParticle(p); }
148 if(kinEnergy < lowEnergyLimit) {
return cross; }
150 G4Exception(
"G4WentzelVIRelModel::ComputeCrossSectionPerAtom",
"em0011",
162 if(cosTetMaxNuc < 1.0) {
191 G4double tlimit = currentMinimalStep;
195 singleScatteringMode =
false;
204 currentRange =
GetRange(particle,preKinEnergy,currentCouple);
207 cosTetMaxNuc = wokvi->
SetupKinematic(preKinEnergy, currentMaterial,
212 if(tlimit > currentRange) { tlimit = currentRange; }
215 if(inside || tlimit < tlimitminfix) {
222 if(currentRange < presafety) {
229 if(stepStatus !=
fGeomBoundary && presafety < tlimitminfix) {
231 if(currentRange < presafety) {
247 0.7*(1.0 - cosTetMaxNuc)*lambdaeff);
250 if(cosThetaMax > cosTetMaxNuc) {
260 if(rcut > rlimit) { rlimit =
std::min(rlimit, rcut*sqrt(rlimit/rcut)); }
262 if(rlimit < tlimit) { tlimit = rlimit; }
264 tlimit =
std::max(tlimit, tlimitminfix);
289 tPathLength = truelength;
290 zPathLength = tPathLength;
292 if(lambdaeff > 0.0 && lambdaeff !=
DBL_MAX) {
293 G4double tau = tPathLength/lambdaeff;
298 zPathLength *= (1.0 - 0.5*tau + tau*tau/6.0);
303 if(currentRange > tPathLength) {
304 e1 =
GetEnergy(particle,currentRange-tPathLength,currentCouple);
306 e1 = 0.5*(e1 + preKinEnergy);
310 zPathLength = lambdaeff*(1.0 -
G4Exp(-tPathLength/lambdaeff));
312 }
else { lambdaeff =
DBL_MAX; }
324 cosThetaMin = cosTetMaxNuc;
330 singleScatteringMode =
true;
331 zPathLength = geomStepLength;
332 tPathLength = geomStepLength;
339 static const G4double singleScatLimit = 1.0e-7;
340 if(geomStepLength < lambdaeff*singleScatLimit*(1.0 - cosTetMaxNuc)) {
341 singleScatteringMode =
true;
343 zPathLength = geomStepLength;
344 tPathLength = geomStepLength;
347 }
else if(geomStepLength != zPathLength) {
350 zPathLength = geomStepLength;
351 G4double tau = geomStepLength/lambdaeff;
352 tPathLength = zPathLength*(1.0 + 0.5*tau + tau*tau/3.0);
357 if(currentRange > tPathLength) {
358 e1 =
GetEnergy(particle,currentRange-tPathLength,currentCouple);
360 e1 = 0.5*(e1 + preKinEnergy);
364 tau = zPathLength/lambdaeff;
366 if(tau < 0.999999) { tPathLength = -lambdaeff*
G4Log(1.0 - tau); }
367 else { tPathLength = currentRange; }
374 if(!singleScatteringMode) { cosThetaMin = 1.0 - 1.5*tPathLength/lambdaeff; }
378 if(cosThetaMin > cosTetMaxNuc) {
381 G4double cross = ComputeXSectionPerVolume();
384 singleScatteringMode =
true;
385 tPathLength = zPathLength;
387 }
else if(xtsec > 0.0) {
389 lambdaeff = 1./cross;
392 tPathLength = zPathLength*(1.0 + 0.5*tau + tau*tau/3.0);
394 else if(tau < 0.999999) { tPathLength = -lambdaeff*
G4Log(1.0 - tau); }
395 else { tPathLength = currentRange; }
397 if(tPathLength > currentRange) { tPathLength = currentRange; }
423 if(preKinEnergy < lowEnergyLimit || tPathLength <= 0.0)
427 if(lambdaeff <
DBL_MAX) { invlambda = 0.5/lambdaeff; }
430 G4double cut = (*currentCuts)[currentMaterialIndex];
441 static const G4double thinlimit = 0.1;
444 if(tPathLength*invlambda > thinlimit && safety > tlimitminfix) {
466 G4double mscfac = zPathLength/tPathLength;
497 for (; i<nelm; ++i) {
if(xsecn[i] >= qsec) {
break; } }
512 if(x2 <= 0.0) { --nMscSteps; }
521 if(!singleScatteringMode) {
530 if(z0 > 0.01) { zzz =
G4Exp(-1.0/z0); }
536 if(cost > 1.0) { cost = 1.0; }
537 else if(cost < -1.0) { cost =-1.0; }
538 sint = sqrt((1.0 - cost)*(1.0 + cost));
545 G4double rms = invsqrt12*sqrt(2*z0);
551 if(d >= 0.0) { dz = sqrt(d) - r; }
552 else { dx = dy = dz = 0.0; }
560 temp.
set(vx1,vy1,cost);
566 }
while (0 < nMscSteps);
596 G4double G4WentzelVIRelModel::ComputeXSectionPerVolume()
600 const G4double* theAtomNumDensityVector =
603 if(nelm > nelments) {
608 G4double cut = (*currentCuts)[currentMaterialIndex];
613 if(cosTetMaxNuc > cosThetaMin) {
return 0.0; }
617 for (
G4int i=0; i<nelm; ++i) {
619 wokvi->
SetupTarget((*theElementVector)[i]->GetZasInt(), cut);
620 G4double density = theAtomNumDensityVector[i];
623 if(costm < cosThetaMin) {
626 if(1.0 > cosThetaMin) {
633 if(nucsec > 0.0) { esec /= nucsec; }
634 xtsec += nucsec*density;
void set(double x, double y, double z)
static G4Pow * GetInstance()
ThreeVector shoot(const G4int Ap, const G4int Af)
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
virtual void StartTracking(G4Track *)
static G4LossTableManager * Instance()
static constexpr double mm
std::vector< G4Element * > G4ElementVector
static constexpr double proton_mass_c2
G4double GetKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetProductionCut(G4int index) const
G4StepStatus GetStepStatus() const
G4double ConvertTrueToGeom(G4double &tLength, G4double &gLength)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double ¤tMinimalStep) override
G4ParticleDefinition * GetDefinition() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual G4double ComputeGeomPathLength(G4double truePathLength) override
virtual void StartTracking(G4Track *) override
G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax)
const G4Step * GetStep() const
G4double ComputeSafety(const G4ThreeVector &position, G4double limit=DBL_MAX)
const G4ElementVector * GetElementVector() const
static G4NistManager * Instance()
static constexpr double twopi
G4StepPoint * GetPreStepPoint() const
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
const G4double * GetVecNbOfAtomsPerVolume() const
const G4MaterialCutsCouple * CurrentCouple() const
const G4ThreeVector & GetPosition() const
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat, G4double cut, G4double tmass)
G4double SetupTarget(G4int Z, G4double cut)
Hep3Vector & rotateUz(const Hep3Vector &)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=DBL_MAX, G4double emax=DBL_MAX) override
static constexpr double eV
static constexpr double amu_c2
G4double ComputeGeomLimit(const G4Track &, G4double &presafety, G4double limit)
G4ThreeVector fDisplacement
G4double GetTransportMeanFreePath(const G4ParticleDefinition *part, G4double kinEnergy)
G4double GetRadlen() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
G4WentzelVIRelModel(G4bool combined=true)
G4double G4Log(G4double x)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=nullptr)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
virtual ~G4WentzelVIRelModel()
G4double PolarAngleLimit() const
G4double GetAtomicMassAmu(const G4String &symb) const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double GetSafety() const
G4MscStepLimitType steppingAlgorithm
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety) override
virtual G4double ComputeTrueStepLength(G4double geomStepLength) override
static constexpr double pi
size_t GetNumberOfElements() const
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio)
G4ProductionCuts * GetProductionCuts() const