Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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//
27// Optical Photon Boundary Process Class Implementation
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//
75
77
78#include "G4ios.hh"
82#include "G4OpProcessSubType.hh"
86#include "G4SystemOfUnits.hh"
89
90//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92 G4ProcessType type)
93 : G4VDiscreteProcess(processName, type)
94{
95 Initialise();
96
97 if(verboseLevel > 0)
98 {
99 G4cout << GetProcessName() << " is created " << G4endl;
100 }
102
104 fModel = glisur;
106 fReflectivity = 1.;
107 fEfficiency = 0.;
108 fTransmittance = 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
128//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
130
131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
133{
134 Initialise();
135}
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139{
143}
144
145//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147 const G4Step& aStep)
148{
152
153 // Get hyperStep from G4ParallelWorldProcess
154 // NOTE: PostSetpDoIt of this process to be invoked after
155 // G4ParallelWorldProcess!
156 const G4Step* pStep = &aStep;
158 if(hStep != nullptr)
159 pStep = hStep;
160
162 {
165 }
166 else
167 {
169 if(verboseLevel > 1)
171 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
172 }
173
176
177 if(verboseLevel > 1)
178 {
179 G4cout << " Photon at Boundary! " << G4endl;
180 if(thePrePV != nullptr)
181 G4cout << " thePrePV: " << thePrePV->GetName() << G4endl;
182 if(thePostPV != nullptr)
183 G4cout << " thePostPV: " << thePostPV->GetName() << G4endl;
184 }
185
186 if(aTrack.GetStepLength() <= fCarTolerance)
187 {
189 if(verboseLevel > 1)
191
192 G4MaterialPropertyVector* groupvel = nullptr;
194 if(aMPT != nullptr)
195 {
196 groupvel = aMPT->GetProperty(kGROUPVEL);
197 }
198 if(groupvel != nullptr)
199 {
202 }
203 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
204 }
205
206 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
207
208 fPhotonMomentum = aParticle->GetTotalMomentum();
209 fOldMomentum = aParticle->GetMomentumDirection();
210 fOldPolarization = aParticle->GetPolarization();
211
212 if(verboseLevel > 1)
213 {
214 G4cout << " Old Momentum Direction: " << fOldMomentum << G4endl
215 << " Old Polarization: " << fOldPolarization << G4endl;
216 }
217
218 G4ThreeVector theGlobalPoint = pStep->GetPostStepPoint()->GetPosition();
219 G4bool valid;
220
221 // ID of Navigator which limits step
225 fGlobalNormal = (iNav[hNavId])->GetGlobalExitNormal(theGlobalPoint, &valid);
226
227 if(valid)
228 {
230 }
231 else
232 {
234 ed << " G4OpBoundaryProcess/PostStepDoIt(): "
235 << " The Navigator reports that it returned an invalid normal" << G4endl;
237 "G4OpBoundaryProcess::PostStepDoIt", "OpBoun01", EventMustBeAborted, ed,
238 "Invalid Surface Normal - Geometry must return valid surface normal");
239 }
240
241 if(fOldMomentum * fGlobalNormal > 0.0)
242 {
243#ifdef G4OPTICAL_DEBUG
245 ed << " G4OpBoundaryProcess/PostStepDoIt(): fGlobalNormal points in a "
246 "wrong direction. "
247 << G4endl
248 << " The momentum of the photon arriving at interface (oldMomentum)"
249 << " must exit the volume cross in the step. " << G4endl
250 << " So it MUST have dot < 0 with the normal that Exits the new "
251 "volume (globalNormal)."
252 << G4endl << " >> The dot product of oldMomentum and global Normal is "
254 << " Old Momentum (during step) = " << fOldMomentum << G4endl
255 << " Global Normal (Exiting New Vol) = " << fGlobalNormal << G4endl
256 << G4endl;
257 G4Exception("G4OpBoundaryProcess::PostStepDoIt", "OpBoun02",
258 EventMustBeAborted, // Or JustWarning to see if it happens
259 // repeatedly on one ray
260 ed,
261 "Invalid Surface Normal - Geometry must return valid surface "
262 "normal pointing in the right direction");
263#else
265#endif
266 }
267
268 G4MaterialPropertyVector* rIndexMPV = nullptr;
270 if(MPT != nullptr)
271 {
272 rIndexMPV = MPT->GetProperty(kRINDEX);
273 }
274 if(rIndexMPV != nullptr)
275 {
277 }
278 else
279 {
281 if(verboseLevel > 1)
285 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
286 }
287
288 fReflectivity = 1.;
289 fEfficiency = 0.;
290 fTransmittance = 0.;
292 fModel = glisur;
295
296 rIndexMPV = nullptr;
297 fOpticalSurface = nullptr;
298
299 G4LogicalSurface* surface =
300 G4LogicalBorderSurface::GetSurface(thePrePV, thePostPV);
301 if(surface == nullptr)
302 {
303 if(thePostPV->GetMotherLogical() == thePrePV->GetLogicalVolume())
304 {
306 if(surface == nullptr)
307 {
308 surface =
310 }
311 }
312 else
313 {
315 if(surface == nullptr)
316 {
317 surface =
319 }
320 }
321 }
322
323 if(surface != nullptr)
324 {
326 dynamic_cast<G4OpticalSurface*>(surface->GetSurfaceProperty());
327 }
328 if(fOpticalSurface != nullptr)
329 {
330 type = fOpticalSurface->GetType();
333
336 if(sMPT != nullptr)
337 {
339 {
340 rIndexMPV = sMPT->GetProperty(kRINDEX);
341 if(rIndexMPV != nullptr)
342 {
344 }
345 else
346 {
348 if(verboseLevel > 1)
352 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
353 }
354 }
355
358 f_iTE = f_iTM = 1;
359
361 if((pp = sMPT->GetProperty(kREFLECTIVITY)))
362 {
364 }
366 {
368 }
369
370 if((pp = sMPT->GetProperty(kEFFICIENCY)))
371 {
373 }
374 if((pp = sMPT->GetProperty(kTRANSMITTANCE)))
375 {
377 }
379 {
381 }
382
383 if(fModel == unified)
384 {
386 ? pp->Value(fPhotonMomentum, idx_lobe)
387 : 0.;
389 ? pp->Value(fPhotonMomentum, idx_spike)
390 : 0.;
392 ? pp->Value(fPhotonMomentum, idx_back)
393 : 0.;
394 }
395 } // end of if(sMPT)
397 {
400 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
401 }
402 } // end of if(fOpticalSurface)
403
404 // DIELECTRIC-DIELECTRIC
405 if(type == dielectric_dielectric)
406 {
407 if(fFinish == polished || fFinish == ground)
408 {
410 {
412 if(verboseLevel > 1)
414 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
415 }
417 rIndexMPV = nullptr;
418 if(MPT != nullptr)
419 {
420 rIndexMPV = MPT->GetProperty(kRINDEX);
421 }
422 if(rIndexMPV != nullptr)
423 {
425 }
426 else
427 {
429 if(verboseLevel > 1)
433 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
434 }
435 }
437 {
439 }
440 else
441 {
442 G4double rand = G4UniformRand();
443 if(rand > fReflectivity + fTransmittance)
444 {
445 DoAbsorption();
446 }
447 else if(rand > fReflectivity)
448 {
452 }
453 else
454 {
456 {
457 DoReflection();
458 }
459 else if(fFinish == groundfrontpainted)
460 {
462 DoReflection();
463 }
464 else
465 {
467 }
468 }
469 }
470 }
471 else if(type == dielectric_metal)
472 {
474 }
475 else if(type == dielectric_LUT)
476 {
478 }
479 else if(type == dielectric_LUTDAVIS)
480 {
482 }
483 else if(type == dielectric_dichroic)
484 {
486 }
487 else
488 {
490 ed << " PostStepDoIt(): Illegal boundary type." << G4endl;
491 G4Exception("G4OpBoundaryProcess", "OpBoun04", JustWarning, ed);
492 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
493 }
494
497
498 if(verboseLevel > 1)
499 {
500 G4cout << " New Momentum Direction: " << fNewMomentum << G4endl
501 << " New Polarization: " << fNewPolarization << G4endl;
503 }
504
507
509 {
510 // not all surface types check that fMaterial2 has an MPT
512 G4MaterialPropertyVector* groupvel = nullptr;
513 if(aMPT != nullptr)
514 {
515 groupvel = aMPT->GetProperty(kGROUPVEL);
516 }
517 if(groupvel != nullptr)
518 {
521 }
522 }
523
524 if(fStatus == Detection && fInvokeSD)
525 InvokeSD(pStep);
526 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
527}
528
529//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
531{
532 G4cout << " *** ";
533 if(fStatus == Undefined)
534 G4cout << "Undefined";
535 else if(fStatus == Transmission)
536 G4cout << "Transmission";
537 else if(fStatus == FresnelRefraction)
538 G4cout << "FresnelRefraction";
539 else if(fStatus == FresnelReflection)
540 G4cout << "FresnelReflection";
542 G4cout << "TotalInternalReflection";
543 else if(fStatus == LambertianReflection)
544 G4cout << "LambertianReflection";
545 else if(fStatus == LobeReflection)
546 G4cout << "LobeReflection";
547 else if(fStatus == SpikeReflection)
548 G4cout << "SpikeReflection";
549 else if(fStatus == BackScattering)
550 G4cout << "BackScattering";
552 G4cout << "PolishedLumirrorAirReflection";
554 G4cout << "PolishedLumirrorGlueReflection";
556 G4cout << "PolishedAirReflection";
558 G4cout << "PolishedTeflonAirReflection";
560 G4cout << "PolishedTiOAirReflection";
562 G4cout << "PolishedTyvekAirReflection";
564 G4cout << "PolishedVM2000AirReflection";
566 G4cout << "PolishedVM2000GlueReflection";
568 G4cout << "EtchedLumirrorAirReflection";
570 G4cout << "EtchedLumirrorGlueReflection";
571 else if(fStatus == EtchedAirReflection)
572 G4cout << "EtchedAirReflection";
574 G4cout << "EtchedTeflonAirReflection";
576 G4cout << "EtchedTiOAirReflection";
578 G4cout << "EtchedTyvekAirReflection";
580 G4cout << "EtchedVM2000AirReflection";
582 G4cout << "EtchedVM2000GlueReflection";
584 G4cout << "GroundLumirrorAirReflection";
586 G4cout << "GroundLumirrorGlueReflection";
587 else if(fStatus == GroundAirReflection)
588 G4cout << "GroundAirReflection";
590 G4cout << "GroundTeflonAirReflection";
592 G4cout << "GroundTiOAirReflection";
594 G4cout << "GroundTyvekAirReflection";
596 G4cout << "GroundVM2000AirReflection";
598 G4cout << "GroundVM2000GlueReflection";
599 else if(fStatus == Absorption)
600 G4cout << "Absorption";
601 else if(fStatus == Detection)
602 G4cout << "Detection";
603 else if(fStatus == NotAtBoundary)
604 G4cout << "NotAtBoundary";
605 else if(fStatus == SameMaterial)
606 G4cout << "SameMaterial";
607 else if(fStatus == StepTooSmall)
608 G4cout << "StepTooSmall";
609 else if(fStatus == NoRINDEX)
610 G4cout << "NoRINDEX";
611 else if(fStatus == Dichroic)
612 G4cout << "Dichroic Transmission";
613 G4cout << " ***" << G4endl;
614}
615
616//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
618 const G4ThreeVector& momentum, const G4ThreeVector& normal) const
619{
620 G4ThreeVector facetNormal;
621 if(fModel == unified || fModel == LUT || fModel == DAVIS)
622 {
623 /* This function codes alpha to a random value taken from the
624 distribution p(alpha) = g(alpha; 0, sigma_alpha)*std::sin(alpha),
625 for alpha > 0 and alpha < 90, where g(alpha; 0, sigma_alpha) is a
626 gaussian distribution with mean 0 and standard deviation sigma_alpha. */
627
628 G4double sigma_alpha = 0.0;
630 sigma_alpha = fOpticalSurface->GetSigmaAlpha();
631 if(sigma_alpha == 0.0)
632 {
633 return normal;
634 }
635
636 G4double f_max = std::min(1.0, 4. * sigma_alpha);
637 G4double alpha, phi, sinAlpha;
638
639 do
640 { // Loop checking, 13-Aug-2015, Peter Gumplinger
641 do
642 { // Loop checking, 13-Aug-2015, Peter Gumplinger
643 alpha = G4RandGauss::shoot(0.0, sigma_alpha);
644 sinAlpha = std::sin(alpha);
645 } while(G4UniformRand() * f_max > sinAlpha || alpha >= halfpi);
646
647 phi = G4UniformRand() * twopi;
648 facetNormal.set(sinAlpha * std::cos(phi), sinAlpha * std::sin(phi),
649 std::cos(alpha));
650 facetNormal.rotateUz(normal);
651 } while(momentum * facetNormal >= 0.0);
652 }
653 else
654 {
655 G4double polish = 1.0;
657 polish = fOpticalSurface->GetPolish();
658 if(polish < 1.0)
659 {
660 do
661 { // Loop checking, 13-Aug-2015, Peter Gumplinger
662 G4ThreeVector smear;
663 do
664 { // Loop checking, 13-Aug-2015, Peter Gumplinger
665 smear.setX(2. * G4UniformRand() - 1.);
666 smear.setY(2. * G4UniformRand() - 1.);
667 smear.setZ(2. * G4UniformRand() - 1.);
668 } while(smear.mag2() > 1.0);
669 facetNormal = normal + (1. - polish) * smear;
670 } while(momentum * facetNormal >= 0.0);
671 facetNormal = facetNormal.unit();
672 }
673 else
674 {
675 facetNormal = normal;
676 }
677 }
678 return facetNormal;
679}
680
681//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
683{
684 G4int n = 0;
685 G4double rand;
686 G4ThreeVector A_trans;
687
688 do
689 {
690 ++n;
691 rand = G4UniformRand();
692 if(rand > fReflectivity && n == 1)
693 {
694 if(rand > fReflectivity + fTransmittance)
695 {
696 DoAbsorption();
697 }
698 else
699 {
703 }
704 break;
705 }
706 else
707 {
709 {
710 if(n > 1)
711 {
714 {
715 DoAbsorption();
716 break;
717 }
718 }
719 }
720 if(fModel == glisur || fFinish == polished)
721 {
722 DoReflection();
723 }
724 else
725 {
726 if(n == 1)
729 {
730 DoReflection();
731 }
732 else if(fStatus == BackScattering)
733 {
736 }
737 else
738 {
740 {
742 {
744 }
745 //else
746 // case of complex rindex needs to be implemented
747 }
750
751 if(f_iTE > 0 && f_iTM > 0)
752 {
756 }
757 else if(f_iTE > 0)
758 {
759 A_trans = (fSint1 > 0.0) ? fOldMomentum.cross(fFacetNormal).unit()
761 fNewPolarization = -A_trans;
762 }
763 else if(f_iTM > 0)
764 {
766 -fNewMomentum.cross(A_trans).unit(); // = -A_paral
767 }
768 }
769 }
772 }
773 // Loop checking, 13-Aug-2015, Peter Gumplinger
774 } while(fNewMomentum * fGlobalNormal < 0.0);
775}
776
777//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
779{
780 G4int thetaIndex, phiIndex;
781 G4double angularDistVal, thetaRad, phiRad;
782 G4ThreeVector perpVectorTheta, perpVectorPhi;
783
786
787 G4int thetaIndexMax = fOpticalSurface->GetThetaIndexMax();
788 G4int phiIndexMax = fOpticalSurface->GetPhiIndexMax();
789
790 G4double rand;
791
792 do
793 {
794 rand = G4UniformRand();
795 if(rand > fReflectivity)
796 {
797 if(rand > fReflectivity + fTransmittance)
798 {
799 DoAbsorption();
800 }
801 else
802 {
806 }
807 break;
808 }
809 else
810 {
811 // Calculate Angle between Normal and Photon Momentum
812 G4double anglePhotonToNormal = fOldMomentum.angle(-fGlobalNormal);
813 // Round to closest integer: LBNL model array has 91 values
814 G4int angleIncident = std::lrint(anglePhotonToNormal / CLHEP::deg);
815
816 // Take random angles THETA and PHI,
817 // and see if below Probability - if not - Redo
818 do
819 {
820 thetaIndex = G4RandFlat::shootInt(thetaIndexMax - 1);
821 phiIndex = G4RandFlat::shootInt(phiIndexMax - 1);
822 // Find probability with the new indeces from LUT
824 angleIncident, thetaIndex, phiIndex);
825 // Loop checking, 13-Aug-2015, Peter Gumplinger
826 } while(!G4BooleanRand(angularDistVal));
827
828 thetaRad = G4double(-90 + 4 * thetaIndex) * pi / 180.;
829 phiRad = G4double(-90 + 5 * phiIndex) * pi / 180.;
830 // Rotate Photon Momentum in Theta, then in Phi
832
833 perpVectorTheta = fNewMomentum.cross(fGlobalNormal);
834 if(perpVectorTheta.mag() < fCarTolerance)
835 {
836 perpVectorTheta = fNewMomentum.orthogonal();
837 }
839 fNewMomentum.rotate(anglePhotonToNormal - thetaRad, perpVectorTheta);
840 perpVectorPhi = perpVectorTheta.cross(fNewMomentum);
841 fNewMomentum = fNewMomentum.rotate(-phiRad, perpVectorPhi);
842
843 // Rotate Polarization too:
847 }
848 // Loop checking, 13-Aug-2015, Peter Gumplinger
849 } while(fNewMomentum * fGlobalNormal <= 0.0);
850}
851
852//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
854{
855 G4int angindex, random, angleIncident;
856 G4double reflectivityValue, elevation, azimuth;
857 G4double anglePhotonToNormal;
858
859 G4int lutbin = fOpticalSurface->GetLUTbins();
860 G4double rand = G4UniformRand();
861
862 G4double sinEl;
863 G4ThreeVector u, vNorm, w;
864
865 do
866 {
867 anglePhotonToNormal = fOldMomentum.angle(-fGlobalNormal);
868
869 // Davis model has 90 reflection bins: round down
870 // don't allow angleIncident to be 90 for anglePhotonToNormal close to 90
871 angleIncident = std::min(static_cast<G4int>(
872 std::floor(anglePhotonToNormal / CLHEP::deg)), 89);
873 reflectivityValue = fOpticalSurface->GetReflectivityLUTValue(angleIncident);
874
875 if(rand > reflectivityValue)
876 {
877 if(fEfficiency > 0.)
878 {
879 DoAbsorption();
880 break;
881 }
882 else
883 {
885
886 if(angleIncident <= 0.01)
887 {
889 break;
890 }
891
892 do
893 {
894 random = G4RandFlat::shootInt(1, lutbin + 1);
895 angindex =
896 (((random * 2) - 1)) + angleIncident * lutbin * 2 + 3640000;
897
898 azimuth =
900 elevation = fOpticalSurface->GetAngularDistributionValueLUT(angindex);
901 } while(elevation == 0. && azimuth == 0.);
902
903 sinEl = std::sin(elevation);
904 vNorm = (fGlobalNormal.cross(fOldMomentum)).unit();
905 u = vNorm.cross(fGlobalNormal) * (sinEl * std::cos(azimuth));
906 vNorm *= (sinEl * std::sin(azimuth));
907 // fGlobalNormal shouldn't be modified here
908 w = (fGlobalNormal *= std::cos(elevation));
909 fNewMomentum = u + vNorm + w;
910
911 // Rotate Polarization too:
915 }
916 }
917 else
918 {
920
921 if(angleIncident == 0)
922 {
924 break;
925 }
926
927 do
928 {
929 random = G4RandFlat::shootInt(1, lutbin + 1);
930 angindex = (((random * 2) - 1)) + (angleIncident - 1) * lutbin * 2;
931
932 azimuth = fOpticalSurface->GetAngularDistributionValueLUT(angindex - 1);
933 elevation = fOpticalSurface->GetAngularDistributionValueLUT(angindex);
934 } while(elevation == 0. && azimuth == 0.);
935
936 sinEl = std::sin(elevation);
937 vNorm = (fGlobalNormal.cross(fOldMomentum)).unit();
938 u = vNorm.cross(fGlobalNormal) * (sinEl * std::cos(azimuth));
939 vNorm *= (sinEl * std::sin(azimuth));
940 // fGlobalNormal shouldn't be modified here
941 w = (fGlobalNormal *= std::cos(elevation));
942
943 fNewMomentum = u + vNorm + w;
944
945 // Rotate Polarization too: (needs revision)
947 }
948 } while(fNewMomentum * fGlobalNormal <= 0.0);
949}
950
951//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
953{
954 // Calculate Angle between Normal and Photon Momentum
955 G4double anglePhotonToNormal = fOldMomentum.angle(-fGlobalNormal);
956
957 // Round it to closest integer
958 G4double angleIncident = std::floor(180. / pi * anglePhotonToNormal + 0.5);
959
960 if(!fDichroicVector)
961 {
964 }
965
967 {
968 G4double wavelength = h_Planck * c_light / fPhotonMomentum;
969 fTransmittance = fDichroicVector->Value(wavelength / nm, angleIncident,
971 perCent;
972 // G4cout << "wavelength: " << std::floor(wavelength/nm)
973 // << "nm" << G4endl;
974 // G4cout << "Incident angle: " << angleIncident << "deg" << G4endl;
975 // G4cout << "Transmittance: "
976 // << std::floor(fTransmittance/perCent) << "%" << G4endl;
977 }
978 else
979 {
981 ed << " G4OpBoundaryProcess/DielectricDichroic(): "
982 << " The dichroic surface has no G4Physics2DVector" << G4endl;
983 G4Exception("G4OpBoundaryProcess::DielectricDichroic", "OpBoun03",
984 FatalException, ed,
985 "A dichroic surface must have an associated G4Physics2DVector");
986 }
987
989 { // Not transmitted, so reflect
990 if(fModel == glisur || fFinish == polished)
991 {
992 DoReflection();
993 }
994 else
995 {
998 {
999 DoReflection();
1000 }
1001 else if(fStatus == BackScattering)
1002 {
1005 }
1006 else
1007 {
1008 G4double PdotN, EdotN;
1009 do
1010 {
1011 if(fStatus == LobeReflection)
1012 {
1014 }
1015 PdotN = fOldMomentum * fFacetNormal;
1016 fNewMomentum = fOldMomentum - (2. * PdotN) * fFacetNormal;
1017 // Loop checking, 13-Aug-2015, Peter Gumplinger
1018 } while(fNewMomentum * fGlobalNormal <= 0.0);
1019
1022 }
1023 }
1024 }
1025 else
1026 {
1027 fStatus = Dichroic;
1030 }
1031}
1032
1033//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1035{
1036 G4bool inside = false;
1037 G4bool swap = false;
1038
1039 if(fFinish == polished)
1040 {
1042 }
1043 else
1044 {
1046 }
1048 G4double cost2 = 0.;
1049 G4double sint2 = 0.;
1050
1051 G4bool surfaceRoughnessCriterionPass = true;
1052 if(fSurfaceRoughness != 0. && fRindex1 > fRindex2)
1053 {
1054 G4double wavelength = h_Planck * c_light / fPhotonMomentum;
1055 G4double surfaceRoughnessCriterion = std::exp(-std::pow(
1056 (4. * pi * fSurfaceRoughness * fRindex1 * cost1 / wavelength), 2));
1057 surfaceRoughnessCriterionPass = G4BooleanRand(surfaceRoughnessCriterion);
1058 }
1059
1060leap:
1061
1062 G4bool through = false;
1063 G4bool done = false;
1064
1065 G4ThreeVector A_trans, A_paral, E1pp, E1pl;
1066 G4double E1_perp, E1_parl;
1067 G4double s1, s2, E2_perp, E2_parl, E2_total, transCoeff;
1068 G4double E2_abs, C_parl, C_perp;
1070
1071 do
1072 {
1073 if(through)
1074 {
1075 swap = !swap;
1076 through = false;
1080 }
1081
1082 if(fFinish == polished)
1083 {
1085 }
1086 else
1087 {
1089 }
1090
1091 cost1 = -fOldMomentum * fFacetNormal;
1092 if(std::abs(cost1) < 1.0 - fCarTolerance)
1093 {
1094 fSint1 = std::sqrt(1. - cost1 * cost1);
1095 sint2 = fSint1 * fRindex1 / fRindex2; // *** Snell's Law ***
1096 // this isn't a sine as we might expect from the name; can be > 1
1097 }
1098 else
1099 {
1100 fSint1 = 0.0;
1101 sint2 = 0.0;
1102 }
1103
1104 // TOTAL INTERNAL REFLECTION
1105 if(sint2 >= 1.0)
1106 {
1107 swap = false;
1108
1110 if(!surfaceRoughnessCriterionPass)
1112 if(fModel == unified && fFinish != polished)
1115 {
1116 DoReflection();
1117 }
1118 else if(fStatus == BackScattering)
1119 {
1122 }
1123 else
1124 {
1125 fNewMomentum =
1129 }
1130 }
1131 // NOT TIR
1132 else if(sint2 < 1.0)
1133 {
1134 // Calculate amplitude for transmission (Q = P x N)
1135 if(cost1 > 0.0)
1136 {
1137 cost2 = std::sqrt(1. - sint2 * sint2);
1138 }
1139 else
1140 {
1141 cost2 = -std::sqrt(1. - sint2 * sint2);
1142 }
1143
1144 if(fSint1 > 0.0)
1145 {
1146 A_trans = (fOldMomentum.cross(fFacetNormal)).unit();
1147 E1_perp = fOldPolarization * A_trans;
1148 E1pp = E1_perp * A_trans;
1149 E1pl = fOldPolarization - E1pp;
1150 E1_parl = E1pl.mag();
1151 }
1152 else
1153 {
1154 A_trans = fOldPolarization;
1155 // Here we Follow Jackson's conventions and set the parallel
1156 // component = 1 in case of a ray perpendicular to the surface
1157 E1_perp = 0.0;
1158 E1_parl = 1.0;
1159 }
1160
1161 s1 = fRindex1 * cost1;
1162 E2_perp = 2. * s1 * E1_perp / (fRindex1 * cost1 + fRindex2 * cost2);
1163 E2_parl = 2. * s1 * E1_parl / (fRindex2 * cost1 + fRindex1 * cost2);
1164 E2_total = E2_perp * E2_perp + E2_parl * E2_parl;
1165 s2 = fRindex2 * cost2 * E2_total;
1166
1167 if(fTransmittance > 0.)
1168 transCoeff = fTransmittance;
1169 else if(cost1 != 0.0)
1170 transCoeff = s2 / s1;
1171 else
1172 transCoeff = 0.0;
1173
1174 // NOT TIR: REFLECTION
1175 if(!G4BooleanRand(transCoeff))
1176 {
1177 swap = false;
1179
1180 if(!surfaceRoughnessCriterionPass)
1182 if(fModel == unified && fFinish != polished)
1185 {
1186 DoReflection();
1187 }
1188 else if(fStatus == BackScattering)
1189 {
1192 }
1193 else
1194 {
1195 fNewMomentum =
1197 if(fSint1 > 0.0)
1198 { // incident ray oblique
1199 E2_parl = fRindex2 * E2_parl / fRindex1 - E1_parl;
1200 E2_perp = E2_perp - E1_perp;
1201 E2_total = E2_perp * E2_perp + E2_parl * E2_parl;
1202 A_paral = (fNewMomentum.cross(A_trans)).unit();
1203 E2_abs = std::sqrt(E2_total);
1204 C_parl = E2_parl / E2_abs;
1205 C_perp = E2_perp / E2_abs;
1206
1207 fNewPolarization = C_parl * A_paral + C_perp * A_trans;
1208 }
1209 else
1210 { // incident ray perpendicular
1211 if(fRindex2 > fRindex1)
1212 {
1214 }
1215 else
1216 {
1218 }
1219 }
1220 }
1221 }
1222 // NOT TIR: TRANSMISSION
1223 else
1224 {
1225 inside = !inside;
1226 through = true;
1228
1229 if(fSint1 > 0.0)
1230 { // incident ray oblique
1231 alpha = cost1 - cost2 * (fRindex2 / fRindex1);
1233 A_paral = (fNewMomentum.cross(A_trans)).unit();
1234 E2_abs = std::sqrt(E2_total);
1235 C_parl = E2_parl / E2_abs;
1236 C_perp = E2_perp / E2_abs;
1237
1238 fNewPolarization = C_parl * A_paral + C_perp * A_trans;
1239 }
1240 else
1241 { // incident ray perpendicular
1244 }
1245 }
1246 }
1247
1250
1252 {
1253 done = (fNewMomentum * fGlobalNormal <= 0.0);
1254 }
1255 else
1256 {
1258 }
1259 // Loop checking, 13-Aug-2015, Peter Gumplinger
1260 } while(!done);
1261
1262 if(inside && !swap)
1263 {
1265 {
1266 G4double rand = G4UniformRand();
1267 if(rand > fReflectivity + fTransmittance)
1268 {
1269 DoAbsorption();
1270 }
1271 else if(rand > fReflectivity)
1272 {
1276 }
1277 else
1278 {
1280 {
1282 }
1283 else
1284 {
1285 swap = !swap;
1288 }
1291
1292 DoReflection();
1293
1296
1297 goto leap;
1298 }
1299 }
1300 }
1301}
1302
1303//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1306{
1307 *condition = Forced;
1308 return DBL_MAX;
1309}
1310
1311//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1313{
1314 return pi - std::acos(fOldMomentum * fFacetNormal /
1316}
1317
1318//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1320 G4double E1_parl,
1321 G4double incidentangle,
1322 G4double realRindex,
1323 G4double imaginaryRindex)
1324{
1325 G4complex reflectivity, reflectivity_TE, reflectivity_TM;
1326 G4complex N1(fRindex1, 0.), N2(realRindex, imaginaryRindex);
1327 G4complex cosPhi;
1328
1329 G4complex u(1., 0.); // unit number 1
1330
1331 G4complex numeratorTE; // E1_perp=1 E1_parl=0 -> TE polarization
1332 G4complex numeratorTM; // E1_parl=1 E1_perp=0 -> TM polarization
1333 G4complex denominatorTE, denominatorTM;
1334 G4complex rTM, rTE;
1335
1339 if(ppR && ppI)
1340 {
1341 G4double rRindex = ppR->Value(fPhotonMomentum, idx_rrindex);
1342 G4double iRindex = ppI->Value(fPhotonMomentum, idx_irindex);
1343 N1 = G4complex(rRindex, iRindex);
1344 }
1345
1346 // Following two equations, rTM and rTE, are from: "Introduction To Modern
1347 // Optics" written by Fowles
1348 cosPhi = std::sqrt(u - ((std::sin(incidentangle) * std::sin(incidentangle)) *
1349 (N1 * N1) / (N2 * N2)));
1350
1351 numeratorTE = N1 * std::cos(incidentangle) - N2 * cosPhi;
1352 denominatorTE = N1 * std::cos(incidentangle) + N2 * cosPhi;
1353 rTE = numeratorTE / denominatorTE;
1354
1355 numeratorTM = N2 * std::cos(incidentangle) - N1 * cosPhi;
1356 denominatorTM = N2 * std::cos(incidentangle) + N1 * cosPhi;
1357 rTM = numeratorTM / denominatorTM;
1358
1359 // This is my (PG) calculaton for reflectivity on a metallic surface
1360 // depending on the fraction of TE and TM polarization
1361 // when TE polarization, E1_parl=0 and E1_perp=1, R=abs(rTE)^2 and
1362 // when TM polarization, E1_parl=1 and E1_perp=0, R=abs(rTM)^2
1363
1364 reflectivity_TE = (rTE * conj(rTE)) * (E1_perp * E1_perp) /
1365 (E1_perp * E1_perp + E1_parl * E1_parl);
1366 reflectivity_TM = (rTM * conj(rTM)) * (E1_parl * E1_parl) /
1367 (E1_perp * E1_perp + E1_parl * E1_parl);
1368 reflectivity = reflectivity_TE + reflectivity_TM;
1369
1370 do
1371 {
1372 if(G4UniformRand() * real(reflectivity) > real(reflectivity_TE))
1373 {
1374 f_iTE = -1;
1375 }
1376 else
1377 {
1378 f_iTE = 1;
1379 }
1380 if(G4UniformRand() * real(reflectivity) > real(reflectivity_TM))
1381 {
1382 f_iTM = -1;
1383 }
1384 else
1385 {
1386 f_iTM = 1;
1387 }
1388 // Loop checking, 13-Aug-2015, Peter Gumplinger
1389 } while(f_iTE < 0 && f_iTM < 0);
1390
1391 return real(reflectivity);
1392}
1393
1394//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1395
1397{
1399 G4double imaginaryRindex =
1401
1402 // calculate FacetNormal
1403 if(fFinish == ground)
1404 {
1406 }
1407 else
1408 {
1410 }
1411
1413 if(std::abs(cost1) < 1.0 - fCarTolerance)
1414 {
1415 fSint1 = std::sqrt(1. - cost1 * cost1);
1416 }
1417 else
1418 {
1419 fSint1 = 0.0;
1420 }
1421
1422 G4ThreeVector A_trans, A_paral, E1pp, E1pl;
1423 G4double E1_perp, E1_parl;
1424
1425 if(fSint1 > 0.0)
1426 {
1427 A_trans = (fOldMomentum.cross(fFacetNormal)).unit();
1428 E1_perp = fOldPolarization * A_trans;
1429 E1pp = E1_perp * A_trans;
1430 E1pl = fOldPolarization - E1pp;
1431 E1_parl = E1pl.mag();
1432 }
1433 else
1434 {
1435 A_trans = fOldPolarization;
1436 // Here we Follow Jackson's conventions and we set the parallel
1437 // component = 1 in case of a ray perpendicular to the surface
1438 E1_perp = 0.0;
1439 E1_parl = 1.0;
1440 }
1441
1442 G4double incidentangle = GetIncidentAngle();
1443
1444 // calculate the reflectivity depending on incident angle,
1445 // polarization and complex refractive
1446 fReflectivity = GetReflectivity(E1_perp, E1_parl, incidentangle, realRindex,
1447 imaginaryRindex);
1448}
1449
1450//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1452{
1453 G4Step aStep = *pStep;
1455
1457 if(sd != nullptr)
1458 return sd->Hit(&aStep);
1459 else
1460 return false;
1461}
1462
1463//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1465{
1466 fInvokeSD = flag;
1468}
1469
1470//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1472{
1473 verboseLevel = verbose;
1475}
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
@ FatalException
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
@ Forced
static const G4double alpha
@ kBACKSCATTERCONSTANT
@ kSPECULARLOBECONSTANT
@ kSPECULARSPIKECONSTANT
G4OpBoundaryProcessStatus
@ PolishedTiOAirReflection
@ GroundTeflonAirReflection
@ EtchedVM2000AirReflection
@ Transmission
@ EtchedVM2000GlueReflection
@ GroundLumirrorGlueReflection
@ LobeReflection
@ GroundTiOAirReflection
@ GroundTyvekAirReflection
@ FresnelReflection
@ NoRINDEX
@ GroundAirReflection
@ EtchedAirReflection
@ PolishedVM2000GlueReflection
@ PolishedAirReflection
@ PolishedTeflonAirReflection
@ EtchedTyvekAirReflection
@ SpikeReflection
@ EtchedLumirrorGlueReflection
@ GroundVM2000AirReflection
@ PolishedTyvekAirReflection
@ PolishedVM2000AirReflection
@ EtchedTiOAirReflection
@ EtchedTeflonAirReflection
@ GroundVM2000GlueReflection
@ BackScattering
@ PolishedLumirrorGlueReflection
@ TotalInternalReflection
@ PolishedLumirrorAirReflection
@ EtchedLumirrorAirReflection
@ GroundLumirrorAirReflection
@ FresnelRefraction
@ Detection
@ fOpBoundary
@ unified
@ DAVIS
@ glisur
@ groundfrontpainted
@ polishedbackpainted
@ polished
@ ground
@ polishedfrontpainted
@ groundbackpainted
G4ProcessType
static constexpr double perCent
Definition: G4SIunits.hh:325
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double nm
Definition: G4SIunits.hh:92
static constexpr double pi
Definition: G4SIunits.hh:55
static constexpr double halfpi
Definition: G4SIunits.hh:57
@ fGeomBoundary
Definition: G4StepStatus.hh:43
G4SurfaceType
@ dielectric_metal
@ dielectric_LUT
@ 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
@ NotAtBoundary
@ SameMaterial
@ Absorption
@ StepTooSmall
@ LambertianReflection
#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:252
G4OpticalSurfaceFinish fFinish
G4OpBoundaryProcessStatus fStatus
G4ThreeVector fNewPolarization
G4double GetReflectivity(G4double E1_perp, G4double E1_parl, G4double incidentangle, G4double RealRindex, G4double ImaginaryRindex)
void BoundaryProcessVerbose(void) const
G4OpBoundaryProcess(const G4String &processName="OpBoundary", G4ProcessType type=fOptical)
G4OpticalSurface * fOpticalSurface
G4OpticalSurfaceModel fModel
G4ThreeVector GetFacetNormal(const G4ThreeVector &Momentum, const G4ThreeVector &Normal) const
G4Physics2DVector * fDichroicVector
G4MaterialPropertyVector * fRealRIndexMPV
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4MaterialPropertyVector * fImagRIndexMPV
virtual void SetInvokeSD(G4bool)
virtual G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *condition) override
G4bool G4BooleanRand(const G4double prob) const
virtual void PreparePhysicsTable(const G4ParticleDefinition &) override
G4bool InvokeSD(const G4Step *step)
G4ThreeVector fOldPolarization
G4bool GetBoundaryInvokeSD() const
void SetBoundaryInvokeSD(G4bool)
void SetBoundaryVerboseLevel(G4int)
G4int GetBoundaryVerboseLevel() const
static G4OpticalParameters * Instance()
G4double GetAngularDistributionValueLUT(G4int)
G4OpticalSurfaceModel GetModel() const
G4double GetSigmaAlpha() const
G4OpticalSurfaceFinish GetFinish() const
G4int GetThetaIndexMax(void) const
G4int GetPhiIndexMax(void) const
G4double GetPolish() const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
G4Physics2DVector * GetDichroicVector()
G4double GetReflectivityLUTValue(G4int)
G4int GetLUTbins(void) const
G4double GetAngularDistributionValue(G4int, G4int, G4int)
static const G4Step * GetHyperStep()
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
void ProposeVelocity(G4double finalVelocity)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
virtual void Initialize(const G4Track &)
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:327
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
G4bool Hit(G4Step *aStep)
static constexpr double deg
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:79
ThreeVector shoot(const G4int Ap, const G4int Af)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
float c_light
Definition: hepunit.py:256
float h_Planck
Definition: hepunit.py:262
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