OpBoundaryProcessSimple.cxx
Go to the documentation of this file.
1 // Adapted for LArSoft by Ben Jones, MIT, March 2010
2 //
3 // This class invokes a simplified model of optical reflections at
4 // boundaries between different materials. The relevant reflectivities
5 // are read in by an implementation of the MaterialPropertiesLoader.
6 //
7 // The required parameters are total reflectance and ratio of diffuse
8 // to specular reflectance. Each photon crossing a boundary with a
9 // defined reflectance is randomly either reflected or absorbed and killed
10 // according to the supplied probability.
11 //
12 // Every reflected photon with a defined diffuse reflection fraction
13 // is then randomly either diffusely or specularly reflected according
14 // to the supplied probability. All materials with no defined
15 // reflectance are assumed to be black and
16 // absorb all incident photons.
17 //
18 // This physics process is loaded in the OpticalPhysics physics constructor.
19 //
20 // This class is based on the G4OpBoundaryProcess class in Geant4
21 
22 //
23 // ********************************************************************
24 // * License and Disclaimer *
25 // * *
26 // * The Geant4 software is copyright of the Copyright Holders of *
27 // * the Geant4 Collaboration. It is provided under the terms and *
28 // * conditions of the Geant4 Software License, included in the file *
29 // * LICENSE and available at http://cern.ch/geant4/license . These *
30 // * include a list of copyright holders. *
31 // * *
32 // * Neither the authors of this software system, nor their employing *
33 // * institutes,nor the agencies providing financial support for this *
34 // * work make any representation or warranty, express or implied, *
35 // * regarding this software system or assume any liability for its *
36 // * use. Please see the license in the file LICENSE and URL above *
37 // * for the full disclaimer and the limitation of liability. *
38 // * *
39 // * This code implementation is the result of the scientific and *
40 // * technical work of the GEANT4 collaboration. *
41 // * By using, copying, modifying or distributing the software (or *
42 // * any work based on the software) you agree to acknowledge its *
43 // * use in resulting scientific publications, and indicate your *
44 // * acceptance of all terms of the Geant4 Software license. *
45 // ********************************************************************
46 //
47 ////////////////////////////////////////////////////////////////////////
48 // Optical Photon Boundary Process Class Implementation
49 ////////////////////////////////////////////////////////////////////////
50 
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"
57 
60 
61 //#define G4DEBUG_OPTICAL
62 
63 namespace larg4 {
64 
65  OpBoundaryProcessSimple::OpBoundaryProcessSimple(const G4String& processName, G4ProcessType type)
66  : G4VDiscreteProcess(processName, type)
67  {
68  G4cout << GetProcessName() << " is created " << G4endl;
69 
70  SetProcessSubType(fOpBoundary);
71 
73 
74  fCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
75  }
76 
77  //Action to take after making each step of an optical photon - described in file header.
78 
79  // PostStepDoIt
80  G4VParticleChange*
81  OpBoundaryProcessSimple::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
82  {
83 
84  // Note - these have to be loaded here, since this object
85  // is constructed before LArG4, and hence before the
86  // LArG4 parameters are read from config
87 
89  aParticleChange.Initialize(aTrack);
90 
91  G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
92  G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
93 
94  if (pPostStepPoint->GetStepStatus() != fGeomBoundary) {
96  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
97  }
98  if (aTrack.GetStepLength() <= fCarTolerance / 2) {
100  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
101  }
102 
103  G4Material* Material1 = pPreStepPoint->GetMaterial();
104  G4Material* Material2 = pPostStepPoint->GetMaterial();
105 
106  if (Material1 != Material2) {
107 
109 
110  fVerbosity = lgp->OpVerbosity();
111 
112  G4ThreeVector NewMomentum;
113  G4ThreeVector NewPolarization;
114 
115  G4ThreeVector theGlobalNormal;
116  G4ThreeVector theFacetNormal;
117 
118  std::string Material1Name = Material1->GetName();
119  std::string Material2Name = Material2->GetName();
120 
121  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
122  G4double thePhotonMomentum = aParticle->GetTotalMomentum();
123  G4ThreeVector OldMomentum = aParticle->GetMomentumDirection();
124  G4ThreeVector OldPolarization = aParticle->GetPolarization();
125 
126  if (fVerbosity > 9)
127  std::cout << "OpBoundaryProcessSimple Debug: Photon " << aTrack.GetTrackID() << " momentum "
128  << thePhotonMomentum << " Material1 " << Material1Name.c_str() << " Material 2 "
129  << Material2Name.c_str() << std::endl;
130 
131  G4ThreeVector theGlobalPoint = pPostStepPoint->GetPosition();
132 
133  G4Navigator* theNavigator =
134  G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking();
135 
136  G4ThreeVector theLocalPoint =
137  theNavigator->GetGlobalToLocalTransform().TransformPoint(theGlobalPoint);
138 
139  // Normal points back into volume
140  G4bool valid;
141  G4ThreeVector theLocalNormal = theNavigator->GetLocalExitNormal(&valid);
142 
143  if (valid) { theLocalNormal = -theLocalNormal; }
144  else {
145  G4cerr << " OpBoundaryProcessSimple/PostStepDoIt(): "
146  << " The Navigator reports that it returned an invalid normal" << G4endl;
147  }
148 
149  theGlobalNormal = theNavigator->GetLocalToGlobalTransform().TransformAxis(theLocalNormal);
150 
151  if (OldMomentum * theGlobalNormal > 0.0) {
152 #ifdef G4DEBUG_OPTICAL
153  G4cerr << " OpBoundaryProcessSimple/PostStepDoIt(): "
154  << " theGlobalNormal points the wrong direction " << G4endl;
155 #endif
156  theGlobalNormal = -theGlobalNormal;
157  }
158 
159  G4MaterialPropertiesTable* aMaterialPropertiesTable;
160  G4MaterialPropertyVector* Reflectance;
161  G4MaterialPropertyVector* DiffuseReflectanceFraction;
162 
163  aMaterialPropertiesTable = Material1->GetMaterialPropertiesTable();
164 
165  if (aMaterialPropertiesTable) {
166 
167  std::stringstream PropertyName;
168  PropertyName << "REFLECTANCE_" << Material2Name.c_str();
169  Reflectance = aMaterialPropertiesTable->GetProperty(PropertyName.str().c_str());
170 
171  PropertyName.str("");
172  PropertyName << "DIFFUSE_REFLECTANCE_FRACTION_" << Material2Name.c_str();
173  DiffuseReflectanceFraction =
174  aMaterialPropertiesTable->GetProperty(PropertyName.str().c_str());
175 
176  if (Reflectance) {
177  double theReflectance = Reflectance->Value(thePhotonMomentum);
178  if (G4BooleanRand(theReflectance)) {
179  if (DiffuseReflectanceFraction) {
180  double theDiffuseReflectanceFraction =
181  DiffuseReflectanceFraction->Value(thePhotonMomentum);
182  if (G4BooleanRand(theDiffuseReflectanceFraction)) { fTheStatus = SimpleDiffuse; }
183  else {
185  }
186  }
187  else {
189  }
190  }
191  else {
193  }
194  }
195  else {
197  }
198  }
199 
200  // Take action according to track status
201 
202  if (fTheStatus == SimpleDiffuse) {
203  NewMomentum = G4LambertianRand(theGlobalNormal);
204  theFacetNormal = (NewMomentum - OldMomentum).unit();
205  }
206 
207  else if (fTheStatus == SimpleSpecular) {
208  theFacetNormal = theGlobalNormal;
209  G4double PdotN = OldMomentum * theFacetNormal;
210  NewMomentum = OldMomentum - (2. * PdotN) * theFacetNormal;
211  }
212 
214  aParticleChange.ProposeTrackStatus(fStopAndKill);
215  }
216 
217  else if ((fTheStatus == NotAtBoundary) || (fTheStatus == StepTooSmall)) {
218  NewMomentum = OldMomentum;
219  }
220  else
221  aParticleChange.ProposeTrackStatus(fStopAndKill);
222 
223  NewMomentum = NewMomentum.unit();
224  aParticleChange.ProposeMomentumDirection(NewMomentum);
225  }
226  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
227  }
228 
229  // Compulsary method for G4DiscreteProcesses. Serves no actual function here.
230 
231  // GetMeanFreePath
232  G4double
233  OpBoundaryProcessSimple::GetMeanFreePath(const G4Track&, G4double, G4ForceCondition* condition)
234  {
235  *condition = Forced;
236 
237  return DBL_MAX;
238  }
239 }
Store parameters for running LArG4.
std::string string
Definition: nybbler.cc:12
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
Geant4 interface.
OpBoundaryProcessSimpleStatus fTheStatus
int OpVerbosity() const
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)