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

Last change on this file since 989 was 963, checked in by garnier, 17 years ago

update processes

File size: 10.2 KB
RevLine 
[819]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//
[963]27// $Id: G4OpBoundaryProcess.hh,v 1.18 2008/11/07 17:59:37 gum Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
[819]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//
52// Author: Peter Gumplinger
53// adopted from work by Werner Keil - April 2/96
54// mail: gum@triumf.ca
55//
56// CVS version tag:
57////////////////////////////////////////////////////////////////////////
58
59#ifndef G4OpBoundaryProcess_h
60#define G4OpBoundaryProcess_h 1
61
62/////////////
63// Includes
64/////////////
65
66#include "globals.hh"
67#include "templates.hh"
68#include "geomdefs.hh"
69#include "Randomize.hh"
[963]70
71#include "G4RandomTools.hh"
72#include "G4RandomDirection.hh"
73
[819]74#include "G4Step.hh"
75#include "G4VDiscreteProcess.hh"
76#include "G4DynamicParticle.hh"
77#include "G4Material.hh"
78#include "G4LogicalBorderSurface.hh"
79#include "G4LogicalSkinSurface.hh"
80#include "G4OpticalSurface.hh"
81#include "G4OpticalPhoton.hh"
82#include "G4TransportationManager.hh"
83
84// Class Description:
85// Discrete Process -- reflection/refraction at optical interfaces.
86// Class inherits publicly from G4VDiscreteProcess.
87// Class Description - End:
88
89/////////////////////
90// Class Definition
91/////////////////////
92
93enum G4OpBoundaryProcessStatus { Undefined,
94 FresnelRefraction, FresnelReflection,
95 TotalInternalReflection,
96 LambertianReflection, LobeReflection,
97 SpikeReflection, BackScattering,
98 Absorption, Detection, NotAtBoundary,
99 SameMaterial, StepTooSmall, NoRINDEX };
100
101class G4OpBoundaryProcess : public G4VDiscreteProcess
102{
103
104private:
105
106 //////////////
107 // Operators
108 //////////////
109
110 // G4OpBoundaryProcess& operator=(const G4OpBoundaryProcess &right);
111
112 // G4OpBoundaryProcess(const G4OpBoundaryProcess &right);
113
114public: // Without description
115
116 ////////////////////////////////
117 // Constructors and Destructor
118 ////////////////////////////////
119
120 G4OpBoundaryProcess(const G4String& processName = "OpBoundary",
121 G4ProcessType type = fOptical);
122
123 ~G4OpBoundaryProcess();
124
125 ////////////
126 // Methods
127 ////////////
128
129public: // With description
130
131 G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
132 // Returns true -> 'is applicable' only for an optical photon.
133
134 G4double GetMeanFreePath(const G4Track& ,
135 G4double ,
136 G4ForceCondition* condition);
137 // Returns infinity; i. e. the process does not limit the step,
138 // but sets the 'Forced' condition for the DoIt to be invoked at
139 // every step. However, only at a boundary will any action be
140 // taken.
141
142 G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
143 const G4Step& aStep);
144 // This is the method implementing boundary processes.
145
146 G4OpticalSurfaceModel GetModel() const;
147 // Returns the optical surface mode.
148
149 G4OpBoundaryProcessStatus GetStatus() const;
150 // Returns the current status.
151
152 G4double GetIncidentAngle();
153 // Returns the incident angle of optical photon
154
155 G4double GetReflectivity(G4double E1_perp,
156 G4double E1_parl,
157 G4double incidentangle,
158 G4double RealRindex,
159 G4double ImaginaryRindex);
160 // Returns the Reflectivity on a metalic surface
161
162 void SetModel(G4OpticalSurfaceModel model);
163 // Set the optical surface model to be followed
164 // (glisur || unified).
165
166private:
167
168 G4bool G4BooleanRand(const G4double prob) const;
169
170 G4ThreeVector GetFacetNormal(const G4ThreeVector& Momentum,
171 const G4ThreeVector& Normal) const;
172
173 void DielectricMetal();
174 void DielectricDielectric();
175
176 void ChooseReflection();
177 void DoAbsorption();
178 void DoReflection();
179
180private:
181
182 G4double thePhotonMomentum;
183
184 G4ThreeVector OldMomentum;
185 G4ThreeVector OldPolarization;
186
187 G4ThreeVector NewMomentum;
188 G4ThreeVector NewPolarization;
189
190 G4ThreeVector theGlobalNormal;
191 G4ThreeVector theFacetNormal;
192
193 G4Material* Material1;
194 G4Material* Material2;
195
196 G4OpticalSurface* OpticalSurface;
197
198 G4double Rindex1;
199 G4double Rindex2;
200
201 G4double cost1, cost2, sint1, sint2;
202
203 G4OpBoundaryProcessStatus theStatus;
204
205 G4OpticalSurfaceModel theModel;
206
207 G4OpticalSurfaceFinish theFinish;
208
209 G4double theReflectivity;
210 G4double theEfficiency;
211 G4double prob_sl, prob_ss, prob_bs;
212
213 G4int iTE, iTM;
214
215 G4double kCarTolerance;
216};
217
218////////////////////
219// Inline methods
220////////////////////
221
222inline
223G4bool G4OpBoundaryProcess::G4BooleanRand(const G4double prob) const
224{
225 /* Returns a random boolean variable with the specified probability */
226
227 return (G4UniformRand() < prob);
228}
229
230inline
231G4bool G4OpBoundaryProcess::IsApplicable(const G4ParticleDefinition&
232 aParticleType)
233{
234 return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
235}
236
237inline
238G4OpticalSurfaceModel G4OpBoundaryProcess::GetModel() const
239{
240 return theModel;
241}
242
243inline
244G4OpBoundaryProcessStatus G4OpBoundaryProcess::GetStatus() const
245{
246 return theStatus;
247}
248
249inline
250void G4OpBoundaryProcess::SetModel(G4OpticalSurfaceModel model)
251{
252 theModel = model;
253}
254
255inline
256void G4OpBoundaryProcess::ChooseReflection()
257{
258 G4double rand = G4UniformRand();
259 if ( rand >= 0.0 && rand < prob_ss ) {
260 theStatus = SpikeReflection;
261 theFacetNormal = theGlobalNormal;
262 }
263 else if ( rand >= prob_ss &&
264 rand <= prob_ss+prob_sl) {
265 theStatus = LobeReflection;
266 }
267 else if ( rand > prob_ss+prob_sl &&
268 rand < prob_ss+prob_sl+prob_bs ) {
269 theStatus = BackScattering;
270 }
271 else {
272 theStatus = LambertianReflection;
273 }
274}
275
276inline
277void G4OpBoundaryProcess::DoAbsorption()
278{
279 theStatus = Absorption;
280
281 if ( G4BooleanRand(theEfficiency) ) {
282
283 // EnergyDeposited =/= 0 means: photon has been detected
284 theStatus = Detection;
285 aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum);
286 }
287 else {
288 aParticleChange.ProposeLocalEnergyDeposit(0.0);
289 }
290
291 NewMomentum = OldMomentum;
292 NewPolarization = OldPolarization;
293
294// aParticleChange.ProposeEnergy(0.0);
295 aParticleChange.ProposeTrackStatus(fStopAndKill);
296}
297
298inline
299void G4OpBoundaryProcess::DoReflection()
300{
301 if ( theStatus == LambertianReflection ) {
302
303 NewMomentum = G4LambertianRand(theGlobalNormal);
304 theFacetNormal = (NewMomentum - OldMomentum).unit();
305
306 }
307 else if ( theFinish == ground ) {
308
309 theStatus = LobeReflection;
310 theFacetNormal = GetFacetNormal(OldMomentum,theGlobalNormal);
311 G4double PdotN = OldMomentum * theFacetNormal;
312 NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
313
314 }
315 else {
316
317 theStatus = SpikeReflection;
318 theFacetNormal = theGlobalNormal;
319 G4double PdotN = OldMomentum * theFacetNormal;
320 NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
321
322 }
323 G4double EdotN = OldPolarization * theFacetNormal;
324 NewPolarization = -OldPolarization + (2.*EdotN)*theFacetNormal;
325}
326
327#endif /* G4OpBoundaryProcess_h */
Note: See TracBrowser for help on using the repository browser.