Geant4  10.01.p01
G4PAIModel.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: G4PAIModel.cc 83400 2014-08-21 14:48:50Z 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 "G4PhysicsLogVector.hh"
58 #include "G4PhysicsFreeVector.hh"
59 #include "G4PhysicsTable.hh"
60 #include "G4ProductionCutsTable.hh"
61 #include "G4MaterialCutsCouple.hh"
62 #include "G4MaterialTable.hh"
63 #include "G4SandiaTable.hh"
64 #include "G4OrderedTable.hh"
65 
66 #include "Randomize.hh"
67 #include "G4Electron.hh"
68 #include "G4Positron.hh"
69 #include "G4Poisson.hh"
70 #include "G4Step.hh"
71 #include "G4Material.hh"
72 #include "G4DynamicParticle.hh"
73 #include "G4ParticleDefinition.hh"
75 #include "G4PAIModelData.hh"
76 #include "G4DeltaAngle.hh"
77 
79 
80 using namespace std;
81 
84  fVerbose(0),
85  fModelData(0),
86  fParticle(0)
87 {
90 
91  fParticleChange = 0;
92 
93  if(p) { SetParticle(p); }
94  else { SetParticle(fElectron); }
95 
96  // default generator
98 
99  isInitialised = false;
100 }
101 
103 
105 {
106  if(IsMaster()) { delete fModelData; }
107 }
108 
110 
112  const G4DataVector& cuts)
113 {
114  if(fVerbose > 0) {
115  G4cout<<"G4PAIModel::Initialise for "<<p->GetParticleName()<<G4endl;
116  }
117 
118  if(isInitialised) { return; }
119  isInitialised = true;
120 
121  SetParticle(p);
123 
124  if(IsMaster()) {
125 
127 
128  if(!fModelData) {
129 
130  G4double tmin = LowEnergyLimit()*fRatio;
132  fModelData = new G4PAIModelData(tmin, tmax, fVerbose);
133  }
134  // Prepare initialization
135  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
136  size_t numOfMat = G4Material::GetNumberOfMaterials();
137  size_t numRegions = fPAIRegionVector.size();
138 
139  if(fVerbose > 0) {
140  G4cout << "G4PAIModel is defined for " << numRegions << " regions "
141  << G4endl;
142  G4cout << " total number of materials " << numOfMat << G4endl;
143  }
144  for(size_t iReg = 0; iReg<numRegions; ++iReg) {
145  const G4Region* curReg = fPAIRegionVector[iReg];
146  G4Region* reg = const_cast<G4Region*>(curReg);
147 
148  for(size_t jMat = 0; jMat<numOfMat; ++jMat) {
149  G4Material* mat = (*theMaterialTable)[jMat];
150  const G4MaterialCutsCouple* cutCouple = reg->FindCouple(mat);
151  size_t n = fMaterialCutsCoupleVector.size();
152  /*
153  G4cout << "Region: " << reg->GetName() << " " << reg
154  << " Couple " << cutCouple
155  << " PAI defined for " << n << " couples"
156  << " jMat= " << jMat << " " << mat->GetName()
157  << G4endl;
158  */
159  if(cutCouple) {
160  if(fVerbose > 0) {
161  G4cout << "Region <" << curReg->GetName() << "> mat <"
162  << mat->GetName() << "> CoupleIndex= "
163  << cutCouple->GetIndex()
164  << " " << p->GetParticleName()
165  << " cutsize= " << cuts.size() << G4endl;
166  }
167  // check if this couple is not already initialized
168  G4bool isnew = true;
169  if(0 < n) {
170  for(size_t i=0; i<n; ++i) {
171  if(cutCouple == fMaterialCutsCoupleVector[i]) {
172  isnew = false;
173  break;
174  }
175  }
176  }
177  // initialise data banks
178  //G4cout << " isNew: " << isnew << " " << cutCouple << G4endl;
179  if(isnew) {
180  fMaterialCutsCoupleVector.push_back(cutCouple);
181  G4double deltaCutInKinEnergy = cuts[cutCouple->GetIndex()];
182  fModelData->Initialise(cutCouple, deltaCutInKinEnergy, this);
183  }
184  }
185  }
186  }
187  }
188 }
189 
191 
193  G4VEmModel* masterModel)
194 {
195  fModelData = static_cast<G4PAIModel*>(masterModel)->GetPAIModelData();
197 }
198 
200 
202  const G4ParticleDefinition* p,
203  G4double kineticEnergy,
204  G4double cutEnergy)
205 {
206  /*
207  G4cout << "===1=== " << CurrentCouple()
208  << " idx= " << CurrentCouple()->GetIndex()
209  << " " << fMaterialCutsCoupleVector[0]
210  << G4endl;
211  */
212  G4int coupleIndex = FindCoupleIndex(CurrentCouple());
213  //G4cout << "===2=== " << coupleIndex << G4endl;
214  if(0 > coupleIndex) { return 0.0; }
215 
216  G4double cut = std::min(MaxSecondaryEnergy(p, kineticEnergy), cutEnergy);
217 
218  G4double scaledTkin = kineticEnergy*fRatio;
219 
220  return fChargeSquare*fModelData->DEDXPerVolume(coupleIndex, scaledTkin, cut);
221 }
222 
224 
226  const G4ParticleDefinition* p,
227  G4double kineticEnergy,
228  G4double cutEnergy,
229  G4double maxEnergy )
230 {
231  G4int coupleIndex = FindCoupleIndex(CurrentCouple());
232  if(0 > coupleIndex) { return 0.0; }
233 
234  G4double tmax = std::min(MaxSecondaryEnergy(p, kineticEnergy), maxEnergy);
235  if(tmax <= cutEnergy) { return 0.0; }
236 
237  G4double scaledTkin = kineticEnergy*fRatio;
238 
239  return fChargeSquare*fModelData->CrossSectionPerVolume(coupleIndex,
240  scaledTkin,
241  cutEnergy, tmax);
242 }
243 
245 //
246 // It is analog of PostStepDoIt in terms of secondary electron.
247 //
248 
249 void G4PAIModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
250  const G4MaterialCutsCouple* matCC,
251  const G4DynamicParticle* dp,
252  G4double tmin,
253  G4double maxEnergy)
254 {
255  G4int coupleIndex = FindCoupleIndex(matCC);
256  if(0 > coupleIndex) { return; }
257 
258  SetParticle(dp->GetDefinition());
259  G4double kineticEnergy = dp->GetKineticEnergy();
260 
261  G4double tmax = MaxSecondaryEnergy(fParticle, kineticEnergy);
262  if(maxEnergy < tmax) { tmax = maxEnergy; }
263  if(tmin >= tmax) { return; }
264 
265  G4ThreeVector direction= dp->GetMomentumDirection();
266  G4double scaledTkin = kineticEnergy*fRatio;
267  G4double totalEnergy = kineticEnergy + fMass;
268  G4double totalMomentum = sqrt(kineticEnergy*(totalEnergy+fMass));
269 
270  G4double deltaTkin =
271  fModelData->SamplePostStepTransfer(coupleIndex, scaledTkin, tmax);
272 
273  //G4cout<<"G4PAIModel::SampleSecondaries; deltaKIn = "<<deltaTkin/keV
274  // <<" keV "<< " Escaled(MeV)= " << scaledTkin << G4endl;
275 
276  if( !(deltaTkin <= 0.) && !(deltaTkin > 0)) {
277  G4cout<<"G4PAIModel::SampleSecondaries; deltaKIn = "<<deltaTkin/keV
278  <<" keV "<< " Escaled(MeV)= " << scaledTkin << G4endl;
279  return;
280  }
281  if( deltaTkin <= 0.) { return; }
282 
283  if( deltaTkin > tmax) { deltaTkin = tmax; }
284 
285  const G4Element* anElement = SelectRandomAtom(matCC,fParticle,kineticEnergy);
286  G4int Z = G4lrint(anElement->GetZ());
287 
289  GetAngularDistribution()->SampleDirection(dp, deltaTkin,
290  Z, matCC->GetMaterial()),
291  deltaTkin);
292 
293  // primary change
294  kineticEnergy -= deltaTkin;
295  G4ThreeVector dir = totalMomentum*direction - deltaRay->GetMomentum();
296  direction = dir.unit();
299 
300  vdp->push_back(deltaRay);
301 }
302 
304 
306  const G4DynamicParticle* aParticle,
307  G4double, G4double step,
308  G4double eloss)
309 {
310  G4int coupleIndex = FindCoupleIndex(matCC);
311  if(0 > coupleIndex) { return eloss; }
312 
313  SetParticle(aParticle->GetDefinition());
314 
315  /*
316  G4cout << "G4PAIModel::SampleFluctuations step(mm)= "<< step/mm
317  << " Eloss(keV)= " << eloss/keV << " in "
318  << matCC->Getmaterial()->GetName() << G4endl;
319  */
320 
321  G4double Tkin = aParticle->GetKineticEnergy();
322  G4double scaledTkin = Tkin*fRatio;
323 
324  G4double loss = fModelData->SampleAlongStepTransfer(coupleIndex, Tkin,
325  scaledTkin,
326  step*fChargeSquare);
327 
328  // G4cout<<"PAIModel AlongStepLoss = "<<loss/keV<<" keV, on step = "
329  //<<step/mm<<" mm"<<G4endl;
330  return loss;
331 
332 }
333 
335 //
336 // Returns the statistical estimation of the energy loss distribution variance
337 //
338 
339 
341  const G4DynamicParticle* aParticle,
342  G4double tmax,
343  G4double step )
344 {
345  G4double particleMass = aParticle->GetMass();
346  G4double electronDensity = material->GetElectronDensity();
347  G4double kineticEnergy = aParticle->GetKineticEnergy();
348  G4double q = aParticle->GetCharge()/eplus;
349  G4double etot = kineticEnergy + particleMass;
350  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*particleMass)/(etot*etot);
351  G4double siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * step
352  * electronDensity * q * q;
353 
354  return siga;
355  /*
356  G4double loss, sumLoss=0., sumLoss2=0., sigma2, meanLoss=0.;
357  for(G4int i = 0; i < fMeanNumber; i++)
358  {
359  loss = SampleFluctuations(material,aParticle,tmax,step,meanLoss);
360  sumLoss += loss;
361  sumLoss2 += loss*loss;
362  }
363  meanLoss = sumLoss/fMeanNumber;
364  sigma2 = meanLoss*meanLoss + (sumLoss2-2*sumLoss*meanLoss)/fMeanNumber;
365  return sigma2;
366  */
367 }
368 
370 
372  G4double kinEnergy)
373 {
374  SetParticle(p);
375  G4double tmax = kinEnergy;
376  if(p == fElectron) { tmax *= 0.5; }
377  else if(p != fPositron) {
378  G4double ratio= electron_mass_c2/fMass;
379  G4double gamma= kinEnergy/fMass + 1.0;
380  tmax = 2.0*electron_mass_c2*(gamma*gamma - 1.) /
381  (1. + 2.0*gamma*ratio + ratio*ratio);
382  }
383  return tmax;
384 }
385 
387 
389 {
390  fPAIRegionVector.push_back(r);
391 }
392 
393 //
394 //
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:622
G4PAIModelData * GetPAIModelData()
Definition: G4PAIModel.hh:153
std::vector< const G4MaterialCutsCouple * > fMaterialCutsCoupleVector
Definition: G4PAIModel.hh:138
const G4String & GetName() const
const G4ParticleDefinition * fPositron
Definition: G4PAIModel.hh:143
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:120
G4int FindCoupleIndex(const G4MaterialCutsCouple *)
Definition: G4PAIModel.hh:168
G4double fMass
Definition: G4PAIModel.hh:146
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:615
std::vector< const G4Region * > fPAIRegionVector
Definition: G4PAIModel.hh:139
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
Definition: G4PAIModel.cc:371
G4double fRatio
Definition: G4PAIModel.hh:147
const G4String & GetName() const
Definition: G4Material.hh:178
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4PAIModel.cc:111
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4PAIModel.cc:192
G4double GetZ() const
Definition: G4Element.hh:131
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:588
std::vector< G4Material * > G4MaterialTable
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:598
G4int fVerbose
Definition: G4PAIModel.hh:134
G4bool IsMaster() const
Definition: G4VEmModel.hh:699
G4ParticleDefinition * GetDefinition() const
G4bool isInitialised
Definition: G4PAIModel.hh:150
int G4int
Definition: G4Types.hh:78
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double, G4double)
Definition: G4PAIModel.cc:340
G4double SampleAlongStepTransfer(G4int coupleIndex, G4double kinEnergy, G4double scaledTkin, G4double stepFactor) const
const G4String & GetParticleName() const
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
Definition: G4PAIModel.cc:225
static const G4double reg
G4double SamplePostStepTransfer(G4int coupleIndex, G4double scaledTkin, G4double scaledTmax) const
void SetParticle(const G4ParticleDefinition *p)
Definition: G4PAIModel.hh:181
const G4ParticleDefinition * fParticle
Definition: G4PAIModel.hh:141
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
Definition: G4PAIModel.cc:249
G4double DEDXPerVolume(G4int coupleIndex, G4double scaledTkin, G4double cut) const
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:433
G4GLOB_DLL std::ostream G4cout
G4double fChargeSquare
Definition: G4PAIModel.hh:148
G4double GetElectronDensity() const
Definition: G4Material.hh:217
G4double GetMass() const
G4PAIModelData * fModelData
Definition: G4PAIModel.hh:136
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
G4double GetCharge() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4MaterialCutsCouple * FindCouple(G4Material *mat)
const G4int n
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:783
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:595
const G4ParticleDefinition * fElectron
Definition: G4PAIModel.hh:142
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double CrossSectionPerVolume(G4int coupleIndex, G4double scaledTkin, G4double tcut, G4double tmax) const
G4ParticleChangeForLoss * fParticleChange
Definition: G4PAIModel.hh:144
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:791
virtual ~G4PAIModel()
Definition: G4PAIModel.cc:104
int G4lrint(double ad)
Definition: templates.hh:163
G4PAIModel(const G4ParticleDefinition *p=0, const G4String &nam="PAI")
Definition: G4PAIModel.cc:82
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:605
void Initialise(const G4MaterialCutsCouple *, G4double cut, G4PAIModel *)
virtual G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double, G4double, G4double)
Definition: G4PAIModel.cc:305
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
Definition: G4PAIModel.cc:201
static const double keV
Definition: G4SIunits.hh:195
double G4double
Definition: G4Types.hh:76
void DefineForRegion(const G4Region *r)
Definition: G4PAIModel.cc:388
static const double eplus
Definition: G4SIunits.hh:178
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:525
const G4Material * GetMaterial() const