Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4PAIModel.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: G4PAIModel.cc 106217 2017-09-21 00:03:23Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class
31 // File name: G4PAIModel.cc
32 //
33 // Author: Vladimir.Grichine@cern.ch on base of V.Ivanchenko model interface
34 //
35 // Creation date: 05.10.2003
36 //
37 // Modifications:
38 //
39 // 17.08.04 V.Grichine, bug fixed for Tkin<=0 in SampleSecondary
40 // 16.08.04 V.Grichine, bug fixed in massRatio for DEDX, CrossSection,
41 // SampleSecondary
42 // 08.04.05 Major optimisation of internal interfaces (V.Ivantchenko)
43 // 26.07.09 Fixed logic to work with several materials (V.Ivantchenko)
44 // 21.11.10 V. Grichine verbose flag for protons and G4PAYySection to
45 // check sandia table
46 // 12.06.13 V. Grichine Bug fixed in SampleSecondaries for scaled Tkin
47 // (fMass -> proton_mass_c2)
48 // 19.08.13 V.Ivanchenko extract data handling to G4PAIModelData class
49 // added sharing of internal data between threads (MT migration)
50 //
51 
52 #include "G4PAIModel.hh"
53 
54 #include "G4SystemOfUnits.hh"
55 #include "G4PhysicalConstants.hh"
56 #include "G4Region.hh"
57 #include "G4MaterialCutsCouple.hh"
58 #include "G4MaterialTable.hh"
59 #include "G4RegionStore.hh"
60 
61 #include "Randomize.hh"
62 #include "G4Electron.hh"
63 #include "G4Positron.hh"
64 #include "G4Poisson.hh"
65 #include "G4Step.hh"
66 #include "G4Material.hh"
67 #include "G4DynamicParticle.hh"
68 #include "G4ParticleDefinition.hh"
70 #include "G4PAIModelData.hh"
71 #include "G4DeltaAngle.hh"
72 
74 
75 using namespace std;
76 
79  fVerbose(0),
80  fModelData(nullptr),
81  fParticle(nullptr)
82 {
85 
86  fParticleChange = nullptr;
87 
88  if(p) { SetParticle(p); }
89  else { SetParticle(fElectron); }
90 
91  // default generator
93  fLowestTcut = 12.5*CLHEP::eV;
94 }
95 
97 
99 {
100  if(IsMaster()) { delete fModelData; }
101 }
102 
104 
106  const G4DataVector& cuts)
107 {
108  if(fVerbose > 0) {
109  G4cout<<"G4PAIModel::Initialise for "<<p->GetParticleName()<<G4endl;
110  }
111  SetParticle(p);
113 
114  if(IsMaster()) {
115 
116  delete fModelData;
117  fMaterialCutsCoupleVector.clear();
118  if(fVerbose > 0) {
119  G4cout << "G4PAIModel instantiates data for " << p->GetParticleName()
120  << G4endl;
121  }
122  G4double tmin = LowEnergyLimit()*fRatio;
124  fModelData = new G4PAIModelData(tmin, tmax, fVerbose);
125 
126  // Prepare initialization
127  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
128  size_t numOfMat = G4Material::GetNumberOfMaterials();
129  size_t numRegions = fPAIRegionVector.size();
130 
131  // protect for unit tests
132  if(0 == numRegions) {
133  G4Exception("G4PAIModel::Initialise()","em0106",JustWarning,
134  "no G4Regions are registered for the PAI model - World is used");
136  ->GetRegion("DefaultRegionForTheWorld", false));
137  numRegions = 1;
138  }
139 
140  if(fVerbose > 0) {
141  G4cout << "G4PAIModel is defined for " << numRegions << " regions "
142  << G4endl;
143  G4cout << " total number of materials " << numOfMat << G4endl;
144  }
145  for(size_t iReg = 0; iReg<numRegions; ++iReg) {
146  const G4Region* curReg = fPAIRegionVector[iReg];
147  G4Region* reg = const_cast<G4Region*>(curReg);
148 
149  for(size_t jMat = 0; jMat<numOfMat; ++jMat) {
150  G4Material* mat = (*theMaterialTable)[jMat];
151  const G4MaterialCutsCouple* cutCouple = reg->FindCouple(mat);
152  size_t n = fMaterialCutsCoupleVector.size();
153  /*
154  G4cout << "Region: " << reg->GetName() << " " << reg
155  << " Couple " << cutCouple
156  << " PAI defined for " << n << " couples"
157  << " jMat= " << jMat << " " << mat->GetName()
158  << G4endl;
159  */
160  if(cutCouple) {
161  if(fVerbose > 0) {
162  G4cout << "Region <" << curReg->GetName() << "> mat <"
163  << mat->GetName() << "> CoupleIndex= "
164  << cutCouple->GetIndex()
165  << " " << p->GetParticleName()
166  << " cutsize= " << cuts.size() << G4endl;
167  }
168  // check if this couple is not already initialized
169  G4bool isnew = true;
170  if(0 < n) {
171  for(size_t i=0; i<n; ++i) {
172  if(cutCouple == fMaterialCutsCoupleVector[i]) {
173  isnew = false;
174  break;
175  }
176  }
177  }
178  // initialise data banks
179  //G4cout << " isNew: " << isnew << " " << cutCouple << G4endl;
180  if(isnew) {
181  fMaterialCutsCoupleVector.push_back(cutCouple);
182  fModelData->Initialise(cutCouple, this);
183  }
184  }
185  }
186  }
188  }
189 }
190 
192 
194  G4VEmModel* masterModel)
195 {
196  SetParticle(p);
197  fModelData = static_cast<G4PAIModel*>(masterModel)->GetPAIModelData();
199  static_cast<G4PAIModel*>(masterModel)->GetVectorOfCouples();
201 }
202 
204 
206  const G4MaterialCutsCouple*)
207 {
208  return fLowestTcut;
209 }
210 
212 
214  const G4ParticleDefinition* p,
216  G4double cutEnergy)
217 {
218  //G4cout << "===1=== " << CurrentCouple()
219  // << " idx= " << CurrentCouple()->GetIndex()
220  // << " " << fMaterialCutsCoupleVector[0]
221  // << G4endl;
222  G4int coupleIndex = FindCoupleIndex(CurrentCouple());
223  //G4cout << "===2=== " << coupleIndex << G4endl;
224  if(0 > coupleIndex) { return 0.0; }
225 
226  G4double cut = std::min(MaxSecondaryEnergy(p, kineticEnergy), cutEnergy);
227 
228  G4double scaledTkin = kineticEnergy*fRatio;
229 
230  return fChargeSquare*fModelData->DEDXPerVolume(coupleIndex, scaledTkin,
231  cut);
232 }
233 
235 
237  const G4ParticleDefinition* p,
239  G4double cutEnergy,
240  G4double maxEnergy )
241 {
242  //G4cout << "===3=== " << CurrentCouple()
243  // << " idx= " << CurrentCouple()->GetIndex()
244  // << " " << fMaterialCutsCoupleVector[0]
245  // << G4endl;
246  G4int coupleIndex = FindCoupleIndex(CurrentCouple());
247  //G4cout << "===4=== " << coupleIndex << G4endl;
248  if(0 > coupleIndex) { return 0.0; }
249 
250  G4double tmax = std::min(MaxSecondaryEnergy(p, kineticEnergy), maxEnergy);
251  if(tmax <= cutEnergy) { return 0.0; }
252 
253  G4double scaledTkin = kineticEnergy*fRatio;
254 
255  return fChargeSquare*fModelData->CrossSectionPerVolume(coupleIndex,
256  scaledTkin,
257  cutEnergy,
258  tmax);
259 }
260 
262 //
263 // It is analog of PostStepDoIt in terms of secondary electron.
264 //
265 
266 void G4PAIModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
267  const G4MaterialCutsCouple* matCC,
268  const G4DynamicParticle* dp,
269  G4double tmin,
270  G4double maxEnergy)
271 {
272  G4int coupleIndex = FindCoupleIndex(matCC);
273 
274  //G4cout << "G4PAIModel::SampleSecondaries: coupleIndex= "<<coupleIndex<<G4endl;
275  if(0 > coupleIndex) { return; }
276 
277  SetParticle(dp->GetDefinition());
279 
280  G4double tmax = MaxSecondaryEnergy(fParticle, kineticEnergy);
281  if(maxEnergy < tmax) { tmax = maxEnergy; }
282  if(tmin >= tmax) { return; }
283 
284  G4ThreeVector direction= dp->GetMomentumDirection();
285  G4double scaledTkin = kineticEnergy*fRatio;
286  G4double totalEnergy = kineticEnergy + fMass;
287  G4double totalMomentum = sqrt(kineticEnergy*(totalEnergy+fMass));
288 
289  G4double deltaTkin =
290  fModelData->SamplePostStepTransfer(coupleIndex, scaledTkin, tmin, tmax);
291 
292  //G4cout<<"G4PAIModel::SampleSecondaries; deltaKIn = "<<deltaTkin/keV
293  // <<" keV "<< " Escaled(MeV)= " << scaledTkin << G4endl;
294 
295  if( !(deltaTkin <= 0.) && !(deltaTkin > 0)) {
296  G4cout<<"G4PAIModel::SampleSecondaries; deltaKIn = "<<deltaTkin/keV
297  <<" keV "<< " Escaled(MeV)= " << scaledTkin << G4endl;
298  return;
299  }
300  if( deltaTkin <= 0.) { return; }
301 
302  if( deltaTkin > tmax) { deltaTkin = tmax; }
303 
304  const G4Element* anElement = SelectRandomAtom(matCC,fParticle,kineticEnergy);
305  G4int Z = G4lrint(anElement->GetZ());
306 
308  GetAngularDistribution()->SampleDirection(dp, deltaTkin,
309  Z, matCC->GetMaterial()),
310  deltaTkin);
311 
312  // primary change
313  kineticEnergy -= deltaTkin;
314  G4ThreeVector dir = totalMomentum*direction - deltaRay->GetMomentum();
315  direction = dir.unit();
318 
319  vdp->push_back(deltaRay);
320 }
321 
323 
325  const G4DynamicParticle* aParticle,
326  G4double tmax, G4double step,
327  G4double eloss)
328 {
329  G4int coupleIndex = FindCoupleIndex(matCC);
330  if(0 > coupleIndex) { return eloss; }
331 
332  SetParticle(aParticle->GetDefinition());
333 
334  /*
335  G4cout << "G4PAIModel::SampleFluctuations step(mm)= "<< step/mm
336  << " Eloss(keV)= " << eloss/keV << " in "
337  << matCC->Getmaterial()->GetName() << G4endl;
338  */
339 
340  G4double Tkin = aParticle->GetKineticEnergy();
341  G4double scaledTkin = Tkin*fRatio;
342 
343  G4double loss = fModelData->SampleAlongStepTransfer(coupleIndex, Tkin,
344  scaledTkin, tmax,
345  step*fChargeSquare);
346 
347  // G4cout<<"PAIModel AlongStepLoss = "<<loss/keV<<" keV, on step = "
348  //<<step/mm<<" mm"<<G4endl;
349  return loss;
350 
351 }
352 
354 //
355 // Returns the statistical estimation of the energy loss distribution variance
356 //
357 
358 
360  const G4DynamicParticle* aParticle,
361  G4double tmax,
362  G4double step )
363 {
364  G4double particleMass = aParticle->GetMass();
365  G4double electronDensity = material->GetElectronDensity();
366  G4double kineticEnergy = aParticle->GetKineticEnergy();
367  G4double q = aParticle->GetCharge()/eplus;
368  G4double etot = kineticEnergy + particleMass;
369  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*particleMass)/(etot*etot);
370  G4double siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * step
371  * electronDensity * q * q;
372 
373  return siga;
374 }
375 
377 
379  G4double kinEnergy)
380 {
381  SetParticle(p);
382  G4double tmax = kinEnergy;
383  if(p == fElectron) { tmax *= 0.5; }
384  else if(p != fPositron) {
386  G4double gamma= kinEnergy/fMass + 1.0;
387  tmax = 2.0*electron_mass_c2*(gamma*gamma - 1.) /
388  (1. + 2.0*gamma*ratio + ratio*ratio);
389  }
390  return tmax;
391 }
392 
394 
396 {
397  fPAIRegionVector.push_back(r);
398 }
399 
401 
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:792
G4int FindCoupleIndex(const G4MaterialCutsCouple *)
Definition: G4PAIModel.hh:175
G4PAIModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="PAI")
Definition: G4PAIModel.cc:77
void DefineForRegion(const G4Region *r) final
Definition: G4PAIModel.cc:395
G4double fLowestTcut
Definition: G4PAIModel.hh:151
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:600
void SetParticle(const G4ParticleDefinition *p)
Definition: G4PAIModel.hh:188
~G4PAIModel() final
Definition: G4PAIModel.cc:98
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:593
G4double fChargeSquare
Definition: G4PAIModel.hh:150
static constexpr double keV
Definition: G4SIunits.hh:216
static constexpr double twopi_mc2_rcl2
G4int fVerbose
Definition: G4PAIModel.hh:136
const G4String & GetName() const
#define G4endl
Definition: G4ios.hh:61
const G4ThreeVector & GetMomentumDirection() const
const char * p
Definition: xmltok.h:285
const G4String & GetParticleName() const
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:552
G4double SamplePostStepTransfer(G4int coupleIndex, G4double scaledTkin, G4double tmin, G4double tmax) const
void SetProposedMomentumDirection(const G4ThreeVector &dir)
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:784
G4ParticleChangeForLoss * fParticleChange
Definition: G4PAIModel.hh:146
const G4ParticleDefinition * fParticle
Definition: G4PAIModel.hh:143
G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double, G4double) final
Definition: G4PAIModel.cc:359
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:616
const G4String & GetName() const
Definition: G4Material.hh:179
Float_t Z
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:116
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4double fRatio
Definition: G4PAIModel.hh:149
G4double DEDXPerVolume(G4int coupleIndex, G4double scaledTkin, G4double cut) const
G4ParticleDefinition * GetDefinition() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) final
Definition: G4PAIModel.cc:193
void Initialise(const G4ParticleDefinition *, const G4DataVector &) final
Definition: G4PAIModel.cc:105
std::vector< const G4MaterialCutsCouple * > fMaterialCutsCoupleVector
Definition: G4PAIModel.hh:140
static constexpr double electron_mass_c2
void Initialise(const G4MaterialCutsCouple *, G4PAIModel *)
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) final
Definition: G4PAIModel.cc:266
G4double GetCharge() const
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
Definition: G4PAIModel.cc:378
G4double fMass
Definition: G4PAIModel.hh:148
G4double CrossSectionPerVolume(G4int coupleIndex, G4double scaledTkin, G4double tcut, G4double tmax) const
static G4Positron * Positron()
Definition: G4Positron.cc:94
static constexpr double eV
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:582
Float_t mat
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) final
Definition: G4PAIModel.cc:236
static G4Electron * Electron()
Definition: G4Electron.cc:94
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:462
Hep3Vector unit() const
G4double SampleAlongStepTransfer(G4int coupleIndex, G4double kinEnergy, G4double scaledTkin, G4double tmax, G4double stepFactor) const
std::vector< G4Material * > G4MaterialTable
static constexpr double eplus
Definition: G4SIunits.hh:199
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) final
Definition: G4PAIModel.cc:213
int G4lrint(double ad)
Definition: templates.hh:151
static G4RegionStore * GetInstance()
int G4int
Definition: G4Types.hh:78
void SetProposedKineticEnergy(G4double proposedKinEnergy)
TDirectory * dir
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) final
Definition: G4PAIModel.cc:205
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
const G4ParticleDefinition * fElectron
Definition: G4PAIModel.hh:144
G4bool IsMaster() const
Definition: G4VEmModel.hh:700
G4MaterialCutsCouple * FindCouple(G4Material *mat)
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:575
G4double GetZ() const
Definition: G4Element.hh:131
Char_t n[5]
const G4Material * GetMaterial() const
G4double GetMass() const
G4PAIModelData * GetPAIModelData()
Definition: G4PAIModel.hh:154
static const G4double reg
const G4ParticleDefinition * fPositron
Definition: G4PAIModel.hh:145
std::vector< const G4Region * > fPAIRegionVector
Definition: G4PAIModel.hh:141
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:609
G4double GetElectronDensity() const
Definition: G4Material.hh:218
const std::vector< const G4MaterialCutsCouple * > & GetVectorOfCouples()
Definition: G4PAIModel.hh:160
G4PAIModelData * fModelData
Definition: G4PAIModel.hh:138
G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double, G4double, G4double) final
Definition: G4PAIModel.cc:324
T min(const T t1, const T t2)
brief Return the smallest of the two arguments