Geant4 10.7.0
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=0, 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=0)
 
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 ()
 

Protected Member Functions

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

Detailed Description

Definition at line 63 of file G4VCrossSectionHandler.hh.

Constructor & Destructor Documentation

◆ G4VCrossSectionHandler() [1/2]

G4VCrossSectionHandler::G4VCrossSectionHandler ( )

Definition at line 87 of file G4VCrossSectionHandler.cc.

88{
89 crossSections = 0;
90 interpolation = 0;
91 Initialise();
93}
void Initialise(G4VDataSetAlgorithm *interpolation=0, 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/2]

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 
)

Definition at line 96 of file G4VCrossSectionHandler.cc.

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

◆ ~G4VCrossSectionHandler()

G4VCrossSectionHandler::~G4VCrossSectionHandler ( )
virtual

Definition at line 111 of file G4VCrossSectionHandler.cc.

112{
113 delete interpolation;
114 interpolation = 0;
115 std::map<G4int,G4VEMDataSet*,std::less<G4int> >::iterator pos;
116
117 for (pos = dataMap.begin(); pos != dataMap.end(); ++pos)
118 {
119 // The following is a workaround for STL ObjectSpace implementation,
120 // which does not support the standard and does not accept
121 // the syntax pos->second
122 // G4VEMDataSet* dataSet = pos->second;
123 G4VEMDataSet* dataSet = (*pos).second;
124 delete dataSet;
125 }
126
127 if (crossSections != 0)
128 {
129 size_t n = crossSections->size();
130 for (size_t i=0; i<n; i++)
131 {
132 delete (*crossSections)[i];
133 }
134 delete crossSections;
135 crossSections = 0;
136 }
137}

Member Function Documentation

◆ ActiveElements()

void G4VCrossSectionHandler::ActiveElements ( )
protected

Definition at line 694 of file G4VCrossSectionHandler.cc.

695{
696 const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
697 if (materialTable == 0)
698 G4Exception("G4VCrossSectionHandler::ActiveElements",
699 "em1001",FatalException,"no MaterialTable found");
700
702
703 for (G4int mLocal2=0; mLocal2<nMaterials; mLocal2++)
704 {
705 const G4Material* material= (*materialTable)[mLocal2];
706 const G4ElementVector* elementVector = material->GetElementVector();
707 const G4int nElements = material->GetNumberOfElements();
708
709 for (G4int iEl=0; iEl<nElements; iEl++)
710 {
711 G4Element* element = (*elementVector)[iEl];
712 G4double Z = element->GetZ();
713 if (!(activeZ.contains(Z)) && Z >= zMin && Z <= zMax)
714 {
715 activeZ.push_back(Z);
716 }
717 }
718 }
719}
std::vector< G4Element * > G4ElementVector
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::vector< G4Material * > G4MaterialTable
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4bool contains(const G4double &) const
G4double GetZ() const
Definition: G4Element.hh:130
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:644
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:637

Referenced by Clear(), and G4VCrossSectionHandler().

◆ BuildCrossSectionsForMaterials()

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

◆ BuildMeanFreePathForMaterials()

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

Definition at line 462 of file G4VCrossSectionHandler.cc.

463{
464 // Builds a CompositeDataSet containing the mean free path for each material
465 // in the material table
466
467 G4DataVector energyVector;
468 G4double dBin = std::log10(eMax/eMin) / nBins;
469
470 for (G4int i=0; i<nBins+1; i++)
471 {
472 energyVector.push_back(std::pow(10., std::log10(eMin)+i*dBin));
473 }
474
475 // Factory method to build cross sections in derived classes,
476 // related to the type of physics process
477
478 if (crossSections != 0)
479 { // Reset the list of cross sections
480 std::vector<G4VEMDataSet*>::iterator mat;
481 if (! crossSections->empty())
482 {
483 for (mat = crossSections->begin(); mat!= crossSections->end(); ++mat)
484 {
485 G4VEMDataSet* set = *mat;
486 delete set;
487 set = 0;
488 }
489 crossSections->clear();
490 delete crossSections;
491 crossSections = 0;
492 }
493 }
494
495 crossSections = BuildCrossSectionsForMaterials(energyVector,energyCuts);
496
497 if (crossSections == 0)
498 {
499 G4Exception("G4VCrossSectionHandler::BuildMeanFreePathForMaterials",
500 "em1010",FatalException,"crossSections = 0");
501 return 0;
502 }
503
505 G4VEMDataSet* materialSet = new G4CompositeEMDataSet(algo);
506 //G4cout << "G4VCrossSectionHandler new dataset " << materialSet << G4endl;
507
508 G4DataVector* energies;
509 G4DataVector* data;
510 G4DataVector* log_energies;
511 G4DataVector* log_data;
512
513
514 const G4ProductionCutsTable* theCoupleTable=
516 size_t numOfCouples = theCoupleTable->GetTableSize();
517
518
519 for (size_t mLocal=0; mLocal<numOfCouples; mLocal++)
520 {
521 energies = new G4DataVector;
522 data = new G4DataVector;
523 log_energies = new G4DataVector;
524 log_data = new G4DataVector;
525 for (G4int bin=0; bin<nBins; bin++)
526 {
527 G4double energy = energyVector[bin];
528 energies->push_back(energy);
529 log_energies->push_back(std::log10(energy));
530 G4VEMDataSet* matCrossSet = (*crossSections)[mLocal];
531 G4double materialCrossSection = 0.0;
532 G4int nElm = matCrossSet->NumberOfComponents();
533 for(G4int j=0; j<nElm; j++) {
534 materialCrossSection += matCrossSet->GetComponent(j)->FindValue(energy);
535 }
536
537 if (materialCrossSection > 0.)
538 {
539 data->push_back(1./materialCrossSection);
540 log_data->push_back(std::log10(1./materialCrossSection));
541 }
542 else
543 {
544 data->push_back(DBL_MAX);
545 log_data->push_back(std::log10(DBL_MAX));
546 }
547 }
549
550 //G4VEMDataSet* dataSet = new G4EMDataSet(m,energies,data,algo,1.,1.);
551
552 G4VEMDataSet* dataSet = new G4EMDataSet(mLocal,energies,data,log_energies,log_data,algoLocal,1.,1.);
553
554 materialSet->AddComponent(dataSet);
555 }
556
557 return materialSet;
558}
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
virtual std::vector< G4VEMDataSet * > * BuildCrossSectionsForMaterials(const G4DataVector &energyVector, const G4DataVector *energyCuts=0)=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 352 of file G4VCrossSectionHandler.cc.

353{
354 // Reset the map of data sets: remove the data sets from the map
355 std::map<G4int,G4VEMDataSet*,std::less<G4int> >::iterator pos;
356
357 if(! dataMap.empty())
358 {
359 for (pos = dataMap.begin(); pos != dataMap.end(); ++pos)
360 {
361 // The following is a workaround for STL ObjectSpace implementation,
362 // which does not support the standard and does not accept
363 // the syntax pos->first or pos->second
364 // G4VEMDataSet* dataSet = pos->second;
365 G4VEMDataSet* dataSet = (*pos).second;
366 delete dataSet;
367 dataSet = 0;
368 G4int i = (*pos).first;
369 dataMap[i] = 0;
370 }
371 dataMap.clear();
372 }
373
374 activeZ.clear();
376}

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

◆ CreateInterpolation()

G4VDataSetAlgorithm * G4VCrossSectionHandler::CreateInterpolation ( )
protectedvirtual

◆ FindValue() [1/2]

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

Definition at line 378 of file G4VCrossSectionHandler.cc.

379{
380 G4double value = 0.;
381
382 std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
383 pos = dataMap.find(Z);
384 if (pos!= dataMap.end())
385 {
386 // The following is a workaround for STL ObjectSpace implementation,
387 // which does not support the standard and does not accept
388 // the syntax pos->first or pos->second
389 // G4VEMDataSet* dataSet = pos->second;
390 G4VEMDataSet* dataSet = (*pos).second;
391 value = dataSet->FindValue(energy);
392 }
393 else
394 {
395 G4cout << "WARNING: G4VCrossSectionHandler::FindValue did not find Z = "
396 << Z << G4endl;
397 }
398 return value;
399}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

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

◆ FindValue() [2/2]

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

Definition at line 401 of file G4VCrossSectionHandler.cc.

403{
404 G4double value = 0.;
405
406 std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
407 pos = dataMap.find(Z);
408 if (pos!= dataMap.end())
409 {
410 // The following is a workaround for STL ObjectSpace implementation,
411 // which does not support the standard and does not accept
412 // the syntax pos->first or pos->second
413 // G4VEMDataSet* dataSet = pos->second;
414 G4VEMDataSet* dataSet = (*pos).second;
415 if (shellIndex >= 0)
416 {
417 G4int nComponents = dataSet->NumberOfComponents();
418 if(shellIndex < nComponents)
419 // - MGP - Why doesn't it use G4VEMDataSet::FindValue directly?
420 value = dataSet->GetComponent(shellIndex)->FindValue(energy);
421 else
422 {
423 G4cout << "WARNING: G4VCrossSectionHandler::FindValue did not find"
424 << " shellIndex= " << shellIndex
425 << " for Z= "
426 << Z << G4endl;
427 }
428 } else {
429 value = dataSet->FindValue(energy);
430 }
431 }
432 else
433 {
434 G4cout << "WARNING: G4VCrossSectionHandler::FindValue did not find Z = "
435 << Z << G4endl;
436 }
437 return value;
438}

◆ GetInterpolation()

const G4VDataSetAlgorithm * G4VCrossSectionHandler::GetInterpolation ( ) const
inlineprotected

Definition at line 125 of file G4VCrossSectionHandler.hh.

125{ return interpolation; }

◆ Initialise()

void G4VCrossSectionHandler::Initialise ( G4VDataSetAlgorithm interpolation = 0,
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 139 of file G4VCrossSectionHandler.cc.

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

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

◆ LoadData()

void G4VCrossSectionHandler::LoadData ( const G4String dataFile)

Definition at line 186 of file G4VCrossSectionHandler.cc.

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

Referenced by G4LivermoreComptonModifiedModel::Initialise().

◆ LoadNonLogData()

void G4VCrossSectionHandler::LoadNonLogData ( const G4String dataFile)

Definition at line 263 of file G4VCrossSectionHandler.cc.

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

◆ LoadShellData()

void G4VCrossSectionHandler::LoadShellData ( const G4String dataFile)

Definition at line 333 of file G4VCrossSectionHandler.cc.

334{
335 size_t nZ = activeZ.size();
336 for (size_t i=0; i<nZ; i++)
337 {
338 G4int Z = (G4int) activeZ[i];
339
340 G4VDataSetAlgorithm* algo = interpolation->Clone();
341 G4VEMDataSet* dataSet = new G4ShellEMDataSet(Z, algo);
342
343 dataSet->LoadData(fileName);
344
345 dataMap[Z] = dataSet;
346 }
347}
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 727 of file G4VCrossSectionHandler.cc.

728{
729 G4int n = 0;
730
731 std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
732 pos = dataMap.find(Z);
733 if (pos!= dataMap.end())
734 {
735 G4VEMDataSet* dataSet = (*pos).second;
736 n = dataSet->NumberOfComponents();
737 }
738 else
739 {
740 G4cout << "WARNING: G4VCrossSectionHandler::NumberOfComponents did not "
741 << "find Z = "
742 << Z << G4endl;
743 }
744 return n;
745}

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

◆ PrintData()

void G4VCrossSectionHandler::PrintData ( ) const

Definition at line 165 of file G4VCrossSectionHandler.cc.

166{
167 std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
168
169 for (pos = dataMap.begin(); pos != dataMap.end(); pos++)
170 {
171 // The following is a workaround for STL ObjectSpace implementation,
172 // which does not support the standard and does not accept
173 // the syntax pos->first or pos->second
174 // G4int z = pos->first;
175 // G4VEMDataSet* dataSet = pos->second;
176 G4int z = (*pos).first;
177 G4VEMDataSet* dataSet = (*pos).second;
178 G4cout << "---- Data set for Z = "
179 << z
180 << G4endl;
181 dataSet->PrintData();
182 G4cout << "--------------------------------------------------" << G4endl;
183 }
184}
virtual void PrintData(void) const =0

Referenced by G4LivermoreIonisationModel::Initialise().

◆ SelectRandomAtom()

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

Definition at line 561 of file G4VCrossSectionHandler.cc.

563{
564 // Select randomly an element within the material, according to the weight
565 // determined by the cross sections in the data set
566
567 const G4Material* material = couple->GetMaterial();
568 G4int nElements = material->GetNumberOfElements();
569
570 // Special case: the material consists of one element
571 if (nElements == 1)
572 {
573 G4int Z = (G4int) material->GetZ();
574 return Z;
575 }
576
577 // Composite material
578
579 const G4ElementVector* elementVector = material->GetElementVector();
580 size_t materialIndex = couple->GetIndex();
581
582 G4VEMDataSet* materialSet = (*crossSections)[materialIndex];
583 G4double materialCrossSection0 = 0.0;
584 G4DataVector cross;
585 cross.clear();
586 for ( G4int i=0; i < nElements; i++ )
587 {
588 G4double cr = materialSet->GetComponent(i)->FindValue(e);
589 materialCrossSection0 += cr;
590 cross.push_back(materialCrossSection0);
591 }
592
593 G4double random = G4UniformRand() * materialCrossSection0;
594
595 for (G4int k=0 ; k < nElements ; k++ )
596 {
597 if (random <= cross[k]) return (G4int) (*elementVector)[k]->GetZ();
598 }
599 // It should never get here
600 return 0;
601}
#define G4UniformRand()
Definition: Randomize.hh:52
const G4Material * GetMaterial() const
G4double GetZ() const
Definition: G4Material.cc:701

Referenced by G4LivermoreIonisationModel::SampleSecondaries().

◆ SelectRandomElement()

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

Definition at line 603 of file G4VCrossSectionHandler.cc.

605{
606 // Select randomly an element within the material, according to the weight determined
607 // by the cross sections in the data set
608
609 const G4Material* material = couple->GetMaterial();
610 G4Element* nullElement = 0;
611 G4int nElements = material->GetNumberOfElements();
612 const G4ElementVector* elementVector = material->GetElementVector();
613
614 // Special case: the material consists of one element
615 if (nElements == 1)
616 {
617 G4Element* element = (*elementVector)[0];
618 return element;
619 }
620 else
621 {
622 // Composite material
623
624 size_t materialIndex = couple->GetIndex();
625
626 G4VEMDataSet* materialSet = (*crossSections)[materialIndex];
627 G4double materialCrossSection0 = 0.0;
628 G4DataVector cross;
629 cross.clear();
630 for (G4int i=0; i<nElements; i++)
631 {
632 G4double cr = materialSet->GetComponent(i)->FindValue(e);
633 materialCrossSection0 += cr;
634 cross.push_back(materialCrossSection0);
635 }
636
637 G4double random = G4UniformRand() * materialCrossSection0;
638
639 for (G4int k=0 ; k < nElements ; k++ )
640 {
641 if (random <= cross[k]) return (*elementVector)[k];
642 }
643 // It should never end up here
644 G4cout << "G4VCrossSectionHandler::SelectRandomElement - no element found" << G4endl;
645 return nullElement;
646 }
647}
size_t GetIndex() const
Definition: G4Element.hh:181

◆ SelectRandomShell()

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

Definition at line 649 of file G4VCrossSectionHandler.cc.

650{
651 // Select randomly a shell, according to the weight determined by the cross sections
652 // in the data set
653
654 // Note for later improvement: it would be useful to add a cache mechanism for already
655 // used shells to improve performance
656
657 G4int shell = 0;
658
659 G4double totCrossSection = FindValue(Z,e);
660 G4double random = G4UniformRand() * totCrossSection;
661 G4double partialSum = 0.;
662
663 G4VEMDataSet* dataSet = 0;
664 std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
665 pos = dataMap.find(Z);
666 // The following is a workaround for STL ObjectSpace implementation,
667 // which does not support the standard and does not accept
668 // the syntax pos->first or pos->second
669 // if (pos != dataMap.end()) dataSet = pos->second;
670 if (pos != dataMap.end())
671 dataSet = (*pos).second;
672 else
673 {
674 G4Exception("G4VCrossSectionHandler::SelectRandomShell",
675 "em1011",FatalException,"unable to load the dataSet");
676 return 0;
677 }
678
679 size_t nShells = dataSet->NumberOfComponents();
680 for (size_t i=0; i<nShells; i++)
681 {
682 const G4VEMDataSet* shellDataSet = dataSet->GetComponent(i);
683 if (shellDataSet != 0)
684 {
685 G4double value = shellDataSet->FindValue(e);
686 partialSum += value;
687 if (random <= partialSum) return i;
688 }
689 }
690 // It should never get here
691 return shell;
692}
G4double FindValue(G4int Z, G4double e) const

Referenced by G4LivermoreIonisationModel::SampleSecondaries().

◆ ValueForMaterial()

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

Definition at line 441 of file G4VCrossSectionHandler.cc.

443{
444 G4double value = 0.;
445
446 const G4ElementVector* elementVector = material->GetElementVector();
447 const G4double* nAtomsPerVolume = material->GetVecNbOfAtomsPerVolume();
448 G4int nElements = material->GetNumberOfElements();
449
450 for (G4int i=0 ; i<nElements ; i++)
451 {
452 G4int Z = (G4int) (*elementVector)[i]->GetZ();
453 G4double elementValue = FindValue(Z,energy);
454 G4double nAtomsVol = nAtomsPerVolume[i];
455 value += nAtomsVol * elementValue;
456 }
457
458 return value;
459}
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204

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