source: trunk/source/processes/optical/include/G4OpBoundaryProcess.hh@ 1199

Last change on this file since 1199 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 12.3 KB
Line 
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: G4OpBoundaryProcess.hh,v 1.22 2009/11/20 01:06:45 gum Exp $
28// GEANT4 tag $Name: geant4-09-03-cand-01 $
29//
30//
31////////////////////////////////////////////////////////////////////////
32// Optical Photon Boundary Process Class Definition
33////////////////////////////////////////////////////////////////////////
34//
35// File: G4OpBoundaryProcess.hh
36// Description: Discrete Process -- reflection/refraction at
37// optical interfaces
38// Version: 1.1
39// Created: 1997-06-18
40// Modified: 2005-07-28 add G4ProcessType to constructor
41// 1999-10-29 add method and class descriptors
42// 1999-10-10 - Fill NewMomentum/NewPolarization in
43// DoAbsorption. These members need to be
44// filled since DoIt calls
45// aParticleChange.SetMomentumChange etc.
46// upon return (thanks to: Clark McGrew)
47// 2006-11-04 - add capability of calculating the reflectivity
48// off a metal surface by way of a complex index
49// of refraction - Thanks to Sehwook Lee and John
50// Hauptman (Dept. of Physics - Iowa State Univ.)
51// 2009-11-10 - add capability of simulating surface reflections
52// with Look-Up-Tables (LUT) containing measured
53// optical reflectance for a variety of surface
54// treatments - Thanks to Martin Janecek and
55// William Moses (Lawrence Berkeley National Lab.)
56//
57// Author: Peter Gumplinger
58// adopted from work by Werner Keil - April 2/96
59// mail: gum@triumf.ca
60//
61////////////////////////////////////////////////////////////////////////
62
63#ifndef G4OpBoundaryProcess_h
64#define G4OpBoundaryProcess_h 1
65
66/////////////
67// Includes
68/////////////
69
70#include "globals.hh"
71#include "templates.hh"
72#include "geomdefs.hh"
73#include "Randomize.hh"
74
75#include "G4RandomTools.hh"
76#include "G4RandomDirection.hh"
77
78#include "G4Step.hh"
79#include "G4VDiscreteProcess.hh"
80#include "G4DynamicParticle.hh"
81#include "G4Material.hh"
82#include "G4LogicalBorderSurface.hh"
83#include "G4LogicalSkinSurface.hh"
84#include "G4OpticalSurface.hh"
85#include "G4OpticalPhoton.hh"
86#include "G4TransportationManager.hh"
87
88// Class Description:
89// Discrete Process -- reflection/refraction at optical interfaces.
90// Class inherits publicly from G4VDiscreteProcess.
91// Class Description - End:
92
93/////////////////////
94// Class Definition
95/////////////////////
96
97enum G4OpBoundaryProcessStatus { Undefined,
98 FresnelRefraction, FresnelReflection,
99 TotalInternalReflection,
100 LambertianReflection, LobeReflection,
101 SpikeReflection, BackScattering,
102 Absorption, Detection, NotAtBoundary,
103 SameMaterial, StepTooSmall, NoRINDEX,
104 PolishedLumirrorAirReflection,
105 PolishedLumirrorGlueReflection,
106 PolishedAirReflection,
107 PolishedTeflonAirReflection,
108 PolishedTiOAirReflection,
109 PolishedTyvekAirReflection,
110 PolishedVM2000AirReflection,
111 PolishedVM2000GlueReflection,
112 EtchedLumirrorAirReflection,
113 EtchedLumirrorGlueReflection,
114 EtchedAirReflection,
115 EtchedTeflonAirReflection,
116 EtchedTiOAirReflection,
117 EtchedTyvekAirReflection,
118 EtchedVM2000AirReflection,
119 EtchedVM2000GlueReflection,
120 GroundLumirrorAirReflection,
121 GroundLumirrorGlueReflection,
122 GroundAirReflection,
123 GroundTeflonAirReflection,
124 GroundTiOAirReflection,
125 GroundTyvekAirReflection,
126 GroundVM2000AirReflection,
127 GroundVM2000GlueReflection };
128
129class G4OpBoundaryProcess : public G4VDiscreteProcess
130{
131
132private:
133
134 //////////////
135 // Operators
136 //////////////
137
138 // G4OpBoundaryProcess& operator=(const G4OpBoundaryProcess &right);
139
140 // G4OpBoundaryProcess(const G4OpBoundaryProcess &right);
141
142public: // Without description
143
144 ////////////////////////////////
145 // Constructors and Destructor
146 ////////////////////////////////
147
148 G4OpBoundaryProcess(const G4String& processName = "OpBoundary",
149 G4ProcessType type = fOptical);
150
151 ~G4OpBoundaryProcess();
152
153 ////////////
154 // Methods
155 ////////////
156
157public: // With description
158
159 G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
160 // Returns true -> 'is applicable' only for an optical photon.
161
162 G4double GetMeanFreePath(const G4Track& ,
163 G4double ,
164 G4ForceCondition* condition);
165 // Returns infinity; i. e. the process does not limit the step,
166 // but sets the 'Forced' condition for the DoIt to be invoked at
167 // every step. However, only at a boundary will any action be
168 // taken.
169
170 G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
171 const G4Step& aStep);
172 // This is the method implementing boundary processes.
173
174 G4OpticalSurfaceModel GetModel() const;
175 // Returns the optical surface mode.
176
177 G4OpBoundaryProcessStatus GetStatus() const;
178 // Returns the current status.
179
180 void SetModel(G4OpticalSurfaceModel model);
181 // Set the optical surface model to be followed
182 // (glisur || unified || LUT).
183
184private:
185
186 G4bool G4BooleanRand(const G4double prob) const;
187
188 G4ThreeVector GetFacetNormal(const G4ThreeVector& Momentum,
189 const G4ThreeVector& Normal) const;
190
191 void DielectricMetal();
192 void DielectricDielectric();
193 void DielectricLUT();
194
195 void ChooseReflection();
196 void DoAbsorption();
197 void DoReflection();
198
199 G4double GetIncidentAngle();
200 // Returns the incident angle of optical photon
201
202 G4double GetReflectivity(G4double E1_perp,
203 G4double E1_parl,
204 G4double incidentangle,
205 G4double RealRindex,
206 G4double ImaginaryRindex);
207 // Returns the Reflectivity on a metalic surface
208
209 void CalculateReflectivity(void);
210
211 void BoundaryProcessVerbose(void) const;
212
213private:
214
215 G4double thePhotonMomentum;
216
217 G4ThreeVector OldMomentum;
218 G4ThreeVector OldPolarization;
219
220 G4ThreeVector NewMomentum;
221 G4ThreeVector NewPolarization;
222
223 G4ThreeVector theGlobalNormal;
224 G4ThreeVector theFacetNormal;
225
226 G4Material* Material1;
227 G4Material* Material2;
228
229 G4OpticalSurface* OpticalSurface;
230
231 G4MaterialPropertyVector* PropertyPointer;
232 G4MaterialPropertyVector* PropertyPointer1;
233 G4MaterialPropertyVector* PropertyPointer2;
234
235 G4double Rindex1;
236 G4double Rindex2;
237
238 G4double cost1, cost2, sint1, sint2;
239
240 G4OpBoundaryProcessStatus theStatus;
241
242 G4OpticalSurfaceModel theModel;
243
244 G4OpticalSurfaceFinish theFinish;
245
246 G4double theReflectivity;
247 G4double theEfficiency;
248 G4double prob_sl, prob_ss, prob_bs;
249
250 G4int iTE, iTM;
251
252 G4double kCarTolerance;
253
254};
255
256////////////////////
257// Inline methods
258////////////////////
259
260inline
261G4bool G4OpBoundaryProcess::G4BooleanRand(const G4double prob) const
262{
263 /* Returns a random boolean variable with the specified probability */
264
265 return (G4UniformRand() < prob);
266}
267
268inline
269G4bool G4OpBoundaryProcess::IsApplicable(const G4ParticleDefinition&
270 aParticleType)
271{
272 return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
273}
274
275inline
276G4OpticalSurfaceModel G4OpBoundaryProcess::GetModel() const
277{
278 return theModel;
279}
280
281inline
282G4OpBoundaryProcessStatus G4OpBoundaryProcess::GetStatus() const
283{
284 return theStatus;
285}
286
287inline
288void G4OpBoundaryProcess::SetModel(G4OpticalSurfaceModel model)
289{
290 theModel = model;
291}
292
293inline
294void G4OpBoundaryProcess::ChooseReflection()
295{
296 G4double rand = G4UniformRand();
297 if ( rand >= 0.0 && rand < prob_ss ) {
298 theStatus = SpikeReflection;
299 theFacetNormal = theGlobalNormal;
300 }
301 else if ( rand >= prob_ss &&
302 rand <= prob_ss+prob_sl) {
303 theStatus = LobeReflection;
304 }
305 else if ( rand > prob_ss+prob_sl &&
306 rand < prob_ss+prob_sl+prob_bs ) {
307 theStatus = BackScattering;
308 }
309 else {
310 theStatus = LambertianReflection;
311 }
312}
313
314inline
315void G4OpBoundaryProcess::DoAbsorption()
316{
317 theStatus = Absorption;
318
319 if ( G4BooleanRand(theEfficiency) ) {
320
321 // EnergyDeposited =/= 0 means: photon has been detected
322 theStatus = Detection;
323 aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum);
324 }
325 else {
326 aParticleChange.ProposeLocalEnergyDeposit(0.0);
327 }
328
329 NewMomentum = OldMomentum;
330 NewPolarization = OldPolarization;
331
332// aParticleChange.ProposeEnergy(0.0);
333 aParticleChange.ProposeTrackStatus(fStopAndKill);
334}
335
336inline
337void G4OpBoundaryProcess::DoReflection()
338{
339 if ( theStatus == LambertianReflection ) {
340
341 NewMomentum = G4LambertianRand(theGlobalNormal);
342 theFacetNormal = (NewMomentum - OldMomentum).unit();
343
344 }
345 else if ( theFinish == ground ) {
346
347 theStatus = LobeReflection;
348 if ( PropertyPointer1 && PropertyPointer2 ){
349 } else {
350 theFacetNormal =
351 GetFacetNormal(OldMomentum,theGlobalNormal);
352 }
353 G4double PdotN = OldMomentum * theFacetNormal;
354 NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
355
356 }
357 else {
358
359 theStatus = SpikeReflection;
360 theFacetNormal = theGlobalNormal;
361 G4double PdotN = OldMomentum * theFacetNormal;
362 NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
363
364 }
365 G4double EdotN = OldPolarization * theFacetNormal;
366 NewPolarization = -OldPolarization + (2.*EdotN)*theFacetNormal;
367}
368
369#endif /* G4OpBoundaryProcess_h */
Note: See TracBrowser for help on using the repository browser.