Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4OpBoundaryProcess.hh
이 파일의 문서화 페이지로 가기
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 104458 2017-05-31 15:54:11Z gcosmo $
28 //
29 //
31 // Optical Photon Boundary Process Class Definition
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 // 2013-06-01 - add the capability of simulating the transmission
56 // of a dichronic filter
57 // 2017-02-24 - add capability of simulating surface reflections
58 // with Look-Up-Tables (LUT) developed in DAVIS
59 //
60 // Author: Peter Gumplinger
61 // adopted from work by Werner Keil - April 2/96
62 // mail: gum@triumf.ca
63 //
65 
66 #ifndef G4OpBoundaryProcess_h
67 #define G4OpBoundaryProcess_h 1
68 
70 // Includes
72 
73 #include "globals.hh"
74 #include "templates.hh"
75 #include "geomdefs.hh"
76 #include "Randomize.hh"
77 
78 #include "G4RandomTools.hh"
79 #include "G4RandomDirection.hh"
80 
81 #include "G4Step.hh"
82 #include "G4VDiscreteProcess.hh"
83 #include "G4DynamicParticle.hh"
84 #include "G4Material.hh"
86 #include "G4LogicalSkinSurface.hh"
87 #include "G4OpticalSurface.hh"
88 #include "G4OpticalPhoton.hh"
90 
91 // Class Description:
92 // Discrete Process -- reflection/refraction at optical interfaces.
93 // Class inherits publicly from G4VDiscreteProcess.
94 // Class Description - End:
95 
97 // Class Definition
99 
132 
134 {
135 
136 public:
137 
139  // Constructors and Destructor
141 
142  G4OpBoundaryProcess(const G4String& processName = "OpBoundary",
143  G4ProcessType type = fOptical);
145 
146 private:
147 
149 
151  // Operators
153 
155 
156 public:
157 
159  // Methods
161 
162  G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
163  // Returns true -> 'is applicable' only for an optical photon.
164 
166  G4double ,
168  // Returns infinity; i. e. the process does not limit the step,
169  // but sets the 'Forced' condition for the DoIt to be invoked at
170  // every step. However, only at a boundary will any action be
171  // taken.
172 
173  G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
174  const G4Step& aStep);
175  // This is the method implementing boundary processes.
176 
178  // Returns the current status.
179 
180  void SetInvokeSD(G4bool );
181  // Set flag for call to InvokeSD method.
182 
183 private:
184 
185  G4bool G4BooleanRand(const G4double prob) const;
186 
188  const G4ThreeVector& Normal) const;
189 
190  void DielectricMetal();
191  void DielectricDielectric();
192 
193  void DielectricLUT();
194  void DielectricLUTDAVIS();
195 
196  void DielectricDichroic();
197 
198  void ChooseReflection();
199  void DoAbsorption();
200  void DoReflection();
201 
203  // Returns the incident angle of optical photon
204 
206  G4double E1_parl,
207  G4double incidentangle,
208  G4double RealRindex,
209  G4double ImaginaryRindex);
210  // Returns the Reflectivity on a metalic surface
211 
212  void CalculateReflectivity(void);
213 
214  void BoundaryProcessVerbose(void) const;
215 
216  // Invoke SD for post step point if the photon is 'detected'
217  G4bool InvokeSD(const G4Step* step);
218 
219 private:
220 
222 
225 
228 
231 
234 
236 
240 
243 
245 
247 
249 
251 
255 
257 
259 
261 
263 
264  size_t idx, idy;
266 
268 };
269 
271 // Inline methods
273 
274 inline
276 {
277  /* Returns a random boolean variable with the specified probability */
278 
279  return (G4UniformRand() < prob);
280 }
281 
282 inline
284  aParticleType)
285 {
286  return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
287 }
288 
289 inline
291 {
292  return theStatus;
293 }
294 
295 inline
297 {
298  fInvokeSD = flag;
299 }
300 
301 inline
303 {
304  G4double rand = G4UniformRand();
305  if ( rand >= 0.0 && rand < prob_ss ) {
308  }
309  else if ( rand >= prob_ss &&
310  rand <= prob_ss+prob_sl) {
312  }
313  else if ( rand > prob_ss+prob_sl &&
314  rand < prob_ss+prob_sl+prob_bs ) {
316  }
317  else {
319  }
320 }
321 
322 inline
324 {
326 
327  if ( G4BooleanRand(theEfficiency) ) {
328 
329  // EnergyDeposited =/= 0 means: photon has been detected
332  }
333  else {
335  }
336 
339 
340 // aParticleChange.ProposeEnergy(0.0);
342 }
343 
344 inline
346 {
347  if ( theStatus == LambertianReflection ) {
348 
351 
352  }
353  else if ( theFinish == ground ) {
354 
357  } else {
360  }
362  NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
363 
364  }
365  else {
366 
370  NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
371 
372  }
374  NewPolarization = -OldPolarization + (2.*EdotN)*theFacetNormal;
375 }
376 
377 #endif /* G4OpBoundaryProcess_h */
G4OpticalSurfaceModel theModel
G4OpBoundaryProcessStatus
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
G4double GetReflectivity(G4double E1_perp, G4double E1_parl, G4double incidentangle, G4double RealRindex, G4double ImaginaryRindex)
G4double condition(const G4ErrorSymMatrix &m)
G4OpBoundaryProcess(const G4String &processName="OpBoundary", G4ProcessType type=fOptical)
G4MaterialPropertyVector * PropertyPointer2
G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *condition)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4Physics2DVector * DichroicVector
G4MaterialPropertyVector * PropertyPointer
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4ThreeVector G4LambertianRand(const G4ThreeVector &normal)
G4OpticalSurfaceFinish theFinish
G4ProcessType
G4OpticalSurface * OpticalSurface
#define G4UniformRand()
Definition: Randomize.hh:53
static G4OpticalPhoton * OpticalPhoton()
Definition: G4Step.hh:76
G4OpBoundaryProcessStatus theStatus
G4OpticalSurfaceFinish
int G4int
Definition: G4Types.hh:78
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
G4ForceCondition
G4bool G4BooleanRand(const G4double prob) const
G4MaterialPropertyVector * PropertyPointer1
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4OpBoundaryProcess & operator=(const G4OpBoundaryProcess &right)
G4ThreeVector GetFacetNormal(const G4ThreeVector &Momentum, const G4ThreeVector &Normal) const
void ProposeTrackStatus(G4TrackStatus status)
G4bool InvokeSD(const G4Step *step)
G4OpBoundaryProcessStatus GetStatus() const
void BoundaryProcessVerbose(void) const
G4OpticalSurfaceModel