Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UCNMaterialPropertiesTable.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//---------------------------------------------------------------------------
28//
29// ClassName: G4UCNMaterialPropertiesTable
30//
31// Class description:
32//
33// A derived class of G4MaterialPropertiesTable in order to save the look-up
34// table for the microroughness probability. The derived class has four
35// pointers to G4double-arrays:
36// (1) integral prob. for reflection
37// (2) maximum probability for reflection (needed for accept-reject method)
38// (3) integral prob. for transmission
39// (4) maximum probability for transmission
40//
41// 12-05-14, adopted from Stefan Heule (PSI) Thesis by P.Gumplinger
42// http://ucn.web.psi.ch/papers/stefanheule_thesis2008.pdf
43// reported in F. Atchison et al., Eur. Phys. J. A 44, 23–29 (2010)
44// Thanks to Geza Zsigmond
45
46#ifndef G4UCNMATERIALPROPERTIESTABLE_HH
47#define G4UCNMATERIALPROPERTIESTABLE_HH 1
48
50
52{
53public:
54
57
58 // returns the pointer to the mr-reflection table
60
61 // returns the pointer to the mr-transmission table
63
64 // Assigns double-array to the table-pointers, currently not used
66
67 // Creates new double arrays and assigns them to the table pointers
69
70 // Reads the MR-parameters from the corresponding fields and starts
71 // the computation of the mr-tables
73
74 // returns the integral prob. value for a theta_i - E pair
76
77 // returns the maximum prob. value for a theta_i - E pair
79
80 // sets the maximum prob. value for a theta_i - E pair
82
83 // returns the mr-prob.
84
85 // arguments:
86 // 1) theta_i
87 // 2) Energy
88 // 3) V_F
89 // 4) theta_o
90 // 5) phi_o
91
93
94 // returns the integral transmission prob. value for a theta_i - E pair
96
97 // returns the maximum transmission prob. for a theta_i - E pair
99
100 // sets the maximum prob. value for a theta_i - E pair
102
103 // returns the mr-transmission-prob.
104
105 // arguments:
106 // 1) theta_i
107 // 2) E
108 // 3) V_F
109 // 4) theta_o
110 // 5) phi_o
111
114
115 // Checks if the validity condition for the microroughness model are
116 // satisfied, cf. Steyerl-paper p. 175
117 G4bool ConditionsValid (G4double E, G4double VFermi, G4double theta_i);
118
119 // Checks if the validity conditions for the transmission of the
120 // microroughness model are satisfied
122
123 // Adds the values for mr-related units to the MaterialPropertiesTable
124
125 // arguments:
126 // 1) w
127 // 2) b
128 // 3) number of angles theta_i in the look-up tables
129 // 4) number of energies in the look-up tables
130 // 5) minimum value of theta_i
131 // 6) maximum value of theta_i
132 // 7) minimum value of E
133 // 8) maximum value of E
134 // 9) number of angles theta_o in the look-up table calculation
135 // 10) number of angles phi_o in the look-up table calculation
136 // 11) angular cut
137
140 G4double);
141
142 // returns b
143 G4double GetRMS() const;
144
145 // returns w
146 G4double GetCorrLen() const;
147
148private:
149
150 // Pointer to the integral reflection probability table
151 G4double* theMicroRoughnessTable;
152
153 // Pointer to the maximum reflection probability table
154 G4double* maxMicroRoughnessTable;
155
156 // Pointer to the integral transmission probability table
157 G4double* theMicroRoughnessTransTable;
158
159 // Pointer to the maximum transmission probability table
160 G4double* maxMicroRoughnessTransTable;
161
162 G4double theta_i_min;
163 G4double theta_i_max;
164 G4double Emin;
165 G4double Emax;
166 G4int no_theta_i;
167 G4int noE;
168 G4double theta_i_step;
169 G4double E_step;
170
171 // RMS roughness and correlation length
172 G4double b, w;
173 G4double AngCut;
174};
175
176// ==========================================================================
177// inline functions
178// ==========================================================================
179
182
183#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4double GetMRMaxTransProbability(G4double, G4double)
void SetMicroRoughnessParameters(G4double, G4double, G4int, G4int, G4double, G4double, G4double, G4double, G4int, G4int, G4double)
G4double GetMRProbability(G4double, G4double, G4double, G4double, G4double)
G4bool ConditionsValid(G4double E, G4double VFermi, G4double theta_i)
void LoadMicroRoughnessTables(G4double *, G4double *, G4double *, G4double *)
G4bool TransConditionsValid(G4double E, G4double VFermi, G4double theta_i)
void SetMRMaxProbability(G4double, G4double, G4double)
void SetMRMaxTransProbability(G4double, G4double, G4double)
G4double GetMRIntTransProbability(G4double, G4double)
G4double GetMRMaxProbability(G4double, G4double)
G4double GetMRTransProbability(G4double, G4double, G4double, G4double, G4double)
G4double GetMRIntProbability(G4double, G4double)