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:394
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 394 of file HFCherenkov.cc.

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

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

394  {
395 
396  wl.clear();
397  wlini.clear();
398  wltrap.clear();
399  wlatten.clear();
400  wlhem.clear();
401  wlqeff.clear();
402  momZ.clear();
403 }
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 349 of file HFCherenkov.cc.

References LogDebug.

Referenced by computeNPE().

349  {
350 
351  double hEMEff = 0;
352  if (wavelength < 400.) {
353  hEMEff = 0.;
354  } else if (wavelength >= 400. && wavelength < 410.) {
355  //hEMEff = .99 * (wavelength - 400.) / 10.;
356  hEMEff = (-1322.453 / wavelength) + 4.2056;
357  } else if (wavelength >= 410.) {
358  hEMEff = 0.99;
359  if (wavelength > 415. && wavelength < 445.) {
360  //abs(wavelength - 430.) < 15.
361  //hEMEff = 0.95;
362  hEMEff = 0.97;
363  } else if (wavelength > 550. && wavelength < 600.) {
364  // abs(wavelength - 575.) < 25.)
365  //hEMEff = 0.96;
366  hEMEff = 0.98;
367  } else if (wavelength > 565. && wavelength <= 635.) { // added later
368  // abs(wavelength - 600.) < 35.)
369  hEMEff = (701.7268 / wavelength) - 0.186;
370  }
371  }
372 #ifdef DebugLog
373  LogDebug("HFShower") << "HFCherenkov::computeHEMEff: wavelength "
374  << wavelength << " hEMEff " << hEMEff;
375 #endif
376  return hEMEff;
377 }
#define LogDebug(id)
int HFCherenkov::computeNbOfPhotons ( double  pBeta,
double  step_length 
)
private

Definition at line 316 of file HFCherenkov.cc.

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

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

316  {
317 
318  double pBeta = beta;
319  double alpha = 0.0073;
320  double step_length = stepL;
321  double theta_C = acos(1./(pBeta*ref_index));
322  double lambdaDiff = (1./lambda1 - 1./lambda2);
323  double cherenPhPerLength = 2 * M_PI * alpha * lambdaDiff*cm;
324  double d_NOfPhotons = cherenPhPerLength * sin(theta_C)*sin(theta_C) * (step_length/cm);
325  int nbOfPhotons = int(d_NOfPhotons);
326 #ifdef DebugLog
327  LogDebug("HFShower") << "HFCherenkov::computeNbOfPhotons: StepLength "
328  << step_length << " theta_C " << theta_C
329  << " lambdaDiff " << lambdaDiff
330  << " cherenPhPerLength " << cherenPhPerLength
331  << " Photons " << d_NOfPhotons << " " << nbOfPhotons;
332 #endif
333  return nbOfPhotons;
334 }
#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 85 of file HFCherenkov.cc.

References aperture, checkSurvive, clearVectors(), computeHEMEff(), computeNbOfPhotons(), computeQEff(), funct::cos(), funct::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().

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

195  {
196  clearVectors();
197  int npe_ = 0;
198  if (!isApplicable(pDef)) {return 0;}
199  if (pBeta < (1/ref_index) || step_length < 0.0001) {
200 #ifdef DebugLog
201  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: pBeta " << pBeta
202  << " 1/mu " << (1/ref_index) << " step_length "
203  << step_length;
204 #endif
205  return 0;
206  }
207 
208  double uv = sqrt(u*u + v*v);
209  int nbOfPhotons = computeNbOfPhotons(pBeta, step_length);
210 #ifdef DebugLog
211  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: pBeta " << pBeta
212  << " u/v/w " << u << "/" << v << "/" << w
213  << " step_length " << step_length
214  << " nbOfPhotons " << nbOfPhotons;
215 #endif
216  if (nbOfPhotons < 0) {
217  return 0;
218  } else if (nbOfPhotons > 0) {
219  double u_ph, v_ph, w_ph=0;
220  for (int i = 0; i < nbOfPhotons; i++) {
221  double rand = G4UniformRand();
222  double theta_C = acos(1./(pBeta*ref_index));
223  double phi_C = 2*M_PI*rand;
224  double sinTheta = sin(theta_C);
225  double cosTheta = cos(theta_C);
226  double sinPhi = sin(phi_C);
227  double cosPhi = cos(phi_C);
228  //photon momentum
229  if (uv < 0.001) { // aligned with z-axis
230  u_ph = sinTheta * cosPhi;
231  v_ph = sinTheta * sinPhi;
232  w_ph = cosTheta;
233  } else { // general case
234  u_ph = u * cosTheta + sinTheta * (v*sinPhi + u*w*cosPhi)/ uv;
235  v_ph = v * cosTheta + sinTheta * (-u*sinPhi + v*w*cosPhi)/ uv;
236  w_ph = w * cosTheta - sinTheta * cosPhi * uv;
237  }
238  double r_lambda = G4UniformRand();
239  double lambda0 = (lambda1 * lambda2) / (lambda2 - r_lambda *
240  (lambda2 - lambda1));
241  double lambda = (lambda0/cm) * pow(double(10),7); // lambda is in nm
242  wlini.push_back(lambda);
243 #ifdef DebugLog
244  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: " <<i <<" lambda "
245  << lambda << " w_ph " << w_ph << " aperture "
246  << aperture;
247 #endif
248  if (w_ph > aperture) { // phton trapped inside PMT glass
249  wltrap.push_back(lambda);
250  rand = G4UniformRand();
251 #ifdef DebugLog
252  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: Random " << rand
253  << " Survive? " << (rand < 1.);
254 #endif
255  if (rand < 1.0) { // survived all the times and sent to photo-cathode
256  wlatten.push_back(lambda);
257  rand = G4UniformRand();
258  double qEffic = computeQEff(lambda);//Quantum efficiency of the PMT
259  rand = G4UniformRand();
260 #ifdef DebugLog
261  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: qEffic "
262  << qEffic << " Random " << rand
263  << " Survive? " <<(rand < qEffic);
264 #endif
265  if (rand < qEffic) { // made photoelectron
266  npe_ += 1;
267  momZ.push_back(w_ph);
268  wl.push_back(lambda);
269  wlqeff.push_back(lambda);
270  } // made pe
271  } // accepted all Cherenkov photons
272  } // end of if(w_ph < w_aperture), trapped inside glass
273  }// end of ++NbOfPhotons
274  } // end of if(NbOfPhotons)}
275 #ifdef DebugLog
276  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: npe " << npe_;
277 #endif
278  return npe_;
279 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
double computeQEff(double wavelength)
Definition: HFCherenkov.cc:336
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< double > wlqeff
Definition: HFCherenkov.h:72
void clearVectors()
Definition: HFCherenkov.cc:394
std::vector< double > wlini
Definition: HFCherenkov.h:68
T sqrt(T t)
Definition: SSEVec.h:28
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:316
Signal rand(Signal arg)
Definition: vlib.cc:442
bool isApplicable(const G4ParticleDefinition *aParticleType)
Definition: HFCherenkov.cc:405
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
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 u_ph, v_ph, 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 sinPhi = sin(phi_C);
65  double cosPhi = cos(phi_C);
66  //photon momentum
67  if (uv < 0.001) { // aligned with z-axis
68  u_ph = sinTheta * cosPhi;
69  v_ph = sinTheta * sinPhi;
70  w_ph = cosTheta;
71  } else { // general case
72  u_ph = u * cosTheta + sinTheta * (v*sinPhi + u*w*cosPhi)/ uv;
73  v_ph = v * cosTheta + sinTheta * (-u*sinPhi + v*w*cosPhi)/ uv;
74  w_ph = w * cosTheta - sinTheta * cosPhi * uv;
75  }
76  if (w_ph > apertureTrap) { // phton trapped inside fiber
77  npe_Dose += 1;
78  }
79  }
80  }
81  int n_photons = npe_Dose;
82  return n_photons;
83 }
int i
Definition: DBlmapReader.cc:9
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T sqrt(T t)
Definition: SSEVec.h:28
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:316
Signal rand(Signal arg)
Definition: vlib.cc:442
mathSSE::Vec4< T > v
double ref_index
Definition: HFCherenkov.h:59
double HFCherenkov::computeQEff ( double  wavelength)
private

Definition at line 336 of file HFCherenkov.cc.

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

Referenced by computeNPE(), and computeNPEinPMT().

336  {
337 
338  double y = (wavelength - 275.) /180.;
339  double func = 1. / (1. + 250.*pow((y/5.),4));
340  double qE_R7525 = 0.77 * y * exp(-y) * func ;
341  double qeff = qE_R7525;
342 #ifdef DebugLog
343  LogDebug("HFShower") << "HFCherenkov::computeQEff: wavelength " << wavelength
344  << " y/func " << y << "/" << func << " qeff " << qeff;
345 #endif
346  return qeff;
347 }
#define LogDebug(id)
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
std::vector< double > HFCherenkov::getMom ( )

Definition at line 311 of file HFCherenkov.cc.

References momZ, and v.

Referenced by HFShower::getHits().

311  {
312  std::vector<double> v = momZ;
313  return v;
314 }
std::vector< double > momZ
Definition: HFCherenkov.h:67
mathSSE::Vec4< T > v
std::vector< double > HFCherenkov::getWL ( )

Definition at line 306 of file HFCherenkov.cc.

References v, and wl.

Referenced by HFShower::getHits().

306  {
307  std::vector<double> v = wl;
308  return v;
309 }
mathSSE::Vec4< T > v
std::vector< double > wl
Definition: HFCherenkov.h:66
std::vector< double > HFCherenkov::getWLAtten ( )

Definition at line 291 of file HFCherenkov.cc.

References v, and wlatten.

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

Definition at line 296 of file HFCherenkov.cc.

References v, and wlhem.

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

Definition at line 281 of file HFCherenkov.cc.

References v, and wlini.

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

Definition at line 301 of file HFCherenkov.cc.

References v, and wlqeff.

301  {
302  std::vector<double> v = wlqeff;
303  return v;
304 }
std::vector< double > wlqeff
Definition: HFCherenkov.h:72
mathSSE::Vec4< T > v
std::vector< double > HFCherenkov::getWLTrap ( )

Definition at line 286 of file HFCherenkov.cc.

References v, and wltrap.

286  {
287  std::vector<double> v = wltrap;
288  return v;
289 }
std::vector< double > wltrap
Definition: HFCherenkov.h:69
mathSSE::Vec4< T > v
bool HFCherenkov::isApplicable ( const G4ParticleDefinition *  aParticleType)
private

Definition at line 405 of file HFCherenkov.cc.

References LogDebug, and tmp.

Referenced by computeNPE(), and computeNPEinPMT().

405  {
406  bool tmp = (aParticleType->GetPDGCharge() != 0);
407 #ifdef DebugLog
408  LogDebug("HFShower") << "HFCherenkov::isApplicable: aParticleType "
409  << aParticleType->GetParticleName() << " PDGCharge "
410  << aParticleType->GetPDGCharge() << " Result " << tmp;
411 #endif
412  return tmp;
413 }
#define LogDebug(id)
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
double HFCherenkov::smearNPE ( G4int  Npe)

Definition at line 379 of file HFCherenkov.cc.

References gain, i, and LogDebug.

379  {
380 
381  double pe = 0.;
382  if (npe > 0) {
383  for (int i = 0; i < npe; ++i) {
384  double val = G4Poisson(gain);
385  pe += (val/gain) + 0.001*G4UniformRand();
386  }
387  }
388 #ifdef DebugLog
389  LogDebug("HFShower") << "HFCherenkov::smearNPE: npe " << npe << " pe " << pe;
390 #endif
391  return pe;
392 }
#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().