Geant4  10.01.p01
G4EmElementSelector.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: G4EmElementSelector.cc 83007 2014-07-24 14:46:57Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4EmElementSelector
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 29.05.2008
38 //
39 // Modifications:
40 //
41 // Class Description:
42 //
43 // Generic helper class for the random selection of an element
44 
45 // -------------------------------------------------------------------
46 //
47 
48 #include "G4EmElementSelector.hh"
49 #include "G4VEmModel.hh"
50 #include "G4SystemOfUnits.hh"
51 
52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54 
56  const G4Material* mat,
57  G4int bins,
58  G4double emin,
59  G4double emax,
60  G4bool):
61  model(mod), material(mat), nbins(bins), cutEnergy(-1.0),
62  lowEnergy(emin), highEnergy(emax)
63 {
65  nElmMinusOne = n - 1;
67  // element = (*theElementVector)[0];
68  if(nElmMinusOne > 0) {
69  xSections.reserve(n);
71  xSections.push_back(v0);
72  v0->SetSpline(false);
73  for(G4int i=1; i<n; ++i) {
75  xSections.push_back(v);
76  }
77  }
78  /*
79  G4cout << "G4EmElementSelector for " << mat->GetName() << " n= " << n
80  << " nbins= " << nbins << " Emin= " << lowEnergy
81  << " Emax= " << highEnergy << G4endl;
82  */
83 }
84 
85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
86 
88 {
89  if(nElmMinusOne > 0) {
90  for(G4int i=0; i<=nElmMinusOne; ++i) { delete xSections[i]; }
91  }
92 }
93 
94 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
95 
97  G4double cut)
98 {
99  //G4cout << "G4EmElementSelector initialise for " << material->GetName()
100  // << G4endl;
101  if(0 == nElmMinusOne || cut == cutEnergy) { return; }
102 
103  cutEnergy = cut;
104  //G4cout << "cut(keV)= " << cut/keV << G4endl;
105  G4double cross;
106 
107  const G4double* theAtomNumDensityVector =
109 
110  // loop over bins
111  for(G4int j=0; j<=nbins; ++j) {
112  G4double e = (xSections[0])->Energy(j);
113  model->SetupForMaterial(part, material, e);
114  cross = 0.0;
115  //G4cout << "j= " << j << " e(MeV)= " << e/MeV << G4endl;
116  for (G4int i=0; i<=nElmMinusOne; ++i) {
117  cross += theAtomNumDensityVector[i]*
119  cutEnergy, e);
120  xSections[i]->PutValue(j, cross);
121  }
122  }
123 
124  // xSections start from null, so use probabilities from the next bin
125  if(0.0 == (*xSections[nElmMinusOne])[0]) {
126  for (G4int i=0; i<=nElmMinusOne; ++i) {
127  xSections[i]->PutValue(0, (*xSections[i])[1]);
128  }
129  }
130  // xSections ends with null, so use probabilities from the previous bin
131  if(0.0 == (*xSections[nElmMinusOne])[nbins]) {
132  for (G4int i=0; i<=nElmMinusOne; ++i) {
133  xSections[i]->PutValue(nbins, (*xSections[i])[nbins-1]);
134  }
135  }
136  // perform normalization
137  for(G4int j=0; j<=nbins; ++j) {
138  cross = (*xSections[nElmMinusOne])[j];
139  // only for positive X-section
140  if(cross > 0.0) {
141  for (G4int i=0; i<nElmMinusOne; ++i) {
142  G4double x = (*xSections[i])[j]/cross;
143  xSections[i]->PutValue(j, x);
144  }
145  }
146  }
147  //G4cout << "======== G4EmElementSelector for the " << model->GetName()
148  // << G4endl;
149 }
150 
151 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
152 
154 {
155  G4cout << "======== G4EmElementSelector for the " << model->GetName();
156  if(part) G4cout << " and " << part->GetParticleName();
157  G4cout << " for " << material->GetName() << " ========" << G4endl;
158  if(0 < nElmMinusOne) {
159  for(G4int i=0; i<nElmMinusOne; i++) {
160  G4cout << " " << (*theElementVector)[i]->GetName() << " : " << G4endl;
161  G4cout << *(xSections[i]) << G4endl;
162  }
163  }
164  G4cout << "Last Element in element vector "
165  << (*theElementVector)[nElmMinusOne]->GetName()
166  << G4endl;
167  G4cout << G4endl;
168 }
169 
170 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
171 
172 
const G4Material * material
const G4String & GetName() const
Definition: G4Material.hh:178
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:395
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
void SetSpline(G4bool)
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
const G4ElementVector * theElementVector
const G4int n
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:315
void Dump(const G4ParticleDefinition *p=0)
#define G4endl
Definition: G4ios.hh:61
const G4String & GetName() const
Definition: G4VEmModel.hh:776
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4EmElementSelector(G4VEmModel *, const G4Material *, G4int bins, G4double emin, G4double emax, G4bool spline=true)
double G4double
Definition: G4Types.hh:76
std::vector< G4PhysicsLogVector * > xSections
void Initialise(const G4ParticleDefinition *, G4double cut=0.0)