Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VCrossSectionHandler Class Referenceabstract

#include <G4VCrossSectionHandler.hh>

+ Inheritance diagram for G4VCrossSectionHandler:

Public Member Functions

 G4VCrossSectionHandler ()
 
 G4VCrossSectionHandler (G4VDataSetAlgorithm *interpolation, G4double minE=250 *CLHEP::eV, G4double maxE=100 *CLHEP::GeV, G4int nBins=200, G4double unitE=CLHEP::MeV, G4double unitData=CLHEP::barn, G4int minZ=1, G4int maxZ=99)
 
virtual ~G4VCrossSectionHandler ()
 
void Initialise (G4VDataSetAlgorithm *interpolation=nullptr, G4double minE=250 *CLHEP::eV, G4double maxE=100 *CLHEP::GeV, G4int numberOfBins=200, G4double unitE=CLHEP::MeV, G4double unitData=CLHEP::barn, G4int minZ=1, G4int maxZ=99)
 
G4int SelectRandomAtom (const G4MaterialCutsCouple *couple, G4double e) const
 
const G4ElementSelectRandomElement (const G4MaterialCutsCouple *material, G4double e) const
 
G4int SelectRandomShell (G4int Z, G4double e) const
 
G4VEMDataSetBuildMeanFreePathForMaterials (const G4DataVector *energyCuts=nullptr)
 
G4double FindValue (G4int Z, G4double e) const
 
G4double FindValue (G4int Z, G4double e, G4int shellIndex) const
 
G4double ValueForMaterial (const G4Material *material, G4double e) const
 
void LoadData (const G4String &dataFile)
 
void LoadNonLogData (const G4String &dataFile)
 
void LoadShellData (const G4String &dataFile)
 
void PrintData () const
 
void Clear ()
 
 G4VCrossSectionHandler (const G4VCrossSectionHandler &)=delete
 
G4VCrossSectionHandleroperator= (const G4VCrossSectionHandler &right)=delete
 

Protected Member Functions

G4int NumberOfComponents (G4int Z) const
 
void ActiveElements ()
 
virtual std::vector< G4VEMDataSet * > * BuildCrossSectionsForMaterials (const G4DataVector &energyVector, const G4DataVector *energyCuts=nullptr)=0
 
virtual G4VDataSetAlgorithmCreateInterpolation ()
 
const G4VDataSetAlgorithmGetInterpolation () const
 

Detailed Description

Definition at line 63 of file G4VCrossSectionHandler.hh.

Constructor & Destructor Documentation

◆ G4VCrossSectionHandler() [1/3]

G4VCrossSectionHandler::G4VCrossSectionHandler ( )
explicit

Definition at line 88 of file G4VCrossSectionHandler.cc.

89{
90 crossSections = 0;
91 interpolation = 0;
92 Initialise();
94}
void Initialise(G4VDataSetAlgorithm *interpolation=nullptr, G4double minE=250 *CLHEP::eV, G4double maxE=100 *CLHEP::GeV, G4int numberOfBins=200, G4double unitE=CLHEP::MeV, G4double unitData=CLHEP::barn, G4int minZ=1, G4int maxZ=99)

◆ G4VCrossSectionHandler() [2/3]

G4VCrossSectionHandler::G4VCrossSectionHandler ( G4VDataSetAlgorithm interpolation,
G4double  minE = 250*CLHEP::eV,
G4double  maxE = 100*CLHEP::GeV,
G4int  nBins = 200,
G4double  unitE = CLHEP::MeV,
G4double  unitData = CLHEP::barn,
G4int  minZ = 1,
G4int  maxZ = 99 
)
explicit

Definition at line 98 of file G4VCrossSectionHandler.cc.

106 : interpolation(algorithm), eMin(minE), eMax(maxE),
107 unit1(unitE), unit2(unitData), zMin(minZ), zMax(maxZ),
108 nBins(bins)
109{
110 crossSections = nullptr;
112}

◆ ~G4VCrossSectionHandler()

G4VCrossSectionHandler::~G4VCrossSectionHandler ( )
virtual

Definition at line 116 of file G4VCrossSectionHandler.cc.

117{
118 delete interpolation;
119 interpolation = nullptr;
120
121 for (auto & pos : dataMap)
122 {
123 G4VEMDataSet* dataSet = pos.second;
124 delete dataSet;
125 }
126
127 if (crossSections != nullptr)
128 {
129 std::size_t n = crossSections->size();
130 for (std::size_t i=0; i<n; i++)
131 {
132 delete (*crossSections)[i];
133 }
134 delete crossSections;
135 crossSections = nullptr;
136 }
137}

◆ G4VCrossSectionHandler() [3/3]

G4VCrossSectionHandler::G4VCrossSectionHandler ( const G4VCrossSectionHandler )
delete

Member Function Documentation

◆ ActiveElements()

void G4VCrossSectionHandler::ActiveElements ( )
protected

Definition at line 658 of file G4VCrossSectionHandler.cc.

659{
660 const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
661 if (materialTable == nullptr)
662 G4Exception("G4VCrossSectionHandler::ActiveElements",
663 "em1001",FatalException,"no MaterialTable found");
664
665 std::size_t nMaterials = G4Material::GetNumberOfMaterials();
666
667 for (std::size_t mLocal2=0; mLocal2<nMaterials; ++mLocal2)
668 {
669 const G4Material* material= (*materialTable)[mLocal2];
670 const G4ElementVector* elementVector = material->GetElementVector();
671 const std::size_t nElements = material->GetNumberOfElements();
672
673 for (std::size_t iEl=0; iEl<nElements; ++iEl)
674 {
675 G4double Z = (*elementVector)[iEl]->GetZ();
676 if (!(activeZ.contains(Z)) && Z >= zMin && Z <= zMax)
677 {
678 activeZ.push_back(Z);
679 }
680 }
681 }
682}
std::vector< const G4Element * > G4ElementVector
@ FatalException
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
const G4int Z[17]
G4bool contains(const G4double &) const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:684
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:677

Referenced by Clear(), and G4VCrossSectionHandler().

◆ BuildCrossSectionsForMaterials()

virtual std::vector< G4VEMDataSet * > * G4VCrossSectionHandler::BuildCrossSectionsForMaterials ( const G4DataVector energyVector,
const G4DataVector energyCuts = nullptr 
)
protectedpure virtual

◆ BuildMeanFreePathForMaterials()

G4VEMDataSet * G4VCrossSectionHandler::BuildMeanFreePathForMaterials ( const G4DataVector energyCuts = nullptr)

Definition at line 441 of file G4VCrossSectionHandler.cc.

442{
443 // Builds a CompositeDataSet containing the mean free path for each material
444 // in the material table
445 G4DataVector energyVector;
446 G4double dBin = std::log10(eMax/eMin) / nBins;
447
448 for (G4int i=0; i<nBins+1; i++)
449 {
450 energyVector.push_back(std::pow(10., std::log10(eMin)+i*dBin));
451 }
452
453 // Factory method to build cross sections in derived classes,
454 // related to the type of physics process
455
456 if (crossSections != nullptr)
457 { // Reset the list of cross sections
458 if (! crossSections->empty())
459 {
460 for (auto mat=crossSections->begin(); mat != crossSections->end(); ++mat)
461 {
462 G4VEMDataSet* set = *mat;
463 delete set;
464 set = nullptr;
465 }
466 crossSections->clear();
467 delete crossSections;
468 crossSections = nullptr;
469 }
470 }
471
472 crossSections = BuildCrossSectionsForMaterials(energyVector,energyCuts);
473
474 if (crossSections == nullptr)
475 {
476 G4Exception("G4VCrossSectionHandler::BuildMeanFreePathForMaterials",
477 "em1010",FatalException,"crossSections = 0");
478 return 0;
479 }
480
482 G4VEMDataSet* materialSet = new G4CompositeEMDataSet(algo);
483
484 G4DataVector* energies;
485 G4DataVector* data;
486 G4DataVector* log_energies;
487 G4DataVector* log_data;
488
489 const G4ProductionCutsTable* theCoupleTable=
491 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
492
493 for (G4int mLocal=0; mLocal<numOfCouples; ++mLocal)
494 {
495 energies = new G4DataVector;
496 data = new G4DataVector;
497 log_energies = new G4DataVector;
498 log_data = new G4DataVector;
499 for (G4int bin=0; bin<nBins; bin++)
500 {
501 G4double energy = energyVector[bin];
502 energies->push_back(energy);
503 log_energies->push_back(std::log10(energy));
504 G4VEMDataSet* matCrossSet = (*crossSections)[mLocal];
505 G4double materialCrossSection = 0.0;
506 G4int nElm = (G4int)matCrossSet->NumberOfComponents();
507 for(G4int j=0; j<nElm; ++j) {
508 materialCrossSection += matCrossSet->GetComponent(j)->FindValue(energy);
509 }
510
511 if (materialCrossSection > 0.)
512 {
513 data->push_back(1./materialCrossSection);
514 log_data->push_back(std::log10(1./materialCrossSection));
515 }
516 else
517 {
518 data->push_back(DBL_MAX);
519 log_data->push_back(std::log10(DBL_MAX));
520 }
521 }
523 G4VEMDataSet* dataSet = new G4EMDataSet(mLocal,energies,data,log_energies,log_data,algoLocal,1.,1.);
524 materialSet->AddComponent(dataSet);
525 }
526 return materialSet;
527}
int G4int
Definition: G4Types.hh:85
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
virtual std::vector< G4VEMDataSet * > * BuildCrossSectionsForMaterials(const G4DataVector &energyVector, const G4DataVector *energyCuts=nullptr)=0
virtual G4VDataSetAlgorithm * CreateInterpolation()
virtual const G4VEMDataSet * GetComponent(G4int componentId) const =0
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
virtual void AddComponent(G4VEMDataSet *dataSet)=0
virtual size_t NumberOfComponents(void) const =0
G4double energy(const ThreeVector &p, const G4double m)
#define DBL_MAX
Definition: templates.hh:62

Referenced by G4LivermoreIonisationModel::Initialise().

◆ Clear()

void G4VCrossSectionHandler::Clear ( )

Definition at line 345 of file G4VCrossSectionHandler.cc.

346{
347 // Reset the map of data sets: remove the data sets from the map
348 if(! dataMap.empty())
349 {
350 for (auto & pos : dataMap)
351 {
352 G4VEMDataSet* dataSet = pos.second;
353 delete dataSet;
354 dataSet = nullptr;
355 G4int i = pos.first;
356 dataMap[i] = nullptr;
357 }
358 dataMap.clear();
359 }
360 activeZ.clear();
362}

Referenced by G4LivermoreIonisationCrossSection::Initialise(), and G4LivermoreIonisationModel::Initialise().

◆ CreateInterpolation()

G4VDataSetAlgorithm * G4VCrossSectionHandler::CreateInterpolation ( )
protectedvirtual

◆ FindValue() [1/2]

G4double G4VCrossSectionHandler::FindValue ( G4int  Z,
G4double  e 
) const

Definition at line 366 of file G4VCrossSectionHandler.cc.

367{
368 G4double value = 0.;
369
370 auto pos = dataMap.find(Z);
371 if (pos!= dataMap.end())
372 {
373 G4VEMDataSet* dataSet = (*pos).second;
374 value = dataSet->FindValue(energy);
375 }
376 else
377 {
378 G4cout << "WARNING: G4VCrossSectionHandler::FindValue did not find Z = "
379 << Z << G4endl;
380 }
381 return value;
382}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

Referenced by G4eIonisationCrossSectionHandler::BuildCrossSectionsForMaterials(), G4CrossSectionHandler::BuildCrossSectionsForMaterials(), G4LivermoreIonisationModel::ComputeDEDXPerVolume(), G4LivermoreIonisationCrossSection::CrossSection(), G4eIonisationCrossSectionHandler::GetCrossSectionAboveThresholdForElement(), SelectRandomShell(), and ValueForMaterial().

◆ FindValue() [2/2]

G4double G4VCrossSectionHandler::FindValue ( G4int  Z,
G4double  e,
G4int  shellIndex 
) const

Definition at line 386 of file G4VCrossSectionHandler.cc.

388{
389 G4double value = 0.;
390 auto pos = dataMap.find(Z);
391 if (pos!= dataMap.cend())
392 {
393 G4VEMDataSet* dataSet = (*pos).second;
394 if (shellIndex >= 0)
395 {
396 G4int nComponents = (G4int)dataSet->NumberOfComponents();
397 if(shellIndex < nComponents)
398 value = dataSet->GetComponent(shellIndex)->FindValue(energy);
399 else
400 {
401 G4cout << "WARNING: G4VCrossSectionHandler::FindValue did not find"
402 << " shellIndex= " << shellIndex
403 << " for Z= "
404 << Z << G4endl;
405 }
406 } else {
407 value = dataSet->FindValue(energy);
408 }
409 }
410 else
411 {
412 G4cout << "WARNING: G4VCrossSectionHandler::FindValue did not find Z = "
413 << Z << G4endl;
414 }
415 return value;
416}

◆ GetInterpolation()

const G4VDataSetAlgorithm * G4VCrossSectionHandler::GetInterpolation ( ) const
inlineprotected

Definition at line 114 of file G4VCrossSectionHandler.hh.

114{ return interpolation; }

◆ Initialise()

void G4VCrossSectionHandler::Initialise ( G4VDataSetAlgorithm interpolation = nullptr,
G4double  minE = 250*CLHEP::eV,
G4double  maxE = 100*CLHEP::GeV,
G4int  numberOfBins = 200,
G4double  unitE = CLHEP::MeV,
G4double  unitData = CLHEP::barn,
G4int  minZ = 1,
G4int  maxZ = 99 
)

Definition at line 141 of file G4VCrossSectionHandler.cc.

146{
147 if (algorithm != nullptr)
148 {
149 delete interpolation;
150 interpolation = algorithm;
151 }
152 else
153 {
154 delete interpolation;
155 interpolation = CreateInterpolation();
156 }
157
158 eMin = minE;
159 eMax = maxE;
160 nBins = numberOfBins;
161 unit1 = unitE;
162 unit2 = unitData;
163 zMin = minZ;
164 zMax = maxZ;
165}

Referenced by G4eCrossSectionHandler::G4eCrossSectionHandler(), G4eIonisationCrossSectionHandler::G4eIonisationCrossSectionHandler(), and G4VCrossSectionHandler().

◆ LoadData()

void G4VCrossSectionHandler::LoadData ( const G4String dataFile)

Definition at line 185 of file G4VCrossSectionHandler.cc.

186{
187 std::size_t nZ = activeZ.size();
188 for (std::size_t i=0; i<nZ; ++i)
189 {
190 G4int Z = G4int(activeZ[i]);
191
192 // Build the complete string identifying the file with the data set
193 const char* path = G4FindDataDir("G4LEDATA");
194 if (!path)
195 {
196 G4Exception("G4VCrossSectionHandler::LoadData",
197 "em0006",FatalException,"G4LEDATA environment variable not set");
198 return;
199 }
200
201 std::ostringstream ost;
202 ost << path << '/' << fileName << Z << ".dat";
203 std::ifstream file(ost.str().c_str());
204 std::filebuf* lsdp = file.rdbuf();
205
206 if (! (lsdp->is_open()) )
207 {
208 G4String excep = "data file: " + ost.str() + " not found";
209 G4Exception("G4VCrossSectionHandler::LoadData",
210 "em0003",FatalException,excep);
211 }
212 G4double a = 0;
213 G4int k = 0;
214 G4int nColumns = 2;
215
216 G4DataVector* orig_reg_energies = new G4DataVector;
217 G4DataVector* orig_reg_data = new G4DataVector;
218 G4DataVector* log_reg_energies = new G4DataVector;
219 G4DataVector* log_reg_data = new G4DataVector;
220
221 do
222 {
223 file >> a;
224
225 if (a==0.) a=1e-300;
226
227 // The file is organized into four columns:
228 // 1st column contains the values of energy
229 // 2nd column contains the corresponding data value
230 // The file terminates with the pattern: -1 -1
231 // -2 -2
232 //
233 if (a != -1 && a != -2)
234 {
235 if (k%nColumns == 0)
236 {
237 orig_reg_energies->push_back(a*unit1);
238 log_reg_energies->push_back(std::log10(a)+std::log10(unit1));
239 }
240 else if (k%nColumns == 1)
241 {
242 orig_reg_data->push_back(a*unit2);
243 log_reg_data->push_back(std::log10(a)+std::log10(unit2));
244 }
245 k++;
246 }
247 }
248 while (a != -2); // End of File
249
250 file.close();
251 G4VDataSetAlgorithm* algo = interpolation->Clone();
252 G4VEMDataSet* dataSet = new G4EMDataSet(Z,orig_reg_energies,orig_reg_data,
253 log_reg_energies,log_reg_data,algo);
254 dataMap[Z] = dataSet;
255 }
256}
const char * G4FindDataDir(const char *)
virtual G4VDataSetAlgorithm * Clone() const =0

◆ LoadNonLogData()

void G4VCrossSectionHandler::LoadNonLogData ( const G4String dataFile)

Definition at line 260 of file G4VCrossSectionHandler.cc.

261{
262 std::size_t nZ = activeZ.size();
263 for (std::size_t i=0; i<nZ; ++i)
264 {
265 G4int Z = G4int(activeZ[i]);
266
267 // Build the complete string identifying the file with the data set
268 const char* path = G4FindDataDir("G4LEDATA");
269 if (!path)
270 {
271 G4Exception("G4VCrossSectionHandler::LoadNonLogData",
272 "em0006",FatalException,"G4LEDATA environment variable not set");
273 return;
274 }
275
276 std::ostringstream ost;
277 ost << path << '/' << fileName << Z << ".dat";
278 std::ifstream file(ost.str().c_str());
279 std::filebuf* lsdp = file.rdbuf();
280
281 if (! (lsdp->is_open()) )
282 {
283 G4String excep = "data file: " + ost.str() + " not found";
284 G4Exception("G4VCrossSectionHandler::LoadNonLogData",
285 "em0003",FatalException,excep);
286 }
287 G4double a = 0;
288 G4int k = 0;
289 G4int nColumns = 2;
290
291 G4DataVector* orig_reg_energies = new G4DataVector;
292 G4DataVector* orig_reg_data = new G4DataVector;
293
294 do
295 {
296 file >> a;
297
298 // The file is organized into four columns:
299 // 1st column contains the values of energy
300 // 2nd column contains the corresponding data value
301 // The file terminates with the pattern: -1 -1
302 // -2 -2
303 //
304 if (a != -1 && a != -2)
305 {
306 if (k%nColumns == 0)
307 {
308 orig_reg_energies->push_back(a*unit1);
309 }
310 else if (k%nColumns == 1)
311 {
312 orig_reg_data->push_back(a*unit2);
313 }
314 k++;
315 }
316 }
317 while (a != -2); // End of File
318
319 file.close();
320 G4VDataSetAlgorithm* algo = interpolation->Clone();
321
322 G4VEMDataSet* dataSet = new G4EMDataSet(Z,orig_reg_energies,orig_reg_data,algo);
323 dataMap[Z] = dataSet;
324 }
325}

◆ LoadShellData()

void G4VCrossSectionHandler::LoadShellData ( const G4String dataFile)

Definition at line 329 of file G4VCrossSectionHandler.cc.

330{
331 std::size_t nZ = activeZ.size();
332 for (std::size_t i=0; i<nZ; ++i)
333 {
334 G4int Z = G4int(activeZ[i]);
335
336 G4VDataSetAlgorithm* algo = interpolation->Clone();
337 G4VEMDataSet* dataSet = new G4ShellEMDataSet(Z, algo);
338 dataSet->LoadData(fileName);
339 dataMap[Z] = dataSet;
340 }
341}
virtual G4bool LoadData(const G4String &fileName)=0

Referenced by G4LivermoreIonisationCrossSection::Initialise(), and G4LivermoreIonisationModel::Initialise().

◆ NumberOfComponents()

G4int G4VCrossSectionHandler::NumberOfComponents ( G4int  Z) const
protected

Definition at line 694 of file G4VCrossSectionHandler.cc.

695{
696 G4int n = 0;
697
698 auto pos = dataMap.find(Z);
699 if (pos!= dataMap.end())
700 {
701 G4VEMDataSet* dataSet = (*pos).second;
702 n = (G4int)dataSet->NumberOfComponents();
703 }
704 else
705 {
706 G4cout << "WARNING: G4VCrossSectionHandler::NumberOfComponents did not "
707 << "find Z = "
708 << Z << G4endl;
709 }
710 return n;
711}

Referenced by G4eIonisationCrossSectionHandler::BuildCrossSectionsForMaterials(), and G4eIonisationCrossSectionHandler::GetCrossSectionAboveThresholdForElement().

◆ operator=()

G4VCrossSectionHandler & G4VCrossSectionHandler::operator= ( const G4VCrossSectionHandler right)
delete

◆ PrintData()

void G4VCrossSectionHandler::PrintData ( ) const

Definition at line 169 of file G4VCrossSectionHandler.cc.

170{
171 for (auto & pos : dataMap)
172 {
173 G4int z = pos.first;
174 G4VEMDataSet* dataSet = pos.second;
175 G4cout << "---- Data set for Z = "
176 << z
177 << G4endl;
178 dataSet->PrintData();
179 G4cout << "--------------------------------------------------" << G4endl;
180 }
181}
virtual void PrintData(void) const =0

Referenced by G4LivermoreIonisationModel::Initialise().

◆ SelectRandomAtom()

G4int G4VCrossSectionHandler::SelectRandomAtom ( const G4MaterialCutsCouple couple,
G4double  e 
) const

Definition at line 531 of file G4VCrossSectionHandler.cc.

533{
534 // Select randomly an element within the material, according to the weight
535 // determined by the cross sections in the data set
536 const G4Material* material = couple->GetMaterial();
537 G4int nElements = (G4int)material->GetNumberOfElements();
538
539 // Special case: the material consists of one element
540 if (nElements == 1)
541 {
542 G4int Z = (G4int) material->GetZ();
543 return Z;
544 }
545
546 // Composite material
547 const G4ElementVector* elementVector = material->GetElementVector();
548 std::size_t materialIndex = couple->GetIndex();
549
550 G4VEMDataSet* materialSet = (*crossSections)[materialIndex];
551 G4double materialCrossSection0 = 0.0;
552 G4DataVector cross;
553 cross.clear();
554 for ( G4int i=0; i < nElements; i++ )
555 {
556 G4double cr = materialSet->GetComponent(i)->FindValue(e);
557 materialCrossSection0 += cr;
558 cross.push_back(materialCrossSection0);
559 }
560
561 G4double random = G4UniformRand() * materialCrossSection0;
562 for (G4int k=0 ; k < nElements ; k++ )
563 {
564 if (random <= cross[k]) return (G4int) (*elementVector)[k]->GetZ();
565 }
566 // It should never get here
567 return 0;
568}
#define G4UniformRand()
Definition: Randomize.hh:52
const G4Material * GetMaterial() const
G4double GetZ() const
Definition: G4Material.cc:745

Referenced by G4LivermoreIonisationModel::SampleSecondaries().

◆ SelectRandomElement()

const G4Element * G4VCrossSectionHandler::SelectRandomElement ( const G4MaterialCutsCouple material,
G4double  e 
) const

Definition at line 572 of file G4VCrossSectionHandler.cc.

574{
575 // Select randomly an element within the material, according to the weight determined
576 // by the cross sections in the data set
577 const G4Material* material = couple->GetMaterial();
578 G4Element* nullElement = 0;
579 G4int nElements = (G4int)material->GetNumberOfElements();
580 const G4ElementVector* elementVector = material->GetElementVector();
581
582 // Special case: the material consists of one element
583 if (nElements == 1)
584 {
585 return (*elementVector)[0];
586 }
587 else
588 {
589 // Composite material
590
591 std::size_t materialIndex = couple->GetIndex();
592
593 G4VEMDataSet* materialSet = (*crossSections)[materialIndex];
594 G4double materialCrossSection0 = 0.0;
595 G4DataVector cross;
596 cross.clear();
597 for (G4int i=0; i<nElements; ++i)
598 {
599 G4double cr = materialSet->GetComponent(i)->FindValue(e);
600 materialCrossSection0 += cr;
601 cross.push_back(materialCrossSection0);
602 }
603
604 G4double random = G4UniformRand() * materialCrossSection0;
605
606 for (G4int k=0 ; k < nElements ; ++k )
607 {
608 if (random <= cross[k]) return (*elementVector)[k];
609 }
610 // It should never end up here
611 G4cout << "G4VCrossSectionHandler::SelectRandomElement - no element found" << G4endl;
612 return nullElement;
613 }
614}

◆ SelectRandomShell()

G4int G4VCrossSectionHandler::SelectRandomShell ( G4int  Z,
G4double  e 
) const

Definition at line 618 of file G4VCrossSectionHandler.cc.

619{
620 // Select randomly a shell, according to the weight determined by the cross sections
621 // in the data set
622 // Note for later improvement: it would be useful to add a cache mechanism for already
623 // used shells to improve performance
624 G4int shell = 0;
625
626 G4double totCrossSection = FindValue(Z,e);
627 G4double random = G4UniformRand() * totCrossSection;
628 G4double partialSum = 0.;
629
630 G4VEMDataSet* dataSet = nullptr;
631 auto pos = dataMap.find(Z);
632 if (pos != dataMap.end())
633 dataSet = (*pos).second;
634 else
635 {
636 G4Exception("G4VCrossSectionHandler::SelectRandomShell",
637 "em1011",FatalException,"unable to load the dataSet");
638 return 0;
639 }
640
641 G4int nShells = (G4int)dataSet->NumberOfComponents();
642 for (G4int i=0; i<nShells; ++i)
643 {
644 const G4VEMDataSet* shellDataSet = dataSet->GetComponent(i);
645 if (shellDataSet != nullptr)
646 {
647 G4double value = shellDataSet->FindValue(e);
648 partialSum += value;
649 if (random <= partialSum) return i;
650 }
651 }
652 // It should never get here
653 return shell;
654}
G4double FindValue(G4int Z, G4double e) const

Referenced by G4LivermoreIonisationModel::SampleSecondaries().

◆ ValueForMaterial()

G4double G4VCrossSectionHandler::ValueForMaterial ( const G4Material material,
G4double  e 
) const

Definition at line 420 of file G4VCrossSectionHandler.cc.

422{
423 G4double value = 0.;
424 const G4ElementVector* elementVector = material->GetElementVector();
425 const G4double* nAtomsPerVolume = material->GetVecNbOfAtomsPerVolume();
426 std::size_t nElements = material->GetNumberOfElements();
427
428 for (std::size_t i=0 ; i<nElements ; ++i)
429 {
430 G4int Z = (G4int) (*elementVector)[i]->GetZ();
431 G4double elementValue = FindValue(Z,energy);
432 G4double nAtomsVol = nAtomsPerVolume[i];
433 value += nAtomsVol * elementValue;
434 }
435
436 return value;
437}
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:201

The documentation for this class was generated from the following files: