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: G4BremsstrahlungCrossSectionHandler.cc,v 1.9 2006/06/29 19:38:42 gunter Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-02-ref-02 $ |
---|
28 | // |
---|
29 | // ------------------------------------------------------------------- |
---|
30 | // |
---|
31 | // GEANT4 Class file |
---|
32 | // |
---|
33 | // |
---|
34 | // File name: G4BremsstrahlungCrossSectionHandler |
---|
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 | |
---|
46 | #include "G4BremsstrahlungCrossSectionHandler.hh" |
---|
47 | #include "G4eBremsstrahlungSpectrum.hh" |
---|
48 | #include "G4DataVector.hh" |
---|
49 | #include "G4CompositeEMDataSet.hh" |
---|
50 | #include "G4VDataSetAlgorithm.hh" |
---|
51 | #include "G4SemiLogInterpolation.hh" |
---|
52 | #include "G4VEMDataSet.hh" |
---|
53 | #include "G4EMDataSet.hh" |
---|
54 | #include "G4Material.hh" |
---|
55 | #include "G4ProductionCutsTable.hh" |
---|
56 | |
---|
57 | G4BremsstrahlungCrossSectionHandler::G4BremsstrahlungCrossSectionHandler(const G4VEnergySpectrum* spec, |
---|
58 | G4VDataSetAlgorithm* ) |
---|
59 | : theBR(spec) |
---|
60 | { |
---|
61 | interp = new G4SemiLogInterpolation(); |
---|
62 | } |
---|
63 | |
---|
64 | |
---|
65 | G4BremsstrahlungCrossSectionHandler::~G4BremsstrahlungCrossSectionHandler() |
---|
66 | { |
---|
67 | delete interp; |
---|
68 | } |
---|
69 | |
---|
70 | |
---|
71 | std::vector<G4VEMDataSet*>* |
---|
72 | G4BremsstrahlungCrossSectionHandler::BuildCrossSectionsForMaterials(const G4DataVector& energyVector, |
---|
73 | const G4DataVector* energyCuts) |
---|
74 | { |
---|
75 | std::vector<G4VEMDataSet*>* set = new std::vector<G4VEMDataSet*>; |
---|
76 | |
---|
77 | G4DataVector* energies; |
---|
78 | G4DataVector* cs; |
---|
79 | G4int nOfBins = energyVector.size(); |
---|
80 | |
---|
81 | const G4ProductionCutsTable* theCoupleTable= |
---|
82 | G4ProductionCutsTable::GetProductionCutsTable(); |
---|
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 | G4VDataSetAlgorithm* algo = interp->Clone(); |
---|
96 | G4VEMDataSet* setForMat = new G4CompositeEMDataSet(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 | G4VDataSetAlgorithm* algol = interp->Clone(); |
---|
121 | G4VEMDataSet* elSet = new G4EMDataSet(i,energies,cs,algol,1.,1.); |
---|
122 | setForMat->AddComponent(elSet); |
---|
123 | } |
---|
124 | set->push_back(setForMat); |
---|
125 | } |
---|
126 | |
---|
127 | return set; |
---|
128 | } |
---|