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

Last change on this file since 1071 was 1055, checked in by garnier, 15 years ago

maj sur la beta de geant 4.9.3

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