G4StationSimulator/G4TankOpBoundaryProcess.h
Go to the documentation of this file.
1 // ********************************************************************
2 // * License and Disclaimer *
3 // * *
4 // * The Geant4 software is copyright of the Copyright Holders of *
5 // * the Geant4 Collaboration. It is provided under the terms and *
6 // * conditions of the Geant4 Software License, included in the file *
7 // * LICENSE and available at http://cern.ch/geant4/license . These *
8 // * include a list of copyright holders. *
9 // * *
10 // * Neither the authors of this software system, nor their employing *
11 // * institutes,nor the agencies providing financial support for this *
12 // * work make any representation or warranty, express or implied, *
13 // * regarding this software system or assume any liability for its *
14 // * use. Please see the license in the file LICENSE and URL above *
15 // * for the full disclaimer and the limitation of liability. *
16 // * *
17 // * This code implementation is the result of the scientific and *
18 // * technical work of the GEANT4 collaboration. *
19 // * By using, copying, modifying or distributing the software (or *
20 // * any work based on the software) you agree to acknowledge its *
21 // * use in resulting scientific publications, and indicate your *
22 // * acceptance of all terms of the Geant4 Software license. *
23 // ********************************************************************
24 //
25 //
26 // GEANT4 tag $Name: geant4-09-00-patch-01 $
27 //
28 //
30 // Optical Photon Boundary Process Class Definition
32 //
33 // File: G4TankOpBoundaryProcess.hh
34 // Description: Discrete Process -- reflection/refraction at
35 // optical interfaces
36 // Version: 1.1
37 // Created: 1997-06-18
38 // Modified: 2005-07-28 add G4ProcessType to constructor
39 // 1999-10-29 add method and class descriptors
40 // 1999-10-10 - Fill NewMomentum/NewPolarization in
41 // DoAbsorption. These members need to be
42 // filled since DoIt calls
43 // aParticleChange.SetMomentumChange etc.
44 // upon return (thanks to: Clark McGrew)
45 //
46 // Author: Peter Gumplinger
47 // adopted from work by Werner Keil - April 2/96
48 // mail: gum@triumf.ca
49 //
50 // CVS version tag:
52 
53 #ifndef _G4StationSimulatorOG_G4TankOpBoundaryProcess_h_
54 #define _G4StationSimulatorOG_G4TankOpBoundaryProcess_h_
55 
56 #include "globals.hh"
57 #include "templates.hh"
58 #include "geomdefs.hh"
59 #include "Randomize.hh"
60 #include "G4Step.hh"
61 #include "G4VDiscreteProcess.hh"
62 #include "G4DynamicParticle.hh"
63 #include "G4Material.hh"
64 #include "G4LogicalBorderSurface.hh"
65 #include "G4LogicalSkinSurface.hh"
66 #include "G4OpticalSurface.hh"
67 #include "G4OpticalPhoton.hh"
68 #include "G4TransportationManager.hh"
69 #include "G4PhysicalConstants.hh"
70 
71 // Discrete Process -- reflection/refraction at optical interfaces.
72 // Class inherits publicly from G4VDiscreteProcess.
73 
74 
75 namespace G4StationSimulatorOG {
76 
92  };
93 
94 
95  class G4TankOpBoundaryProcess : public G4VDiscreteProcess {
96 
97  public:
98  G4TankOpBoundaryProcess(const G4String& processName = "OpBoundary",
99  const G4ProcessType type = fOptical);
100 
102 
103  public:
104  // Returns true -> 'is applicable' only for an optical photon.
105  virtual G4bool IsApplicable(const G4ParticleDefinition& particleType) override
106  { return &particleType == G4OpticalPhoton::OpticalPhoton(); }
107 
108  // Returns infinity; i. e. the process does not limit the step, but sets
109  // the 'Forced' condition for the DoIt to be invoked at every step.
110  // However, only at a boundary will any action be taken.
111  virtual G4double GetMeanFreePath(const G4Track& /*track*/,
112  const G4double /*previousStepSize*/,
113  G4ForceCondition* const condition) override
114  { *condition = Forced; return DBL_MAX; }
115 
116  // This is the method implementing boundary processes.
117  virtual G4VParticleChange* PostStepDoIt(const G4Track& track, const G4Step& step) override;
118 
119  // Set the optical surface model to be followed (glisur || unified).
120  void SetModel(const G4OpticalSurfaceModel model) { fModel = model; }
121 
122  private:
123  G4bool G4BooleanRand(const G4double prob) const { return G4UniformRand() < prob; }
124  // Returns a random boolean variable with the specified probability
125 
126  G4ThreeVector G4IsotropicRand() const;
127 
128  G4ThreeVector G4LambertianRand(const G4ThreeVector& normal) const;
129 
130  G4ThreeVector G4PlaneVectorRand(const G4ThreeVector& normal) const;
131 
132  G4ThreeVector GetFacetNormal(const G4ThreeVector& momentum, const G4ThreeVector& normal) const;
133 
134  void DielectricMetal();
135 
136  void DielectricDielectric();
137 
138  void ChooseReflection();
139 
140  void DoAbsorption();
141 
142  void DoReflection();
143 
144  private:
145  G4double fPhotonMomentum = 0;
146 
147  G4ThreeVector fOldMomentum;
148  G4ThreeVector fOldPolarization;
149 
150  G4ThreeVector fNewMomentum;
151  G4ThreeVector fNewPolarization;
152 
153  G4ThreeVector fGlobalNormal;
154  G4ThreeVector fFacetNormal;
155 
156  G4Material* fMaterial1 = nullptr;
157  G4Material* fMaterial2 = nullptr;
158 
159  G4OpticalSurface* fOpticalSurface = nullptr;
160 
161  G4double fRIndex1 = 0;
162  G4double fRIndex2 = 0;
163 
164  /*G4double cost1 = 0;
165  G4double cost2 = 0;
166  G4double sint1 = 0;
167  G4double sint2 = 0;*/
168 
170  G4OpticalSurfaceModel fModel = glisur;
171  G4OpticalSurfaceFinish fFinish = polished;
172 
173  G4double fReflectivity = 0;
174  G4double fEfficiency = 0;
175  G4double fProb_sl = 0;
176  G4double fProb_ss = 0;
177  G4double fProb_bs = 0;
178  G4double kCarTolerance = 0;
179 
180  };
181 
182 }
183 
184 
185 #endif
G4ThreeVector G4LambertianRand(const G4ThreeVector &normal) const
G4ThreeVector GetFacetNormal(const G4ThreeVector &momentum, const G4ThreeVector &normal) const
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step) override
G4TankOpBoundaryProcess(const G4String &processName="OpBoundary", const G4ProcessType type=fOptical)
virtual G4bool IsApplicable(const G4ParticleDefinition &particleType) override
G4ThreeVector G4PlaneVectorRand(const G4ThreeVector &normal) const
virtual G4double GetMeanFreePath(const G4Track &, const G4double, G4ForceCondition *const condition) override

, generated on Tue Sep 26 2023.