Geant4 11.1.1
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
90//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92 G4ProcessType ptype)
93 : G4VDiscreteProcess(processName, ptype)
94{
95 Initialise();
96
97 if(verboseLevel > 0)
98 {
99 G4cout << GetProcessName() << " is created " << G4endl;
100 }
102
103 fStatus = Undefined;
104 fModel = glisur;
105 fFinish = polished;
106 fReflectivity = 1.;
107 fEfficiency = 0.;
108 fTransmittance = 0.;
109 fSurfaceRoughness = 0.;
110 fProb_sl = 0.;
111 fProb_ss = 0.;
112 fProb_bs = 0.;
113
114 fRealRIndexMPV = nullptr;
115 fImagRIndexMPV = nullptr;
116 fMaterial1 = nullptr;
117 fMaterial2 = nullptr;
118 fOpticalSurface = nullptr;
120
121 f_iTE = f_iTM = 0;
122 fPhotonMomentum = 0.;
123 fRindex1 = fRindex2 = 1.;
124 fSint1 = 0.;
125 fDichroicVector = nullptr;
126
127 fNumWarnings = 0;
128}
129
130//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
135{
136 Initialise();
137}
138
139//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141{
145}
146
147//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
149 const G4Step& aStep)
150{
151 fStatus = Undefined;
154
155 // Get hyperStep from G4ParallelWorldProcess
156 // NOTE: PostSetpDoIt of this process to be invoked after
157 // G4ParallelWorldProcess!
158 const G4Step* pStep = &aStep;
160 if(hStep != nullptr)
161 pStep = hStep;
162
164 {
165 fMaterial1 = pStep->GetPreStepPoint()->GetMaterial();
166 fMaterial2 = pStep->GetPostStepPoint()->GetMaterial();
167 }
168 else
169 {
170 fStatus = NotAtBoundary;
171 if(verboseLevel > 1)
172 BoundaryProcessVerbose();
173 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
174 }
175
178
179 if(verboseLevel > 1)
180 {
181 G4cout << " Photon at Boundary! " << G4endl;
182 if(thePrePV != nullptr)
183 G4cout << " thePrePV: " << thePrePV->GetName() << G4endl;
184 if(thePostPV != nullptr)
185 G4cout << " thePostPV: " << thePostPV->GetName() << G4endl;
186 }
187
188 G4double stepLength = aTrack.GetStepLength();
189 if(stepLength <= fCarTolerance)
190 {
191 fStatus = StepTooSmall;
192 if(verboseLevel > 1)
193 BoundaryProcessVerbose();
194
195 G4MaterialPropertyVector* groupvel = nullptr;
197 if(aMPT != nullptr)
198 {
199 groupvel = aMPT->GetProperty(kGROUPVEL);
200 }
201
202 if(groupvel != nullptr)
203 {
205 groupvel->Value(fPhotonMomentum, idx_groupvel));
206 }
207 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
208 }
209 else if (stepLength <= 10.*fCarTolerance && fNumWarnings < 10)
210 { // see bug 2510
211 ++fNumWarnings;
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(fNumWarnings == 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 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;
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();
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(fFinish == polishedbackpainted || fFinish == groundbackpainted)
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();
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(fFinish == polishedbackpainted || fFinish == groundbackpainted)
418 {
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();
454 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
455 }
456 }
457 if(fFinish == polishedbackpainted || fFinish == groundbackpainted)
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 fStatus = Transmission;
471 fNewMomentum = fOldMomentum;
472 fNewPolarization = fOldPolarization;
473 }
474 else
475 {
476 if(fFinish == polishedfrontpainted)
477 {
478 DoReflection();
479 }
480 else if(fFinish == groundfrontpainted)
481 {
482 fStatus = LambertianReflection;
483 DoReflection();
484 }
485 else
486 {
487 DielectricDielectric();
488 }
489 }
490 }
491 }
492 else if(type == dielectric_metal)
493 {
494 DielectricMetal();
495 }
496 else if(type == dielectric_LUT)
497 {
498 DielectricLUT();
499 }
500 else if(type == dielectric_LUTDAVIS)
501 {
502 DielectricLUTDAVIS();
503 }
504 else if(type == dielectric_dichroic)
505 {
506 DielectricDichroic();
507 }
508 else if(type == coated)
509 {
510 CoatedDielectricDielectric();
511 }
512 else
513 {
515 ed << " PostStepDoIt(): Illegal boundary type." << G4endl;
516 G4Exception("G4OpBoundaryProcess", "OpBoun04", JustWarning, ed);
517 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
518 }
519
520 fNewMomentum = fNewMomentum.unit();
521 fNewPolarization = fNewPolarization.unit();
522
523 if(verboseLevel > 1)
524 {
525 G4cout << " New Momentum Direction: " << fNewMomentum << G4endl
526 << " New Polarization: " << fNewPolarization << G4endl;
527 BoundaryProcessVerbose();
528 }
529
531 aParticleChange.ProposePolarization(fNewPolarization);
532
533 if(fStatus == FresnelRefraction || fStatus == Transmission)
534 {
535 // not all surface types check that fMaterial2 has an MPT
537 G4MaterialPropertyVector* groupvel = nullptr;
538 if(aMPT != nullptr)
539 {
540 groupvel = aMPT->GetProperty(kGROUPVEL);
541 }
542 if(groupvel != nullptr)
543 {
545 groupvel->Value(fPhotonMomentum, idx_groupvel));
546 }
547 }
548
549 if(fStatus == Detection && fInvokeSD)
550 InvokeSD(pStep);
551 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
552}
553
554//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
555void G4OpBoundaryProcess::BoundaryProcessVerbose() const
556{
557 G4cout << " *** ";
558 if(fStatus == Undefined)
559 G4cout << "Undefined";
560 else if(fStatus == Transmission)
561 G4cout << "Transmission";
562 else if(fStatus == FresnelRefraction)
563 G4cout << "FresnelRefraction";
564 else if(fStatus == FresnelReflection)
565 G4cout << "FresnelReflection";
566 else if(fStatus == TotalInternalReflection)
567 G4cout << "TotalInternalReflection";
568 else if(fStatus == LambertianReflection)
569 G4cout << "LambertianReflection";
570 else if(fStatus == LobeReflection)
571 G4cout << "LobeReflection";
572 else if(fStatus == SpikeReflection)
573 G4cout << "SpikeReflection";
574 else if(fStatus == BackScattering)
575 G4cout << "BackScattering";
576 else if(fStatus == PolishedLumirrorAirReflection)
577 G4cout << "PolishedLumirrorAirReflection";
578 else if(fStatus == PolishedLumirrorGlueReflection)
579 G4cout << "PolishedLumirrorGlueReflection";
580 else if(fStatus == PolishedAirReflection)
581 G4cout << "PolishedAirReflection";
582 else if(fStatus == PolishedTeflonAirReflection)
583 G4cout << "PolishedTeflonAirReflection";
584 else if(fStatus == PolishedTiOAirReflection)
585 G4cout << "PolishedTiOAirReflection";
586 else if(fStatus == PolishedTyvekAirReflection)
587 G4cout << "PolishedTyvekAirReflection";
588 else if(fStatus == PolishedVM2000AirReflection)
589 G4cout << "PolishedVM2000AirReflection";
590 else if(fStatus == PolishedVM2000GlueReflection)
591 G4cout << "PolishedVM2000GlueReflection";
592 else if(fStatus == EtchedLumirrorAirReflection)
593 G4cout << "EtchedLumirrorAirReflection";
594 else if(fStatus == EtchedLumirrorGlueReflection)
595 G4cout << "EtchedLumirrorGlueReflection";
596 else if(fStatus == EtchedAirReflection)
597 G4cout << "EtchedAirReflection";
598 else if(fStatus == EtchedTeflonAirReflection)
599 G4cout << "EtchedTeflonAirReflection";
600 else if(fStatus == EtchedTiOAirReflection)
601 G4cout << "EtchedTiOAirReflection";
602 else if(fStatus == EtchedTyvekAirReflection)
603 G4cout << "EtchedTyvekAirReflection";
604 else if(fStatus == EtchedVM2000AirReflection)
605 G4cout << "EtchedVM2000AirReflection";
606 else if(fStatus == EtchedVM2000GlueReflection)
607 G4cout << "EtchedVM2000GlueReflection";
608 else if(fStatus == GroundLumirrorAirReflection)
609 G4cout << "GroundLumirrorAirReflection";
610 else if(fStatus == GroundLumirrorGlueReflection)
611 G4cout << "GroundLumirrorGlueReflection";
612 else if(fStatus == GroundAirReflection)
613 G4cout << "GroundAirReflection";
614 else if(fStatus == GroundTeflonAirReflection)
615 G4cout << "GroundTeflonAirReflection";
616 else if(fStatus == GroundTiOAirReflection)
617 G4cout << "GroundTiOAirReflection";
618 else if(fStatus == GroundTyvekAirReflection)
619 G4cout << "GroundTyvekAirReflection";
620 else if(fStatus == GroundVM2000AirReflection)
621 G4cout << "GroundVM2000AirReflection";
622 else if(fStatus == GroundVM2000GlueReflection)
623 G4cout << "GroundVM2000GlueReflection";
624 else if(fStatus == Absorption)
625 G4cout << "Absorption";
626 else if(fStatus == Detection)
627 G4cout << "Detection";
628 else if(fStatus == NotAtBoundary)
629 G4cout << "NotAtBoundary";
630 else if(fStatus == SameMaterial)
631 G4cout << "SameMaterial";
632 else if(fStatus == StepTooSmall)
633 G4cout << "StepTooSmall";
634 else if(fStatus == NoRINDEX)
635 G4cout << "NoRINDEX";
636 else if(fStatus == Dichroic)
637 G4cout << "Dichroic Transmission";
638 else if(fStatus == CoatedDielectricReflection)
639 G4cout << "Coated Dielectric Reflection";
640 else if(fStatus == CoatedDielectricRefraction)
641 G4cout << "Coated Dielectric Refraction";
643 G4cout << "Coated Dielectric Frustrated Transmission";
644
645 G4cout << " ***" << G4endl;
646}
647
648//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
649G4ThreeVector G4OpBoundaryProcess::GetFacetNormal(
650 const G4ThreeVector& momentum, const G4ThreeVector& normal) const
651{
652 G4ThreeVector facetNormal;
653 if(fModel == unified || fModel == LUT || fModel == DAVIS)
654 {
655 /* This function codes alpha to a random value taken from the
656 distribution p(alpha) = g(alpha; 0, sigma_alpha)*std::sin(alpha),
657 for alpha > 0 and alpha < 90, where g(alpha; 0, sigma_alpha) is a
658 gaussian distribution with mean 0 and standard deviation sigma_alpha. */
659
660 G4double sigma_alpha = 0.0;
661 if(fOpticalSurface)
662 sigma_alpha = fOpticalSurface->GetSigmaAlpha();
663 if(sigma_alpha == 0.0)
664 {
665 return normal;
666 }
667
668 G4double f_max = std::min(1.0, 4. * sigma_alpha);
669 G4double alpha, phi, sinAlpha;
670
671 do
672 { // Loop checking, 13-Aug-2015, Peter Gumplinger
673 do
674 { // Loop checking, 13-Aug-2015, Peter Gumplinger
675 alpha = G4RandGauss::shoot(0.0, sigma_alpha);
676 sinAlpha = std::sin(alpha);
677 } while(G4UniformRand() * f_max > sinAlpha || alpha >= halfpi);
678
679 phi = G4UniformRand() * twopi;
680 facetNormal.set(sinAlpha * std::cos(phi), sinAlpha * std::sin(phi),
681 std::cos(alpha));
682 facetNormal.rotateUz(normal);
683 } while(momentum * facetNormal >= 0.0);
684 }
685 else
686 {
687 G4double polish = 1.0;
688 if(fOpticalSurface)
689 polish = fOpticalSurface->GetPolish();
690 if(polish < 1.0)
691 {
692 do
693 { // Loop checking, 13-Aug-2015, Peter Gumplinger
694 G4ThreeVector smear;
695 do
696 { // Loop checking, 13-Aug-2015, Peter Gumplinger
697 smear.setX(2. * G4UniformRand() - 1.);
698 smear.setY(2. * G4UniformRand() - 1.);
699 smear.setZ(2. * G4UniformRand() - 1.);
700 } while(smear.mag2() > 1.0);
701 facetNormal = normal + (1. - polish) * smear;
702 } while(momentum * facetNormal >= 0.0);
703 facetNormal = facetNormal.unit();
704 }
705 else
706 {
707 facetNormal = normal;
708 }
709 }
710 return facetNormal;
711}
712
713//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
714void G4OpBoundaryProcess::DielectricMetal()
715{
716 G4int n = 0;
717 G4double rand;
718 G4ThreeVector A_trans;
719
720 do
721 {
722 ++n;
723 rand = G4UniformRand();
724 if(rand > fReflectivity && n == 1)
725 {
726 if(rand > fReflectivity + fTransmittance)
727 {
728 DoAbsorption();
729 }
730 else
731 {
732 fStatus = Transmission;
733 fNewMomentum = fOldMomentum;
734 fNewPolarization = fOldPolarization;
735 }
736 break;
737 }
738 else
739 {
740 if(fRealRIndexMPV && fImagRIndexMPV)
741 {
742 if(n > 1)
743 {
744 CalculateReflectivity();
745 if(!G4BooleanRand(fReflectivity))
746 {
747 DoAbsorption();
748 break;
749 }
750 }
751 }
752 if(fModel == glisur || fFinish == polished)
753 {
754 DoReflection();
755 }
756 else
757 {
758 if(n == 1)
759 ChooseReflection();
760 if(fStatus == LambertianReflection)
761 {
762 DoReflection();
763 }
764 else if(fStatus == BackScattering)
765 {
766 fNewMomentum = -fOldMomentum;
767 fNewPolarization = -fOldPolarization;
768 }
769 else
770 {
771 if(fStatus == LobeReflection)
772 {
773 if(!fRealRIndexMPV || !fImagRIndexMPV)
774 {
775 fFacetNormal = GetFacetNormal(fOldMomentum, fGlobalNormal);
776 }
777 // else
778 // case of complex rindex needs to be implemented
779 }
780 fNewMomentum =
781 fOldMomentum - 2. * fOldMomentum * fFacetNormal * fFacetNormal;
782
783 if(f_iTE > 0 && f_iTM > 0)
784 {
785 fNewPolarization =
786 -fOldPolarization +
787 (2. * fOldPolarization * fFacetNormal * fFacetNormal);
788 }
789 else if(f_iTE > 0)
790 {
791 A_trans = (fSint1 > 0.0) ? fOldMomentum.cross(fFacetNormal).unit()
792 : fOldPolarization;
793 fNewPolarization = -A_trans;
794 }
795 else if(f_iTM > 0)
796 {
797 fNewPolarization =
798 -fNewMomentum.cross(A_trans).unit(); // = -A_paral
799 }
800 }
801 }
802 fOldMomentum = fNewMomentum;
803 fOldPolarization = fNewPolarization;
804 }
805 // Loop checking, 13-Aug-2015, Peter Gumplinger
806 } while(fNewMomentum * fGlobalNormal < 0.0);
807}
808
809//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
810void G4OpBoundaryProcess::DielectricLUT()
811{
812 G4int thetaIndex, phiIndex;
813 G4double angularDistVal, thetaRad, phiRad;
814 G4ThreeVector perpVectorTheta, perpVectorPhi;
815
817 G4int(fFinish) + (G4int(NoRINDEX) - G4int(groundbackpainted)));
818
819 G4int thetaIndexMax = fOpticalSurface->GetThetaIndexMax();
820 G4int phiIndexMax = fOpticalSurface->GetPhiIndexMax();
821
822 G4double rand;
823
824 do
825 {
826 rand = G4UniformRand();
827 if(rand > fReflectivity)
828 {
829 if(rand > fReflectivity + fTransmittance)
830 {
831 DoAbsorption();
832 }
833 else
834 {
835 fStatus = Transmission;
836 fNewMomentum = fOldMomentum;
837 fNewPolarization = fOldPolarization;
838 }
839 break;
840 }
841 else
842 {
843 // Calculate Angle between Normal and Photon Momentum
844 G4double anglePhotonToNormal = fOldMomentum.angle(-fGlobalNormal);
845 // Round to closest integer: LBNL model array has 91 values
846 G4int angleIncident = (G4int)std::lrint(anglePhotonToNormal / CLHEP::deg);
847
848 // Take random angles THETA and PHI,
849 // and see if below Probability - if not - Redo
850 do
851 {
852 thetaIndex = (G4int)G4RandFlat::shootInt(thetaIndexMax - 1);
853 phiIndex = (G4int)G4RandFlat::shootInt(phiIndexMax - 1);
854 // Find probability with the new indeces from LUT
855 angularDistVal = fOpticalSurface->GetAngularDistributionValue(
856 angleIncident, thetaIndex, phiIndex);
857 // Loop checking, 13-Aug-2015, Peter Gumplinger
858 } while(!G4BooleanRand(angularDistVal));
859
860 thetaRad = G4double(-90 + 4 * thetaIndex) * pi / 180.;
861 phiRad = G4double(-90 + 5 * phiIndex) * pi / 180.;
862 // Rotate Photon Momentum in Theta, then in Phi
863 fNewMomentum = -fOldMomentum;
864
865 perpVectorTheta = fNewMomentum.cross(fGlobalNormal);
866 if(perpVectorTheta.mag() < fCarTolerance)
867 {
868 perpVectorTheta = fNewMomentum.orthogonal();
869 }
870 fNewMomentum =
871 fNewMomentum.rotate(anglePhotonToNormal - thetaRad, perpVectorTheta);
872 perpVectorPhi = perpVectorTheta.cross(fNewMomentum);
873 fNewMomentum = fNewMomentum.rotate(-phiRad, perpVectorPhi);
874
875 // Rotate Polarization too:
876 fFacetNormal = (fNewMomentum - fOldMomentum).unit();
877 fNewPolarization = -fOldPolarization +
878 (2. * fOldPolarization * fFacetNormal * fFacetNormal);
879 }
880 // Loop checking, 13-Aug-2015, Peter Gumplinger
881 } while(fNewMomentum * fGlobalNormal <= 0.0);
882}
883
884//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
885void G4OpBoundaryProcess::DielectricLUTDAVIS()
886{
887 G4int angindex, random, angleIncident;
888 G4double reflectivityValue, elevation, azimuth;
889 G4double anglePhotonToNormal;
890
891 G4int lutbin = fOpticalSurface->GetLUTbins();
892 G4double rand = G4UniformRand();
893
894 G4double sinEl;
895 G4ThreeVector u, vNorm, w;
896
897 do
898 {
899 anglePhotonToNormal = fOldMomentum.angle(-fGlobalNormal);
900
901 // Davis model has 90 reflection bins: round down
902 // don't allow angleIncident to be 90 for anglePhotonToNormal close to 90
903 angleIncident = std::min(
904 static_cast<G4int>(std::floor(anglePhotonToNormal / CLHEP::deg)), 89);
905 reflectivityValue = fOpticalSurface->GetReflectivityLUTValue(angleIncident);
906
907 if(rand > reflectivityValue)
908 {
909 if(fEfficiency > 0.)
910 {
911 DoAbsorption();
912 break;
913 }
914 else
915 {
916 fStatus = Transmission;
917
918 if(angleIncident <= 0.01)
919 {
920 fNewMomentum = fOldMomentum;
921 break;
922 }
923
924 do
925 {
926 random = (G4int)G4RandFlat::shootInt(1, lutbin + 1);
927 angindex =
928 (((random * 2) - 1)) + angleIncident * lutbin * 2 + 3640000;
929
930 azimuth =
931 fOpticalSurface->GetAngularDistributionValueLUT(angindex - 1);
932 elevation = fOpticalSurface->GetAngularDistributionValueLUT(angindex);
933 } while(elevation == 0. && azimuth == 0.);
934
935 sinEl = std::sin(elevation);
936 vNorm = (fGlobalNormal.cross(fOldMomentum)).unit();
937 u = vNorm.cross(fGlobalNormal) * (sinEl * std::cos(azimuth));
938 vNorm *= (sinEl * std::sin(azimuth));
939 // fGlobalNormal shouldn't be modified here
940 w = (fGlobalNormal *= std::cos(elevation));
941 fNewMomentum = u + vNorm + w;
942
943 // Rotate Polarization too:
944 fFacetNormal = (fNewMomentum - fOldMomentum).unit();
945 fNewPolarization = -fOldPolarization + (2. * fOldPolarization *
946 fFacetNormal * fFacetNormal);
947 }
948 }
949 else
950 {
951 fStatus = LobeReflection;
952
953 if(angleIncident == 0)
954 {
955 fNewMomentum = -fOldMomentum;
956 break;
957 }
958
959 do
960 {
961 random = (G4int)G4RandFlat::shootInt(1, lutbin + 1);
962 angindex = (((random * 2) - 1)) + (angleIncident - 1) * lutbin * 2;
963
964 azimuth = fOpticalSurface->GetAngularDistributionValueLUT(angindex - 1);
965 elevation = fOpticalSurface->GetAngularDistributionValueLUT(angindex);
966 } while(elevation == 0. && azimuth == 0.);
967
968 sinEl = std::sin(elevation);
969 vNorm = (fGlobalNormal.cross(fOldMomentum)).unit();
970 u = vNorm.cross(fGlobalNormal) * (sinEl * std::cos(azimuth));
971 vNorm *= (sinEl * std::sin(azimuth));
972 // fGlobalNormal shouldn't be modified here
973 w = (fGlobalNormal *= std::cos(elevation));
974
975 fNewMomentum = u + vNorm + w;
976
977 // Rotate Polarization too: (needs revision)
978 fNewPolarization = fOldPolarization;
979 }
980 } while(fNewMomentum * fGlobalNormal <= 0.0);
981}
982
983//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
984void G4OpBoundaryProcess::DielectricDichroic()
985{
986 // Calculate Angle between Normal and Photon Momentum
987 G4double anglePhotonToNormal = fOldMomentum.angle(-fGlobalNormal);
988
989 // Round it to closest integer
990 G4double angleIncident = std::floor(180. / pi * anglePhotonToNormal + 0.5);
991
992 if(!fDichroicVector)
993 {
994 if(fOpticalSurface)
995 fDichroicVector = fOpticalSurface->GetDichroicVector();
996 }
997
998 if(fDichroicVector)
999 {
1000 G4double wavelength = h_Planck * c_light / fPhotonMomentum;
1001 fTransmittance = fDichroicVector->Value(wavelength / nm, angleIncident,
1002 idx_dichroicX, idx_dichroicY) *
1003 perCent;
1004 // G4cout << "wavelength: " << std::floor(wavelength/nm)
1005 // << "nm" << G4endl;
1006 // G4cout << "Incident angle: " << angleIncident << "deg" << G4endl;
1007 // G4cout << "Transmittance: "
1008 // << std::floor(fTransmittance/perCent) << "%" << G4endl;
1009 }
1010 else
1011 {
1013 ed << " G4OpBoundaryProcess/DielectricDichroic(): "
1014 << " The dichroic surface has no G4Physics2DVector" << G4endl;
1015 G4Exception("G4OpBoundaryProcess::DielectricDichroic", "OpBoun03",
1016 FatalException, ed,
1017 "A dichroic surface must have an associated G4Physics2DVector");
1018 }
1019
1020 if(!G4BooleanRand(fTransmittance))
1021 { // Not transmitted, so reflect
1022 if(fModel == glisur || fFinish == polished)
1023 {
1024 DoReflection();
1025 }
1026 else
1027 {
1028 ChooseReflection();
1029 if(fStatus == LambertianReflection)
1030 {
1031 DoReflection();
1032 }
1033 else if(fStatus == BackScattering)
1034 {
1035 fNewMomentum = -fOldMomentum;
1036 fNewPolarization = -fOldPolarization;
1037 }
1038 else
1039 {
1040 G4double PdotN, EdotN;
1041 do
1042 {
1043 if(fStatus == LobeReflection)
1044 {
1045 fFacetNormal = GetFacetNormal(fOldMomentum, fGlobalNormal);
1046 }
1047 PdotN = fOldMomentum * fFacetNormal;
1048 fNewMomentum = fOldMomentum - (2. * PdotN) * fFacetNormal;
1049 // Loop checking, 13-Aug-2015, Peter Gumplinger
1050 } while(fNewMomentum * fGlobalNormal <= 0.0);
1051
1052 EdotN = fOldPolarization * fFacetNormal;
1053 fNewPolarization = -fOldPolarization + (2. * EdotN) * fFacetNormal;
1054 }
1055 }
1056 }
1057 else
1058 {
1059 fStatus = Dichroic;
1060 fNewMomentum = fOldMomentum;
1061 fNewPolarization = fOldPolarization;
1062 }
1063}
1064
1065//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1066void G4OpBoundaryProcess::DielectricDielectric()
1067{
1068 G4bool inside = false;
1069 G4bool swap = false;
1070
1071 if(fFinish == polished)
1072 {
1073 fFacetNormal = fGlobalNormal;
1074 }
1075 else
1076 {
1077 fFacetNormal = GetFacetNormal(fOldMomentum, fGlobalNormal);
1078 }
1079 G4double cost1 = -fOldMomentum * fFacetNormal;
1080 G4double cost2 = 0.;
1081 G4double sint2 = 0.;
1082
1083 G4bool surfaceRoughnessCriterionPass = true;
1084 if(fSurfaceRoughness != 0. && fRindex1 > fRindex2)
1085 {
1086 G4double wavelength = h_Planck * c_light / fPhotonMomentum;
1087 G4double surfaceRoughnessCriterion = std::exp(-std::pow(
1088 (4. * pi * fSurfaceRoughness * fRindex1 * cost1 / wavelength), 2));
1089 surfaceRoughnessCriterionPass = G4BooleanRand(surfaceRoughnessCriterion);
1090 }
1091
1092leap:
1093
1094 G4bool through = false;
1095 G4bool done = false;
1096
1097 G4ThreeVector A_trans, A_paral, E1pp, E1pl;
1098 G4double E1_perp, E1_parl;
1099 G4double s1, s2, E2_perp, E2_parl, E2_total, transCoeff;
1100 G4double E2_abs, C_parl, C_perp;
1102
1103 do
1104 {
1105 if(through)
1106 {
1107 swap = !swap;
1108 through = false;
1109 fGlobalNormal = -fGlobalNormal;
1110 G4SwapPtr(fMaterial1, fMaterial2);
1111 G4SwapObj(&fRindex1, &fRindex2);
1112 }
1113
1114 if(fFinish == polished)
1115 {
1116 fFacetNormal = fGlobalNormal;
1117 }
1118 else
1119 {
1120 fFacetNormal = GetFacetNormal(fOldMomentum, fGlobalNormal);
1121 }
1122
1123 cost1 = -fOldMomentum * fFacetNormal;
1124 if(std::abs(cost1) < 1.0 - fCarTolerance)
1125 {
1126 fSint1 = std::sqrt(1. - cost1 * cost1);
1127 sint2 = fSint1 * fRindex1 / fRindex2; // *** Snell's Law ***
1128 // this isn't a sine as we might expect from the name; can be > 1
1129 }
1130 else
1131 {
1132 fSint1 = 0.0;
1133 sint2 = 0.0;
1134 }
1135
1136 // TOTAL INTERNAL REFLECTION
1137 if(sint2 >= 1.0)
1138 {
1139 swap = false;
1140
1141 fStatus = TotalInternalReflection;
1142 if(!surfaceRoughnessCriterionPass)
1143 fStatus = LambertianReflection;
1144 if(fModel == unified && fFinish != polished)
1145 ChooseReflection();
1146 if(fStatus == LambertianReflection)
1147 {
1148 DoReflection();
1149 }
1150 else if(fStatus == BackScattering)
1151 {
1152 fNewMomentum = -fOldMomentum;
1153 fNewPolarization = -fOldPolarization;
1154 }
1155 else
1156 {
1157 fNewMomentum =
1158 fOldMomentum - 2. * fOldMomentum * fFacetNormal * fFacetNormal;
1159 fNewPolarization = -fOldPolarization + (2. * fOldPolarization *
1160 fFacetNormal * fFacetNormal);
1161 }
1162 }
1163 // NOT TIR
1164 else if(sint2 < 1.0)
1165 {
1166 // Calculate amplitude for transmission (Q = P x N)
1167 if(cost1 > 0.0)
1168 {
1169 cost2 = std::sqrt(1. - sint2 * sint2);
1170 }
1171 else
1172 {
1173 cost2 = -std::sqrt(1. - sint2 * sint2);
1174 }
1175
1176 if(fSint1 > 0.0)
1177 {
1178 A_trans = (fOldMomentum.cross(fFacetNormal)).unit();
1179 E1_perp = fOldPolarization * A_trans;
1180 E1pp = E1_perp * A_trans;
1181 E1pl = fOldPolarization - E1pp;
1182 E1_parl = E1pl.mag();
1183 }
1184 else
1185 {
1186 A_trans = fOldPolarization;
1187 // Here we Follow Jackson's conventions and set the parallel
1188 // component = 1 in case of a ray perpendicular to the surface
1189 E1_perp = 0.0;
1190 E1_parl = 1.0;
1191 }
1192
1193 s1 = fRindex1 * cost1;
1194 E2_perp = 2. * s1 * E1_perp / (fRindex1 * cost1 + fRindex2 * cost2);
1195 E2_parl = 2. * s1 * E1_parl / (fRindex2 * cost1 + fRindex1 * cost2);
1196 E2_total = E2_perp * E2_perp + E2_parl * E2_parl;
1197 s2 = fRindex2 * cost2 * E2_total;
1198
1199 if(fTransmittance > 0.)
1200 transCoeff = fTransmittance;
1201 else if(cost1 != 0.0)
1202 transCoeff = s2 / s1;
1203 else
1204 transCoeff = 0.0;
1205
1206 // NOT TIR: REFLECTION
1207 if(!G4BooleanRand(transCoeff))
1208 {
1209 swap = false;
1210 fStatus = FresnelReflection;
1211
1212 if(!surfaceRoughnessCriterionPass)
1213 fStatus = LambertianReflection;
1214 if(fModel == unified && fFinish != polished)
1215 ChooseReflection();
1216 if(fStatus == LambertianReflection)
1217 {
1218 DoReflection();
1219 }
1220 else if(fStatus == BackScattering)
1221 {
1222 fNewMomentum = -fOldMomentum;
1223 fNewPolarization = -fOldPolarization;
1224 }
1225 else
1226 {
1227 fNewMomentum =
1228 fOldMomentum - 2. * fOldMomentum * fFacetNormal * fFacetNormal;
1229 if(fSint1 > 0.0)
1230 { // incident ray oblique
1231 E2_parl = fRindex2 * E2_parl / fRindex1 - E1_parl;
1232 E2_perp = E2_perp - E1_perp;
1233 E2_total = E2_perp * E2_perp + E2_parl * E2_parl;
1234 A_paral = (fNewMomentum.cross(A_trans)).unit();
1235 E2_abs = std::sqrt(E2_total);
1236 C_parl = E2_parl / E2_abs;
1237 C_perp = E2_perp / E2_abs;
1238
1239 fNewPolarization = C_parl * A_paral + C_perp * A_trans;
1240 }
1241 else
1242 { // incident ray perpendicular
1243 if(fRindex2 > fRindex1)
1244 {
1245 fNewPolarization = -fOldPolarization;
1246 }
1247 else
1248 {
1249 fNewPolarization = fOldPolarization;
1250 }
1251 }
1252 }
1253 }
1254 // NOT TIR: TRANSMISSION
1255 else
1256 {
1257 inside = !inside;
1258 through = true;
1259 fStatus = FresnelRefraction;
1260
1261 if(fSint1 > 0.0)
1262 { // incident ray oblique
1263 alpha = cost1 - cost2 * (fRindex2 / fRindex1);
1264 fNewMomentum = (fOldMomentum + alpha * fFacetNormal).unit();
1265 A_paral = (fNewMomentum.cross(A_trans)).unit();
1266 E2_abs = std::sqrt(E2_total);
1267 C_parl = E2_parl / E2_abs;
1268 C_perp = E2_perp / E2_abs;
1269
1270 fNewPolarization = C_parl * A_paral + C_perp * A_trans;
1271 }
1272 else
1273 { // incident ray perpendicular
1274 fNewMomentum = fOldMomentum;
1275 fNewPolarization = fOldPolarization;
1276 }
1277 }
1278 }
1279
1280 fOldMomentum = fNewMomentum.unit();
1281 fOldPolarization = fNewPolarization.unit();
1282
1283 if(fStatus == FresnelRefraction)
1284 {
1285 done = (fNewMomentum * fGlobalNormal <= 0.0);
1286 }
1287 else
1288 {
1289 done = (fNewMomentum * fGlobalNormal >= -fCarTolerance);
1290 }
1291 // Loop checking, 13-Aug-2015, Peter Gumplinger
1292 } while(!done);
1293
1294 if(inside && !swap)
1295 {
1296 if(fFinish == polishedbackpainted || fFinish == groundbackpainted)
1297 {
1298 G4double rand = G4UniformRand();
1299 if(rand > fReflectivity + fTransmittance)
1300 {
1301 DoAbsorption();
1302 }
1303 else if(rand > fReflectivity)
1304 {
1305 fStatus = Transmission;
1306 fNewMomentum = fOldMomentum;
1307 fNewPolarization = fOldPolarization;
1308 }
1309 else
1310 {
1311 if(fStatus != FresnelRefraction)
1312 {
1313 fGlobalNormal = -fGlobalNormal;
1314 }
1315 else
1316 {
1317 swap = !swap;
1318 G4SwapPtr(fMaterial1, fMaterial2);
1319 G4SwapObj(&fRindex1, &fRindex2);
1320 }
1321 if(fFinish == groundbackpainted)
1322 fStatus = LambertianReflection;
1323
1324 DoReflection();
1325
1326 fGlobalNormal = -fGlobalNormal;
1327 fOldMomentum = fNewMomentum;
1328
1329 goto leap;
1330 }
1331 }
1332 }
1333}
1334
1335//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1338{
1339 *condition = Forced;
1340 return DBL_MAX;
1341}
1342
1343//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1344G4double G4OpBoundaryProcess::GetIncidentAngle()
1345{
1346 return pi - std::acos(fOldMomentum * fFacetNormal /
1347 (fOldMomentum.mag() * fFacetNormal.mag()));
1348}
1349
1350//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1351G4double G4OpBoundaryProcess::GetReflectivity(G4double E1_perp,
1352 G4double E1_parl,
1353 G4double incidentangle,
1354 G4double realRindex,
1355 G4double imaginaryRindex)
1356{
1357 G4complex reflectivity, reflectivity_TE, reflectivity_TM;
1358 G4complex N1(fRindex1, 0.), N2(realRindex, imaginaryRindex);
1359 G4complex cosPhi;
1360
1361 G4complex u(1., 0.); // unit number 1
1362
1363 G4complex numeratorTE; // E1_perp=1 E1_parl=0 -> TE polarization
1364 G4complex numeratorTM; // E1_parl=1 E1_perp=0 -> TM polarization
1365 G4complex denominatorTE, denominatorTM;
1366 G4complex rTM, rTE;
1367
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 cosPhi = std::sqrt(u - ((std::sin(incidentangle) * std::sin(incidentangle)) *
1381 (N1 * N1) / (N2 * N2)));
1382
1383 numeratorTE = N1 * std::cos(incidentangle) - N2 * cosPhi;
1384 denominatorTE = N1 * std::cos(incidentangle) + N2 * cosPhi;
1385 rTE = numeratorTE / denominatorTE;
1386
1387 numeratorTM = N2 * std::cos(incidentangle) - N1 * cosPhi;
1388 denominatorTM = N2 * std::cos(incidentangle) + N1 * cosPhi;
1389 rTM = numeratorTM / denominatorTM;
1390
1391 // This is my (PG) calculaton for reflectivity on a metallic surface
1392 // depending on the fraction of TE and TM polarization
1393 // when TE polarization, E1_parl=0 and E1_perp=1, R=abs(rTE)^2 and
1394 // when TM polarization, E1_parl=1 and E1_perp=0, R=abs(rTM)^2
1395
1396 reflectivity_TE = (rTE * conj(rTE)) * (E1_perp * E1_perp) /
1397 (E1_perp * E1_perp + E1_parl * E1_parl);
1398 reflectivity_TM = (rTM * conj(rTM)) * (E1_parl * E1_parl) /
1399 (E1_perp * E1_perp + E1_parl * E1_parl);
1400 reflectivity = reflectivity_TE + reflectivity_TM;
1401
1402 do
1403 {
1404 if(G4UniformRand() * real(reflectivity) > real(reflectivity_TE))
1405 {
1406 f_iTE = -1;
1407 }
1408 else
1409 {
1410 f_iTE = 1;
1411 }
1412 if(G4UniformRand() * real(reflectivity) > real(reflectivity_TM))
1413 {
1414 f_iTM = -1;
1415 }
1416 else
1417 {
1418 f_iTM = 1;
1419 }
1420 // Loop checking, 13-Aug-2015, Peter Gumplinger
1421 } while(f_iTE < 0 && f_iTM < 0);
1422
1423 return real(reflectivity);
1424}
1425
1426//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1427void G4OpBoundaryProcess::CalculateReflectivity()
1428{
1429 G4double realRindex = fRealRIndexMPV->Value(fPhotonMomentum, idx_rrindex);
1430 G4double imaginaryRindex =
1431 fImagRIndexMPV->Value(fPhotonMomentum, idx_irindex);
1432
1433 // calculate FacetNormal
1434 if(fFinish == ground)
1435 {
1436 fFacetNormal = GetFacetNormal(fOldMomentum, fGlobalNormal);
1437 }
1438 else
1439 {
1440 fFacetNormal = fGlobalNormal;
1441 }
1442
1443 G4double cost1 = -fOldMomentum * fFacetNormal;
1444 if(std::abs(cost1) < 1.0 - fCarTolerance)
1445 {
1446 fSint1 = std::sqrt(1. - cost1 * cost1);
1447 }
1448 else
1449 {
1450 fSint1 = 0.0;
1451 }
1452
1453 G4ThreeVector A_trans, A_paral, E1pp, E1pl;
1454 G4double E1_perp, E1_parl;
1455
1456 if(fSint1 > 0.0)
1457 {
1458 A_trans = (fOldMomentum.cross(fFacetNormal)).unit();
1459 E1_perp = fOldPolarization * A_trans;
1460 E1pp = E1_perp * A_trans;
1461 E1pl = fOldPolarization - E1pp;
1462 E1_parl = E1pl.mag();
1463 }
1464 else
1465 {
1466 A_trans = fOldPolarization;
1467 // Here we Follow Jackson's conventions and we set the parallel
1468 // component = 1 in case of a ray perpendicular to the surface
1469 E1_perp = 0.0;
1470 E1_parl = 1.0;
1471 }
1472
1473 G4double incidentangle = GetIncidentAngle();
1474
1475 // calculate the reflectivity depending on incident angle,
1476 // polarization and complex refractive
1477 fReflectivity = GetReflectivity(E1_perp, E1_parl, incidentangle, realRindex,
1478 imaginaryRindex);
1479}
1480
1481//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1482G4bool G4OpBoundaryProcess::InvokeSD(const G4Step* pStep)
1483{
1484 G4Step aStep = *pStep;
1485 aStep.AddTotalEnergyDeposit(fPhotonMomentum);
1486
1488 if(sd != nullptr)
1489 return sd->Hit(&aStep);
1490 else
1491 return false;
1492}
1493
1494//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1496{
1497 fInvokeSD = flag;
1499}
1500
1501//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1503{
1504 verboseLevel = verbose;
1506}
1507
1508//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1509void G4OpBoundaryProcess::CoatedDielectricDielectric()
1510{
1511 G4MaterialPropertyVector* pp = nullptr;
1512
1514 if((pp = MPT->GetProperty(kRINDEX)))
1515 {
1516 fRindex2 = pp->Value(fPhotonMomentum, idx_rindex2);
1517 }
1518
1519 MPT = fOpticalSurface->GetMaterialPropertiesTable();
1520 if((pp = MPT->GetProperty(kCOATEDRINDEX)))
1521 {
1522 fCoatedRindex = pp->Value(fPhotonMomentum, idx_coatedrindex);
1523 }
1525 {
1526 fCoatedThickness = MPT->GetConstProperty(kCOATEDTHICKNESS);
1527 }
1529 {
1530 fCoatedFrustratedTransmission =
1532 }
1533
1534 G4double sintTL;
1535 G4double wavelength = h_Planck * c_light / fPhotonMomentum;
1536 G4double PdotN;
1537 G4double E1_perp, E1_parl;
1538 G4double s1, E2_perp, E2_parl, E2_total, transCoeff;
1539 G4double E2_abs, C_parl, C_perp;
1541 G4ThreeVector A_trans, A_paral, E1pp, E1pl;
1542 //G4bool Inside = false;
1543 //G4bool Swap = false;
1544 G4bool through = false;
1545 G4bool done = false;
1546
1547 do {
1548 if (through)
1549 {
1550 //Swap = !Swap;
1551 through = false;
1552 fGlobalNormal = -fGlobalNormal;
1553 G4SwapPtr(fMaterial1, fMaterial2);
1554 G4SwapObj(&fRindex1, &fRindex2);
1555 }
1556
1557 if(fFinish == polished)
1558 {
1559 fFacetNormal = fGlobalNormal;
1560 }
1561 else
1562 {
1563 fFacetNormal = GetFacetNormal(fOldMomentum, fGlobalNormal);
1564 }
1565
1566 PdotN = fOldMomentum * fFacetNormal;
1567 G4double cost1 = -PdotN;
1568 G4double sint2, cost2 = 0.;
1569
1570 if (std::abs(cost1) < 1.0 - fCarTolerance)
1571 {
1572 fSint1 = std::sqrt(1. - cost1 * cost1);
1573 sint2 = fSint1 * fRindex1 / fRindex2;
1574 sintTL = fSint1 * fRindex1 / fCoatedRindex;
1575 } else
1576 {
1577 fSint1 = 0.0;
1578 sint2 = 0.0;
1579 sintTL = 0.0;
1580 }
1581
1582 if (fSint1 > 0.0)
1583 {
1584 A_trans = fOldMomentum.cross(fFacetNormal);
1585 A_trans = A_trans.unit();
1586 E1_perp = fOldPolarization * A_trans;
1587 E1pp = E1_perp * A_trans;
1588 E1pl = fOldPolarization - E1pp;
1589 E1_parl = E1pl.mag();
1590 }
1591 else
1592 {
1593 A_trans = fOldPolarization;
1594 E1_perp = 0.0;
1595 E1_parl = 1.0;
1596 }
1597
1598 s1 = fRindex1 * cost1;
1599
1600 if (cost1 > 0.0)
1601 {
1602 cost2 = std::sqrt(1. - sint2 * sint2);
1603 }
1604 else
1605 {
1606 cost2 = -std::sqrt(1. - sint2 * sint2);
1607 }
1608
1609 transCoeff = 0.0;
1610
1611 if (sintTL >= 1.0)
1612 { // --> Angle > Angle Limit
1613 //Swap = false;
1614 }
1615 E2_perp = 2. * s1 * E1_perp / (fRindex1 * cost1 + fRindex2 * cost2);
1616 E2_parl = 2. * s1 * E1_parl / (fRindex2 * cost1 + fRindex1 * cost2);
1617 E2_total = E2_perp * E2_perp + E2_parl * E2_parl;
1618
1619 transCoeff = 1. - GetReflectivityThroughThinLayer(
1620 sintTL, E1_perp, E1_parl, wavelength, cost1, cost2);
1621 if (!G4BooleanRand(transCoeff))
1622 {
1623 if(verboseLevel > 2)
1624 G4cout << "Reflection from " << fMaterial1->GetName() << " to "
1625 << fMaterial2->GetName() << G4endl;
1626
1627 //Swap = false;
1628
1629 if (sintTL >= 1.0)
1630 {
1631 fStatus = TotalInternalReflection;
1632 }
1633 else
1634 {
1636 }
1637
1638 PdotN = fOldMomentum * fFacetNormal;
1639 fNewMomentum = fOldMomentum - (2. * PdotN) * fFacetNormal;
1640
1641 if (fSint1 > 0.0) { // incident ray oblique
1642
1643 E2_parl = fRindex2 * E2_parl / fRindex1 - E1_parl;
1644 E2_perp = E2_perp - E1_perp;
1645 E2_total = E2_perp * E2_perp + E2_parl * E2_parl;
1646 A_paral = fNewMomentum.cross(A_trans);
1647 A_paral = A_paral.unit();
1648 E2_abs = std::sqrt(E2_total);
1649 C_parl = E2_parl / E2_abs;
1650 C_perp = E2_perp / E2_abs;
1651
1652 fNewPolarization = C_parl * A_paral + C_perp * A_trans;
1653
1654 }
1655 else
1656 { // incident ray perpendicular
1657 if (fRindex2 > fRindex1)
1658 {
1659 fNewPolarization = -fOldPolarization;
1660 }
1661 else
1662 {
1663 fNewPolarization = fOldPolarization;
1664 }
1665 }
1666
1667 } else { // photon gets transmitted
1668 if (verboseLevel > 2)
1669 G4cout << "Transmission from " << fMaterial1->GetName() << " to "
1670 << fMaterial2->GetName() << G4endl;
1671
1672 //Inside = !Inside;
1673 through = true;
1674
1675 if (fEfficiency > 0.)
1676 {
1677 DoAbsorption();
1678 return;
1679 }
1680 else
1681 {
1682 if (sintTL >= 1.0)
1683 {
1685 }
1686 else
1687 {
1689 }
1690
1691 if (fSint1 > 0.0) { // incident ray oblique
1692
1693 alpha = cost1 - cost2 * (fRindex2 / fRindex1);
1694 fNewMomentum = fOldMomentum + alpha * fFacetNormal;
1695 fNewMomentum = fNewMomentum.unit();
1696 A_paral = fNewMomentum.cross(A_trans);
1697 A_paral = A_paral.unit();
1698 E2_abs = std::sqrt(E2_total);
1699 C_parl = E2_parl / E2_abs;
1700 C_perp = E2_perp / E2_abs;
1701
1702 fNewPolarization = C_parl * A_paral + C_perp * A_trans;
1703
1704 }
1705 else
1706 { // incident ray perpendicular
1707 fNewMomentum = fOldMomentum;
1708 fNewPolarization = fOldPolarization;
1709 }
1710 }
1711 }
1712
1713 fOldMomentum = fNewMomentum.unit();
1714 fOldPolarization = fNewPolarization.unit();
1715 if ((fStatus == CoatedDielectricFrustratedTransmission) ||
1716 (fStatus == CoatedDielectricRefraction))
1717 {
1718 done = (fNewMomentum * fGlobalNormal <= 0.0);
1719 }
1720 else
1721 {
1722 done = (fNewMomentum * fGlobalNormal >= -fCarTolerance);
1723 }
1724
1725 } while (!done);
1726}
1727
1728//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1729G4double G4OpBoundaryProcess::GetReflectivityThroughThinLayer(G4double sinTL,
1730 G4double E1_perp,
1731 G4double E1_parl,
1732 G4double wavelength, G4double cost1, G4double cost2) {
1733 G4complex Reflectivity, Reflectivity_TE, Reflectivity_TM;
1734 G4double gammaTL, costTL;
1735
1736 G4complex i(0, 1);
1737 G4complex rTM, rTE;
1738 G4complex r1toTL, rTLto2;
1739 G4double k0 = 2 * pi / wavelength;
1740
1741 // Angle > Angle limit
1742 if (sinTL >= 1.0) {
1743 if (fCoatedFrustratedTransmission) { //Frustrated transmission
1744
1745 if (cost1 > 0.0)
1746 {
1747 gammaTL = std::sqrt(fRindex1 * fRindex1 * fSint1 * fSint1 -
1748 fCoatedRindex * fCoatedRindex);
1749 }
1750 else
1751 {
1752 gammaTL = -std::sqrt(fRindex1 * fRindex1 * fSint1 * fSint1 -
1753 fCoatedRindex * fCoatedRindex);
1754 }
1755
1756 // TE
1757 r1toTL = (fRindex1 * cost1 - i * gammaTL) / (fRindex1 * cost1 + i * gammaTL);
1758 rTLto2 = (i * gammaTL - fRindex2 * cost2) / (i * gammaTL + fRindex2 * cost2);
1759 if (cost1 != 0.0)
1760 {
1761 rTE = (r1toTL + rTLto2 * std::exp(-2 * k0 * fCoatedThickness * gammaTL)) /
1762 (1.0 + r1toTL * rTLto2 * std::exp(-2 * k0 * fCoatedThickness * gammaTL));
1763 }
1764 // TM
1765 r1toTL = (fRindex1 * i * gammaTL - fCoatedRindex * fCoatedRindex * cost1) /
1766 (fRindex1 * i * gammaTL + fCoatedRindex * fCoatedRindex * cost1);
1767 rTLto2 = (fCoatedRindex * fCoatedRindex * cost2 - fRindex2 * i * gammaTL) /
1768 (fCoatedRindex * fCoatedRindex * cost2 + fRindex2 * i * gammaTL);
1769 if (cost1 != 0.0)
1770 {
1771 rTM = (r1toTL + rTLto2 * std::exp(-2 * k0 * fCoatedThickness * gammaTL)) /
1772 (1.0 + r1toTL * rTLto2 * std::exp(-2 * k0 * fCoatedThickness * gammaTL));
1773 }
1774 }
1775 else
1776 { //Total reflection
1777 return(1.);
1778 }
1779 }
1780
1781 // Angle <= Angle limit
1782 else //if (sinTL < 1.0)
1783 {
1784 if (cost1 > 0.0)
1785 {
1786 costTL = std::sqrt(1. - sinTL * sinTL);
1787 }
1788 else
1789 {
1790 costTL = -std::sqrt(1. - sinTL * sinTL);
1791 }
1792 // TE
1793 r1toTL = (fRindex1 * cost1 - fCoatedRindex * costTL) / (fRindex1 * cost1 + fCoatedRindex * costTL);
1794 rTLto2 = (fCoatedRindex * costTL - fRindex2 * cost2) / (fCoatedRindex * costTL + fRindex2 * cost2);
1795 if (cost1 != 0.0)
1796 {
1797 rTE = (r1toTL + rTLto2 * std::exp(2.0 * i * k0 * fCoatedRindex * fCoatedThickness * costTL)) /
1798 (1.0 + r1toTL * rTLto2 * std::exp(2.0 * i * k0 * fCoatedRindex * fCoatedThickness * costTL));
1799 }
1800 // TM
1801 r1toTL = (fRindex1 * costTL - fCoatedRindex * cost1) / (fRindex1 * costTL + fCoatedRindex * cost1);
1802 rTLto2 = (fCoatedRindex * cost2 - fRindex2 * costTL) / (fCoatedRindex * cost2 + fRindex2 * costTL);
1803 if (cost1 != 0.0)
1804 {
1805 rTM = (r1toTL + rTLto2 * std::exp(2.0 * i * k0 * fCoatedRindex * fCoatedThickness * costTL)) /
1806 (1.0 + r1toTL * rTLto2 * std::exp(2.0 * i * k0 * fCoatedRindex * fCoatedThickness * costTL));
1807 }
1808 }
1809
1810 Reflectivity_TE = (rTE * conj(rTE)) * (E1_perp * E1_perp) / (E1_perp * E1_perp + E1_parl * E1_parl);
1811 Reflectivity_TM = (rTM * conj(rTM)) * (E1_parl * E1_parl) / (E1_perp * E1_perp + E1_parl * E1_parl);
1812 Reflectivity = Reflectivity_TE + Reflectivity_TM;
1813
1814 return real(Reflectivity);
1815}
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
@ FatalException
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
@ Forced
@ kBACKSCATTERCONSTANT
@ kSPECULARLOBECONSTANT
@ kSPECULARSPIKECONSTANT
@ kCOATEDFRUSTRATEDTRANSMISSION
G4OpBoundaryProcessStatus
@ PolishedTiOAirReflection
@ GroundTeflonAirReflection
@ EtchedVM2000AirReflection
@ Undefined
@ NotAtBoundary
@ Transmission
@ CoatedDielectricReflection
@ EtchedVM2000GlueReflection
@ GroundLumirrorGlueReflection
@ LobeReflection
@ GroundTiOAirReflection
@ GroundTyvekAirReflection
@ FresnelReflection
@ NoRINDEX
@ GroundAirReflection
@ EtchedAirReflection
@ PolishedVM2000GlueReflection
@ PolishedAirReflection
@ PolishedTeflonAirReflection
@ SameMaterial
@ EtchedTyvekAirReflection
@ SpikeReflection
@ EtchedLumirrorGlueReflection
@ GroundVM2000AirReflection
@ PolishedTyvekAirReflection
@ PolishedVM2000AirReflection
@ EtchedTiOAirReflection
@ EtchedTeflonAirReflection
@ GroundVM2000GlueReflection
@ BackScattering
@ PolishedLumirrorGlueReflection
@ Absorption
@ TotalInternalReflection
@ StepTooSmall
@ CoatedDielectricFrustratedTransmission
@ PolishedLumirrorAirReflection
@ EtchedLumirrorAirReflection
@ GroundLumirrorAirReflection
@ LambertianReflection
@ CoatedDielectricRefraction
@ FresnelRefraction
@ Detection
@ fOpBoundary
@ unified
@ DAVIS
@ glisur
@ groundfrontpainted
@ polishedbackpainted
@ polished
@ ground
@ polishedfrontpainted
@ groundbackpainted
G4ProcessType
@ fGeomBoundary
Definition: G4StepStatus.hh:43
G4SurfaceType
@ dielectric_metal
@ dielectric_LUT
@ coated
@ dielectric_dielectric
@ dielectric_LUTDAVIS
@ dielectric_dichroic
@ 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
#define G4endl
Definition: G4ios.hh:57
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 &)
Definition: ThreeVector.cc:33
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:24
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)
G4SurfaceProperty * GetSurfaceProperty() const
G4bool ConstPropertyExists(const G4String &key) const
G4double GetConstProperty(const G4String &key) const
G4MaterialPropertyVector * GetProperty(const char *key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:251
const G4String & GetName() const
Definition: G4Material.hh:172
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
void SetBoundaryInvokeSD(G4bool)
void SetBoundaryVerboseLevel(G4int)
G4int GetBoundaryVerboseLevel() const
static G4OpticalParameters * Instance()
G4int GetLUTbins() const
G4double GetAngularDistributionValueLUT(G4int)
G4OpticalSurfaceModel GetModel() const
G4double GetSigmaAlpha() const
G4OpticalSurfaceFinish GetFinish() const
G4int GetThetaIndexMax() const
G4double GetPolish() const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
G4Physics2DVector * GetDichroicVector()
G4int GetPhiIndexMax() const
G4double GetReflectivityLUTValue(G4int)
G4double GetAngularDistributionValue(G4int, G4int, G4int)
static const G4Step * GetHyperStep()
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
void Initialize(const G4Track &) override
void ProposeVelocity(G4double finalVelocity)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double Value(G4double x, G4double y, std::size_t &lastidx, std::size_t &lastidy) const
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
Definition: G4Step.hh:62
void AddTotalEnergyDeposit(G4double value)
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
const G4SurfaceType & GetType() 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 &)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4LogicalVolume * GetMotherLogical() const
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331
G4int verboseLevel
Definition: G4VProcess.hh:360
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
G4bool Hit(G4Step *aStep)
const G4double pi
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