Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 73607 2013-09-02 10:04:03Z 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 
78 ////////////////////////////////////////////////////////////////////////
79 
80 using namespace std;
81 
84  fVerbose(0),
85  fModelData(0),
86  fParticle(0)
87 {
88  fElectron = G4Electron::Electron();
89  fPositron = G4Positron::Positron();
90 
91  fParticleChange = 0;
92 
93  if(p) { SetParticle(p); }
94  else { SetParticle(fElectron); }
95 
96  // default generator
98 
99  isInitialised = false;
100 }
101 
102 ////////////////////////////////////////////////////////////////////////////
103 
105 {
106  if(IsMaster()) { delete fModelData; }
107 }
108 
109 ////////////////////////////////////////////////////////////////////////////
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);
122  fParticleChange = GetParticleChangeForLoss();
123 
124  if(IsMaster()) {
125 
127 
128  if(!fModelData) {
129 
130  G4double tmin = LowEnergyLimit()*fRatio;
131  G4double tmax = HighEnergyLimit()*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  for(size_t iReg = 0; iReg < numRegions; ++iReg) {
140  const G4Region* curReg = fPAIRegionVector[iReg];
141  G4Region* reg = const_cast<G4Region*>(curReg);
142 
143  for(size_t jMat = 0; jMat < numOfMat; ++jMat) {
144  G4Material* mat = (*theMaterialTable)[jMat];
145  const G4MaterialCutsCouple* cutCouple = reg->FindCouple(mat);
146  //G4cout << "Couple <" << fCutCouple << G4endl;
147  if(cutCouple) {
148  /*
149  G4cout << "Reg <" <<curReg->GetName() << "> mat <"
150  << fMaterial->GetName() << "> fCouple= "
151  << fCutCouple << " idx= " << fCutCouple->GetIndex()
152  <<" " << p->GetParticleName() <<G4endl;
153  // G4cout << cuts.size() << G4endl;
154  */
155  // check if this couple is not already initialized
156  size_t n = fMaterialCutsCoupleVector.size();
157  if(0 < n) {
158  for(size_t i=0; i<fMaterialCutsCoupleVector.size(); ++i) {
159  if(cutCouple == fMaterialCutsCoupleVector[i]) {
160  break;
161  }
162  }
163  }
164  // initialise data banks
165  fMaterialCutsCoupleVector.push_back(cutCouple);
166  G4double deltaCutInKinEnergy = cuts[cutCouple->GetIndex()];
167  fModelData->Initialise(cutCouple, deltaCutInKinEnergy, this);
168  }
169  }
170  }
171  }
172 }
173 
174 /////////////////////////////////////////////////////////////////////////
175 
177  G4VEmModel* masterModel)
178 {
179  fModelData = static_cast<G4PAIModel*>(masterModel)->GetPAIModelData();
181 }
182 
183 //////////////////////////////////////////////////////////////////////////////
184 
186  const G4ParticleDefinition* p,
187  G4double kineticEnergy,
188  G4double cutEnergy)
189 {
190  G4int coupleIndex = FindCoupleIndex(CurrentCouple());
191  if(0 > coupleIndex) { return 0.0; }
192 
193  G4double cut = std::min(MaxSecondaryEnergy(p, kineticEnergy), cutEnergy);
194 
195  G4double scaledTkin = kineticEnergy*fRatio;
196 
197  return fChargeSquare*fModelData->DEDXPerVolume(coupleIndex, scaledTkin, cut);
198 }
199 
200 /////////////////////////////////////////////////////////////////////////
201 
203  const G4ParticleDefinition* p,
204  G4double kineticEnergy,
205  G4double cutEnergy,
206  G4double maxEnergy )
207 {
208  G4int coupleIndex = FindCoupleIndex(CurrentCouple());
209  if(0 > coupleIndex) { return 0.0; }
210 
211  G4double tmax = std::min(MaxSecondaryEnergy(p, kineticEnergy), maxEnergy);
212  if(tmax <= cutEnergy) { return 0.0; }
213 
214  G4double scaledTkin = kineticEnergy*fRatio;
215 
216  return fChargeSquare*fModelData->CrossSectionPerVolume(coupleIndex,
217  scaledTkin,
218  cutEnergy, tmax);
219 }
220 
221 ///////////////////////////////////////////////////////////////////////////
222 //
223 // It is analog of PostStepDoIt in terms of secondary electron.
224 //
225 
226 void G4PAIModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
227  const G4MaterialCutsCouple* matCC,
228  const G4DynamicParticle* dp,
229  G4double tmin,
230  G4double maxEnergy)
231 {
232  G4int coupleIndex = FindCoupleIndex(matCC);
233  if(0 > coupleIndex) { return; }
234 
235  SetParticle(dp->GetDefinition());
236  G4double kineticEnergy = dp->GetKineticEnergy();
237 
238  G4double tmax = MaxSecondaryEnergy(fParticle, kineticEnergy);
239  if(maxEnergy < tmax) { tmax = maxEnergy; }
240  if(tmin >= tmax) { return; }
241 
242  G4ThreeVector direction= dp->GetMomentumDirection();
243  G4double scaledTkin = kineticEnergy*fRatio;
244  G4double totalEnergy = kineticEnergy + fMass;
245  G4double totalMomentum = sqrt(kineticEnergy*(totalEnergy+fMass));
246 
247  G4double deltaTkin =
248  fModelData->SamplePostStepTransfer(coupleIndex, scaledTkin);
249 
250  // G4cout<<"G4PAIModel::SampleSecondaries; deltaKIn = "<<deltaTkin/keV
251  // <<" keV "<<G4endl;
252 
253  if( deltaTkin <= 0. && fVerbose > 0)
254  {
255  G4cout<<"G4PAIModel::SampleSecondary e- deltaTkin = "<<deltaTkin<<G4endl;
256  }
257  if( deltaTkin <= 0.) { return; }
258 
259  if( deltaTkin > tmax) { deltaTkin = tmax; }
260 
261  const G4Element* anElement = SelectRandomAtom(matCC,fParticle,kineticEnergy);
262  G4int Z = G4lrint(anElement->GetZ());
263 
264  G4DynamicParticle* deltaRay = new G4DynamicParticle(fElectron,
265  GetAngularDistribution()->SampleDirection(dp, deltaTkin,
266  Z, matCC->GetMaterial()),
267  deltaTkin);
268 
269  // primary change
270  kineticEnergy -= deltaTkin;
271  G4ThreeVector dir = totalMomentum*direction - deltaRay->GetMomentum();
272  direction = dir.unit();
273  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
274  fParticleChange->SetProposedMomentumDirection(direction);
275 
276  vdp->push_back(deltaRay);
277 }
278 
279 ///////////////////////////////////////////////////////////////////////
280 
282  const G4DynamicParticle* aParticle,
283  G4double, G4double step,
284  G4double eloss)
285 {
286  G4int coupleIndex = FindCoupleIndex(matCC);
287  if(0 > coupleIndex) { return eloss; }
288 
289  SetParticle(aParticle->GetDefinition());
290 
291  /*
292  G4cout << "G4PAIModel::SampleFluctuations step(mm)= "<< step/mm
293  << " Eloss(keV)= " << eloss/keV << " in "
294  << matCC->Getmaterial()->GetName() << G4endl;
295  */
296 
297  G4double Tkin = aParticle->GetKineticEnergy();
298  G4double scaledTkin = Tkin*fRatio;
299 
300  G4double loss = fModelData->SampleAlongStepTransfer(coupleIndex, Tkin,
301  scaledTkin,
302  step*fChargeSquare);
303 
304  // G4cout<<"PAIModel AlongStepLoss = "<<loss/keV<<" keV, on step = "
305  //<<step/mm<<" mm"<<G4endl;
306  return loss;
307 
308 }
309 
310 //////////////////////////////////////////////////////////////////////
311 //
312 // Returns the statistical estimation of the energy loss distribution variance
313 //
314 
315 
317  const G4DynamicParticle* aParticle,
318  G4double tmax,
319  G4double step )
320 {
321  G4double particleMass = aParticle->GetMass();
322  G4double electronDensity = material->GetElectronDensity();
323  G4double kineticEnergy = aParticle->GetKineticEnergy();
324  G4double q = aParticle->GetCharge()/eplus;
325  G4double etot = kineticEnergy + particleMass;
326  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*particleMass)/(etot*etot);
327  G4double siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * step
328  * electronDensity * q * q;
329 
330  return siga;
331  /*
332  G4double loss, sumLoss=0., sumLoss2=0., sigma2, meanLoss=0.;
333  for(G4int i = 0; i < fMeanNumber; i++)
334  {
335  loss = SampleFluctuations(material,aParticle,tmax,step,meanLoss);
336  sumLoss += loss;
337  sumLoss2 += loss*loss;
338  }
339  meanLoss = sumLoss/fMeanNumber;
340  sigma2 = meanLoss*meanLoss + (sumLoss2-2*sumLoss*meanLoss)/fMeanNumber;
341  return sigma2;
342  */
343 }
344 
345 /////////////////////////////////////////////////////////////////////
346 
348  G4double kinEnergy)
349 {
350  SetParticle(p);
351  G4double tmax = kinEnergy;
352  if(p == fElectron) { tmax *= 0.5; }
353  else if(p != fPositron) {
354  G4double ratio= electron_mass_c2/fMass;
355  G4double gamma= kinEnergy/fMass + 1.0;
356  tmax = 2.0*electron_mass_c2*(gamma*gamma - 1.) /
357  (1. + 2.0*gamma*ratio + ratio*ratio);
358  }
359  return tmax;
360 }
361 
362 ///////////////////////////////////////////////////////////////
363 
365 {
366  fPAIRegionVector.push_back(r);
367 }
368 
369 //
370 //
371 /////////////////////////////////////////////////
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
G4PAIModelData * GetPAIModelData()
Definition: G4PAIModel.hh:153
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
G4double GetKineticEnergy() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
Definition: G4PAIModel.cc:347
const char * p
Definition: xmltok.h:285
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4PAIModel.cc:111
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4PAIModel.cc:176
G4double GetZ() const
Definition: G4Element.hh:131
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:564
std::vector< G4Material * > G4MaterialTable
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:578
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double, G4double)
Definition: G4PAIModel.cc:316
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:202
string material
Definition: eplot.py:19
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
Definition: G4PAIModel.cc:226
G4double DEDXPerVolume(G4int coupleIndex, G4double scaledTkin, G4double cut) const
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:426
G4GLOB_DLL std::ostream G4cout
G4double GetElectronDensity() const
Definition: G4Material.hh:215
G4double GetMass() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetCharge() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
float electron_mass_c2
Definition: hepunit.py:274
G4MaterialCutsCouple * FindCouple(G4Material *mat)
const G4int n
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:571
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double CrossSectionPerVolume(G4int coupleIndex, G4double scaledTkin, G4double tcut, G4double tmax) const
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:768
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
Hep3Vector unit() const
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:585
void Initialise(const G4MaterialCutsCouple *, G4double cut, G4PAIModel *)
virtual G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double, G4double, G4double)
Definition: G4PAIModel.cc:281
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:185
double G4double
Definition: G4Types.hh:76
void DefineForRegion(const G4Region *r)
Definition: G4PAIModel.cc:364
G4double SamplePostStepTransfer(G4int coupleIndex, G4double scaledTkin) const
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:510
const G4Material * GetMaterial() const