71std::vector<G4double>* G4mplIonisationWithDeltaModel::dedx0 =
nullptr;
77 twoln10(
std::log(100.0)),
80 beta2lim(betalim*betalim),
81 bg2lim(beta2lim*(1.0 + beta2lim))
83 nmpl =
G4lrint(std::abs(magCharge) * 2 * fine_structure_const);
84 if(nmpl > 6) { nmpl = 6; }
85 else if(nmpl < 1) { nmpl = 1; }
86 pi_hbarc2_over_mc2 = pi * hbarc * hbarc / electron_mass_c2;
87 chargeSquare = magCharge * magCharge;
88 dedxlim = 45.*nmpl*nmpl*GeV*cm2/g;
89 fParticleChange =
nullptr;
91 G4cout <<
"### Monopole ionisation model with d-electron production, Gmag= "
92 << magCharge/eplus <<
G4endl;
111 std::min(
LowEnergyLimit(),0.1*mass*(1./sqrt(1. - betalow*betalow) - 1.));
127 if(!dedx0) { dedx0 =
new std::vector<G4double>; }
132 if(n < numOfCouples) { dedx0->resize(numOfCouples); }
136 for(
G4int i=0; i<numOfCouples; ++i) {
141 G4double vF2 = 2*electron_Compton_length*g4calc->
A13(3.*pi*pi*eDensity);
142 (*dedx0)[i] = pi_hbarc2_over_mc2*eDensity*nmpl*nmpl*
143 (
G4Log(vF2/fine_structure_const) - 0.5)/vF2;
167 G4double cutEnergy = std::min(tmax, maxEnergy);
169 G4double tau = kineticEnergy / mass;
182 if(beta >= betalim) {
183 dedx = ComputeDEDXAhlen(material, bg2, cutEnergy);
187 G4double dedx2 = ComputeDEDXAhlen(material, bg2lim, cutEnergy);
192 dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
201G4mplIonisationWithDeltaModel::ComputeDEDXAhlen(
const G4Material* material,
210 0.5*(
G4Log(2.0*electron_mass_c2*bg2*cutEnergy/(eexc*eexc)) -1.0);
214 if(nmpl > 1) { k = 0.346; }
217 const G4double B[7] = { 0.0, 0.248, 0.672, 1.022, 1.243, 1.464, 1.685};
219 dedx += 0.5 * k -
B[nmpl];
226 dedx *= pi_hbarc2_over_mc2 * eDensity * nmpl * nmpl;
228 dedx = std::max(dedx, 0.0);
243 G4double maxEnergy = std::min(tmax, maxKinEnergy);
245 G4double cross = (cutEnergy < maxEnergy)
246 ? (0.5/cutEnergy - 0.5/maxEnergy)*pi_hbarc2_over_mc2 * nmpl * nmpl : 0.0;
277 G4double maxKinEnergy = std::min(maxEnergy,tmax);
278 if(minKinEnergy >= maxKinEnergy) {
return; }
284 G4double totEnergy = kineticEnergy + mass;
285 G4double etot2 = totEnergy*totEnergy;
286 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
290 G4double deltaKinEnergy = minKinEnergy*maxKinEnergy
291 /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
294 G4double totMomentum = totEnergy*sqrt(beta2);
296 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
297 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
298 (deltaMomentum * totMomentum);
299 cost = std::min(cost, 1.0);
301 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
305 G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost);
313 vdp->push_back(delta);
316 kineticEnergy -= deltaKinEnergy;
317 G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
318 finalP = finalP.
unit();
336 siga = std::sqrt(siga);
337 G4double twomeanLoss = meanLoss + meanLoss;
339 if(twomeanLoss < siga) {
343 x = (loss - meanLoss)/siga;
348 loss = G4RandGauss::shoot(meanLoss,siga);
350 }
while (0.0 > loss || loss > twomeanLoss);
368 siga = (tmax/(beta*beta) - 0.5*tcut) * twopi_mc2_rcl2 * length
381 return 2.0*electron_mass_c2*tau*(tau + 2.);
G4double B(G4double temperature)
G4double G4Log(G4double x)
G4GLOB_DLL std::ostream G4cout
Hep3Vector & rotateUz(const Hep3Vector &)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
G4double DensityCorrection(G4double x) const
G4double GetMeanExcitationEnergy() const
const G4Material * GetMaterial() const
G4IonisParamMat * GetIonisation() const
G4double GetElectronDensity() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double GetPDGMass() const
static G4Pow * GetInstance()
G4double A13(G4double A) const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetHighEnergyLimit(G4double)
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
void SetLowEnergyLimit(G4double)
const G4MaterialCutsCouple * CurrentCouple() const
G4ParticleChangeForLoss * GetParticleChangeForLoss()
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4mplIonisationWithDeltaModel(G4double mCharge, const G4String &nam="mplIonisationWithDelta")
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
virtual G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length, const G4double meanLoss) override
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
void SetParticle(const G4ParticleDefinition *p)
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length) override
virtual ~G4mplIonisationWithDeltaModel()
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override