Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4OpBoundaryProcess.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// Optical Photon Boundary Process Class Definition
31////////////////////////////////////////////////////////////////////////
32//
33// File: G4OpBoundaryProcess.hh
34// Description: Discrete Process -- reflection/refraction at
35// optical interfaces
36// Version: 1.1
37// Created: 1997-06-18
38// Modified: 2005-07-28 add G4ProcessType to constructor
39// 1999-10-29 add method and class descriptors
40// 1999-10-10 - Fill NewMomentum/NewPolarization in
41// DoAbsorption. These members need to be
42// filled since DoIt calls
43// aParticleChange.SetMomentumChange etc.
44// upon return (thanks to: Clark McGrew)
45// 2006-11-04 - add capability of calculating the reflectivity
46// off a metal surface by way of a complex index
47// of refraction - Thanks to Sehwook Lee and John
48// Hauptman (Dept. of Physics - Iowa State Univ.)
49// 2009-11-10 - add capability of simulating surface reflections
50// with Look-Up-Tables (LUT) containing measured
51// optical reflectance for a variety of surface
52// treatments - Thanks to Martin Janecek and
53// William Moses (Lawrence Berkeley National Lab.)
54// 2013-06-01 - add the capability of simulating the transmission
55// of a dichronic filter
56// 2017-02-24 - add capability of simulating surface reflections
57// with Look-Up-Tables (LUT) developed in DAVIS
58//
59// Author: Peter Gumplinger
60// adopted from work by Werner Keil - April 2/96
61//
62////////////////////////////////////////////////////////////////////////
63
64#ifndef G4OpBoundaryProcess_h
65#define G4OpBoundaryProcess_h 1
66
67#include "G4RandomTools.hh"
68#include "G4VDiscreteProcess.hh"
69#include "G4OpticalSurface.hh"
70#include "G4OpticalPhoton.hh"
71
73{
115
117{
118 public:
119 explicit G4OpBoundaryProcess(const G4String& processName = "OpBoundary",
120 G4ProcessType type = fOptical);
121 virtual ~G4OpBoundaryProcess();
122
123 virtual G4bool IsApplicable(
124 const G4ParticleDefinition& aParticleType) override;
125 // Returns true -> 'is applicable' only for an optical photon.
126
128 G4ForceCondition* condition) override;
129 // Returns infinity; i. e. the process does not limit the step, but sets the
130 // 'Forced' condition for the DoIt to be invoked at every step. However, only
131 // at a boundary will any action be taken.
132
134 const G4Step& aStep) override;
135 // This is the method implementing boundary processes.
136
137 virtual G4OpBoundaryProcessStatus GetStatus() const;
138 // Returns the current status.
139
140 virtual void SetInvokeSD(G4bool);
141 // Set flag for call to InvokeSD method.
142
143 virtual void PreparePhysicsTable(const G4ParticleDefinition&) override;
144
145 virtual void Initialise();
146
147 private:
148 G4OpBoundaryProcess(const G4OpBoundaryProcess& right) = delete;
149 G4OpBoundaryProcess& operator=(const G4OpBoundaryProcess& right) = delete;
150
151 G4bool G4BooleanRand(const G4double prob) const;
152
153 G4ThreeVector GetFacetNormal(const G4ThreeVector& Momentum,
154 const G4ThreeVector& Normal) const;
155
156 void DielectricMetal();
157 void DielectricDielectric();
158
159 void DielectricLUT();
160 void DielectricLUTDAVIS();
161
162 void DielectricDichroic();
163
164 void ChooseReflection();
165 void DoAbsorption();
166 void DoReflection();
167
168 G4double GetIncidentAngle();
169 // Returns the incident angle of optical photon
170
171 G4double GetReflectivity(G4double E1_perp, G4double E1_parl,
172 G4double incidentangle, G4double RealRindex,
173 G4double ImaginaryRindex);
174 // Returns the Reflectivity on a metalic surface
175
176 void CalculateReflectivity(void);
177
178 void BoundaryProcessVerbose(void) const;
179
180 // Invoke SD for post step point if the photon is 'detected'
181 G4bool InvokeSD(const G4Step* step);
182
183 G4double thePhotonMomentum;
184
185 G4ThreeVector OldMomentum;
186 G4ThreeVector OldPolarization;
187
188 G4ThreeVector NewMomentum;
189 G4ThreeVector NewPolarization;
190
191 G4ThreeVector theGlobalNormal;
192 G4ThreeVector theFacetNormal;
193
194 G4Material* Material1;
195 G4Material* Material2;
196
197 G4OpticalSurface* OpticalSurface;
198
199 G4MaterialPropertyVector* fRealRIndexMPV;
200 G4MaterialPropertyVector* fImagRIndexMPV;
201
202 G4double Rindex1;
203 G4double Rindex2;
204
205 G4double cost1, cost2, sint1, sint2;
206
208
209 G4OpticalSurfaceModel theModel;
210
211 G4OpticalSurfaceFinish theFinish;
212
213 G4double theReflectivity;
214 G4double theEfficiency;
215 G4double theTransmittance;
216
217 G4double theSurfaceRoughness;
218
219 G4double prob_sl, prob_ss, prob_bs;
220
221 G4int iTE, iTM;
222
223 G4double kCarTolerance;
224
225 size_t idx, idy;
226 G4Physics2DVector* DichroicVector;
227
228 G4bool fInvokeSD;
229
230 size_t idx_rindex1 = 0;
231 size_t idx_rindex_surface = 0;
232 size_t idx_reflect = 0;
233 size_t idx_eff = 0;
234 size_t idx_trans = 0;
235 size_t idx_lobe = 0;
236 size_t idx_spike = 0;
237 size_t idx_back = 0;
238 size_t idx_rindex2 = 0;
239 size_t idx_groupvel = 0;
240 size_t idx_rrindex = 0;
241 size_t idx_irindex = 0;
242};
243
244////////////////////
245// Inline methods
246////////////////////
247
248inline G4bool G4OpBoundaryProcess::G4BooleanRand(const G4double prob) const
249{
250 /* Returns a random boolean variable with the specified probability */
251 return (G4UniformRand() < prob);
252}
253
255 const G4ParticleDefinition& aParticleType)
256{
257 return (&aParticleType == G4OpticalPhoton::OpticalPhoton());
258}
259
261{
262 return theStatus;
263}
264
265inline void G4OpBoundaryProcess::SetInvokeSD(G4bool flag) { fInvokeSD = flag; }
266
267inline void G4OpBoundaryProcess::ChooseReflection()
268{
269 G4double rand = G4UniformRand();
270 if(rand >= 0.0 && rand < prob_ss)
271 {
272 theStatus = SpikeReflection;
273 theFacetNormal = theGlobalNormal;
274 }
275 else if(rand >= prob_ss && rand <= prob_ss + prob_sl)
276 {
277 theStatus = LobeReflection;
278 }
279 else if(rand > prob_ss + prob_sl && rand < prob_ss + prob_sl + prob_bs)
280 {
281 theStatus = BackScattering;
282 }
283 else
284 {
285 theStatus = LambertianReflection;
286 }
287}
288
289inline void G4OpBoundaryProcess::DoAbsorption()
290{
291 theStatus = Absorption;
292
293 if(G4BooleanRand(theEfficiency))
294 {
295 // EnergyDeposited =/= 0 means: photon has been detected
296 theStatus = Detection;
298 }
299 else
300 {
302 }
303
304 NewMomentum = OldMomentum;
305 NewPolarization = OldPolarization;
306
308}
309
310inline void G4OpBoundaryProcess::DoReflection()
311{
312 if(theStatus == LambertianReflection)
313 {
314 NewMomentum = G4LambertianRand(theGlobalNormal);
315 theFacetNormal = (NewMomentum - OldMomentum).unit();
316 }
317 else if(theFinish == ground)
318 {
319 theStatus = LobeReflection;
320 if(fRealRIndexMPV && fImagRIndexMPV)
321 {
322 //
323 }
324 else
325 {
326 theFacetNormal = GetFacetNormal(OldMomentum, theGlobalNormal);
327 }
328 G4double PdotN = OldMomentum * theFacetNormal;
329 NewMomentum = OldMomentum - (2. * PdotN) * theFacetNormal;
330 }
331 else
332 {
333 theStatus = SpikeReflection;
334 theFacetNormal = theGlobalNormal;
335 G4double PdotN = OldMomentum * theFacetNormal;
336 NewMomentum = OldMomentum - (2. * PdotN) * theFacetNormal;
337 }
338 G4double EdotN = OldPolarization * theFacetNormal;
339 NewPolarization = -OldPolarization + (2. * EdotN) * theFacetNormal;
340}
341
342#endif /* G4OpBoundaryProcess_h */
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4OpBoundaryProcessStatus
@ PolishedTiOAirReflection
@ GroundTeflonAirReflection
@ EtchedVM2000AirReflection
@ Undefined
@ NotAtBoundary
@ Transmission
@ EtchedVM2000GlueReflection
@ GroundLumirrorGlueReflection
@ LobeReflection
@ GroundTiOAirReflection
@ GroundTyvekAirReflection
@ FresnelReflection
@ NoRINDEX
@ GroundAirReflection
@ EtchedAirReflection
@ PolishedVM2000GlueReflection
@ PolishedAirReflection
@ PolishedTeflonAirReflection
@ SameMaterial
@ EtchedTyvekAirReflection
@ SpikeReflection
@ EtchedLumirrorGlueReflection
@ GroundVM2000AirReflection
@ PolishedTyvekAirReflection
@ PolishedVM2000AirReflection
@ EtchedTiOAirReflection
@ EtchedTeflonAirReflection
@ GroundVM2000GlueReflection
@ BackScattering
@ PolishedLumirrorGlueReflection
@ Absorption
@ TotalInternalReflection
@ StepTooSmall
@ PolishedLumirrorAirReflection
@ EtchedLumirrorAirReflection
@ GroundLumirrorAirReflection
@ LambertianReflection
@ FresnelRefraction
@ Detection
G4OpticalSurfaceModel
G4OpticalSurfaceFinish
@ ground
G4ProcessType
@ fOptical
G4ThreeVector G4LambertianRand(const G4ThreeVector &normal)
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
virtual G4OpBoundaryProcessStatus GetStatus() const
virtual G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
virtual void SetInvokeSD(G4bool)
virtual G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *condition) override
virtual void PreparePhysicsTable(const G4ParticleDefinition &) override
static G4OpticalPhoton * OpticalPhoton()
Definition: G4Step.hh:62
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:327