test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HFCherenkov.cc
Go to the documentation of this file.
1 // File: HFCherenkov.cc
3 // Description: Generate Cherenkov photons for HF
5 
7 
8 #include "G4Poisson.hh"
9 #include "G4ParticleDefinition.hh"
10 #include "G4NavigationHistory.hh"
11 #include "TMath.h"
12 
13 #include "Randomize.hh"
14 #include "G4SystemOfUnits.hh"
15 
16 //#define DebugLog
17 
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 }
43 
45 
47  double u, double v, double w,
48  double step_length, double zFiber,
49  double dose, int npe_Dose) {
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 }
81 
82 int HFCherenkov::computeNPE(G4Step * aStep, G4ParticleDefinition* pDef,
83  double pBeta, double u, double v, double w,
84  double step_length, double zFiber,
85  double dose, int npe_Dose) {
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 }
228 
229 int HFCherenkov::computeNPEinPMT(G4ParticleDefinition* pDef, double pBeta,
230  double u, double v, double w,
231  double step_length){
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  rand = G4UniformRand();
289  double qEffic = computeQEff(lambda);//Quantum efficiency of the PMT
290  rand = G4UniformRand();
291 #ifdef DebugLog
292  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: qEffic "
293  << qEffic << " Random " << rand
294  << " Survive? " <<(rand < qEffic);
295 #endif
296  if (rand < qEffic) { // made photoelectron
297  npe_ += 1;
298  momZ.push_back(w_ph);
299  wl.push_back(lambda);
300  wlqeff.push_back(lambda);
301  } // made pe
302  } // accepted all Cherenkov photons
303  } // end of if(w_ph < w_aperture), trapped inside glass
304  }// end of ++NbOfPhotons
305  } // end of if(NbOfPhotons)}
306 #ifdef DebugLog
307  LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: npe " << npe_;
308 #endif
309  return npe_;
310 }
311 
312 std::vector<double> HFCherenkov::getWLIni() {
313  std::vector<double> v = wlini;
314  return v;
315 }
316 
317 std::vector<double> HFCherenkov::getWLTrap() {
318  std::vector<double> v = wltrap;
319  return v;
320 }
321 
322 std::vector<double> HFCherenkov::getWLAtten() {
323  std::vector<double> v = wlatten;
324  return v;
325 }
326 
327 std::vector<double> HFCherenkov::getWLHEM() {
328  std::vector<double> v = wlhem;
329  return v;
330 }
331 
332 std::vector<double> HFCherenkov::getWLQEff() {
333  std::vector<double> v = wlqeff;
334  return v;
335 }
336 
337 std::vector<double> HFCherenkov::getWL() {
338  std::vector<double> v = wl;
339  return v;
340 }
341 
342 std::vector<double> HFCherenkov::getMom() {
343  std::vector<double> v = momZ;
344  return v;
345 }
346 
347 int HFCherenkov::computeNbOfPhotons(double beta, G4double stepL) {
348 
349  double pBeta = beta;
350  double alpha = 0.0073;
351  double step_length = stepL;
352  double theta_C = acos(1./(pBeta*ref_index));
353  double lambdaDiff = (1./lambda1 - 1./lambda2);
354  double cherenPhPerLength = 2 * M_PI * alpha * lambdaDiff*cm;
355  double d_NOfPhotons = cherenPhPerLength * sin(theta_C)*sin(theta_C) * (step_length/cm);
356  int nbOfPhotons = int(d_NOfPhotons);
357 #ifdef DebugLog
358  LogDebug("HFShower") << "HFCherenkov::computeNbOfPhotons: StepLength "
359  << step_length << " theta_C " << theta_C
360  << " lambdaDiff " << lambdaDiff
361  << " cherenPhPerLength " << cherenPhPerLength
362  << " Photons " << d_NOfPhotons << " " << nbOfPhotons;
363 #endif
364  return nbOfPhotons;
365 }
366 
367 double HFCherenkov::computeQEff(double wavelength) {
368 
369  double qeff(0.);
370  if (UseNewPMT) {
371  if (wavelength<=350) {
372  qeff=2.45867*(TMath::Landau(wavelength,353.820,59.1324));
373  } else if (wavelength>350 && wavelength<500) {
374  qeff= 0.441989*exp(-pow((wavelength-358.371),2)/(2*pow((138.277),2)));
375  } else if (wavelength>=500 && wavelength<550) {
376  qeff= 0.271862*exp(-pow((wavelength-491.505),2)/(2*pow((47.0418),2)));
377  } else if (wavelength>=550) {
378  qeff= 0.137297*exp(-pow((wavelength-520.260),2)/(2*pow((75.5023),2)));
379  }
380 #ifdef DebugLog
381  LogDebug("HFShower") << "HFCherenkov:: for new PMT : wavelength === "
382  << wavelength << "\tqeff ===\t" << qeff;
383 #endif
384  } else {
385  double y = (wavelength - 275.) /180.;
386  double func = 1. / (1. + 250.*pow((y/5.),4));
387  double qE_R7525 = 0.77 * y * exp(-y) * func ;
388  qeff = qE_R7525;
389 #ifdef DebugLog
390  LogDebug("HFShower") << "HFCherenkov:: for old PMT : wavelength === "
391  << wavelength << "; qeff = " << qeff;
392 #endif
393  }
394 
395 #ifdef DebugLog
396  LogDebug("HFShower") << "HFCherenkov::computeQEff: wavelength " << wavelength
397  << " y/func " << y << "/" << func << " qeff " << qeff;
398 #endif
399  return qeff;
400 }
401 
402 double HFCherenkov::computeHEMEff(double wavelength) {
403 
404  double hEMEff = 0;
405  if (wavelength < 400.) {
406  hEMEff = 0.;
407  } else if (wavelength >= 400. && wavelength < 410.) {
408  //hEMEff = .99 * (wavelength - 400.) / 10.;
409  hEMEff = (-1322.453 / wavelength) + 4.2056;
410  } else if (wavelength >= 410.) {
411  hEMEff = 0.99;
412  if (wavelength > 415. && wavelength < 445.) {
413  //abs(wavelength - 430.) < 15.
414  //hEMEff = 0.95;
415  hEMEff = 0.97;
416  } else if (wavelength > 550. && wavelength < 600.) {
417  // abs(wavelength - 575.) < 25.)
418  //hEMEff = 0.96;
419  hEMEff = 0.98;
420  } else if (wavelength > 565. && wavelength <= 635.) { // added later
421  // abs(wavelength - 600.) < 35.)
422  hEMEff = (701.7268 / wavelength) - 0.186;
423  }
424  }
425 #ifdef DebugLog
426  LogDebug("HFShower") << "HFCherenkov::computeHEMEff: wavelength "
427  << wavelength << " hEMEff " << hEMEff;
428 #endif
429  return hEMEff;
430 }
431 
432 double HFCherenkov::smearNPE(int npe) {
433 
434  double pe = 0.;
435  if (npe > 0) {
436  for (int i = 0; i < npe; ++i) {
437  double val = G4Poisson(gain);
438  pe += (val/gain) + 0.001*G4UniformRand();
439  }
440  }
441 #ifdef DebugLog
442  LogDebug("HFShower") << "HFCherenkov::smearNPE: npe " << npe << " pe " << pe;
443 #endif
444  return pe;
445 }
446 
448 
449  wl.clear();
450  wlini.clear();
451  wltrap.clear();
452  wlatten.clear();
453  wlhem.clear();
454  wlqeff.clear();
455  momZ.clear();
456 }
457 
458 bool HFCherenkov::isApplicable(const G4ParticleDefinition* aParticleType) {
459  bool tmp = (aParticleType->GetPDGCharge() != 0);
460 #ifdef DebugLog
461  LogDebug("HFShower") << "HFCherenkov::isApplicable: aParticleType "
462  << aParticleType->GetParticleName() << " PDGCharge "
463  << aParticleType->GetPDGCharge() << " Result " << tmp;
464 #endif
465  return tmp;
466 }
#define LogDebug(id)
const double beta
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
double computeQEff(double wavelength)
Definition: HFCherenkov.cc:367
float alpha
Definition: AMPTWrapper.h:95
virtual ~HFCherenkov()
Definition: HFCherenkov.cc:44
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
HFCherenkov(edm::ParameterSet const &p)
Definition: HFCherenkov.cc:18
std::vector< double > wlqeff
Definition: HFCherenkov.h:74
void clearVectors()
Definition: HFCherenkov.cc:447
int computeNPEinPMT(G4ParticleDefinition *pDef, double pBeta, double u, double v, double w, double step_length)
Definition: HFCherenkov.cc:229
bool UseNewPMT
Definition: HFCherenkov.h:65
std::vector< double > wlini
Definition: HFCherenkov.h:70
double computeHEMEff(double wavelength)
Definition: HFCherenkov.cc:402
double smearNPE(G4int Npe)
Definition: HFCherenkov.cc:432
double sinPsimax
Definition: HFCherenkov.h:63
std::vector< double > getWL()
Definition: HFCherenkov.cc:337
std::vector< double > getWLQEff()
Definition: HFCherenkov.cc:332
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
double apertureTrap
Definition: HFCherenkov.h:62
std::vector< double > getMom()
Definition: HFCherenkov.cc:342
#define M_PI
std::vector< double > getWLIni()
Definition: HFCherenkov.cc:312
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
int computeNbOfPhotons(double pBeta, double step_length)
Definition: HFCherenkov.cc:347
std::vector< double > getWLTrap()
Definition: HFCherenkov.cc:317
Signal rand(Signal arg)
Definition: vlib.cc:442
std::vector< double > getWLAtten()
Definition: HFCherenkov.cc:322
bool isApplicable(const G4ParticleDefinition *aParticleType)
Definition: HFCherenkov.cc:458
std::vector< double > wlatten
Definition: HFCherenkov.h:72
double gain
Definition: HFCherenkov.h:63
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
int computeNPhTrapped(double pBeta, double u, double v, double w, double step_length, double zFiber, double Dose, int Npe_Dose)
Definition: HFCherenkov.cc:46
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double ref_index
Definition: HFCherenkov.h:60
int computeNPE(G4Step *step, G4ParticleDefinition *pDef, double pBeta, double u, double v, double w, double step_length, double zFiber, double Dose, int Npe_Dose)
Definition: HFCherenkov.cc:82
std::vector< double > getWLHEM()
Definition: HFCherenkov.cc:327
std::vector< double > wl
Definition: HFCherenkov.h:68