Geant4  10.00.p01
G4VEmModel.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: G4VEmModel.cc 76333 2013-11-08 14:31:50Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4VEmModel
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 25.07.2005
38 //
39 // Modifications:
40 // 25.10.2005 Set default highLimit=100.TeV (V.Ivanchenko)
41 // 06.02.2006 add method ComputeMeanFreePath() (mma)
42 // 16.02.2009 Move implementations of virtual methods to source (VI)
43 //
44 //
45 // Class Description:
46 //
47 // Abstract interface to energy loss models
48 
49 // -------------------------------------------------------------------
50 //
51 
52 #include "G4VEmModel.hh"
53 #include "G4LossTableManager.hh"
54 #include "G4ProductionCutsTable.hh"
57 #include "G4SystemOfUnits.hh"
58 #include "G4Log.hh"
59 
60 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62 
63 const G4double log106 = 6*G4Log(10.);
64 
66  flucModel(0),anglModel(0), name(nam), lowLimit(0.1*CLHEP::keV),
67  highLimit(100.0*CLHEP::TeV),eMinActive(0.0),eMaxActive(DBL_MAX),
68  polarAngleLimit(CLHEP::pi),secondaryThreshold(DBL_MAX),
69  theLPMflag(false),flagDeexcitation(false),flagForceBuildTable(false),
70  isMaster(true),fElementData(0),pParticleChange(0),xSectionTable(0),
71  theDensityFactor(0),theDensityIdx(0),fCurrentCouple(0),fCurrentElement(0),
72  nsec(5)
73 {
74  xsec.resize(nsec);
75  nSelectors = 0;
76  elmSelectors = 0;
77  localElmSelectors = true;
78  localTable = true;
79  useAngularGenerator = false;
80  idxTable = 0;
81 
83  fManager->Register(this);
84 }
85 
86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87 
89 {
90  if(localElmSelectors) {
91  if(nSelectors > 0) {
92  for(G4int i=0; i<nSelectors; ++i) {
93  delete (*elmSelectors)[i];
94  }
95  }
96  delete elmSelectors;
97  }
98  delete anglModel;
99 
100  if(localTable) { delete xSectionTable; }
101 
102  fManager->DeRegister(this);
103 }
104 
105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106 
108 {
110  if (pParticleChange) {
111  p = static_cast<G4ParticleChangeForLoss*>(pParticleChange);
112  } else {
113  p = new G4ParticleChangeForLoss();
114  pParticleChange = p;
115  }
116  return p;
117 }
118 
119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120 
122 {
124  if (pParticleChange) {
125  p = static_cast<G4ParticleChangeForGamma*>(pParticleChange);
126  } else {
127  p = new G4ParticleChangeForGamma();
128  pParticleChange = p;
129  }
130  return p;
131 }
132 
133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134 
136  const G4DataVector& cuts)
137 {
138  // using spline for element selectors should be investigated in details
139  // because small number of points may provide biased results
140  // large number of points requires significant increase of memory
141  //G4bool spline = fManager->SplineFlag();
142  G4bool spline = false;
143 
144  //G4cout << "IES: for " << GetName() << " Emin(MeV)= " << lowLimit/MeV
145  // << " Emax(MeV)= " << highLimit/MeV << G4endl;
146 
147  // two times less bins because probability functon is normalized
148  // so correspondingly is more smooth
149  if(highLimit <= lowLimit) { return; }
150 
151  G4ProductionCutsTable* theCoupleTable=
153  G4int numOfCouples = theCoupleTable->GetTableSize();
154 
155  // prepare vector
156  if(!elmSelectors) {
157  elmSelectors = new std::vector<G4EmElementSelector*>;
158  }
159  if(numOfCouples > nSelectors) {
160  for(G4int i=nSelectors; i<numOfCouples; ++i) {
161  elmSelectors->push_back(0);
162  }
163  nSelectors = numOfCouples;
164  }
165 
166  // initialise vector
167  for(G4int i=0; i<numOfCouples; ++i) {
168 
169  // no need in element selectors for infionite cuts
170  if(cuts[i] == DBL_MAX) { continue; }
171 
172  fCurrentCouple = theCoupleTable->GetMaterialCutsCouple(i);
173  const G4Material* material = fCurrentCouple->GetMaterial();
174 
175  // selector already exist check if should be deleted
176  G4bool create = true;
177  if((*elmSelectors)[i]) {
178  if(material == ((*elmSelectors)[i])->GetMaterial()) { create = false; }
179  else { delete (*elmSelectors)[i]; }
180  }
181  if(create) {
182  G4double emin = std::max(lowLimit,
183  MinPrimaryEnergy(material, part, cuts[i]));
184  G4double emax = std::max(highLimit, 10*emin);
186  *G4Log(emax/emin)/log106);
187  if(nbins < 3) { nbins = 3; }
188 
189  (*elmSelectors)[i] = new G4EmElementSelector(this,material,nbins,
190  emin,emax,spline);
191  }
192  ((*elmSelectors)[i])->Initialise(part, cuts[i]);
193  /*
194  G4cout << "G4VEmModel::InitialiseElmSelectors i= " << i
195  << " idx= " << fCurrentCouple->GetIndex()
196  << " " << part->GetParticleName()
197  << " for " << GetName() << " cut= " << cuts[i]
198  << " " << (*elmSelectors)[i] << G4endl;
199  ((*elmSelectors)[i])->Dump(part);
200  */
201  }
202 }
203 
204 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
205 
207  G4VEmModel*)
208 {}
209 
210 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
211 
213  const G4Material* material)
214 {
215  if(material) {
216  const G4ElementVector* theElementVector = material->GetElementVector();
217  G4int n = material->GetNumberOfElements();
218  for(G4int i=0; i<n; ++i) {
219  G4int Z = G4lrint(((*theElementVector)[i])->GetZ());
220  InitialiseForElement(part, Z);
221  }
222  } else {
223  //G4cout << "G4VEmModel::InitialiseForMaterial for " << GetName();
224  //if(part) { G4cout << " and " << part->GetParticleName(); }
225  //G4cout << " with no material" << G4endl;
226  }
227 }
228 
229 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
230 
232 {}
233 
234 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
235 
237  const G4ParticleDefinition*,
239 {
240  return 0.0;
241 }
242 
243 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
244 
246  const G4ParticleDefinition* p,
247  G4double ekin,
248  G4double emin,
249  G4double emax)
250 {
251  SetupForMaterial(p, material, ekin);
252  G4double cross = 0.0;
253  const G4ElementVector* theElementVector = material->GetElementVector();
254  const G4double* theAtomNumDensityVector =
255  material->GetVecNbOfAtomsPerVolume();
256  G4int nelm = material->GetNumberOfElements();
257  if(nelm > nsec) {
258  xsec.resize(nelm);
259  nsec = nelm;
260  }
261  for (G4int i=0; i<nelm; i++) {
262  cross += theAtomNumDensityVector[i]*
263  ComputeCrossSectionPerAtom(p,(*theElementVector)[i],ekin,emin,emax);
264  xsec[i] = cross;
265  }
266  return cross;
267 }
268 
269 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
270 
272 {}
273 
274 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
275 
277  const G4ParticleDefinition* pd,
278  G4double kinEnergy,
279  G4double tcut,
280  G4double tmax)
281 {
282  const G4ElementVector* theElementVector = material->GetElementVector();
283  G4int n = material->GetNumberOfElements() - 1;
284  fCurrentElement = (*theElementVector)[n];
285  if (n > 0) {
286  G4double x = G4UniformRand()*
287  G4VEmModel::CrossSectionPerVolume(material,pd,kinEnergy,tcut,tmax);
288  for(G4int i=0; i<n; ++i) {
289  if (x <= xsec[i]) {
290  fCurrentElement = (*theElementVector)[i];
291  break;
292  }
293  }
294  }
295  return fCurrentElement;
296 }
297 
298 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
299 
303 {
304  return 0.0;
305 }
306 
307 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
308 
310 {}
311 
312 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
313 
315 {
317  track.GetMaterial(), track.GetKineticEnergy());
318 }
319 
320 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
321 
323  const G4Material*, G4double)
324 {
326  return q*q;
327 }
328 
329 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
330 
332  const G4Material*, G4double)
333 {
334  return p->GetPDGCharge();
335 }
336 
337 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
338 
340  const G4DynamicParticle*,
342 {}
343 
344 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
345 
347  const G4ParticleDefinition* p, G4double e)
348 {
349  SetCurrentCouple(couple);
350  return e*e*CrossSectionPerVolume(couple->GetMaterial(),p,e,0.0,DBL_MAX);
351 }
352 
353 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
354 
356  const G4ParticleDefinition*,
357  G4double)
358 {
359  return 0.0;
360 }
361 
362 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
363 
365  const G4MaterialCutsCouple*)
366 {
367  return 0.0;
368 }
369 
370 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
371 
373  G4double kineticEnergy)
374 {
375  return kineticEnergy;
376 }
377 
378 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
379 
381  const G4Material*, G4double)
382 {}
383 
384 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
385 
386 void
388 {
389  if(p && pParticleChange != p) { pParticleChange = p; }
390  flucModel = f;
391 }
392 
393 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
394 
396 {
397  if(p != xSectionTable) {
398  if(xSectionTable && localTable) {
400  delete xSectionTable;
401  }
402  xSectionTable = p;
403  }
404  localTable = isLocal;
405 }
406 
407 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4bool localTable
Definition: G4VEmModel.hh:389
size_t idxTable
Definition: G4VEmModel.hh:402
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:245
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
Definition: G4VEmModel.cc:231
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:271
static G4LossTableManager * Instance()
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
std::vector< G4Element * > G4ElementVector
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
Definition: G4VEmModel.cc:339
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
void DeRegister(G4VEnergyLossProcess *p)
virtual void DefineForRegion(const G4Region *)
Definition: G4VEmModel.cc:309
G4String name
Definition: TRTMaterials.hh:40
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
Definition: G4VEmModel.cc:355
const G4double pi
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:206
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:380
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
Definition: G4VEmModel.cc:364
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
const G4Element * fCurrentElement
Definition: G4VEmModel.hh:410
G4VEmFluctuationModel * flucModel
Definition: G4VEmModel.hh:372
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
G4int nSelectors
Definition: G4VEmModel.hh:392
G4LossTableManager * fManager
Definition: G4VEmModel.hh:408
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
G4double GetKineticEnergy() const
virtual void InitialiseForMaterial(const G4ParticleDefinition *, const G4Material *)
Definition: G4VEmModel.cc:212
#define G4UniformRand()
Definition: Randomize.hh:87
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=0)
Definition: G4VEmModel.cc:387
const G4MaterialCutsCouple * fCurrentCouple
Definition: G4VEmModel.hh:409
std::vector< G4EmElementSelector * > * elmSelectors
Definition: G4VEmModel.hh:393
bool G4bool
Definition: G4Types.hh:79
const G4ParticleDefinition * GetParticleDefinition() const
virtual ~G4VEmModel()
Definition: G4VEmModel.cc:88
G4VEmAngularDistribution * anglModel
Definition: G4VEmModel.hh:373
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
Definition: G4VEmModel.cc:395
void Register(G4VEnergyLossProcess *p)
const G4int n
G4int nsec
Definition: G4VEmModel.hh:412
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:300
virtual G4double ChargeSquareRatio(const G4Track &)
Definition: G4VEmModel.cc:314
G4Material * GetMaterial() const
G4double G4Log(G4double x)
Definition: G4Log.hh:227
std::vector< G4double > xsec
Definition: G4VEmModel.hh:413
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:322
virtual G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:331
static G4ProductionCutsTable * GetProductionCutsTable()
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:236
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
int G4lrint(double ad)
Definition: templates.hh:163
G4bool localElmSelectors
Definition: G4VEmModel.hh:390
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4int GetNumberOfBinsPerDecade() const
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:419
G4double highLimit
Definition: G4VEmModel.hh:379
G4bool useAngularGenerator
Definition: G4VEmModel.hh:391
G4double lowLimit
Definition: G4VEmModel.hh:378
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:398
const G4double log106
Definition: G4VEmModel.cc:63
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:372
static const double TeV
Definition: G4SIunits.hh:197
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
static const double keV
Definition: G4SIunits.hh:195
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:399
double G4double
Definition: G4Types.hh:76
static const double eplus
Definition: G4SIunits.hh:178
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:510
void clearAndDestroy()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:346
const G4Material * GetMaterial() const