Geant4 Cross Reference |
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 // 28 // 29 //////////////////////////////////////////////////////////////////////// 30 // Optical Photon Boundary Process Class Definition 31 //////////////////////////////////////////////////////////////////////// 32 // 33 // File: G4OpBoundaryProcess.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 // 2006-11-04 - add capability of calculating the reflectivity 46 // off a metal surface by way of a complex index 47 // of refraction - Thanks to Sehwook Lee and John 48 // Hauptman (Dept. of Physics - Iowa State Univ.) 49 // 2009-11-10 - add capability of simulating surface reflections 50 // with Look-Up-Tables (LUT) containing measured 51 // optical reflectance for a variety of surface 52 // treatments - Thanks to Martin Janecek and 53 // William Moses (Lawrence Berkeley National Lab.) 54 // 2013-06-01 - add the capability of simulating the transmission 55 // of a dichronic filter 56 // 2017-02-24 - add capability of simulating surface reflections 57 // with Look-Up-Tables (LUT) developed in DAVIS 58 // 59 // Author: Peter Gumplinger 60 // adopted from work by Werner Keil - April 2/96 61 // 62 //////////////////////////////////////////////////////////////////////// 63 64 #ifndef G4OpBoundaryProcess_h 65 #define G4OpBoundaryProcess_h 1 66 67 #include "G4OpticalPhoton.hh" 68 #include "G4OpticalSurface.hh" 69 #include "G4RandomTools.hh" 70 #include "G4VDiscreteProcess.hh" 71 72 enum G4OpBoundaryProcessStatus 73 { 74 Undefined, 75 Transmission, 76 FresnelRefraction, 77 FresnelReflection, 78 TotalInternalReflection, 79 LambertianReflection, 80 LobeReflection, 81 SpikeReflection, 82 BackScattering, 83 Absorption, 84 Detection, 85 NotAtBoundary, 86 SameMaterial, 87 StepTooSmall, 88 NoRINDEX, 89 PolishedLumirrorAirReflection, 90 PolishedLumirrorGlueReflection, 91 PolishedAirReflection, 92 PolishedTeflonAirReflection, 93 PolishedTiOAirReflection, 94 PolishedTyvekAirReflection, 95 PolishedVM2000AirReflection, 96 PolishedVM2000GlueReflection, 97 EtchedLumirrorAirReflection, 98 EtchedLumirrorGlueReflection, 99 EtchedAirReflection, 100 EtchedTeflonAirReflection, 101 EtchedTiOAirReflection, 102 EtchedTyvekAirReflection, 103 EtchedVM2000AirReflection, 104 EtchedVM2000GlueReflection, 105 GroundLumirrorAirReflection, 106 GroundLumirrorGlueReflection, 107 GroundAirReflection, 108 GroundTeflonAirReflection, 109 GroundTiOAirReflection, 110 GroundTyvekAirReflection, 111 GroundVM2000AirReflection, 112 GroundVM2000GlueReflection, 113 Dichroic, 114 CoatedDielectricReflection, 115 CoatedDielectricRefraction, 116 CoatedDielectricFrustratedTransmission 117 }; 118 119 class G4OpBoundaryProcess : public G4VDiscreteProcess 120 { 121 public: 122 explicit G4OpBoundaryProcess(const G4String& processName = "OpBoundary", 123 G4ProcessType type = fOptical); 124 virtual ~G4OpBoundaryProcess(); 125 126 virtual G4bool IsApplicable( 127 const G4ParticleDefinition& aParticleType) override; 128 // Returns true -> 'is applicable' only for an optical photon. 129 130 virtual G4double GetMeanFreePath(const G4Track&, G4double, 131 G4ForceCondition* condition) override; 132 // Returns infinity; i. e. the process does not limit the step, but sets the 133 // 'Forced' condition for the DoIt to be invoked at every step. However, only 134 // at a boundary will any action be taken. 135 136 G4VParticleChange* PostStepDoIt(const G4Track& aTrack, 137 const G4Step& aStep) override; 138 // This is the method implementing boundary processes. 139 140 virtual G4OpBoundaryProcessStatus GetStatus() const; 141 // Returns the current status. 142 143 virtual void SetInvokeSD(G4bool); 144 // Set flag for call to InvokeSD method. 145 146 virtual void PreparePhysicsTable(const G4ParticleDefinition&) override; 147 148 virtual void Initialise(); 149 150 void SetVerboseLevel(G4int); 151 152 private: 153 G4OpBoundaryProcess(const G4OpBoundaryProcess& right) = delete; 154 G4OpBoundaryProcess& operator=(const G4OpBoundaryProcess& right) = delete; 155 156 G4bool G4BooleanRand(const G4double prob) const; 157 158 G4ThreeVector GetFacetNormal(const G4ThreeVector& Momentum, 159 const G4ThreeVector& Normal) const; 160 161 void DielectricMetal(); 162 void DielectricDielectric(); 163 164 void DielectricLUT(); 165 void DielectricLUTDAVIS(); 166 167 void DielectricDichroic(); 168 void CoatedDielectricDielectric(); 169 170 void ChooseReflection(); 171 void DoAbsorption(); 172 void DoReflection(); 173 174 G4double GetIncidentAngle(); 175 // Returns the incident angle of optical photon 176 177 G4double GetReflectivity(G4double E1_perp, G4double E1_parl, 178 G4double incidentangle, G4double RealRindex, 179 G4double ImaginaryRindex); 180 // Returns the Reflectivity on a metallic surface 181 182 G4double GetReflectivityThroughThinLayer(G4double sinTL, G4double E1_perp, 183 G4double E1_parl, G4double wavelength, 184 G4double cost1, G4double cost2); 185 // Returns the Reflectivity on a coated surface 186 187 void CalculateReflectivity(); 188 189 void BoundaryProcessVerbose() const; 190 191 // Invoke SD for post step point if the photon is 'detected' 192 G4bool InvokeSD(const G4Step* step); 193 194 G4ThreeVector fOldMomentum; 195 G4ThreeVector fOldPolarization; 196 197 G4ThreeVector fNewMomentum; 198 G4ThreeVector fNewPolarization; 199 200 G4ThreeVector fGlobalNormal; 201 G4ThreeVector fFacetNormal; 202 203 const G4Material* fMaterial1; 204 const G4Material* fMaterial2; 205 206 G4OpticalSurface* fOpticalSurface; 207 208 G4MaterialPropertyVector* fRealRIndexMPV; 209 G4MaterialPropertyVector* fImagRIndexMPV; 210 G4Physics2DVector* fDichroicVector; 211 212 G4double fPhotonMomentum; 213 G4double fRindex1; 214 G4double fRindex2; 215 216 G4double fSint1; 217 218 G4double fReflectivity; 219 G4double fEfficiency; 220 G4double fTransmittance; 221 G4double fSurfaceRoughness; 222 223 G4double fProb_sl, fProb_ss, fProb_bs; 224 G4double fCarTolerance; 225 226 // Used by CoatedDielectricDielectric() 227 G4double fCoatedRindex, fCoatedThickness; 228 229 G4OpBoundaryProcessStatus fStatus; 230 G4OpticalSurfaceModel fModel; 231 G4OpticalSurfaceFinish fFinish; 232 233 G4int f_iTE, f_iTM; 234 235 G4int fNumSmallStepWarnings = 0; // number of times small step warning printed 236 G4int fNumBdryTypeWarnings = 0; // number of times boundary type warning printed 237 238 size_t idx_dichroicX = 0; 239 size_t idx_dichroicY = 0; 240 size_t idx_rindex1 = 0; 241 size_t idx_rindex_surface = 0; 242 size_t idx_reflect = 0; 243 size_t idx_eff = 0; 244 size_t idx_trans = 0; 245 size_t idx_lobe = 0; 246 size_t idx_spike = 0; 247 size_t idx_back = 0; 248 size_t idx_rindex2 = 0; 249 size_t idx_groupvel = 0; 250 size_t idx_rrindex = 0; 251 size_t idx_irindex = 0; 252 size_t idx_coatedrindex = 0; 253 254 // Used by CoatedDielectricDielectric() 255 G4bool fCoatedFrustratedTransmission = true; 256 257 G4bool fInvokeSD; 258 }; 259 260 //////////////////// 261 // Inline methods 262 //////////////////// 263 264 inline G4bool G4OpBoundaryProcess::G4BooleanRand(const G4double prob) const 265 { 266 // Returns a random boolean variable with the specified probability 267 return (G4UniformRand() < prob); 268 } 269 270 inline G4bool G4OpBoundaryProcess::IsApplicable( 271 const G4ParticleDefinition& aParticleType) 272 { 273 return (&aParticleType == G4OpticalPhoton::OpticalPhoton()); 274 } 275 276 inline G4OpBoundaryProcessStatus G4OpBoundaryProcess::GetStatus() const 277 { 278 return fStatus; 279 } 280 281 inline void G4OpBoundaryProcess::ChooseReflection() 282 { 283 G4double rand = G4UniformRand(); 284 if(rand < fProb_ss) 285 { 286 fStatus = SpikeReflection; 287 fFacetNormal = fGlobalNormal; 288 } 289 else if(rand < fProb_ss + fProb_sl) 290 { 291 fStatus = LobeReflection; 292 } 293 else if(rand < fProb_ss + fProb_sl + fProb_bs) 294 { 295 fStatus = BackScattering; 296 } 297 else 298 { 299 fStatus = LambertianReflection; 300 } 301 } 302 303 inline void G4OpBoundaryProcess::DoAbsorption() 304 { 305 fStatus = Absorption; 306 307 if(G4BooleanRand(fEfficiency)) 308 { 309 // EnergyDeposited =/= 0 means: photon has been detected 310 fStatus = Detection; 311 aParticleChange.ProposeLocalEnergyDeposit(fPhotonMomentum); 312 } 313 else 314 { 315 aParticleChange.ProposeLocalEnergyDeposit(0.0); 316 } 317 318 fNewMomentum = fOldMomentum; 319 fNewPolarization = fOldPolarization; 320 321 aParticleChange.ProposeTrackStatus(fStopAndKill); 322 } 323 324 inline void G4OpBoundaryProcess::DoReflection() 325 { 326 if(fStatus == LambertianReflection) 327 { 328 fNewMomentum = G4LambertianRand(fGlobalNormal); 329 fFacetNormal = (fNewMomentum - fOldMomentum).unit(); 330 } 331 else if(fFinish == ground) 332 { 333 fStatus = LobeReflection; 334 if(!fRealRIndexMPV || !fImagRIndexMPV) 335 { 336 fFacetNormal = GetFacetNormal(fOldMomentum, fGlobalNormal); 337 } 338 // else 339 // complex ref. index to be implemented 340 fNewMomentum = 341 fOldMomentum - (2. * fOldMomentum * fFacetNormal * fFacetNormal); 342 } 343 else 344 { 345 fStatus = SpikeReflection; 346 fFacetNormal = fGlobalNormal; 347 fNewMomentum = 348 fOldMomentum - (2. * fOldMomentum * fFacetNormal * fFacetNormal); 349 } 350 fNewPolarization = 351 -fOldPolarization + (2. * fOldPolarization * fFacetNormal * fFacetNormal); 352 } 353 354 #endif /* G4OpBoundaryProcess_h */ 355