Geant4 11.2.2
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// ClassName: G4UCNMaterialPropertiesTable
29//
30// Class description:
31//
32// A derived class of G4MaterialPropertiesTable in order to save the look-up
33// table for the microroughness probability. The derived class has four
34// pointers to G4double-arrays:
35// (1) integral prob. for reflection
36// (2) maximum probability for reflection (needed for accept-reject method)
37// (3) integral prob. for transmission
38// (4) maximum probability for transmission
39//
40// 12-05-14, adopted from Stefan Heule (PSI) Thesis by P.Gumplinger
41// http://ucn.web.psi.ch/papers/stefanheule_thesis2008.pdf
42// reported in F. Atchison et al., Eur. Phys. J. A 44, 23–29 (2010)
43// Thanks to Geza Zsigmond
44
45#ifndef G4UCNMATERIALPROPERTIESTABLE_HH
46#define G4UCNMATERIALPROPERTIESTABLE_HH 1
47
49
51{
52 public:
55
56 // returns the pointer to the mr-reflection table
58
59 // returns the pointer to the mr-transmission table
61
62 // Assigns double-array to the table-pointers, currently not used
64
65 // Creates new double arrays and assigns them to the table pointers
67
68 // Reads the MR-parameters from the corresponding fields and starts
69 // the computation of the mr-tables
71
72 // returns the integral prob. value for a theta_i - E pair
74
75 // returns the maximum prob. value for a theta_i - E pair
77
78 // sets the maximum prob. value for a theta_i - E pair
80
81 // returns the mr-prob.
82
83 // arguments:
84 // 1) theta_i
85 // 2) Energy
86 // 3) V_F
87 // 4) theta_o
88 // 5) phi_o
90
91 // returns the integral transmission prob. value for a theta_i - E pair
93
94 // returns the maximum transmission prob. for a theta_i - E pair
96
97 // sets the maximum prob. value for a theta_i - E pair
99
100 // returns the mr-transmission-prob.
101 // arguments:
102 // 1) theta_i
103 // 2) E
104 // 3) V_F
105 // 4) theta_o
106 // 5) phi_o
108
109 // Checks if the validity condition for the microroughness model are
110 // satisfied, cf. Steyerl-paper p. 175
111 G4bool ConditionsValid(G4double E, G4double VFermi, G4double theta_i);
112
113 // Checks if the validity conditions for the transmission of the
114 // microroughness model are satisfied
116
117 // Adds the values for mr-related units to the MaterialPropertiesTable
118 // arguments:
119 // 1) w
120 // 2) b
121 // 3) number of angles theta_i in the look-up tables
122 // 4) number of energies in the look-up tables
123 // 5) minimum value of theta_i
124 // 6) maximum value of theta_i
125 // 7) minimum value of E
126 // 8) maximum value of E
127 // 9) number of angles theta_o in the look-up table calculation
128 // 10) number of angles phi_o in the look-up table calculation
129 // 11) angular cut
132
133 // returns b
134 G4double GetRMS() const;
135
136 // returns w
137 G4double GetCorrLen() const;
138
139 private:
140 // Pointer to the integral reflection probability table
141 G4double* theMicroRoughnessTable;
142
143 // Pointer to the maximum reflection probability table
144 G4double* maxMicroRoughnessTable;
145
146 // Pointer to the integral transmission probability table
147 G4double* theMicroRoughnessTransTable;
148
149 // Pointer to the maximum transmission probability table
150 G4double* maxMicroRoughnessTransTable;
151
152 G4double theta_i_min;
153 G4double theta_i_max;
154 G4double Emin;
155 G4double Emax;
156 G4int no_theta_i;
157 G4int noE;
158 G4double theta_i_step;
159 G4double E_step;
160
161 // RMS roughness and correlation length
162 G4double b, w;
163 G4double AngCut;
164};
165
166// ==========================================================================
167// inline functions
168// ==========================================================================
169
172
173#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)