Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4RDBremsstrahlungCrossSectionHandler.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 // GEANT4 tag $Name: geant4-09-01-ref-00 $
28 //
29 // -------------------------------------------------------------------
30 //
31 // GEANT4 Class file
32 //
33 //
34 // File name: G4RDBremsstrahlungCrossSectionHandler
35 //
36 // Author: V.Ivanchenko (Vladimir.Ivanchenko@cern.ch)
37 //
38 // Creation date: 25 September 2001
39 //
40 // Modifications:
41 // 10.10.2001 MGP Revision to improve code quality and consistency with design
42 // 21.01.2003 VI cut per region
43 //
44 // -------------------------------------------------------------------
45 
48 #include "G4DataVector.hh"
50 #include "G4RDVDataSetAlgorithm.hh"
52 #include "G4RDVEMDataSet.hh"
53 #include "G4RDEMDataSet.hh"
54 #include "G4Material.hh"
55 #include "G4ProductionCutsTable.hh"
56 
59  : theBR(spec)
60 {
61  interp = new G4RDSemiLogInterpolation();
62 }
63 
64 
66 {
67  delete interp;
68 }
69 
70 
71 std::vector<G4RDVEMDataSet*>*
73  const G4DataVector* energyCuts)
74 {
75  std::vector<G4RDVEMDataSet*>* set = new std::vector<G4RDVEMDataSet*>;
76 
77  G4DataVector* energies;
78  G4DataVector* cs;
79  G4int nOfBins = energyVector.size();
80 
81  const G4ProductionCutsTable* theCoupleTable=
83  size_t numOfCouples = theCoupleTable->GetTableSize();
84 
85  for (size_t m=0; m<numOfCouples; m++) {
86 
87  const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(m);
88  const G4Material* material= couple->GetMaterial();
89  const G4ElementVector* elementVector = material->GetElementVector();
90  const G4double* nAtomsPerVolume = material->GetVecNbOfAtomsPerVolume();
91  G4int nElements = material->GetNumberOfElements();
92 
93  G4double tcut = (*energyCuts)[m];
94 
95  G4RDVDataSetAlgorithm* algo = interp->Clone();
96  G4RDVEMDataSet* setForMat = new G4RDCompositeEMDataSet(algo,1.,1.);
97 
98  for (G4int i=0; i<nElements; i++) {
99 
100  G4int Z = (G4int) ((*elementVector)[i]->GetZ());
101  energies = new G4DataVector;
102  cs = new G4DataVector;
103  G4double density = nAtomsPerVolume[i];
104 
105  for (G4int bin=0; bin<nOfBins; bin++) {
106 
107  G4double e = energyVector[bin];
108  energies->push_back(e);
109  G4double value = 0.0;
110 
111  if(e > tcut) {
112  G4double elemCs = FindValue(Z, e);
113 
114  value = theBR->Probability(Z, tcut, e, e);
115 
116  value *= elemCs*density;
117  }
118  cs->push_back(value);
119  }
120  G4RDVDataSetAlgorithm* algol = interp->Clone();
121  G4RDVEMDataSet* elSet = new G4RDEMDataSet(i,energies,cs,algol,1.,1.);
122  setForMat->AddComponent(elSet);
123  }
124  set->push_back(setForMat);
125  }
126 
127  return set;
128 }
std::vector< G4RDVEMDataSet * > * BuildCrossSectionsForMaterials(const G4DataVector &energyVector, const G4DataVector *energyCuts)
virtual G4RDVDataSetAlgorithm * Clone() const =0
std::vector< G4Element * > G4ElementVector
tuple bin
Definition: plottest35.py:22
G4double FindValue(G4int Z, G4double e) const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
string material
Definition: eplot.py:19
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
static constexpr double m
Definition: G4SIunits.hh:129
const XML_Char int const XML_Char * value
Definition: expat.h:331
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4RDBremsstrahlungCrossSectionHandler(const G4RDVEnergySpectrum *spectrum, G4RDVDataSetAlgorithm *interpolation)
double G4double
Definition: G4Types.hh:76
virtual void AddComponent(G4RDVEMDataSet *dataSet)=0
const G4Material * GetMaterial() const
virtual G4double Probability(G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const =0