Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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//
27// $Id$
28
29// class description
30//
31// This class is an interface to G4StaticSandiaData.
32// it provides - Sandia coeff for an element, given its Z
33// - sandia coeff for a material, given a pointer to it
34//
35
36//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
37//
38// History:
39//
40// 10.06.97 created. V. Grichine
41// 18.11.98 simplified public interface; new methods for materials. mma
42// 30.01.01 major bug in the computation of AoverAvo and in the units (/g!)
43// in GetSandiaCofPerAtom(). mma
44// 03.04.01 fnulcof[4] added; returned if energy < emin
45// 05.03.04 V.Grichine, new methods for old sorting algorithm for PAI model
46//
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
48
49#ifndef G4SANDIATABLE_HH
50#define G4SANDIATABLE_HH
51
52#include <assert.h>
53
54#include "G4OrderedTable.hh"
55#include "G4ios.hh"
56#include "globals.hh"
57
59
60class G4Material;
61
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
63
65{
66public: // with description
67
69
71
72 //main computation per atom:
73 static G4double* GetSandiaCofPerAtom(G4int Z, G4double energy);
74 static G4double GetZtoA(G4int Z);
75
76 //per volume of a material:
81
86
87 void SetVerbose(G4int ver){fVerbose = ver;};
88
89public: // without description
90
91 G4SandiaTable(__void__&);
92 // Fake default constructor for usage restricted to direct object
93 // persistency for clients requiring preallocation of memory for
94 // persistifiable objects.
95
96private:
97
98 void ComputeMatSandiaMatrix();
99 void ComputeMatSandiaMatrixPAI();
100
101 // methods per atom
102 G4int GetNbOfIntervals (G4int Z);
104 G4double GetIonizationPot (G4int Z);
105
106 // static members of the class
107 static const G4int fNumberOfElements;
108 static const G4int fIntervalLimit;
109 static const G4int fNumberOfIntervals;
110
111 static const G4double fSandiaTable[981][5];
112 static const G4int fNbOfIntervals[101];
113 static const G4double fZtoAratio[101];
114 static const G4double fIonizationPotentials[101];
115 static const G4double funitc[5];
116
117 static G4int fCumulInterval[101];
118 static G4double fSandiaCofPerAtom[4];
119
120 // members of the class
121 G4Material* fMaterial;
122 G4int fMatNbOfIntervals;
123 G4OrderedTable* fMatSandiaMatrix;
124 G4OrderedTable* fMatSandiaMatrixPAI;
125
126 G4double fnulcof[4];
127
128/////////////////////////////////////////////////////////////////////
129//
130// Methods for old implementation of PAI model
131// Will be removed for the next major release
132
133public: // without description
134
136
137 inline void SandiaSwap( G4double** da,
138 G4int i,
139 G4int j );
140
141 void SandiaSort( G4double** da,
142 G4int sz );
143
145 G4int el );
146
148 const G4double* fractionW,
149 G4int el,
150 G4int mi );
151
152 inline G4double GetPhotoAbsorpCof(G4int i , G4int j) const;
153
154 inline G4int GetMaxInterval() const;
155
156 inline G4double** GetPointerToCof();
157
158private:
159
160 void ComputeMatTable();
161
162 // copy constructor and hide assignment operator
164 G4SandiaTable & operator=(const G4SandiaTable &right);
165
166//////////////////////////////////////////////////////////////////////////
167//
168// data members for PAI model
169
170private:
171
172 G4double** fPhotoAbsorptionCof ; // SandiaTable for mixture
173
174 G4int fMaxInterval ;
175 G4int fVerbose;
176
177
178//
179//
180//////////////////////////////////////////////////////////////////////////
181
182};
183
184///////////////////////////////////////////////////////////////////////
185//
186// Inline methods for PAI model, will be removed in next major release
187
188inline G4int
190 return fMaxInterval;
191}
192
193inline G4double**
195{
196 if(!fPhotoAbsorptionCof) { ComputeMatTable(); }
197 return fPhotoAbsorptionCof;
198}
199
200inline void
202 G4int i,
203 G4int j )
204{
205 G4double tmp = da[i][0] ;
206 da[i][0] = da[j][0] ;
207 da[j][0] = tmp ;
208}
209
210inline
212{
213 return fPhotoAbsorptionCof[i][j]*funitc[j];
214}
215
216//
217//
218////////////////////////////////////////////////////////////////////////////
219
220
221#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
void SetVerbose(G4int ver)
G4double GetSandiaMatTablePAI(G4int, G4int)
void SandiaSwap(G4double **da, G4int i, G4int j)
G4int GetMaxInterval() const
static G4double GetZtoA(G4int Z)
G4double GetSandiaCofForMaterialPAI(G4int, G4int)
static G4double * GetSandiaCofPerAtom(G4int Z, G4double energy)
G4int GetMatNbOfIntervals()
G4int SandiaMixing(G4int Z[], const G4double *fractionW, G4int el, G4int mi)
G4double GetSandiaMatTable(G4int, G4int)
G4double GetSandiaCofForMaterial(G4int, G4int)
void SandiaSort(G4double **da, G4int sz)
G4OrderedTable * GetSandiaMatrixPAI()
G4int SandiaIntervals(G4int Z[], G4int el)
G4double GetPhotoAbsorpCof(G4int i, G4int j) const
G4double ** GetPointerToCof()