CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/ElectroWeakAnalysis/ZEE/interface/SCEnergyCorrections.h

Go to the documentation of this file.
00001 #include "DataFormats/EgammaReco/interface/BasicClusterFwd.h"
00002 
00003 reco::CaloClusterPtrVector CaloClusterVectorCopier(const reco::SuperCluster& sc);
00004 
00005 reco::SuperCluster fEtaScCorr(const reco::SuperCluster& sc)
00006 {
00007         reco::CaloClusterPtrVector bcs = CaloClusterVectorCopier(sc);           
00008         double ieta = fabs(sc.eta())*(5/0.087);
00009         double p0 = 40.2198;
00010         double p1 = -3.03103e-6;
00011         double newE;
00012 //                      std::cout << "Corrected E = Raw E * (1+ p1*(ieta - p0)*(ieta - p0))"<< std::endl;
00013         if ( ieta < p0 ) newE = sc.rawEnergy();
00014         else newE = sc.rawEnergy()/(1 + p1*(ieta - p0)*(ieta - p0));
00015 
00016         reco::SuperCluster corrSc(newE, sc.position(), sc.seed(), bcs, sc.preshowerEnergy(), 0., 0.);
00017         return corrSc;
00018 }
00019 
00020 reco::SuperCluster fBremScCorr(const reco::SuperCluster& sc, const edm::ParameterSet& ps)
00021 {
00022         std::vector<double> fBrem = ps.getParameter<std::vector<double> >("fBremVec");
00023         double bremFrLowThr = ps.getParameter<double>("brLinearLowThr");
00024         double bremFrHighThr = ps.getParameter<double>("brLinearHighThr");
00025         
00026         reco::CaloClusterPtrVector bcs = CaloClusterVectorCopier(sc);           
00027         double bremFrac = sc.phiWidth()/sc.etaWidth();  
00028         double newE = sc.energy();
00029         if(fabs(sc.eta()) < 1.479)
00030         {
00031                 reco::SuperCluster fEtaSC = fEtaScCorr(sc);
00032                 reco::CaloClusterPtrVector bcs = CaloClusterVectorCopier(sc);           
00033                 newE = fEtaSC.energy();
00034         }
00035 
00036         if(bremFrac < bremFrLowThr)  bremFrac = bremFrLowThr;
00037         if(bremFrac < bremFrHighThr)  bremFrac = bremFrHighThr;
00038         
00039         double p0 = fBrem[0]; 
00040         double p1 = fBrem[1]; 
00041         double p2 = fBrem[2]; 
00042         double p3 = fBrem[3]; 
00043         double p4 = fBrem[4]; 
00044         //
00045         double threshold = p4; 
00046           
00047         double y = p0*threshold*threshold + p1*threshold + p2;
00048         double yprime = 2*p0*threshold + p1;
00049         double a = p3;
00050         double b = yprime - 2*a*threshold;
00051         double c = y - a*threshold*threshold - b*threshold;
00052          
00053         double fCorr = 1;
00054         if( bremFrac < threshold ) 
00055             fCorr = p0*bremFrac*bremFrac + p1*bremFrac + p2;
00056         else 
00057             fCorr = a*bremFrac*bremFrac + b*bremFrac + c;
00058 
00059         newE /= fCorr;
00060         reco::SuperCluster corrSc(newE, sc.position(), sc.seed(), bcs, sc.preshowerEnergy(), 0., 0.);
00061         return corrSc;
00062 }
00063 
00064 reco::SuperCluster fEtEtaCorr(const reco::SuperCluster& sc, const edm::ParameterSet& ps)
00065 {
00066   // et -- Et of the SuperCluster (with respect to (0,0,0))
00067   // eta -- eta of the SuperCluster
00068         std::vector<double> fEtEtaParams = ps.getParameter<std::vector<double> >("fEtEtaParamsVec");
00069 
00070         reco::SuperCluster fBremSC = fBremScCorr(sc, ps);
00071         reco::CaloClusterPtrVector bcs = CaloClusterVectorCopier(sc);           
00072 
00073         double eta = sc.eta();
00074         double et = fBremSC.energy()/cosh(eta);
00075         double fCorr = 0.;
00076 
00077         double p0 = fEtEtaParams[0] + fEtEtaParams[1]/(et + fEtEtaParams[ 2]) + fEtEtaParams[ 3]/(et*et);
00078         double p1 = fEtEtaParams[4] + fEtEtaParams[5]/(et + fEtEtaParams[ 6]) + fEtEtaParams[ 7]/(et*et);
00079         double p2 = fEtEtaParams[8] + fEtEtaParams[9]/(et + fEtEtaParams[10]) + fEtEtaParams[11]/(et*et);
00080 
00081         fCorr = 
00082         p0 + 
00083         p1 * atan(fEtEtaParams[12]*(fEtEtaParams[13]-fabs(eta))) + fEtEtaParams[14] * fabs(eta) + 
00084         p1 * fEtEtaParams[15] * fabs(eta) +
00085         p2 * fEtEtaParams[16] * eta * eta; 
00086 
00087         if ( fCorr < 0.5 ) fCorr = 0.5;
00088 
00089         double newE = et/(fCorr*cosh(eta));
00090         reco::SuperCluster corrSc(newE, sc.position(), sc.seed(), bcs, sc.preshowerEnergy(), 0., 0.);
00091         return corrSc;
00092 }
00093 
00094 reco::SuperCluster fEAddScCorr(const reco::SuperCluster& sc, double Ecorr)
00095 {
00096         reco::CaloClusterPtrVector bcs = CaloClusterVectorCopier(sc);           
00097         
00098         double newE = sc.rawEnergy()+Ecorr; 
00099         reco::SuperCluster corrSc(newE, sc.position(), sc.seed(), bcs, sc.preshowerEnergy(), 0., 0.);
00100         return corrSc;
00101 }
00102         
00103 
00104 reco::CaloClusterPtrVector CaloClusterVectorCopier(const reco::SuperCluster& sc)
00105 {
00106         reco::CaloClusterPtrVector clusters_v;
00107 
00108         for(reco::CaloCluster_iterator cluster = sc.clustersBegin(); cluster != sc.clustersEnd(); cluster ++)
00109         {
00110                 clusters_v.push_back(*cluster);
00111         }
00112 
00113         return clusters_v;
00114 }