LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
larg4::OpBoundaryProcessSimple Class Reference

Discrete process for reflection and diffusion at optical interfaces. More...

#include "OpBoundaryProcessSimple.hh"

Inheritance diagram for larg4::OpBoundaryProcessSimple:

Public Member Functions

 OpBoundaryProcessSimple (const G4String &processName="OpBoundary", G4ProcessType type=fOptical)
 
 ~OpBoundaryProcessSimple ()
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *condition)
 
G4VParticleChange * PostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
OpBoundaryProcessSimpleStatus GetStatus () const
 

Private Member Functions

G4bool G4BooleanRand (const G4double prob) const
 

Private Attributes

OpBoundaryProcessSimpleStatus fTheStatus
 
G4double fCarTolerance
 
int fVerbosity
 

Detailed Description

Discrete process for reflection and diffusion at optical interfaces.

See also
G4VDiscreteProcess

This class invokes a simplified model of optical reflections at boundaries between different materials. The relevant reflectivities are ultimately read from detinfo::LArProperties via larg4::MaterialPropertiesLoader.

The required parameters are total reflectance (detinfo::LArProperties::SurfaceReflectances()) and ratio of diffuse to specular reflectance (detinfo::LArProperties::SurfaceReflectanceDiffuseFractions()). Each photon crossing a boundary with a defined reflectance is randomly either reflected or absorbed and killed according to the supplied probability.

Every reflected photon with a defined diffuse reflection fraction is then randomly either diffusely or specularly reflected according to the supplied probability. All materials with no defined reflectance are assumed to be black and absorb all incident photons.

This physics process is loaded in larg4::OpticalPhysics physics constructor.

This class is based on the G4OpBoundaryProcess class in Geant4 and was adapted for LArSoft by Ben Jones, MIT, March 2010.

Definition at line 108 of file OpBoundaryProcessSimple.hh.

Constructor & Destructor Documentation

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

Definition at line 73 of file OpBoundaryProcessSimple.cxx.

References fCarTolerance, fTheStatus, and larg4::Undefined.

75  : G4VDiscreteProcess(processName, type)
76  {
77  if ( 1 > 0) {
78  G4cout << GetProcessName() << " is created " << G4endl;
79  }
80 
81 
82  SetProcessSubType(fOpBoundary);
83 
85 
86  fCarTolerance = G4GeometryTolerance::GetInstance()
87  ->GetSurfaceTolerance();
88 
89  }
OpBoundaryProcessSimpleStatus fTheStatus
larg4::OpBoundaryProcessSimple::~OpBoundaryProcessSimple ( )

Definition at line 98 of file OpBoundaryProcessSimple.cxx.

99  {
100  }

Member Function Documentation

G4bool larg4::OpBoundaryProcessSimple::G4BooleanRand ( const G4double  prob) const
inlineprivate

Definition at line 163 of file OpBoundaryProcessSimple.hh.

Referenced by PostStepDoIt().

164  {
165  /* Returns a random boolean variable with the specified probability */
166  return (G4UniformRand() < prob);
167  }
G4double larg4::OpBoundaryProcessSimple::GetMeanFreePath ( const G4Track &  ,
G4double  ,
G4ForceCondition *  condition 
)

Definition at line 297 of file OpBoundaryProcessSimple.cxx.

300  {
301  *condition = Forced;
302 
303  return DBL_MAX;
304  }
OpBoundaryProcessSimpleStatus larg4::OpBoundaryProcessSimple::GetStatus ( ) const
inline

Definition at line 177 of file OpBoundaryProcessSimple.hh.

References fTheStatus.

178  {
179  return fTheStatus;
180  }
OpBoundaryProcessSimpleStatus fTheStatus
G4bool larg4::OpBoundaryProcessSimple::IsApplicable ( const G4ParticleDefinition &  aParticleType)
inline

Definition at line 170 of file OpBoundaryProcessSimple.hh.

172  {
173  return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
174  }
G4VParticleChange * larg4::OpBoundaryProcessSimple::PostStepDoIt ( const G4Track &  aTrack,
const G4Step &  aStep 
)

Definition at line 109 of file OpBoundaryProcessSimple.cxx.

References fCarTolerance, fTheStatus, fVerbosity, G4BooleanRand(), larg4::NotAtBoundary, sim::LArG4Parameters::OpVerbosity(), larg4::SimpleAbsorbed, larg4::SimpleAbsorbedNoRefl, larg4::SimpleDiffuse, larg4::SimpleSpecular, larg4::StepTooSmall, and larg4::Undefined.

110  {
111 
112  // Note - these have to be loaded here, since this object
113  // is constructed before LArG4, and hence before the
114  // LArG4 parameters are read from config
115 
117  aParticleChange.Initialize(aTrack);
118 
119 
120  G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
121  G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
122 
123 
124  if (pPostStepPoint->GetStepStatus() != fGeomBoundary){
126  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
127  }
128  if (aTrack.GetStepLength()<=fCarTolerance/2){
130  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
131  }
132 
133  G4Material* Material1 = pPreStepPoint -> GetMaterial();
134  G4Material* Material2 = pPostStepPoint -> GetMaterial();
135 
136 
137  if(Material1 != Material2)
138  {
139 
141 
142  fVerbosity = lgp->OpVerbosity();
143 
144 
145 
146  G4double thePhotonMomentum;
147 
148  G4ThreeVector OldMomentum;
149  G4ThreeVector OldPolarization;
150 
151  G4ThreeVector NewMomentum;
152  G4ThreeVector NewPolarization;
153 
154  G4ThreeVector theGlobalNormal;
155  G4ThreeVector theFacetNormal;
156 
157  std::string Material1Name = Material1->GetName();
158  std::string Material2Name = Material2->GetName();
159 
160  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
161 
162 
163  thePhotonMomentum = aParticle->GetTotalMomentum();
164  OldMomentum = aParticle->GetMomentumDirection();
165  OldPolarization = aParticle->GetPolarization();
166 
167 
168 
169  if (fVerbosity>9) std::cout<<"OpBoundaryProcessSimple Debug: Photon " <<aTrack.GetTrackID() <<" momentum "<<thePhotonMomentum << " Material1 " << Material1Name.c_str() << " Material 2 " << Material2Name.c_str() << std::endl;
170 
171  G4ThreeVector theGlobalPoint = pPostStepPoint->GetPosition();
172 
173  G4Navigator* theNavigator =
174  G4TransportationManager::GetTransportationManager()->
175  GetNavigatorForTracking();
176 
177  G4ThreeVector theLocalPoint = theNavigator->
178  GetGlobalToLocalTransform().
179  TransformPoint(theGlobalPoint);
180 
181  G4ThreeVector theLocalNormal; // Normal points back into volume
182 
183  G4bool valid;
184  theLocalNormal = theNavigator->GetLocalExitNormal(&valid);
185 
186  if (valid) {
187  theLocalNormal = -theLocalNormal;
188  }
189  else {
190  G4cerr << " OpBoundaryProcessSimple/PostStepDoIt(): "
191  << " The Navigator reports that it returned an invalid normal"
192  << G4endl;
193  }
194 
195  theGlobalNormal = theNavigator->GetLocalToGlobalTransform().
196  TransformAxis(theLocalNormal);
197 
198  if (OldMomentum * theGlobalNormal > 0.0) {
199 #ifdef G4DEBUG_OPTICAL
200  G4cerr << " OpBoundaryProcessSimple/PostStepDoIt(): "
201  << " theGlobalNormal points the wrong direction "
202  << G4endl;
203 #endif
204  theGlobalNormal = -theGlobalNormal;
205  }
206 
207  G4MaterialPropertiesTable* aMaterialPropertiesTable;
208  G4MaterialPropertyVector* Reflectance;
209  G4MaterialPropertyVector* DiffuseReflectanceFraction;
210 
211 
212  aMaterialPropertiesTable = Material1->GetMaterialPropertiesTable();
213 
214 
215  if (aMaterialPropertiesTable) {
216 
217  std::stringstream PropertyName;
218  PropertyName<<"REFLECTANCE_"<<Material2Name.c_str();
219  Reflectance = aMaterialPropertiesTable->GetProperty(PropertyName.str().c_str());
220 
221  PropertyName.str("");
222  PropertyName<<"DIFFUSE_REFLECTANCE_FRACTION_"<<Material2Name.c_str();
223  DiffuseReflectanceFraction=aMaterialPropertiesTable->GetProperty(PropertyName.str().c_str());
224 
225  if (Reflectance)
226  {
227  double theReflectance = Reflectance->Value(thePhotonMomentum);
228  if( G4BooleanRand(theReflectance))
229  {
230  if(DiffuseReflectanceFraction)
231  {
232  double theDiffuseReflectanceFraction = DiffuseReflectanceFraction->Value(thePhotonMomentum);
233  if(G4BooleanRand(theDiffuseReflectanceFraction))
234  {
236  }
237  else
238  {
240  }
241  }
242  else
243  {
245  }
246  }
247  else
248  {
250  }
251 
252  }
253  else
254  {
256  }
257  }
258 
259 
260 
261 
262  // Take action according to track status
263 
264  if ( fTheStatus == SimpleDiffuse ) {
265  NewMomentum = G4LambertianRand(theGlobalNormal);
266  theFacetNormal = (NewMomentum - OldMomentum).unit();
267  }
268 
269  else if (fTheStatus==SimpleSpecular) {
270  theFacetNormal = theGlobalNormal;
271  G4double PdotN = OldMomentum * theFacetNormal;
272  NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
273  }
274 
276  aParticleChange.ProposeTrackStatus(fStopAndKill);
277  }
278 
280  {
281  NewMomentum=OldMomentum;
282  }
283  else aParticleChange.ProposeTrackStatus(fStopAndKill);
284 
285  NewMomentum = NewMomentum.unit();
286  aParticleChange.ProposeMomentumDirection(NewMomentum);
287  }
288  return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
289 
290  }
OpBoundaryProcessSimpleStatus fTheStatus
int OpVerbosity() const
G4bool G4BooleanRand(const G4double prob) const

Member Data Documentation

G4double larg4::OpBoundaryProcessSimple::fCarTolerance
private

Definition at line 153 of file OpBoundaryProcessSimple.hh.

Referenced by OpBoundaryProcessSimple(), and PostStepDoIt().

OpBoundaryProcessSimpleStatus larg4::OpBoundaryProcessSimple::fTheStatus
private

Definition at line 152 of file OpBoundaryProcessSimple.hh.

Referenced by GetStatus(), OpBoundaryProcessSimple(), and PostStepDoIt().

int larg4::OpBoundaryProcessSimple::fVerbosity
private

Definition at line 155 of file OpBoundaryProcessSimple.hh.

Referenced by PostStepDoIt().


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