CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
G4mplIonisationWithDeltaModel.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4mplIonisationWithDeltaModel.cc,v 1.1 2010/07/29 23:05:19 sunanda Exp $
27 // GEANT4 tag $Name: CMSSW_5_3_11_patch3 $
28 //
29 // -------------------------------------------------------------------
30 //
31 // GEANT4 Class header file
32 //
33 //
34 // File name: G4mplIonisationWithDeltaModel
35 //
36 // Author: Vladimir Ivanchenko
37 //
38 // Creation date: 06.09.2005
39 //
40 // Modifications:
41 // 12.08.2007 Changing low energy approximation and extrapolation.
42 // Small bug fixing and refactoring (M. Vladymyrov)
43 // 13.11.2007 Use low-energy asymptotic from [3] (V.Ivanchenko)
44 //
45 //
46 // -------------------------------------------------------------------
47 // References
48 // [1] Steven P. Ahlen: Energy loss of relativistic heavy ionizing particles,
49 // S.P. Ahlen, Rev. Mod. Phys 52(1980), p121
50 // [2] K.A. Milton arXiv:hep-ex/0602040
51 // [3] S.P. Ahlen and K. Kinoshita, Phys. Rev. D26 (1982) 2347
52 
53 
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
56 
57 #include "G4mplIonisationWithDeltaModel.hh"
58 #include "Randomize.hh"
59 #include "G4LossTableManager.hh"
60 #include "G4ParticleChangeForLoss.hh"
61 #include "G4Electron.hh"
62 #include "G4DynamicParticle.hh"
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65 
66 using namespace std;
67 
68 G4mplIonisationWithDeltaModel::G4mplIonisationWithDeltaModel(G4double mCharge, const G4String& nam)
69  : G4VEmModel(nam),G4VEmFluctuationModel(nam),
70  magCharge(mCharge),
71  twoln10(log(100.0)),
72  betalow(0.01),
73  betalim(0.1),
74  beta2lim(betalim*betalim),
75  bg2lim(beta2lim*(1.0 + beta2lim))
76 {
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;
83  fParticleChange = 0;
84  theElectron = G4Electron::Electron();
85  G4cout << "### Monopole ionisation model with d-electron production, Gmag= "
86  << magCharge/eplus << G4endl;
87 }
88 
89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
90 
91 G4mplIonisationWithDeltaModel::~G4mplIonisationWithDeltaModel()
92 {}
93 
94 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
95 
96 void
97 G4mplIonisationWithDeltaModel::Initialise(const G4ParticleDefinition* p,
98  const G4DataVector&)
99 {
100  monopole = p;
101  mass = monopole->GetPDGMass();
102  if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
103 }
104 
105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
106 
107 G4double
108 G4mplIonisationWithDeltaModel::ComputeDEDXPerVolume(const G4Material* material,
109  const G4ParticleDefinition* p,
110  G4double kineticEnergy,
111  G4double maxEnergy)
112 {
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);
119  G4double beta = sqrt(beta2);
120 
121  // low-energy asymptotic formula
122  G4double dedx = dedxlim*beta*material->GetDensity();
123 
124  // above asymptotic
125  if(beta > betalow) {
126 
127  // high energy
128  if(beta >= betalim) {
129  dedx = ComputeDEDXAhlen(material, bg2, cutEnergy);
130 
131  } else {
132 
133  G4double dedx1 = dedxlim*betalow*material->GetDensity();
134  G4double dedx2 = ComputeDEDXAhlen(material, bg2lim, cutEnergy);
135 
136  // extrapolation between two formula
137  G4double kapa2 = beta - betalow;
138  G4double kapa1 = betalim - beta;
139  dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
140  }
141  }
142  return dedx;
143 }
144 
145 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
146 
147 G4double
148 G4mplIonisationWithDeltaModel::ComputeDEDXAhlen(const G4Material* material,
149  G4double bg2,
150  G4double cutEnergy)
151 {
152  G4double eDensity = material->GetElectronDensity();
153  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
154 
155  // Ahlen's formula for nonconductors, [1]p157, f(5.7)
156  G4double dedx =
157  0.5*(log(2.0 * electron_mass_c2 * bg2*cutEnergy / (eexc*eexc)) - 1.0);
158 
159  // Kazama et al. cross-section correction
160  G4double k = 0.406;
161  if(nmpl > 1) { k = 0.346; }
162 
163  // Bloch correction
164  const G4double B[7] = { 0.0, 0.248, 0.672, 1.022, 1.243, 1.464, 1.685};
165 
166  dedx += 0.5 * k - B[nmpl];
167 
168  // density effect correction
169  G4double x = log(bg2)/twoln10;
170  dedx -= material->GetIonisation()->DensityCorrection(x);
171 
172  // now compute the total ionization loss
173  dedx *= pi_hbarc2_over_mc2 * eDensity * nmpl * nmpl;
174 
175  if (dedx < 0.0) { dedx = 0; }
176  return dedx;
177 }
178 
179 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
180 
181 G4double
182 G4mplIonisationWithDeltaModel::ComputeCrossSectionPerElectron(
183  const G4ParticleDefinition* p,
184  G4double kineticEnergy,
185  G4double cutEnergy,
186  G4double maxKinEnergy)
187 {
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;
193  }
194  return cross;
195 }
196 
197 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
198 
199 G4double
200 G4mplIonisationWithDeltaModel::ComputeCrossSectionPerAtom(
201  const G4ParticleDefinition* p,
202  G4double kineticEnergy,
203  G4double Z, G4double,
204  G4double cutEnergy,
205  G4double maxEnergy)
206 {
207  G4double cross =
208  Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
209  return cross;
210 }
211 
212 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
213 
214 void
215 G4mplIonisationWithDeltaModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
216  const G4MaterialCutsCouple*,
217  const G4DynamicParticle* dp,
218  G4double minKinEnergy,
219  G4double maxEnergy)
220 {
221  G4double kineticEnergy = dp->GetKineticEnergy();
222  G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
223 
224  G4double maxKinEnergy = std::min(maxEnergy,tmax);
225  if(minKinEnergy >= maxKinEnergy) { return; }
226 
227  //G4cout << "G4mplIonisationWithDeltaModel::SampleSecondaries: E(GeV)= "
228  // << kineticEnergy/GeV << " M(GeV)= " << mass/GeV
229  // << " tmin(MeV)= " << minKinEnergy/MeV << G4endl;
230 
231  G4double totEnergy = kineticEnergy + mass;
232  G4double etot2 = totEnergy*totEnergy;
233  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
234 
235  // sampling without nuclear size effect
236  G4double q = G4UniformRand();
237  G4double deltaKinEnergy = minKinEnergy*maxKinEnergy
238  /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
239 
240  // delta-electron is produced
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; }
247 
248  G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
249 
250  G4double phi = twopi * G4UniformRand() ;
251 
252  G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost);
253  G4ThreeVector direction = dp->GetMomentumDirection();
254  deltaDirection.rotateUz(direction);
255 
256  // create G4DynamicParticle object for delta ray
257  G4DynamicParticle* delta =
258  new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
259 
260  vdp->push_back(delta);
261 
262  // Change kinematics of primary particle
263  kineticEnergy -= deltaKinEnergy;
264  G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
265  finalP = finalP.unit();
266 
267  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
268  fParticleChange->SetProposedMomentumDirection(finalP);
269 }
270 
271 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
272 
273 G4double G4mplIonisationWithDeltaModel::SampleFluctuations(
274  const G4Material* material,
275  const G4DynamicParticle* dp,
276  G4double& tmax,
277  G4double& length,
278  G4double& meanLoss)
279 {
280  G4double siga = Dispersion(material,dp,tmax,length);
281  G4double loss = meanLoss;
282  siga = sqrt(siga);
283  G4double twomeanLoss = meanLoss + meanLoss;
284 
285  if(twomeanLoss < siga) {
286  G4double x;
287  do {
288  loss = twomeanLoss*G4UniformRand();
289  x = (loss - meanLoss)/siga;
290  } while (1.0 - 0.5*x*x < G4UniformRand());
291  } else {
292  do {
293  loss = G4RandGauss::shoot(meanLoss,siga);
294  } while (0.0 > loss || loss > twomeanLoss);
295  }
296  return loss;
297 }
298 
299 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
300 
301 G4double
302 G4mplIonisationWithDeltaModel::Dispersion(const G4Material* material,
303  const G4DynamicParticle* dp,
304  G4double& tmax,
305  G4double& length)
306 {
307  G4double siga = 0.0;
308  G4double tau = dp->GetKineticEnergy()/mass;
309  if(tau > 0.0) {
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;
315  }
316  return siga;
317 }
318 
319 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
320 
321 G4double
322 G4mplIonisationWithDeltaModel::MaxSecondaryEnergy(const G4ParticleDefinition*,
323  G4double kinEnergy)
324 {
325  G4double tau = kinEnergy/mass;
326  return 2.0*electron_mass_c2*tau*(tau + 2.);
327 }
328 
329 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
const double Z[kNumberCalorimeter]
const double beta
dbl * delta
Definition: mlp_gen.cc:36
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define abs(x)
Definition: mlp_lapack.h:159
#define min(a, b)
Definition: mlp_lapack.h:161
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
Definition: Activities.doc:4
const Double_t pi
T sqrt(T t)
Definition: SSEVec.h:46
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
int k[5][pyjets_maxn]
static const double tmax[3]
tuple mass
Definition: scaleCards.py:27
x
Definition: VDTMath.h:216
Basic3DVector cross(const Basic3DVector &v) const
Vector product, or &quot;cross&quot; product, with a vector of same type.
Definition: DDAxes.h:10