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