Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4UCNBoundaryProcess.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: G4UCNBoundaryProcess.hh 71487 2013-06-17 08:19:40Z gcosmo $
28 //
29 //
31 // Ultra Cold Neutron (UCN) Boundary Process Class Definition
33 //
34 // File: G4UCNBoundaryProcess.hh
35 // Description: Discrete Process -- Boundary Process of UCN
36 // Version: 1.0
37 // Created: 2014-06-04
38 // Author: Peter Gumplinger
39 // Adopted from: UCNMaterialBoundary by Peter Fierlinger 4.9.2004
40 // Updated:
41 // mail: gum@triumf.ca
42 //
44 
45 #ifndef G4UCNBOUNDARYPROCESS_HH
46 #define G4UCNBOUNDARYPROCESS_HH 1
47 
49 // Includes
51 
52 #include "G4VDiscreteProcess.hh"
53 
54 #include "G4Neutron.hh"
55 
57 
59 
60 // Class Description:
61 // Discrete Process -- Boundary Process of Ultra Cold Neutrons.
62 // Reflects/Absorpts UCN at boundaries.
63 // Class inherits publicly from G4VDiscreteProcess.
64 // Class Description - End:
65 
67 // Class Definition
69 
70 
81  };
82 
84 {
85 
86 public:
87 
89  // Constructors and Destructor
91 
92  G4UCNBoundaryProcess(const G4String& processName = "UCNBoundaryProcess",
93  G4ProcessType type = fUCN);
94  virtual ~G4UCNBoundaryProcess();
95 
96 private:
97 
99 
101  // Operators
103 
104  G4UCNBoundaryProcess& operator=(const G4UCNBoundaryProcess &right);
105 
106 public:
107 
109  // Methods
111 
112  G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
113  // Returns true -> 'is applicable' only for an UCN.
114 
115  G4double GetMeanFreePath(const G4Track& aTrack,
116  G4double ,
118 
119  G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
120  const G4Step& aStep);
121 
122 private:
123 
124  G4UCNBoundaryProcessMessenger* fMessenger;
125 
126  G4double neV;
127 
128  G4double kCarTolerance;
129 
130  G4UCNBoundaryProcessStatus theStatus;
131 
132  G4Material* Material1;
133  G4Material* Material2;
134 
135  // the G4UCNMaterialPropertiesTable of PreStepPoint
136  G4UCNMaterialPropertiesTable* aMaterialPropertiesTable1;
137  // the G4UCNMaterialPropertiesTable of PostStepPoint
138  G4UCNMaterialPropertiesTable* aMaterialPropertiesTable2;
139 
140  G4bool UseMicroRoughnessReflection;
141  G4bool DoMicroRoughnessReflection;
142 
144  // Methods
146 
147  G4bool High(G4double , G4double );
148 
149  G4bool Loss(G4double , G4double , G4double );
150 
151  G4bool SpinFlip(G4double );
152 
153  G4double Reflectivity(G4double , G4double );
154 
156 
157  G4double Transmit(G4double, G4double );
158 
159  G4ThreeVector LDiffRefl(G4ThreeVector );
160 
161 public:
162 
165  G4double , G4double );
168  G4double , G4double , G4double& );
169 
170 private:
171 
174  G4ThreeVector MRDiffTrans(G4ThreeVector , G4double , G4double ,
176 
177  G4RotationMatrix GetCoordinateTransformMatrix(G4ThreeVector ,
178  G4ThreeVector );
179 
180  void BoundaryProcessVerbose() const;
181 
182  // Invoke SD for post step point if the photon is 'detected'
183  G4bool InvokeSD(const G4Step* step);
184 
185 private:
186 
187  G4int nNoMPT, nNoMRT, nNoMRCondition;
188  G4int nAbsorption, nEzero, nFlip;
189  G4int aSpecularReflection, bSpecularReflection;
190  G4int bLambertianReflection;
191  G4int aMRDiffuseReflection, bMRDiffuseReflection;
192  G4int nSnellTransmit, mSnellTransmit;
193  G4int aMRDiffuseTransmit;
194 
195  G4double ftheta_o, fphi_o;
196 
197 public:
198 
199  void SetMicroRoughness(G4bool );
201 
203 
204  void BoundaryProcessSummary() const;
205 
207  {aMaterialPropertiesTable1 = MPT;}
208 
210  {aMaterialPropertiesTable2 = MPT;}
211 
212  G4double GetTheta_o() {return ftheta_o;};
213  G4double GetPhi_o() {return fphi_o;};
214 
215 };
216 
218 // Inline methods
220 
221 inline G4bool
223 {
224  return ( &aParticleType == G4Neutron::NeutronDefinition() );
225 }
226 
227 inline
229 {
230  return theStatus;
231 }
232 
233 inline
234 G4bool G4UCNBoundaryProcess::High(G4double Energy, G4double FermiPotDiff)
235 {
236  // Returns true for Energy > Fermi Potential Difference
237 
238  return (Energy > FermiPotDiff);
239 }
240 
241 inline
243 {
244  UseMicroRoughnessReflection = active;
245 }
246 
247 inline
249 {
250  return UseMicroRoughnessReflection;
251 }
252 
253 #endif /* G4UCNBOUNDARYPROCESS_HH */
G4double condition(const G4ErrorSymMatrix &m)
G4UCNBoundaryProcessStatus GetStatus() const
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
int G4int
Definition: G4Types.hh:78
G4UCNBoundaryProcess(const G4String &processName="UCNBoundaryProcess", G4ProcessType type=fUCN)
bool G4bool
Definition: G4Types.hh:79
Definition: G4Step.hh:76
void BoundaryProcessSummary() const
G4ThreeVector MRreflect(G4double, G4ThreeVector, G4ThreeVector, G4double, G4double)
void SetMaterialPropertiesTable1(G4UCNMaterialPropertiesTable *MPT)
void SetMaterialPropertiesTable2(G4UCNMaterialPropertiesTable *MPT)
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *condition)
G4ThreeVector MRreflectHigh(G4double, G4double, G4double, G4ThreeVector, G4ThreeVector, G4double, G4double, G4double &)
double G4double
Definition: G4Types.hh:76
G4ForceCondition
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:99
G4UCNBoundaryProcessStatus
G4ProcessType