/home/docs/checkouts/readthedocs.org/user_builds/ratpac/checkouts/latest/src/physics/include/RAT/GLG4PMTOpticalModel.hh Source File

Ratpac-two: /home/docs/checkouts/readthedocs.org/user_builds/ratpac/checkouts/latest/src/physics/include/RAT/GLG4PMTOpticalModel.hh Source File
Ratpac-two
GLG4PMTOpticalModel.hh
Go to the documentation of this file.
1 
14 #ifndef __GLG4PMTOpticalModel_hh__
15 #define __GLG4PMTOpticalModel_hh__
16 
17 #include <RAT/Log.hh>
18 #include <map>
19 #include <utility>
20 #include <vector>
21 
22 #include "G4LogicalBorderSurface.hh"
23 #include "G4LogicalVolume.hh"
24 #include "G4MaterialPropertyVector.hh"
25 #include "G4OpticalSurface.hh"
26 #include "G4UImessenger.hh"
27 #include "G4VFastSimulationModel.hh"
28 
29 class G4UIcommand;
30 class G4UIdirectory;
31 class G4Region;
32 
33 G4complex carcsin(G4complex theta); // complex sin^-1
34 G4complex gfunc(G4complex ni, G4complex nj, G4complex ti, G4complex tj);
35 G4complex rfunc(G4complex ni, G4complex nj, G4complex ti, G4complex tj);
36 G4complex trfunc(G4complex ni, G4complex nj, G4complex ti, G4complex tj, G4complex tk);
37 
38 class GLG4PMTOpticalModel : public G4VFastSimulationModel, public G4UImessenger {
39  public:
40  //-------------------------
41  // Constructor, destructor
42  //-------------------------
43  // 28-Jul-2006 WGS: Must define a G4Region for Fast Simulations
44  // (change from Geant 4.7 to Geant 4.8).
45  GLG4PMTOpticalModel(G4String, G4Region *, G4LogicalVolume *body, G4OpticalSurface *op_surface,
46  double efficiency_correction = 1.0, double dynodeTop = 0.0, double dynodeRadius = 0.0,
47  double prepulseProb = 0.0, double photocathode_MINrho = 0.0, double photocathode_MAXrho = 0.0);
48  // Note: There is no GLG4PMTOpticalModel(G4String) constructor.
50 
51  //------------------------------
52  // Virtual methods of the base
53  // class to be coded by the user
54  //------------------------------
55 
56  G4bool IsApplicable(const G4ParticleDefinition &);
57  G4bool ModelTrigger(const G4FastTrack &);
58  void DoIt(const G4FastTrack &, G4FastStep &);
59 
60  // following two methods are for G4UImessenger
61  void SetNewValue(G4UIcommand *command, G4String newValues);
62  G4String GetCurrentValue(G4UIcommand *command);
63 
64  void SetEfficiencyCorrection(std::map<int, double> _EffiCorr) {
65  EfficiencyCorrection = _EffiCorr;
66  RAT::info << GetName() << ": Individual efficiency correction table set" << newline;
67  }
68  void DumpEfficiencyCorrectionTable() {
69  RAT::info << "Individual correction table for the PMT efficiencies of " << GetName() << ":\nPMT ID corr. factor"
70  << newline;
71  for (std::map<int, double>::iterator iter = EfficiencyCorrection.begin(); iter != EfficiencyCorrection.end();
72  iter++) {
73  RAT::info << iter->first << "," << iter->second << newline;
74  }
75  }
76  void fillPMTVector(double code, double A, double An, double T, double R, double collection_eff, double N_pe, double x,
77  double y, double z, double gx, double gy, double gz, double wavelength, double time) {
78  pmtHitVectorIndex.push_back(code);
79  pmtHitVectorIndex.push_back(A);
80  pmtHitVectorIndex.push_back(An);
81  pmtHitVectorIndex.push_back(T);
82  pmtHitVectorIndex.push_back(R);
83  pmtHitVectorIndex.push_back(collection_eff);
84  pmtHitVectorIndex.push_back(N_pe);
85  pmtHitVectorIndex.push_back(x);
86  pmtHitVectorIndex.push_back(y);
87  pmtHitVectorIndex.push_back(z);
88  pmtHitVectorIndex.push_back(gx);
89  pmtHitVectorIndex.push_back(gy);
90  pmtHitVectorIndex.push_back(gz);
91  pmtHitVectorIndex.push_back(wavelength);
92  pmtHitVectorIndex.push_back(time);
93  pmtHitVector.push_back(pmtHitVectorIndex);
94  pmtHitVectorIndex.resize(0);
95  }
96  std::vector<std::vector<double>> GetPMTVector() { return pmtHitVector; }
97  static std::vector<std::vector<double>> pmtHitVector;
98 
99  private:
100  // material property vector pointers, initialized in constructor,
101  // so we don't have to look them up every time DoIt is called.
102  G4MaterialPropertyVector *_rindex_glass; // function of photon energy
103  G4MaterialPropertyVector *_rindex_photocathode; // function of photon energy
104  G4MaterialPropertyVector *_kindex_photocathode; // n~ = n+ik, as in M.D.Lay
105  G4MaterialPropertyVector *_efficiency_photocathode; // necessary?
106  G4MaterialPropertyVector *_thickness_photocathode; // function of Z (mm)
107  // interior solid (vacuum), also initialized in constructor,
108  // so we don't have to look it up each time DoIt is called.
109  // Note it is implicitly assumed everywhere in the code that this solid
110  // is vacuum-filled and placed in the body with no offset or rotation.
111  G4VSolid *_inner1_solid;
112  G4VPhysicalVolume *_inner1_phys, *_inner2_phys, *_central_gap_phys;
113 
114  // "luxury level" -- how fancy should the optical model be?
115  G4int _luxlevel;
116  // MFB
117  G4double _rho;
118  G4double _rhoAvg;
119  G4double _rhoDif;
120  G4double _erfProb;
121  G4double _photocathode_MINrho;
122  G4double _photocathode_MAXrho;
123  G4double _efficiency_correction; // global efficiency correction, default to 1.0
124  G4double _dynodeTop;
125  G4double _dynodeRadius;
126  G4double _prepulseProb;
127 
128  G4bool _applyCorrection;
129 
130  // verbose level -- how verbose to be (diagnostics and such)
131  G4int _verbosity;
132 
133  // directory for commands
134  static G4UIdirectory *fgCmdDir;
135 
136  // "current values" of many parameters, for efficiency
137  // [I claim it is quicker to access these than to
138  // push them on the stack when calling CalculateCoefficients, Reflect, etc.]
139  // The following are set by DoIt() prior to any CalculateCoefficients() call.
140  G4double _photon_energy; // energy of current photon
141  G4double _wavelength; // wavelength of current photon
142  G4double _n1; // index of refraction of curr. medium at wavelength
143  G4double _n2, _k2; // index of refraction of photocathode at wavelength
144  G4double _n3; // index of refraction of far side at wavelength
145  G4double _efficiency; // efficiency of photocathode at wavelength (?)
146  G4double _thickness; // thickness of photocathode at current position
147  G4double _cos_theta1; // cosine of angle of incidence
148  // The following are set by CalculateCoefficients()
149  // and used by DoIt(), Refract(), and Reflect():
150  G4double _sin_theta1; // sine of angle of incidence
151  G4double _sin_theta3; // sine of angle of refraction
152  G4double _cos_theta3; // cosine of angle of refraction
153  G4double fR_s; // reflection coefficient for s-polarized light
154  G4double fT_s; // transmission coefficient for s-polarized light
155  G4double fR_p; // reflection coefficient for p-polarized light
156  G4double fT_p; // transmission coefficient for p-polarized light
157  G4double fR_n; // reference for fR_s/p at normal incidence
158  G4double fT_n; // reference for fT_s/p at normal incidence
159 
160  // private methods
161  void CalculateCoefficients(); // calculate and set fR_s, etc.
162  void Reflect(G4ThreeVector &dir, G4ThreeVector &pol, G4ThreeVector &norm);
163  void Refract(G4ThreeVector &dir, G4ThreeVector &pol, G4ThreeVector &norm);
164  int GetPMTID(const G4FastTrack &fastTrack);
165 
166  static double surfaceTolerance;
167 
168  std::map<int, double> EfficiencyCorrection;
169  std::vector<double> pmtHitVectorIndex;
170 };
171 
172 #endif
Definition: GLG4PMTOpticalModel.hh:38