57 #include "G4mplIonisationWithDeltaModel.hh"
58 #include "Randomize.hh"
59 #include "G4LossTableManager.hh"
60 #include "G4ParticleChangeForLoss.hh"
61 #include "G4Electron.hh"
62 #include "G4DynamicParticle.hh"
68 G4mplIonisationWithDeltaModel::G4mplIonisationWithDeltaModel(G4double mCharge,
const G4String& nam)
69 : G4VEmModel(nam),G4VEmFluctuationModel(nam),
74 beta2lim(betalim*betalim),
75 bg2lim(beta2lim*(1.0 + beta2lim))
77 nmpl = G4int(
abs(magCharge) * 2 * fine_structure_const + 0.5);
78 if(nmpl > 6) { nmpl = 6; }
79 else if(nmpl < 1) { nmpl = 1; }
80 pi_hbarc2_over_mc2 =
pi * hbarc * hbarc / electron_mass_c2;
81 chargeSquare = magCharge * magCharge;
82 dedxlim = 45.*nmpl*nmpl*GeV*cm2/
g;
84 theElectron = G4Electron::Electron();
85 G4cout <<
"### Monopole ionisation model with d-electron production, Gmag= "
86 << magCharge/eplus << G4endl;
91 G4mplIonisationWithDeltaModel::~G4mplIonisationWithDeltaModel()
97 G4mplIonisationWithDeltaModel::Initialise(
const G4ParticleDefinition*
p,
101 mass = monopole->GetPDGMass();
102 if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
108 G4mplIonisationWithDeltaModel::ComputeDEDXPerVolume(
const G4Material* material,
109 const G4ParticleDefinition* p,
110 G4double kineticEnergy,
113 G4double
tmax = MaxSecondaryEnergy(p,kineticEnergy);
114 G4double cutEnergy =
std::min(tmax, maxEnergy);
115 G4double
tau = kineticEnergy /
mass;
116 G4double
gam = tau + 1.0;
117 G4double bg2 = tau * (tau + 2.0);
118 G4double beta2 = bg2 / (
gam *
gam);
122 G4double dedx = dedxlim*beta*material->GetDensity();
128 if(beta >= betalim) {
129 dedx = ComputeDEDXAhlen(material, bg2, cutEnergy);
133 G4double dedx1 = dedxlim*betalow*material->GetDensity();
134 G4double dedx2 = ComputeDEDXAhlen(material, bg2lim, cutEnergy);
137 G4double kapa2 = beta - betalow;
138 G4double kapa1 = betalim -
beta;
139 dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
148 G4mplIonisationWithDeltaModel::ComputeDEDXAhlen(
const G4Material* material,
152 G4double eDensity = material->GetElectronDensity();
153 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
157 0.5*(
log(2.0 * electron_mass_c2 * bg2*cutEnergy / (eexc*eexc)) - 1.0);
161 if(nmpl > 1) { k = 0.346; }
164 const G4double
B[7] = { 0.0, 0.248, 0.672, 1.022, 1.243, 1.464, 1.685};
166 dedx += 0.5 * k - B[nmpl];
169 G4double
x =
log(bg2)/twoln10;
170 dedx -= material->GetIonisation()->DensityCorrection(x);
173 dedx *= pi_hbarc2_over_mc2 * eDensity * nmpl * nmpl;
175 if (dedx < 0.0) { dedx = 0; }
182 G4mplIonisationWithDeltaModel::ComputeCrossSectionPerElectron(
183 const G4ParticleDefinition* p,
184 G4double kineticEnergy,
186 G4double maxKinEnergy)
188 G4double
cross = 0.0;
189 G4double
tmax = MaxSecondaryEnergy(p, kineticEnergy);
190 G4double maxEnergy =
min(tmax,maxKinEnergy);
191 if(cutEnergy < maxEnergy) {
192 cross = (1.0/cutEnergy - 1.0/maxEnergy)*twopi_mc2_rcl2*chargeSquare;
200 G4mplIonisationWithDeltaModel::ComputeCrossSectionPerAtom(
201 const G4ParticleDefinition* p,
202 G4double kineticEnergy,
203 G4double
Z, G4double,
208 Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
215 G4mplIonisationWithDeltaModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
216 const G4MaterialCutsCouple*,
217 const G4DynamicParticle* dp,
218 G4double minKinEnergy,
221 G4double kineticEnergy = dp->GetKineticEnergy();
222 G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
224 G4double maxKinEnergy =
std::min(maxEnergy,tmax);
225 if(minKinEnergy >= maxKinEnergy) {
return; }
231 G4double totEnergy = kineticEnergy +
mass;
232 G4double etot2 = totEnergy*totEnergy;
233 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*
mass)/etot2;
236 G4double
q = G4UniformRand();
237 G4double deltaKinEnergy = minKinEnergy*maxKinEnergy
238 /(minKinEnergy*(1.0 -
q) + maxKinEnergy*q);
241 G4double totMomentum = totEnergy*
sqrt(beta2);
242 G4double deltaMomentum =
243 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
244 G4double
cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
245 (deltaMomentum * totMomentum);
246 if(cost > 1.0) { cost = 1.0; }
248 G4double sint =
sqrt((1.0 - cost)*(1.0 + cost));
250 G4double
phi = twopi * G4UniformRand() ;
252 G4ThreeVector deltaDirection(sint*
cos(phi),sint*
sin(phi), cost);
253 G4ThreeVector direction = dp->GetMomentumDirection();
254 deltaDirection.rotateUz(direction);
257 G4DynamicParticle*
delta =
258 new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
260 vdp->push_back(delta);
263 kineticEnergy -= deltaKinEnergy;
264 G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
265 finalP = finalP.unit();
267 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
268 fParticleChange->SetProposedMomentumDirection(finalP);
273 G4double G4mplIonisationWithDeltaModel::SampleFluctuations(
274 const G4Material* material,
275 const G4DynamicParticle* dp,
280 G4double siga = Dispersion(material,dp,tmax,length);
281 G4double loss = meanLoss;
283 G4double twomeanLoss = meanLoss + meanLoss;
285 if(twomeanLoss < siga) {
288 loss = twomeanLoss*G4UniformRand();
289 x = (loss - meanLoss)/siga;
290 }
while (1.0 - 0.5*x*x < G4UniformRand());
293 loss = G4RandGauss::shoot(meanLoss,siga);
294 }
while (0.0 > loss || loss > twomeanLoss);
302 G4mplIonisationWithDeltaModel::Dispersion(
const G4Material* material,
303 const G4DynamicParticle* dp,
308 G4double
tau = dp->GetKineticEnergy()/
mass;
310 G4double electronDensity = material->GetElectronDensity();
311 G4double
gam = tau + 1.0;
312 G4double invbeta2 = (gam*
gam)/(tau * (tau+2.0));
313 siga = (invbeta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
314 * electronDensity * chargeSquare;
322 G4mplIonisationWithDeltaModel::MaxSecondaryEnergy(
const G4ParticleDefinition*,
325 G4double tau = kinEnergy/
mass;
326 return 2.0*electron_mass_c2*tau*(tau + 2.);
const double Z[kNumberCalorimeter]
Sin< T >::type sin(const T &t)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Cos< T >::type cos(const T &t)
static const double tmax[3]
Basic3DVector cross(const Basic3DVector &v) const
Vector product, or "cross" product, with a vector of same type.