Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
28//
29////////////////////////////////////////////////////////////////////////
30// Ultra Cold Neutron (UCN) Boundary Process Class Definition
31////////////////////////////////////////////////////////////////////////
32//
33// File: G4UCNBoundaryProcess.hh
34// Description: Discrete Process -- Boundary Process of UCN
35// Version: 1.0
36// Created: 2014-06-04
37// Author: Peter Gumplinger
38// Adopted from: UCNMaterialBoundary by Peter Fierlinger 4.9.2004
39// Updated:
40// mail: [email protected]
41//
42////////////////////////////////////////////////////////////////////////
43
44#ifndef G4UCNBOUNDARYPROCESS_HH
45#define G4UCNBOUNDARYPROCESS_HH 1
46
47/////////////
48// Includes
49/////////////
50
51#include "G4VDiscreteProcess.hh"
52
53#include "G4Neutron.hh"
54
56
58
59// Class Description:
60// Discrete Process -- Boundary Process of Ultra Cold Neutrons.
61// Reflects/Absorpts UCN at boundaries.
62// Class inherits publicly from G4VDiscreteProcess.
63// Class Description - End:
64
65/////////////////////
66// Class Definition
67/////////////////////
68
69
80 };
81
83{
84
85public:
86
87 ////////////////////////////////
88 // Constructors and Destructor
89 ////////////////////////////////
90
91 G4UCNBoundaryProcess(const G4String& processName = "UCNBoundaryProcess",
92 G4ProcessType type = fUCN);
93 virtual ~G4UCNBoundaryProcess();
94
95private:
96
98
99 //////////////
100 // Operators
101 //////////////
102
103 G4UCNBoundaryProcess& operator=(const G4UCNBoundaryProcess &right);
104
105public:
106
107 ////////////
108 // Methods
109 ///////////
110
111 G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
112 // Returns true -> 'is applicable' only for an UCN.
113
114 G4double GetMeanFreePath(const G4Track& aTrack,
115 G4double ,
117
119 const G4Step& aStep);
120
121private:
122
124
125 G4double neV;
126
127 G4double kCarTolerance;
128
130
131 G4Material* Material1;
132 G4Material* Material2;
133
134 // the G4UCNMaterialPropertiesTable of PreStepPoint
135 G4UCNMaterialPropertiesTable* aMaterialPropertiesTable1;
136 // the G4UCNMaterialPropertiesTable of PostStepPoint
137 G4UCNMaterialPropertiesTable* aMaterialPropertiesTable2;
138
139 G4bool UseMicroRoughnessReflection;
140 G4bool DoMicroRoughnessReflection;
141
142 ////////////
143 // Methods
144 ///////////
145
146 G4bool High(G4double , G4double );
147
148 G4bool Loss(G4double , G4double , G4double );
149
150 G4bool SpinFlip(G4double );
151
152 G4double Reflectivity(G4double , G4double );
153
155
156 G4double Transmit(G4double, G4double );
157
158 G4ThreeVector LDiffRefl(G4ThreeVector );
159
160public:
161
164 G4double , G4double );
168
169private:
170
175
176 G4RotationMatrix GetCoordinateTransformMatrix(G4ThreeVector ,
178
179 void BoundaryProcessVerbose() const;
180
181 // Invoke SD for post step point if the photon is 'detected'
182 G4bool InvokeSD(const G4Step* step);
183
184private:
185
186 G4int nNoMPT, nNoMRT, nNoMRCondition;
187 G4int nAbsorption, nEzero, nFlip;
188 G4int aSpecularReflection, bSpecularReflection;
189 G4int bLambertianReflection;
190 G4int aMRDiffuseReflection, bMRDiffuseReflection;
191 G4int nSnellTransmit, mSnellTransmit;
192 G4int aMRDiffuseTransmit;
193
194 G4double ftheta_o, fphi_o;
195
196public:
197
200
202
203 void BoundaryProcessSummary() const;
204
206 {aMaterialPropertiesTable1 = MPT;}
207
209 {aMaterialPropertiesTable2 = MPT;}
210
211 G4double GetTheta_o() {return ftheta_o;};
212 G4double GetPhi_o() {return fphi_o;};
213
214};
215
216////////////////////
217// Inline methods
218////////////////////
219
220inline G4bool
222{
223 return ( &aParticleType == G4Neutron::NeutronDefinition() );
224}
225
226inline
228{
229 return theStatus;
230}
231
232inline
233G4bool G4UCNBoundaryProcess::High(G4double Energy, G4double FermiPotDiff)
234{
235 // Returns true for Energy > Fermi Potential Difference
236
237 return (Energy > FermiPotDiff);
238}
239
240inline
242{
243 UseMicroRoughnessReflection = active;
244}
245
246inline
248{
249 return UseMicroRoughnessReflection;
250}
251
252#endif /* G4UCNBOUNDARYPROCESS_HH */
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4ProcessType
@ fUCN
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4UCNBoundaryProcessStatus
@ NotAtBoundary
@ NoMRCondition
@ SnellTransmit
@ MRDiffuseReflection
@ SameMaterial
@ SpecularReflection
@ Absorption
@ StepTooSmall
@ MRDiffuseTransmit
@ LambertianReflection
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:98
Definition: G4Step.hh:62
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *condition)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
G4ThreeVector MRreflectHigh(G4double, G4double, G4double, G4ThreeVector, G4ThreeVector, G4double, G4double, G4double &)
void SetMaterialPropertiesTable2(G4UCNMaterialPropertiesTable *MPT)
G4UCNBoundaryProcessStatus GetStatus() const
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
void SetMaterialPropertiesTable1(G4UCNMaterialPropertiesTable *MPT)
G4ThreeVector MRreflect(G4double, G4ThreeVector, G4ThreeVector, G4double, G4double)