CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoEcal/EgammaClusterAlgos/src/EgammaSCEnergyCorrectionAlgo.cc

Go to the documentation of this file.
00001 //
00002 // $Id: EgammaSCEnergyCorrectionAlgo.cc,v 1.45 2011/02/17 22:47:03 argiro Exp $
00003 // Author: David Evans, Bristol
00004 //
00005 #include "RecoEcal/EgammaClusterAlgos/interface/EgammaSCEnergyCorrectionAlgo.h"
00006 #include "RecoEcal/EgammaCoreTools/interface/SuperClusterShapeAlgo.h"
00007 #include <iostream>
00008 #include <string>
00009 #include <vector>
00010 
00011 EgammaSCEnergyCorrectionAlgo::EgammaSCEnergyCorrectionAlgo(float noise, 
00012                                                            reco::CaloCluster::AlgoId theAlgo,
00013                                                            const edm::ParameterSet& pset,
00014                                                            EgammaSCEnergyCorrectionAlgo::VerbosityLevel verbosity
00015                                                            )
00016 
00017 {
00018   sigmaElectronicNoise_ = noise;
00019   verbosity_ = verbosity;  
00020 }
00021 
00022 
00023 reco::SuperCluster EgammaSCEnergyCorrectionAlgo::applyCorrection(const reco::SuperCluster &cl, 
00024                                                                  const EcalRecHitCollection &rhc, reco::CaloCluster::AlgoId theAlgo, 
00025                                                                  const CaloSubdetectorGeometry* geometry,
00026                                                                  EcalClusterFunctionBaseClass* energyCorrectionFunction) {      
00027 
00028         
00029   // A little bit of trivial info to be sure all is well
00030 
00031   if (verbosity_ <= pINFO)
00032   {
00033     std::cout << "   EgammaSCEnergyCorrectionAlgo::applyCorrection" << std::endl;
00034     std::cout << "   SC has energy " << cl.energy() << std::endl;
00035     std::cout << "   Will correct now.... " << std::endl;
00036   }
00037 
00038   // Get the seed cluster       
00039   reco::CaloClusterPtr seedC = cl.seed();
00040 
00041   if (verbosity_ <= pINFO)
00042   {
00043     std::cout << "   Seed cluster energy... " << seedC->energy() << std::endl;
00044   }
00045 
00046 
00047   // Find the algorithm used to construct the basic clusters making up the supercluster 
00048   if (verbosity_ <= pINFO) 
00049   {
00050     std::cout << "   The seed cluster used algo " << theAlgo;  
00051   }
00052  
00053   // Find the detector region of the supercluster
00054   // where is the seed cluster?
00055   std::vector<std::pair<DetId, float> > seedHits = seedC->hitsAndFractions();  
00056   EcalSubdetector theBase = EcalSubdetector(seedHits.at(0).first.subdetId());
00057 
00058   if (verbosity_ <= pINFO)
00059   {
00060     std::cout << "   seed cluster location == " << theBase << std::endl;
00061   }
00062 
00063   // Get number of crystals 2sigma above noise in seed basiccluster      
00064   int nCryGT2Sigma = nCrystalsGT2Sigma(*seedC,rhc);
00065   if (verbosity_ <= pINFO)
00066   {
00067     std::cout << "   nCryGT2Sigma " << nCryGT2Sigma << std::endl;
00068   }
00069 
00070   // Supercluster enery - seed basiccluster energy
00071   float bremsEnergy = cl.energy() - seedC->energy();
00072   if (verbosity_ <= pINFO)
00073   {
00074     std::cout << "   bremsEnergy " << bremsEnergy << std::endl;
00075   }
00076 
00077   //Create the pointer ot class SuperClusterShapeAlgo
00078   //which calculates phiWidth and etaWidth
00079   SuperClusterShapeAlgo  SCShape(&rhc, geometry);
00080 
00081   double phiWidth = 0.;
00082   double etaWidth = 0.;
00083   //Calculate phiWidth & etaWidth for SuperClusters
00084   SCShape.Calculate_Covariances(cl);
00085   phiWidth = SCShape.phiWidth();
00086   etaWidth = SCShape.etaWidth();
00087 
00088   // Calculate the new supercluster energy 
00089   //as a function of number of crystals in the seed basiccluster for Endcap 
00090   //or apply new Enegry SCale correction
00091   float newEnergy = 0;
00092 
00093   reco::SuperCluster tmp = cl;
00094   tmp.setPhiWidth(phiWidth); 
00095   tmp.setEtaWidth(etaWidth); 
00096     
00097   if ( theAlgo == reco::CaloCluster::hybrid || theAlgo == reco::CaloCluster::dynamicHybrid ) {
00098     newEnergy = tmp.rawEnergy() + energyCorrectionFunction->getValue(tmp, 3);
00099 
00100   } else if  ( theAlgo == reco::CaloCluster::multi5x5 ) {     
00101     newEnergy = tmp.rawEnergy() + tmp.preshowerEnergy() + energyCorrectionFunction->getValue(tmp, 5);
00102 
00103   } else {  
00104     //Apply f(nCry) correction on island algo and fixedMatrix algo 
00105     newEnergy = seedC->energy()/fNCrystals(nCryGT2Sigma, theAlgo, theBase)+bremsEnergy;
00106   } 
00107 
00108   // Create a new supercluster with the corrected energy 
00109   if (verbosity_ <= pINFO)
00110     {
00111       std::cout << "   UNCORRECTED SC has energy... " << cl.energy() << std::endl;
00112       std::cout << "   CORRECTED SC has energy... " << newEnergy << std::endl;
00113     }
00114 
00115   reco::SuperCluster corrCl =cl;
00116   
00117   corrCl.setEnergy(newEnergy);
00118   corrCl.setPhiWidth(phiWidth);
00119   corrCl.setEtaWidth(etaWidth);
00120   
00121   return corrCl;
00122 }
00123 
00124 float EgammaSCEnergyCorrectionAlgo::fNCrystals(int nCry, reco::CaloCluster::AlgoId theAlgo, EcalSubdetector theBase) const {
00125 
00126   float p0 = 0, p1 = 0, p2 = 0, p3 = 0, p4 = 0;
00127   float x  =  nCry;
00128   float result =1.f;
00129  
00130   if((theBase == EcalBarrel) && (theAlgo == reco::CaloCluster::hybrid))  {
00131     if (nCry<=10) 
00132       {
00133         p0 =  6.32879e-01f; 
00134         p1 =  1.14893e-01f; 
00135         p2 = -2.45705e-02f; 
00136         p3 =  2.53074e-03f; 
00137         p4 = -9.29654e-05f; 
00138       } 
00139     else if (nCry>10 && nCry<=30) 
00140       {
00141         p0 =  6.93196e-01f; 
00142         p1 =  4.44034e-02f; 
00143         p2 = -2.82229e-03f; 
00144         p3 =  8.19495e-05f; 
00145         p4 = -8.96645e-07f; 
00146       } 
00147     else 
00148       {
00149         p0 =  5.65474e+00f; 
00150         p1 = -6.31640e-01f; 
00151         p2 =  3.14218e-02f; 
00152         p3 = -6.84256e-04f; 
00153         p4 =  5.50659e-06f; 
00154       }
00155     if (x > 40.f) x = 40.f;
00156   }
00157   
00158   else if((theBase == EcalEndcap) && (theAlgo == reco::CaloCluster::hybrid)) {
00159     if (verbosity_ <= pERROR)
00160       {
00161         std::cout << "ERROR! HybridEFRYsc called" << std::endl;
00162       } 
00163     return 1.f;
00164   }
00165   
00166   else if((theBase == EcalBarrel) && (theAlgo == reco::CaloCluster::island)) { 
00167     p0 = 4.69976e-01f;     // extracted from fit to all endcap classes with Ebrem = 0.
00168     p1 = 1.45900e-01f;
00169     p2 = -1.61359e-02f;
00170     p3 = 7.99423e-04f;
00171     p4 = -1.47873e-05f;
00172     if (x > 16.f) x = 16.f;
00173   }
00174   
00175   else if((theBase == EcalEndcap) && (theAlgo == reco::CaloCluster::island)) {    
00176     p0 = 4.69976e-01f;     // extracted from fit to all endcap classes with Ebrem = 0.
00177     p1 = 1.45900e-01f;
00178     p2 = -1.61359e-02f;
00179     p3 = 7.99423e-04f;
00180     p4 = -1.47873e-05f;
00181     if (x > 16.f) x = 16.f;
00182   }
00183   
00184   else {
00185     if (verbosity_ <= pINFO)
00186       {
00187         std::cout << "trying to correct unknown cluster!!!" << std::endl;
00188       }
00189     return 1.f;
00190   }
00191   result = p0 + x*(p1 + x*(p2 + x*(p3 + x*p4)));
00192   
00193   //Rescale energy scale correction to take into account change in calibrated
00194   //RecHit definition introduced in CMSSW_1_5_0
00195   float const ebfact = 1.f/0.965f; 
00196   float const eefact = 1.f/0.975f; 
00197     
00198   if(theBase == EcalBarrel) {
00199     result*=ebfact;
00200   } else {
00201     result*=eefact;
00202   }
00203   
00204   return result;  
00205 }
00206 
00207 int EgammaSCEnergyCorrectionAlgo::nCrystalsGT2Sigma(reco::BasicCluster const & seed, EcalRecHitCollection const & rhc) const
00208 {
00209   // return number of crystals 2Sigma above
00210   // electronic noise
00211   
00212   std::vector<std::pair<DetId,float > > const & hits = seed.hitsAndFractions();
00213 
00214   if (verbosity_ <= pINFO)
00215   {
00216     std::cout << "      EgammaSCEnergyCorrectionAlgo::nCrystalsGT2Sigma" << std::endl;
00217     std::cout << "      Will calculate number of crystals above 2sigma noise" << std::endl;
00218     std::cout << "      sigmaElectronicNoise = " << sigmaElectronicNoise_ << std::endl;
00219     std::cout << "      There are " << hits.size() << " recHits" << std::endl;
00220   }
00221 
00222   int nCry = 0;
00223   std::vector<std::pair<DetId,float > >::const_iterator hit;
00224   EcalRecHitCollection::const_iterator aHit;
00225   for(hit = hits.begin(); hit != hits.end(); hit++) {
00226       // need to get hit by DetID in order to get energy
00227       aHit = rhc.find((*hit).first);
00228       if((*aHit).energy()>2.f*sigmaElectronicNoise_) nCry++;
00229     }
00230 
00231   if (verbosity_ <= pINFO)
00232   {
00233     std::cout << "         " << nCry << " of these above 2sigma noise" << std::endl;  
00234   }
00235  
00236   return nCry;
00237 }
00238 
00239 // apply crack correction
00240 
00241 reco::SuperCluster 
00242 EgammaSCEnergyCorrectionAlgo::applyCrackCorrection(const reco::SuperCluster &cl,
00243                                                    EcalClusterFunctionBaseClass* crackCorrectionFunction){
00244 
00245 
00246   double crackcor = 1.; 
00247 
00248   for(reco::CaloCluster_iterator cIt = cl.clustersBegin(); cIt != cl.clustersEnd(); ++cIt) {
00249 
00250     const reco::CaloClusterPtr cc = *cIt; 
00251     crackcor *= ( (cl.rawEnergy() +
00252                    cc->energy()*(crackCorrectionFunction->getValue(*cc)-1.)) / 
00253                    cl.rawEnergy() );   
00254   }// loop on BCs
00255   
00256 
00257   reco::SuperCluster corrCl=cl;
00258   corrCl.setEnergy(cl.energy()*crackcor);
00259   
00260 
00261   return corrCl;
00262 }
00263