Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4BremsstrahlungCrossSectionHandler.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// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32//
33// File name: G4BremsstrahlungCrossSectionHandler
34//
35// Author: V.Ivanchenko (Vladimir.Ivanchenko@cern.ch)
36//
37// Creation date: 25 September 2001
38//
39// Modifications:
40//
41// 10.10.2001 MGP Revision to improve code quality and consistency with design
42// 21.01.2003 VI cut per region
43// 03.03.2009 LP Added public method to make a easier migration of
44// G4LowEnergyBremsstrahlung to G4LivermoreBremsstrahlungModel
45//
46// 15 Jul 2009 Nicolas A. Karakatsanis
47//
48// - BuildCrossSectionForMaterials method was revised in order to calculate the
49// logarithmic values of the loaded data.
50// It retrieves the data values from the G4EMLOW data files but, then, calculates the
51// respective log values and loads them to seperate data structures.
52// The EM data sets, initialized this way, contain both non-log and log values.
53// These initialized data sets can enhance the computing performance of data interpolation
54// operations
55//
56//
57//
58//
59// -------------------------------------------------------------------
60
63#include "G4DataVector.hh"
67#include "G4VEMDataSet.hh"
68#include "G4EMDataSet.hh"
69#include "G4Material.hh"
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73
76 : theBR(spec)
77{
78 interp = new G4SemiLogInterpolation();
79}
80
81//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
82
84{
85 delete interp;
86}
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89
90std::vector<G4VEMDataSet*>*
92 const G4DataVector* energyCuts)
93{
94 std::vector<G4VEMDataSet*>* set = new std::vector<G4VEMDataSet*>;
95
96 G4DataVector* energies;
97 G4DataVector* cs;
98
99 G4DataVector* log_energies;
100 G4DataVector* log_cs;
101
102 G4int nOfBins = energyVector.size();
103
104 const G4ProductionCutsTable* theCoupleTable=
106 size_t numOfCouples = theCoupleTable->GetTableSize();
107
108 for (size_t mLocal=0; mLocal<numOfCouples; mLocal++) {
109
110 const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(mLocal);
111 const G4Material* material= couple->GetMaterial();
112 const G4ElementVector* elementVector = material->GetElementVector();
113 const G4double* nAtomsPerVolume = material->GetVecNbOfAtomsPerVolume();
114 G4int nElements = material->GetNumberOfElements();
115
116 G4double tcut = (*energyCuts)[mLocal];
117
118 G4VDataSetAlgorithm* algo = interp->Clone();
119 G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo,1.,1.);
120
121 for (G4int i=0; i<nElements; i++) {
122
123 G4int Z = (G4int) ((*elementVector)[i]->GetZ());
124
125 energies = new G4DataVector;
126 cs = new G4DataVector;
127
128 log_energies = new G4DataVector;
129 log_cs = new G4DataVector;
130
131 G4double density = nAtomsPerVolume[i];
132
133 for (G4int bin=0; bin<nOfBins; bin++) {
134
135 G4double e = energyVector[bin];
136 energies->push_back(e);
137 if (e==0.) e=1e-300;
138 log_energies->push_back(std::log10(e));
139 G4double value = 0.0;
140
141 if(e > tcut) {
142 G4double elemCs = FindValue(Z, e);
143
144 value = theBR->Probability(Z, tcut, e, e);
145
146 value *= elemCs*density;
147 }
148 cs->push_back(value);
149
150 if (value==0.) value=1e-300;
151 log_cs->push_back(std::log10(value));
152 }
153 G4VDataSetAlgorithm* algol = interp->Clone();
154
155 //G4VEMDataSet* elSet = new G4EMDataSet(i,energies,cs,algol,1.,1.);
156
157 G4VEMDataSet* elSet = new G4EMDataSet(i,energies,cs,log_energies,log_cs,algol,1.,1.);
158
159 setForMat->AddComponent(elSet);
160 }
161 set->push_back(setForMat);
162 }
163
164 return set;
165}
166
167//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
168
170 G4double cutEnergy,
171 G4int Z)
172{
173 G4double value = 0.;
174 if(energy > cutEnergy)
175 {
176 G4double elemCs = FindValue(Z, energy);
177 value = theBR->Probability(Z,cutEnergy, energy, energy);
178 value *= elemCs;
179 }
180 return value;
181}
std::vector< G4Element * > G4ElementVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double GetCrossSectionAboveThresholdForElement(G4double energy, G4double cutEnergy, G4int Z)
G4BremsstrahlungCrossSectionHandler(const G4VEnergySpectrum *spectrum, G4VDataSetAlgorithm *interpolation)
std::vector< G4VEMDataSet * > * BuildCrossSectionsForMaterials(const G4DataVector &energyVector, const G4DataVector *energyCuts)
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:205
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) 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