source: trunk/source/processes/electromagnetic/lowenergy/src/G4PenelopeCrossSectionHandler.cc

Last change on this file was 1347, checked in by garnier, 15 years ago

geant4 tag 9.4

File size: 5.3 KB
Line 
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// GEANT4 Class file
27//
28//
29// File name: G4PenelopeCrossSectionHandler
30//
31// Author: Luciano Pandola
32//
33// Creation date: 04 Jul 2003
34//
35// -------------------------------------------------------------------
36#include "G4PenelopeCrossSectionHandler.hh"
37#include "G4PenelopeIonisation.hh"
38#include "G4DataVector.hh"
39#include "G4CompositeEMDataSet.hh"
40#include "G4VDataSetAlgorithm.hh"
41#include "G4LinLogLogInterpolation.hh"
42#include "G4SemiLogInterpolation.hh"
43#include "G4VEMDataSet.hh"
44#include "G4EMDataSet.hh"
45#include "G4Material.hh"
46#include "G4ProductionCutsTable.hh"
47#include "G4ParticleDefinition.hh"
48
49G4PenelopeCrossSectionHandler::G4PenelopeCrossSectionHandler(G4PenelopeIonisation* theProcess,
50 const G4ParticleDefinition& aParticleType,
51 G4VDataSetAlgorithm* alg,
52 G4double emin,
53 G4double emax,
54 G4int nbin)
55 : G4VCrossSectionHandler(),
56 process(theProcess),
57 particle(aParticleType)
58{
59 G4VCrossSectionHandler::Initialise(alg, emin, emax, nbin);
60 interp = new G4LinLogLogInterpolation();
61}
62
63
64G4PenelopeCrossSectionHandler::~G4PenelopeCrossSectionHandler()
65{
66 delete interp;
67}
68
69
70std::vector<G4VEMDataSet*>* G4PenelopeCrossSectionHandler::BuildCrossSectionsForMaterials(
71 const G4DataVector& energyVector,
72 const G4DataVector* energyCuts)
73{
74 //G4int verbose = 0;
75 std::vector<G4VEMDataSet*>* set = new std::vector<G4VEMDataSet*>;
76
77 G4DataVector* energies;
78 G4DataVector* cs;
79 G4int nOfBins = energyVector.size();
80
81
82 const G4ProductionCutsTable* theCoupleTable=
83 G4ProductionCutsTable::GetProductionCutsTable();
84 size_t numOfCouples = theCoupleTable->GetTableSize();
85
86 for (size_t m=0; m<numOfCouples; m++) {
87
88 const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(m);
89 const G4Material* material= couple->GetMaterial();
90 const G4ElementVector* elementVector = material->GetElementVector();
91 const G4double* nAtomsPerVolume = material->GetAtomicNumDensityVector();
92 G4double electronVolumeDensity =
93 material->GetTotNbOfElectPerVolume(); //electron density
94 G4int nElements = material->GetNumberOfElements();
95
96 /*
97 if(verbose > 0) {
98 G4cout << "Penelope CS for " << m << "th material "
99 << material->GetName()
100 << " eEl= " << nElements << G4endl;
101 }
102 */
103
104 G4double tcut = (*energyCuts)[m];
105
106 G4VDataSetAlgorithm* algo = interp->Clone();
107 G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo,1.,1.);
108
109 for (G4int i=0; i<nElements; i++) {
110
111 G4int Z = (G4int) (*elementVector)[i]->GetZ();
112 energies = new G4DataVector;
113 cs = new G4DataVector;
114 G4double density = nAtomsPerVolume[i];
115
116 for (G4int bin=0; bin<nOfBins; bin++) {
117
118 G4double e = energyVector[bin];
119 energies->push_back(e);
120 G4double value = 0.0;
121
122 if(e > tcut) {
123 G4double cross = FindValue(Z, e);
124 G4double p = process->CalculateCrossSectionsRatio(e,tcut,Z,electronVolumeDensity,
125 particle);
126 value += cross * p * density;
127
128 /*
129 if(verbose>0 && m == 0 && e>=1. && e<=0.) {
130 G4cout << "G4PenIonCrossSH: e(MeV)= " << e/MeV
131 << " cross= " << cross
132 << " p= " << p
133 << " value= " << value
134 << " tcut(MeV)= " << tcut/MeV
135 << " rho= " << density
136 << " Z= " << Z
137 << G4endl;
138 }
139 */
140
141 }
142 cs->push_back(value);
143 }
144
145 G4VDataSetAlgorithm* algo = interp->Clone();
146 G4VEMDataSet* elSet = new G4EMDataSet(i,energies,cs,algo,1.,1.);
147 setForMat->AddComponent(elSet);
148 }
149 set->push_back(setForMat);
150 }
151 return set;
152}
153
154
Note: See TracBrowser for help on using the repository browser.