Geant4  10.02.p03
G4OpBoundaryProcess Class Reference

#include <G4OpBoundaryProcess.hh>

Inheritance diagram for G4OpBoundaryProcess:
Collaboration diagram for G4OpBoundaryProcess:

Public Member Functions

 G4OpBoundaryProcess (const G4String &processName="OpBoundary", G4ProcessType type=fOptical)
 
 ~G4OpBoundaryProcess ()
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *condition)
 
G4VParticleChange * PostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
G4OpBoundaryProcessStatus GetStatus () const
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChange * AtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChange * AlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Private Member Functions

 G4OpBoundaryProcess (const G4OpBoundaryProcess &right)
 
G4OpBoundaryProcessoperator= (const G4OpBoundaryProcess &right)
 
G4bool G4BooleanRand (const G4double prob) const
 
G4ThreeVector GetFacetNormal (const G4ThreeVector &Momentum, const G4ThreeVector &Normal) const
 
void DielectricMetal ()
 
void DielectricDielectric ()
 
void DielectricLUT ()
 
void DielectricDichroic ()
 
void ChooseReflection ()
 
void DoAbsorption ()
 
void DoReflection ()
 
G4double GetIncidentAngle ()
 
G4double GetReflectivity (G4double E1_perp, G4double E1_parl, G4double incidentangle, G4double RealRindex, G4double ImaginaryRindex)
 
void CalculateReflectivity (void)
 
void BoundaryProcessVerbose (void) const
 
G4bool InvokeSD (const G4Step *step)
 

Private Attributes

G4double thePhotonMomentum
 
G4ThreeVector OldMomentum
 
G4ThreeVector OldPolarization
 
G4ThreeVector NewMomentum
 
G4ThreeVector NewPolarization
 
G4ThreeVector theGlobalNormal
 
G4ThreeVector theFacetNormal
 
G4MaterialMaterial1
 
G4MaterialMaterial2
 
G4OpticalSurfaceOpticalSurface
 
G4MaterialPropertyVectorPropertyPointer
 
G4MaterialPropertyVectorPropertyPointer1
 
G4MaterialPropertyVectorPropertyPointer2
 
G4double Rindex1
 
G4double Rindex2
 
G4double cost1
 
G4double cost2
 
G4double sint1
 
G4double sint2
 
G4OpBoundaryProcessStatus theStatus
 
G4OpticalSurfaceModel theModel
 
G4OpticalSurfaceFinish theFinish
 
G4double theReflectivity
 
G4double theEfficiency
 
G4double theTransmittance
 
G4double theSurfaceRoughness
 
G4double prob_sl
 
G4double prob_ss
 
G4double prob_bs
 
G4int iTE
 
G4int iTM
 
G4double kCarTolerance
 
size_t idx
 
size_t idy
 
G4Physics2DVectorDichroicVector
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChange * pParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 131 of file G4OpBoundaryProcess.hh.

Constructor & Destructor Documentation

◆ G4OpBoundaryProcess() [1/2]

G4OpBoundaryProcess::G4OpBoundaryProcess ( const G4String processName = "OpBoundary",
G4ProcessType  type = fOptical 
)

Definition at line 103 of file G4OpBoundaryProcess.cc.

105  : G4VDiscreteProcess(processName, type)
106 {
107  if ( verboseLevel > 0) {
108  G4cout << GetProcessName() << " is created " << G4endl;
109  }
110 
112 
114  theModel = glisur;
116  theReflectivity = 1.;
117  theEfficiency = 0.;
118  theTransmittance = 0.;
119 
120  theSurfaceRoughness = 0.;
121 
122  prob_sl = 0.;
123  prob_ss = 0.;
124  prob_bs = 0.;
125 
126  PropertyPointer = NULL;
127  PropertyPointer1 = NULL;
128  PropertyPointer2 = NULL;
129 
130  Material1 = NULL;
131  Material2 = NULL;
132 
133  OpticalSurface = NULL;
134 
137 
138  iTE = iTM = 0;
139  thePhotonMomentum = 0.;
140  Rindex1 = Rindex2 = 1.;
141  cost1 = cost2 = sint1 = sint2 = 0.;
142 
143  idx = idy = 0;
144  DichroicVector = NULL;
145 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4Physics2DVector * DichroicVector
G4OpticalSurface * OpticalSurface
G4double GetSurfaceTolerance() const
G4MaterialPropertyVector * PropertyPointer1
G4OpticalSurfaceFinish theFinish
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4OpBoundaryProcessStatus theStatus
G4GLOB_DLL std::ostream G4cout
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
G4OpticalSurfaceModel theModel
G4MaterialPropertyVector * PropertyPointer
#define G4endl
Definition: G4ios.hh:61
G4MaterialPropertyVector * PropertyPointer2
static G4GeometryTolerance * GetInstance()
Here is the call graph for this function:

◆ ~G4OpBoundaryProcess()

G4OpBoundaryProcess::~G4OpBoundaryProcess ( )

Definition at line 155 of file G4OpBoundaryProcess.cc.

155 {}

◆ G4OpBoundaryProcess() [2/2]

G4OpBoundaryProcess::G4OpBoundaryProcess ( const G4OpBoundaryProcess right)
private

Member Function Documentation

◆ BoundaryProcessVerbose()

void G4OpBoundaryProcess::BoundaryProcessVerbose ( void  ) const
private

Definition at line 545 of file G4OpBoundaryProcess.cc.

546 {
547  if ( theStatus == Undefined )
548  G4cout << " *** Undefined *** " << G4endl;
549  if ( theStatus == Transmission )
550  G4cout << " *** Transmission *** " << G4endl;
551  if ( theStatus == FresnelRefraction )
552  G4cout << " *** FresnelRefraction *** " << G4endl;
553  if ( theStatus == FresnelReflection )
554  G4cout << " *** FresnelReflection *** " << G4endl;
556  G4cout << " *** TotalInternalReflection *** " << G4endl;
558  G4cout << " *** LambertianReflection *** " << G4endl;
559  if ( theStatus == LobeReflection )
560  G4cout << " *** LobeReflection *** " << G4endl;
561  if ( theStatus == SpikeReflection )
562  G4cout << " *** SpikeReflection *** " << G4endl;
563  if ( theStatus == BackScattering )
564  G4cout << " *** BackScattering *** " << G4endl;
566  G4cout << " *** PolishedLumirrorAirReflection *** " << G4endl;
568  G4cout << " *** PolishedLumirrorGlueReflection *** " << G4endl;
570  G4cout << " *** PolishedAirReflection *** " << G4endl;
572  G4cout << " *** PolishedTeflonAirReflection *** " << G4endl;
574  G4cout << " *** PolishedTiOAirReflection *** " << G4endl;
576  G4cout << " *** PolishedTyvekAirReflection *** " << G4endl;
578  G4cout << " *** PolishedVM2000AirReflection *** " << G4endl;
580  G4cout << " *** PolishedVM2000GlueReflection *** " << G4endl;
582  G4cout << " *** EtchedLumirrorAirReflection *** " << G4endl;
584  G4cout << " *** EtchedLumirrorGlueReflection *** " << G4endl;
586  G4cout << " *** EtchedAirReflection *** " << G4endl;
588  G4cout << " *** EtchedTeflonAirReflection *** " << G4endl;
590  G4cout << " *** EtchedTiOAirReflection *** " << G4endl;
592  G4cout << " *** EtchedTyvekAirReflection *** " << G4endl;
594  G4cout << " *** EtchedVM2000AirReflection *** " << G4endl;
596  G4cout << " *** EtchedVM2000GlueReflection *** " << G4endl;
598  G4cout << " *** GroundLumirrorAirReflection *** " << G4endl;
600  G4cout << " *** GroundLumirrorGlueReflection *** " << G4endl;
602  G4cout << " *** GroundAirReflection *** " << G4endl;
604  G4cout << " *** GroundTeflonAirReflection *** " << G4endl;
606  G4cout << " *** GroundTiOAirReflection *** " << G4endl;
608  G4cout << " *** GroundTyvekAirReflection *** " << G4endl;
610  G4cout << " *** GroundVM2000AirReflection *** " << G4endl;
612  G4cout << " *** GroundVM2000GlueReflection *** " << G4endl;
613  if ( theStatus == Absorption )
614  G4cout << " *** Absorption *** " << G4endl;
615  if ( theStatus == Detection )
616  G4cout << " *** Detection *** " << G4endl;
617  if ( theStatus == NotAtBoundary )
618  G4cout << " *** NotAtBoundary *** " << G4endl;
619  if ( theStatus == SameMaterial )
620  G4cout << " *** SameMaterial *** " << G4endl;
621  if ( theStatus == StepTooSmall )
622  G4cout << " *** StepTooSmall *** " << G4endl;
623  if ( theStatus == NoRINDEX )
624  G4cout << " *** NoRINDEX *** " << G4endl;
625  if ( theStatus == Dichroic )
626  G4cout << " *** Dichroic Transmission *** " << G4endl;
627 }
G4OpBoundaryProcessStatus theStatus
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
Here is the caller graph for this function:

◆ CalculateReflectivity()

void G4OpBoundaryProcess::CalculateReflectivity ( void  )
private

Definition at line 1286 of file G4OpBoundaryProcess.cc.

1287 {
1288  G4double RealRindex =
1290  G4double ImaginaryRindex =
1292 
1293  // calculate FacetNormal
1294  if ( theFinish == ground ) {
1295  theFacetNormal =
1297  } else {
1299  }
1300 
1302  cost1 = -PdotN;
1303 
1304  if (std::abs(cost1) < 1.0 - kCarTolerance) {
1305  sint1 = std::sqrt(1. - cost1*cost1);
1306  } else {
1307  sint1 = 0.0;
1308  }
1309 
1310  G4ThreeVector A_trans, A_paral, E1pp, E1pl;
1311  G4double E1_perp, E1_parl;
1312 
1313  if (sint1 > 0.0 ) {
1314  A_trans = OldMomentum.cross(theFacetNormal);
1315  A_trans = A_trans.unit();
1316  E1_perp = OldPolarization * A_trans;
1317  E1pp = E1_perp * A_trans;
1318  E1pl = OldPolarization - E1pp;
1319  E1_parl = E1pl.mag();
1320  }
1321  else {
1322  A_trans = OldPolarization;
1323  // Here we Follow Jackson's conventions and we set the
1324  // parallel component = 1 in case of a ray perpendicular
1325  // to the surface
1326  E1_perp = 0.0;
1327  E1_parl = 1.0;
1328  }
1329 
1330  //calculate incident angle
1331  G4double incidentangle = GetIncidentAngle();
1332 
1333  //calculate the reflectivity depending on incident angle,
1334  //polarization and complex refractive
1335 
1336  theReflectivity =
1337  GetReflectivity(E1_perp, E1_parl, incidentangle,
1338  RealRindex, ImaginaryRindex);
1339 }
G4double GetReflectivity(G4double E1_perp, G4double E1_parl, G4double incidentangle, G4double RealRindex, G4double ImaginaryRindex)
G4MaterialPropertyVector * PropertyPointer1
G4OpticalSurfaceFinish theFinish
Hep3Vector cross(const Hep3Vector &) const
double mag() const
Hep3Vector unit() const
G4double Value(G4double theEnergy, size_t &lastidx) const
G4ThreeVector GetFacetNormal(const G4ThreeVector &Momentum, const G4ThreeVector &Normal) const
G4MaterialPropertyVector * PropertyPointer2
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ChooseReflection()

void G4OpBoundaryProcess::ChooseReflection ( )
inlineprivate

Definition at line 286 of file G4OpBoundaryProcess.hh.

287 {
288  G4double rand = G4UniformRand();
289  if ( rand >= 0.0 && rand < prob_ss ) {
292  }
293  else if ( rand >= prob_ss &&
294  rand <= prob_ss+prob_sl) {
296  }
297  else if ( rand > prob_ss+prob_sl &&
298  rand < prob_ss+prob_sl+prob_bs ) {
300  }
301  else {
303  }
304 }
G4OpBoundaryProcessStatus theStatus
#define G4UniformRand()
Definition: Randomize.hh:97
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ DielectricDichroic()

void G4OpBoundaryProcess::DielectricDichroic ( )
private

Definition at line 872 of file G4OpBoundaryProcess.cc.

873 {
874  // Calculate Angle between Normal and Photon Momentum
875  G4double anglePhotonToNormal = OldMomentum.angle(-theGlobalNormal);
876 
877  // Round it to closest integer
878  G4double angleIncident = std::floor(180/pi*anglePhotonToNormal+0.5);
879 
880  if (!DichroicVector) {
882  }
883 
884 
885  if (DichroicVector) {
888  DichroicVector->Value(wavelength/nm,angleIncident,idx,idy)*perCent;
889 // G4cout << "wavelength: " << std::floor(wavelength/nm)
890 // << "nm" << G4endl;
891 // G4cout << "Incident angle: " << angleIncident << "deg" << G4endl;
892 // G4cout << "Transmittance: "
893 // << std::floor(theTransmittance/perCent) << "%" << G4endl;
894  } else {
896  ed << " G4OpBoundaryProcess/DielectricDichroic(): "
897  << " The dichroic surface has no G4Physics2DVector"
898  << G4endl;
899  G4Exception("G4OpBoundaryProcess::DielectricDichroic", "OpBoun03",
900  FatalException,ed,
901  "A dichroic surface must have an associated G4Physics2DVector");
902  }
903 
904  if ( !G4BooleanRand(theTransmittance) ) { // Not transmitted, so reflect
905 
906  if ( theModel == glisur || theFinish == polished ) {
907  DoReflection();
908  } else {
910  if ( theStatus == LambertianReflection ) {
911  DoReflection();
912  } else if ( theStatus == BackScattering ) {
915  } else {
916  G4double PdotN, EdotN;
917  do {
920  PdotN = OldMomentum * theFacetNormal;
921  NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
922  // Loop checking, 13-Aug-2015, Peter Gumplinger
923  } while (NewMomentum * theGlobalNormal <= 0.0);
925  NewPolarization = -OldPolarization + (2.*EdotN)*theFacetNormal;
926  }
927  }
928 
929  } else {
930 
934 
935  }
936 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4Physics2DVector * DichroicVector
float h_Planck
Definition: hepunit.py:263
G4OpticalSurface * OpticalSurface
G4OpticalSurfaceFinish theFinish
double angle(const Hep3Vector &) const
G4OpBoundaryProcessStatus theStatus
G4bool G4BooleanRand(const G4double prob) const
static const double nm
Definition: G4SIunits.hh:111
static const double perCent
Definition: G4SIunits.hh:329
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4OpticalSurfaceModel theModel
G4ThreeVector GetFacetNormal(const G4ThreeVector &Momentum, const G4ThreeVector &Normal) const
static const double pi
Definition: G4SIunits.hh:74
G4Physics2DVector * GetDichroicVector()
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double Value(G4double x, G4double y, size_t &lastidx, size_t &lastidy) const
float c_light
Definition: hepunit.py:257
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DielectricDielectric()

void G4OpBoundaryProcess::DielectricDielectric ( )
private

Definition at line 938 of file G4OpBoundaryProcess.cc.

939 {
940  G4bool Inside = false;
941  G4bool Swap = false;
942 
943  G4bool SurfaceRoughnessCriterionPass = 1;
944  if (theSurfaceRoughness != 0. && Rindex1 > Rindex2) {
946  G4double SurfaceRoughnessCriterion =
947  std::exp(-std::pow((4*pi*theSurfaceRoughness*Rindex1*cost1/wavelength),2));
948  SurfaceRoughnessCriterionPass =
949  G4BooleanRand(SurfaceRoughnessCriterion);
950  }
951 
952  leap:
953 
954  G4bool Through = false;
955  G4bool Done = false;
956 
957  G4double PdotN, EdotN;
958 
959  G4ThreeVector A_trans, A_paral, E1pp, E1pl;
960  G4double E1_perp, E1_parl;
961  G4double s1, s2, E2_perp, E2_parl, E2_total, TransCoeff;
962  G4double E2_abs, C_parl, C_perp;
963  G4double alpha;
964 
965  do {
966 
967  if (Through) {
968  Swap = !Swap;
969  Through = false;
973  }
974 
975  if ( theFinish == polished ) {
977  }
978  else {
981  }
982 
983  PdotN = OldMomentum * theFacetNormal;
985 
986  cost1 = - PdotN;
987  if (std::abs(cost1) < 1.0-kCarTolerance){
988  sint1 = std::sqrt(1.-cost1*cost1);
989  sint2 = sint1*Rindex1/Rindex2; // *** Snell's Law ***
990  }
991  else {
992  sint1 = 0.0;
993  sint2 = 0.0;
994  }
995 
996  if (sint2 >= 1.0) {
997 
998  // Simulate total internal reflection
999 
1000  if (Swap) Swap = !Swap;
1001 
1003 
1004  if ( !SurfaceRoughnessCriterionPass ) theStatus =
1006 
1007  if ( theModel == unified && theFinish != polished )
1008  ChooseReflection();
1009 
1010  if ( theStatus == LambertianReflection ) {
1011  DoReflection();
1012  }
1013  else if ( theStatus == BackScattering ) {
1016  }
1017  else {
1018 
1019  PdotN = OldMomentum * theFacetNormal;
1020  NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
1021  EdotN = OldPolarization * theFacetNormal;
1022  NewPolarization = -OldPolarization + (2.*EdotN)*theFacetNormal;
1023 
1024  }
1025  }
1026  else if (sint2 < 1.0) {
1027 
1028  // Calculate amplitude for transmission (Q = P x N)
1029 
1030  if (cost1 > 0.0) {
1031  cost2 = std::sqrt(1.-sint2*sint2);
1032  }
1033  else {
1034  cost2 = -std::sqrt(1.-sint2*sint2);
1035  }
1036 
1037  if (sint1 > 0.0) {
1038  A_trans = OldMomentum.cross(theFacetNormal);
1039  A_trans = A_trans.unit();
1040  E1_perp = OldPolarization * A_trans;
1041  E1pp = E1_perp * A_trans;
1042  E1pl = OldPolarization - E1pp;
1043  E1_parl = E1pl.mag();
1044  }
1045  else {
1046  A_trans = OldPolarization;
1047  // Here we Follow Jackson's conventions and we set the
1048  // parallel component = 1 in case of a ray perpendicular
1049  // to the surface
1050  E1_perp = 0.0;
1051  E1_parl = 1.0;
1052  }
1053 
1054  s1 = Rindex1*cost1;
1055  E2_perp = 2.*s1*E1_perp/(Rindex1*cost1+Rindex2*cost2);
1056  E2_parl = 2.*s1*E1_parl/(Rindex2*cost1+Rindex1*cost2);
1057  E2_total = E2_perp*E2_perp + E2_parl*E2_parl;
1058  s2 = Rindex2*cost2*E2_total;
1059 
1060  if (theTransmittance > 0) TransCoeff = theTransmittance;
1061  else if (cost1 != 0.0) TransCoeff = s2/s1;
1062  else TransCoeff = 0.0;
1063 
1064  if ( !G4BooleanRand(TransCoeff) ) {
1065 
1066  // Simulate reflection
1067 
1068  if (Swap) Swap = !Swap;
1069 
1071 
1072  if ( !SurfaceRoughnessCriterionPass ) theStatus =
1074 
1075  if ( theModel == unified && theFinish != polished )
1076  ChooseReflection();
1077 
1078  if ( theStatus == LambertianReflection ) {
1079  DoReflection();
1080  }
1081  else if ( theStatus == BackScattering ) {
1084  }
1085  else {
1086 
1087  PdotN = OldMomentum * theFacetNormal;
1088  NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
1089 
1090  if (sint1 > 0.0) { // incident ray oblique
1091 
1092  E2_parl = Rindex2*E2_parl/Rindex1 - E1_parl;
1093  E2_perp = E2_perp - E1_perp;
1094  E2_total = E2_perp*E2_perp + E2_parl*E2_parl;
1095  A_paral = NewMomentum.cross(A_trans);
1096  A_paral = A_paral.unit();
1097  E2_abs = std::sqrt(E2_total);
1098  C_parl = E2_parl/E2_abs;
1099  C_perp = E2_perp/E2_abs;
1100 
1101  NewPolarization = C_parl*A_paral + C_perp*A_trans;
1102 
1103  }
1104 
1105  else { // incident ray perpendicular
1106 
1107  if (Rindex2 > Rindex1) {
1109  }
1110  else {
1112  }
1113 
1114  }
1115  }
1116  }
1117  else { // photon gets transmitted
1118 
1119  // Simulate transmission/refraction
1120 
1121  Inside = !Inside;
1122  Through = true;
1124 
1125  if (sint1 > 0.0) { // incident ray oblique
1126 
1127  alpha = cost1 - cost2*(Rindex2/Rindex1);
1129  NewMomentum = NewMomentum.unit();
1130 // PdotN = -cost2;
1131  A_paral = NewMomentum.cross(A_trans);
1132  A_paral = A_paral.unit();
1133  E2_abs = std::sqrt(E2_total);
1134  C_parl = E2_parl/E2_abs;
1135  C_perp = E2_perp/E2_abs;
1136 
1137  NewPolarization = C_parl*A_paral + C_perp*A_trans;
1138 
1139  }
1140  else { // incident ray perpendicular
1141 
1142  NewMomentum = OldMomentum;
1144 
1145  }
1146  }
1147  }
1148 
1149  OldMomentum = NewMomentum.unit();
1151 
1152  if (theStatus == FresnelRefraction) {
1153  Done = (NewMomentum * theGlobalNormal <= 0.0);
1154  }
1155  else {
1156  Done = (NewMomentum * theGlobalNormal >= -kCarTolerance);
1157  }
1158 
1159  // Loop checking, 13-Aug-2015, Peter Gumplinger
1160  } while (!Done);
1161 
1162  if (Inside && !Swap) {
1163  if( theFinish == polishedbackpainted ||
1165 
1166  G4double rand = G4UniformRand();
1167  if ( rand > theReflectivity ) {
1168  if (rand > theReflectivity + theTransmittance) {
1169  DoAbsorption();
1170  } else {
1172  NewMomentum = OldMomentum;
1174  }
1175  }
1176  else {
1177  if (theStatus != FresnelRefraction ) {
1179  }
1180  else {
1181  Swap = !Swap;
1184  }
1185  if ( theFinish == groundbackpainted )
1187 
1188  DoReflection();
1189 
1192 
1193  goto leap;
1194  }
1195  }
1196  }
1197 }
void G4SwapObj(T *a, T *b)
Definition: templates.hh:129
float h_Planck
Definition: hepunit.py:263
G4OpticalSurfaceFinish theFinish
G4OpBoundaryProcessStatus theStatus
#define G4UniformRand()
Definition: Randomize.hh:97
Hep3Vector cross(const Hep3Vector &) const
G4bool G4BooleanRand(const G4double prob) const
bool G4bool
Definition: G4Types.hh:79
double mag() const
Hep3Vector unit() const
G4OpticalSurfaceModel theModel
G4ThreeVector GetFacetNormal(const G4ThreeVector &Momentum, const G4ThreeVector &Normal) const
static const double pi
Definition: G4SIunits.hh:74
void G4SwapPtr(T *&a, T *&b)
Definition: templates.hh:121
double G4double
Definition: G4Types.hh:76
static const G4double alpha
float c_light
Definition: hepunit.py:257
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DielectricLUT()

void G4OpBoundaryProcess::DielectricLUT ( )
private

Definition at line 802 of file G4OpBoundaryProcess.cc.

803 {
804  G4int thetaIndex, phiIndex;
805  G4double AngularDistributionValue, thetaRad, phiRad, EdotN;
806  G4ThreeVector PerpendicularVectorTheta, PerpendicularVectorPhi;
807 
810 
811  G4int thetaIndexMax = OpticalSurface->GetThetaIndexMax();
812  G4int phiIndexMax = OpticalSurface->GetPhiIndexMax();
813 
814  G4double rand;
815 
816  do {
817  rand = G4UniformRand();
818  if ( rand > theReflectivity ) {
819  if (rand > theReflectivity + theTransmittance) {
820  DoAbsorption();
821  } else {
825  }
826  break;
827  }
828  else {
829  // Calculate Angle between Normal and Photon Momentum
830  G4double anglePhotonToNormal =
832  // Round it to closest integer
833  G4int angleIncident = G4int(std::floor(180/pi*anglePhotonToNormal+0.5));
834 
835  // Take random angles THETA and PHI,
836  // and see if below Probability - if not - Redo
837  do {
838  thetaIndex = G4RandFlat::shootInt(thetaIndexMax-1);
839  phiIndex = G4RandFlat::shootInt(phiIndexMax-1);
840  // Find probability with the new indeces from LUT
841  AngularDistributionValue = OpticalSurface ->
842  GetAngularDistributionValue(angleIncident,
843  thetaIndex,
844  phiIndex);
845  // Loop checking, 13-Aug-2015, Peter Gumplinger
846  } while ( !G4BooleanRand(AngularDistributionValue) );
847 
848  thetaRad = (-90 + 4*thetaIndex)*pi/180;
849  phiRad = (-90 + 5*phiIndex)*pi/180;
850  // Rotate Photon Momentum in Theta, then in Phi
852 
853  PerpendicularVectorTheta = NewMomentum.cross(theGlobalNormal);
854  if (PerpendicularVectorTheta.mag() < kCarTolerance )
855  PerpendicularVectorTheta = NewMomentum.orthogonal();
856  NewMomentum =
857  NewMomentum.rotate(anglePhotonToNormal-thetaRad,
858  PerpendicularVectorTheta);
859  PerpendicularVectorPhi =
860  PerpendicularVectorTheta.cross(NewMomentum);
861  NewMomentum = NewMomentum.rotate(-phiRad,PerpendicularVectorPhi);
862 
863  // Rotate Polarization too:
866  NewPolarization = -OldPolarization + (2.*EdotN)*theFacetNormal;
867  }
868  // Loop checking, 13-Aug-2015, Peter Gumplinger
869  } while (NewMomentum * theGlobalNormal <= 0.0);
870 }
G4OpticalSurface * OpticalSurface
G4OpBoundaryProcessStatus
int G4int
Definition: G4Types.hh:78
G4OpticalSurfaceFinish theFinish
double angle(const Hep3Vector &) const
G4OpBoundaryProcessStatus theStatus
#define G4UniformRand()
Definition: Randomize.hh:97
Hep3Vector cross(const Hep3Vector &) const
G4bool G4BooleanRand(const G4double prob) const
double mag() const
Hep3Vector orthogonal() const
static const double pi
Definition: G4SIunits.hh:74
G4int GetThetaIndexMax(void) const
double G4double
Definition: G4Types.hh:76
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:28
G4int GetPhiIndexMax(void) const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DielectricMetal()

void G4OpBoundaryProcess::DielectricMetal ( )
private

Definition at line 709 of file G4OpBoundaryProcess.cc.

710 {
711  G4int n = 0;
712  G4double rand, PdotN, EdotN;
713  G4ThreeVector A_trans, A_paral;
714 
715  do {
716 
717  n++;
718 
719  rand = G4UniformRand();
720  if ( rand > theReflectivity && n == 1 ) {
721  if (rand > theReflectivity + theTransmittance) {
722  DoAbsorption();
723  } else {
727  }
728  break;
729  }
730  else {
731 
733  if ( n > 1 ) {
735  if ( !G4BooleanRand(theReflectivity) ) {
736  DoAbsorption();
737  break;
738  }
739  }
740  }
741 
742  if ( theModel == glisur || theFinish == polished ) {
743 
744  DoReflection();
745 
746  } else {
747 
748  if ( n == 1 ) ChooseReflection();
749 
750  if ( theStatus == LambertianReflection ) {
751  DoReflection();
752  }
753  else if ( theStatus == BackScattering ) {
756  }
757  else {
758 
761  } else {
764  }
765  }
766 
767  PdotN = OldMomentum * theFacetNormal;
768  NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
770 
771  if (sint1 > 0.0 ) {
772  A_trans = OldMomentum.cross(theFacetNormal);
773  A_trans = A_trans.unit();
774  } else {
775  A_trans = OldPolarization;
776  }
777  A_paral = NewMomentum.cross(A_trans);
778  A_paral = A_paral.unit();
779 
780  if(iTE>0&&iTM>0) {
781  NewPolarization =
782  -OldPolarization + (2.*EdotN)*theFacetNormal;
783  } else if (iTE>0) {
784  NewPolarization = -A_trans;
785  } else if (iTM>0) {
786  NewPolarization = -A_paral;
787  }
788 
789  }
790 
791  }
792 
795 
796  }
797 
798  // Loop checking, 13-Aug-2015, Peter Gumplinger
799  } while (NewMomentum * theGlobalNormal < 0.0);
800 }
G4MaterialPropertyVector * PropertyPointer1
int G4int
Definition: G4Types.hh:78
G4OpticalSurfaceFinish theFinish
Char_t n[5]
G4OpBoundaryProcessStatus theStatus
#define G4UniformRand()
Definition: Randomize.hh:97
Hep3Vector cross(const Hep3Vector &) const
G4bool G4BooleanRand(const G4double prob) const
Hep3Vector unit() const
G4OpticalSurfaceModel theModel
G4ThreeVector GetFacetNormal(const G4ThreeVector &Momentum, const G4ThreeVector &Normal) const
G4MaterialPropertyVector * PropertyPointer2
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DoAbsorption()

void G4OpBoundaryProcess::DoAbsorption ( )
inlineprivate

Definition at line 307 of file G4OpBoundaryProcess.hh.

308 {
310 
311  if ( G4BooleanRand(theEfficiency) ) {
312 
313  // EnergyDeposited =/= 0 means: photon has been detected
315  aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum);
316  }
317  else {
318  aParticleChange.ProposeLocalEnergyDeposit(0.0);
319  }
320 
323 
324 // aParticleChange.ProposeEnergy(0.0);
325  aParticleChange.ProposeTrackStatus(fStopAndKill);
326 }
G4OpBoundaryProcessStatus theStatus
G4bool G4BooleanRand(const G4double prob) const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DoReflection()

void G4OpBoundaryProcess::DoReflection ( )
inlineprivate

Definition at line 329 of file G4OpBoundaryProcess.hh.

330 {
331  if ( theStatus == LambertianReflection ) {
332 
335 
336  }
337  else if ( theFinish == ground ) {
338 
341  } else {
344  }
346  NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
347 
348  }
349  else {
350 
354  NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
355 
356  }
358  NewPolarization = -OldPolarization + (2.*EdotN)*theFacetNormal;
359 }
G4ThreeVector G4LambertianRand(const G4ThreeVector &normal)
G4MaterialPropertyVector * PropertyPointer1
G4OpticalSurfaceFinish theFinish
G4OpBoundaryProcessStatus theStatus
G4ThreeVector GetFacetNormal(const G4ThreeVector &Momentum, const G4ThreeVector &Normal) const
G4MaterialPropertyVector * PropertyPointer2
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ G4BooleanRand()

G4bool G4OpBoundaryProcess::G4BooleanRand ( const G4double  prob) const
inlineprivate

Definition at line 265 of file G4OpBoundaryProcess.hh.

266 {
267  /* Returns a random boolean variable with the specified probability */
268 
269  return (G4UniformRand() < prob);
270 }
#define G4UniformRand()
Definition: Randomize.hh:97
Here is the caller graph for this function:

◆ GetFacetNormal()

G4ThreeVector G4OpBoundaryProcess::GetFacetNormal ( const G4ThreeVector Momentum,
const G4ThreeVector Normal 
) const
private

Definition at line 630 of file G4OpBoundaryProcess.cc.

632 {
633  G4ThreeVector FacetNormal;
634 
635  if (theModel == unified || theModel == LUT) {
636 
637  /* This function code alpha to a random value taken from the
638  distribution p(alpha) = g(alpha; 0, sigma_alpha)*std::sin(alpha),
639  for alpha > 0 and alpha < 90, where g(alpha; 0, sigma_alpha)
640  is a gaussian distribution with mean 0 and standard deviation
641  sigma_alpha. */
642 
643  G4double alpha;
644 
645  G4double sigma_alpha = 0.0;
646  if (OpticalSurface) sigma_alpha = OpticalSurface->GetSigmaAlpha();
647 
648  if (sigma_alpha == 0.0) return FacetNormal = Normal;
649 
650  G4double f_max = std::min(1.0,4.*sigma_alpha);
651 
652  G4double phi, SinAlpha, CosAlpha, SinPhi, CosPhi, unit_x, unit_y, unit_z;
653  G4ThreeVector tmpNormal;
654 
655  do {
656  do {
657  alpha = G4RandGauss::shoot(0.0,sigma_alpha);
658  // Loop checking, 13-Aug-2015, Peter Gumplinger
659  } while (G4UniformRand()*f_max > std::sin(alpha) || alpha >= halfpi );
660 
661  phi = G4UniformRand()*twopi;
662 
663  SinAlpha = std::sin(alpha);
664  CosAlpha = std::cos(alpha);
665  SinPhi = std::sin(phi);
666  CosPhi = std::cos(phi);
667 
668  unit_x = SinAlpha * CosPhi;
669  unit_y = SinAlpha * SinPhi;
670  unit_z = CosAlpha;
671 
672  FacetNormal.setX(unit_x);
673  FacetNormal.setY(unit_y);
674  FacetNormal.setZ(unit_z);
675 
676  tmpNormal = Normal;
677 
678  FacetNormal.rotateUz(tmpNormal);
679  // Loop checking, 13-Aug-2015, Peter Gumplinger
680  } while (Momentum * FacetNormal >= 0.0);
681  }
682  else {
683 
684  G4double polish = 1.0;
685  if (OpticalSurface) polish = OpticalSurface->GetPolish();
686 
687  if (polish < 1.0) {
688  do {
689  G4ThreeVector smear;
690  do {
691  smear.setX(2.*G4UniformRand()-1.0);
692  smear.setY(2.*G4UniformRand()-1.0);
693  smear.setZ(2.*G4UniformRand()-1.0);
694  // Loop checking, 13-Aug-2015, Peter Gumplinger
695  } while (smear.mag()>1.0);
696  smear = (1.-polish) * smear;
697  FacetNormal = Normal + smear;
698  // Loop checking, 13-Aug-2015, Peter Gumplinger
699  } while (Momentum * FacetNormal >= 0.0);
700  FacetNormal = FacetNormal.unit();
701  }
702  else {
703  FacetNormal = Normal;
704  }
705  }
706  return FacetNormal;
707 }
ThreeVector shoot(const G4int Ap, const G4int Af)
static const double halfpi
Definition: G4SIunits.hh:76
G4OpticalSurface * OpticalSurface
void setY(double)
void setZ(double)
void setX(double)
G4double GetPolish() const
#define G4UniformRand()
Definition: Randomize.hh:97
double mag() const
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
static const double twopi
Definition: G4SIunits.hh:75
G4double GetSigmaAlpha() const
G4OpticalSurfaceModel theModel
double G4double
Definition: G4Types.hh:76
static const G4double alpha
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetIncidentAngle()

G4double G4OpBoundaryProcess::GetIncidentAngle ( )
private

Definition at line 1211 of file G4OpBoundaryProcess.cc.

1212 {
1214  G4double magP= OldMomentum.mag();
1215  G4double magN= theFacetNormal.mag();
1216  G4double incidentangle = pi - std::acos(PdotN/(magP*magN));
1217 
1218  return incidentangle;
1219 }
double mag() const
static const double pi
Definition: G4SIunits.hh:74
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetMeanFreePath()

G4double G4OpBoundaryProcess::GetMeanFreePath ( const G4Track &  ,
G4double  ,
G4ForceCondition *  condition 
)
virtual

Implements G4VDiscreteProcess.

Definition at line 1202 of file G4OpBoundaryProcess.cc.

1205 {
1206  *condition = Forced;
1207 
1208  return DBL_MAX;
1209 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83

◆ GetReflectivity()

G4double G4OpBoundaryProcess::GetReflectivity ( G4double  E1_perp,
G4double  E1_parl,
G4double  incidentangle,
G4double  RealRindex,
G4double  ImaginaryRindex 
)
private

Definition at line 1221 of file G4OpBoundaryProcess.cc.

1226 {
1227  G4complex Reflectivity, Reflectivity_TE, Reflectivity_TM;
1228  G4complex N1(Rindex1, 0), N2(RealRindex, ImaginaryRindex);
1229  G4complex CosPhi;
1230 
1231  G4complex u(1,0); //unit number 1
1232 
1233  G4complex numeratorTE; // E1_perp=1 E1_parl=0 -> TE polarization
1234  G4complex numeratorTM; // E1_parl=1 E1_perp=0 -> TM polarization
1235  G4complex denominatorTE, denominatorTM;
1236  G4complex rTM, rTE;
1237 
1238  G4MaterialPropertiesTable* aMaterialPropertiesTable =
1240  G4MaterialPropertyVector* aPropertyPointerR =
1241  aMaterialPropertiesTable->GetProperty("REALRINDEX");
1242  G4MaterialPropertyVector* aPropertyPointerI =
1243  aMaterialPropertiesTable->GetProperty("IMAGINARYRINDEX");
1244  if (aPropertyPointerR && aPropertyPointerI) {
1245  G4double RRindex = aPropertyPointerR->Value(thePhotonMomentum);
1246  G4double IRindex = aPropertyPointerI->Value(thePhotonMomentum);
1247  N1 = G4complex(RRindex,IRindex);
1248  }
1249 
1250  // Following two equations, rTM and rTE, are from: "Introduction To Modern
1251  // Optics" written by Fowles
1252 
1253  CosPhi=std::sqrt(u-((std::sin(incidentangle)*std::sin(incidentangle))*(N1*N1)/(N2*N2)));
1254 
1255  numeratorTE = N1*std::cos(incidentangle) - N2*CosPhi;
1256  denominatorTE = N1*std::cos(incidentangle) + N2*CosPhi;
1257  rTE = numeratorTE/denominatorTE;
1258 
1259  numeratorTM = N2*std::cos(incidentangle) - N1*CosPhi;
1260  denominatorTM = N2*std::cos(incidentangle) + N1*CosPhi;
1261  rTM = numeratorTM/denominatorTM;
1262 
1263  // This is my calculaton for reflectivity on a metalic surface
1264  // depending on the fraction of TE and TM polarization
1265  // when TE polarization, E1_parl=0 and E1_perp=1, R=abs(rTE)^2 and
1266  // when TM polarization, E1_parl=1 and E1_perp=0, R=abs(rTM)^2
1267 
1268  Reflectivity_TE = (rTE*conj(rTE))*(E1_perp*E1_perp)
1269  / (E1_perp*E1_perp + E1_parl*E1_parl);
1270  Reflectivity_TM = (rTM*conj(rTM))*(E1_parl*E1_parl)
1271  / (E1_perp*E1_perp + E1_parl*E1_parl);
1272  Reflectivity = Reflectivity_TE + Reflectivity_TM;
1273 
1274  do {
1275  if(G4UniformRand()*real(Reflectivity) > real(Reflectivity_TE))
1276  {iTE = -1;}else{iTE = 1;}
1277  if(G4UniformRand()*real(Reflectivity) > real(Reflectivity_TM))
1278  {iTM = -1;}else{iTM = 1;}
1279  // Loop checking, 13-Aug-2015, Peter Gumplinger
1280  } while(iTE<0&&iTM<0);
1281 
1282  return real(Reflectivity);
1283 
1284 }
G4MaterialPropertyVector * GetProperty(const char *key)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
std::complex< G4double > G4complex
Definition: G4Types.hh:81
#define G4UniformRand()
Definition: Randomize.hh:97
G4double Value(G4double theEnergy, size_t &lastidx) const
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetStatus()

G4OpBoundaryProcessStatus G4OpBoundaryProcess::GetStatus ( ) const
inline

Definition at line 280 of file G4OpBoundaryProcess.hh.

281 {
282  return theStatus;
283 }
G4OpBoundaryProcessStatus theStatus
Here is the caller graph for this function:

◆ InvokeSD()

G4bool G4OpBoundaryProcess::InvokeSD ( const G4Step *  step)
private

Definition at line 1341 of file G4OpBoundaryProcess.cc.

1342 {
1343  G4Step aStep = *pStep;
1344 
1345  aStep.AddTotalEnergyDeposit(thePhotonMomentum);
1346 
1347  G4VSensitiveDetector* sd = aStep.GetPostStepPoint()->GetSensitiveDetector();
1348  if (sd) return sd->Hit(&aStep);
1349  else return false;
1350 }
G4bool Hit(G4Step *aStep)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IsApplicable()

G4bool G4OpBoundaryProcess::IsApplicable ( const G4ParticleDefinition aParticleType)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 273 of file G4OpBoundaryProcess.hh.

275 {
276  return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
277 }
static G4OpticalPhoton * OpticalPhoton()
Here is the call graph for this function:

◆ operator=()

G4OpBoundaryProcess& G4OpBoundaryProcess::operator= ( const G4OpBoundaryProcess right)
private

◆ PostStepDoIt()

G4VParticleChange * G4OpBoundaryProcess::PostStepDoIt ( const G4Track &  aTrack,
const G4Step &  aStep 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 166 of file G4OpBoundaryProcess.cc.

167 {
169 
170  aParticleChange.Initialize(aTrack);
171  aParticleChange.ProposeVelocity(aTrack.GetVelocity());
172 
173  // Get hyperStep from G4ParallelWorldProcess
174  // NOTE: PostSetpDoIt of this process should be
175  // invoked after G4ParallelWorldProcess!
176 
177  const G4Step* pStep = &aStep;
178 
179  const G4Step* hStep = G4ParallelWorldProcess::GetHyperStep();
180 
181  if (hStep) pStep = hStep;
182 
183  G4bool isOnBoundary =
184  (pStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary);
185 
186  if (isOnBoundary) {
187  Material1 = pStep->GetPreStepPoint()->GetMaterial();
188  Material2 = pStep->GetPostStepPoint()->GetMaterial();
189  } else {
192  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
193  }
194 
195  G4VPhysicalVolume* thePrePV =
196  pStep->GetPreStepPoint() ->GetPhysicalVolume();
197  G4VPhysicalVolume* thePostPV =
198  pStep->GetPostStepPoint()->GetPhysicalVolume();
199 
200  if ( verboseLevel > 0 ) {
201  G4cout << " Photon at Boundary! " << G4endl;
202  if (thePrePV) G4cout << " thePrePV: " << thePrePV->GetName() << G4endl;
203  if (thePostPV) G4cout << " thePostPV: " << thePostPV->GetName() << G4endl;
204  }
205 
206  if (aTrack.GetStepLength()<=kCarTolerance/2){
209  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
210  }
211 
212  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
213 
214  thePhotonMomentum = aParticle->GetTotalMomentum();
215  OldMomentum = aParticle->GetMomentumDirection();
216  OldPolarization = aParticle->GetPolarization();
217 
218  if ( verboseLevel > 0 ) {
219  G4cout << " Old Momentum Direction: " << OldMomentum << G4endl;
220  G4cout << " Old Polarization: " << OldPolarization << G4endl;
221  }
222 
223  G4ThreeVector theGlobalPoint = pStep->GetPostStepPoint()->GetPosition();
224 
225  G4bool valid;
226  // Use the new method for Exit Normal in global coordinates,
227  // which provides the normal more reliably.
228 
229  // ID of Navigator which limits step
230 
232  std::vector<G4Navigator*>::iterator iNav =
234  GetActiveNavigatorsIterator();
236  (iNav[hNavId])->GetGlobalExitNormal(theGlobalPoint,&valid);
237 
238  if (valid) {
240  }
241  else
242  {
244  ed << " G4OpBoundaryProcess/PostStepDoIt(): "
245  << " The Navigator reports that it returned an invalid normal"
246  << G4endl;
247  G4Exception("G4OpBoundaryProcess::PostStepDoIt", "OpBoun01",
249  "Invalid Surface Normal - Geometry must return valid surface normal");
250  }
251 
252  if (OldMomentum * theGlobalNormal > 0.0) {
253 #ifdef G4OPTICAL_DEBUG
255  ed << " G4OpBoundaryProcess/PostStepDoIt(): "
256  << " theGlobalNormal points in a wrong direction. "
257  << G4endl;
258  ed << " The momentum of the photon arriving at interface (oldMomentum)"
259  << " must exit the volume cross in the step. " << G4endl;
260  ed << " So it MUST have dot < 0 with the normal that Exits the new volume (globalNormal)." << G4endl;
261  ed << " >> The dot product of oldMomentum and global Normal is " << OldMomentum*theGlobalNormal << G4endl;
262  ed << " Old Momentum (during step) = " << OldMomentum << G4endl;
263  ed << " Global Normal (Exiting New Vol) = " << theGlobalNormal << G4endl;
264  ed << G4endl;
265  G4Exception("G4OpBoundaryProcess::PostStepDoIt", "OpBoun02",
266  EventMustBeAborted, // Or JustWarning to see if it happens repeatedbly on one ray
267  ed,
268  "Invalid Surface Normal - Geometry must return valid surface normal pointing in the right direction");
269 #else
271 #endif
272  }
273 
274  G4MaterialPropertiesTable* aMaterialPropertiesTable;
275  G4MaterialPropertyVector* Rindex;
276 
277  aMaterialPropertiesTable = Material1->GetMaterialPropertiesTable();
278  if (aMaterialPropertiesTable) {
279  Rindex = aMaterialPropertiesTable->GetProperty("RINDEX");
280  }
281  else {
284  aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum);
285  aParticleChange.ProposeTrackStatus(fStopAndKill);
286  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
287  }
288 
289  if (Rindex) {
290  Rindex1 = Rindex->Value(thePhotonMomentum);
291  }
292  else {
295  aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum);
296  aParticleChange.ProposeTrackStatus(fStopAndKill);
297  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
298  }
299 
300  theReflectivity = 1.;
301  theEfficiency = 0.;
302  theTransmittance = 0.;
303 
304  theSurfaceRoughness = 0.;
305 
306  theModel = glisur;
308 
310 
311  Rindex = NULL;
312  OpticalSurface = NULL;
313 
314  G4LogicalSurface* Surface = NULL;
315 
316  Surface = G4LogicalBorderSurface::GetSurface(thePrePV, thePostPV);
317 
318  if (Surface == NULL){
319  G4bool enteredDaughter= (thePostPV->GetMotherLogical() ==
320  thePrePV ->GetLogicalVolume());
321  if(enteredDaughter){
322  Surface =
324  if(Surface == NULL)
325  Surface =
327  }
328  else {
329  Surface =
331  if(Surface == NULL)
332  Surface =
334  }
335  }
336 
337  if (Surface) OpticalSurface =
338  dynamic_cast <G4OpticalSurface*> (Surface->GetSurfaceProperty());
339 
340  if (OpticalSurface) {
341 
342  type = OpticalSurface->GetType();
345 
346  aMaterialPropertiesTable = OpticalSurface->
347  GetMaterialPropertiesTable();
348 
349  if (aMaterialPropertiesTable) {
350 
353  Rindex = aMaterialPropertiesTable->GetProperty("RINDEX");
354  if (Rindex) {
355  Rindex2 = Rindex->Value(thePhotonMomentum);
356  }
357  else {
360  aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum);
361  aParticleChange.ProposeTrackStatus(fStopAndKill);
362  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
363  }
364  }
365 
367  aMaterialPropertiesTable->GetProperty("REFLECTIVITY");
369  aMaterialPropertiesTable->GetProperty("REALRINDEX");
371  aMaterialPropertiesTable->GetProperty("IMAGINARYRINDEX");
372 
373  iTE = 1;
374  iTM = 1;
375 
376  if (PropertyPointer) {
377 
380 
381  } else if (PropertyPointer1 && PropertyPointer2) {
382 
384 
385  }
386 
388  aMaterialPropertiesTable->GetProperty("EFFICIENCY");
389  if (PropertyPointer) {
390  theEfficiency =
392  }
393 
395  aMaterialPropertiesTable->GetProperty("TRANSMITTANCE");
396  if (PropertyPointer) {
399  }
400 
401  if (aMaterialPropertiesTable->
402  ConstPropertyExists("SURFACEROUGHNESS"))
403  theSurfaceRoughness = aMaterialPropertiesTable->
404  GetConstProperty("SURFACEROUGHNESS");
405 
406  if ( theModel == unified ) {
408  aMaterialPropertiesTable->GetProperty("SPECULARLOBECONSTANT");
409  if (PropertyPointer) {
410  prob_sl =
412  } else {
413  prob_sl = 0.0;
414  }
415 
417  aMaterialPropertiesTable->GetProperty("SPECULARSPIKECONSTANT");
418  if (PropertyPointer) {
419  prob_ss =
421  } else {
422  prob_ss = 0.0;
423  }
424 
426  aMaterialPropertiesTable->GetProperty("BACKSCATTERCONSTANT");
427  if (PropertyPointer) {
428  prob_bs =
430  } else {
431  prob_bs = 0.0;
432  }
433  }
434  }
435  else if (theFinish == polishedbackpainted ||
437  aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum);
438  aParticleChange.ProposeTrackStatus(fStopAndKill);
439  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
440  }
441  }
442 
443  if (type == dielectric_dielectric ) {
444  if (theFinish == polished || theFinish == ground ) {
445 
446  if (Material1 == Material2){
449  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
450  }
451  aMaterialPropertiesTable =
453  if (aMaterialPropertiesTable)
454  Rindex = aMaterialPropertiesTable->GetProperty("RINDEX");
455  if (Rindex) {
456  Rindex2 = Rindex->Value(thePhotonMomentum);
457  }
458  else {
461  aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum);
462  aParticleChange.ProposeTrackStatus(fStopAndKill);
463  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
464  }
465  }
466  }
467 
468  if (type == dielectric_metal) {
469 
470  DielectricMetal();
471 
472  }
473  else if (type == dielectric_LUT) {
474 
475  DielectricLUT();
476 
477  }
478  else if (type == dielectric_dichroic) {
479 
481 
482  }
483  else if (type == dielectric_dielectric) {
484 
485  if ( theFinish == polishedbackpainted ||
488  }
489  else {
490  G4double rand = G4UniformRand();
491  if ( rand > theReflectivity ) {
492  if (rand > theReflectivity + theTransmittance) {
493  DoAbsorption();
494  } else {
498  }
499  }
500  else {
501  if ( theFinish == polishedfrontpainted ) {
502  DoReflection();
503  }
504  else if ( theFinish == groundfrontpainted ) {
506  DoReflection();
507  }
508  else {
510  }
511  }
512  }
513  }
514  else {
515 
516  G4cerr << " Error: G4BoundaryProcess: illegal boundary type " << G4endl;
517  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
518 
519  }
520 
523 
524  if ( verboseLevel > 0) {
525  G4cout << " New Momentum Direction: " << NewMomentum << G4endl;
526  G4cout << " New Polarization: " << NewPolarization << G4endl;
528  }
529 
530  aParticleChange.ProposeMomentumDirection(NewMomentum);
531  aParticleChange.ProposePolarization(NewPolarization);
532 
534  G4MaterialPropertyVector* groupvel =
536  G4double finalVelocity = groupvel->Value(thePhotonMomentum);
537  aParticleChange.ProposeVelocity(finalVelocity);
538  }
539 
540  if ( theStatus == Detection ) InvokeSD(pStep);
541 
542  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
543 }
const G4SurfaceType & GetType() const
G4MaterialPropertyVector * GetProperty(const char *key)
G4int verboseLevel
Definition: G4VProcess.hh:368
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void BoundaryProcessVerbose(void) const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
G4OpticalSurface * OpticalSurface
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
G4double GetTotalMomentum() const
G4OpticalSurfaceModel GetModel() const
G4SurfaceType
G4MaterialPropertyVector * PropertyPointer1
int G4int
Definition: G4Types.hh:78
G4OpticalSurfaceFinish theFinish
static G4LogicalBorderSurface * GetSurface(const G4VPhysicalVolume *vol1, const G4VPhysicalVolume *vol2)
G4OpBoundaryProcessStatus theStatus
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
Hep3Vector unit() const
G4double Value(G4double theEnergy, size_t &lastidx) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4TransportationManager * GetTransportationManager()
G4OpticalSurfaceModel theModel
G4bool InvokeSD(const G4Step *step)
const G4String & GetName() const
G4MaterialPropertyVector * PropertyPointer
const G4ThreeVector & GetMomentumDirection() const
G4OpticalSurfaceFinish GetFinish() const
const G4ThreeVector & GetPolarization() const
G4LogicalVolume * GetMotherLogical() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
#define G4endl
Definition: G4ios.hh:61
G4SurfaceProperty * GetSurfaceProperty() const
G4MaterialPropertyVector * PropertyPointer2
double G4double
Definition: G4Types.hh:76
G4LogicalVolume * GetLogicalVolume() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
static const G4Step * GetHyperStep()
G4GLOB_DLL std::ostream G4cerr
static G4LogicalSkinSurface * GetSurface(const G4LogicalVolume *vol)
Here is the call graph for this function:

Member Data Documentation

◆ cost1

G4double G4OpBoundaryProcess::cost1
private

Definition at line 236 of file G4OpBoundaryProcess.hh.

◆ cost2

G4double G4OpBoundaryProcess::cost2
private

Definition at line 236 of file G4OpBoundaryProcess.hh.

◆ DichroicVector

G4Physics2DVector* G4OpBoundaryProcess::DichroicVector
private

Definition at line 257 of file G4OpBoundaryProcess.hh.

◆ idx

size_t G4OpBoundaryProcess::idx
private

Definition at line 256 of file G4OpBoundaryProcess.hh.

◆ idy

size_t G4OpBoundaryProcess::idy
private

Definition at line 256 of file G4OpBoundaryProcess.hh.

◆ iTE

G4int G4OpBoundaryProcess::iTE
private

Definition at line 252 of file G4OpBoundaryProcess.hh.

◆ iTM

G4int G4OpBoundaryProcess::iTM
private

Definition at line 252 of file G4OpBoundaryProcess.hh.

◆ kCarTolerance

G4double G4OpBoundaryProcess::kCarTolerance
private

Definition at line 254 of file G4OpBoundaryProcess.hh.

◆ Material1

G4Material* G4OpBoundaryProcess::Material1
private

Definition at line 224 of file G4OpBoundaryProcess.hh.

◆ Material2

G4Material* G4OpBoundaryProcess::Material2
private

Definition at line 225 of file G4OpBoundaryProcess.hh.

◆ NewMomentum

G4ThreeVector G4OpBoundaryProcess::NewMomentum
private

Definition at line 218 of file G4OpBoundaryProcess.hh.

◆ NewPolarization

G4ThreeVector G4OpBoundaryProcess::NewPolarization
private

Definition at line 219 of file G4OpBoundaryProcess.hh.

◆ OldMomentum

G4ThreeVector G4OpBoundaryProcess::OldMomentum
private

Definition at line 215 of file G4OpBoundaryProcess.hh.

◆ OldPolarization

G4ThreeVector G4OpBoundaryProcess::OldPolarization
private

Definition at line 216 of file G4OpBoundaryProcess.hh.

◆ OpticalSurface

G4OpticalSurface* G4OpBoundaryProcess::OpticalSurface
private

Definition at line 227 of file G4OpBoundaryProcess.hh.

◆ prob_bs

G4double G4OpBoundaryProcess::prob_bs
private

Definition at line 250 of file G4OpBoundaryProcess.hh.

◆ prob_sl

G4double G4OpBoundaryProcess::prob_sl
private

Definition at line 250 of file G4OpBoundaryProcess.hh.

◆ prob_ss

G4double G4OpBoundaryProcess::prob_ss
private

Definition at line 250 of file G4OpBoundaryProcess.hh.

◆ PropertyPointer

G4MaterialPropertyVector* G4OpBoundaryProcess::PropertyPointer
private

Definition at line 229 of file G4OpBoundaryProcess.hh.

◆ PropertyPointer1

G4MaterialPropertyVector* G4OpBoundaryProcess::PropertyPointer1
private

Definition at line 230 of file G4OpBoundaryProcess.hh.

◆ PropertyPointer2

G4MaterialPropertyVector* G4OpBoundaryProcess::PropertyPointer2
private

Definition at line 231 of file G4OpBoundaryProcess.hh.

◆ Rindex1

G4double G4OpBoundaryProcess::Rindex1
private

Definition at line 233 of file G4OpBoundaryProcess.hh.

◆ Rindex2

G4double G4OpBoundaryProcess::Rindex2
private

Definition at line 234 of file G4OpBoundaryProcess.hh.

◆ sint1

G4double G4OpBoundaryProcess::sint1
private

Definition at line 236 of file G4OpBoundaryProcess.hh.

◆ sint2

G4double G4OpBoundaryProcess::sint2
private

Definition at line 236 of file G4OpBoundaryProcess.hh.

◆ theEfficiency

G4double G4OpBoundaryProcess::theEfficiency
private

Definition at line 245 of file G4OpBoundaryProcess.hh.

◆ theFacetNormal

G4ThreeVector G4OpBoundaryProcess::theFacetNormal
private

Definition at line 222 of file G4OpBoundaryProcess.hh.

◆ theFinish

G4OpticalSurfaceFinish G4OpBoundaryProcess::theFinish
private

Definition at line 242 of file G4OpBoundaryProcess.hh.

◆ theGlobalNormal

G4ThreeVector G4OpBoundaryProcess::theGlobalNormal
private

Definition at line 221 of file G4OpBoundaryProcess.hh.

◆ theModel

G4OpticalSurfaceModel G4OpBoundaryProcess::theModel
private

Definition at line 240 of file G4OpBoundaryProcess.hh.

◆ thePhotonMomentum

G4double G4OpBoundaryProcess::thePhotonMomentum
private

Definition at line 213 of file G4OpBoundaryProcess.hh.

◆ theReflectivity

G4double G4OpBoundaryProcess::theReflectivity
private

Definition at line 244 of file G4OpBoundaryProcess.hh.

◆ theStatus

G4OpBoundaryProcessStatus G4OpBoundaryProcess::theStatus
private

Definition at line 238 of file G4OpBoundaryProcess.hh.

◆ theSurfaceRoughness

G4double G4OpBoundaryProcess::theSurfaceRoughness
private

Definition at line 248 of file G4OpBoundaryProcess.hh.

◆ theTransmittance

G4double G4OpBoundaryProcess::theTransmittance
private

Definition at line 246 of file G4OpBoundaryProcess.hh.


The documentation for this class was generated from the following files: