Geant4-11
G4SandiaTable.hh
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// class description
27//
28// This class is an interface to G4StaticSandiaData.
29// it provides - Sandia coeff for an element, given its Z
30// - sandia coeff for a material, given a pointer to it
31//
32
33//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
34//
35// History:
36//
37// 10.06.97 created. V. Grichine
38// 18.11.98 simplified public interface; new methods for materials. mma
39// 30.01.01 major bug in the computation of AoverAvo and in the units (/g!)
40// in GetSandiaCofPerAtom(). mma
41// 03.04.01 fnulcof[4] added; returned if energy < emin
42// 05.03.04 V.Grichine, new methods for old sorting algorithm for PAI model
43// 21.21.13 V.Ivanchenko, changed signature of methods, reduced number of
44// static variables, methods
45//
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
47
48#ifndef G4SANDIATABLE_HH
49#define G4SANDIATABLE_HH
50
51#include "G4OrderedTable.hh"
52#include "G4ios.hh"
53#include "globals.hh"
54#include <assert.h>
55#include <vector>
56
58
59class G4Material;
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
62
64{
65public: // with description
66
68
70
71 //main computation per atom:
73 std::vector<G4double>& coeff) const;
74
76 std::vector<G4double>& coeff) const;
77
80
81 static G4double GetZtoA(G4int Z);
82
83 //per volume of a material:
88
91
92 inline void SetVerbose(G4int ver) { fVerbose = ver; };
93
94public: // without description
95
96 G4SandiaTable(__void__&);
97 // Fake default constructor for usage restricted to direct object
98 // persistency for clients requiring preallocation of memory for
99 // persistifiable objects.
100
101private:
102
105
106 // methods per atom
108
109#ifdef G4VERBOSE
110 static G4int PrintErrorZ(G4int Z, const G4String&);
111 static void PrintErrorV(const G4String&);
112#endif
113
114 // computed once
116 static const G4double funitc[5];
117
118 // used at initialisation
119 std::vector<G4double> fSandiaCofPerAtom;
120
121 // members of the class
126
128//
129// Methods for implementation of PAI model
130//
132
133public: // without description
134
135 G4SandiaTable(G4int matIndex);
136
138
139 void Initialize(const G4Material*);
140
142
143 G4int SandiaMixing(G4int Z[], const G4double* fractionW,
144 G4int el, G4int mi);
145
147
148 G4int GetMaxInterval() const;
149
150 inline G4bool GetLowerI1() {return fLowerI1;};
151 inline void SetLowerI1(G4bool flag) {fLowerI1=flag;};
152
153 // operators
154 G4bool operator==(const G4SandiaTable&) const = delete;
155 G4bool operator!=(const G4SandiaTable&) const = delete;
157 G4SandiaTable & operator=(const G4SandiaTable &right) = delete;
158
159private:
160
161 void ComputeMatTable();
162
163 void SandiaSwap(G4double** da, G4int i, G4int j);
164
165 void SandiaSort(G4double** da, G4int sz);
166
168
169 static const G4double fSandiaTable[981][5];
171 static const G4int fIntervalLimit;
173 static const G4int fH2OlowerInt;
174
175 // data members for PAI model
176 G4double** fPhotoAbsorptionCof; // SandiaTable for mixture
177
181};
182
183#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
static const G4int fIntervalLimit
G4int fMatNbOfIntervals
void SetVerbose(G4int ver)
void Initialize(const G4Material *)
void ComputeMatTable()
void SandiaSwap(G4double **da, G4int i, G4int j)
G4int GetMaxInterval() const
G4SandiaTable & operator=(const G4SandiaTable &right)=delete
void ComputeMatSandiaMatrix()
G4bool operator==(const G4SandiaTable &) const =delete
G4double GetWaterEnergyLimit() const
G4int GetMatNbOfIntervals() const
static G4double GetZtoA(G4int Z)
G4double GetWaterCofForMaterial(G4int, G4int) const
G4SandiaTable(G4SandiaTable &)=delete
G4double GetSandiaCofForMaterial(G4int, G4int) const
G4OrderedTable * fMatSandiaMatrix
static G4int fCumulInterval[101]
const G4Material * fMaterial
const G4double * GetSandiaCofForMaterialPAI(G4double energy) const
G4double GetSandiaPerAtom(G4int Z, G4int, G4int) const
void GetSandiaCofWater(G4double energy, std::vector< G4double > &coeff) const
void SetLowerI1(G4bool flag)
static const G4int fNumberOfElements
static const G4double fSandiaTable[981][5]
static const G4int fNumberOfIntervals
static const G4double funitc[5]
G4bool operator!=(const G4SandiaTable &) const =delete
G4double GetSandiaMatTablePAI(G4int, G4int) const
G4int SandiaMixing(G4int Z[], const G4double *fractionW, G4int el, G4int mi)
void ComputeMatSandiaMatrixPAI()
static const G4int fH2OlowerInt
G4double ** fPhotoAbsorptionCof
G4OrderedTable * fMatSandiaMatrixPAI
std::vector< G4double > fSandiaCofPerAtom
void GetSandiaCofPerAtom(G4int Z, G4double energy, std::vector< G4double > &coeff) const
void SandiaSort(G4double **da, G4int sz)
G4int SandiaIntervals(G4int Z[], G4int el)
G4double GetSandiaMatTable(G4int, G4int) const
G4bool GetLowerI1()
G4double GetPhotoAbsorpCof(G4int i, G4int j) const
G4double ** GetPointerToCof()
G4double energy(const ThreeVector &p, const G4double m)