OpBoundaryProcessSimple.hh
Go to the documentation of this file.
1 /**
2  * @file larsim/LArG4/OpBoundaryProcessSimple.hh
3  * @author Ben Jones (bjpjones@mit.edu)
4  * @date March 2010
5  * @see larsim/LArG4/OpBoundaryProcessSimple.cxx
6  */
7 
8 // ********************************************************************
9 // * License and Disclaimer *
10 // * *
11 // * The Geant4 software is copyright of the Copyright Holders of *
12 // * the Geant4 Collaboration. It is provided under the terms and *
13 // * conditions of the Geant4 Software License, included in the file *
14 // * LICENSE and available at http://cern.ch/geant4/license . These *
15 // * include a list of copyright holders. *
16 // * *
17 // * Neither the authors of this software system, nor their employing *
18 // * institutes,nor the agencies providing financial support for this *
19 // * work make any representation or warranty, express or implied, *
20 // * regarding this software system or assume any liability for its *
21 // * use. Please see the license in the file LICENSE and URL above *
22 // * for the full disclaimer and the limitation of liability. *
23 // * *
24 // * This code implementation is the result of the scientific and *
25 // * technical work of the GEANT4 collaboration. *
26 // * By using, copying, modifying or distributing the software (or *
27 // * any work based on the software) you agree to acknowledge its *
28 // * use in resulting scientific publications, and indicate your *
29 // * acceptance of all terms of the Geant4 Software license. *
30 // ********************************************************************
31 //
32 //
33 // GEANT4 tag $Name: $
34 //
35 //
36 ////////////////////////////////////////////////////////////////////////
37 // Optical Photon Boundary Process Class Definition
38 ////////////////////////////////////////////////////////////////////////
39 
40 #ifndef OpBoundaryProcessSimple_h
41 #define OpBoundaryProcessSimple_h 1
42 
43 #include "Geant4/G4ForceCondition.hh"
44 #include "Geant4/G4OpticalPhoton.hh"
45 #include "Geant4/G4ProcessType.hh"
46 #include "Geant4/G4String.hh"
47 #include "Geant4/G4Types.hh"
48 #include "Geant4/G4VDiscreteProcess.hh"
49 #include "Geant4/Randomize.hh"
50 
51 class G4ParticleDefinition;
52 class G4Step;
53 class G4Track;
54 class G4VParticleChange;
55 
56 /////////////////////
57 // Class Definition
58 /////////////////////
59 
60 namespace larg4 {
61 
62  // Possible statuses of a particle after each step.
72  };
73 
74  /**
75  * @brief Discrete process for reflection and diffusion at optical interfaces.
76  * @see `G4VDiscreteProcess`
77  *
78  * This class invokes a simplified model of optical reflections at
79  * boundaries between different materials. The relevant reflectivities
80  * are ultimately read from `detinfo::LArProperties` via
81  * `larg4::MaterialPropertiesLoader`.
82  *
83  * The required parameters are total reflectance
84  * (`detinfo::LArProperties::SurfaceReflectances()`)
85  * and ratio of diffuse to specular reflectance
86  * (`detinfo::LArProperties::SurfaceReflectanceDiffuseFractions()`). Each
87  * photon crossing a boundary with a defined reflectance is randomly either
88  * reflected or absorbed and killed according to the supplied probability.
89  *
90  * Every reflected photon with a defined diffuse reflection fraction
91  * is then randomly either diffusely or specularly reflected according
92  * to the supplied probability. All materials with no defined
93  * reflectance are assumed to be black and absorb all incident photons.
94  *
95  * This physics process is loaded in `larg4::OpticalPhysics` physics
96  * constructor.
97  *
98  * This class is based on the `G4OpBoundaryProcess` class in Geant4 and was
99  * adapted for LArSoft by Ben Jones, MIT, March 2010.
100  */
101 
102  class OpBoundaryProcessSimple : public G4VDiscreteProcess {
103  public:
104  OpBoundaryProcessSimple(const G4String& processName = "OpBoundary",
105  G4ProcessType type = fOptical);
106 
107  // Returns true -> 'is applicable' only for an optical photon.
108  G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
109 
110  // Returns infinity; i. e. the process does not limit the step,
111  // but sets the 'Forced' condition for the DoIt to be invoked at
112  // every step. However, only at a boundary will any action be
113  // taken.
114  G4double GetMeanFreePath(const G4Track&, G4double, G4ForceCondition* condition);
115 
116  // This is the method implementing boundary processes.
117  G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep);
118 
119  // Returns the current status.
121 
122  private:
123  // Generate a random bool to decide which process to execute
124  G4bool G4BooleanRand(const G4double prob) const;
125 
127  G4double fCarTolerance;
128 
130  };
131 
132  inline G4bool
133  OpBoundaryProcessSimple::G4BooleanRand(const G4double prob) const
134  {
135  /* Returns a random boolean variable with the specified probability */
136  return (G4UniformRand() < prob);
137  }
138 
139  inline G4bool
140  OpBoundaryProcessSimple::IsApplicable(const G4ParticleDefinition& aParticleType)
141  {
142  return (&aParticleType == G4OpticalPhoton::OpticalPhoton());
143  }
144 
147  {
148  return fTheStatus;
149  }
150 
151 }
152 
153 #endif /* OpBoundaryProcessSimple_h */
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
Geant4 interface.
OpBoundaryProcessSimpleStatus GetStatus() const
OpBoundaryProcessSimpleStatus fTheStatus
Discrete process for reflection and diffusion at optical interfaces.
G4bool G4BooleanRand(const G4double prob) const
G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *condition)
OpBoundaryProcessSimple(const G4String &processName="OpBoundary", G4ProcessType type=fOptical)