CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
HFCherenkov Class Reference

#include <HFCherenkov.h>

Public Member Functions

void clearVectors ()
 
int computeNPE (G4ParticleDefinition *pDef, double pBeta, double u, double v, double w, double step_length, double zFiber, double Dose, int Npe_Dose)
 
int computeNPEinPMT (G4ParticleDefinition *pDef, double pBeta, double u, double v, double w, double step_length)
 
int computeNPhTrapped (double pBeta, double u, double v, double w, double step_length, double zFiber, double Dose, int Npe_Dose)
 
std::vector< double > getMom ()
 
std::vector< double > getWL ()
 
std::vector< double > getWLAtten ()
 
std::vector< double > getWLHEM ()
 
std::vector< double > getWLIni ()
 
std::vector< double > getWLQEff ()
 
std::vector< double > getWLTrap ()
 
 HFCherenkov (edm::ParameterSet const &p)
 
double smearNPE (G4int Npe)
 
virtual ~HFCherenkov ()
 

Private Member Functions

double computeHEMEff (double wavelength)
 
int computeNbOfPhotons (double pBeta, double step_length)
 
double computeQEff (double wavelength)
 
bool isApplicable (const G4ParticleDefinition *aParticleType)
 

Private Attributes

double aperture
 
double apertureTrap
 
bool checkSurvive
 
double gain
 
double lambda1
 
double lambda2
 
std::vector< double > momZ
 
G4ThreeVector phMom
 
double ref_index
 
std::vector< double > wl
 
std::vector< double > wlatten
 
std::vector< double > wlhem
 
std::vector< double > wlini
 
std::vector< double > wlqeff
 
std::vector< double > wltrap
 

Detailed Description

Definition at line 18 of file HFCherenkov.h.

Constructor & Destructor Documentation

HFCherenkov::HFCherenkov ( edm::ParameterSet const &  p)

Definition at line 13 of file HFCherenkov.cc.

References aperture, apertureTrap, checkSurvive, clearVectors(), funct::cos(), gain, edm::ParameterSet::getParameter(), lambda1, lambda2, funct::pow(), and ref_index.

13  {
14 
15  //Simple configurables
16  //static SimpleConfigurable<double> p1(1.459, "HFCherenkov:RefIndex");
17  //static SimpleConfigurable<double> p2(280.0, "HFCherenkov:Lambda1");
18  //static SimpleConfigurable<double> p3(700.0, "HFCherenkov:Lambda2");
19  //static SimpleConfigurable<double> p4(0.33, "HFCherenkov:Aperture");
20  //static SimpleConfigurable<double> p5(0.22, "HFCherenkov:ApertureTrapped");
21  //static SimpleConfigurable<double> p6(0.33, "HFCherenkov:Gain");
22  //static SimpleConfigurable<bool> p7(false, "HFCherenkov:CheckSurvive");
23 
24  ref_index = m_HF.getParameter<double>("RefIndex");
25  lambda1 = ((m_HF.getParameter<double>("Lambda1"))/pow(double(10),7))*cm;
26  lambda2 = ((m_HF.getParameter<double>("Lambda2"))/pow(double(10),7))*cm;
27  aperture = cos(asin(m_HF.getParameter<double>("Aperture")));
28  apertureTrap = cos(asin(m_HF.getParameter<double>("ApertureTrapped")));
29  gain = m_HF.getParameter<double>("Gain");
30  checkSurvive = m_HF.getParameter<bool>("CheckSurvive");
31 
32  edm::LogInfo("HFShower") << "HFCherenkov:: initialised with ref_index "
33  << ref_index << " lambda1/lambda2 (cm) "
34  << lambda1/cm << "/" << lambda2/cm
35  << " aperture(total/trapped) " << aperture << "/"
36  << apertureTrap << " Check photon survival in HF "
37  << checkSurvive << " Gain " << gain;
38 
39  clearVectors();
40 }
void clearVectors()
Definition: HFCherenkov.cc:379
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double aperture
Definition: HFCherenkov.h:61
double apertureTrap
Definition: HFCherenkov.h:61
double gain
Definition: HFCherenkov.h:62
double lambda1
Definition: HFCherenkov.h:60
double lambda2
Definition: HFCherenkov.h:60
bool checkSurvive
Definition: HFCherenkov.h:63
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double ref_index
Definition: HFCherenkov.h:59
HFCherenkov::~HFCherenkov ( )
virtual

Definition at line 42 of file HFCherenkov.cc.

42 {}

Member Function Documentation

void HFCherenkov::clearVectors ( )

Definition at line 379 of file HFCherenkov.cc.

References momZ, wl, wlatten, wlhem, wlini, wlqeff, and wltrap.

Referenced by computeNPE(), computeNPEinPMT(), and HFCherenkov().

379  {
380 
381  wl.clear();
382  wlini.clear();
383  wltrap.clear();
384  wlatten.clear();
385  wlhem.clear();
386  wlqeff.clear();
387  momZ.clear();
388 }
std::vector< double > wlqeff
Definition: HFCherenkov.h:72
std::vector< double > wlini
Definition: HFCherenkov.h:68
std::vector< double > wltrap
Definition: HFCherenkov.h:69
std::vector< double > wlhem
Definition: HFCherenkov.h:71
std::vector< double > wlatten
Definition: HFCherenkov.h:70
std::vector< double > momZ
Definition: HFCherenkov.h:67
std::vector< double > wl
Definition: HFCherenkov.h:66
double HFCherenkov::computeHEMEff ( double  wavelength)
private

Definition at line 334 of file HFCherenkov.cc.

References LogDebug.

Referenced by computeNPE().

334  {
335 
336  double hEMEff = 0;
337  if (wavelength < 400.) {
338  hEMEff = 0.;
339  } else if (wavelength >= 400. && wavelength < 410.) {
340  //hEMEff = .99 * (wavelength - 400.) / 10.;
341  hEMEff = (-1322.453 / wavelength) + 4.2056;
342  } else if (wavelength >= 410.) {
343  hEMEff = 0.99;
344  if (wavelength > 415. && wavelength < 445.) {
345  //abs(wavelength - 430.) < 15.
346  //hEMEff = 0.95;
347  hEMEff = 0.97;
348  } else if (wavelength > 550. && wavelength < 600.) {
349  // abs(wavelength - 575.) < 25.)
350  //hEMEff = 0.96;
351  hEMEff = 0.98;
352  } else if (wavelength > 565. && wavelength <= 635.) { // added later
353  // abs(wavelength - 600.) < 35.)
354  hEMEff = (701.7268 / wavelength) - 0.186;
355  }
356  }
357 #ifdef DebugLog
358  LogDebug("HFShower") << "HFCherenkov::computeHEMEff: wavelength "
359  << wavelength << " hEMEff " << hEMEff;
360 #endif
361  return hEMEff;
362 }
#define LogDebug(id)
int HFCherenkov::computeNbOfPhotons ( double  pBeta,
double  step_length 
)
private

Definition at line 301 of file HFCherenkov.cc.

References alpha, beta, lambda1, lambda2, LogDebug, M_PI, ref_index, and funct::sin().

Referenced by computeNPE(), computeNPEinPMT(), and computeNPhTrapped().

301  {
302 
303  double pBeta = beta;
304  double alpha = 0.0073;
305  double step_length = stepL;
306  double theta_C = acos(1./(pBeta*ref_index));
307  double lambdaDiff = (1./lambda1 - 1./lambda2);
308  double cherenPhPerLength = 2 * M_PI * alpha * lambdaDiff*cm;
309  double d_NOfPhotons = cherenPhPerLength * sin(theta_C)*sin(theta_C) * (step_length/cm);
310  int nbOfPhotons = int(d_NOfPhotons);
311 #ifdef DebugLog
312  LogDebug("HFShower") << "HFCherenkov::computeNbOfPhotons: StepLength "
313  << step_length << " theta_C " << theta_C
314  << " lambdaDiff " << lambdaDiff
315  << " cherenPhPerLength " << cherenPhPerLength
316  << " Photons " << d_NOfPhotons << " " << nbOfPhotons;
317 #endif
318  return nbOfPhotons;
319 }
#define LogDebug(id)
const double beta
float alpha
Definition: AMPTWrapper.h:95
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define M_PI
Definition: BFit3D.cc:3
double lambda1
Definition: HFCherenkov.h:60
double lambda2
Definition: HFCherenkov.h:60
double ref_index
Definition: HFCherenkov.h:59
int HFCherenkov::computeNPE ( G4ParticleDefinition *  pDef,
double  pBeta,
double  u,
double  v,
double  w,
double  step_length,
double  zFiber,
double  Dose,
int  Npe_Dose 
)

Definition at line 80 of file HFCherenkov.cc.

References aperture, checkSurvive, clearVectors(), computeHEMEff(), computeNbOfPhotons(), computeQEff(), funct::cos(), create_public_lumi_plots::exp, i, isApplicable(), lambda1, lambda2, LogDebug, M_PI, momZ, funct::pow(), rand(), ref_index, funct::sin(), mathSSE::sqrt(), wl, wlatten, wlhem, wlini, wlqeff, and wltrap.

Referenced by HFShower::getHits().

83  {
84 
85  clearVectors();
86  if (!isApplicable(pDef)) {return 0;}
87  if (pBeta < (1/ref_index) || step_length < 0.0001) {
88 #ifdef DebugLog
89  LogDebug("HFShower") << "HFCherenkov::computeNPE: pBeta " << pBeta
90  << " 1/mu " << (1/ref_index) << " step_length "
91  << step_length;
92 #endif
93  return 0;
94  }
95 
96  double uv = sqrt(u*u + v*v);
97  int nbOfPhotons = computeNbOfPhotons(pBeta, step_length);
98 #ifdef DebugLog
99  LogDebug("HFShower") << "HFCherenkov::computeNPE: pBeta " << pBeta
100  << " u/v/w " << u << "/" << v << "/" << w
101  << " step_length " << step_length << " zFib " << zFiber
102  << " nbOfPhotons " << nbOfPhotons;
103 #endif
104  if (nbOfPhotons < 0) {
105  return 0;
106  } else if (nbOfPhotons > 0) {
107  double w_ph=0;
108  for (int i = 0; i < nbOfPhotons; i++) {
109  double rand = G4UniformRand();
110  double theta_C = acos(1./(pBeta*ref_index));
111  double phi_C = 2*M_PI*rand;
112  double sinTheta = sin(theta_C);
113  double cosTheta = cos(theta_C);
114  double cosPhi = cos(phi_C);
115  //photon momentum
116  if (uv < 0.001) { // aligned with z-axis
117  w_ph = cosTheta;
118  } else { // general case
119  w_ph = w * cosTheta - sinTheta * cosPhi * uv;
120  }
121  double r_lambda = G4UniformRand();
122  double lambda0 = (lambda1 * lambda2) / (lambda2 - r_lambda *
123  (lambda2 - lambda1));
124  double lambda = (lambda0/cm) * pow(double(10),7); // lambda is in nm
125  wlini.push_back(lambda);
126 #ifdef DebugLog
127  LogDebug("HFShower") << "HFCherenkov::computeNPE: " << i << " lambda "
128  << lambda << " w_ph " << w_ph << " aperture "
129  << aperture;
130 #endif
131  if (w_ph > aperture) { // phton trapped inside fiber
132  wltrap.push_back(lambda);
133  double prob_HF = 1.0; //photon survived in HF
134  double a0_inv = 0.1234; //meter^-1
135  double prob_MX = exp( - 0.5 * a0_inv ); //light mixer
136  if (checkSurvive) {
137  double a_inv = a0_inv + 0.14 * pow(dose,0.30);
138  double z_meters = zFiber;
139  prob_HF = exp(-z_meters * a_inv ); //photon survived in HF
140  }
141  rand = G4UniformRand();
142 #ifdef DebugLog
143  LogDebug("HFShower") << "HFCherenkov::computeNPE: probHF " << prob_HF
144  << " prob_MX " << prob_MX << " Random " << rand
145  << " Survive? " << (rand < (prob_HF * prob_MX));
146 #endif
147  if (rand < (prob_HF * prob_MX)) { // survived and sent to light mixer
148  wlatten.push_back(lambda);
149  rand = G4UniformRand();
150  double effHEM = computeHEMEff(lambda);
151 #ifdef DebugLog
152  LogDebug("HFShower") << "HFCherenkov::computeNPE: w_ph " << w_ph
153  << " effHEM " << effHEM << " Random " << rand
154  << " Survive? " << (w_ph>0.997||(rand<effHEM));
155 #endif
156  if (w_ph>0.997 || (rand<effHEM)) { // survived HEM
157  wlhem.push_back(lambda);
158  double qEffic = computeQEff(lambda);
159  rand = G4UniformRand();
160 #ifdef DebugLog
161  LogDebug("HFShower") << "HFCherenkov::computeNPE: qEffic "
162  << qEffic << " Random " << rand
163  << " Survive? " <<(rand < qEffic);
164 #endif
165  if (rand < qEffic) { // made photoelectron
166  npe_Dose += 1;
167  momZ.push_back(w_ph);
168  wl.push_back(lambda);
169  wlqeff.push_back(lambda);
170  } // made pe
171  } // passed HEM
172  } // passed fiber
173  } // end of if(w_ph < w_aperture), trapped inside fiber
174  }// end of ++NbOfPhotons
175  } // end of if(NbOfPhotons)}
176  int npe = npe_Dose; // Nb of photoelectrons
177 #ifdef DebugLog
178  LogDebug("HFShower") << "HFCherenkov::computeNPE: npe " << npe;
179 #endif
180  return npe;
181 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
double computeQEff(double wavelength)
Definition: HFCherenkov.cc:321
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< double > wlqeff
Definition: HFCherenkov.h:72
void clearVectors()
Definition: HFCherenkov.cc:379
std::vector< double > wlini
Definition: HFCherenkov.h:68
double computeHEMEff(double wavelength)
Definition: HFCherenkov.cc:334
T sqrt(T t)
Definition: SSEVec.h:46
std::vector< double > wltrap
Definition: HFCherenkov.h:69
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double aperture
Definition: HFCherenkov.h:61
std::vector< double > wlhem
Definition: HFCherenkov.h:71
#define M_PI
Definition: BFit3D.cc:3
int computeNbOfPhotons(double pBeta, double step_length)
Definition: HFCherenkov.cc:301
Signal rand(Signal arg)
Definition: vlib.cc:442
bool isApplicable(const G4ParticleDefinition *aParticleType)
Definition: HFCherenkov.cc:390
std::vector< double > wlatten
Definition: HFCherenkov.h:70
double lambda1
Definition: HFCherenkov.h:60
double lambda2
Definition: HFCherenkov.h:60
bool checkSurvive
Definition: HFCherenkov.h:63
std::vector< double > momZ
Definition: HFCherenkov.h:67
mathSSE::Vec4< T > v
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double ref_index
Definition: HFCherenkov.h:59
T w() const
std::vector< double > wl
Definition: HFCherenkov.h:66
int HFCherenkov::computeNPEinPMT ( G4ParticleDefinition *  pDef,
double  pBeta,
double  u,
double  v,
double  w,
double  step_length 
)

Definition at line 183 of file HFCherenkov.cc.

References aperture, clearVectors(), computeNbOfPhotons(), computeQEff(), funct::cos(), i, isApplicable(), lambda1, lambda2, LogDebug, M_PI, momZ, funct::pow(), rand(), ref_index, funct::sin(), mathSSE::sqrt(), wl, wlatten, wlini, wlqeff, and wltrap.

Referenced by HFShowerFibreBundle::getHits(), and HFShowerPMT::getHits().

185  {
186  clearVectors();
187  int npe_ = 0;
188  if (!isApplicable(pDef)) {return 0;}
189  if (pBeta < (1/ref_index) || step_length < 0.0001) {
190 #ifdef DebugLog
191  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: pBeta " << pBeta
192  << " 1/mu " << (1/ref_index) << " step_length "
193  << step_length;
194 #endif
195  return 0;
196  }
197 
198  double uv = sqrt(u*u + v*v);
199  int nbOfPhotons = computeNbOfPhotons(pBeta, step_length);
200 #ifdef DebugLog
201  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: pBeta " << pBeta
202  << " u/v/w " << u << "/" << v << "/" << w
203  << " step_length " << step_length
204  << " nbOfPhotons " << nbOfPhotons;
205 #endif
206  if (nbOfPhotons < 0) {
207  return 0;
208  } else if (nbOfPhotons > 0) {
209  double w_ph=0;
210  for (int i = 0; i < nbOfPhotons; i++) {
211  double rand = G4UniformRand();
212  double theta_C = acos(1./(pBeta*ref_index));
213  double phi_C = 2*M_PI*rand;
214  double sinTheta = sin(theta_C);
215  double cosTheta = cos(theta_C);
216  double cosPhi = cos(phi_C);
217  //photon momentum
218  if (uv < 0.001) { // aligned with z-axis
219  w_ph = cosTheta;
220  } else { // general case
221  w_ph = w * cosTheta - sinTheta * cosPhi * uv;
222  }
223  double r_lambda = G4UniformRand();
224  double lambda0 = (lambda1 * lambda2) / (lambda2 - r_lambda *
225  (lambda2 - lambda1));
226  double lambda = (lambda0/cm) * pow(double(10),7); // lambda is in nm
227  wlini.push_back(lambda);
228 #ifdef DebugLog
229  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: " <<i <<" lambda "
230  << lambda << " w_ph " << w_ph << " aperture "
231  << aperture;
232 #endif
233  if (w_ph > aperture) { // phton trapped inside PMT glass
234  wltrap.push_back(lambda);
235  rand = G4UniformRand();
236 #ifdef DebugLog
237  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: Random " << rand
238  << " Survive? " << (rand < 1.);
239 #endif
240  if (rand < 1.0) { // survived all the times and sent to photo-cathode
241  wlatten.push_back(lambda);
242  rand = G4UniformRand();
243  double qEffic = computeQEff(lambda);//Quantum efficiency of the PMT
244  rand = G4UniformRand();
245 #ifdef DebugLog
246  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: qEffic "
247  << qEffic << " Random " << rand
248  << " Survive? " <<(rand < qEffic);
249 #endif
250  if (rand < qEffic) { // made photoelectron
251  npe_ += 1;
252  momZ.push_back(w_ph);
253  wl.push_back(lambda);
254  wlqeff.push_back(lambda);
255  } // made pe
256  } // accepted all Cherenkov photons
257  } // end of if(w_ph < w_aperture), trapped inside glass
258  }// end of ++NbOfPhotons
259  } // end of if(NbOfPhotons)}
260 #ifdef DebugLog
261  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: npe " << npe_;
262 #endif
263  return npe_;
264 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
double computeQEff(double wavelength)
Definition: HFCherenkov.cc:321
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< double > wlqeff
Definition: HFCherenkov.h:72
void clearVectors()
Definition: HFCherenkov.cc:379
std::vector< double > wlini
Definition: HFCherenkov.h:68
T sqrt(T t)
Definition: SSEVec.h:46
std::vector< double > wltrap
Definition: HFCherenkov.h:69
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double aperture
Definition: HFCherenkov.h:61
#define M_PI
Definition: BFit3D.cc:3
int computeNbOfPhotons(double pBeta, double step_length)
Definition: HFCherenkov.cc:301
Signal rand(Signal arg)
Definition: vlib.cc:442
bool isApplicable(const G4ParticleDefinition *aParticleType)
Definition: HFCherenkov.cc:390
std::vector< double > wlatten
Definition: HFCherenkov.h:70
double lambda1
Definition: HFCherenkov.h:60
double lambda2
Definition: HFCherenkov.h:60
std::vector< double > momZ
Definition: HFCherenkov.h:67
mathSSE::Vec4< T > v
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double ref_index
Definition: HFCherenkov.h:59
T w() const
std::vector< double > wl
Definition: HFCherenkov.h:66
int HFCherenkov::computeNPhTrapped ( double  pBeta,
double  u,
double  v,
double  w,
double  step_length,
double  zFiber,
double  Dose,
int  Npe_Dose 
)

Definition at line 44 of file HFCherenkov.cc.

References apertureTrap, computeNbOfPhotons(), funct::cos(), i, M_PI, rand(), ref_index, funct::sin(), and mathSSE::sqrt().

47  {
48 
49  if (pBeta < (1/ref_index) || step_length < 0.0001) {return 0;}
50 
51  double uv = sqrt(u*u + v*v);
52  int nbOfPhotons = computeNbOfPhotons(pBeta, step_length);
53 
54  if (nbOfPhotons < 0) {
55  return 0;
56  } else if (nbOfPhotons > 0) {
57  double w_ph=0;
58  for (int i = 0; i < nbOfPhotons; i++) {
59  double rand = G4UniformRand();
60  double theta_C = acos(1./(pBeta*ref_index));
61  double phi_C = 2*M_PI*rand;
62  double sinTheta = sin(theta_C);
63  double cosTheta = cos(theta_C);
64  double cosPhi = cos(phi_C);
65  //photon momentum
66  if (uv < 0.001) { // aligned with z-axis
67  w_ph = cosTheta;
68  } else { // general case
69  w_ph = w * cosTheta - sinTheta * cosPhi * uv;
70  }
71  if (w_ph > apertureTrap) { // phton trapped inside fiber
72  npe_Dose += 1;
73  }
74  }
75  }
76  int n_photons = npe_Dose;
77  return n_photons;
78 }
int i
Definition: DBlmapReader.cc:9
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T sqrt(T t)
Definition: SSEVec.h:46
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double apertureTrap
Definition: HFCherenkov.h:61
#define M_PI
Definition: BFit3D.cc:3
int computeNbOfPhotons(double pBeta, double step_length)
Definition: HFCherenkov.cc:301
Signal rand(Signal arg)
Definition: vlib.cc:442
mathSSE::Vec4< T > v
double ref_index
Definition: HFCherenkov.h:59
T w() const
double HFCherenkov::computeQEff ( double  wavelength)
private

Definition at line 321 of file HFCherenkov.cc.

References create_public_lumi_plots::exp, LogDebug, funct::pow(), and detailsBasic3DVector::y.

Referenced by computeNPE(), and computeNPEinPMT().

321  {
322 
323  double y = (wavelength - 275.) /180.;
324  double func = 1. / (1. + 250.*pow((y/5.),4));
325  double qE_R7525 = 0.77 * y * exp(-y) * func ;
326  double qeff = qE_R7525;
327 #ifdef DebugLog
328  LogDebug("HFShower") << "HFCherenkov::computeQEff: wavelength " << wavelength
329  << " y/func " << y << "/" << func << " qeff " << qeff;
330 #endif
331  return qeff;
332 }
#define LogDebug(id)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
std::vector< double > HFCherenkov::getMom ( )

Definition at line 296 of file HFCherenkov.cc.

References momZ, and v.

Referenced by HFShower::getHits().

296  {
297  std::vector<double> v = momZ;
298  return v;
299 }
std::vector< double > momZ
Definition: HFCherenkov.h:67
mathSSE::Vec4< T > v
std::vector< double > HFCherenkov::getWL ( )

Definition at line 291 of file HFCherenkov.cc.

References v, and wl.

Referenced by HFShower::getHits().

291  {
292  std::vector<double> v = wl;
293  return v;
294 }
mathSSE::Vec4< T > v
std::vector< double > wl
Definition: HFCherenkov.h:66
std::vector< double > HFCherenkov::getWLAtten ( )

Definition at line 276 of file HFCherenkov.cc.

References v, and wlatten.

276  {
277  std::vector<double> v = wlatten;
278  return v;
279 }
std::vector< double > wlatten
Definition: HFCherenkov.h:70
mathSSE::Vec4< T > v
std::vector< double > HFCherenkov::getWLHEM ( )

Definition at line 281 of file HFCherenkov.cc.

References v, and wlhem.

281  {
282  std::vector<double> v = wlhem;
283  return v;
284 }
std::vector< double > wlhem
Definition: HFCherenkov.h:71
mathSSE::Vec4< T > v
std::vector< double > HFCherenkov::getWLIni ( )

Definition at line 266 of file HFCherenkov.cc.

References v, and wlini.

266  {
267  std::vector<double> v = wlini;
268  return v;
269 }
std::vector< double > wlini
Definition: HFCherenkov.h:68
mathSSE::Vec4< T > v
std::vector< double > HFCherenkov::getWLQEff ( )

Definition at line 286 of file HFCherenkov.cc.

References v, and wlqeff.

286  {
287  std::vector<double> v = wlqeff;
288  return v;
289 }
std::vector< double > wlqeff
Definition: HFCherenkov.h:72
mathSSE::Vec4< T > v
std::vector< double > HFCherenkov::getWLTrap ( )

Definition at line 271 of file HFCherenkov.cc.

References v, and wltrap.

271  {
272  std::vector<double> v = wltrap;
273  return v;
274 }
std::vector< double > wltrap
Definition: HFCherenkov.h:69
mathSSE::Vec4< T > v
bool HFCherenkov::isApplicable ( const G4ParticleDefinition *  aParticleType)
private

Definition at line 390 of file HFCherenkov.cc.

References LogDebug, and tmp.

Referenced by computeNPE(), and computeNPEinPMT().

390  {
391  bool tmp = (aParticleType->GetPDGCharge() != 0);
392 #ifdef DebugLog
393  LogDebug("HFShower") << "HFCherenkov::isApplicable: aParticleType "
394  << aParticleType->GetParticleName() << " PDGCharge "
395  << aParticleType->GetPDGCharge() << " Result " << tmp;
396 #endif
397  return tmp;
398 }
#define LogDebug(id)
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
double HFCherenkov::smearNPE ( G4int  Npe)

Definition at line 364 of file HFCherenkov.cc.

References gain, i, and LogDebug.

364  {
365 
366  double pe = 0.;
367  if (npe > 0) {
368  for (int i = 0; i < npe; ++i) {
369  double val = G4Poisson(gain);
370  pe += (val/gain) + 0.001*G4UniformRand();
371  }
372  }
373 #ifdef DebugLog
374  LogDebug("HFShower") << "HFCherenkov::smearNPE: npe " << npe << " pe " << pe;
375 #endif
376  return pe;
377 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
double gain
Definition: HFCherenkov.h:62

Member Data Documentation

double HFCherenkov::aperture
private

Definition at line 61 of file HFCherenkov.h.

Referenced by computeNPE(), computeNPEinPMT(), and HFCherenkov().

double HFCherenkov::apertureTrap
private

Definition at line 61 of file HFCherenkov.h.

Referenced by computeNPhTrapped(), and HFCherenkov().

bool HFCherenkov::checkSurvive
private

Definition at line 63 of file HFCherenkov.h.

Referenced by computeNPE(), and HFCherenkov().

double HFCherenkov::gain
private

Definition at line 62 of file HFCherenkov.h.

Referenced by HFCherenkov(), and smearNPE().

double HFCherenkov::lambda1
private

Definition at line 60 of file HFCherenkov.h.

Referenced by computeNbOfPhotons(), computeNPE(), computeNPEinPMT(), and HFCherenkov().

double HFCherenkov::lambda2
private

Definition at line 60 of file HFCherenkov.h.

Referenced by computeNbOfPhotons(), computeNPE(), computeNPEinPMT(), and HFCherenkov().

std::vector<double> HFCherenkov::momZ
private

Definition at line 67 of file HFCherenkov.h.

Referenced by clearVectors(), computeNPE(), computeNPEinPMT(), and getMom().

G4ThreeVector HFCherenkov::phMom
private

Definition at line 65 of file HFCherenkov.h.

double HFCherenkov::ref_index
private
std::vector<double> HFCherenkov::wl
private

Definition at line 66 of file HFCherenkov.h.

Referenced by clearVectors(), computeNPE(), computeNPEinPMT(), and getWL().

std::vector<double> HFCherenkov::wlatten
private

Definition at line 70 of file HFCherenkov.h.

Referenced by clearVectors(), computeNPE(), computeNPEinPMT(), and getWLAtten().

std::vector<double> HFCherenkov::wlhem
private

Definition at line 71 of file HFCherenkov.h.

Referenced by clearVectors(), computeNPE(), and getWLHEM().

std::vector<double> HFCherenkov::wlini
private

Definition at line 68 of file HFCherenkov.h.

Referenced by clearVectors(), computeNPE(), computeNPEinPMT(), and getWLIni().

std::vector<double> HFCherenkov::wlqeff
private

Definition at line 72 of file HFCherenkov.h.

Referenced by clearVectors(), computeNPE(), computeNPEinPMT(), and getWLQEff().

std::vector<double> HFCherenkov::wltrap
private

Definition at line 69 of file HFCherenkov.h.

Referenced by clearVectors(), computeNPE(), computeNPEinPMT(), and getWLTrap().