Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PenelopeBremsstrahlungFS.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
32// 02 May 2011 L.Pandola Remove dependency on CLHEP::HepMatrix
33// 24 May 2011 L.Pandola Renamed (make v2008 as default Penelope)
34// 03 Oct 2013 L.Pandola Migration to MT
35// 30 Oct 2013 L.Pandola Use G4Cache to avoid new/delete of the
36// data vector on the fly in SampleGammaEnergy()
37//
38
39//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
42#include "G4PhysicsTable.hh"
43#include "G4Material.hh"
44#include "Randomize.hh"
45#include "G4AutoDelete.hh"
46#include "G4Exp.hh"
48#include "G4SystemOfUnits.hh"
49
50//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
51
53 fReducedXSTable(nullptr),fEffectiveZSq(nullptr),fSamplingTable(nullptr),
54 fPBcut(nullptr),fVerbosity(verbosity)
55{
56 fCache.Put(0);
57 G4double tempvector[fNBinsX] =
58 {1.0e-12,0.025e0,0.05e0,0.075e0,0.1e0,0.15e0,0.2e0,0.25e0,
59 0.3e0,0.35e0,0.4e0,0.45e0,0.5e0,0.55e0,0.6e0,0.65e0,0.7e0,
60 0.75e0,0.8e0,0.85e0,0.9e0,0.925e0,0.95e0,0.97e0,0.99e0,
61 0.995e0,0.999e0,0.9995e0,0.9999e0,0.99995e0,0.99999e0,1.0e0};
62
63 for (std::size_t ix=0;ix<fNBinsX;++ix)
64 theXGrid[ix] = tempvector[ix];
65
66 for (std::size_t i=0;i<fNBinsE;++i)
67 theEGrid[i] = 0.;
68
69 fElementData = new std::map<G4int,G4DataVector*>;
70}
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73
75{
77
78 //The G4Physics*Vector pointers contained in the fCache are automatically deleted by
79 //the G4AutoDelete so there is no need to take care of them manually
80
81 //Clear manually fElementData
82 if (fElementData)
83 {
84 for (auto& item : (*fElementData))
85 delete item.second;
86 delete fElementData;
87 fElementData = nullptr;
88 }
89}
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...
92
93
95{
96 //Just to check
97 if (!isMaster)
98 G4Exception("G4PenelopeBremsstrahlungFS::ClearTables()",
99 "em0100",FatalException,"Worker thread in this method");
100
101 if (fReducedXSTable)
102 {
103 for (auto& item : (*fReducedXSTable))
104 {
105 G4PhysicsTable* tab = item.second;
106 tab->clearAndDestroy();
107 delete tab;
108 }
109 fReducedXSTable->clear();
110 delete fReducedXSTable;
111 fReducedXSTable = nullptr;
112 }
113
114 if (fSamplingTable)
115 {
116 for (auto& item : (*fSamplingTable))
117 {
118 G4PhysicsTable* tab = item.second;
119 tab->clearAndDestroy();
120 delete tab;
121 }
122 fSamplingTable->clear();
123 delete fSamplingTable;
124 fSamplingTable = nullptr;
125 }
126 if (fPBcut)
127 {
128 /*
129 std::map< std::pair<const G4Material*,G4double> ,G4PhysicsFreeVector*>::iterator kk;
130 for (kk=fPBcut->begin(); kk != fPBcut->end(); kk++)
131 delete kk->second;
132 */
133 delete fPBcut;
134 fPBcut = nullptr;
135 }
136
137 if (fEffectiveZSq)
138 {
139 delete fEffectiveZSq;
140 fEffectiveZSq = nullptr;
141 }
142
143 return;
144}
145
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
147
149{
150 if (!fEffectiveZSq)
151 {
153 ed << "The container for the <Z^2> values is not initialized" << G4endl;
154 G4Exception("G4PenelopeBremsstrahlungFS::GetEffectiveZSquared()",
155 "em2007",FatalException,ed);
156 return 0;
157 }
158 //found in the table: return it
159 if (fEffectiveZSq->count(material))
160 return fEffectiveZSq->find(material)->second;
161 else
162 {
164 ed << "The value of <Z^2> is not properly set for material " <<
165 material->GetName() << G4endl;
166 //requires running of BuildScaledXSTable()
167 G4Exception("G4PenelopeBremsstrahlungFS::GetEffectiveZSquared()",
168 "em2008",FatalException,ed);
169 }
170 return 0;
171}
172
173//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
174
176 G4double cut,G4bool isMaster)
177{
178 //Corresponds to subroutines EBRaW and EBRaR of PENELOPE
179 /*
180 This method generates the table of the scaled energy-loss cross section from
181 bremsstrahlung emission for the given material. Original data are read from
182 file. The table is normalized according to the Berger-Seltzer cross section.
183 */
184
185 //Just to check
186 if (!isMaster)
187 G4Exception("G4PenelopeBremsstrahlungFS::BuildScaledXSTable()",
188 "em0100",FatalException,"Worker thread in this method");
189
190 if (fVerbosity > 2)
191 {
192 G4cout << "Entering in G4PenelopeBremsstrahlungFS::BuildScaledXSTable for " <<
193 material->GetName() << G4endl;
194 G4cout << "Threshold = " << cut/keV << " keV, isMaster= " << isMaster <<
195 G4endl;
196 }
197
198 //This method should be accessed by the master only
199 if (!fSamplingTable)
200 fSamplingTable =
201 new std::map< std::pair<const G4Material*,G4double> , G4PhysicsTable*>;
202 if (!fPBcut)
203 fPBcut =
204 new std::map< std::pair<const G4Material*,G4double> , G4PhysicsFreeVector* >;
205
206 //check if the container exists (if not, create it)
207 if (!fReducedXSTable)
208 fReducedXSTable = new std::map< std::pair<const G4Material*,G4double> ,
210 if (!fEffectiveZSq)
211 fEffectiveZSq = new std::map<const G4Material*,G4double>;
212
213 //*********************************************************************
214 //Determine the equivalent atomic number <Z^2>
215 //*********************************************************************
216 std::vector<G4double> *StechiometricFactors = new std::vector<G4double>;
217 std::size_t nElements = material->GetNumberOfElements();
218 const G4ElementVector* elementVector = material->GetElementVector();
219 const G4double* fractionVector = material->GetFractionVector();
220 for (std::size_t i=0;i<nElements;i++)
221 {
222 G4double fraction = fractionVector[i];
223 G4double atomicWeigth = (*elementVector)[i]->GetA()/(g/mole);
224 StechiometricFactors->push_back(fraction/atomicWeigth);
225 }
226 //Find max
227 G4double MaxStechiometricFactor = 0.;
228 for (std::size_t i=0;i<nElements;i++)
229 {
230 if ((*StechiometricFactors)[i] > MaxStechiometricFactor)
231 MaxStechiometricFactor = (*StechiometricFactors)[i];
232 }
233 //Normalize
234 for (std::size_t i=0;i<nElements;i++)
235 if (MaxStechiometricFactor > 0.)
236 (*StechiometricFactors)[i] /= MaxStechiometricFactor;
237
238 G4double sumz2 = 0;
239 G4double sums = 0;
240 for (std::size_t i=0;i<nElements;i++)
241 {
242 G4double Z = (*elementVector)[i]->GetZ();
243 sumz2 += (*StechiometricFactors)[i]*Z*Z;
244 sums += (*StechiometricFactors)[i];
245 }
246 G4double ZBR2 = sumz2/sums;
247
248 fEffectiveZSq->insert(std::make_pair(material,ZBR2));
249
250 //*********************************************************************
251 // loop on elements and read data files
252 //*********************************************************************
253 G4DataVector* tempData = new G4DataVector(fNBinsE);
254 G4DataVector* tempMatrix = new G4DataVector(fNBinsE*fNBinsX,0.);
255
256 for (std::size_t iel=0;iel<nElements;iel++)
257 {
258 G4double Z = (*elementVector)[iel]->GetZ();
259 G4int iZ = (G4int) Z;
260 G4double wgt = (*StechiometricFactors)[iel]*Z*Z/ZBR2;
261
262 //the element is not already loaded
263 if (!fElementData->count(iZ))
264 {
265 ReadDataFile(iZ);
266 if (!fElementData->count(iZ))
267 {
269 ed << "Error in G4PenelopeBremsstrahlungFS::BuildScaledXSTable" << G4endl;
270 ed << "Unable to retrieve data for element " << iZ << G4endl;
271 G4Exception("G4PenelopeBremsstrahlungFS::BuildScaledXSTable()",
272 "em2009",FatalException,ed);
273 }
274 }
275
276 G4DataVector* atomData = fElementData->find(iZ)->second;
277
278 for (std::size_t ie=0;ie<fNBinsE;++ie)
279 {
280 (*tempData)[ie] += wgt*(*atomData)[ie*(fNBinsX+1)+fNBinsX]; //last column contains total XS
281 for (std::size_t ix=0;ix<fNBinsX;++ix)
282 (*tempMatrix)[ie*fNBinsX+ix] += wgt*(*atomData)[ie*(fNBinsX+1)+ix];
283 }
284 }
285
286 //*********************************************************************
287 // the total energy loss spectrum is re-normalized to reproduce the total
288 // scaled cross section of Berger and Seltzer
289 //*********************************************************************
290 for (std::size_t ie=0;ie<fNBinsE;++ie)
291 {
292 //for each energy, calculate integral of dSigma/dx over dx
293 G4double* tempData2 = new G4double[fNBinsX];
294 for (std::size_t ix=0;ix<fNBinsX;++ix)
295 tempData2[ix] = (*tempMatrix)[ie*fNBinsX+ix];
296 G4double rsum = GetMomentumIntegral(tempData2,1.0,0);
297 delete[] tempData2;
298 G4double fact = millibarn*(theEGrid[ie]+electron_mass_c2)*(1./fine_structure_const)/
299 (classic_electr_radius*classic_electr_radius*(theEGrid[ie]+2.0*electron_mass_c2));
300 G4double fnorm = (*tempData)[ie]/(rsum*fact);
301 G4double TST = 100.*std::fabs(fnorm-1.0);
302 if (TST > 1.0)
303 {
305 ed << "G4PenelopeBremsstrahlungFS. Corrupted data files?" << G4endl;
306 G4cout << "TST= " << TST << "; fnorm = " << fnorm << G4endl;
307 G4cout << "rsum = " << rsum << G4endl;
308 G4cout << "fact = " << fact << G4endl;
309 G4cout << ie << " " << theEGrid[ie]/keV << " " << (*tempData)[ie]/barn << G4endl;
310 G4Exception("G4PenelopeBremsstrahlungFS::BuildScaledXSTable()",
311 "em2010",FatalException,ed);
312 }
313 for (std::size_t ix=0;ix<fNBinsX;++ix)
314 (*tempMatrix)[ie*fNBinsX+ix] *= fnorm;
315 }
316
317 //*********************************************************************
318 // create and fill the tables
319 //*********************************************************************
320 G4PhysicsTable* thePhysicsTable = new G4PhysicsTable();
321 // the table will contain 32 G4PhysicsFreeVectors with different
322 // values of x. Each of the G4PhysicsFreeVectors has a profile of
323 // log(XS) vs. log(E)
324
325 //reserve space of the vectors. Everything is log-log
326 //I add one extra "fake" point at low energy, since the Penelope
327 //table starts at 1 keV
328 for (std::size_t i=0;i<fNBinsX;++i)
329 thePhysicsTable->push_back(new G4PhysicsFreeVector(fNBinsE+1));
330
331 for (std::size_t ix=0;ix<fNBinsX;++ix)
332 {
333 G4PhysicsFreeVector* theVec =
334 (G4PhysicsFreeVector*) ((*thePhysicsTable)[ix]);
335 for (std::size_t ie=0;ie<fNBinsE;++ie)
336 {
337 G4double logene = G4Log(theEGrid[ie]);
338 G4double aValue = (*tempMatrix)[ie*fNBinsX+ix];
339 if (aValue < 1e-20*millibarn) //protection against log(0)
340 aValue = 1e-20*millibarn;
341 theVec->PutValues(ie+1,logene,G4Log(aValue));
342 }
343 //Add fake point at 1 eV using an extrapolation with the derivative
344 //at the first valid point (Penelope approach)
345 G4double derivative = ((*theVec)[2]-(*theVec)[1])/(theVec->Energy(2) - theVec->Energy(1));
346 G4double log1eV = G4Log(1*eV);
347 G4double val1eV = (*theVec)[1]+derivative*(log1eV-theVec->Energy(1));
348 //fake point at very low energy
349 theVec->PutValues(0,log1eV,val1eV);
350 }
351 std::pair<const G4Material*,G4double> theKey = std::make_pair(material,cut);
352 fReducedXSTable->insert(std::make_pair(theKey,thePhysicsTable));
353
354 delete StechiometricFactors;
355 delete tempData;
356 delete tempMatrix;
357
358 //Do here also the initialization of the energy sampling
359 if (!(fSamplingTable->count(theKey)))
360 InitializeEnergySampling(material,cut);
361
362 return;
363}
364
365//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
366
367void G4PenelopeBremsstrahlungFS::ReadDataFile(G4int Z)
368{
369 const char* path = G4FindDataDir("G4LEDATA");
370 if (!path)
371 {
372 G4String excep = "G4PenelopeBremsstrahlungFS - G4LEDATA environment variable not set!";
373 G4Exception("G4PenelopeBremsstrahlungFS::ReadDataFile()",
374 "em0006",FatalException,excep);
375 return;
376 }
377 /*
378 Read the cross section file
379 */
380 std::ostringstream ost;
381 if (Z>9)
382 ost << path << "/penelope/bremsstrahlung/pdebr" << Z << ".p08";
383 else
384 ost << path << "/penelope/bremsstrahlung/pdebr0" << Z << ".p08";
385 std::ifstream file(ost.str().c_str());
386 if (!file.is_open())
387 {
388 G4String excep = "G4PenelopeBremsstrahlungFS - data file " +
389 G4String(ost.str()) + " not found!";
390 G4Exception("G4PenelopeBremsstrahlungFS::ReadDataFile()",
391 "em0003",FatalException,excep);
392 return;
393 }
394
395 G4int readZ =0;
396 file >> readZ;
397
398 //check the right file is opened.
399 if (readZ != Z)
400 {
402 ed << "Corrupted data file for Z=" << Z << G4endl;
403 G4Exception("G4PenelopeBremsstrahlungFS::ReadDataFile()",
404 "em0005",FatalException,ed);
405 return;
406 }
407
408 G4DataVector* theMatrix = new G4DataVector(fNBinsE*(fNBinsX+1),0.); //initialized with zeros
409
410 for (std::size_t ie=0;ie<fNBinsE;++ie)
411 {
412 G4double myDouble = 0;
413 file >> myDouble; //energy (eV)
414 if (!theEGrid[ie]) //fill only the first time
415 theEGrid[ie] = myDouble*eV;
416 //
417 for (std::size_t ix=0;ix<fNBinsX;++ix)
418 {
419 file >> myDouble;
420 (*theMatrix)[ie*(fNBinsX+1)+ix] = myDouble*millibarn;
421 }
422 file >> myDouble; //total cross section
423 (*theMatrix)[ie*(fNBinsX+1)+fNBinsX] = myDouble*millibarn;
424 }
425
426 if (fElementData)
427 fElementData->insert(std::make_pair(Z,theMatrix));
428 else
429 delete theMatrix;
430 file.close();
431 return;
432}
433
434//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
435
437 G4double xup,G4int momOrder) const
438//x is always the gridX
439{
440 //Corresponds to the function RLMOM of Penelope
441 //This method performs the calculation of the integral of (x^momOrder)*y over the interval
442 //from x[0] to xup, obtained by linear interpolation on a table of y.
443 //The independent variable is assumed to take positive values only.
444 //
445 std::size_t size = fNBinsX;
446 const G4double eps = 1e-35;
447
448 //Check that the call is valid
449 if (momOrder<-1 || size<2 || theXGrid[0]<0)
450 {
451 G4Exception("G4PenelopeBremsstrahlungFS::GetMomentumIntegral()",
452 "em2011",FatalException,"Invalid call");
453 }
454
455 for (std::size_t i=1;i<size;++i)
456 {
457 if (theXGrid[i]<0 || theXGrid[i]<theXGrid[i-1])
458 {
460 ed << "Invalid call for bin " << i << G4endl;
461 G4Exception("G4PenelopeBremsstrahlungFS::GetMomentumIntegral()",
462 "em2012",FatalException,ed);
463 }
464 }
465
466 //Compute the integral
467 G4double result = 0;
468 if (xup < theXGrid[0])
469 return result;
470 G4bool loopAgain = true;
471 G4double xt = std::min(xup,theXGrid[size-1]);
472 G4double xtc = 0;
473 for (std::size_t i=0;i<size-1;++i)
474 {
475 G4double x1 = std::max(theXGrid[i],eps);
476 G4double y1 = y[i];
477 G4double x2 = std::max(theXGrid[i+1],eps);
478 G4double y2 = y[i+1];
479 if (xt < x2)
480 {
481 xtc = xt;
482 loopAgain = false;
483 }
484 else
485 xtc = x2;
486 G4double dx = x2-x1;
487 G4double dy = y2-y1;
488 G4double ds = 0;
489 if (std::fabs(dx)>1e-14*std::fabs(dy))
490 {
491 G4double b=dy/dx;
492 G4double a=y1-b*x1;
493 if (momOrder == -1)
494 ds = a*G4Log(xtc/x1)+b*(xtc-x1);
495 else if (momOrder == 0) //speed it up, not using pow()
496 ds = a*(xtc-x1) + 0.5*b*(xtc*xtc-x1*x1);
497 else
498 ds = a*(std::pow(xtc,momOrder+1)-std::pow(x1,momOrder+1))/((G4double) (momOrder + 1))
499 + b*(std::pow(xtc,momOrder+2)-std::pow(x1,momOrder+2))/((G4double) (momOrder + 2));
500 }
501 else
502 ds = 0.5*(y1+y2)*(xtc-x1)*std::pow(xtc,momOrder);
503 result += ds;
504 if (!loopAgain)
505 return result;
506 }
507 return result;
508}
509
510//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
511
513 const G4double cut) const
514{
515 //check if it already contains the entry
516 std::pair<const G4Material*,G4double> theKey = std::make_pair(mat,cut);
517
518 if (!(fReducedXSTable->count(theKey)))
519 {
520 G4Exception("G4PenelopeBremsstrahlungFS::GetScaledXSTable()",
521 "em2013",FatalException,"Unable to retrieve the cross section table");
522 }
523
524 return fReducedXSTable->find(theKey)->second;
525}
526
527//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
528
529void G4PenelopeBremsstrahlungFS::InitializeEnergySampling(const G4Material* material,
530 G4double cut)
531{
532 if (fVerbosity > 2)
533 G4cout << "Entering in G4PenelopeBremsstrahlungFS::InitializeEnergySampling() for " <<
534 material->GetName() << G4endl;
535
536 //This method should be accessed by the master only
537 std::pair<const G4Material*,G4double> theKey = std::make_pair(material,cut);
538
539 G4PhysicsTable* thePhysicsTable = new G4PhysicsTable();
540 // the table will contain 57 G4PhysicsFreeVectors with different
541 // values of E.
542 G4PhysicsFreeVector* thePBvec = new G4PhysicsFreeVector(fNBinsE);
543
544 //I reserve space of the vectors.
545 for (std::size_t i=0;i<fNBinsE;++i)
546 thePhysicsTable->push_back(new G4PhysicsFreeVector(fNBinsX));
547
548 //Retrieve the table. Must already exist at this point, because this
549 //method is invoked by GetScaledXSTable()
550 if (!(fReducedXSTable->count(theKey)))
551 G4Exception("G4PenelopeBremsstrahlungFS::InitializeEnergySampling()",
552 "em2013",FatalException,"Unable to retrieve the cross section table");
553 G4PhysicsTable* theTableReduced = fReducedXSTable->find(theKey)->second;
554
555 for (std::size_t ie=0;ie<fNBinsE;++ie)
556 {
557 G4PhysicsFreeVector* theVec =
558 (G4PhysicsFreeVector*) ((*thePhysicsTable)[ie]);
559 //Fill the table
560 G4double value = 0; //first value
561 theVec->PutValues(0,theXGrid[0],value);
562 for (std::size_t ix=1;ix<fNBinsX;++ix)
563 {
564 //Here calculate the cumulative distribution
565 // int_{0}^{x} dSigma(x',E)/dx' (1/x') dx'
566 G4PhysicsFreeVector* v1 = (G4PhysicsFreeVector*) (*theTableReduced)[ix-1];
567 G4PhysicsFreeVector* v2 = (G4PhysicsFreeVector*) (*theTableReduced)[ix];
568
569 G4double x1=std::max(theXGrid[ix-1],1.0e-35);
570 //Remember: the table fReducedXSTable has a fake first point in energy
571 //so, it contains one more bin than fNBinsE.
572 G4double y1=G4Exp((*v1)[ie+1]);
573 G4double x2=std::max(theXGrid[ix],1.0e-35);
574 G4double y2=G4Exp((*v2)[ie+1]);
575 G4double B = (y2-y1)/(x2-x1);
576 G4double A = y1-B*x1;
577 G4double dS = A*G4Log(x2/x1)+B*(x2-x1);
578 value += dS;
579 theVec->PutValues(ix,theXGrid[ix],value);
580 }
581 //fill the PB vector
582 G4double xc = cut/theEGrid[ie];
583 //Fill a temp data vector
584 G4double* tempData = new G4double[fNBinsX];
585 for (std::size_t ix=0;ix<fNBinsX;++ix)
586 {
587 G4PhysicsFreeVector* vv = (G4PhysicsFreeVector*) (*theTableReduced)[ix];
588 tempData[ix] = G4Exp((*vv)[ie+1]);
589 }
590 G4double pbval = (xc<=1) ?
591 GetMomentumIntegral(tempData,xc,-1) :
592 GetMomentumIntegral(tempData,1.0,-1);
593 thePBvec->PutValues(ie,theEGrid[ie],pbval);
594 delete[] tempData;
595 }
596
597 fSamplingTable->insert(std::make_pair(theKey,thePhysicsTable));
598 fPBcut->insert(std::make_pair(theKey,thePBvec));
599 return;
600}
601
602//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
603
605 const G4double cut) const
606{
607 std::pair<const G4Material*,G4double> theKey = std::make_pair(mat,cut);
608 if (!(fSamplingTable->count(theKey)) || !(fPBcut->count(theKey)) ||
609 !(fReducedXSTable->count(theKey)))
610 {
612 ed << "Unable to retrieve the SamplingTable for " << mat->GetName() << G4endl;
613 G4Exception("G4PenelopeBremsstrahlungFS::SampleGammaEnergy()",
614 "em2014",FatalException,ed);
615 return 0.;
616 }
617 const G4PhysicsTable* theTableInte = fSamplingTable->find(theKey)->second;
618 const G4PhysicsTable* theTableRed = fReducedXSTable->find(theKey)->second;
619
620 //Find the energy bin using bi-partition
621 std::size_t eBin = 0;
622 G4bool firstOrLastBin = false;
623
624 if (energy < theEGrid[0]) //below first bin
625 {
626 eBin = 0;
627 firstOrLastBin = true;
628 }
629 else if (energy > theEGrid[fNBinsE-1]) //after last bin
630 {
631 eBin = fNBinsE-1;
632 firstOrLastBin = true;
633 }
634 else
635 {
636 std::size_t i=0;
637 std::size_t j=fNBinsE-1;
638 while ((j-i)>1)
639 {
640 std::size_t k = (i+j)/2;
641 if (energy > theEGrid[k])
642 i = k;
643 else
644 j = k;
645 }
646 eBin = i;
647 }
648
649 //Get the appropriate physics vector
650 const G4PhysicsFreeVector* theVec1 = (G4PhysicsFreeVector*) (*theTableInte)[eBin];
651
652 //Use a "temporary" vector which contains the linear interpolation of the x spectra
653 //in energy. The temporary vector is thread-local, so that there is no conflict.
654 //This is achieved via G4Cache. The theTempVect is allocated only once per thread
655 //(member variable), but it is overwritten at every call of this method
656 //(because the interpolation factors change!)
657 G4PhysicsFreeVector* theTempVec = fCache.Get();
658 if (!theTempVec) //First time this thread gets the cache
659 {
660 theTempVec = new G4PhysicsFreeVector(fNBinsX);
661 fCache.Put(theTempVec);
662 // The G4AutoDelete takes care here to clean up the vectors
663 G4AutoDelete::Register(theTempVec);
664 if (fVerbosity > 4)
665 G4cout << "Creating new instance of G4PhysicsFreeVector() on the worker" << G4endl;
666 }
667
668 //theTempVect is allocated only once (member variable), but it is overwritten at
669 //every call of this method (because the interpolation factors change!)
670 if (!firstOrLastBin)
671 {
672 const G4PhysicsFreeVector* theVec2 = (G4PhysicsFreeVector*) (*theTableInte)[eBin+1];
673 for (std::size_t iloop=0;iloop<fNBinsX;++iloop)
674 {
675 G4double val = (*theVec1)[iloop]+(((*theVec2)[iloop]-(*theVec1)[iloop]))*
676 (energy-theEGrid[eBin])/(theEGrid[eBin+1]-theEGrid[eBin]);
677 theTempVec->PutValues(iloop,theXGrid[iloop],val);
678 }
679 }
680 else //first or last bin, no interpolation
681 {
682 for (std::size_t iloop=0;iloop<fNBinsX;++iloop)
683 theTempVec->PutValues(iloop,theXGrid[iloop],(*theVec1)[iloop]);
684 }
685
686 //Start the game
687 G4double pbcut = (*(fPBcut->find(theKey)->second))[eBin];
688
689 if (!firstOrLastBin) //linear interpolation on pbcut as well
690 {
691 pbcut = (*(fPBcut->find(theKey)->second))[eBin] +
692 ((*(fPBcut->find(theKey)->second))[eBin+1]-(*(fPBcut->find(theKey)->second))[eBin])*
693 (energy-theEGrid[eBin])/(theEGrid[eBin+1]-theEGrid[eBin]);
694 }
695
696 G4double pCumulative = (*theTempVec)[fNBinsX-1]; //last value
697
698 G4double eGamma = 0;
699 G4int nIterations = 0;
700 do
701 {
702 G4double pt = pbcut + G4UniformRand()*(pCumulative - pbcut);
703 nIterations++;
704
705 //find where it is
706 std::size_t ibin = 0;
707 if (pt < (*theTempVec)[0])
708 ibin = 0;
709 else if (pt > (*theTempVec)[fNBinsX-1])
710 {
711 //We observed problems due to numerical rounding here (STT).
712 //delta here is a tiny positive number
713 G4double delta = pt-(*theTempVec)[fNBinsX-1];
714 if (delta < pt*1e-10) // very small! Numerical rounding only
715 {
716 ibin = fNBinsX-2;
718 ed << "Found that (pt > (*theTempVec)[fNBinsX-1]) with pt = " << pt <<
719 " , (*theTempVec)[fNBinsX-1] = " << (*theTempVec)[fNBinsX-1] << " and delta = " <<
720 (pt-(*theTempVec)[fNBinsX-1]) << G4endl;
721 ed << "Possible symptom of problem with numerical precision" << G4endl;
722 G4Exception("G4PenelopeBremsstrahlungFS::SampleGammaEnergy()",
723 "em2015",JustWarning,ed);
724 }
725 else //real problem
726 {
728 ed << "Crash at (pt > (*theTempVec)[fNBinsX-1]) with pt = " << pt <<
729 " , (*theTempVec)[fNBinsX-1]=" << (*theTempVec)[fNBinsX-1] << " and fNBinsX = " <<
730 fNBinsX << G4endl;
731 ed << "Material: " << mat->GetName() << ", energy = " << energy/keV << " keV" <<
732 G4endl;
733 G4Exception("G4PenelopeBremsstrahlungFS::SampleGammaEnergy()",
734 "em2015",FatalException,ed);
735 }
736 }
737 else
738 {
739 std::size_t i=0;
740 std::size_t j=fNBinsX-1;
741 while ((j-i)>1)
742 {
743 std::size_t k = (i+j)/2;
744 if (pt > (*theTempVec)[k])
745 i = k;
746 else
747 j = k;
748 }
749 ibin = i;
750 }
751
752 G4double w1 = theXGrid[ibin];
753 G4double w2 = theXGrid[ibin+1];
754
755 const G4PhysicsFreeVector* v1 = (G4PhysicsFreeVector*) (*theTableRed)[ibin];
756 const G4PhysicsFreeVector* v2 = (G4PhysicsFreeVector*) (*theTableRed)[ibin+1];
757 //Remember: the table fReducedXSTable has a fake first point in energy
758 //so, it contains one more bin than fNBinsE.
759 G4double pdf1 = G4Exp((*v1)[eBin+1]);
760 G4double pdf2 = G4Exp((*v2)[eBin+1]);
761 G4double deltaW = w2-w1;
762 G4double dpdfb = pdf2-pdf1;
763 G4double B = dpdfb/deltaW;
764 G4double A = pdf1-B*w1;
765 //I already made an interpolation in energy, so I can use the actual value for the
766 //calculation of the wbcut, instead of the grid values (except for the last bin)
767 G4double wbcut = (cut < energy) ? cut/energy : 1.0;
768 if (firstOrLastBin) //this is an particular case: no interpolation available
769 wbcut = (cut < theEGrid[eBin]) ? cut/theEGrid[eBin] : 1.0;
770
771 if (w1 < wbcut)
772 w1 = wbcut;
773 if (w2 < w1)
774 {
775 //This configuration can happen if initially wbcut > w2 > w1. Due to the previous
776 //statement, (w1 = wbcut), it becomes wbcut = w1 > w2. In this case, it is not a
777 //real problem. It becomes a problem if w2 < w1 before the w1 = wbcut statement. Issue
778 //a warning only in this specific case.
779 if (w2 > wbcut)
780 {
782 ed << "Warning in G4PenelopeBremsstrahlungFS::SampleX()" << G4endl;
783 ed << "Conflicting end-point values: w1=" << w1 << "; w2 = " << w2 << G4endl;
784 ed << "wbcut = " << wbcut << " energy= " << energy/keV << " keV" << G4endl;
785 ed << "cut = " << cut/keV << " keV" << G4endl;
786 G4Exception("G4PenelopeBremsstrahlungFS::SampleGammaEnergy()","em2015",
787 JustWarning,ed);
788 }
789 return w1*energy;
790 }
791
792 G4double pmax = std::max(A+B*w1,A+B*w2);
793 G4bool loopAgain = false;
794 do
795 {
796 loopAgain = false;
797 eGamma = w1* std::pow((w2/w1),G4UniformRand());
798 if (G4UniformRand()*pmax > (A+B*eGamma))
799 loopAgain = true;
800 }while(loopAgain);
801 eGamma *= energy;
802 if (nIterations > 100) //protection against infinite loops
803 return eGamma;
804 }while(eGamma < cut); //repeat if sampled sub-cut!
805
806 return eGamma;
807}
G4double B(G4double temperature)
std::vector< const G4Element * > G4ElementVector
const char * G4FindDataDir(const char *)
@ 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
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
const G4ElementVector * GetElementVector() const
const G4double * GetFractionVector() const
std::size_t GetNumberOfElements() const
const G4String & GetName() const
void BuildScaledXSTable(const G4Material *material, G4double cut, G4bool isMaster)
G4double GetEffectiveZSquared(const G4Material *mat) const
const G4PhysicsTable * GetScaledXSTable(const G4Material *, const G4double cut) const
void ClearTables(G4bool isMaster=true)
Reserved for the master model: they build and handle tables.
G4double SampleGammaEnergy(G4double energy, const G4Material *, const G4double cut) const
G4double GetMomentumIntegral(G4double *y, G4double up, G4int momOrder) const
G4PenelopeBremsstrahlungFS(G4int verbosity=0)
Only master models are supposed to create instances.
void PutValues(const std::size_t index, const G4double energy, const G4double value)
void push_back(G4PhysicsVector *)
void clearAndDestroy()
G4double Energy(const std::size_t index) const
void Register(T *inst)