CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/SimG4CMS/Calo/src/HFCherenkov.cc

Go to the documentation of this file.
00001 
00002 // File:  HFCherenkov.cc
00003 // Description: Generate Cherenkov photons for HF
00005 
00006 #include "SimG4CMS/Calo/interface/HFCherenkov.h"
00007 
00008 #include "G4Poisson.hh"
00009 #include "G4ParticleDefinition.hh"
00010 
00011 //#define DebugLog
00012 
00013 HFCherenkov::HFCherenkov(edm::ParameterSet const & m_HF) {
00014 
00015   //Simple configurables
00016   //static SimpleConfigurable<double> p1(1.459, "HFCherenkov:RefIndex");
00017   //static SimpleConfigurable<double> p2(280.0, "HFCherenkov:Lambda1");
00018   //static SimpleConfigurable<double> p3(700.0, "HFCherenkov:Lambda2");
00019   //static SimpleConfigurable<double> p4(0.33,  "HFCherenkov:Aperture");
00020   //static SimpleConfigurable<double> p5(0.22,  "HFCherenkov:ApertureTrapped");
00021   //static SimpleConfigurable<double> p6(0.33,  "HFCherenkov:Gain");
00022   //static SimpleConfigurable<bool>   p7(false, "HFCherenkov:CheckSurvive");
00023 
00024   ref_index    = m_HF.getParameter<double>("RefIndex");
00025   lambda1      = ((m_HF.getParameter<double>("Lambda1"))/pow(double(10),7))*cm;
00026   lambda2      = ((m_HF.getParameter<double>("Lambda2"))/pow(double(10),7))*cm;
00027   aperture     = cos(asin(m_HF.getParameter<double>("Aperture")));
00028   apertureTrap = cos(asin(m_HF.getParameter<double>("ApertureTrapped")));
00029   gain         = m_HF.getParameter<double>("Gain");
00030   checkSurvive = m_HF.getParameter<bool>("CheckSurvive");
00031 
00032   edm::LogInfo("HFShower") << "HFCherenkov:: initialised with ref_index " 
00033                            << ref_index << " lambda1/lambda2 (cm) " 
00034                            << lambda1/cm << "/" << lambda2/cm
00035                            << " aperture(total/trapped) " << aperture << "/"
00036                            << apertureTrap << " Check photon survival in HF " 
00037                            << checkSurvive << " Gain " << gain;
00038 
00039   clearVectors();
00040 }
00041 
00042 HFCherenkov::~HFCherenkov() {}
00043 
00044 int HFCherenkov::computeNPhTrapped(double pBeta, 
00045                                    double u, double v, double w, 
00046                                    double step_length, double zFiber,
00047                                    double dose, int npe_Dose) {
00048 
00049   if (pBeta < (1/ref_index) || step_length < 0.0001) {return 0;}
00050 
00051   double uv = sqrt(u*u + v*v);
00052   int nbOfPhotons = computeNbOfPhotons(pBeta, step_length);
00053 
00054   if (nbOfPhotons < 0) {
00055     return 0;
00056   } else if (nbOfPhotons > 0) {
00057     double w_ph=0;
00058     for (int i = 0; i < nbOfPhotons; i++) {
00059       double rand     = G4UniformRand();
00060       double theta_C  = acos(1./(pBeta*ref_index));
00061       double phi_C    = 2*M_PI*rand;
00062       double sinTheta = sin(theta_C);
00063       double cosTheta = cos(theta_C);
00064       double cosPhi   = cos(phi_C);
00065       //photon momentum
00066       if (uv < 0.001) { // aligned with z-axis
00067         w_ph = cosTheta;
00068       } else { // general case
00069         w_ph = w * cosTheta  - sinTheta * cosPhi * uv;
00070       }
00071       if (w_ph > apertureTrap) { // phton trapped inside fiber
00072         npe_Dose += 1; 
00073       }
00074     }
00075   }
00076   int n_photons = npe_Dose;
00077   return n_photons;
00078 }
00079 
00080 int HFCherenkov::computeNPE(G4ParticleDefinition* pDef, double pBeta,
00081                             double u, double v, double w,
00082                             double step_length, double zFiber,
00083                             double dose, int npe_Dose) {
00084 
00085   clearVectors();
00086   if (!isApplicable(pDef)) {return 0;}
00087   if (pBeta < (1/ref_index) || step_length < 0.0001) {
00088 #ifdef DebugLog
00089     LogDebug("HFShower") << "HFCherenkov::computeNPE: pBeta " << pBeta 
00090                          << " 1/mu " << (1/ref_index) << " step_length " 
00091                          << step_length;
00092 #endif
00093     return 0;
00094   }
00095    
00096   double uv = sqrt(u*u + v*v);
00097   int nbOfPhotons = computeNbOfPhotons(pBeta, step_length);
00098 #ifdef DebugLog
00099   LogDebug("HFShower") << "HFCherenkov::computeNPE: pBeta " << pBeta 
00100                        << " u/v/w " << u << "/" << v << "/" << w 
00101                        << " step_length " << step_length << " zFib " << zFiber 
00102                        << " nbOfPhotons " << nbOfPhotons;
00103 #endif
00104   if (nbOfPhotons < 0) {
00105     return 0;
00106   } else if (nbOfPhotons > 0) {
00107     double w_ph=0;
00108     for (int i = 0; i < nbOfPhotons; i++) {
00109       double rand     = G4UniformRand();
00110       double theta_C  = acos(1./(pBeta*ref_index));
00111       double phi_C    = 2*M_PI*rand;
00112       double sinTheta = sin(theta_C);
00113       double cosTheta = cos(theta_C);
00114       double cosPhi   = cos(phi_C);
00115       //photon momentum
00116       if (uv < 0.001) { // aligned with z-axis
00117         w_ph = cosTheta;
00118       } else { // general case
00119         w_ph = w * cosTheta - sinTheta * cosPhi * uv;
00120       }
00121       double r_lambda = G4UniformRand();
00122       double lambda0 = (lambda1 * lambda2) / (lambda2 - r_lambda *
00123                                               (lambda2 - lambda1));
00124       double lambda  = (lambda0/cm) * pow(double(10),7); // lambda is in nm
00125       wlini.push_back(lambda);
00126 #ifdef DebugLog
00127       LogDebug("HFShower") << "HFCherenkov::computeNPE: " << i << " lambda "
00128                            << lambda << " w_ph " << w_ph << " aperture "
00129                            << aperture;
00130 #endif
00131       if (w_ph > aperture) { // phton trapped inside fiber
00132         wltrap.push_back(lambda);
00133         double prob_HF  = 1.0; //photon survived in HF
00134         double a0_inv   = 0.1234;  //meter^-1
00135         double prob_MX  = exp( - 0.5 * a0_inv ); //light mixer
00136         if (checkSurvive) {
00137           double a_inv = a0_inv + 0.14 * pow(dose,0.30);
00138           double z_meters = zFiber;
00139           prob_HF  = exp(-z_meters * a_inv ); //photon survived in HF
00140         }
00141         rand = G4UniformRand();
00142 #ifdef DebugLog
00143         LogDebug("HFShower") << "HFCherenkov::computeNPE: probHF " << prob_HF
00144                              << " prob_MX " << prob_MX << " Random " << rand 
00145                              << " Survive? " << (rand < (prob_HF * prob_MX));
00146 #endif
00147         if (rand < (prob_HF * prob_MX)) { // survived and sent to light mixer
00148           wlatten.push_back(lambda);
00149           rand = G4UniformRand();
00150           double effHEM = computeHEMEff(lambda);
00151 #ifdef DebugLog
00152           LogDebug("HFShower") << "HFCherenkov::computeNPE: w_ph " << w_ph 
00153                                << " effHEM " << effHEM << " Random " << rand 
00154                                << " Survive? " << (w_ph>0.997||(rand<effHEM));
00155 #endif
00156           if (w_ph>0.997 || (rand<effHEM)) { // survived HEM
00157             wlhem.push_back(lambda);
00158             double qEffic = computeQEff(lambda);
00159             rand = G4UniformRand();
00160 #ifdef DebugLog
00161             LogDebug("HFShower") << "HFCherenkov::computeNPE: qEffic "
00162                                  << qEffic << " Random " << rand
00163                                  << " Survive? " <<(rand < qEffic);
00164 #endif
00165             if (rand < qEffic) { // made photoelectron
00166               npe_Dose += 1;
00167               momZ.push_back(w_ph);
00168               wl.push_back(lambda);
00169               wlqeff.push_back(lambda);
00170             } // made pe
00171           } // passed HEM
00172         } // passed fiber
00173       } // end of  if(w_ph < w_aperture), trapped inside fiber
00174     }// end of ++NbOfPhotons
00175   } // end of if(NbOfPhotons)}
00176   int npe =  npe_Dose; // Nb of photoelectrons
00177 #ifdef DebugLog
00178   LogDebug("HFShower") << "HFCherenkov::computeNPE: npe " << npe;
00179 #endif
00180   return npe;
00181 }
00182 
00183 int HFCherenkov::computeNPEinPMT(G4ParticleDefinition* pDef, double pBeta, 
00184                                  double u, double v, double w, 
00185                                  double step_length){
00186   clearVectors();
00187   int npe_ = 0;
00188   if (!isApplicable(pDef)) {return 0;}
00189   if (pBeta < (1/ref_index) || step_length < 0.0001) {
00190 #ifdef DebugLog
00191     LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: pBeta " << pBeta 
00192                          << " 1/mu " << (1/ref_index) << " step_length " 
00193                          << step_length;
00194 #endif
00195     return 0;
00196   }
00197    
00198   double uv = sqrt(u*u + v*v);
00199   int nbOfPhotons = computeNbOfPhotons(pBeta, step_length);
00200 #ifdef DebugLog
00201   LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: pBeta " << pBeta 
00202                        << " u/v/w " << u << "/" << v << "/" << w 
00203                        << " step_length " << step_length  
00204                        << " nbOfPhotons " << nbOfPhotons;
00205 #endif   
00206   if (nbOfPhotons < 0) {
00207     return 0;
00208   } else if (nbOfPhotons > 0) {
00209     double w_ph=0;
00210     for (int i = 0; i < nbOfPhotons; i++) {
00211       double rand     = G4UniformRand();
00212       double theta_C  = acos(1./(pBeta*ref_index));
00213       double phi_C    = 2*M_PI*rand;
00214       double sinTheta = sin(theta_C);
00215       double cosTheta = cos(theta_C);
00216       double cosPhi   = cos(phi_C); 
00217       //photon momentum
00218       if (uv < 0.001) { // aligned with z-axis
00219         w_ph = cosTheta;
00220       } else { // general case
00221         w_ph = w * cosTheta - sinTheta * cosPhi * uv;
00222       }
00223       double r_lambda = G4UniformRand();
00224       double lambda0 = (lambda1 * lambda2) / (lambda2 - r_lambda *
00225                                               (lambda2 - lambda1));
00226       double lambda  = (lambda0/cm) * pow(double(10),7); // lambda is in nm
00227       wlini.push_back(lambda);
00228 #ifdef DebugLog
00229       LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: " <<i <<" lambda "
00230                            << lambda << " w_ph " << w_ph << " aperture " 
00231                            << aperture;
00232 #endif
00233       if (w_ph > aperture) { // phton trapped inside PMT glass
00234         wltrap.push_back(lambda);
00235         rand = G4UniformRand();
00236 #ifdef DebugLog
00237         LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: Random " << rand
00238                              << " Survive? " << (rand < 1.);
00239 #endif
00240         if (rand < 1.0) { // survived all the times and sent to photo-cathode
00241           wlatten.push_back(lambda);
00242           rand = G4UniformRand();
00243           double qEffic = computeQEff(lambda);//Quantum efficiency of the PMT
00244           rand = G4UniformRand();
00245 #ifdef DebugLog
00246           LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: qEffic " 
00247                                << qEffic << " Random " << rand 
00248                                << " Survive? " <<(rand < qEffic);
00249 #endif
00250           if (rand < qEffic) { // made photoelectron
00251             npe_ += 1;
00252             momZ.push_back(w_ph);
00253             wl.push_back(lambda);
00254             wlqeff.push_back(lambda);
00255           } // made pe
00256         } // accepted all Cherenkov photons
00257       } // end of  if(w_ph < w_aperture), trapped inside glass
00258     }// end of ++NbOfPhotons
00259   } // end of if(NbOfPhotons)}
00260 #ifdef DebugLog
00261   LogDebug("HFShower") << "HFCherenkov::computeNPEinPMT: npe " << npe_;
00262 #endif
00263   return npe_;
00264 }
00265 
00266 std::vector<double>  HFCherenkov::getWLIni() {
00267   std::vector<double> v = wlini;
00268   return v;
00269 }
00270 
00271 std::vector<double>  HFCherenkov::getWLTrap() {
00272   std::vector<double> v = wltrap;
00273   return v;
00274 }
00275 
00276 std::vector<double>  HFCherenkov::getWLAtten() {
00277   std::vector<double> v = wlatten;
00278   return v;
00279 }
00280 
00281 std::vector<double>  HFCherenkov::getWLHEM() {
00282   std::vector<double> v  = wlhem;
00283   return v;
00284 }
00285 
00286 std::vector<double>  HFCherenkov::getWLQEff() {
00287   std::vector<double> v = wlqeff;
00288   return v;
00289 }
00290 
00291 std::vector<double>  HFCherenkov::getWL() {
00292   std::vector<double> v = wl;
00293   return v;
00294 }
00295 
00296 std::vector<double>  HFCherenkov::getMom() {
00297   std::vector<double> v = momZ;
00298   return v;
00299 }
00300 
00301 int HFCherenkov::computeNbOfPhotons(double beta, G4double stepL) {
00302 
00303   double pBeta = beta;
00304   double alpha = 0.0073;
00305   double step_length = stepL;
00306   double theta_C = acos(1./(pBeta*ref_index));
00307   double lambdaDiff = (1./lambda1 - 1./lambda2);
00308   double cherenPhPerLength = 2 * M_PI * alpha * lambdaDiff*cm;
00309   double d_NOfPhotons = cherenPhPerLength * sin(theta_C)*sin(theta_C) *  (step_length/cm);
00310   int nbOfPhotons = int(d_NOfPhotons);
00311 #ifdef DebugLog
00312   LogDebug("HFShower") << "HFCherenkov::computeNbOfPhotons: StepLength " 
00313                        << step_length << " theta_C " << theta_C 
00314                        << " lambdaDiff " << lambdaDiff
00315                        << " cherenPhPerLength " << cherenPhPerLength 
00316                        << " Photons " << d_NOfPhotons << " " << nbOfPhotons;
00317 #endif
00318   return nbOfPhotons;
00319 }
00320 
00321 double HFCherenkov::computeQEff(double wavelength) {
00322 
00323   double y        = (wavelength - 275.) /180.;
00324   double func     = 1. / (1. + 250.*pow((y/5.),4));
00325   double qE_R7525 = 0.77 * y * exp(-y) * func ;
00326   double qeff     = qE_R7525;
00327 #ifdef DebugLog
00328   LogDebug("HFShower") << "HFCherenkov::computeQEff: wavelength " << wavelength
00329                        << " y/func " << y << "/" << func << " qeff " << qeff;
00330 #endif
00331   return qeff;
00332 }
00333 
00334 double HFCherenkov::computeHEMEff(double wavelength) {
00335 
00336   double hEMEff = 0;
00337   if (wavelength < 400.) {
00338     hEMEff = 0.;
00339   } else if (wavelength >= 400. && wavelength < 410.) {
00340     //hEMEff = .99 * (wavelength - 400.) / 10.;
00341     hEMEff = (-1322.453 / wavelength) + 4.2056;
00342   } else if (wavelength >= 410.) {
00343     hEMEff = 0.99;
00344     if (wavelength > 415. && wavelength < 445.) {
00345       //abs(wavelength - 430.) < 15.
00346       //hEMEff = 0.95;
00347       hEMEff = 0.97;
00348     } else if (wavelength > 550. && wavelength < 600.) {
00349       // abs(wavelength - 575.) < 25.)
00350       //hEMEff = 0.96;
00351       hEMEff = 0.98;
00352     } else if (wavelength > 565. && wavelength <= 635.) { // added later
00353       // abs(wavelength - 600.) < 35.)
00354       hEMEff = (701.7268 / wavelength) - 0.186;
00355     }
00356   }
00357 #ifdef DebugLog
00358   LogDebug("HFShower") << "HFCherenkov::computeHEMEff: wavelength "
00359                        << wavelength << " hEMEff " << hEMEff;
00360 #endif
00361   return hEMEff;
00362 }
00363 
00364 double HFCherenkov::smearNPE(int npe) {
00365 
00366   double pe = 0.;
00367   if (npe > 0) {
00368     for (int i = 0; i < npe; ++i) {
00369       double val =  G4Poisson(gain);
00370       pe += (val/gain) + 0.001*G4UniformRand();
00371     }
00372   }
00373 #ifdef DebugLog
00374   LogDebug("HFShower") << "HFCherenkov::smearNPE: npe " << npe << " pe " << pe;
00375 #endif
00376   return pe;
00377 }
00378 
00379 void HFCherenkov::clearVectors() {
00380 
00381   wl.clear();
00382   wlini.clear();
00383   wltrap.clear();
00384   wlatten.clear();
00385   wlhem.clear();
00386   wlqeff.clear();
00387   momZ.clear();
00388 }
00389 
00390 bool HFCherenkov::isApplicable(const G4ParticleDefinition* aParticleType) {
00391   bool tmp = (aParticleType->GetPDGCharge() != 0);
00392 #ifdef DebugLog
00393   LogDebug("HFShower") << "HFCherenkov::isApplicable: aParticleType " 
00394                        << aParticleType->GetParticleName() << " PDGCharge " 
00395                        << aParticleType->GetPDGCharge() << " Result " << tmp;
00396 #endif
00397   return tmp;
00398 }