CMS 3D CMS Logo

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 HFCherenkov::HFCherenkov(edm::ParameterSet const & p) {
00012 
00013   edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFCherenkov");
00014   //Simple configurables
00015   //static SimpleConfigurable<double> p1(1.459, "HFCherenkov:RefIndex");
00016   //static SimpleConfigurable<double> p2(280.0, "HFCherenkov:Lambda1");
00017   //static SimpleConfigurable<double> p3(700.0, "HFCherenkov:Lambda2");
00018   //static SimpleConfigurable<double> p4(0.33,  "HFCherenkov:Aperture");
00019   //static SimpleConfigurable<double> p5(0.22,  "HFCherenkov:ApertureTrapped");
00020   //static SimpleConfigurable<double> p6(0.33,  "HFCherenkov:Gain");
00021   //static SimpleConfigurable<bool>   p7(false, "HFCherenkov:CheckSurvive");
00022 
00023   ref_index    = m_HF.getParameter<double>("RefIndex");
00024   lambda1      = ((m_HF.getParameter<double>("Lambda1"))/pow(double(10),7))*cm;
00025   lambda2      = ((m_HF.getParameter<double>("Lambda2"))/pow(double(10),7))*cm;
00026   aperture     = cos(asin(m_HF.getParameter<double>("Aperture")));
00027   apertureTrap = cos(asin(m_HF.getParameter<double>("ApertureTrapped")));
00028   gain         = m_HF.getParameter<double>("Gain");
00029   checkSurvive = m_HF.getParameter<bool>("CheckSurvive");
00030 
00031   edm::LogInfo("HFShower") << "HFCherenkov:: initialised with ref_index " 
00032                            << ref_index << " lambda1/lambda2 (cm) " 
00033                            << lambda1/cm << "/" << lambda2/cm
00034                            << " aperture(total/trapped) " << aperture << "/"
00035                            << apertureTrap << " Check photon survival in HF " 
00036                            << checkSurvive << " Gain " << gain;
00037 
00038   clearVectors();
00039 }
00040 
00041 HFCherenkov::~HFCherenkov() {}
00042         
00043 int HFCherenkov::computeNPhTrapped(double pBeta, 
00044                                    double u, double v, double w, 
00045                                    double step_length, double zFiber, 
00046                                    double dose, int npe_Dose) {
00047    
00048   if (pBeta < (1/ref_index) || step_length < 0.0001) {return 0;}
00049    
00050   double uv = sqrt(u*u + v*v);
00051   int nbOfPhotons = computeNbOfPhotons(pBeta, step_length);
00052 
00053   if (nbOfPhotons < 0) {
00054     return 0;
00055   } else if (nbOfPhotons > 0) {
00056     double u_ph, v_ph, w_ph=0;
00057     for (int i = 0; i < nbOfPhotons; i++) {
00058       double rand     = G4UniformRand();
00059       double theta_C  = acos(1./(pBeta*ref_index));
00060       double phi_C    = 2*M_PI*rand;
00061       double sinTheta = sin(theta_C);
00062       double cosTheta = cos(theta_C);
00063       double sinPhi   = sin(phi_C);
00064       double cosPhi   = cos(phi_C); 
00065       //photon momentum
00066       if (uv < 0.001) { // aligned with z-axis
00067         u_ph = sinTheta * cosPhi;
00068         v_ph = sinTheta * sinPhi;
00069         w_ph = cosTheta;
00070       } else { // general case
00071         u_ph = u * cosTheta  + sinTheta * (v*sinPhi + u*w*cosPhi)/ uv;
00072         v_ph = v * cosTheta  + sinTheta * (-u*sinPhi + v*w*cosPhi)/ uv;
00073         w_ph = w * cosTheta  - sinTheta * cosPhi * uv;
00074       }
00075       if (w_ph > apertureTrap) { // phton trapped inside fiber
00076         npe_Dose += 1; 
00077       }
00078     }
00079   }
00080   int n_photons = npe_Dose;
00081   return n_photons;
00082 }
00083 
00084 
00085 int HFCherenkov::computeNPE(G4ParticleDefinition* pDef, double pBeta, 
00086                             double u, double v, double w, 
00087                             double step_length, double zFiber, 
00088                             double dose, int npe_Dose) {
00089 
00090   clearVectors();
00091   if (!isApplicable(pDef)) {return 0;}
00092   if (pBeta < (1/ref_index) || step_length < 0.0001) {
00093     LogDebug("HFShower") << "HFCherenkov::computeNPE: pBeta " << pBeta 
00094                          << " 1/mu " << (1/ref_index) << " step_length " 
00095                          << step_length;
00096     return 0;
00097   }
00098    
00099   double uv = sqrt(u*u + v*v);
00100   int nbOfPhotons = computeNbOfPhotons(pBeta, step_length);
00101   LogDebug("HFShower") << "HFCherenkov::computeNPE: pBeta " << pBeta 
00102                        << " u/v/w " << u << "/" << v << "/" << w 
00103                        << " step_length " << step_length << " zFib " << zFiber 
00104                        << " nbOfPhotons " << nbOfPhotons;
00105    
00106   if (nbOfPhotons < 0) {
00107     return 0;
00108   } else if (nbOfPhotons > 0) {
00109     double u_ph, v_ph, w_ph=0;
00110     for (int i = 0; i < nbOfPhotons; i++) {
00111       double rand     = G4UniformRand();
00112       double theta_C  = acos(1./(pBeta*ref_index));
00113       double phi_C    = 2*M_PI*rand;
00114       double sinTheta = sin(theta_C);
00115       double cosTheta = cos(theta_C);
00116       double sinPhi   = sin(phi_C);
00117       double cosPhi   = cos(phi_C); 
00118       //photon momentum
00119       if (uv < 0.001) { // aligned with z-axis
00120         u_ph = sinTheta * cosPhi;
00121         v_ph = sinTheta * sinPhi;
00122         w_ph = cosTheta;
00123       } else { // general case
00124         u_ph = u * cosTheta  + sinTheta * (v*sinPhi + u*w*cosPhi)/ uv;
00125         v_ph = v * cosTheta  + sinTheta * (-u*sinPhi + v*w*cosPhi)/ uv;
00126         w_ph = w * cosTheta - sinTheta * cosPhi * uv;
00127       }
00128       double r_lambda = G4UniformRand();
00129       double lambda0 = (lambda1 * lambda2) / (lambda2 - r_lambda *
00130                                               (lambda2 - lambda1));
00131       double lambda  = (lambda0/cm) * pow(double(10),7); // lambda is in nm
00132       wlini.push_back(lambda);
00133       LogDebug("HFShower") << "HFCherenkov::computeNPE: " << i << " lambda " 
00134                            << lambda << " w_ph " << w_ph << " aperture " 
00135                            << aperture;
00136       if (w_ph > aperture) { // phton trapped inside fiber
00137         wltrap.push_back(lambda);
00138         double prob_HF  = 1.0; //photon survived in HF
00139         double a0_inv   = 0.1234;  //meter^-1
00140         double prob_MX  = exp( - 0.5 * a0_inv ); //light mixer
00141         if (checkSurvive) {
00142           double a_inv = a0_inv + 0.14 * pow(dose,0.30);
00143           double z_meters = zFiber;
00144           prob_HF  = exp(-z_meters * a_inv ); //photon survived in HF
00145         }
00146         rand = G4UniformRand();
00147         LogDebug("HFShower") << "HFCherenkov::computeNPE: probHF " << prob_HF
00148                              << " prob_MX " << prob_MX << " Random " << rand 
00149                              << " Survive? " << (rand < (prob_HF * prob_MX));
00150         if (rand < (prob_HF * prob_MX)) { // survived and sent to light mixer
00151           wlatten.push_back(lambda);
00152           rand = G4UniformRand();
00153           double effHEM = computeHEMEff(lambda);
00154           LogDebug("HFShower") << "HFCherenkov::computeNPE: w_ph " << w_ph 
00155                                << " effHEM " << effHEM << " Random " << rand 
00156                                << " Survive? " << (w_ph>0.997||(rand<effHEM));
00157           if (w_ph>0.997 || (rand<effHEM)) { // survived HEM
00158             wlhem.push_back(lambda);
00159             double qEffic = computeQEff(lambda);
00160             rand = G4UniformRand();
00161             LogDebug("HFShower") << "HFCherenkov::computeNPE: qEffic " 
00162                                  << qEffic << " Random " << rand 
00163                                  << " Survive? " <<(rand < qEffic);
00164             if (rand < qEffic) { // made photoelectron
00165               npe_Dose += 1;
00166               momZ.push_back(w_ph);
00167               wl.push_back(lambda);
00168               wlqeff.push_back(lambda);
00169             } // made pe
00170           } // passed HEM
00171         } // passed fiber
00172       } // end of  if(w_ph < w_aperture), trapped inside fiber
00173     }// end of ++NbOfPhotons
00174   } // end of if(NbOfPhotons)}
00175   int npe =  npe_Dose; // Nb of photoelectrons
00176   LogDebug("HFShower") << "HFCherenkov::computeNPE: npe " << npe;
00177   return npe;
00178 }
00179 
00180 std::vector<double>  HFCherenkov::getWLIni() {
00181   std::vector<double> v = wlini;
00182   return v;
00183 }
00184 
00185 std::vector<double>  HFCherenkov::getWLTrap() {
00186   std::vector<double> v = wltrap;
00187   return v;
00188 }
00189 
00190 std::vector<double>  HFCherenkov::getWLAtten() {
00191   std::vector<double> v = wlatten;
00192   return v;
00193 }
00194 
00195 std::vector<double>  HFCherenkov::getWLHEM() {
00196   std::vector<double> v  = wlhem;
00197   return v;
00198 }
00199 
00200 std::vector<double>  HFCherenkov::getWLQEff() {
00201   std::vector<double> v = wlqeff;
00202   return v;
00203 }
00204 
00205 std::vector<double>  HFCherenkov::getWL() {
00206   std::vector<double> v = wl;
00207   return v;
00208 }
00209 
00210 std::vector<double>  HFCherenkov::getMom() {
00211   std::vector<double> v = momZ;
00212   return v;
00213 }
00214 
00215 int HFCherenkov::computeNbOfPhotons(double beta, G4double stepL) {
00216 
00217   double pBeta = beta;
00218   double alpha = 0.0073;
00219   double step_length = stepL;
00220   double theta_C = acos(1./(pBeta*ref_index));
00221   double lambdaDiff = (1./lambda1 - 1./lambda2);
00222   double cherenPhPerLength = 2 * M_PI * alpha * lambdaDiff*cm;
00223   double d_NOfPhotons = cherenPhPerLength * sin(theta_C)*sin(theta_C) *  (step_length/cm);
00224   int nbOfPhotons = int(d_NOfPhotons);
00225   LogDebug("HFShower") << "HFCherenkov::computeNbOfPhotons: StepLength " 
00226                        << step_length << " theta_C " << theta_C 
00227                        << " lambdaDiff " << lambdaDiff
00228                        << " cherenPhPerLength " << cherenPhPerLength 
00229                        << " Photons " << d_NOfPhotons << " " << nbOfPhotons;
00230   return nbOfPhotons;
00231 }
00232 
00233 double HFCherenkov::computeQEff(double wavelength) {
00234 
00235   double y        = (wavelength - 275.) /180.;
00236   double func     = 1. / (1. + 250.*pow((y/5.),4));
00237   double qE_R7525 = 0.77 * y * exp(-y) * func ;
00238   double qeff     = qE_R7525;
00239   LogDebug("HFShower") << "HFCherenkov::computeQEff: wavelength " << wavelength
00240                        << " y/func " << y << "/" << func << " qeff " << qeff;
00241   return qeff;
00242 }
00243 
00244 double HFCherenkov::computeHEMEff(double wavelength) {
00245 
00246   double hEMEff = 0;
00247   if (wavelength < 400.) {
00248     hEMEff = 0.;
00249   } else if (wavelength >= 400. && wavelength < 410.) {
00250     //hEMEff = .99 * (wavelength - 400.) / 10.;
00251     hEMEff = (-1322.453 / wavelength) + 4.2056;
00252   } else if (wavelength >= 410.) {
00253     hEMEff = 0.99;
00254     if (wavelength > 415. && wavelength < 445.) {
00255       //abs(wavelength - 430.) < 15.
00256       //hEMEff = 0.95;
00257       hEMEff = 0.97;
00258     } else if (wavelength > 550. && wavelength < 600.) { 
00259       // abs(wavelength - 575.) < 25.) 
00260       //hEMEff = 0.96;
00261       hEMEff = 0.98;
00262     } else if (wavelength > 565. && wavelength <= 635.) { // added later 
00263       // abs(wavelength - 600.) < 35.) 
00264       hEMEff = (701.7268 / wavelength) - 0.186;
00265     }
00266   }
00267   LogDebug("HFShower") << "HFCherenkov::computeHEMEff: wavelength " 
00268                        << wavelength << " hEMEff " << hEMEff;
00269   return hEMEff;
00270 }
00271 
00272 double HFCherenkov::smearNPE(int npe) {
00273 
00274   double pe = 0.;
00275   if (npe > 0) {
00276     for (int i = 0; i < npe; ++i) {
00277       double val =  G4Poisson(gain);
00278       pe += (val/gain) + 0.001*G4UniformRand();
00279     }
00280   }
00281   LogDebug("HFShower") << "HFCherenkov::smearNPE: npe " << npe << " pe " << pe;
00282   return pe;
00283 }
00284 
00285 void HFCherenkov::clearVectors() {
00286 
00287   wl.clear();
00288   wlini.clear();
00289   wltrap.clear();
00290   wlatten.clear();
00291   wlhem.clear();
00292   wlqeff.clear();
00293   momZ.clear();
00294 }
00295 
00296 bool HFCherenkov::isApplicable(const G4ParticleDefinition* aParticleType) {
00297   bool tmp = (aParticleType->GetPDGCharge() != 0);
00298   LogDebug("HFShower") << "HFCherenkov::isApplicable: aParticleType " 
00299                        << aParticleType->GetParticleName() << " PDGCharge " 
00300                        << aParticleType->GetPDGCharge() << " Result " << tmp;
00301   return tmp;
00302 }
00303    
00304  

Generated on Tue Jun 9 17:46:49 2009 for CMSSW by  doxygen 1.5.4