Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4VEmModel.cc
이 파일의 문서화 페이지로 가기
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 110572 2018-05-30 13:08:12Z 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 "G4ElementData.hh"
54 #include "G4LossTableManager.hh"
55 #include "G4ProductionCutsTable.hh"
58 #include "G4SystemOfUnits.hh"
59 #include "G4Log.hh"
60 #include "Randomize.hh"
61 #include <iostream>
62 
63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65 
67  flucModel(nullptr),anglModel(nullptr), name(nam), lowLimit(0.1*CLHEP::keV),
68  highLimit(100.0*CLHEP::TeV),eMinActive(0.0),eMaxActive(DBL_MAX),
69  polarAngleLimit(CLHEP::pi),secondaryThreshold(DBL_MAX),
70  theLPMflag(false),flagDeexcitation(false),flagForceBuildTable(false),
71  isMaster(true),fElementData(nullptr),pParticleChange(nullptr),
72  xSectionTable(nullptr),theDensityFactor(nullptr),theDensityIdx(nullptr),
73  lossFlucFlag(true),inveplus(1.0/CLHEP::eplus),fCurrentCouple(nullptr),
74  fCurrentElement(nullptr),fCurrentIsotope(nullptr),
75  fTripletModel(nullptr),nsec(5)
76 {
77  xsec.resize(nsec);
78  nSelectors = 0;
79  elmSelectors = nullptr;
80  localElmSelectors = true;
81  localTable = true;
82  useAngularGenerator = false;
83  isLocked = false;
84  idxTable = 0;
85 
87  fEmManager->Register(this);
88 }
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
91 
93 {
94  if(localElmSelectors) {
95  for(G4int i=0; i<nSelectors; ++i) {
96  delete (*elmSelectors)[i];
97  }
98  delete elmSelectors;
99  }
100  delete anglModel;
101 
102  if(localTable && xSectionTable) {
104  delete xSectionTable;
105  xSectionTable = nullptr;
106  }
107  if(isMaster && fElementData) {
108  delete fElementData;
109  fElementData = nullptr;
110  }
111  fEmManager->DeRegister(this);
112 }
113 
114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
115 
117 {
118  G4ParticleChangeForLoss* p = nullptr;
119  if (pParticleChange) {
120  p = static_cast<G4ParticleChangeForLoss*>(pParticleChange);
121  } else {
122  p = new G4ParticleChangeForLoss();
123  pParticleChange = p;
124  }
126  return p;
127 }
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
130 
132 {
133  G4ParticleChangeForGamma* p = nullptr;
134  if (pParticleChange) {
135  p = static_cast<G4ParticleChangeForGamma*>(pParticleChange);
136  } else {
137  p = new G4ParticleChangeForGamma();
138  pParticleChange = p;
139  }
141  return p;
142 }
143 
144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145 
147  const G4DataVector& cuts)
148 {
149  // using spline for element selectors should be investigated in details
150  // because small number of points may provide biased results
151  // large number of points requires significant increase of memory
152  G4bool spline = false;
153 
154  //G4cout << "IES: for " << GetName() << " Emin(MeV)= " << lowLimit/MeV
155  // << " Emax(MeV)= " << highLimit/MeV << G4endl;
156 
157  // two times less bins because probability functon is normalized
158  // so correspondingly is more smooth
159  if(highLimit <= lowLimit) { return; }
160 
162 
163  G4ProductionCutsTable* theCoupleTable=
165  G4int numOfCouples = theCoupleTable->GetTableSize();
166 
167  // prepare vector
168  if(!elmSelectors) {
169  elmSelectors = new std::vector<G4EmElementSelector*>;
170  }
171  if(numOfCouples > nSelectors) {
172  for(G4int i=nSelectors; i<numOfCouples; ++i) {
173  elmSelectors->push_back(nullptr);
174  }
175  nSelectors = numOfCouples;
176  }
177 
178  // initialise vector
179  for(G4int i=0; i<numOfCouples; ++i) {
180 
181  // no need in element selectors for infionite cuts
182  if(cuts[i] == DBL_MAX) { continue; }
183 
184  fCurrentCouple = theCoupleTable->GetMaterialCutsCouple(i);
185  const G4Material* material = fCurrentCouple->GetMaterial();
186 
187  // selector already exist check if should be deleted
188  G4bool create = true;
189  if((*elmSelectors)[i]) {
190  if(material == ((*elmSelectors)[i])->GetMaterial()) { create = false; }
191  else { delete (*elmSelectors)[i]; }
192  }
193  if(create) {
194  G4double emin = std::max(lowLimit,
195  MinPrimaryEnergy(material, part, cuts[i]));
196  G4double emax = std::max(highLimit, 10*emin);
197  static const G4double invlog106 = 1.0/(6*G4Log(10.));
198  G4int nbins = (G4int)(nbinsPerDec*G4Log(emax/emin)*invlog106);
199  nbins = std::max(nbins, 3);
200 
201  (*elmSelectors)[i] = new G4EmElementSelector(this,material,nbins,
202  emin,emax,spline);
203  }
204  ((*elmSelectors)[i])->Initialise(part, cuts[i]);
205  /*
206  G4cout << "G4VEmModel::InitialiseElmSelectors i= " << i
207  << " idx= " << fCurrentCouple->GetIndex()
208  << " " << part->GetParticleName()
209  << " for " << GetName() << " cut= " << cuts[i]
210  << " " << (*elmSelectors)[i] << G4endl;
211  ((*elmSelectors)[i])->Dump(part);
212  */
213  }
214 }
215 
216 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
217 
219  G4VEmModel*)
220 {}
221 
222 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
223 
225  const G4Material* material)
226 {
227  if(material) {
228  size_t n = material->GetNumberOfElements();
229  for(size_t i=0; i<n; ++i) {
230  G4int Z = material->GetElement(i)->GetZasInt();
231  InitialiseForElement(part, Z);
232  }
233  }
234 }
235 
236 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
237 
239 {}
240 
241 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
242 
244  const G4ParticleDefinition*,
246 {
247  return 0.0;
248 }
249 
250 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
251 
253  const G4ParticleDefinition* p,
254  G4double ekin,
255  G4double emin,
256  G4double emax)
257 {
258  SetupForMaterial(p, material, ekin);
259  G4double cross = 0.0;
260  const G4double* theAtomNumDensityVector =
261  material->GetVecNbOfAtomsPerVolume();
262  G4int nelm = material->GetNumberOfElements();
263  if(nelm > nsec) {
264  xsec.resize(nelm);
265  nsec = nelm;
266  }
267  for (G4int i=0; i<nelm; ++i) {
268  cross += theAtomNumDensityVector[i]*
269  ComputeCrossSectionPerAtom(p,material->GetElement(i),ekin,emin,emax);
270  xsec[i] = cross;
271  }
272  return cross;
273 }
274 
275 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
276 
278  const G4ParticleDefinition*,
279  G4double)
280 {
281  return 0.0;
282 }
283 
284 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
285 
287 {}
288 
289 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
290 
292  const G4ParticleDefinition* pd,
293  G4double kinEnergy,
294  G4double tcut,
295  G4double tmax)
296 {
297  size_t n = material->GetNumberOfElements();
298  fCurrentElement = material->GetElement(0);
299  if (n > 1) {
301  G4VEmModel::CrossSectionPerVolume(material,pd,kinEnergy,tcut,tmax);
302  for(size_t i=0; i<n; ++i) {
303  if (x <= xsec[i]) {
304  fCurrentElement = material->GetElement(i);
305  break;
306  }
307  }
308  }
309  return fCurrentElement;
310 }
311 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
312 
314 {
315  // this algorith assumes that cross section is proportional to
316  // number electrons multiplied by number of atoms
317  size_t nn = mat->GetNumberOfElements();
318  fCurrentElement = mat->GetElement(0);
319  if(1 < nn) {
320  const G4double* at = mat->GetVecNbOfAtomsPerVolume();
322  for(size_t i=0; i<nn; ++i) {
323  tot -= at[i];
324  if(tot <= 0.0) {
325  fCurrentElement = mat->GetElement(i);
326  break;
327  }
328  }
329  }
330  return fCurrentElement->GetZasInt();
331 }
332 
333 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
334 
336 {
337  SetCurrentElement(elm);
338  size_t ni = elm->GetNumberOfIsotopes();
339  fCurrentIsotope = elm->GetIsotope(0);
340  size_t idx = 0;
341  if(ni > 1) {
342  const G4double* ab = elm->GetRelativeAbundanceVector();
344  for(; idx<ni; ++idx) {
345  x -= ab[idx];
346  if (x <= 0.0) {
347  fCurrentIsotope = elm->GetIsotope(idx);
348  break;
349  }
350  }
351  }
352  return fCurrentIsotope->GetN();
353 }
354 
355 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
356 
360 {
361  return 0.0;
362 }
363 
364 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
365 
366 G4double
368  G4int, G4int,
370 {
371  return 0.0;
372 }
373 
374 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
375 
377 {}
378 
379 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
380 
382 {
384  track.GetMaterial(), track.GetKineticEnergy());
385 }
386 
387 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
388 
390  const G4Material*, G4double)
391 {
392  G4double q = p->GetPDGCharge()*inveplus;
393  return q*q;
394 }
395 
396 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
397 
399  const G4Material*, G4double)
400 {
401  return p->GetPDGCharge();
402 }
403 
404 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
405 
407  const G4DynamicParticle*,
409 {}
410 
411 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
412 
415 {
416  SetCurrentCouple(couple);
417  return e*e*CrossSectionPerVolume(couple->GetMaterial(),p,e,0.0,DBL_MAX);
418 }
419 
420 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
421 
423  const G4ParticleDefinition*,
424  G4double)
425 {
426  return 0.0;
427 }
428 
429 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
430 
432  const G4MaterialCutsCouple*)
433 {
434  return 0.0;
435 }
436 
437 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
438 
441 {
442  return kineticEnergy;
443 }
444 
445 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
446 
448  const G4Material*, G4double)
449 {}
450 
451 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
452 
453 void
455 {
456  if(p && pParticleChange != p) { pParticleChange = p; }
457  if(flucModel != f) { flucModel = f; }
458 }
459 
460 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
461 
463 {
464  if(p != xSectionTable) {
465  if(xSectionTable && localTable) {
467  delete xSectionTable;
468  }
469  xSectionTable = p;
470  }
471  localTable = isLocal;
472 }
473 
474 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
475 
476 void G4VEmModel::ModelDescription(std::ostream& outFile) const
477 {
478  ModelDescription(outFile, G4String("\n"));
479 }
480 
481 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
482 
483 void G4VEmModel::ModelDescription(std::ostream& outFile,
484  G4String endOfLine) const
485 {
486  outFile << "The description for this model has not been written yet."
487  << endOfLine;
488 }
489 
490 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Float_t x
Definition: compare.C:6
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
Definition: G4VEmModel.cc:462
G4double GetKineticEnergy() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const XML_Char * name
Definition: expat.h:151
virtual G4double GetPartialCrossSection(const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:277
G4double lowLimit
Definition: G4VEmModel.hh:409
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:447
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.cc:313
G4double highLimit
Definition: G4VEmModel.hh:410
G4int NumberOfBinsPerDecade() const
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
G4VEmModel * fTripletModel
Definition: G4VEmModel.hh:446
const G4MaterialCutsCouple * fCurrentCouple
Definition: G4VEmModel.hh:443
const G4Isotope * fCurrentIsotope
Definition: G4VEmModel.hh:445
static constexpr double keV
Definition: G4SIunits.hh:216
G4VEmAngularDistribution * anglModel
Definition: G4VEmModel.hh:404
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:389
G4bool isMaster
Definition: G4VEmModel.hh:418
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.cc:335
const char * p
Definition: xmltok.h:285
virtual ~G4VEmModel()
Definition: G4VEmModel.cc:92
G4double GetPDGCharge() const
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:207
G4ElementData * fElementData
Definition: G4VEmModel.hh:430
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:552
static const G4double emax
const G4double inveplus
const G4ParticleDefinition * GetParticleDefinition() const
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:455
virtual void DefineForRegion(const G4Region *)
Definition: G4VEmModel.cc:376
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:469
G4double G4Log(G4double x)
Definition: G4Log.hh:230
Float_t Z
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:116
static const G4double ab
static constexpr double TeV
Definition: G4SIunits.hh:218
void Register(G4VEnergyLossProcess *p)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
TString part[npart]
Definition: Style.C:32
std::vector< G4double > xsec
Definition: G4VEmModel.hh:449
G4double inveplus
Definition: G4VEmModel.hh:437
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
Definition: G4VEmModel.cc:238
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
G4VEmFluctuationModel * flucModel
Definition: G4VEmModel.hh:403
#define G4UniformRand()
Definition: Randomize.hh:53
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
Definition: G4VEmModel.cc:406
G4bool localElmSelectors
Definition: G4VEmModel.hh:421
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:286
void DeRegister(G4VEnergyLossProcess *p)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:439
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:357
Float_t mat
virtual G4double ComputeCrossSectionPerShell(const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:367
G4Material * GetMaterial() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
virtual void InitialiseForMaterial(const G4ParticleDefinition *, const G4Material *)
Definition: G4VEmModel.cc:224
G4int GetZasInt() const
Definition: G4Element.hh:132
static G4ProductionCutsTable * GetProductionCutsTable()
virtual G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:398
static constexpr double eplus
Definition: G4SIunits.hh:199
G4int GetN() const
Definition: G4Isotope.hh:94
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:252
int G4int
Definition: G4Types.hh:78
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:431
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:454
G4int nSelectors
Definition: G4VEmModel.hh:424
G4bool useAngularGenerator
Definition: G4VEmModel.hh:422
G4LossTableManager * fEmManager
Definition: G4VEmModel.hh:426
static G4LossTableManager * Instance()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
void clearAndDestroy()
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
const G4Element * fCurrentElement
Definition: G4VEmModel.hh:444
G4int nsec
Definition: G4VEmModel.hh:448
Char_t n[5]
const G4Material * GetMaterial() const
static constexpr double pi
Definition: G4SIunits.hh:75
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:131
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
Definition: G4VEmModel.cc:422
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:210
virtual void ModelDescription(std::ostream &outFile) const
Definition: G4VEmModel.cc:476
virtual G4double ChargeSquareRatio(const G4Track &)
Definition: G4VEmModel.cc:381
G4bool localTable
Definition: G4VEmModel.hh:420
size_t idxTable
Definition: G4VEmModel.hh:435
G4bool isLocked
Definition: G4VEmModel.hh:423
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:432
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:66
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:243
#define DBL_MAX
Definition: templates.hh:83
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:203
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
Definition: G4VEmModel.cc:431
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:413
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:218
static G4EmParameters * Instance()
size_t GetNumberOfElements() const
Definition: G4Material.hh:187
std::vector< G4EmElementSelector * > * elmSelectors
Definition: G4VEmModel.hh:425
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159