51 #include "Geant4/G4GeometryTolerance.hh" 52 #include "Geant4/G4Navigator.hh" 53 #include "Geant4/G4OpProcessSubType.hh" 54 #include "Geant4/G4RandomTools.hh" 55 #include "Geant4/G4TransportationManager.hh" 56 #include "Geant4/G4ios.hh" 66 : G4VDiscreteProcess(processName, type)
68 G4cout << GetProcessName() <<
" is created " << G4endl;
70 SetProcessSubType(fOpBoundary);
74 fCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
89 aParticleChange.Initialize(aTrack);
91 G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
92 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
94 if (pPostStepPoint->GetStepStatus() != fGeomBoundary) {
96 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
100 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
103 G4Material* Material1 = pPreStepPoint->GetMaterial();
104 G4Material* Material2 = pPostStepPoint->GetMaterial();
106 if (Material1 != Material2) {
112 G4ThreeVector NewMomentum;
113 G4ThreeVector NewPolarization;
115 G4ThreeVector theGlobalNormal;
116 G4ThreeVector theFacetNormal;
121 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
122 G4double thePhotonMomentum = aParticle->GetTotalMomentum();
123 G4ThreeVector OldMomentum = aParticle->GetMomentumDirection();
124 G4ThreeVector OldPolarization = aParticle->GetPolarization();
127 std::cout <<
"OpBoundaryProcessSimple Debug: Photon " << aTrack.GetTrackID() <<
" momentum " 128 << thePhotonMomentum <<
" Material1 " << Material1Name.c_str() <<
" Material 2 " 131 G4ThreeVector theGlobalPoint = pPostStepPoint->GetPosition();
133 G4Navigator* theNavigator =
134 G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking();
136 G4ThreeVector theLocalPoint =
137 theNavigator->GetGlobalToLocalTransform().TransformPoint(theGlobalPoint);
141 G4ThreeVector theLocalNormal = theNavigator->GetLocalExitNormal(&valid);
143 if (valid) { theLocalNormal = -theLocalNormal; }
145 G4cerr <<
" OpBoundaryProcessSimple/PostStepDoIt(): " 146 <<
" The Navigator reports that it returned an invalid normal" << G4endl;
149 theGlobalNormal = theNavigator->GetLocalToGlobalTransform().TransformAxis(theLocalNormal);
151 if (OldMomentum * theGlobalNormal > 0.0) {
152 #ifdef G4DEBUG_OPTICAL 153 G4cerr <<
" OpBoundaryProcessSimple/PostStepDoIt(): " 154 <<
" theGlobalNormal points the wrong direction " << G4endl;
156 theGlobalNormal = -theGlobalNormal;
159 G4MaterialPropertiesTable* aMaterialPropertiesTable;
160 G4MaterialPropertyVector* Reflectance;
161 G4MaterialPropertyVector* DiffuseReflectanceFraction;
163 aMaterialPropertiesTable = Material1->GetMaterialPropertiesTable();
165 if (aMaterialPropertiesTable) {
167 std::stringstream PropertyName;
168 PropertyName <<
"REFLECTANCE_" << Material2Name.c_str();
169 Reflectance = aMaterialPropertiesTable->GetProperty(PropertyName.str().c_str());
171 PropertyName.str(
"");
172 PropertyName <<
"DIFFUSE_REFLECTANCE_FRACTION_" << Material2Name.c_str();
173 DiffuseReflectanceFraction =
174 aMaterialPropertiesTable->GetProperty(PropertyName.str().c_str());
177 double theReflectance = Reflectance->Value(thePhotonMomentum);
179 if (DiffuseReflectanceFraction) {
180 double theDiffuseReflectanceFraction =
181 DiffuseReflectanceFraction->Value(thePhotonMomentum);
203 NewMomentum = G4LambertianRand(theGlobalNormal);
204 theFacetNormal = (NewMomentum - OldMomentum).unit();
208 theFacetNormal = theGlobalNormal;
209 G4double PdotN = OldMomentum * theFacetNormal;
210 NewMomentum = OldMomentum - (2. * PdotN) * theFacetNormal;
214 aParticleChange.ProposeTrackStatus(fStopAndKill);
218 NewMomentum = OldMomentum;
221 aParticleChange.ProposeTrackStatus(fStopAndKill);
223 NewMomentum = NewMomentum.unit();
224 aParticleChange.ProposeMomentumDirection(NewMomentum);
226 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
Store parameters for running LArG4.
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
OpBoundaryProcessSimpleStatus fTheStatus
G4bool G4BooleanRand(const G4double prob) const
G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *condition)
QTextStream & endl(QTextStream &s)
OpBoundaryProcessSimple(const G4String &processName="OpBoundary", G4ProcessType type=fOptical)