Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PenelopeBremsstrahlungModel.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// Author: Luciano Pandola
28//
29// History:
30// --------
31// 23 Nov 2010 L Pandola First complete implementation, Penelope v2008
32// 24 May 2011 L. Pandola Renamed (make default Penelope)
33// 13 Mar 2012 L. Pandola Updated the interface for the angular generator
34// 18 Jul 2012 L. Pandola Migrate to the new interface of the angular generator, which
35// now provides the G4ThreeVector and takes care of rotation
36// 02 Oct 2013 L. Pandola Migrated to MT
37// 17 Oct 2013 L. Pandola Partially revert the MT migration: the angular generator is
38// kept as thread-local, and created/managed by the workers.
39//
40
43#include "G4SystemOfUnits.hh"
49#include "G4DynamicParticle.hh"
50#include "G4Gamma.hh"
51#include "G4Electron.hh"
52#include "G4Positron.hh"
56#include "G4PhysicsLogVector.hh"
57#include "G4PhysicsTable.hh"
58#include "G4Exp.hh"
59
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
61namespace { G4Mutex PenelopeBremsstrahlungModelMutex = G4MUTEX_INITIALIZER; }
62
64 const G4String& nam)
65 :G4VEmModel(nam),fParticleChange(nullptr),fParticle(nullptr),
66 fPenelopeFSHelper(nullptr),fPenelopeAngular(nullptr),fEnergyGrid(nullptr),
67 fXSTableElectron(nullptr),fXSTablePositron(nullptr),
68 fIsInitialised(false),fLocalTable(false)
69
70{
71 fIntrinsicLowEnergyLimit = 100.0*eV;
72 fIntrinsicHighEnergyLimit = 100.0*GeV;
73 nBins = 200;
74
75 if (part)
76 SetParticle(part);
77
78 SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
79 //
81 //
82 fVerboseLevel= 0;
83 // Verbosity scale:
84 // 0 = nothing
85 // 1 = warning for energy non-conservation
86 // 2 = details of energy budget
87 // 3 = calculation of cross sections, file openings, sampling of atoms
88 // 4 = entering in methods
89
90 // Atomic deexcitation model activated by default
92
93}
94
95//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
96
98{
99 if (IsMaster() || fLocalTable)
100 {
101 ClearTables();
102 if (fPenelopeFSHelper)
103 delete fPenelopeFSHelper;
104 }
105 // This is thread-local at the moment
106 if (fPenelopeAngular)
107 delete fPenelopeAngular;
108}
109
110//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
111
113 const G4DataVector& theCuts)
114{
115 if (fVerboseLevel > 3)
116 G4cout << "Calling G4PenelopeBremsstrahlungModel::Initialise()" << G4endl;
117
118 SetParticle(part);
119
120 if (IsMaster() && part == fParticle)
121 {
122 if (!fPenelopeFSHelper)
123 fPenelopeFSHelper = new G4PenelopeBremsstrahlungFS(fVerboseLevel);
124 if (!fPenelopeAngular)
125 fPenelopeAngular = new G4PenelopeBremsstrahlungAngular();
126 //Clear and re-build the tables
127 ClearTables();
128
129 //forces the cleaning of tables, in this specific case
130 if (fPenelopeAngular)
131 fPenelopeAngular->Initialize();
132
133 //Set the number of bins for the tables. 20 points per decade
134 nBins = (std::size_t) (20*std::log10(HighEnergyLimit()/LowEnergyLimit()));
135 nBins = std::max(nBins,(std::size_t)100);
136 fEnergyGrid = new G4PhysicsLogVector(LowEnergyLimit(),
138 nBins-1); //one hidden bin is added
139
140 fXSTableElectron = new
141 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*>;
142 fXSTablePositron = new
143 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*>;
144
145 G4ProductionCutsTable* theCoupleTable =
147
148 //Build tables for all materials
149 for (G4int i=0;i<(G4int)theCoupleTable->GetTableSize();++i)
150 {
151 const G4Material* theMat =
152 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
153 //Forces the building of the helper tables
154 fPenelopeFSHelper->BuildScaledXSTable(theMat,theCuts.at(i),IsMaster());
155 fPenelopeAngular->PrepareTables(theMat,IsMaster());
156 BuildXSTable(theMat,theCuts.at(i));
157
158 }
159
160 if (fVerboseLevel > 2) {
161 G4cout << "Penelope Bremsstrahlung model v2008 is initialized " << G4endl
162 << "Energy range: "
163 << LowEnergyLimit() / keV << " keV - "
164 << HighEnergyLimit() / GeV << " GeV."
165 << G4endl;
166 }
167 }
168
169 if(fIsInitialised) return;
171 fIsInitialised = true;
172}
173
174//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
175
177 G4VEmModel *masterModel)
178{
179 if (fVerboseLevel > 3)
180 G4cout << "Calling G4PenelopeBremsstrahlungModel::InitialiseLocal()" << G4endl;
181 //
182 //Check that particle matches: one might have multiple master models (e.g.
183 //for e+ and e-).
184 //
185 if (part == fParticle)
186 {
187 //Get the const table pointers from the master to the workers
188 const G4PenelopeBremsstrahlungModel* theModel =
189 static_cast<G4PenelopeBremsstrahlungModel*> (masterModel);
190
191 //Copy pointers to the data tables
192 fEnergyGrid = theModel->fEnergyGrid;
193 fXSTableElectron = theModel->fXSTableElectron;
194 fXSTablePositron = theModel->fXSTablePositron;
195 fPenelopeFSHelper = theModel->fPenelopeFSHelper;
196
197 //created in each thread and initialized.
198 if (!fPenelopeAngular)
199 fPenelopeAngular = new G4PenelopeBremsstrahlungAngular();
200 //forces the cleaning of tables, in this specific case
201 if (fPenelopeAngular)
202 fPenelopeAngular->Initialize();
203
204 G4ProductionCutsTable* theCoupleTable =
206 //Build tables for all materials
207 for (G4int i=0;i<(G4int)theCoupleTable->GetTableSize();++i)
208 {
209 const G4Material* theMat =
210 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
211 fPenelopeAngular->PrepareTables(theMat,IsMaster());
212 }
213
214 //copy the data
215 nBins = theModel->nBins;
216
217 //Same verbosity for all workers, as the master
218 fVerboseLevel = theModel->fVerboseLevel;
219 }
220 return;
221}
222
223//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
224
226 const G4ParticleDefinition* theParticle,
227 G4double energy,
228 G4double cutEnergy,
229 G4double)
230{
231 //
232 if (fVerboseLevel > 3)
233 G4cout << "Calling CrossSectionPerVolume() of G4PenelopeBremsstrahlungModel" << G4endl;
234
235 SetupForMaterial(theParticle, material, energy);
236 G4double crossPerMolecule = 0.;
237
238 G4PenelopeCrossSection* theXS = GetCrossSectionTableForCouple(theParticle,material,
239 cutEnergy);
240 if (theXS)
241 crossPerMolecule = theXS->GetHardCrossSection(energy);
242
243 G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
244 G4double atPerMol = fOscManager->GetAtomsPerMolecule(material);
245
246 if (fVerboseLevel > 3)
247 G4cout << "Material " << material->GetName() << " has " << atPerMol <<
248 "atoms per molecule" << G4endl;
249
250 G4double moleculeDensity = 0.;
251 if (atPerMol)
252 moleculeDensity = atomDensity/atPerMol;
253
254 G4double crossPerVolume = crossPerMolecule*moleculeDensity;
255
256 if (fVerboseLevel > 2)
257 {
258 G4cout << "G4PenelopeBremsstrahlungModel " << G4endl;
259 G4cout << "Mean free path for gamma emission > " << cutEnergy/keV << " keV at " <<
260 energy/keV << " keV = " <<
261 (crossPerVolume? (1./crossPerVolume)/mm : DBL_MAX) << " mm" << G4endl;
262 }
263
264 return crossPerVolume;
265}
266
267//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
268
269//This is a dummy method. Never inkoved by the tracking, it just issues
270//a warning if one tries to get Cross Sections per Atom via the
271//G4EmCalculator.
273 G4double,
274 G4double,
275 G4double,
276 G4double,
277 G4double)
278{
279 G4cout << "*** G4PenelopeBremsstrahlungModel -- WARNING ***" << G4endl;
280 G4cout << "Penelope Bremsstrahlung model v2008 does not calculate cross section _per atom_ " << G4endl;
281 G4cout << "so the result is always zero. For physics values, please invoke " << G4endl;
282 G4cout << "GetCrossSectionPerVolume() or GetMeanFreePath() via the G4EmCalculator" << G4endl;
283 return 0;
284}
285
286//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
287
289 const G4ParticleDefinition* theParticle,
290 G4double kineticEnergy,
291 G4double cutEnergy)
292{
293 if (fVerboseLevel > 3)
294 G4cout << "Calling ComputeDEDX() of G4PenelopeBremsstrahlungModel" << G4endl;
295
296 G4PenelopeCrossSection* theXS = GetCrossSectionTableForCouple(theParticle,material,
297 cutEnergy);
298 G4double sPowerPerMolecule = 0.0;
299 if (theXS)
300 sPowerPerMolecule = theXS->GetSoftStoppingPower(kineticEnergy);
301
302 G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
303 G4double atPerMol = fOscManager->GetAtomsPerMolecule(material);
304
305 G4double moleculeDensity = 0.;
306 if (atPerMol)
307 moleculeDensity = atomDensity/atPerMol;
308
309 G4double sPowerPerVolume = sPowerPerMolecule*moleculeDensity;
310
311 if (fVerboseLevel > 2)
312 {
313 G4cout << "G4PenelopeBremsstrahlungModel " << G4endl;
314 G4cout << "Stopping power < " << cutEnergy/keV << " keV at " <<
315 kineticEnergy/keV << " keV = " <<
316 sPowerPerVolume/(keV/mm) << " keV/mm" << G4endl;
317 }
318 return sPowerPerVolume;
319}
320
321//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
322
323void G4PenelopeBremsstrahlungModel::SampleSecondaries(std::vector<G4DynamicParticle*>*fvect,
324 const G4MaterialCutsCouple* couple,
325 const G4DynamicParticle* aDynamicParticle,
326 G4double cutG,
327 G4double)
328{
329 if (fVerboseLevel > 3)
330 G4cout << "Calling SampleSecondaries() of G4PenelopeBremsstrahlungModel" << G4endl;
331
332 G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
333 const G4Material* material = couple->GetMaterial();
334
335 if (kineticEnergy <= fIntrinsicLowEnergyLimit)
336 {
337 fParticleChange->SetProposedKineticEnergy(0.);
338 fParticleChange->ProposeLocalEnergyDeposit(kineticEnergy);
339 return ;
340 }
341
342 G4ParticleMomentum particleDirection0 = aDynamicParticle->GetMomentumDirection();
343 //This is the momentum
344 G4ThreeVector initialMomentum = aDynamicParticle->GetMomentum();
345
346 //Not enough energy to produce a secondary! Return with nothing happened
347 if (kineticEnergy < cutG)
348 return;
349
350 if (fVerboseLevel > 3)
351 G4cout << "Going to sample gamma energy for: " <<material->GetName() << " " <<
352 "energy = " << kineticEnergy/keV << ", cut = " << cutG/keV << G4endl;
353
354 //Sample gamma's energy according to the spectrum
355 G4double gammaEnergy =
356 fPenelopeFSHelper->SampleGammaEnergy(kineticEnergy,material,cutG);
357
358 if (fVerboseLevel > 3)
359 G4cout << "Sampled gamma energy: " << gammaEnergy/keV << " keV" << G4endl;
360
361 //Now sample the direction for the Gamma. Notice that the rotation is already done
362 //Z is unused here, I plug 0. The information is in the material pointer
363 G4ThreeVector gammaDirection1 =
364 fPenelopeAngular->SampleDirection(aDynamicParticle,gammaEnergy,0,material);
365
366 if (fVerboseLevel > 3)
367 G4cout << "Sampled cosTheta for e-: " << gammaDirection1.cosTheta() << G4endl;
368
369 G4double residualPrimaryEnergy = kineticEnergy-gammaEnergy;
370 if (residualPrimaryEnergy < 0)
371 {
372 //Ok we have a problem, all energy goes with the gamma
373 gammaEnergy += residualPrimaryEnergy;
374 residualPrimaryEnergy = 0.0;
375 }
376
377 //Produce final state according to momentum conservation
378 G4ThreeVector particleDirection1 = initialMomentum - gammaEnergy*gammaDirection1;
379 particleDirection1 = particleDirection1.unit(); //normalize
380
381 //Update the primary particle
382 if (residualPrimaryEnergy > 0.)
383 {
384 fParticleChange->ProposeMomentumDirection(particleDirection1);
385 fParticleChange->SetProposedKineticEnergy(residualPrimaryEnergy);
386 }
387 else
388 {
389 fParticleChange->SetProposedKineticEnergy(0.);
390 }
391
392 //Now produce the photon
394 gammaDirection1,
395 gammaEnergy);
396 fvect->push_back(theGamma);
397
398 if (fVerboseLevel > 1)
399 {
400 G4cout << "-----------------------------------------------------------" << G4endl;
401 G4cout << "Energy balance from G4PenelopeBremsstrahlung" << G4endl;
402 G4cout << "Incoming primary energy: " << kineticEnergy/keV << " keV" << G4endl;
403 G4cout << "-----------------------------------------------------------" << G4endl;
404 G4cout << "Outgoing primary energy: " << residualPrimaryEnergy/keV << " keV" << G4endl;
405 G4cout << "Bremsstrahlung photon " << gammaEnergy/keV << " keV" << G4endl;
406 G4cout << "Total final state: " << (residualPrimaryEnergy+gammaEnergy)/keV
407 << " keV" << G4endl;
408 G4cout << "-----------------------------------------------------------" << G4endl;
409 }
410
411 if (fVerboseLevel > 0)
412 {
413 G4double energyDiff = std::fabs(residualPrimaryEnergy+gammaEnergy-kineticEnergy);
414 if (energyDiff > 0.05*keV)
415 G4cout << "Warning from G4PenelopeBremsstrahlung: problem with energy conservation: "
416 <<
417 (residualPrimaryEnergy+gammaEnergy)/keV <<
418 " keV (final) vs. " <<
419 kineticEnergy/keV << " keV (initial)" << G4endl;
420 }
421 return;
422}
423
424//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
425
426void G4PenelopeBremsstrahlungModel::ClearTables()
427{
428 if (!IsMaster() && !fLocalTable)
429 //Should not be here!
430 G4Exception("G4PenelopeBremsstrahlungModel::ClearTables()",
431 "em0100",FatalException,"Worker thread in this method");
432
433 if (fXSTableElectron)
434 {
435 for (auto& item : (*fXSTableElectron))
436 delete item.second;
437 delete fXSTableElectron;
438 fXSTableElectron = nullptr;
439 }
440 if (fXSTablePositron)
441 {
442 for (auto& item : (*fXSTablePositron))
443 delete item.second;
444 delete fXSTablePositron;
445 fXSTablePositron = nullptr;
446 }
447 /*
448 if (fEnergyGrid)
449 delete fEnergyGrid;
450 */
451 if (fPenelopeFSHelper)
452 fPenelopeFSHelper->ClearTables(IsMaster());
453
454 if (fVerboseLevel > 2)
455 G4cout << "G4PenelopeBremsstrahlungModel: cleared tables" << G4endl;
456
457 return;
458}
459
460//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
461
464{
465 return fIntrinsicLowEnergyLimit;
466}
467
468//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
469
470void G4PenelopeBremsstrahlungModel::BuildXSTable(const G4Material* mat,G4double cut)
471{
472 if (!IsMaster() && !fLocalTable)
473 //Should not be here!
474 G4Exception("G4PenelopeBremsstrahlungModel::BuildXSTable()",
475 "em0100",FatalException,"Worker thread in this method");
476
477 //The key of the map
478 std::pair<const G4Material*,G4double> theKey = std::make_pair(mat,cut);
479
480 //The table already exists
481 if (fXSTableElectron->count(theKey) && fXSTablePositron->count(theKey))
482 return;
483
484 //
485 //This method fills the G4PenelopeCrossSection containers for electrons or positrons
486 //and for the given material/cut couple.
487 //Equivalent of subroutines EBRaT and PINaT of Penelope
488 //
489 if (fVerboseLevel > 2)
490 {
491 G4cout << "G4PenelopeBremsstrahlungModel: going to build cross section table " << G4endl;
492 G4cout << "for e+/e- in " << mat->GetName() << " for Ecut(gamma)= " <<
493 cut/keV << " keV " << G4endl;
494 }
495
496 //Tables have been already created (checked by GetCrossSectionTableForCouple)
497 if (fEnergyGrid->GetVectorLength() != nBins)
498 {
500 ed << "Energy Grid looks not initialized" << G4endl;
501 ed << nBins << " " << fEnergyGrid->GetVectorLength() << G4endl;
502 G4Exception("G4PenelopeBremsstrahlungModel::BuildXSTable()",
503 "em2016",FatalException,ed);
504 }
505
506 G4PenelopeCrossSection* XSEntry = new G4PenelopeCrossSection(nBins);
507 G4PenelopeCrossSection* XSEntry2 = new G4PenelopeCrossSection(nBins);
508 const G4PhysicsTable* table = fPenelopeFSHelper->GetScaledXSTable(mat,cut);
509
510 //loop on the energy grid
511 for (std::size_t bin=0;bin<nBins;++bin)
512 {
513 G4double energy = fEnergyGrid->GetLowEdgeEnergy(bin);
514 G4double XH0=0, XH1=0, XH2=0;
515 G4double XS0=0, XS1=0, XS2=0;
516
517 //Global xs factor
518 G4double fact = fPenelopeFSHelper->GetEffectiveZSquared(mat)*
519 ((energy+electron_mass_c2)*(energy+electron_mass_c2)/
520 (energy*(energy+2.0*electron_mass_c2)));
521
522 G4double restrictedCut = cut/energy;
523
524 //Now I need the dSigma/dX profile - interpolated on energy - for
525 //the 32-point x grid. Interpolation is log-log
526 std::size_t nBinsX = fPenelopeFSHelper->GetNBinsX();
527 G4double* tempData = new G4double[nBinsX];
528 G4double logene = G4Log(energy);
529 for (std::size_t ix=0;ix<nBinsX;++ix)
530 {
531 //find dSigma/dx for the given E. X belongs to the 32-point grid.
532 G4double val = (*table)[ix]->Value(logene);
533 tempData[ix] = G4Exp(val); //back to the real value!
534 }
535
536 G4double XH0A = 0.;
537 if (restrictedCut <= 1) //calculate only if we are above threshold!
538 XH0A = fPenelopeFSHelper->GetMomentumIntegral(tempData,1.0,-1) -
539 fPenelopeFSHelper->GetMomentumIntegral(tempData,restrictedCut,-1);
540 G4double XS1A = fPenelopeFSHelper->GetMomentumIntegral(tempData,
541 restrictedCut,0);
542 G4double XS2A = fPenelopeFSHelper->GetMomentumIntegral(tempData,
543 restrictedCut,1);
544 G4double XH1A=0, XH2A=0;
545 if (restrictedCut <=1)
546 {
547 XH1A = fPenelopeFSHelper->GetMomentumIntegral(tempData,1.0,0) -
548 XS1A;
549 XH2A = fPenelopeFSHelper->GetMomentumIntegral(tempData,1.0,1) -
550 XS2A;
551 }
552 delete[] tempData;
553
554 XH0 = XH0A*fact;
555 XS1 = XS1A*fact*energy;
556 XH1 = XH1A*fact*energy;
557 XS2 = XS2A*fact*energy*energy;
558 XH2 = XH2A*fact*energy*energy;
559
560 XSEntry->AddCrossSectionPoint(bin,energy,XH0,XH1,XH2,XS0,XS1,XS2);
561
562 //take care of positrons
563 G4double posCorrection = GetPositronXSCorrection(mat,energy);
564 XSEntry2->AddCrossSectionPoint(bin,energy,XH0*posCorrection,
565 XH1*posCorrection,
566 XH2*posCorrection,
567 XS0,
568 XS1*posCorrection,
569 XS2*posCorrection);
570 }
571
572 //Insert in the appropriate table
573 fXSTableElectron->insert(std::make_pair(theKey,XSEntry));
574 fXSTablePositron->insert(std::make_pair(theKey,XSEntry2));
575
576 return;
577}
578
579
580//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
581
583G4PenelopeBremsstrahlungModel::GetCrossSectionTableForCouple(const G4ParticleDefinition* part,
584 const G4Material* mat,
585 G4double cut)
586{
587 if (part != G4Electron::Electron() && part != G4Positron::Positron())
588 {
590 ed << "Invalid particle: " << part->GetParticleName() << G4endl;
591 G4Exception("G4PenelopeBremsstrahlungModel::GetCrossSectionTableForCouple()",
592 "em0001",FatalException,ed);
593 return nullptr;
594 }
595
596 if (part == G4Electron::Electron())
597 {
598 //Either Initialize() was not called, or we are in a slave and InitializeLocal() was
599 //not invoked
600 if (!fXSTableElectron)
601 {
602 //create a **thread-local** version of the table. Used only for G4EmCalculator and
603 //Unit Tests
604 G4String excep = "The Cross Section Table for e- was not initialized correctly!";
605 G4Exception("G4PenelopeBremsstrahlungModel::GetCrossSectionTableForCouple()",
606 "em2013",JustWarning,excep);
607 fLocalTable = true;
608 fXSTableElectron = new
609 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*>;
610 if (!fEnergyGrid)
611 fEnergyGrid = new G4PhysicsLogVector(LowEnergyLimit(),
613 nBins-1); //one hidden bin is added
614 if (!fPenelopeFSHelper)
615 fPenelopeFSHelper = new G4PenelopeBremsstrahlungFS(fVerboseLevel);
616 }
617 std::pair<const G4Material*,G4double> theKey = std::make_pair(mat,cut);
618 if (fXSTableElectron->count(theKey)) //table already built
619 return fXSTableElectron->find(theKey)->second;
620 else
621 {
622 //If we are here, it means that Initialize() was inkoved, but the MaterialTable was
623 //not filled up. This can happen in a UnitTest or via G4EmCalculator
624 if (fVerboseLevel > 0)
625 {
626 //G4Exception (warning) is issued only in verbose mode
628 ed << "Unable to find e- table for " << mat->GetName() << " at Ecut(gamma)= "
629 << cut/keV << " keV " << G4endl;
630 ed << "This can happen only in Unit Tests or via G4EmCalculator" << G4endl;
631 G4Exception("G4PenelopeBremsstrahlungModel::GetCrossSectionTableForCouple()",
632 "em2009",JustWarning,ed);
633 }
634 //protect file reading via autolock
635 G4AutoLock lock(&PenelopeBremsstrahlungModelMutex);
636 fPenelopeFSHelper->BuildScaledXSTable(mat,cut,true); //pretend to be a master
637 BuildXSTable(mat,cut);
638 lock.unlock();
639 //now it should be ok
640 return fXSTableElectron->find(theKey)->second;
641 }
642 }
643 if (part == G4Positron::Positron())
644 {
645 //Either Initialize() was not called, or we are in a slave and InitializeLocal() was
646 //not invoked
647 if (!fXSTablePositron)
648 {
649 G4String excep = "The Cross Section Table for e+ was not initialized correctly!";
650 G4Exception("G4PenelopeBremsstrahlungModel::GetCrossSectionTableForCouple()",
651 "em2013",JustWarning,excep);
652 fLocalTable = true;
653 fXSTablePositron = new
654 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*>;
655 if (!fEnergyGrid)
656 fEnergyGrid = new G4PhysicsLogVector(LowEnergyLimit(),
658 nBins-1); //one hidden bin is added
659 if (!fPenelopeFSHelper)
660 fPenelopeFSHelper = new G4PenelopeBremsstrahlungFS(fVerboseLevel);
661 }
662 std::pair<const G4Material*,G4double> theKey = std::make_pair(mat,cut);
663 if (fXSTablePositron->count(theKey)) //table already built
664 return fXSTablePositron->find(theKey)->second;
665 else
666 {
667 //If we are here, it means that Initialize() was inkoved, but the MaterialTable was
668 //not filled up. This can happen in a UnitTest or via G4EmCalculator
669 if (fVerboseLevel > 0)
670 {
671 //Issue a G4Exception (warning) only in verbose mode
673 ed << "Unable to find e+ table for " << mat->GetName() << " at Ecut(gamma)= "
674 << cut/keV << " keV " << G4endl;
675 ed << "This can happen only in Unit Tests or via G4EmCalculator" << G4endl;
676 G4Exception("G4PenelopeBremsstrahlungModel::GetCrossSectionTableForCouple()",
677 "em2009",JustWarning,ed);
678 }
679 //protect file reading via autolock
680 G4AutoLock lock(&PenelopeBremsstrahlungModelMutex);
681 fPenelopeFSHelper->BuildScaledXSTable(mat,cut,true); //pretend to be a master
682 BuildXSTable(mat,cut);
683 lock.unlock();
684 //now it should be ok
685 return fXSTablePositron->find(theKey)->second;
686 }
687 }
688 return nullptr;
689}
690
691//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
692
693G4double G4PenelopeBremsstrahlungModel::GetPositronXSCorrection(const G4Material* mat,
694 G4double energy)
695{
696 //The electron-to-positron correction factor is set equal to the ratio of the
697 //radiative stopping powers for positrons and electrons, which has been calculated
698 //by Kim et al. (1986) (cf. Berger and Seltzer, 1982). Here, it is used an
699 //analytical approximation which reproduces the tabulated values with 0.5%
700 //accuracy
701 G4double t=G4Log(1.0+1e6*energy/
702 (electron_mass_c2*fPenelopeFSHelper->GetEffectiveZSquared(mat)));
703 G4double corr = 1.0-G4Exp(-t*(1.2359e-1-t*(6.1274e-2-t*
704 (3.1516e-2-t*(7.7446e-3-t*(1.0595e-3-t*
705 (7.0568e-5-t*
706 1.8080e-6)))))));
707 return corr;
708}
709
710//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...
711
712void G4PenelopeBremsstrahlungModel::SetParticle(const G4ParticleDefinition* p)
713{
714 if(!fParticle) {
715 fParticle = p;
716 }
717}
G4TemplateAutoLock< G4Mutex > G4AutoLock
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:132
G4double G4Log(G4double x)
Definition G4Log.hh:169
G4ThreeVector G4ParticleMomentum
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
double cosTheta() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
const G4Material * GetMaterial() const
G4double GetTotNbOfAtomsPerVolume() const
const G4String & GetName() const
const G4String & GetParticleName() const
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *theParticle, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override
G4PenelopeBremsstrahlungModel(const G4ParticleDefinition *p=nullptr, const G4String &processName="PenBrem")
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *theParticle, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
G4double GetSoftStoppingPower(G4double energy) const
Returns the total stopping power due to soft collisions.
G4double GetHardCrossSection(G4double energy) const
Returns hard cross section at the given energy.
void AddCrossSectionPoint(size_t binNumber, G4double energy, G4double XH0, G4double XH1, G4double XH2, G4double XS0, G4double XS1, G4double XS2)
static G4PenelopeOscillatorManager * GetOscillatorManager()
static G4Positron * Positron()
Definition G4Positron.cc:90
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetHighEnergyLimit(G4double)
G4double LowEnergyLimit() const
G4bool IsMaster() const
G4double HighEnergyLimit() const
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
void SetDeexcitationFlag(G4bool val)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
G4ParticleChangeForLoss * GetParticleChangeForLoss()
G4double energy(const ThreeVector &p, const G4double m)
#define DBL_MAX
Definition templates.hh:62