Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4eIonisationCrossSectionHandler Class Reference

#include <G4eIonisationCrossSectionHandler.hh>

Inheritance diagram for G4eIonisationCrossSectionHandler:
Collaboration diagram for G4eIonisationCrossSectionHandler:

Public Member Functions

 G4eIonisationCrossSectionHandler (const G4VEnergySpectrum *spec, G4VDataSetAlgorithm *alg, G4double emin, G4double emax, G4int nbin)
 
 ~G4eIonisationCrossSectionHandler ()
 
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 G4VDataSetAlgorithmCreateInterpolation ()
 
const G4VDataSetAlgorithmGetInterpolation () const
 

Detailed Description

Definition at line 63 of file G4eIonisationCrossSectionHandler.hh.

Constructor & Destructor Documentation

G4eIonisationCrossSectionHandler::G4eIonisationCrossSectionHandler ( const G4VEnergySpectrum spec,
G4VDataSetAlgorithm alg,
G4double  emin,
G4double  emax,
G4int  nbin 
)

Definition at line 73 of file G4eIonisationCrossSectionHandler.cc.

77  theParam(spec),verbose(0)
78 {
79  G4VCrossSectionHandler::Initialise(alg, emin, emax, nbin);
80  interp = new G4LinLogLogInterpolation();
81 }
static const G4double emax
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)

Here is the call graph for this function:

G4eIonisationCrossSectionHandler::~G4eIonisationCrossSectionHandler ( )

Definition at line 84 of file G4eIonisationCrossSectionHandler.cc.

85 {
86  delete interp;
87 }

Member Function Documentation

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

Implements G4VCrossSectionHandler.

Definition at line 90 of file G4eIonisationCrossSectionHandler.cc.

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->GetAtomicNumDensityVector();
114  G4int nElements = material->GetNumberOfElements();
115 
116  if(verbose > 0)
117  {
118  G4cout << "eIonisation CS for " << mLocal << "th material "
119  << material->GetName()
120  << " eEl= " << nElements << G4endl;
121  }
122 
123  G4double tcut = (*energyCuts)[mLocal];
124 
125  G4VDataSetAlgorithm* algo = interp->Clone();
126  G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo,1.,1.);
127 
128  for (G4int i=0; i<nElements; i++) {
129 
130  G4int Z = (G4int) (*elementVector)[i]->GetZ();
131  G4int nShells = NumberOfComponents(Z);
132 
133  energies = new G4DataVector;
134  cs = new G4DataVector;
135 
136  log_energies = new G4DataVector;
137  log_cs = new G4DataVector;
138 
139  G4double density = nAtomsPerVolume[i];
140 
141  for (G4int bin=0; bin<nOfBins; bin++) {
142 
143  G4double e = energyVector[bin];
144  energies->push_back(e);
145  log_energies->push_back(std::log10(e));
146  G4double value = 0.0;
147  G4double log_value = -300;
148 
149  if(e > tcut) {
150  for (G4int n=0; n<nShells; n++) {
151  G4double cross = FindValue(Z, e, n);
152  G4double p = theParam->Probability(Z, tcut, e, e, n);
153  value += cross * p * density;
154 
155  if(verbose>0 && mLocal == 0 && e>=1. && e<=0.)
156  {
157  G4cout << "G4eIonCrossSH: e(MeV)= " << e/MeV
158  << " n= " << n
159  << " cross= " << cross
160  << " p= " << p
161  << " value= " << value
162  << " tcut(MeV)= " << tcut/MeV
163  << " rho= " << density
164  << " Z= " << Z
165  << G4endl;
166  }
167 
168  }
169  if (value == 0.) value = 1e-300;
170  log_value = std::log10(value);
171  }
172  cs->push_back(value);
173  log_cs->push_back(log_value);
174  }
175  G4VDataSetAlgorithm* algoLocal = interp->Clone();
176 
177  //G4VEMDataSet* elSet = new G4EMDataSet(i,energies,cs,algoLocal,1.,1.);
178 
179  G4VEMDataSet* elSet = new G4EMDataSet(i,energies,cs,log_energies,log_cs,algoLocal,1.,1.);
180 
181  setForMat->AddComponent(elSet);
182  }
183  set->push_back(setForMat);
184  }
185 
186  return set;
187 }
std::vector< G4Element * > G4ElementVector
tuple bin
Definition: plottest35.py:22
const char * p
Definition: xmltok.h:285
const G4String & GetName() const
Definition: G4Material.hh:178
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
virtual G4VDataSetAlgorithm * Clone() const =0
string material
Definition: eplot.py:19
G4int NumberOfComponents(G4int Z) const
G4double FindValue(G4int Z, G4double e) const
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
virtual void AddComponent(G4VEMDataSet *dataSet)=0
const G4int n
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:216
virtual G4double Probability(G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const =0
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76
const G4Material * GetMaterial() const

Here is the call graph for this function:

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

Definition at line 189 of file G4eIonisationCrossSectionHandler.cc.

192 {
193  G4int nShells = NumberOfComponents(Z);
194  G4double value = 0.;
195  if(energy > cutEnergy)
196  {
197  for (G4int n=0; n<nShells; n++) {
198  G4double cross = FindValue(Z, energy, n);
199  G4double p = theParam->Probability(Z, cutEnergy, energy, energy, n);
200  value += cross * p;
201  }
202  }
203  return value;
204 }
const char * p
Definition: xmltok.h:285
int G4int
Definition: G4Types.hh:78
G4int NumberOfComponents(G4int Z) const
G4double FindValue(G4int Z, G4double e) const
const XML_Char int const XML_Char * value
Definition: expat.h:331
const G4int n
virtual G4double Probability(G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const =0
G4double energy(const ThreeVector &p, const G4double m)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:


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