G4StationSimulator/G4TankFastCerenkov.h
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * DISCLAIMER *
4 // * *
5 // * The following disclaimer summarizes all the specific disclaimers *
6 // * of contributors to this software. The specific disclaimers,which *
7 // * govern, are listed with their locations in: *
8 // * http://cern.ch/geant4/license *
9 // * *
10 // * Neither the authors of this software system, nor their employing *
11 // * institutes,nor the agencies providing financial support for this *
12 // * work make any representation or warranty, express or implied, *
13 // * regarding this software system or assume any liability for its *
14 // * use. *
15 // * *
16 // * This code implementation is the intellectual property of the *
17 // * GEANT4 collaboration. *
18 // * By copying, distributing or modifying the Program (or any work *
19 // * based on the Program) you indicate your acceptance of this *
20 // * statement, and all its terms. *
21 // ********************************************************************
22 
23 #ifndef _G4StationSimulatorOG_G4TankFastCerenkov_h_
24 #define _G4StationSimulatorOG_G4TankFastCerenkov_h_
25 
26 #include <sdet/PMT.h>
27 #include <utl/Point.h>
28 #include <sevt/Station.h>
29 
30 #include <globals.hh>
31 #include <templates.hh>
32 #include <Randomize.hh>
33 #include <G4ThreeVector.hh>
34 #include <G4ParticleMomentum.hh>
35 #include <G4Step.hh>
36 #include <G4VContinuousProcess.hh>
37 #include <G4OpticalPhoton.hh>
38 #include <G4DynamicParticle.hh>
39 #include <G4Material.hh>
40 #include <G4PhysicsTable.hh>
41 #include <G4MaterialPropertiesTable.hh>
42 #include <G4PhysicsOrderedFreeVector.hh>
43 #include <G4ios.hh>
44 #include <G4Poisson.hh>
45 
46 
47 namespace G4StationSimulatorOG {
48 
49  class G4StationSimulator;
50 
51 
52  class G4TankFastCerenkov : public G4VContinuousProcess {
53 
54  public:
55  G4TankFastCerenkov(const G4String& processName = "Cerenkov");
56  virtual ~G4TankFastCerenkov();
57 
58  static void GetDataFromConstruction();
59 
60  G4bool IsApplicable(const G4ParticleDefinition& particle) override
61  { return particle.GetPDGCharge(); }
62 
63  void SetTrackSecondariesFirst(const G4bool state)
64  { fTrackSecondariesFirst = state; }
65 
66  void SetMaxNumPhotonsPerStep(const G4int numPhotons)
67  { fMaxPhotons = numPhotons; }
68 
69  void DumpPhysicsTable() const;
70 
71  G4VParticleChange* AlongStepDoIt(const G4Track& track, const G4Step& step) override;
72 
73  G4double GetContinuousStepLimit(const G4Track&, G4double, G4double, G4double&) override;
74 
75  private:
76  void BuildThePhysicsTable();
77  G4double GetAverageNumberOfPhotons(const G4DynamicParticle*,
78  const G4Material*,
79  G4MaterialPropertyVector*) const;
80  void TrackPhotonInTank();
81  void PropagateInTank(G4int flag);
82 
83  template<int pmtId>
84  G4bool TransitionToDome(G4int flag);
85 
86  template<int pmtId>
87  G4bool PropagateInDome(G4int& flag);
88 
89  template<int pmtId>
90  G4bool TransitionToInterface(G4int& flag);
91 
92  template<int pmtId>
93  G4bool PropagateInInterface(G4int& flag);
94 
95  G4double GetSphereIntersect(const G4ThreeVector&, const G4double) const;
96 
97  G4double GetEllipsoidIntersect(const G4ThreeVector&,
98  const G4double, const G4double) const;
99 
100  void DoRayleighScatter();
101  G4bool ScatterOffRoof();
102  G4bool ScatterOffWall();
103  G4bool ScatterOffFloor();
104  void BackScatter();
105  void DiffuseScatterHorizontal(const G4ThreeVector& normal);
106  void DiffuseScatterVertical(const G4ThreeVector& normal);
107  void SpikeScatter(const G4ThreeVector& normal);
108  void LobeScatterHorizontal(const G4ThreeVector& normal);
109  void LobeScatterVertical(const G4ThreeVector& normal);
110 
111  G4bool fTrackSecondariesFirst = false;
112  G4int fMaxPhotons = 0;
113  G4ThreeVector fPhotonMomentum;
114  G4ThreeVector fPhotonPolarization;
115  G4ThreeVector fPhotonPosition;
116  G4double fPhotonTime;
118 
119  G4double fInvMFP;
120  G4double fRayleighFrac;
121  G4double fAbsorbFrac;
122 
123  static double fTankHeight;
124  static double fTankHalfHeight;
125  static double fTankThickness;
126  static double fTankRadius;
127  static double fTankRadiusSq;
128  static double fDomeRadius;
129  static double fInterfaceRadius;
130  static double fPMTRadius;
131  static double fHeightz;
132  static double fDomeRadiusz;
133  static double fInterfaceRadiusz;
134  static double fPMTRadiusz;
135  static double fDomeRadiusSq;
136  static double fInterfaceRadiusSq;
137  static double fPMTRadiusSq;
138  static double fDomeRadiuszSq;
139  static double fInterfaceRadiuszSq;
140  static double fPMTRadiuszSq;
141 
142  static double fSigmaAlpha;
143 
145 
149 
156 
157  // PMT positions (the zeroth position is dummy)
158  static G4ThreeVector fPMTPos[4];
159  // The height of the roof, referred to each PMT
160  static G4double fRoofPos[4];
161 
162  G4double fSampledRI;
163  G4double fRayleighConst;
164 
166 
167  protected:
168  G4PhysicsTable* fPhysicsTable = nullptr;
169 
170  };
171 
172 }
173 
174 
175 #endif
Class to hold collection (x,y) points and provide interpolation between them.
G4double GetSphereIntersect(const G4ThreeVector &, const G4double) const
G4bool IsApplicable(const G4ParticleDefinition &particle) override
class that handles Geant4 SD Station simulation adopted from G4TankSimulator
G4TankFastCerenkov(const G4String &processName="Cerenkov")
G4double GetAverageNumberOfPhotons(const G4DynamicParticle *, const G4Material *, G4MaterialPropertyVector *) const
struct particle_info particle[80]
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step) override
G4double GetContinuousStepLimit(const G4Track &, G4double, G4double, G4double &) override
G4double GetEllipsoidIntersect(const G4ThreeVector &, const G4double, const G4double) const

, generated on Tue Sep 26 2023.