Geant4-11
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//
30// Optical Photon Boundary Process Class Definition
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//
63
64#ifndef G4OpBoundaryProcess_h
65#define G4OpBoundaryProcess_h 1
66
67#include "G4OpticalPhoton.hh"
68#include "G4OpticalSurface.hh"
69#include "G4RandomTools.hh"
70#include "G4VDiscreteProcess.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
148
149 private:
152
153 G4bool G4BooleanRand(const G4double prob) const;
154
156 const G4ThreeVector& Normal) const;
157
158 void DielectricMetal();
160
161 void DielectricLUT();
162 void DielectricLUTDAVIS();
163
164 void DielectricDichroic();
165
166 void ChooseReflection();
167 void DoAbsorption();
168 void DoReflection();
169
171 // Returns the incident angle of optical photon
172
173 G4double GetReflectivity(G4double E1_perp, G4double E1_parl,
174 G4double incidentangle, G4double RealRindex,
175 G4double ImaginaryRindex);
176 // Returns the Reflectivity on a metalic surface
177
178 void CalculateReflectivity(void);
179
180 void BoundaryProcessVerbose(void) const;
181
182 // Invoke SD for post step point if the photon is 'detected'
183 G4bool InvokeSD(const G4Step* step);
184
187
190
193
196
198
202
206
208
213
216
220
222
223 size_t idx_dichroicX = 0;
224 size_t idx_dichroicY = 0;
225 size_t idx_rindex1 = 0;
227 size_t idx_reflect = 0;
228 size_t idx_eff = 0;
229 size_t idx_trans = 0;
230 size_t idx_lobe = 0;
231 size_t idx_spike = 0;
232 size_t idx_back = 0;
233 size_t idx_rindex2 = 0;
234 size_t idx_groupvel = 0;
235 size_t idx_rrindex = 0;
236 size_t idx_irindex = 0;
237
239};
240
242// Inline methods
244
246{
247 // Returns a random boolean variable with the specified probability
248 return (G4UniformRand() < prob);
249}
250
252 const G4ParticleDefinition& aParticleType)
253{
254 return (&aParticleType == G4OpticalPhoton::OpticalPhoton());
255}
256
258{
259 return fStatus;
260}
261
263{
264 G4double rand = G4UniformRand();
265 if(rand < fProb_ss)
266 {
269 }
270 else if(rand < fProb_ss + fProb_sl)
271 {
273 }
274 else if(rand < fProb_ss + fProb_sl + fProb_bs)
275 {
277 }
278 else
279 {
281 }
282}
283
285{
287
289 {
290 // EnergyDeposited =/= 0 means: photon has been detected
293 }
294 else
295 {
297 }
298
301
303}
304
306{
308 {
311 }
312 else if(fFinish == ground)
313 {
316 {
318 }
319 // else
320 // complex ref. index to be implemented
323 }
324 else
325 {
330 }
333}
334
335#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
G4OpticalSurfaceFinish fFinish
G4OpBoundaryProcessStatus fStatus
G4OpBoundaryProcess & operator=(const G4OpBoundaryProcess &right)=delete
G4ThreeVector fNewPolarization
G4double GetReflectivity(G4double E1_perp, G4double E1_parl, G4double incidentangle, G4double RealRindex, G4double ImaginaryRindex)
void BoundaryProcessVerbose(void) const
G4OpBoundaryProcess(const G4String &processName="OpBoundary", G4ProcessType type=fOptical)
virtual G4OpBoundaryProcessStatus GetStatus() const
G4OpticalSurface * fOpticalSurface
virtual G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
G4OpticalSurfaceModel fModel
G4ThreeVector GetFacetNormal(const G4ThreeVector &Momentum, const G4ThreeVector &Normal) const
G4OpBoundaryProcess(const G4OpBoundaryProcess &right)=delete
G4Physics2DVector * fDichroicVector
G4MaterialPropertyVector * fRealRIndexMPV
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4MaterialPropertyVector * fImagRIndexMPV
virtual void SetInvokeSD(G4bool)
virtual G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *condition) override
G4bool G4BooleanRand(const G4double prob) const
virtual void PreparePhysicsTable(const G4ParticleDefinition &) override
G4bool InvokeSD(const G4Step *step)
G4ThreeVector fOldPolarization
static G4OpticalPhoton * OpticalPhoton()
Definition: G4Step.hh:62
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:327