Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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// $Id$
28//
29//
30////////////////////////////////////////////////////////////////////////
31// Optical Photon Boundary Process Class Definition
32////////////////////////////////////////////////////////////////////////
33//
34// File: G4OpBoundaryProcess.hh
35// Description: Discrete Process -- reflection/refraction at
36// optical interfaces
37// Version: 1.1
38// Created: 1997-06-18
39// Modified: 2005-07-28 add G4ProcessType to constructor
40// 1999-10-29 add method and class descriptors
41// 1999-10-10 - Fill NewMomentum/NewPolarization in
42// DoAbsorption. These members need to be
43// filled since DoIt calls
44// aParticleChange.SetMomentumChange etc.
45// upon return (thanks to: Clark McGrew)
46// 2006-11-04 - add capability of calculating the reflectivity
47// off a metal surface by way of a complex index
48// of refraction - Thanks to Sehwook Lee and John
49// Hauptman (Dept. of Physics - Iowa State Univ.)
50// 2009-11-10 - add capability of simulating surface reflections
51// with Look-Up-Tables (LUT) containing measured
52// optical reflectance for a variety of surface
53// treatments - Thanks to Martin Janecek and
54// William Moses (Lawrence Berkeley National Lab.)
55//
56// Author: Peter Gumplinger
57// adopted from work by Werner Keil - April 2/96
58// mail: gum@triumf.ca
59//
60////////////////////////////////////////////////////////////////////////
61
62#ifndef G4OpBoundaryProcess_h
63#define G4OpBoundaryProcess_h 1
64
65/////////////
66// Includes
67/////////////
68
69#include "globals.hh"
70#include "templates.hh"
71#include "geomdefs.hh"
72#include "Randomize.hh"
73
74#include "G4RandomTools.hh"
75#include "G4RandomDirection.hh"
76
77#include "G4Step.hh"
78#include "G4VDiscreteProcess.hh"
79#include "G4DynamicParticle.hh"
80#include "G4Material.hh"
83#include "G4OpticalSurface.hh"
84#include "G4OpticalPhoton.hh"
86
87// Class Description:
88// Discrete Process -- reflection/refraction at optical interfaces.
89// Class inherits publicly from G4VDiscreteProcess.
90// Class Description - End:
91
92/////////////////////
93// Class Definition
94/////////////////////
95
127
129{
130
131public:
132
133 ////////////////////////////////
134 // Constructors and Destructor
135 ////////////////////////////////
136
137 G4OpBoundaryProcess(const G4String& processName = "OpBoundary",
138 G4ProcessType type = fOptical);
140
141private:
142
144
145 //////////////
146 // Operators
147 //////////////
148
149 G4OpBoundaryProcess& operator=(const G4OpBoundaryProcess &right);
150
151public:
152
153 ////////////
154 // Methods
155 ////////////
156
157 G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
158 // Returns true -> 'is applicable' only for an optical photon.
159
161 G4double ,
163 // Returns infinity; i. e. the process does not limit the step,
164 // but sets the 'Forced' condition for the DoIt to be invoked at
165 // every step. However, only at a boundary will any action be
166 // taken.
167
169 const G4Step& aStep);
170 // This is the method implementing boundary processes.
171
173 // Returns the current status.
174
175private:
176
177 G4bool G4BooleanRand(const G4double prob) const;
178
179 G4ThreeVector GetFacetNormal(const G4ThreeVector& Momentum,
180 const G4ThreeVector& Normal) const;
181
182 void DielectricMetal();
183 void DielectricDielectric();
184 void DielectricLUT();
185
186 void ChooseReflection();
187 void DoAbsorption();
188 void DoReflection();
189
190 G4double GetIncidentAngle();
191 // Returns the incident angle of optical photon
192
193 G4double GetReflectivity(G4double E1_perp,
194 G4double E1_parl,
195 G4double incidentangle,
196 G4double RealRindex,
197 G4double ImaginaryRindex);
198 // Returns the Reflectivity on a metalic surface
199
200 void CalculateReflectivity(void);
201
202 void BoundaryProcessVerbose(void) const;
203
204private:
205
206 G4double thePhotonMomentum;
207
208 G4ThreeVector OldMomentum;
209 G4ThreeVector OldPolarization;
210
211 G4ThreeVector NewMomentum;
212 G4ThreeVector NewPolarization;
213
214 G4ThreeVector theGlobalNormal;
215 G4ThreeVector theFacetNormal;
216
217 G4Material* Material1;
218 G4Material* Material2;
219
220 G4OpticalSurface* OpticalSurface;
221
222 G4MaterialPropertyVector* PropertyPointer;
223 G4MaterialPropertyVector* PropertyPointer1;
224 G4MaterialPropertyVector* PropertyPointer2;
225
226 G4double Rindex1;
227 G4double Rindex2;
228
229 G4double cost1, cost2, sint1, sint2;
230
232
233 G4OpticalSurfaceModel theModel;
234
235 G4OpticalSurfaceFinish theFinish;
236
237 G4double theReflectivity;
238 G4double theEfficiency;
239 G4double theTransmittance;
240 G4double prob_sl, prob_ss, prob_bs;
241
242 G4int iTE, iTM;
243
244 G4double kCarTolerance;
245
246};
247
248////////////////////
249// Inline methods
250////////////////////
251
252inline
253G4bool G4OpBoundaryProcess::G4BooleanRand(const G4double prob) const
254{
255 /* Returns a random boolean variable with the specified probability */
256
257 return (G4UniformRand() < prob);
258}
259
260inline
262 aParticleType)
263{
264 return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
265}
266
267inline
269{
270 return theStatus;
271}
272
273inline
274void G4OpBoundaryProcess::ChooseReflection()
275{
276 G4double rand = G4UniformRand();
277 if ( rand >= 0.0 && rand < prob_ss ) {
278 theStatus = SpikeReflection;
279 theFacetNormal = theGlobalNormal;
280 }
281 else if ( rand >= prob_ss &&
282 rand <= prob_ss+prob_sl) {
283 theStatus = LobeReflection;
284 }
285 else if ( rand > prob_ss+prob_sl &&
286 rand < prob_ss+prob_sl+prob_bs ) {
287 theStatus = BackScattering;
288 }
289 else {
290 theStatus = LambertianReflection;
291 }
292}
293
294inline
295void G4OpBoundaryProcess::DoAbsorption()
296{
297 theStatus = Absorption;
298
299 if ( G4BooleanRand(theEfficiency) ) {
300
301 // EnergyDeposited =/= 0 means: photon has been detected
302 theStatus = Detection;
304 }
305 else {
307 }
308
309 NewMomentum = OldMomentum;
310 NewPolarization = OldPolarization;
311
312// aParticleChange.ProposeEnergy(0.0);
314}
315
316inline
317void G4OpBoundaryProcess::DoReflection()
318{
319 if ( theStatus == LambertianReflection ) {
320
321 NewMomentum = G4LambertianRand(theGlobalNormal);
322 theFacetNormal = (NewMomentum - OldMomentum).unit();
323
324 }
325 else if ( theFinish == ground ) {
326
327 theStatus = LobeReflection;
328 if ( PropertyPointer1 && PropertyPointer2 ){
329 } else {
330 theFacetNormal =
331 GetFacetNormal(OldMomentum,theGlobalNormal);
332 }
333 G4double PdotN = OldMomentum * theFacetNormal;
334 NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
335
336 }
337 else {
338
339 theStatus = SpikeReflection;
340 theFacetNormal = theGlobalNormal;
341 G4double PdotN = OldMomentum * theFacetNormal;
342 NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
343
344 }
345 G4double EdotN = OldPolarization * theFacetNormal;
346 NewPolarization = -OldPolarization + (2.*EdotN)*theFacetNormal;
347}
348
349#endif /* G4OpBoundaryProcess_h */
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4OpBoundaryProcessStatus
@ PolishedTiOAirReflection
@ GroundTeflonAirReflection
@ EtchedVM2000AirReflection
@ Undefined
@ NotAtBoundary
@ EtchedVM2000GlueReflection
@ GroundLumirrorGlueReflection
@ LobeReflection
@ GroundTiOAirReflection
@ GroundTyvekAirReflection
@ FresnelReflection
@ 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
G4OpticalSurfaceModel
G4OpticalSurfaceFinish
@ ground
G4ProcessType
@ fOptical
G4ThreeVector G4LambertianRand(const G4ThreeVector &normal)
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4UniformRand()
Definition: Randomize.hh:53
G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *condition)
G4OpBoundaryProcessStatus GetStatus() const
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
static G4OpticalPhoton * OpticalPhoton()
Definition: G4Step.hh:78
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289