Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4OpBoundaryProcess.cc
Go to the documentation of this file.
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// Optical Photon Boundary Process Class Implementation
28////////////////////////////////////////////////////////////////////////
29//
30// File: G4OpBoundaryProcess.cc
31// Description: Discrete Process -- reflection/refraction at
32// optical interfaces
33// Version: 1.1
34// Created: 1997-06-18
35// Modified: 1998-05-25 - Correct parallel component of polarization
36// (thanks to: Stefano Magni + Giovanni Pieri)
37// 1998-05-28 - NULL Rindex pointer before reuse
38// (thanks to: Stefano Magni)
39// 1998-06-11 - delete *sint1 in oblique reflection
40// (thanks to: Giovanni Pieri)
41// 1998-06-19 - move from GetLocalExitNormal() to the new
42// method: GetLocalExitNormal(&valid) to get
43// the surface normal in all cases
44// 1998-11-07 - NULL OpticalSurface pointer before use
45// comparison not sharp for: std::abs(cost1) < 1.0
46// remove sin1, sin2 in lines 556,567
47// (thanks to Stefano Magni)
48// 1999-10-10 - Accommodate changes done in DoAbsorption by
49// changing logic in DielectricMetal
50// 2001-10-18 - avoid Linux (gcc-2.95.2) warning about variables
51// might be used uninitialized in this function
52// moved E2_perp, E2_parl and E2_total out of 'if'
53// 2003-11-27 - Modified line 168-9 to reflect changes made to
54// G4OpticalSurface class ( by Fan Lei)
55// 2004-02-02 - Set theStatus = Undefined at start of DoIt
56// 2005-07-28 - add G4ProcessType to constructor
57// 2006-11-04 - add capability of calculating the reflectivity
58// off a metal surface by way of a complex index
59// of refraction - Thanks to Sehwook Lee and John
60// Hauptman (Dept. of Physics - Iowa State Univ.)
61// 2009-11-10 - add capability of simulating surface reflections
62// with Look-Up-Tables (LUT) containing measured
63// optical reflectance for a variety of surface
64// treatments - Thanks to Martin Janecek and
65// William Moses (Lawrence Berkeley National Lab.)
66// 2013-06-01 - add the capability of simulating the transmission
67// of a dichronic filter
68// 2017-02-24 - add capability of simulating surface reflections
69// with Look-Up-Tables (LUT) developed in DAVIS
70//
71// Author: Peter Gumplinger
72// adopted from work by Werner Keil - April 2/96
73//
74////////////////////////////////////////////////////////////////////////
75
77
78#include "G4ios.hh"
82#include "G4OpProcessSubType.hh"
86#include "G4SystemOfUnits.hh"
89
90namespace {
91 constexpr G4double hc = CLHEP::h_Planck * CLHEP::c_light;
92}
93
94//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96 G4ProcessType ptype)
97 : G4VDiscreteProcess(processName, ptype)
98{
99 Initialise();
100
101 if(verboseLevel > 0)
102 {
103 G4cout << GetProcessName() << " is created " << G4endl;
104 }
106
107 fStatus = Undefined;
108 fModel = glisur;
109 fFinish = polished;
110 fReflectivity = 1.;
111 fEfficiency = 0.;
112 fTransmittance = 0.;
113 fSurfaceRoughness = 0.;
114 fProb_sl = 0.;
115 fProb_ss = 0.;
116 fProb_bs = 0.;
117
118 fRealRIndexMPV = nullptr;
119 fImagRIndexMPV = nullptr;
120 fMaterial1 = nullptr;
121 fMaterial2 = nullptr;
122 fOpticalSurface = nullptr;
124
125 f_iTE = f_iTM = 0;
126 fPhotonMomentum = 0.;
127 fRindex1 = fRindex2 = 1.;
128 fSint1 = 0.;
129 fDichroicVector = nullptr;
130}
131
132//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134
135//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148
149//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
151 const G4Step& aStep)
152{
153 fStatus = Undefined;
154 aParticleChange.Initialize(aTrack);
155 aParticleChange.ProposeVelocity(aTrack.GetVelocity());
156
157 // Get hyperStep from G4ParallelWorldProcess - NOTE: PostSetpDoIt of this
158 // process to be invoked after G4ParallelWorldProcess!
160 const G4Step* pStep = (hStep != nullptr) ? hStep : &aStep;
161
163 {
164 fMaterial1 = pStep->GetPreStepPoint()->GetMaterial();
165 fMaterial2 = pStep->GetPostStepPoint()->GetMaterial();
166 }
167 else
168 {
169 fStatus = NotAtBoundary;
170 if(verboseLevel > 1)
171 BoundaryProcessVerbose();
172 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
173 }
174
177
178 if(verboseLevel > 1)
179 {
180 G4cout << " Photon at Boundary! " << G4endl;
181 if(thePrePV != nullptr)
182 G4cout << " thePrePV: " << thePrePV->GetName() << G4endl;
183 if(thePostPV != nullptr)
184 G4cout << " thePostPV: " << thePostPV->GetName() << G4endl;
185 }
186
187 G4double stepLength = aTrack.GetStepLength();
188 if(stepLength <= fCarTolerance)
189 {
190 fStatus = StepTooSmall;
191 if(verboseLevel > 1)
192 BoundaryProcessVerbose();
193
194 G4MaterialPropertyVector* groupvel = nullptr;
195 G4MaterialPropertiesTable* aMPT = fMaterial2->GetMaterialPropertiesTable();
196 if(aMPT != nullptr)
197 {
198 groupvel = aMPT->GetProperty(kGROUPVEL);
199 }
200
201 if(groupvel != nullptr)
202 {
203 aParticleChange.ProposeVelocity(
204 groupvel->Value(fPhotonMomentum, idx_groupvel));
205 }
206 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
207 }
208 else if (stepLength <= 10. * fCarTolerance && fNumSmallStepWarnings < 10)
209 { // see bug 2510
210 ++fNumSmallStepWarnings;
211 if(verboseLevel > 0)
212 {
214 ed << "G4OpBoundaryProcess: "
215 << "Opticalphoton step length: " << stepLength/mm << " mm." << G4endl
216 << "This is larger than the threshold " << fCarTolerance/mm << " mm "
217 "to set status StepTooSmall." << G4endl
218 << "Boundary scattering may be incorrect. ";
219 if(fNumSmallStepWarnings == 10)
220 {
221 ed << G4endl << "*** Step size warnings stopped.";
222 }
223 G4Exception("G4OpBoundaryProcess", "OpBoun06", JustWarning, ed, "");
224 }
225 }
226
227 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
228
229 fPhotonMomentum = aParticle->GetTotalMomentum();
230 fOldMomentum = aParticle->GetMomentumDirection();
231 fOldPolarization = aParticle->GetPolarization();
232
233 if(verboseLevel > 1)
234 {
235 G4cout << " Old Momentum Direction: " << fOldMomentum << G4endl
236 << " Old Polarization: " << fOldPolarization << G4endl;
237 }
238
239 G4ThreeVector theGlobalPoint = pStep->GetPostStepPoint()->GetPosition();
240 G4bool valid;
241
242 // ID of Navigator which limits step
246 fGlobalNormal = (iNav[hNavId])->GetGlobalExitNormal(theGlobalPoint, &valid);
247
248 if(valid)
249 {
250 fGlobalNormal = -fGlobalNormal;
251 }
252 else
253 {
255 ed << " G4OpBoundaryProcess/PostStepDoIt(): "
256 << " The Navigator reports that it returned an invalid normal" << G4endl;
258 "G4OpBoundaryProcess::PostStepDoIt", "OpBoun01", EventMustBeAborted, ed,
259 "Invalid Surface Normal - Geometry must return a valid surface normal");
260 }
261
262 if(fOldMomentum * fGlobalNormal > 0.0)
263 {
264#ifdef G4OPTICAL_DEBUG
266 ed << " G4OpBoundaryProcess/PostStepDoIt(): fGlobalNormal points in a "
267 "wrong direction. "
268 << G4endl
269 << " The momentum of the photon arriving at interface (oldMomentum)"
270 << " must exit the volume cross in the step. " << G4endl
271 << " So it MUST have dot < 0 with the normal that Exits the new "
272 "volume (globalNormal)."
273 << G4endl << " >> The dot product of oldMomentum and global Normal is "
274 << fOldMomentum * fGlobalNormal << G4endl
275 << " Old Momentum (during step) = " << fOldMomentum << G4endl
276 << " Global Normal (Exiting New Vol) = " << fGlobalNormal << G4endl
277 << G4endl;
278 G4Exception("G4OpBoundaryProcess::PostStepDoIt", "OpBoun02",
279 EventMustBeAborted, // Or JustWarning to see if it happens
280 // repeatedly on one ray
281 ed,
282 "Invalid Surface Normal - Geometry must return valid surface "
283 "normal pointing in the right direction");
284#else
285 fGlobalNormal = -fGlobalNormal;
286#endif
287 }
288
289 G4MaterialPropertyVector* rIndexMPV = nullptr;
290 G4MaterialPropertiesTable* MPT = fMaterial1->GetMaterialPropertiesTable();
291 if(MPT != nullptr)
292 {
293 rIndexMPV = MPT->GetProperty(kRINDEX);
294 }
295 if(rIndexMPV != nullptr)
296 {
297 fRindex1 = rIndexMPV->Value(fPhotonMomentum, idx_rindex1);
298 }
299 else
300 {
301 fStatus = NoRINDEX;
302 if(verboseLevel > 1)
303 BoundaryProcessVerbose();
304 aParticleChange.ProposeLocalEnergyDeposit(fPhotonMomentum);
305 aParticleChange.ProposeTrackStatus(fStopAndKill);
306 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
307 }
308
309 fReflectivity = 1.;
310 fEfficiency = 0.;
311 fTransmittance = 0.;
312 fSurfaceRoughness = 0.;
313 fModel = glisur;
314 fFinish = polished;
316
317 rIndexMPV = nullptr;
318 fOpticalSurface = nullptr;
319
320 G4LogicalSurface* surface =
321 G4LogicalBorderSurface::GetSurface(thePrePV, thePostPV);
322 if(surface == nullptr)
323 {
324 if(thePostPV->GetMotherLogical() == thePrePV->GetLogicalVolume())
325 {
327 if(surface == nullptr)
328 {
329 surface =
331 }
332 }
333 else
334 {
336 if(surface == nullptr)
337 {
338 surface =
340 }
341 }
342 }
343
344 if(surface != nullptr)
345 {
346 fOpticalSurface =
347 dynamic_cast<G4OpticalSurface*>(surface->GetSurfaceProperty());
348 }
349 if(fOpticalSurface != nullptr)
350 {
351 type = fOpticalSurface->GetType();
352 fModel = fOpticalSurface->GetModel();
353 fFinish = fOpticalSurface->GetFinish();
354
356 fOpticalSurface->GetMaterialPropertiesTable();
357 if(sMPT != nullptr)
358 {
359 if(IsBackpainted(fFinish))
360 {
361 rIndexMPV = sMPT->GetProperty(kRINDEX);
362 if(rIndexMPV != nullptr)
363 {
364 fRindex2 = rIndexMPV->Value(fPhotonMomentum, idx_rindex_surface);
365 }
366 else
367 {
368 fStatus = NoRINDEX;
369 if(verboseLevel > 1)
370 BoundaryProcessVerbose();
371 aParticleChange.ProposeLocalEnergyDeposit(fPhotonMomentum);
372 aParticleChange.ProposeTrackStatus(fStopAndKill);
373 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
374 }
375 }
376
377 fRealRIndexMPV = sMPT->GetProperty(kREALRINDEX);
378 fImagRIndexMPV = sMPT->GetProperty(kIMAGINARYRINDEX);
379 f_iTE = f_iTM = 1;
380
382 if((pp = sMPT->GetProperty(kREFLECTIVITY)))
383 {
384 fReflectivity = pp->Value(fPhotonMomentum, idx_reflect);
385 }
386 else if(fRealRIndexMPV && fImagRIndexMPV)
387 {
388 CalculateReflectivity();
389 }
390
391 if((pp = sMPT->GetProperty(kEFFICIENCY)))
392 {
393 fEfficiency = pp->Value(fPhotonMomentum, idx_eff);
394 }
395 if((pp = sMPT->GetProperty(kTRANSMITTANCE)))
396 {
397 fTransmittance = pp->Value(fPhotonMomentum, idx_trans);
398 }
400 {
401 fSurfaceRoughness = sMPT->GetConstProperty(kSURFACEROUGHNESS);
402 }
403
404 if(fModel == unified)
405 {
406 fProb_sl = (pp = sMPT->GetProperty(kSPECULARLOBECONSTANT))
407 ? pp->Value(fPhotonMomentum, idx_lobe)
408 : 0.;
409 fProb_ss = (pp = sMPT->GetProperty(kSPECULARSPIKECONSTANT))
410 ? pp->Value(fPhotonMomentum, idx_spike)
411 : 0.;
412 fProb_bs = (pp = sMPT->GetProperty(kBACKSCATTERCONSTANT))
413 ? pp->Value(fPhotonMomentum, idx_back)
414 : 0.;
415 }
416 } // end of if(sMPT)
417 else if(IsBackpainted(fFinish))
418 {
419 aParticleChange.ProposeLocalEnergyDeposit(fPhotonMomentum);
420 aParticleChange.ProposeTrackStatus(fStopAndKill);
421 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
422 }
423 } // end of if(fOpticalSurface)
424
425 // DIELECTRIC-DIELECTRIC
426 if(type == dielectric_dielectric)
427 {
428 if(fFinish == polished || fFinish == ground)
429 {
430 if(fMaterial1 == fMaterial2)
431 {
432 fStatus = SameMaterial;
433 if(verboseLevel > 1)
434 BoundaryProcessVerbose();
435 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
436 }
437 MPT = fMaterial2->GetMaterialPropertiesTable();
438 rIndexMPV = nullptr;
439 if(MPT != nullptr)
440 {
441 rIndexMPV = MPT->GetProperty(kRINDEX);
442 }
443 if(rIndexMPV != nullptr)
444 {
445 fRindex2 = rIndexMPV->Value(fPhotonMomentum, idx_rindex2);
446 }
447 else
448 {
449 fStatus = NoRINDEX;
450 if(verboseLevel > 1)
451 BoundaryProcessVerbose();
452 aParticleChange.ProposeLocalEnergyDeposit(fPhotonMomentum);
453 aParticleChange.ProposeTrackStatus(fStopAndKill);
454 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
455 }
456 }
457 if(IsBackpainted(fFinish))
458 {
459 DielectricDielectric();
460 }
461 else
462 {
463 G4double rand = G4UniformRand();
464 if(rand > fReflectivity + fTransmittance)
465 {
466 DoAbsorption();
467 }
468 else if(rand > fReflectivity)
469 {
470 DoTransmission();
471 }
472 else
473 {
474 if(fFinish == polishedfrontpainted)
475 {
476 DoReflection();
477 }
478 else if(fFinish == groundfrontpainted)
479 {
480 fStatus = LambertianReflection;
481 DoReflection();
482 }
483 else
484 {
485 DielectricDielectric();
486 }
487 }
488 }
489 }
490 else if(type == dielectric_metal)
491 {
492 DielectricMetal();
493 }
494 else if(type == dielectric_LUT)
495 {
496 DielectricLUT();
497 }
498 else if(type == dielectric_LUTDAVIS)
499 {
500 DielectricLUTDAVIS();
501 }
502 else if(type == dielectric_dichroic)
503 {
504 DielectricDichroic();
505 }
506 else if(type == coated)
507 {
508 CoatedDielectricDielectric();
509 }
510 else
511 {
512 if(fNumBdryTypeWarnings <= 10)
513 {
514 ++fNumBdryTypeWarnings;
515 if(verboseLevel > 0)
516 {
518 ed << " PostStepDoIt(): Illegal boundary type." << G4endl;
519 if(fNumBdryTypeWarnings == 10)
520 {
521 ed << "** Boundary type warnings stopped." << G4endl;
522 }
523 G4Exception("G4OpBoundaryProcess", "OpBoun04", JustWarning, ed);
524 }
525 }
526 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
527 }
528
529 fNewMomentum = fNewMomentum.unit();
530 fNewPolarization = fNewPolarization.unit();
531
532 if(verboseLevel > 1)
533 {
534 G4cout << " New Momentum Direction: " << fNewMomentum << G4endl
535 << " New Polarization: " << fNewPolarization << G4endl;
536 BoundaryProcessVerbose();
537 }
538
539 aParticleChange.ProposeMomentumDirection(fNewMomentum);
540 aParticleChange.ProposePolarization(fNewPolarization);
541
542 if(fStatus == FresnelRefraction || fStatus == Transmission)
543 {
544 // not all surface types check that fMaterial2 has an MPT
545 G4MaterialPropertiesTable* aMPT = fMaterial2->GetMaterialPropertiesTable();
546 G4MaterialPropertyVector* groupvel = nullptr;
547 if(aMPT != nullptr)
548 {
549 groupvel = aMPT->GetProperty(kGROUPVEL);
550 }
551 if(groupvel != nullptr)
552 {
553 aParticleChange.ProposeVelocity(
554 groupvel->Value(fPhotonMomentum, idx_groupvel));
555 }
556 }
557
558 if(fStatus == Detection && fInvokeSD)
559 InvokeSD(pStep);
560 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
561}
562
563//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
564void G4OpBoundaryProcess::BoundaryProcessVerbose() const
565{
566 G4cout << " *** ";
567 if(fStatus == Undefined)
568 G4cout << "Undefined";
569 else if(fStatus == Transmission)
570 G4cout << "Transmission";
571 else if(fStatus == FresnelRefraction)
572 G4cout << "FresnelRefraction";
573 else if(fStatus == FresnelReflection)
574 G4cout << "FresnelReflection";
575 else if(fStatus == TotalInternalReflection)
576 G4cout << "TotalInternalReflection";
577 else if(fStatus == LambertianReflection)
578 G4cout << "LambertianReflection";
579 else if(fStatus == LobeReflection)
580 G4cout << "LobeReflection";
581 else if(fStatus == SpikeReflection)
582 G4cout << "SpikeReflection";
583 else if(fStatus == BackScattering)
584 G4cout << "BackScattering";
585 else if(fStatus == PolishedLumirrorAirReflection)
586 G4cout << "PolishedLumirrorAirReflection";
587 else if(fStatus == PolishedLumirrorGlueReflection)
588 G4cout << "PolishedLumirrorGlueReflection";
589 else if(fStatus == PolishedAirReflection)
590 G4cout << "PolishedAirReflection";
591 else if(fStatus == PolishedTeflonAirReflection)
592 G4cout << "PolishedTeflonAirReflection";
593 else if(fStatus == PolishedTiOAirReflection)
594 G4cout << "PolishedTiOAirReflection";
595 else if(fStatus == PolishedTyvekAirReflection)
596 G4cout << "PolishedTyvekAirReflection";
597 else if(fStatus == PolishedVM2000AirReflection)
598 G4cout << "PolishedVM2000AirReflection";
599 else if(fStatus == PolishedVM2000GlueReflection)
600 G4cout << "PolishedVM2000GlueReflection";
601 else if(fStatus == EtchedLumirrorAirReflection)
602 G4cout << "EtchedLumirrorAirReflection";
603 else if(fStatus == EtchedLumirrorGlueReflection)
604 G4cout << "EtchedLumirrorGlueReflection";
605 else if(fStatus == EtchedAirReflection)
606 G4cout << "EtchedAirReflection";
607 else if(fStatus == EtchedTeflonAirReflection)
608 G4cout << "EtchedTeflonAirReflection";
609 else if(fStatus == EtchedTiOAirReflection)
610 G4cout << "EtchedTiOAirReflection";
611 else if(fStatus == EtchedTyvekAirReflection)
612 G4cout << "EtchedTyvekAirReflection";
613 else if(fStatus == EtchedVM2000AirReflection)
614 G4cout << "EtchedVM2000AirReflection";
615 else if(fStatus == EtchedVM2000GlueReflection)
616 G4cout << "EtchedVM2000GlueReflection";
617 else if(fStatus == GroundLumirrorAirReflection)
618 G4cout << "GroundLumirrorAirReflection";
619 else if(fStatus == GroundLumirrorGlueReflection)
620 G4cout << "GroundLumirrorGlueReflection";
621 else if(fStatus == GroundAirReflection)
622 G4cout << "GroundAirReflection";
623 else if(fStatus == GroundTeflonAirReflection)
624 G4cout << "GroundTeflonAirReflection";
625 else if(fStatus == GroundTiOAirReflection)
626 G4cout << "GroundTiOAirReflection";
627 else if(fStatus == GroundTyvekAirReflection)
628 G4cout << "GroundTyvekAirReflection";
629 else if(fStatus == GroundVM2000AirReflection)
630 G4cout << "GroundVM2000AirReflection";
631 else if(fStatus == GroundVM2000GlueReflection)
632 G4cout << "GroundVM2000GlueReflection";
633 else if(fStatus == Absorption)
634 G4cout << "Absorption";
635 else if(fStatus == Detection)
636 G4cout << "Detection";
637 else if(fStatus == NotAtBoundary)
638 G4cout << "NotAtBoundary";
639 else if(fStatus == SameMaterial)
640 G4cout << "SameMaterial";
641 else if(fStatus == StepTooSmall)
642 G4cout << "StepTooSmall";
643 else if(fStatus == NoRINDEX)
644 G4cout << "NoRINDEX";
645 else if(fStatus == Dichroic)
646 G4cout << "Dichroic Transmission";
647 else if(fStatus == CoatedDielectricReflection)
648 G4cout << "Coated Dielectric Reflection";
649 else if(fStatus == CoatedDielectricRefraction)
650 G4cout << "Coated Dielectric Refraction";
652 G4cout << "Coated Dielectric Frustrated Transmission";
653
654 G4cout << " ***" << G4endl;
655}
656
657//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
658G4ThreeVector G4OpBoundaryProcess::GetFacetNormal(
659 const G4ThreeVector& momentum, const G4ThreeVector& normal) const
660{
661 G4ThreeVector facetNormal;
662 if(fModel == unified || fModel == LUT || fModel == DAVIS)
663 {
664 /* This function codes alpha to a random value taken from the
665 distribution p(alpha) = g(alpha; 0, sigma_alpha)*std::sin(alpha),
666 for alpha > 0 and alpha < 90, where g(alpha; 0, sigma_alpha) is a
667 gaussian distribution with mean 0 and standard deviation sigma_alpha. */
668
669 G4double sigma_alpha = 0.0;
670 if(fOpticalSurface)
671 sigma_alpha = fOpticalSurface->GetSigmaAlpha();
672 if(sigma_alpha == 0.0)
673 {
674 return normal;
675 }
676
677 G4double f_max = std::min(1.0, 4. * sigma_alpha);
678 G4double alpha, phi, sinAlpha;
679
680 do
681 { // Loop checking, 13-Aug-2015, Peter Gumplinger
682 do
683 { // Loop checking, 13-Aug-2015, Peter Gumplinger
684 alpha = G4RandGauss::shoot(0.0, sigma_alpha);
685 sinAlpha = std::sin(alpha);
686 } while(G4UniformRand() * f_max > sinAlpha || alpha >= halfpi);
687
688 phi = G4UniformRand() * twopi;
689 facetNormal.set(sinAlpha * std::cos(phi), sinAlpha * std::sin(phi),
690 std::cos(alpha));
691 facetNormal.rotateUz(normal);
692 } while(momentum * facetNormal >= 0.0);
693 }
694 else
695 {
696 G4double polish = 1.0;
697 if(fOpticalSurface)
698 polish = fOpticalSurface->GetPolish();
699 if(polish < 1.0)
700 {
701 do
702 { // Loop checking, 13-Aug-2015, Peter Gumplinger
703 G4ThreeVector smear;
704 do
705 { // Loop checking, 13-Aug-2015, Peter Gumplinger
706 smear.setX(2. * G4UniformRand() - 1.);
707 smear.setY(2. * G4UniformRand() - 1.);
708 smear.setZ(2. * G4UniformRand() - 1.);
709 } while(smear.mag2() > 1.0);
710 facetNormal = normal + (1. - polish) * smear;
711 } while(momentum * facetNormal >= 0.0);
712 facetNormal = facetNormal.unit();
713 }
714 else
715 {
716 facetNormal = normal;
717 }
718 }
719 return facetNormal;
720}
721
722//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
723void G4OpBoundaryProcess::DielectricMetal()
724{
725 G4int n = 0;
726 G4double rand;
727 G4ThreeVector A_trans;
728
729 do
730 {
731 ++n;
732 rand = G4UniformRand();
733 if(rand > fReflectivity && n == 1)
734 {
735 if(rand > fReflectivity + fTransmittance)
736 {
737 DoAbsorption();
738 }
739 else
740 {
741 DoTransmission();
742 }
743 break;
744 }
745 else
746 {
747 if(fRealRIndexMPV && fImagRIndexMPV)
748 {
749 if(n > 1)
750 {
751 CalculateReflectivity();
752 if(!G4BooleanRand(fReflectivity))
753 {
754 DoAbsorption();
755 break;
756 }
757 }
758 }
759 if(fModel == glisur || fFinish == polished)
760 {
761 DoReflection();
762 }
763 else
764 {
765 if(n == 1)
766 ChooseReflection();
767 if(fStatus == LambertianReflection)
768 {
769 DoReflection();
770 }
771 else if(fStatus == BackScattering)
772 {
773 fNewMomentum = -fOldMomentum;
774 fNewPolarization = -fOldPolarization;
775 }
776 else
777 {
778 if(fStatus == LobeReflection)
779 {
780 if(!fRealRIndexMPV || !fImagRIndexMPV)
781 {
782 fFacetNormal = GetFacetNormal(fOldMomentum, fGlobalNormal);
783 }
784 // else
785 // case of complex rindex needs to be implemented
786 }
787 fNewMomentum = fOldMomentum -
788 2. * (fOldMomentum * fFacetNormal) * fFacetNormal;
789
790 if(f_iTE > 0 && f_iTM > 0)
791 {
792 fNewPolarization = -fOldPolarization +
793 2. * (fOldPolarization * fFacetNormal) * fFacetNormal;
794 }
795 else if(f_iTE > 0)
796 {
797 A_trans = (fSint1 > 0.0) ? fOldMomentum.cross(fFacetNormal).unit()
798 : fOldPolarization;
799 fNewPolarization = -A_trans;
800 }
801 else if(f_iTM > 0)
802 {
803 fNewPolarization =
804 -fNewMomentum.cross(A_trans).unit(); // = -A_paral
805 }
806 }
807 }
808 fOldMomentum = fNewMomentum;
809 fOldPolarization = fNewPolarization;
810 }
811 // Loop checking, 13-Aug-2015, Peter Gumplinger
812 } while(fNewMomentum * fGlobalNormal < 0.0);
813}
814
815//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
816void G4OpBoundaryProcess::DielectricLUT()
817{
818 G4int thetaIndex, phiIndex;
819 G4double angularDistVal, thetaRad, phiRad;
820 G4ThreeVector perpVectorTheta, perpVectorPhi;
821
823 G4int(fFinish) + (G4int(NoRINDEX) - G4int(groundbackpainted)));
824
825 G4int thetaIndexMax = fOpticalSurface->GetThetaIndexMax();
826 G4int phiIndexMax = fOpticalSurface->GetPhiIndexMax();
827
828 G4double rand;
829
830 do
831 {
832 rand = G4UniformRand();
833 if(rand > fReflectivity)
834 {
835 if(rand > fReflectivity + fTransmittance)
836 {
837 DoAbsorption();
838 }
839 else
840 {
841 DoTransmission();
842 }
843 break;
844 }
845 else
846 {
847 // Calculate Angle between Normal and Photon Momentum
848 G4double anglePhotonToNormal = fOldMomentum.angle(-fGlobalNormal);
849 // Round to closest integer: LBNL model array has 91 values
850 G4int angleIncident = (G4int)std::lrint(anglePhotonToNormal / CLHEP::deg);
851
852 // Take random angles THETA and PHI,
853 // and see if below Probability - if not - Redo
854 do
855 {
856 thetaIndex = (G4int)G4RandFlat::shootInt(thetaIndexMax - 1);
857 phiIndex = (G4int)G4RandFlat::shootInt(phiIndexMax - 1);
858 // Find probability with the new indeces from LUT
859 angularDistVal = fOpticalSurface->GetAngularDistributionValue(
860 angleIncident, thetaIndex, phiIndex);
861 // Loop checking, 13-Aug-2015, Peter Gumplinger
862 } while(!G4BooleanRand(angularDistVal));
863
864 thetaRad = G4double(-90 + 4 * thetaIndex) * pi / 180.;
865 phiRad = G4double(-90 + 5 * phiIndex) * pi / 180.;
866 // Rotate Photon Momentum in Theta, then in Phi
867 fNewMomentum = -fOldMomentum;
868
869 perpVectorTheta = fNewMomentum.cross(fGlobalNormal);
870 if(perpVectorTheta.mag() < fCarTolerance)
871 {
872 perpVectorTheta = fNewMomentum.orthogonal();
873 }
874 fNewMomentum =
875 fNewMomentum.rotate(anglePhotonToNormal - thetaRad, perpVectorTheta);
876 perpVectorPhi = perpVectorTheta.cross(fNewMomentum);
877 fNewMomentum = fNewMomentum.rotate(-phiRad, perpVectorPhi);
878
879 // Rotate Polarization too:
880 fFacetNormal = (fNewMomentum - fOldMomentum).unit();
881 fNewPolarization = -fOldPolarization +
882 2. * (fOldPolarization * fFacetNormal) * fFacetNormal;
883 }
884 // Loop checking, 13-Aug-2015, Peter Gumplinger
885 } while(fNewMomentum * fGlobalNormal <= 0.0);
886}
887
888//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
889void G4OpBoundaryProcess::DielectricLUTDAVIS()
890{
891 G4int angindex, random, angleIncident;
892 G4double reflectivityValue, elevation, azimuth;
893 G4double anglePhotonToNormal;
894
895 G4int lutbin = fOpticalSurface->GetLUTbins();
896 G4double rand = G4UniformRand();
897
898 G4double sinEl;
899 G4ThreeVector u, vNorm, w;
900
901 do
902 {
903 anglePhotonToNormal = fOldMomentum.angle(-fGlobalNormal);
904
905 // Davis model has 90 reflection bins: round down
906 // don't allow angleIncident to be 90 for anglePhotonToNormal close to 90
907 angleIncident = std::min(
908 static_cast<G4int>(std::floor(anglePhotonToNormal / CLHEP::deg)), 89);
909 reflectivityValue = fOpticalSurface->GetReflectivityLUTValue(angleIncident);
910
911 if(rand > reflectivityValue)
912 {
913 if(fEfficiency > 0.)
914 {
915 DoAbsorption();
916 break;
917 }
918 else
919 {
920 fStatus = Transmission;
921
922 if(angleIncident <= 0.01)
923 {
924 fNewMomentum = fOldMomentum;
925 break;
926 }
927
928 do
929 {
930 random = (G4int)G4RandFlat::shootInt(1, lutbin + 1);
931 angindex =
932 (((random * 2) - 1)) + angleIncident * lutbin * 2 + 3640000;
933
934 azimuth =
935 fOpticalSurface->GetAngularDistributionValueLUT(angindex - 1);
936 elevation = fOpticalSurface->GetAngularDistributionValueLUT(angindex);
937 } while(elevation == 0. && azimuth == 0.);
938
939 sinEl = std::sin(elevation);
940 vNorm = (fGlobalNormal.cross(fOldMomentum)).unit();
941 u = vNorm.cross(fGlobalNormal) * (sinEl * std::cos(azimuth));
942 vNorm *= (sinEl * std::sin(azimuth));
943 // fGlobalNormal shouldn't be modified here
944 w = (fGlobalNormal *= std::cos(elevation));
945 fNewMomentum = u + vNorm + w;
946
947 // Rotate Polarization too:
948 fFacetNormal = (fNewMomentum - fOldMomentum).unit();
949 fNewPolarization = -fOldPolarization +
950 2. * (fOldPolarization * fFacetNormal) * fFacetNormal;
951 }
952 }
953 else
954 {
955 fStatus = LobeReflection;
956
957 if(angleIncident == 0)
958 {
959 fNewMomentum = -fOldMomentum;
960 break;
961 }
962
963 do
964 {
965 random = (G4int)G4RandFlat::shootInt(1, lutbin + 1);
966 angindex = (((random * 2) - 1)) + (angleIncident - 1) * lutbin * 2;
967
968 azimuth = fOpticalSurface->GetAngularDistributionValueLUT(angindex - 1);
969 elevation = fOpticalSurface->GetAngularDistributionValueLUT(angindex);
970 } while(elevation == 0. && azimuth == 0.);
971
972 sinEl = std::sin(elevation);
973 vNorm = (fGlobalNormal.cross(fOldMomentum)).unit();
974 u = vNorm.cross(fGlobalNormal) * (sinEl * std::cos(azimuth));
975 vNorm *= (sinEl * std::sin(azimuth));
976 // fGlobalNormal shouldn't be modified here
977 w = (fGlobalNormal *= std::cos(elevation));
978
979 fNewMomentum = u + vNorm + w;
980
981 // Rotate Polarization too: (needs revision)
982 fNewPolarization = fOldPolarization;
983 }
984 } while(fNewMomentum * fGlobalNormal <= 0.0);
985}
986
987//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
988void G4OpBoundaryProcess::DielectricDichroic()
989{
990 // Calculate Angle between Normal and Photon Momentum
991 G4double anglePhotonToNormal = fOldMomentum.angle(-fGlobalNormal);
992
993 // Round it to closest integer
994 G4double angleIncident = std::floor(180. / pi * anglePhotonToNormal + 0.5);
995
996 if(!fDichroicVector && fOpticalSurface)
997 {
998 fDichroicVector = fOpticalSurface->GetDichroicVector();
999 }
1000
1001 if(fDichroicVector)
1002 {
1003 G4double wavelength = hc / fPhotonMomentum;
1004 fTransmittance = fDichroicVector->Value(wavelength / nm, angleIncident,
1005 idx_dichroicX, idx_dichroicY) *
1006 perCent;
1007 // G4cout << "wavelength: " << std::floor(wavelength/nm)
1008 // << "nm" << G4endl;
1009 // G4cout << "Incident angle: " << angleIncident << "deg" << G4endl;
1010 // G4cout << "Transmittance: "
1011 // << std::floor(fTransmittance/perCent) << "%" << G4endl;
1012 }
1013 else
1014 {
1016 ed << " G4OpBoundaryProcess/DielectricDichroic(): "
1017 << " The dichroic surface has no G4Physics2DVector" << G4endl;
1018 G4Exception("G4OpBoundaryProcess::DielectricDichroic", "OpBoun03",
1019 FatalException, ed,
1020 "A dichroic surface must have an associated G4Physics2DVector");
1021 }
1022
1023 if(!G4BooleanRand(fTransmittance))
1024 { // Not transmitted, so reflect
1025 if(fModel == glisur || fFinish == polished)
1026 {
1027 DoReflection();
1028 }
1029 else
1030 {
1031 ChooseReflection();
1032 if(fStatus == LambertianReflection)
1033 {
1034 DoReflection();
1035 }
1036 else if(fStatus == BackScattering)
1037 {
1038 fNewMomentum = -fOldMomentum;
1039 fNewPolarization = -fOldPolarization;
1040 }
1041 else
1042 {
1043 G4double PdotN, EdotN;
1044 do
1045 {
1046 if(fStatus == LobeReflection)
1047 {
1048 fFacetNormal = GetFacetNormal(fOldMomentum, fGlobalNormal);
1049 }
1050 PdotN = fOldMomentum * fFacetNormal;
1051 fNewMomentum = fOldMomentum - (2. * PdotN) * fFacetNormal;
1052 // Loop checking, 13-Aug-2015, Peter Gumplinger
1053 } while(fNewMomentum * fGlobalNormal <= 0.0);
1054
1055 EdotN = fOldPolarization * fFacetNormal;
1056 fNewPolarization = -fOldPolarization + (2. * EdotN) * fFacetNormal;
1057 }
1058 }
1059 }
1060 else
1061 {
1062 fStatus = Dichroic;
1063 fNewMomentum = fOldMomentum;
1064 fNewPolarization = fOldPolarization;
1065 }
1066}
1067
1068//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1069void G4OpBoundaryProcess::DielectricDielectric()
1070{
1071 fFacetNormal = (fFinish == polished) ? fGlobalNormal :
1072 GetFacetNormal(fOldMomentum, fGlobalNormal);
1073
1074 G4double cost1 = -fOldMomentum * fFacetNormal;
1075
1076 G4bool surfaceRoughnessCriterionPass = true;
1077 if(fSurfaceRoughness != 0. && fRindex1 > fRindex2)
1078 {
1079 G4double wavelength = hc / fPhotonMomentum;
1080 G4double surfaceRoughnessCriterion = std::exp(-std::pow(
1081 (4. * pi * fSurfaceRoughness * fRindex1 * cost1 / wavelength), 2));
1082 surfaceRoughnessCriterionPass = G4BooleanRand(surfaceRoughnessCriterion);
1083 }
1084
1085 G4bool inside = false; // flag whether the photon passes through the boundary
1086 G4bool swap = false; // flag for swapping of materials and r-indices
1087
1088 const G4int maxBoundaryActions = 100;
1089 G4int numBoundaryActions = 0;
1090
1091 do
1092 {
1093 // Apply the boundary action at a dielectric-dielectric interface
1094 ApplyDielectricBoundaryTransition(cost1, surfaceRoughnessCriterionPass,
1095 inside, swap);
1096
1097 // For backpainted surfaces, allow multiple internal bounces if refracted
1098 if((inside && !swap) && IsBackpainted(fFinish))
1099 {
1100 G4double rand = G4UniformRand();
1101 if(rand > fReflectivity)
1102 {
1103 if(rand > fReflectivity + fTransmittance)
1104 {
1105 DoAbsorption();
1106 }
1107 else
1108 {
1109 DoTransmission();
1110 }
1111 break;
1112 }
1113 else
1114 {
1115 if(fStatus != FresnelRefraction)
1116 {
1117 fGlobalNormal = -fGlobalNormal;
1118 }
1119 else
1120 {
1121 swap = !swap;
1122 G4SwapPtr(fMaterial1, fMaterial2);
1123 G4SwapObj(&fRindex1, &fRindex2);
1124 }
1125 if(fFinish == groundbackpainted)
1126 fStatus = LambertianReflection;
1127
1128 DoReflection();
1129
1130 fGlobalNormal = -fGlobalNormal;
1131 fOldMomentum = fNewMomentum;
1132
1133 ++numBoundaryActions;
1134 }
1135 }
1136 else
1137 {
1138 // No more action, stop
1139 break;
1140 }
1141 } while(numBoundaryActions < maxBoundaryActions);
1142
1143 if(numBoundaryActions == maxBoundaryActions)
1144 {
1146 ed << " DielectricDielectric(): numBoundaryActions exceeds"
1147 << " maxBoundaryActions = " << maxBoundaryActions << G4endl;
1148 G4Exception("G4OpBoundaryProcess", "OpBoun05", JustWarning, ed);
1149 }
1150}
1151
1152//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1153void G4OpBoundaryProcess::ApplyDielectricBoundaryTransition(G4double cost1,
1154 G4bool roughnessCriterionPass,
1155 G4bool& inside, G4bool& swap)
1156{
1157 // Apply boundary actions at a dielectric–dielectric interface using the
1158 // incident angle (cost1) and surface roughness criterion
1159
1160 G4bool through = false; // Flag indicating whether the photon is refracted
1161 G4bool done = false;
1162
1163 G4ThreeVector A_trans, A_paral, E1pp, E1pl;
1164 G4double E1_perp, E1_parl;
1165 G4double s1, s2, E2_perp, E2_parl, E2_total;
1166 G4double E2_abs;
1168
1169 G4double cost2 = 0.;
1170 G4double sint2 = 0.;
1171
1172 do
1173 {
1174 if(through)
1175 {
1176 swap = !swap;
1177 through = false;
1178 fGlobalNormal = -fGlobalNormal;
1179 G4SwapPtr(fMaterial1, fMaterial2);
1180 G4SwapObj(&fRindex1, &fRindex2);
1181 }
1182
1183 fFacetNormal = (fFinish == polished) ? fGlobalNormal :
1184 GetFacetNormal(fOldMomentum, fGlobalNormal);
1185
1186 cost1 = -fOldMomentum * fFacetNormal;
1187 if(std::abs(cost1) < 1.0 - fCarTolerance)
1188 {
1189 fSint1 = std::sqrt(1. - cost1 * cost1);
1190 sint2 = fSint1 * fRindex1 / fRindex2; // *** Snell's Law ***
1191 // this isn't a sine as we might expect from the name; can be > 1
1192 }
1193 else
1194 {
1195 fSint1 = 0.0;
1196 sint2 = 0.0;
1197 }
1198
1199 // TOTAL INTERNAL REFLECTION
1200 if(sint2 >= 1.0)
1201 {
1202 swap = false;
1203
1204 fStatus = TotalInternalReflection;
1205 if(!roughnessCriterionPass)
1206 fStatus = LambertianReflection;
1207 if(fModel == unified && fFinish != polished)
1208 ChooseReflection();
1209 if(fStatus == LambertianReflection)
1210 {
1211 DoReflection();
1212 }
1213 else if(fStatus == BackScattering)
1214 {
1215 fNewMomentum = -fOldMomentum;
1216 fNewPolarization = -fOldPolarization;
1217 }
1218 else
1219 {
1220 fNewMomentum = fOldMomentum -
1221 2. * (fOldMomentum * fFacetNormal) * fFacetNormal;
1222 fNewPolarization = -fOldPolarization +
1223 2. * (fOldPolarization * fFacetNormal) * fFacetNormal;
1224 }
1225 }
1226 // NOT TIR
1227 else // if(sint2 < 1.0)
1228 {
1229 // Calculate amplitude for transmission (Q = P x N)
1230 cost2 = (cost1 > 0.0 ? 1.0 : -1.0) * std::sqrt(1. - sint2 * sint2);
1231
1232 if(fSint1 > 0.0)
1233 {
1234 A_trans = (fOldMomentum.cross(fFacetNormal)).unit();
1235 E1_perp = fOldPolarization * A_trans;
1236 E1pp = E1_perp * A_trans;
1237 E1pl = fOldPolarization - E1pp;
1238 E1_parl = E1pl.mag();
1239 }
1240 else
1241 {
1242 A_trans = fOldPolarization;
1243 // Here we Follow Jackson's conventions and set the parallel
1244 // component = 1 in case of a ray perpendicular to the surface
1245 E1_perp = 0.0;
1246 E1_parl = 1.0;
1247 }
1248
1249 s1 = fRindex1 * cost1;
1250 E2_perp = 2. * s1 * E1_perp / (fRindex1 * cost1 + fRindex2 * cost2);
1251 E2_parl = 2. * s1 * E1_parl / (fRindex2 * cost1 + fRindex1 * cost2);
1252 E2_total = E2_perp * E2_perp + E2_parl * E2_parl;
1253 s2 = fRindex2 * cost2 * E2_total;
1254
1255 // D.Sawkey, 24 May 24
1256 // Transmittance has already been taken into account in PostStepDoIt.
1257 // For e.g. specular surfaces, the ratio of Fresnel refraction to
1258 // reflection should be given by the math, not material property
1259 // TRANSMITTANCE
1260 //if(fTransmittance > 0.)
1261 // transCoeff = fTransmittance;
1262 //else if(cost1 != 0.0)
1263 G4double transCoeff = (cost1 != 0.0) ? s2 / s1 : 0.0;
1264
1265 // NOT TIR: REFLECTION
1266 if(!G4BooleanRand(transCoeff))
1267 {
1268 swap = false;
1269 fStatus = FresnelReflection;
1270
1271 if(!roughnessCriterionPass)
1272 fStatus = LambertianReflection;
1273 if(fModel == unified && fFinish != polished)
1274 ChooseReflection();
1275 if(fStatus == LambertianReflection)
1276 {
1277 DoReflection();
1278 }
1279 else if(fStatus == BackScattering)
1280 {
1281 fNewMomentum = -fOldMomentum;
1282 fNewPolarization = -fOldPolarization;
1283 }
1284 else
1285 {
1286 fNewMomentum =
1287 fOldMomentum - 2. * fOldMomentum * fFacetNormal * fFacetNormal;
1288 if(fSint1 > 0.0)
1289 { // incident ray oblique
1290 E2_parl = fRindex2 * E2_parl / fRindex1 - E1_parl;
1291 E2_perp = E2_perp - E1_perp;
1292 E2_total = E2_perp * E2_perp + E2_parl * E2_parl;
1293 A_paral = (fNewMomentum.cross(A_trans)).unit();
1294 E2_abs = std::sqrt(E2_total);
1295 fNewPolarization = (E2_parl * A_paral + E2_perp * A_trans) / E2_abs;
1296 }
1297 else
1298 { // incident ray perpendicular
1299 fNewPolarization =
1300 (fRindex2 > fRindex1 ? -1.0 : 1.0) * fOldPolarization;
1301 }
1302 }
1303 }
1304 // NOT TIR: TRANSMISSION
1305 else
1306 {
1307 inside = !inside;
1308 through = true;
1309 fStatus = FresnelRefraction;
1310
1311 if(fSint1 > 0.0)
1312 { // incident ray oblique
1313 alpha = cost1 - cost2 * (fRindex2 / fRindex1);
1314 fNewMomentum = (fOldMomentum + alpha * fFacetNormal).unit();
1315 A_paral = (fNewMomentum.cross(A_trans)).unit();
1316 E2_abs = std::sqrt(E2_total);
1317
1318 fNewPolarization = (E2_parl * A_paral + E2_perp * A_trans) / E2_abs;
1319 }
1320 else
1321 { // incident ray perpendicular
1322 fNewMomentum = fOldMomentum;
1323 fNewPolarization = fOldPolarization;
1324 }
1325 }
1326 }
1327
1328 fOldMomentum = fNewMomentum.unit();
1329 fOldPolarization = fNewPolarization.unit();
1330
1331 if(fStatus == FresnelRefraction)
1332 {
1333 done = (fNewMomentum * fGlobalNormal <= 0.0);
1334 }
1335 else
1336 {
1337 done = (fNewMomentum * fGlobalNormal >= -fCarTolerance);
1338 }
1339 // Loop checking, 13-Aug-2015, Peter Gumplinger
1340 } while(!done);
1341}
1342
1343//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1350
1351//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1352G4double G4OpBoundaryProcess::GetIncidentAngle()
1353{
1354 return pi - std::acos(fOldMomentum * fFacetNormal /
1355 (fOldMomentum.mag() * fFacetNormal.mag()));
1356}
1357
1358//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1359G4double G4OpBoundaryProcess::GetReflectivity(G4double E1_perp,
1360 G4double E1_parl,
1361 G4double incidentangle,
1362 G4double realRindex,
1363 G4double imaginaryRindex)
1364{
1365 G4complex N1(fRindex1, 0.);
1366 G4complex N2(realRindex, imaginaryRindex);
1367
1368 G4MaterialPropertiesTable* MPT = fMaterial1->GetMaterialPropertiesTable();
1371 if(ppR && ppI)
1372 {
1373 G4double rRindex = ppR->Value(fPhotonMomentum, idx_rrindex);
1374 G4double iRindex = ppI->Value(fPhotonMomentum, idx_irindex);
1375 N1 = G4complex(rRindex, iRindex);
1376 }
1377
1378 // Following two equations, rTM and rTE, are from: "Introduction To Modern
1379 // Optics" written by Fowles
1380 G4double sint = std::sin(incidentangle);
1381 G4double cost = std::cos(incidentangle);
1382 G4complex u(1., 0.); // unit number 1
1383 G4complex cosPhi = std::sqrt(u - ((sint * sint) * (N1 * N1) / (N2 * N2)));
1384
1385 // TE polarization <- E1_perp=1 E1_parl=0
1386 // TM polarization <- E1_parl=1 E1_perp=0
1387 G4complex rTE = (N1 * cost - N2 * cosPhi) / (N1 * cost + N2 * cosPhi);
1388 G4complex rTM = (N2 * cost - N1 * cosPhi) / (N2 * cost + N1 * cosPhi);
1389
1390 // This is my (PG) calculaton for reflectivity on a metallic surface
1391 // depending on the fraction of TE and TM polarization
1392 // when TE polarization, E1_parl=0 and E1_perp=1, R=abs(rTE)^2 and
1393 // when TM polarization, E1_parl=1 and E1_perp=0, R=abs(rTM)^2
1394 G4double E1_perp_sq = E1_perp * E1_perp;
1395 G4double E1_parl_sq = E1_parl * E1_parl;
1396 G4double E1_sq = E1_perp_sq + E1_parl_sq;
1397
1398 G4complex reflectivity_TE = (rTE * conj(rTE)) * E1_perp_sq / E1_sq;
1399 G4complex reflectivity_TM = (rTM * conj(rTM)) * E1_parl_sq / E1_sq;
1400 G4complex reflectivity = reflectivity_TE + reflectivity_TM;
1401
1402 // Real component of the reflectivity
1403 G4double rr = real(reflectivity);
1404 do
1405 {
1406 f_iTE = (G4UniformRand() * rr > real(reflectivity_TE)) ? -1 : 1;
1407 f_iTM = (G4UniformRand() * rr > real(reflectivity_TM)) ? -1 : 1;
1408 // Loop checking, 13-Aug-2015, Peter Gumplinger
1409 } while(f_iTE < 0 && f_iTM < 0);
1410
1411 return rr;
1412}
1413
1414//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1415void G4OpBoundaryProcess::CalculateReflectivity()
1416{
1417 G4double realRindex = fRealRIndexMPV->Value(fPhotonMomentum, idx_rrindex);
1418 G4double imagRindex = fImagRIndexMPV->Value(fPhotonMomentum, idx_irindex);
1419
1420 // calculate FacetNormal
1421 fFacetNormal = (fFinish == ground) ?
1422 GetFacetNormal(fOldMomentum, fGlobalNormal) : fGlobalNormal;
1423
1424 G4double cost1 = -fOldMomentum * fFacetNormal;
1425 fSint1 = (std::abs(cost1) < 1.0 - fCarTolerance) ?
1426 std::sqrt(1. - cost1 * cost1) : 0.0;
1427
1428 G4ThreeVector A_trans, A_paral, E1pp, E1pl;
1429 G4double E1_perp, E1_parl;
1430
1431 if(fSint1 > 0.0)
1432 {
1433 A_trans = (fOldMomentum.cross(fFacetNormal)).unit();
1434 E1_perp = fOldPolarization * A_trans;
1435 E1pp = E1_perp * A_trans;
1436 E1pl = fOldPolarization - E1pp;
1437 E1_parl = E1pl.mag();
1438 }
1439 else
1440 {
1441 A_trans = fOldPolarization;
1442 // Here we Follow Jackson's conventions and we set the parallel
1443 // component = 1 in case of a ray perpendicular to the surface
1444 E1_perp = 0.0;
1445 E1_parl = 1.0;
1446 }
1447
1448 G4double incidentangle = GetIncidentAngle();
1449
1450 // calculate the reflectivity depending on incident angle,
1451 // polarization and complex refractive
1452 fReflectivity = GetReflectivity(E1_perp, E1_parl, incidentangle, realRindex,
1453 imagRindex);
1454}
1455
1456//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1457G4bool G4OpBoundaryProcess::InvokeSD(const G4Step* pStep)
1458{
1459 G4Step aStep = *pStep;
1460 aStep.AddTotalEnergyDeposit(fPhotonMomentum);
1461
1462 G4VSensitiveDetector* sd = aStep.GetPostStepPoint()->GetSensitiveDetector();
1463 if(sd != nullptr)
1464 return sd->Hit(&aStep);
1465 else
1466 return false;
1467}
1468
1469//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1471{
1472 fInvokeSD = flag;
1474}
1475
1476//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1482
1483//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1484void G4OpBoundaryProcess::CoatedDielectricDielectric()
1485{
1486 G4MaterialPropertyVector* pp = nullptr;
1487
1489 if((pp = MPT->GetProperty(kRINDEX)))
1490 {
1491 fRindex2 = pp->Value(fPhotonMomentum, idx_rindex2);
1492 }
1493
1494 MPT = fOpticalSurface->GetMaterialPropertiesTable();
1495 if((pp = MPT->GetProperty(kCOATEDRINDEX)))
1496 {
1497 fCoatedRindex = pp->Value(fPhotonMomentum, idx_coatedrindex);
1498 }
1500 {
1501 fCoatedThickness = MPT->GetConstProperty(kCOATEDTHICKNESS);
1502 }
1504 {
1505 fCoatedFrustratedTransmission =
1507 }
1508
1509 G4double sintTL;
1510 G4double wavelength = hc / fPhotonMomentum;
1511 G4double PdotN;
1512 G4double E1_perp, E1_parl;
1513 G4double s1, E2_perp, E2_parl, E2_total;
1514 G4double E2_abs;
1516 G4ThreeVector A_trans, A_paral, E1pp, E1pl;
1517 //G4bool Inside = false;
1518 //G4bool Swap = false;
1519 G4bool through = false;
1520 G4bool done = false;
1521
1522 do
1523 {
1524 if (through)
1525 {
1526 //Swap = !Swap;
1527 through = false;
1528 fGlobalNormal = -fGlobalNormal;
1529 G4SwapPtr(fMaterial1, fMaterial2);
1530 G4SwapObj(&fRindex1, &fRindex2);
1531 }
1532
1533 fFacetNormal = (fFinish == polished) ? fGlobalNormal :
1534 GetFacetNormal(fOldMomentum, fGlobalNormal);
1535
1536 PdotN = fOldMomentum * fFacetNormal;
1537 G4double cost1 = -PdotN;
1538 G4double sint2, cost2 = 0.;
1539
1540 if (std::abs(cost1) < 1.0 - fCarTolerance)
1541 {
1542 fSint1 = std::sqrt(1. - cost1 * cost1);
1543 sint2 = fSint1 * fRindex1 / fRindex2;
1544 sintTL = fSint1 * fRindex1 / fCoatedRindex;
1545 }
1546 else
1547 {
1548 fSint1 = 0.0;
1549 sint2 = 0.0;
1550 sintTL = 0.0;
1551 }
1552
1553 if (fSint1 > 0.0)
1554 {
1555 A_trans = fOldMomentum.cross(fFacetNormal);
1556 A_trans = A_trans.unit();
1557 E1_perp = fOldPolarization * A_trans;
1558 E1pp = E1_perp * A_trans;
1559 E1pl = fOldPolarization - E1pp;
1560 E1_parl = E1pl.mag();
1561 }
1562 else
1563 {
1564 A_trans = fOldPolarization;
1565 E1_perp = 0.0;
1566 E1_parl = 1.0;
1567 }
1568
1569 s1 = fRindex1 * cost1;
1570
1571 cost2 = (cost1 > 0.0 ? 1.0 : -1.0) * std::sqrt(1. - sint2 * sint2);
1572
1573 E2_perp = 2. * s1 * E1_perp / (fRindex1 * cost1 + fRindex2 * cost2);
1574 E2_parl = 2. * s1 * E1_parl / (fRindex2 * cost1 + fRindex1 * cost2);
1575 E2_total = E2_perp * E2_perp + E2_parl * E2_parl;
1576
1577 G4double transCoeff = 1. - GetReflectivityThroughThinLayer(
1578 sintTL, E1_perp, E1_parl, wavelength, cost1, cost2);
1579 if (!G4BooleanRand(transCoeff))
1580 {
1581 if(verboseLevel > 2)
1582 G4cout << "Reflection from " << fMaterial1->GetName() << " to "
1583 << fMaterial2->GetName() << G4endl;
1584
1585 //Swap = false;
1586
1587 fStatus = (sintTL >= 1.0) ?
1589
1590 PdotN = fOldMomentum * fFacetNormal;
1591 fNewMomentum = fOldMomentum - (2. * PdotN) * fFacetNormal;
1592
1593 if (fSint1 > 0.0)
1594 { // incident ray oblique
1595 E2_parl = fRindex2 * E2_parl / fRindex1 - E1_parl;
1596 E2_perp = E2_perp - E1_perp;
1597 E2_total = E2_perp * E2_perp + E2_parl * E2_parl;
1598 A_paral = fNewMomentum.cross(A_trans);
1599 A_paral = A_paral.unit();
1600 E2_abs = std::sqrt(E2_total);
1601
1602 fNewPolarization = (E2_parl * A_paral + E2_perp * A_trans) / E2_abs;
1603 }
1604 else
1605 { // incident ray perpendicular
1606 fNewPolarization = (fRindex2 > fRindex1 ? -1. : 1.) * fOldPolarization;
1607 }
1608 }
1609 else
1610 { // photon gets transmitted
1611 if (verboseLevel > 2)
1612 G4cout << "Transmission from " << fMaterial1->GetName() << " to "
1613 << fMaterial2->GetName() << G4endl;
1614
1615 //Inside = !Inside;
1616 through = true;
1617
1618 if (fEfficiency > 0.)
1619 {
1620 DoAbsorption();
1621 return;
1622 }
1623 else
1624 {
1625 fStatus = (sintTL >= 1.0) ? CoatedDielectricFrustratedTransmission :
1627
1628 if (fSint1 > 0.0)
1629 { // incident ray oblique
1630 alpha = cost1 - cost2 * (fRindex2 / fRindex1);
1631 fNewMomentum = fOldMomentum + alpha * fFacetNormal;
1632 fNewMomentum = fNewMomentum.unit();
1633 A_paral = fNewMomentum.cross(A_trans);
1634 A_paral = A_paral.unit();
1635 E2_abs = std::sqrt(E2_total);
1636
1637 fNewPolarization = (E2_parl * A_paral + E2_perp * A_trans) / E2_abs;
1638 }
1639 else
1640 { // incident ray perpendicular
1641 fNewMomentum = fOldMomentum;
1642 fNewPolarization = fOldPolarization;
1643 }
1644 }
1645 }
1646
1647 fOldMomentum = fNewMomentum.unit();
1648 fOldPolarization = fNewPolarization.unit();
1649 if ((fStatus == CoatedDielectricFrustratedTransmission) ||
1650 (fStatus == CoatedDielectricRefraction))
1651 {
1652 done = (fNewMomentum * fGlobalNormal <= 0.0);
1653 }
1654 else
1655 {
1656 done = (fNewMomentum * fGlobalNormal >= -fCarTolerance);
1657 }
1658 } while (!done);
1659}
1660
1661//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1662G4double G4OpBoundaryProcess::GetReflectivityThroughThinLayer(G4double sinTL,
1663 G4double E1_perp,
1664 G4double E1_parl,
1665 G4double wavelength, G4double cost1, G4double cost2)
1666{
1667 G4double gammaTL, costTL;
1668
1669 G4complex i(0, 1);
1670 G4complex rTM, rTE;
1671 G4complex r1toTL, rTLto2;
1672 G4double k0 = 2 * pi / wavelength;
1673
1674 // Angle > Angle limit
1675 if (sinTL >= 1.0)
1676 {
1677 if (fCoatedFrustratedTransmission)
1678 { //Frustrated transmission
1679 G4double coatedRindex_sq = fCoatedRindex * fCoatedRindex;
1680 G4double gamma = fRindex1 * fRindex1 * fSint1 * fSint1 - coatedRindex_sq;
1681 gammaTL = (cost1 > 0.0 ? 1.0 : -1.0) * std::sqrt(gamma);
1682
1683 G4double afact = std::exp(-2 * k0 * fCoatedThickness * gammaTL);
1684 G4double r1 = fRindex1 * cost1;
1685 G4double r2 = fRindex2 * cost2;
1686 G4complex rc = i * gammaTL;
1687 // TE
1688 r1toTL = (r1 - rc) / (r1 + rc);
1689 rTLto2 = (rc - r2) / (rc + r2);
1690 if (cost1 != 0.0)
1691 {
1692 rTE = (r1toTL + rTLto2 * afact) / (1.0 + r1toTL * rTLto2 * afact);
1693 }
1694 // TM
1695 r1toTL = (fRindex1 * rc - coatedRindex_sq * cost1) /
1696 (fRindex1 * rc + coatedRindex_sq * cost1);
1697 rTLto2 = (coatedRindex_sq * cost2 - fRindex2 * rc) /
1698 (coatedRindex_sq * cost2 + fRindex2 * rc);
1699 if (cost1 != 0.0)
1700 {
1701 rTM = (r1toTL + rTLto2 * afact) / (1.0 + r1toTL * rTLto2 * afact);
1702 }
1703 }
1704 else
1705 { //Total reflection
1706 return(1.);
1707 }
1708 }
1709
1710 // Angle <= Angle limit
1711 else //if (sinTL < 1.0)
1712 {
1713 costTL = (cost1 > 0.0 ? 1.0 : -1.0) * std::sqrt(1. - sinTL * sinTL);
1714 G4double rc = fCoatedRindex * costTL;
1715 G4complex afact = std::exp(2.0 * i * k0 * fCoatedThickness * rc);
1716
1717 // TE
1718 G4double r1 = fRindex1 * cost1;
1719 G4double r2 = fRindex2 * cost2;
1720
1721 r1toTL = (r1 - rc) / (r1 + rc);
1722 rTLto2 = (rc - r2) / (r2 + rc);
1723 if (cost1 != 0.0)
1724 {
1725 rTE = (r1toTL + rTLto2 * afact) / (1.0 + r1toTL * rTLto2 * afact);
1726 }
1727 // TM
1728 r1 = fRindex1 * costTL;
1729 r2 = fRindex2 * costTL;
1730 rc = fCoatedRindex * cost1;
1731 r1toTL = (r1 - rc) / (r1 + rc);
1732
1733 rc = fCoatedRindex * cost2;
1734 rTLto2 = (rc - r2) / (rc + r2);
1735 if (cost1 != 0.0)
1736 {
1737 rTM = (r1toTL + rTLto2 * afact) / (1.0 + r1toTL * rTLto2 * afact);
1738 }
1739 }
1740
1741 G4double E1_perp_sq = E1_perp * E1_perp;
1742 G4double E1_parl_sq = E1_parl * E1_parl;
1743 G4double E1_sq = E1_perp_sq + E1_parl_sq;
1744
1745 G4complex Reflectivity_TE = (rTE * conj(rTE)) * E1_perp_sq / E1_sq;
1746 G4complex Reflectivity_TM = (rTM * conj(rTM)) * E1_parl_sq / E1_sq;
1747 G4complex Reflectivity = Reflectivity_TE + Reflectivity_TM;
1748
1749 return real(Reflectivity);
1750}
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
@ FatalException
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4ForceCondition
@ Forced
@ kCOATEDFRUSTRATEDTRANSMISSION
G4PhysicsFreeVector G4MaterialPropertyVector
G4OpBoundaryProcessStatus
@ PolishedTiOAirReflection
@ GroundTeflonAirReflection
@ EtchedVM2000AirReflection
@ NotAtBoundary
@ Transmission
@ CoatedDielectricReflection
@ EtchedVM2000GlueReflection
@ GroundLumirrorGlueReflection
@ LobeReflection
@ GroundTiOAirReflection
@ GroundTyvekAirReflection
@ FresnelReflection
@ GroundAirReflection
@ EtchedAirReflection
@ PolishedVM2000GlueReflection
@ PolishedAirReflection
@ PolishedTeflonAirReflection
@ SameMaterial
@ EtchedTyvekAirReflection
@ SpikeReflection
@ EtchedLumirrorGlueReflection
@ GroundVM2000AirReflection
@ PolishedTyvekAirReflection
@ PolishedVM2000AirReflection
@ EtchedTiOAirReflection
@ EtchedTeflonAirReflection
@ GroundVM2000GlueReflection
@ BackScattering
@ PolishedLumirrorGlueReflection
@ TotalInternalReflection
@ StepTooSmall
@ CoatedDielectricFrustratedTransmission
@ PolishedLumirrorAirReflection
@ EtchedLumirrorAirReflection
@ GroundLumirrorAirReflection
@ LambertianReflection
@ CoatedDielectricRefraction
@ FresnelRefraction
@ fOpBoundary
@ unified
@ groundfrontpainted
@ polished
@ ground
@ polishedfrontpainted
@ groundbackpainted
G4ProcessType
@ fGeomBoundary
G4SurfaceType
@ dielectric_metal
@ dielectric_LUT
@ dielectric_dielectric
@ dielectric_LUTDAVIS
@ dielectric_dichroic
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
std::complex< G4double > G4complex
Definition G4Types.hh:88
int G4int
Definition G4Types.hh:85
@ NotAtBoundary
@ LambertianReflection
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
Hep3Vector orthogonal() const
void setY(double)
double mag2() const
Hep3Vector cross(const Hep3Vector &) const
double angle(const Hep3Vector &) const
void setZ(double)
double mag() const
void set(double x, double y, double z)
void setX(double)
Hep3Vector & rotateUz(const Hep3Vector &)
Hep3Vector & rotate(double, const Hep3Vector &)
const G4ThreeVector & GetMomentumDirection() const
G4double GetTotalMomentum() const
const G4ThreeVector & GetPolarization() const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
static G4LogicalBorderSurface * GetSurface(const G4VPhysicalVolume *vol1, const G4VPhysicalVolume *vol2)
static G4LogicalSkinSurface * GetSurface(const G4LogicalVolume *vol)
G4LogicalSurface is an abstraction of a geometrical surface, it is an abstract base class for differe...
G4SurfaceProperty * GetSurfaceProperty() const
G4bool ConstPropertyExists(const G4String &key) const
G4double GetConstProperty(const G4String &key) const
G4MaterialPropertyVector * GetProperty(const char *key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
virtual ~G4OpBoundaryProcess()
G4OpBoundaryProcess(const G4String &processName="OpBoundary", G4ProcessType type=fOptical)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
virtual void SetInvokeSD(G4bool)
virtual G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *condition) override
virtual void PreparePhysicsTable(const G4ParticleDefinition &) override
G4bool GetBoundaryInvokeSD() const
G4int GetBoundaryVerboseLevel() const
static G4OpticalParameters * Instance()
static const G4Step * GetHyperStep()
G4double Value(const G4double energy, std::size_t &lastidx) const
G4StepStatus GetStepStatus() const
G4Material * GetMaterial() const
const G4ThreeVector & GetPosition() const
G4VSensitiveDetector * GetSensitiveDetector() const
G4VPhysicalVolume * GetPhysicalVolume() const
void AddTotalEnergyDeposit(G4double value)
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4double GetVelocity() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetStepLength() const
std::vector< G4Navigator * >::iterator GetActiveNavigatorsIterator()
static G4TransportationManager * GetTransportationManager()
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4VDiscreteProcess(const G4String &aName, G4ProcessType aType=fNotDefined)
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
G4LogicalVolume * GetMotherLogical() const
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
G4ParticleChange aParticleChange
G4int verboseLevel
void SetProcessSubType(G4int)
const G4String & GetProcessName() const
G4bool Hit(G4Step *aStep)
const G4double pi
const G4double hc
[MeV*fm]
void swap(T &lhs, T &rhs)
Definition pugixml.cc:7568
void G4SwapObj(T *a, T *b)
Definition templates.hh:112
void G4SwapPtr(T *&a, T *&b)
Definition templates.hh:104
#define DBL_MAX
Definition templates.hh:62