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 (G4Step *step, 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
 
double aperturetrapped
 
bool checkSurvive
 
double fibreR
 
double gain
 
double lambda1
 
double lambda2
 
std::vector< double > momZ
 
G4ThreeVector phMom
 
double ref_index
 
double sinPsimax
 
bool UseNewPMT
 
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 19 of file HFCherenkov.h.

Constructor & Destructor Documentation

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

Definition at line 18 of file HFCherenkov.cc.

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

18  {
19 
20  ref_index = m_HF.getParameter<double>("RefIndex");
21  lambda1 = ((m_HF.getParameter<double>("Lambda1"))/pow(double(10),7))*cm;
22  lambda2 = ((m_HF.getParameter<double>("Lambda2"))/pow(double(10),7))*cm;
23  aperture = cos(asin(m_HF.getParameter<double>("Aperture")));
24  apertureTrap = cos(asin(m_HF.getParameter<double>("ApertureTrapped")));
25  aperturetrapped = m_HF.getParameter<double>("CosApertureTrapped");
26  gain = m_HF.getParameter<double>("Gain");
27  checkSurvive = m_HF.getParameter<bool>("CheckSurvive");
28  UseNewPMT = m_HF.getParameter<bool>("UseR7600UPMT");
29  sinPsimax = m_HF.getUntrackedParameter<double>("SinPsiMax",0.5);
30  fibreR = m_HF.getUntrackedParameter<double>("FibreR",0.3)*mm;
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 << "|" << aperturetrapped
37  << " Check photon survival in HF " << checkSurvive
38  << " Gain " << gain << " useNewPMT " << UseNewPMT
39  << " FibreR " << fibreR;
40 
41  clearVectors();
42 }
double fibreR
Definition: HFCherenkov.h:63
double aperturetrapped
Definition: HFCherenkov.h:62
void clearVectors()
Definition: HFCherenkov.cc:446
bool UseNewPMT
Definition: HFCherenkov.h:65
double sinPsimax
Definition: HFCherenkov.h:63
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double aperture
Definition: HFCherenkov.h:62
double apertureTrap
Definition: HFCherenkov.h:62
double gain
Definition: HFCherenkov.h:63
double lambda1
Definition: HFCherenkov.h:61
double lambda2
Definition: HFCherenkov.h:61
bool checkSurvive
Definition: HFCherenkov.h:64
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double ref_index
Definition: HFCherenkov.h:60
HFCherenkov::~HFCherenkov ( )
virtual

Definition at line 44 of file HFCherenkov.cc.

44 {}

Member Function Documentation

void HFCherenkov::clearVectors ( )

Definition at line 446 of file HFCherenkov.cc.

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

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

446  {
447 
448  wl.clear();
449  wlini.clear();
450  wltrap.clear();
451  wlatten.clear();
452  wlhem.clear();
453  wlqeff.clear();
454  momZ.clear();
455 }
std::vector< double > wlqeff
Definition: HFCherenkov.h:74
std::vector< double > wlini
Definition: HFCherenkov.h:70
std::vector< double > wltrap
Definition: HFCherenkov.h:71
std::vector< double > wlhem
Definition: HFCherenkov.h:73
std::vector< double > wlatten
Definition: HFCherenkov.h:72
std::vector< double > momZ
Definition: HFCherenkov.h:69
std::vector< double > wl
Definition: HFCherenkov.h:68
double HFCherenkov::computeHEMEff ( double  wavelength)
private

Definition at line 401 of file HFCherenkov.cc.

References LogDebug.

Referenced by computeNPE().

401  {
402 
403  double hEMEff = 0;
404  if (wavelength < 400.) {
405  hEMEff = 0.;
406  } else if (wavelength >= 400. && wavelength < 410.) {
407  //hEMEff = .99 * (wavelength - 400.) / 10.;
408  hEMEff = (-1322.453 / wavelength) + 4.2056;
409  } else if (wavelength >= 410.) {
410  hEMEff = 0.99;
411  if (wavelength > 415. && wavelength < 445.) {
412  //abs(wavelength - 430.) < 15.
413  //hEMEff = 0.95;
414  hEMEff = 0.97;
415  } else if (wavelength > 550. && wavelength < 600.) {
416  // abs(wavelength - 575.) < 25.)
417  //hEMEff = 0.96;
418  hEMEff = 0.98;
419  } else if (wavelength > 565. && wavelength <= 635.) { // added later
420  // abs(wavelength - 600.) < 35.)
421  hEMEff = (701.7268 / wavelength) - 0.186;
422  }
423  }
424 #ifdef DebugLog
425  LogDebug("HFShower") << "HFCherenkov::computeHEMEff: wavelength "
426  << wavelength << " hEMEff " << hEMEff;
427 #endif
428  return hEMEff;
429 }
#define LogDebug(id)
int HFCherenkov::computeNbOfPhotons ( double  pBeta,
double  step_length 
)
private

Definition at line 346 of file HFCherenkov.cc.

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

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

346  {
347 
348  double pBeta = beta;
349  double alpha = 0.0073;
350  double step_length = stepL;
351  double theta_C = acos(1./(pBeta*ref_index));
352  double lambdaDiff = (1./lambda1 - 1./lambda2);
353  double cherenPhPerLength = 2 * M_PI * alpha * lambdaDiff*cm;
354  double d_NOfPhotons = cherenPhPerLength * sin(theta_C)*sin(theta_C) * (step_length/cm);
355  int nbOfPhotons = int(d_NOfPhotons);
356 #ifdef DebugLog
357  LogDebug("HFShower") << "HFCherenkov::computeNbOfPhotons: StepLength "
358  << step_length << " theta_C " << theta_C
359  << " lambdaDiff " << lambdaDiff
360  << " cherenPhPerLength " << cherenPhPerLength
361  << " Photons " << d_NOfPhotons << " " << nbOfPhotons;
362 #endif
363  return nbOfPhotons;
364 }
#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
double lambda1
Definition: HFCherenkov.h:61
double lambda2
Definition: HFCherenkov.h:61
double ref_index
Definition: HFCherenkov.h:60
int HFCherenkov::computeNPE ( G4Step *  step,
G4ParticleDefinition *  pDef,
double  pBeta,
double  u,
double  v,
double  w,
double  step_length,
double  zFiber,
double  Dose,
int  Npe_Dose 
)

Definition at line 82 of file HFCherenkov.cc.

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

Referenced by HFShower::getHits().

85  {
86 
87  clearVectors();
88  if (!isApplicable(pDef)) {return 0;}
89  if (pBeta < (1/ref_index) || step_length < 0.0001) {
90 #ifdef DebugLog
91  LogDebug("HFShower") << "HFCherenkov::computeNPE: pBeta " << pBeta
92  << " 1/mu " << (1/ref_index) << " step_length "
93  << step_length;
94 #endif
95  return 0;
96  }
97 
98  double uv = sqrt(u*u + v*v);
99  int nbOfPhotons = computeNbOfPhotons(pBeta, step_length)
100  *aStep->GetTrack()->GetWeight();
101 #ifdef DebugLog
102  LogDebug("HFShower") << "HFCherenkov::computeNPE: pBeta " << pBeta
103  << " u/v/w " << u << "/" << v << "/" << w
104  << " step_length " << step_length << " zFib " << zFiber
105  << " nbOfPhotons " << nbOfPhotons;
106 #endif
107  if (nbOfPhotons < 0) {
108  return 0;
109  } else if (nbOfPhotons > 0) {
110  G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
111  G4TouchableHandle theTouchable = preStepPoint->GetTouchableHandle();
112  G4ThreeVector localprepos = theTouchable->GetHistory()->
113  GetTopTransform().TransformPoint(aStep->GetPreStepPoint()->GetPosition());
114  G4ThreeVector localpostpos = theTouchable->GetHistory()->
115  GetTopTransform().TransformPoint(aStep->GetPostStepPoint()->GetPosition());
116 
117  double length=sqrt((localpostpos.x()-localprepos.x())*(localpostpos.x()-localprepos.x())
118  +(localpostpos.y()-localprepos.y())*(localpostpos.y()-localprepos.y()));
119  double yemit = std::sqrt(fibreR*fibreR-length*length/4.);
120 
121  double u_ph=0,v_ph=0, w_ph=0;
122  for (int i = 0; i < nbOfPhotons; i++) {
123  double rand = G4UniformRand();
124  double theta_C = acos(1./(pBeta*ref_index));
125  double phi_C = 2*M_PI*rand;
126  double sinTheta = sin(theta_C);
127  double cosTheta = cos(theta_C);
128  double cosPhi = cos(phi_C);
129  double sinPhi = sin(phi_C);
130  //photon momentum
131  if (uv < 0.001) { // aligned with z-axis
132  u_ph = sinTheta * cosPhi ;
133  v_ph = sinTheta * sinPhi;
134  w_ph = cosTheta;
135  } else { // general case
136  u_ph = uv * cosTheta + sinTheta * cosPhi * w;
137  v_ph = sinTheta * sinPhi;
138  w_ph = w * cosTheta - sinTheta * cosPhi * uv;
139  }
140  double r_lambda = G4UniformRand();
141  double lambda0 = (lambda1 * lambda2) / (lambda2 - r_lambda *
142  (lambda2 - lambda1));
143  double lambda = (lambda0/cm) * pow(double(10),7); // lambda is in nm
144  wlini.push_back(lambda);
145 #ifdef DebugLog
146  LogDebug("HFShower") << "HFCherenkov::computeNPE: " << i << " lambda "
147  << lambda << " w_ph " << w_ph << " aperture "
148  << aperture;
149 #endif
150 // --------------
151  double xemit=length*(G4UniformRand()-0.5);
152  double gam=atan2(yemit,xemit);
153  double eps=atan2(v_ph,u_ph);
154  double sinBeta=sin(gam-eps);
155  double rho=sqrt(xemit*xemit+yemit*yemit);
156  double sinEta=rho/fibreR*sinBeta;
157  double cosEta=sqrt(1.-sinEta*sinEta);
158  double sinPsi=sqrt(1.-w_ph*w_ph);
159  double cosKsi=cosEta*sinPsi;
160 #ifdef DebugLog
161  if (cosKsi < aperturetrapped && w_ph>0.) {
162  LogDebug("HFShower") << "HFCherenkov::Trapped photon : " << u_ph << " "
163  << v_ph << " " << w_ph << " " << xemit << " "
164  << gam << " " << eps << " " << sinBeta << " "
165  << rho << " " << sinEta << " " << cosEta << " "
166  << " " << sinPsi << " " << cosKsi;
167  } else {
168  LogDebug("HFShower") << "HFCherenkov::Rejected photon : " << u_ph <<" "
169  << v_ph << " " << w_ph << " " << xemit << " "
170  << gam << " " << eps << " " << sinBeta << " "
171  << rho << " " << sinEta << " " << cosEta << " "
172  << " " << sinPsi << " " << cosKsi;
173  }
174 #endif
175  if (cosKsi < aperturetrapped // photon trapped inside fiber
176  && w_ph>0. // and moves to PMT
177  && sinPsi < sinPsimax) { // and is not reflected at fiber end
178  wltrap.push_back(lambda);
179  double prob_HF = 1.0; //photon survived in HF
180  double a0_inv = 0.1234; //meter^-1
181  double prob_MX = exp( - 0.5 * a0_inv ); //light mixer
182  if (checkSurvive) {
183  double a_inv = a0_inv + 0.14 * pow(dose,0.30);
184  double z_meters = zFiber;
185  prob_HF = exp(-z_meters * a_inv ); //photon survived in HF
186  }
187  rand = G4UniformRand();
188 #ifdef DebugLog
189  LogDebug("HFShower") << "HFCherenkov::computeNPE: probHF " << prob_HF
190  << " prob_MX " << prob_MX << " Random " << rand
191  << " Survive? " << (rand < (prob_HF * prob_MX));
192 #endif
193  if (rand < (prob_HF * prob_MX)) { // survived and sent to light mixer
194  wlatten.push_back(lambda);
195  rand = G4UniformRand();
196  double effHEM = computeHEMEff(lambda);
197 #ifdef DebugLog
198  LogDebug("HFShower") << "HFCherenkov::computeNPE: w_ph " << w_ph
199  << " effHEM " << effHEM << " Random " << rand
200  << " Survive? " << (w_ph>0.997||(rand<effHEM));
201 #endif
202  if (w_ph>0.997 || (rand<effHEM)) { // survived HEM
203  wlhem.push_back(lambda);
204  double qEffic = computeQEff(lambda);
205  rand = G4UniformRand();
206 #ifdef DebugLog
207  LogDebug("HFShower") << "HFCherenkov::computeNPE: qEffic "
208  << qEffic << " Random " << rand
209  << " Survive? " <<(rand < qEffic);
210 #endif
211  if (rand < qEffic) { // made photoelectron
212  npe_Dose += 1;
213  momZ.push_back(w_ph);
214  wl.push_back(lambda);
215  wlqeff.push_back(lambda);
216  } // made pe
217  } // passed HEM
218  } // passed fiber
219  } // end of if(w_ph < w_aperture), trapped inside fiber
220  }// end of ++NbOfPhotons
221  } // end of if(NbOfPhotons)}
222  int npe = npe_Dose; // Nb of photoelectrons
223 #ifdef DebugLog
224  LogDebug("HFShower") << "HFCherenkov::computeNPE: npe " << npe;
225 #endif
226  return npe;
227 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
double computeQEff(double wavelength)
Definition: HFCherenkov.cc:366
double fibreR
Definition: HFCherenkov.h:63
const double w
Definition: UKUtility.cc:23
double aperturetrapped
Definition: HFCherenkov.h:62
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< double > wlqeff
Definition: HFCherenkov.h:74
void clearVectors()
Definition: HFCherenkov.cc:446
std::vector< double > wlini
Definition: HFCherenkov.h:70
double computeHEMEff(double wavelength)
Definition: HFCherenkov.cc:401
double sinPsimax
Definition: HFCherenkov.h:63
T sqrt(T t)
Definition: SSEVec.h:18
std::vector< double > wltrap
Definition: HFCherenkov.h:71
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double aperture
Definition: HFCherenkov.h:62
std::vector< double > wlhem
Definition: HFCherenkov.h:73
#define M_PI
int computeNbOfPhotons(double pBeta, double step_length)
Definition: HFCherenkov.cc:346
Signal rand(Signal arg)
Definition: vlib.cc:442
bool isApplicable(const G4ParticleDefinition *aParticleType)
Definition: HFCherenkov.cc:457
std::vector< double > wlatten
Definition: HFCherenkov.h:72
double lambda1
Definition: HFCherenkov.h:61
double lambda2
Definition: HFCherenkov.h:61
bool checkSurvive
Definition: HFCherenkov.h:64
std::vector< double > momZ
Definition: HFCherenkov.h:69
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double ref_index
Definition: HFCherenkov.h:60
std::vector< double > wl
Definition: HFCherenkov.h:68
int HFCherenkov::computeNPEinPMT ( G4ParticleDefinition *  pDef,
double  pBeta,
double  u,
double  v,
double  w,
double  step_length 
)

Definition at line 229 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().

231  {
232  clearVectors();
233  int npe_ = 0;
234  if (!isApplicable(pDef)) {return 0;}
235  if (pBeta < (1/ref_index) || step_length < 0.0001) {
236 #ifdef DebugLog
237  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: pBeta " << pBeta
238  << " 1/mu " << (1/ref_index) << " step_length "
239  << step_length;
240 #endif
241  return 0;
242  }
243 
244  double uv = sqrt(u*u + v*v);
245  int nbOfPhotons = computeNbOfPhotons(pBeta, step_length);
246 #ifdef DebugLog
247  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: pBeta " << pBeta
248  << " u/v/w " << u << "/" << v << "/" << w
249  << " step_length " << step_length
250  << " nbOfPhotons " << nbOfPhotons;
251 #endif
252  if (nbOfPhotons < 0) {
253  return 0;
254  } else if (nbOfPhotons > 0) {
255  double w_ph=0;
256  for (int i = 0; i < nbOfPhotons; i++) {
257  double rand = G4UniformRand();
258  double theta_C = acos(1./(pBeta*ref_index));
259  double phi_C = 2*M_PI*rand;
260  double sinTheta = sin(theta_C);
261  double cosTheta = cos(theta_C);
262  double cosPhi = cos(phi_C);
263  //photon momentum
264  if (uv < 0.001) { // aligned with z-axis
265  w_ph = cosTheta;
266  } else { // general case
267  w_ph = w * cosTheta - sinTheta * cosPhi * uv;
268  }
269  double r_lambda = G4UniformRand();
270  double lambda0 = (lambda1 * lambda2) / (lambda2 - r_lambda *
271  (lambda2 - lambda1));
272  double lambda = (lambda0/cm) * pow(double(10),7); // lambda is in nm
273  wlini.push_back(lambda);
274 #ifdef DebugLog
275  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: " <<i <<" lambda "
276  << lambda << " w_ph " << w_ph << " aperture "
277  << aperture;
278 #endif
279  if (w_ph > aperture) { // phton trapped inside PMT glass
280  wltrap.push_back(lambda);
281  rand = G4UniformRand();
282 #ifdef DebugLog
283  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: Random " << rand
284  << " Survive? " << (rand < 1.);
285 #endif
286  if (rand < 1.0) { // survived all the times and sent to photo-cathode
287  wlatten.push_back(lambda);
288  double qEffic = computeQEff(lambda);//Quantum efficiency of the PMT
289  rand = G4UniformRand();
290 #ifdef DebugLog
291  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: qEffic "
292  << qEffic << " Random " << rand
293  << " Survive? " <<(rand < qEffic);
294 #endif
295  if (rand < qEffic) { // made photoelectron
296  npe_ += 1;
297  momZ.push_back(w_ph);
298  wl.push_back(lambda);
299  wlqeff.push_back(lambda);
300  } // made pe
301  } // accepted all Cherenkov photons
302  } // end of if(w_ph < w_aperture), trapped inside glass
303  }// end of ++NbOfPhotons
304  } // end of if(NbOfPhotons)}
305 #ifdef DebugLog
306  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: npe " << npe_;
307 #endif
308  return npe_;
309 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
double computeQEff(double wavelength)
Definition: HFCherenkov.cc:366
const double w
Definition: UKUtility.cc:23
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< double > wlqeff
Definition: HFCherenkov.h:74
void clearVectors()
Definition: HFCherenkov.cc:446
std::vector< double > wlini
Definition: HFCherenkov.h:70
T sqrt(T t)
Definition: SSEVec.h:18
std::vector< double > wltrap
Definition: HFCherenkov.h:71
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double aperture
Definition: HFCherenkov.h:62
#define M_PI
int computeNbOfPhotons(double pBeta, double step_length)
Definition: HFCherenkov.cc:346
Signal rand(Signal arg)
Definition: vlib.cc:442
bool isApplicable(const G4ParticleDefinition *aParticleType)
Definition: HFCherenkov.cc:457
std::vector< double > wlatten
Definition: HFCherenkov.h:72
double lambda1
Definition: HFCherenkov.h:61
double lambda2
Definition: HFCherenkov.h:61
std::vector< double > momZ
Definition: HFCherenkov.h:69
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double ref_index
Definition: HFCherenkov.h:60
std::vector< double > wl
Definition: HFCherenkov.h:68
int HFCherenkov::computeNPhTrapped ( double  pBeta,
double  u,
double  v,
double  w,
double  step_length,
double  zFiber,
double  Dose,
int  Npe_Dose 
)

Definition at line 46 of file HFCherenkov.cc.

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

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

Definition at line 366 of file HFCherenkov.cc.

References create_public_lumi_plots::exp, RecoJets_EventContent_cff::func, LogDebug, funct::pow(), UseNewPMT, and y.

Referenced by computeNPE(), and computeNPEinPMT().

366  {
367 
368  double qeff(0.);
369  if (UseNewPMT) {
370  if (wavelength<=350) {
371  qeff=2.45867*(TMath::Landau(wavelength,353.820,59.1324));
372  } else if (wavelength>350 && wavelength<500) {
373  qeff= 0.441989*exp(-pow((wavelength-358.371),2)/(2*pow((138.277),2)));
374  } else if (wavelength>=500 && wavelength<550) {
375  qeff= 0.271862*exp(-pow((wavelength-491.505),2)/(2*pow((47.0418),2)));
376  } else if (wavelength>=550) {
377  qeff= 0.137297*exp(-pow((wavelength-520.260),2)/(2*pow((75.5023),2)));
378  }
379 #ifdef DebugLog
380  LogDebug("HFShower") << "HFCherenkov:: for new PMT : wavelength === "
381  << wavelength << "\tqeff ===\t" << qeff;
382 #endif
383  } else {
384  double y = (wavelength - 275.) /180.;
385  double func = 1. / (1. + 250.*pow((y/5.),4));
386  double qE_R7525 = 0.77 * y * exp(-y) * func ;
387  qeff = qE_R7525;
388 #ifdef DebugLog
389  LogDebug("HFShower") << "HFCherenkov:: for old PMT : wavelength === "
390  << wavelength << "; qeff = " << qeff;
391 #endif
392  }
393 
394 #ifdef DebugLog
395  LogDebug("HFShower") << "HFCherenkov::computeQEff: wavelength " << wavelength
396  << " y/func " << y << "/" << func << " qeff " << qeff;
397 #endif
398  return qeff;
399 }
#define LogDebug(id)
bool UseNewPMT
Definition: HFCherenkov.h:65
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
std::vector< double > HFCherenkov::getMom ( )

Definition at line 341 of file HFCherenkov.cc.

References momZ, and findQualityFiles::v.

Referenced by HFShower::getHits().

341  {
342  std::vector<double> v = momZ;
343  return v;
344 }
std::vector< double > momZ
Definition: HFCherenkov.h:69
std::vector< double > HFCherenkov::getWL ( )

Definition at line 336 of file HFCherenkov.cc.

References findQualityFiles::v, and wl.

Referenced by HFShower::getHits().

336  {
337  std::vector<double> v = wl;
338  return v;
339 }
std::vector< double > wl
Definition: HFCherenkov.h:68
std::vector< double > HFCherenkov::getWLAtten ( )

Definition at line 321 of file HFCherenkov.cc.

References findQualityFiles::v, and wlatten.

321  {
322  std::vector<double> v = wlatten;
323  return v;
324 }
std::vector< double > wlatten
Definition: HFCherenkov.h:72
std::vector< double > HFCherenkov::getWLHEM ( )

Definition at line 326 of file HFCherenkov.cc.

References findQualityFiles::v, and wlhem.

326  {
327  std::vector<double> v = wlhem;
328  return v;
329 }
std::vector< double > wlhem
Definition: HFCherenkov.h:73
std::vector< double > HFCherenkov::getWLIni ( )

Definition at line 311 of file HFCherenkov.cc.

References findQualityFiles::v, and wlini.

311  {
312  std::vector<double> v = wlini;
313  return v;
314 }
std::vector< double > wlini
Definition: HFCherenkov.h:70
std::vector< double > HFCherenkov::getWLQEff ( )

Definition at line 331 of file HFCherenkov.cc.

References findQualityFiles::v, and wlqeff.

331  {
332  std::vector<double> v = wlqeff;
333  return v;
334 }
std::vector< double > wlqeff
Definition: HFCherenkov.h:74
std::vector< double > HFCherenkov::getWLTrap ( )

Definition at line 316 of file HFCherenkov.cc.

References findQualityFiles::v, and wltrap.

316  {
317  std::vector<double> v = wltrap;
318  return v;
319 }
std::vector< double > wltrap
Definition: HFCherenkov.h:71
bool HFCherenkov::isApplicable ( const G4ParticleDefinition *  aParticleType)
private

Definition at line 457 of file HFCherenkov.cc.

References LogDebug, and tmp.

Referenced by computeNPE(), and computeNPEinPMT().

457  {
458  bool tmp = (aParticleType->GetPDGCharge() != 0);
459 #ifdef DebugLog
460  LogDebug("HFShower") << "HFCherenkov::isApplicable: aParticleType "
461  << aParticleType->GetParticleName() << " PDGCharge "
462  << aParticleType->GetPDGCharge() << " Result " << tmp;
463 #endif
464  return tmp;
465 }
#define LogDebug(id)
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
double HFCherenkov::smearNPE ( G4int  Npe)

Definition at line 431 of file HFCherenkov.cc.

References gain, i, LogDebug, and SpecificationBuilder_cfi::val().

431  {
432 
433  double pe = 0.;
434  if (npe > 0) {
435  for (int i = 0; i < npe; ++i) {
436  double val = G4Poisson(gain);
437  pe += (val/gain) + 0.001*G4UniformRand();
438  }
439  }
440 #ifdef DebugLog
441  LogDebug("HFShower") << "HFCherenkov::smearNPE: npe " << npe << " pe " << pe;
442 #endif
443  return pe;
444 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
double gain
Definition: HFCherenkov.h:63

Member Data Documentation

double HFCherenkov::aperture
private

Definition at line 62 of file HFCherenkov.h.

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

double HFCherenkov::apertureTrap
private

Definition at line 62 of file HFCherenkov.h.

Referenced by computeNPhTrapped(), and HFCherenkov().

double HFCherenkov::aperturetrapped
private

Definition at line 62 of file HFCherenkov.h.

Referenced by computeNPE(), and HFCherenkov().

bool HFCherenkov::checkSurvive
private

Definition at line 64 of file HFCherenkov.h.

Referenced by computeNPE(), and HFCherenkov().

double HFCherenkov::fibreR
private

Definition at line 63 of file HFCherenkov.h.

Referenced by computeNPE(), and HFCherenkov().

double HFCherenkov::gain
private

Definition at line 63 of file HFCherenkov.h.

Referenced by HFCherenkov(), and smearNPE().

double HFCherenkov::lambda1
private

Definition at line 61 of file HFCherenkov.h.

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

double HFCherenkov::lambda2
private

Definition at line 61 of file HFCherenkov.h.

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

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

Definition at line 69 of file HFCherenkov.h.

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

G4ThreeVector HFCherenkov::phMom
private

Definition at line 67 of file HFCherenkov.h.

double HFCherenkov::ref_index
private
double HFCherenkov::sinPsimax
private

Definition at line 63 of file HFCherenkov.h.

Referenced by computeNPE(), and HFCherenkov().

bool HFCherenkov::UseNewPMT
private

Definition at line 65 of file HFCherenkov.h.

Referenced by computeQEff(), and HFCherenkov().

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

Definition at line 68 of file HFCherenkov.h.

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

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

Definition at line 72 of file HFCherenkov.h.

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

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

Definition at line 73 of file HFCherenkov.h.

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

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

Definition at line 70 of file HFCherenkov.h.

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

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

Definition at line 74 of file HFCherenkov.h.

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

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

Definition at line 71 of file HFCherenkov.h.

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