LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
OpBoundaryProcessSimple.hh
Go to the documentation of this file.
1 
9 // ********************************************************************
10 // * License and Disclaimer *
11 // * *
12 // * The Geant4 software is copyright of the Copyright Holders of *
13 // * the Geant4 Collaboration. It is provided under the terms and *
14 // * conditions of the Geant4 Software License, included in the file *
15 // * LICENSE and available at http://cern.ch/geant4/license . These *
16 // * include a list of copyright holders. *
17 // * *
18 // * Neither the authors of this software system, nor their employing *
19 // * institutes,nor the agencies providing financial support for this *
20 // * work make any representation or warranty, express or implied, *
21 // * regarding this software system or assume any liability for its *
22 // * use. Please see the license in the file LICENSE and URL above *
23 // * for the full disclaimer and the limitation of liability. *
24 // * *
25 // * This code implementation is the result of the scientific and *
26 // * technical work of the GEANT4 collaboration. *
27 // * By using, copying, modifying or distributing the software (or *
28 // * any work based on the software) you agree to acknowledge its *
29 // * use in resulting scientific publications, and indicate your *
30 // * acceptance of all terms of the Geant4 Software license. *
31 // ********************************************************************
32 //
33 //
34 // GEANT4 tag $Name: $
35 //
36 //
38 // Optical Photon Boundary Process Class Definition
40 
41 
42 
43 #ifndef OpBoundaryProcessSimple_h
44 #define OpBoundaryProcessSimple_h 1
45 
46 
47 #include "Geant4/globals.hh"
48 #include "Geant4/templates.hh"
49 #include "Geant4/geomdefs.hh"
50 #include "Geant4/Randomize.hh"
51 
52 #include "Geant4/G4RandomTools.hh"
53 #include "Geant4/G4RandomDirection.hh"
54 
55 #include "Geant4/G4Step.hh"
56 #include "Geant4/G4VDiscreteProcess.hh"
57 #include "Geant4/G4DynamicParticle.hh"
58 #include "Geant4/G4Material.hh"
59 #include "Geant4/G4LogicalBorderSurface.hh"
60 #include "Geant4/G4LogicalSkinSurface.hh"
61 #include "Geant4/G4OpticalSurface.hh"
62 #include "Geant4/G4OpticalPhoton.hh"
63 #include "Geant4/G4TransportationManager.hh"
64 
65 
66 #include "TH1.h"
67 #include "TTree.h"
68 
69 
70 
72 // Class Definition
74 
75 namespace larg4 {
76 
77  // Possible statuses of a particle after each step.
80 
108 class OpBoundaryProcessSimple : public G4VDiscreteProcess
109  {
110 
111 
112  public:
113 
114 
115  // Constructors and Destructor
116 
117  OpBoundaryProcessSimple(const G4String& processName = "OpBoundary",
118  G4ProcessType type = fOptical);
119 
121 
122 
123 
124  public:
125 
126  G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
127  // Returns true -> 'is applicable' only for an optical photon.
128 
129  G4double GetMeanFreePath(const G4Track& ,
130  G4double ,
131  G4ForceCondition* condition);
132  // Returns infinity; i. e. the process does not limit the step,
133  // but sets the 'Forced' condition for the DoIt to be invoked at
134  // every step. However, only at a boundary will any action be
135  // taken.
136 
137  G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
138  const G4Step& aStep);
139  // This is the method implementing boundary processes.
140 
142  // Returns the current status.
143 
144 
145  private:
146 
147  G4bool G4BooleanRand(const G4double prob) const;
148  // Generate a random bool to decide which process to execute
149 
150  private:
151 
153  G4double fCarTolerance;
154 
156 
157  };
158 
159 
160  // Inline methods
161 
162  inline
163  G4bool OpBoundaryProcessSimple::G4BooleanRand(const G4double prob) const
164  {
165  /* Returns a random boolean variable with the specified probability */
166  return (G4UniformRand() < prob);
167  }
168 
169  inline
170  G4bool OpBoundaryProcessSimple::IsApplicable(const G4ParticleDefinition&
171  aParticleType)
172  {
173  return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
174  }
175 
176  inline
178  {
179  return fTheStatus;
180  }
181 
182 }
183 
184 #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)