Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4MuBremsstrahlungModel.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$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32//
33// File name: G4MuBremsstrahlungModel
34//
35// Author: Vladimir Ivanchenko on base of Laszlo Urban code
36//
37// Creation date: 24.06.2002
38//
39// Modifications:
40//
41// 04-12-02 Change G4DynamicParticle constructor in PostStepDoIt (V.Ivanchenko)
42// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
43// 24-01-03 Fix for compounds (V.Ivanchenko)
44// 27-01-03 Make models region aware (V.Ivanchenko)
45// 13-02-03 Add name (V.Ivanchenko)
46// 10-02-04 Add lowestKinEnergy (V.Ivanchenko)
47// 08-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
48// 03-08-05 Angular correlations according to PRM (V.Ivanchenko)
49// 13-02-06 add ComputeCrossSectionPerAtom (mma)
50// 21-03-06 Fix problem of initialisation in case when cuts are not defined (VI)
51// 07-11-07 Improve sampling of final state (A.Bogdanov)
52// 28-02-08 Use precomputed Z^1/3 and Log(A) (V.Ivanchenko)
53//
54
55//
56// Class Description:
57//
58//
59// -------------------------------------------------------------------
60//
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
63
66#include "G4SystemOfUnits.hh"
67#include "G4Gamma.hh"
68#include "G4MuonMinus.hh"
69#include "G4MuonPlus.hh"
70#include "Randomize.hh"
71#include "G4Material.hh"
72#include "G4Element.hh"
73#include "G4ElementVector.hh"
76
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
79
80using namespace std;
81
83 const G4String& nam)
84 : G4VEmModel(nam),
85 particle(0),
86 sqrte(sqrt(exp(1.))),
87 bh(202.4),
88 bh1(446.),
89 btf(183.),
90 btf1(1429.),
91 fParticleChange(0),
92 lowestKinEnergy(1.0*GeV),
93 minThreshold(1.0*keV)
94{
95 theGamma = G4Gamma::Gamma();
97
98 mass = rmass = cc = coeff = 1.0;
99
100 fDN[0] = 0.0;
101 for(G4int i=1; i<93; ++i) {
102 G4double dn = 1.54*nist->GetA27(i);
103 fDN[i] = dn;
104 if(1 < i) {
105 fDN[i] /= std::pow(dn, 1./G4double(i));
106 }
107 }
108
109 if(p) { SetParticle(p); }
110}
111
112//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
113
115{
116 size_t n = partialSumSigma.size();
117 if(n > 0) {
118 for(size_t i=0; i<n; i++) {
119 delete partialSumSigma[i];
120 }
121 }
122}
123
124//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125
128{
129 return minThreshold;
130}
131
132//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
133
135 const G4DataVector& cuts)
136{
137 if(p) { SetParticle(p); }
138
139 // partial cross section is computed for fixed energy
140 G4double fixedEnergy = 0.5*HighEnergyLimit();
141
142 const G4ProductionCutsTable* theCoupleTable=
144 if(theCoupleTable) {
145 G4int numOfCouples = theCoupleTable->GetTableSize();
146
147 G4int nn = partialSumSigma.size();
148 G4int nc = cuts.size();
149
150 // do we need to perform initialisation?
151 if(nn == numOfCouples) { return; }
152
153 // clear old data
154 if(nn > 0) {
155 for (G4int ii=0; ii<nn; ii++){
156 G4DataVector* a = partialSumSigma[ii];
157 if ( a ) { delete a; }
158 }
159 partialSumSigma.clear();
160 }
161 // fill new data
162 if (numOfCouples>0) {
163 for (G4int i=0; i<numOfCouples; i++) {
164 G4double cute = DBL_MAX;
165
166 // protection for usage with extrapolator
167 if(i < nc) { cute = cuts[i]; }
168
169 const G4MaterialCutsCouple* couple =
170 theCoupleTable->GetMaterialCutsCouple(i);
171 const G4Material* material = couple->GetMaterial();
172 G4DataVector* dv = ComputePartialSumSigma(material,fixedEnergy,cute);
173 partialSumSigma.push_back(dv);
174 }
175 }
176 }
177
178 // define pointer to G4ParticleChange
179 if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
180}
181
182//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
183
185 const G4Material* material,
187 G4double kineticEnergy,
188 G4double cutEnergy)
189{
190 G4double dedx = 0.0;
191 if (kineticEnergy <= lowestKinEnergy) return dedx;
192
193 G4double tmax = kineticEnergy;
194 G4double cut = std::min(cutEnergy,tmax);
195 if(cut < minThreshold) cut = minThreshold;
196
197 const G4ElementVector* theElementVector = material->GetElementVector();
198 const G4double* theAtomicNumDensityVector =
199 material->GetAtomicNumDensityVector();
200
201 // loop for elements in the material
202 for (size_t i=0; i<material->GetNumberOfElements(); i++) {
203
204 G4double loss =
205 ComputMuBremLoss((*theElementVector)[i]->GetZ(), kineticEnergy, cut);
206
207 dedx += loss*theAtomicNumDensityVector[i];
208 }
209 // G4cout << "BR e= " << kineticEnergy << " dedx= " << dedx << G4endl;
210 if(dedx < 0.) dedx = 0.;
211 return dedx;
212}
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
215
217 G4double tkin, G4double cut)
218{
219 G4double totalEnergy = mass + tkin;
220 G4double ak1 = 0.05;
221 G4int k2=5;
222 G4double xgi[]={0.03377,0.16940,0.38069,0.61931,0.83060,0.96623};
223 G4double wgi[]={0.08566,0.18038,0.23396,0.23396,0.18038,0.08566};
224 G4double loss = 0.;
225
226 G4double vcut = cut/totalEnergy;
227 G4double vmax = tkin/totalEnergy;
228
229 G4double aaa = 0.;
230 G4double bbb = vcut;
231 if(vcut>vmax) bbb=vmax ;
232 G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ;
233 G4double hhh=(bbb-aaa)/float(kkk) ;
234
235 G4double aa = aaa;
236 for(G4int l=0; l<kkk; l++)
237 {
238 for(G4int i=0; i<6; i++)
239 {
240 G4double ep = (aa + xgi[i]*hhh)*totalEnergy;
241 loss += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
242 }
243 aa += hhh;
244 }
245
246 loss *=hhh*totalEnergy ;
247
248 return loss;
249}
250
251//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
252
254 G4double tkin,
255 G4double Z,
256 G4double cut)
257{
258 G4double totalEnergy = tkin + mass;
259 G4double ak1 = 2.3;
260 G4int k2 = 4;
261 G4double xgi[]={0.03377,0.16940,0.38069,0.61931,0.83060,0.96623};
262 G4double wgi[]={0.08566,0.18038,0.23396,0.23396,0.18038,0.08566};
263 G4double cross = 0.;
264
265 if(cut >= tkin) return cross;
266
267 G4double vcut = cut/totalEnergy;
268 G4double vmax = tkin/totalEnergy;
269
270 G4double aaa = log(vcut);
271 G4double bbb = log(vmax);
272 G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ;
273 G4double hhh = (bbb-aaa)/G4double(kkk);
274
275 G4double aa = aaa;
276
277 for(G4int l=0; l<kkk; l++)
278 {
279 for(G4int i=0; i<6; i++)
280 {
281 G4double ep = exp(aa + xgi[i]*hhh)*totalEnergy;
282 cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
283 }
284 aa += hhh;
285 }
286
287 cross *=hhh;
288
289 //G4cout << "BR e= " << tkin<< " cross= " << cross/barn << G4endl;
290
291 return cross;
292}
293
294//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
295
297 G4double tkin,
298 G4double Z,
299 G4double gammaEnergy)
300// differential cross section
301{
302 G4double dxsection = 0.;
303
304 if( gammaEnergy > tkin) return dxsection ;
305
306 G4double E = tkin + mass ;
307 G4double v = gammaEnergy/E ;
308 G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ;
309 G4double rab0=delta*sqrte ;
310
311 G4int iz = G4int(Z);
312 if(iz < 1) iz = 1;
313 else if(iz > 92) iz = 92;
314
315 G4double z13 = 1.0/nist->GetZ13(iz);
316 G4double dnstar = fDN[iz];
317
318 G4double b,b1;
319
320 if(1 == iz) {
321 b = bh;
322 b1 = bh1;
323 } else {
324 b = btf;
325 b1 = btf1;
326 }
327
328 // nucleus contribution logarithm
329 G4double rab1=b*z13;
330 G4double fn=log(rab1/(dnstar*(electron_mass_c2+rab0*rab1))*
331 (mass+delta*(dnstar*sqrte-2.))) ;
332 if(fn <0.) fn = 0. ;
333 // electron contribution logarithm
334 G4double epmax1=E/(1.+0.5*mass*rmass/E) ;
335 G4double fe=0.;
336 if(gammaEnergy<epmax1)
337 {
338 G4double rab2=b1*z13*z13 ;
339 fe=log(rab2*mass/((1.+delta*rmass/(electron_mass_c2*sqrte))*
340 (electron_mass_c2+rab0*rab2))) ;
341 if(fe<0.) fe=0. ;
342 }
343
344 dxsection = coeff*(1.-v*(1. - 0.75*v))*Z*(fn*Z + fe)/gammaEnergy;
345
346 return dxsection;
347}
348
349//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
350
353 G4double kineticEnergy,
355 G4double cutEnergy,
356 G4double maxEnergy)
357{
358 G4double cross = 0.0;
359 if (kineticEnergy <= lowestKinEnergy) return cross;
360 G4double tmax = std::min(maxEnergy, kineticEnergy);
361 G4double cut = std::min(cutEnergy, kineticEnergy);
362 if(cut < minThreshold) cut = minThreshold;
363 if (cut >= tmax) return cross;
364
365 cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut);
366 if(tmax < kineticEnergy) {
367 cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
368 }
369 return cross;
370}
371
372//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
373
374G4DataVector* G4MuBremsstrahlungModel::ComputePartialSumSigma(
375 const G4Material* material,
376 G4double kineticEnergy,
377 G4double cut)
378
379// Build the table of cross section per element.
380// The table is built for material
381// This table is used to select randomly an element in the material.
382{
383 G4int nElements = material->GetNumberOfElements();
384 const G4ElementVector* theElementVector = material->GetElementVector();
385 const G4double* theAtomNumDensityVector =
386 material->GetAtomicNumDensityVector();
387
388 G4DataVector* dv = new G4DataVector();
389
390 G4double cross = 0.0;
391
392 for (G4int i=0; i<nElements; i++ ) {
393 cross += theAtomNumDensityVector[i]
394 * ComputeMicroscopicCrossSection(kineticEnergy,
395 (*theElementVector)[i]->GetZ(), cut);
396 dv->push_back(cross);
397 }
398 return dv;
399}
400
401//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
402
404 std::vector<G4DynamicParticle*>* vdp,
405 const G4MaterialCutsCouple* couple,
406 const G4DynamicParticle* dp,
407 G4double minEnergy,
408 G4double maxEnergy)
409{
410 G4double kineticEnergy = dp->GetKineticEnergy();
411 // check against insufficient energy
412 G4double tmax = std::min(kineticEnergy, maxEnergy);
413 G4double tmin = std::min(kineticEnergy, minEnergy);
414 if(tmin < minThreshold) tmin = minThreshold;
415 if(tmin >= tmax) return;
416
417 // ===== sampling of energy transfer ======
418
419 G4ParticleMomentum partDirection = dp->GetMomentumDirection();
420
421 // select randomly one element constituing the material
422 const G4Element* anElement = SelectRandomAtom(couple);
423 G4double Z = anElement->GetZ();
424
425 G4double totalEnergy = kineticEnergy + mass;
426 G4double totalMomentum = sqrt(kineticEnergy*(kineticEnergy + 2.0*mass));
427
428 G4double func1 = tmin*
429 ComputeDMicroscopicCrossSection(kineticEnergy,Z,tmin);
430
431 G4double lnepksi, epksi;
432 G4double func2;
433
434 G4double xmin = log(tmin/MeV);
435 G4double xmax = log(kineticEnergy/tmin);
436
437 do {
438 lnepksi = xmin + G4UniformRand()*xmax;
439 epksi = MeV*exp(lnepksi);
440 func2 = epksi*ComputeDMicroscopicCrossSection(kineticEnergy,Z,epksi);
441
442 } while(func2 < func1*G4UniformRand());
443
444 G4double gEnergy = epksi;
445
446 // ===== sample angle =====
447
448 G4double gam = totalEnergy/mass;
449 G4double rmax = gam*std::min(1.0, totalEnergy/gEnergy - 1.0);
450 G4double rmax2= rmax*rmax;
451 G4double x = G4UniformRand()*rmax2/(1.0 + rmax2);
452
453 G4double theta = sqrt(x/(1.0 - x))/gam;
454 G4double sint = sin(theta);
455 G4double phi = twopi * G4UniformRand() ;
456 G4double dirx = sint*cos(phi), diry = sint*sin(phi), dirz = cos(theta) ;
457
458 G4ThreeVector gDirection(dirx, diry, dirz);
459 gDirection.rotateUz(partDirection);
460
461 partDirection *= totalMomentum;
462 partDirection -= gEnergy*gDirection;
463 partDirection = partDirection.unit();
464
465 // primary change
466 kineticEnergy -= gEnergy;
467 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
468 fParticleChange->SetProposedMomentumDirection(partDirection);
469
470 // save secondary
471 G4DynamicParticle* aGamma =
472 new G4DynamicParticle(theGamma,gDirection,gEnergy);
473 vdp->push_back(aGamma);
474}
475
476//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
477
478const G4Element* G4MuBremsstrahlungModel::SelectRandomAtom(
479 const G4MaterialCutsCouple* couple) const
480{
481 // select randomly 1 element within the material
482
483 const G4Material* material = couple->GetMaterial();
484 G4int nElements = material->GetNumberOfElements();
485 const G4ElementVector* theElementVector = material->GetElementVector();
486 if(1 == nElements) { return (*theElementVector)[0]; }
487 else if(1 > nElements) { return 0; }
488
489 G4DataVector* dv = partialSumSigma[couple->GetIndex()];
490 G4double rval = G4UniformRand()*((*dv)[nElements-1]);
491 for (G4int i=0; i<nElements; i++) {
492 if (rval <= (*dv)[i]) { return (*theElementVector)[i]; }
493 }
494 return (*theElementVector)[nElements-1];
495}
496
497//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
std::vector< G4Element * > G4ElementVector
G4fissionEvent * fe
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:215
G4MuBremsstrahlungModel(const G4ParticleDefinition *p=0, const G4String &nam="MuBrem")
void SetParticle(const G4ParticleDefinition *)
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double gammaEnergy)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
G4double ComputMuBremLoss(G4double Z, G4double tkin, G4double cut)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4double GetA27(G4int Z)
G4double GetZ13(G4double Z)
static G4NistManager * Instance()
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:95
#define DBL_MAX
Definition: templates.hh:83