Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions
G4BremsstrahlungCrossSectionHandler Class Reference

#include <G4BremsstrahlungCrossSectionHandler.hh>

Inheritance diagram for G4BremsstrahlungCrossSectionHandler:
G4VCrossSectionHandler

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 G4VDataSetAlgorithmCreateInterpolation ()
 
const G4VDataSetAlgorithmGetInterpolation () const
 

Detailed Description

Definition at line 62 of file G4BremsstrahlungCrossSectionHandler.hh.

Constructor & Destructor Documentation

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

Definition at line 74 of file G4BremsstrahlungCrossSectionHandler.cc.

76  : theBR(spec)
77 {
78  interp = new G4SemiLogInterpolation();
79 }
G4BremsstrahlungCrossSectionHandler::~G4BremsstrahlungCrossSectionHandler ( )

Definition at line 83 of file G4BremsstrahlungCrossSectionHandler.cc.

84 {
85  delete interp;
86 }

Member Function Documentation

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

Implements G4VCrossSectionHandler.

Definition at line 91 of file G4BremsstrahlungCrossSectionHandler.cc.

References G4VEMDataSet::AddComponent(), plottest35::bin, G4VDataSetAlgorithm::Clone(), density, G4VCrossSectionHandler::FindValue(), G4Material::GetElementVector(), G4MaterialCutsCouple::GetMaterial(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4Material::GetNumberOfElements(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetTableSize(), G4Material::GetVecNbOfAtomsPerVolume(), eplot::material, and G4VEnergySpectrum::Probability().

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 }
std::vector< G4Element * > G4ElementVector
tuple bin
Definition: plottest35.py:22
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
virtual G4VDataSetAlgorithm * Clone() const =0
string material
Definition: eplot.py:19
G4double density
Definition: TRTMaterials.hh:39
G4double FindValue(G4int Z, G4double e) const
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
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
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
const XML_Char int const XML_Char * value
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
double G4double
Definition: G4Types.hh:76
const G4Material * GetMaterial() const
G4double G4BremsstrahlungCrossSectionHandler::GetCrossSectionAboveThresholdForElement ( G4double  energy,
G4double  cutEnergy,
G4int  Z 
)

Definition at line 169 of file G4BremsstrahlungCrossSectionHandler.cc.

References G4VCrossSectionHandler::FindValue(), and G4VEnergySpectrum::Probability().

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 }
G4double FindValue(G4int Z, G4double e) const
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
virtual G4double Probability(G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const =0
const XML_Char int const XML_Char * value
double G4double
Definition: G4Types.hh:76

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