Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BremsstrahlungCrossSectionHandler Class Reference

#include <G4BremsstrahlungCrossSectionHandler.hh>

+ Inheritance diagram for G4BremsstrahlungCrossSectionHandler:

Public Member Functions

 G4BremsstrahlungCrossSectionHandler (const G4VEnergySpectrum *spectrum, G4VDataSetAlgorithm *interpolation)
 
 ~G4BremsstrahlungCrossSectionHandler ()
 
G4double GetCrossSectionAboveThresholdForElement (G4double energy, G4double cutEnergy, G4int Z)
 
- Public Member Functions inherited from G4VCrossSectionHandler
 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

std::vector< G4VEMDataSet * > * BuildCrossSectionsForMaterials (const G4DataVector &energyVector, const G4DataVector *energyCuts)
 
- Protected Member Functions inherited from G4VCrossSectionHandler
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 62 of file G4BremsstrahlungCrossSectionHandler.hh.

Constructor & Destructor Documentation

◆ G4BremsstrahlungCrossSectionHandler()

G4BremsstrahlungCrossSectionHandler::G4BremsstrahlungCrossSectionHandler ( const G4VEnergySpectrum spectrum,
G4VDataSetAlgorithm interpolation 
)

Definition at line 73 of file G4BremsstrahlungCrossSectionHandler.cc.

75 : theBR(spec)
76{
77 interp = new G4SemiLogInterpolation();
78}

◆ ~G4BremsstrahlungCrossSectionHandler()

G4BremsstrahlungCrossSectionHandler::~G4BremsstrahlungCrossSectionHandler ( )

Definition at line 82 of file G4BremsstrahlungCrossSectionHandler.cc.

83{
84 delete interp;
85}

Member Function Documentation

◆ BuildCrossSectionsForMaterials()

std::vector< G4VEMDataSet * > * G4BremsstrahlungCrossSectionHandler::BuildCrossSectionsForMaterials ( const G4DataVector energyVector,
const G4DataVector energyCuts 
)
protectedvirtual

Implements G4VCrossSectionHandler.

Definition at line 90 of file G4BremsstrahlungCrossSectionHandler.cc.

92{
93 std::vector<G4VEMDataSet*>* set = new std::vector<G4VEMDataSet*>;
94
95 G4DataVector* energies;
96 G4DataVector* cs;
97
98 G4DataVector* log_energies;
99 G4DataVector* log_cs;
100
101 G4int nOfBins = energyVector.size();
102
103 const G4ProductionCutsTable* theCoupleTable=
105 size_t numOfCouples = theCoupleTable->GetTableSize();
106
107 for (size_t mLocal=0; mLocal<numOfCouples; mLocal++) {
108
109 const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(mLocal);
110 const G4Material* material= couple->GetMaterial();
111 const G4ElementVector* elementVector = material->GetElementVector();
112 const G4double* nAtomsPerVolume = material->GetVecNbOfAtomsPerVolume();
113 G4int nElements = material->GetNumberOfElements();
114
115 G4double tcut = (*energyCuts)[mLocal];
116
117 G4VDataSetAlgorithm* algo = interp->Clone();
118 G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo,1.,1.);
119
120 for (G4int i=0; i<nElements; i++) {
121
122 G4int Z = (G4int) ((*elementVector)[i]->GetZ());
123
124 energies = new G4DataVector;
125 cs = new G4DataVector;
126
127 log_energies = new G4DataVector;
128 log_cs = new G4DataVector;
129
130 G4double density = nAtomsPerVolume[i];
131
132 for (G4int bin=0; bin<nOfBins; bin++) {
133
134 G4double e = energyVector[bin];
135 energies->push_back(e);
136 if (e==0.) e=1e-300;
137 log_energies->push_back(std::log10(e));
138 G4double value = 0.0;
139
140 if(e > tcut) {
141 G4double elemCs = FindValue(Z, e);
142
143 value = theBR->Probability(Z, tcut, e, e);
144
145 value *= elemCs*density;
146 }
147 cs->push_back(value);
148
149 if (value==0.) value=1e-300;
150 log_cs->push_back(std::log10(value));
151 }
152 G4VDataSetAlgorithm* algol = interp->Clone();
153
154 //G4VEMDataSet* elSet = new G4EMDataSet(i,energies,cs,algol,1.,1.);
155
156 G4VEMDataSet* elSet = new G4EMDataSet(i,energies,cs,log_energies,log_cs,algol,1.,1.);
157
158 setForMat->AddComponent(elSet);
159 }
160 set->push_back(setForMat);
161 }
162
163 return set;
164}
std::vector< G4Element * > G4ElementVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double FindValue(G4int Z, G4double e) const
virtual G4VDataSetAlgorithm * Clone() const =0
virtual void AddComponent(G4VEMDataSet *dataSet)=0
virtual G4double Probability(G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const =0

◆ GetCrossSectionAboveThresholdForElement()

G4double G4BremsstrahlungCrossSectionHandler::GetCrossSectionAboveThresholdForElement ( G4double  energy,
G4double  cutEnergy,
G4int  Z 
)

Definition at line 168 of file G4BremsstrahlungCrossSectionHandler.cc.

171{
172 G4double value = 0.;
173 if(energy > cutEnergy)
174 {
175 G4double elemCs = FindValue(Z, energy);
176 value = theBR->Probability(Z,cutEnergy, energy, energy);
177 value *= elemCs;
178 }
179 return value;
180}

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