Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VCrossSectionHandler.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//
28// Author: Maria Grazia Pia (Maria.Grazia.Pia@cern.ch)
29//
30// History:
31// -----------
32// 16 Sep 2001 MGP Created
33// 26 Sep 2001 V.Ivanchenko Hide copy constructor and assignement operator
34// 18 Apr 2002 V.Ivanchenko Move member function ValueForMaterial to public
35// 21 Jan 2003 V.Ivanchenko Cut per region
36// 15 Jul 2009 N.A.Karakatsanis New methods added for loading logarithmic data
37// to enhance computing performance of interpolation
38//
39// -------------------------------------------------------------------
40
41// Class description:
42// Low Energy Electromagnetic Physics
43// Base class for cross section manager for an electromagnetic physics process
44
45// -------------------------------------------------------------------
46
47#ifndef G4VCROSSSECTIONHANDLER_HH
48#define G4VCROSSSECTIONHANDLER_HH 1
49
50#include <map>
51#include <vector>
53
54#include "globals.hh"
55#include "G4DataVector.hh"
57
59class G4VEMDataSet;
60class G4Material;
61class G4Element;
62
64
65public:
66 explicit G4VCrossSectionHandler();
67 explicit G4VCrossSectionHandler(G4VDataSetAlgorithm* interpolation,
68 G4double minE = 250*CLHEP::eV,
69 G4double maxE = 100*CLHEP::GeV,
70 G4int nBins = 200,
71 G4double unitE = CLHEP::MeV,
72 G4double unitData = CLHEP::barn,
73 G4int minZ = 1, G4int maxZ = 99);
74
76
77 void Initialise(G4VDataSetAlgorithm* interpolation = nullptr,
78 G4double minE = 250*CLHEP::eV,
79 G4double maxE = 100*CLHEP::GeV,
80 G4int numberOfBins = 200,
81 G4double unitE = CLHEP::MeV,
82 G4double unitData = CLHEP::barn,
83 G4int minZ = 1, G4int maxZ = 99);
84
87 G4double e) const;
88
90 G4VEMDataSet* BuildMeanFreePathForMaterials(const G4DataVector* energyCuts = nullptr);
92 G4double FindValue(G4int Z, G4double e, G4int shellIndex) const;
93 G4double ValueForMaterial(const G4Material* material, G4double e) const;
94 void LoadData(const G4String& dataFile);
95 void LoadNonLogData(const G4String& dataFile);
96 void LoadShellData(const G4String& dataFile);
97 void PrintData() const;
98 void Clear();
99
102
103protected:
104
106 void ActiveElements();
107
108 // Factory method
109 virtual std::vector<G4VEMDataSet*>* BuildCrossSectionsForMaterials(const G4DataVector& energyVector,
110 const G4DataVector* energyCuts = nullptr) = 0;
111
112 // Factory method
114 const G4VDataSetAlgorithm* GetInterpolation() const { return interpolation; }
115
116private:
117 G4VDataSetAlgorithm* interpolation;
118 G4DataVector activeZ;
119
120 std::map<G4int,G4VEMDataSet*,std::less<G4int> > dataMap;
121 std::vector<G4VEMDataSet*>* crossSections;
122
123 G4double eMin;
124 G4double eMax;
125 G4double unit1;
126 G4double unit2;
127
128 G4int zMin;
129 G4int zMax;
130 G4int nBins;
131};
132
133#endif
134
135
136
137
138
139
140
141
142
143
144
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4VDataSetAlgorithm * GetInterpolation() const
G4double ValueForMaterial(const G4Material *material, G4double e) const
void LoadShellData(const G4String &dataFile)
G4VCrossSectionHandler & operator=(const G4VCrossSectionHandler &right)=delete
G4VEMDataSet * BuildMeanFreePathForMaterials(const G4DataVector *energyCuts=nullptr)
G4double FindValue(G4int Z, G4double e) const
G4int NumberOfComponents(G4int Z) const
virtual std::vector< G4VEMDataSet * > * BuildCrossSectionsForMaterials(const G4DataVector &energyVector, const G4DataVector *energyCuts=nullptr)=0
void LoadData(const G4String &dataFile)
void LoadNonLogData(const G4String &dataFile)
G4int SelectRandomAtom(const G4MaterialCutsCouple *couple, G4double e) const
G4VCrossSectionHandler(const G4VCrossSectionHandler &)=delete
G4int SelectRandomShell(G4int Z, G4double e) const
const G4Element * SelectRandomElement(const G4MaterialCutsCouple *material, G4double e) const
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)
virtual G4VDataSetAlgorithm * CreateInterpolation()