CMS 3D CMS Logo

ElectronLikelihood.cc

Go to the documentation of this file.
00001 #include "FWCore/Utilities/interface/Exception.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
00004 #include "RecoEgamma/ElectronIdentification/interface/ElectronLikelihood.h"
00005 #include <iostream>
00006 
00007 
00008 ElectronLikelihood::ElectronLikelihood (const ElectronLikelihoodCalibration *calibration,
00009                                         const std::vector<double> & fisherEBLt15,
00010                                         const std::vector<double> & fisherEBGt15,
00011                                         const std::vector<double> & fisherEELt15,
00012                                         const std::vector<double> & fisherEEGt15,
00013                                         const std::vector<double> & eleFracsEBlt15,
00014                                         const std::vector<double> & piFracsEBlt15,
00015                                         const std::vector<double> & eleFracsEElt15,
00016                                         const std::vector<double> & piFracsEElt15,
00017                                         const std::vector<double> & eleFracsEBgt15,
00018                                         const std::vector<double> & piFracsEBgt15,
00019                                         const std::vector<double> & eleFracsEEgt15,
00020                                         const std::vector<double> & piFracsEEgt15,
00021                                         double eleWeight,
00022                                         double piWeight,
00023                                         LikelihoodSwitches eleIDSwitches,
00024                                         std::string signalWeightSplitting,
00025                                         std::string backgroundWeightSplitting,
00026                                         bool splitSignalPdfs,
00027                                         bool splitBackgroundPdfs) :
00028   _EBlt15lh (new LikelihoodPdfProduct ("electronID_EB_ptLt15_likelihood",0,0)) ,
00029   _EElt15lh (new LikelihoodPdfProduct ("electronID_EE_ptLt15_likelihood",1,0)) ,
00030   _EBgt15lh (new LikelihoodPdfProduct ("electronID_EB_ptGt15_likelihood",0,1)) ,
00031   _EEgt15lh (new LikelihoodPdfProduct ("electronID_EE_ptGt15_likelihood",1,1)) ,
00032   m_eleIDSwitches (eleIDSwitches) ,
00033   m_signalWeightSplitting (signalWeightSplitting), 
00034   m_backgroundWeightSplitting (backgroundWeightSplitting),
00035   m_splitSignalPdfs (splitSignalPdfs), 
00036   m_splitBackgroundPdfs (splitBackgroundPdfs)  
00037 {
00038   Setup (calibration,
00039          fisherEBLt15, fisherEBGt15,
00040          fisherEELt15, fisherEEGt15,
00041          eleFracsEBlt15, piFracsEBlt15,
00042          eleFracsEElt15, piFracsEElt15,
00043          eleFracsEBgt15, piFracsEBgt15,
00044          eleFracsEEgt15, piFracsEEgt15,
00045          eleWeight, piWeight,
00046          signalWeightSplitting, backgroundWeightSplitting,
00047          splitSignalPdfs, splitBackgroundPdfs) ;
00048 }
00049 
00050 
00051 
00052 // --------------------------------------------------------
00053 
00054 
00055 
00056 ElectronLikelihood::~ElectronLikelihood () {
00057   delete _EBlt15lh ;
00058   delete _EElt15lh ;
00059   delete _EBgt15lh ;
00060   delete _EEgt15lh ;
00061 }
00062 
00063 
00064 
00065 // --------------------------------------------------------
00066 
00067 
00068 void 
00069 ElectronLikelihood::Setup (const ElectronLikelihoodCalibration *calibration,
00070                            const std::vector<double> & fisherEBLt15,
00071                            const std::vector<double> & fisherEBGt15,
00072                            const std::vector<double> & fisherEELt15,
00073                            const std::vector<double> & fisherEEGt15,
00074                            const std::vector<double> & eleFracsEBlt15,
00075                            const std::vector<double> & piFracsEBlt15,
00076                            const std::vector<double> & eleFracsEElt15,
00077                            const std::vector<double> & piFracsEElt15,
00078                            const std::vector<double> & eleFracsEBgt15,
00079                            const std::vector<double> & piFracsEBgt15,
00080                            const std::vector<double> & eleFracsEEgt15,
00081                            const std::vector<double> & piFracsEEgt15,
00082                            double eleWeight,
00083                            double piWeight,
00084                            std::string signalWeightSplitting,
00085                            std::string backgroundWeightSplitting,
00086                            bool splitSignalPdfs,
00087                            bool splitBackgroundPdfs) 
00088 {
00089 
00090   // ECAL BARREL LIKELIHOOD - Pt < 15 GeV region
00091   _EBlt15lh->initFromDB (calibration) ;
00092 
00093   _EBlt15lh->addSpecies ("electrons",eleWeight) ;
00094   _EBlt15lh->addSpecies ("hadrons",  piWeight) ;
00095 
00096   if(signalWeightSplitting.compare("fullclass")==0) {
00097     _EBlt15lh->setSplitFrac ("electrons", "fullclass0", eleFracsEBlt15[0]) ;
00098     _EBlt15lh->setSplitFrac ("electrons", "fullclass1", eleFracsEBlt15[1]) ;
00099     _EBlt15lh->setSplitFrac ("electrons", "fullclass2", eleFracsEBlt15[2]) ;
00100     _EBlt15lh->setSplitFrac ("electrons", "fullclass3", eleFracsEBlt15[3]) ;
00101   }
00102   else if(signalWeightSplitting.compare("class")==0) {
00103     _EBlt15lh->setSplitFrac ("electrons", "class0", eleFracsEBlt15[0]+eleFracsEBlt15[1]+eleFracsEBlt15[2]) ;
00104     _EBlt15lh->setSplitFrac ("electrons", "class1", eleFracsEBlt15[3]) ;
00105   }
00106   else {
00107     throw cms::Exception("BadConfig") << "Only class (non-showering / showering)"
00108                                       << " and fullclass (golden / bigbrem / narrow / showering)" 
00109                                       << " splitting is implemented right now";
00110   }
00111 
00112   if (m_eleIDSwitches.m_useDeltaPhiIn)   _EBlt15lh->addPdf ("electrons", "dPhiVtx",     splitSignalPdfs) ;   
00113   if (m_eleIDSwitches.m_useDeltaEtaCalo) _EBlt15lh->addPdf ("electrons", "dEtaCalo",    splitSignalPdfs) ;
00114   if (m_eleIDSwitches.m_useEoverPOut)    _EBlt15lh->addPdf ("electrons", "EoPout",      splitSignalPdfs) ;
00115   if (m_eleIDSwitches.m_useHoverE)       _EBlt15lh->addPdf ("electrons", "HoE",         splitSignalPdfs) ;
00116   if (m_eleIDSwitches.m_useShapeFisher)  _EBlt15lh->addPdf ("electrons", "shapeFisher", splitSignalPdfs) ;
00117   if (m_eleIDSwitches.m_useSigmaEtaEta)  _EBlt15lh->addPdf ("electrons", "sigmaEtaEta", splitSignalPdfs) ;
00118   if (m_eleIDSwitches.m_useE9overE25)    _EBlt15lh->addPdf ("electrons", "s9s25",       splitSignalPdfs) ;
00119 
00120   if(backgroundWeightSplitting.compare("fullclass")==0) {
00121     _EBlt15lh->setSplitFrac ("hadrons", "fullclass0", piFracsEBlt15[0]) ;
00122     _EBlt15lh->setSplitFrac ("hadrons", "fullclass1", piFracsEBlt15[1]) ;
00123     _EBlt15lh->setSplitFrac ("hadrons", "fullclass2", piFracsEBlt15[2]) ;
00124     _EBlt15lh->setSplitFrac ("hadrons", "fullclass3", piFracsEBlt15[3]) ;
00125   }
00126   else if(backgroundWeightSplitting.compare("class")==0) {
00127     _EBlt15lh->setSplitFrac ("hadrons", "class0", piFracsEBlt15[0]+piFracsEBlt15[1]+piFracsEBlt15[2]) ;
00128     _EBlt15lh->setSplitFrac ("hadrons", "class1", piFracsEBlt15[3]) ;
00129   }
00130   else {
00131     throw cms::Exception("BadConfig") << "Only class (non-showering / showering)"
00132                                       << " and fullclass (golden / bigbrem / narrow / showering)" 
00133                                       << " splitting is implemented right now";
00134   }
00135 
00136   if (m_eleIDSwitches.m_useDeltaPhiIn)   _EBlt15lh->addPdf ("hadrons", "dPhiVtx",       splitBackgroundPdfs) ;  
00137   if (m_eleIDSwitches.m_useDeltaEtaCalo) _EBlt15lh->addPdf ("hadrons", "dEtaCalo",      splitBackgroundPdfs) ;
00138   if (m_eleIDSwitches.m_useEoverPOut)    _EBlt15lh->addPdf ("hadrons", "EoPout",        splitBackgroundPdfs) ;
00139   if (m_eleIDSwitches.m_useHoverE)       _EBlt15lh->addPdf ("hadrons", "HoE",           splitBackgroundPdfs) ;
00140   if (m_eleIDSwitches.m_useShapeFisher)  _EBlt15lh->addPdf ("hadrons", "shapeFisher",   splitBackgroundPdfs) ;
00141   if (m_eleIDSwitches.m_useSigmaEtaEta)  _EBlt15lh->addPdf ("hadrons", "sigmaEtaEta",   splitBackgroundPdfs) ;
00142   if (m_eleIDSwitches.m_useE9overE25)    _EBlt15lh->addPdf ("hadrons", "s9s25",         splitBackgroundPdfs) ;
00143 
00144 
00145   // ECAL BARREL LIKELIHOOD - Pt >= 15 GeV region
00146   _EBgt15lh->initFromDB (calibration) ;
00147 
00148   _EBgt15lh->addSpecies ("electrons",eleWeight) ;  
00149   _EBgt15lh->addSpecies ("hadrons",piWeight) ;
00150 
00151   if(signalWeightSplitting.compare("fullclass")==0) {
00152     _EBgt15lh->setSplitFrac ("electrons", "fullclass0", eleFracsEBgt15[0]) ;
00153     _EBgt15lh->setSplitFrac ("electrons", "fullclass1", eleFracsEBgt15[1]) ;
00154     _EBgt15lh->setSplitFrac ("electrons", "fullclass2", eleFracsEBgt15[2]) ;
00155     _EBgt15lh->setSplitFrac ("electrons", "fullclass3", eleFracsEBgt15[3]) ;
00156   }
00157   else if(signalWeightSplitting.compare("class")==0) {
00158     _EBgt15lh->setSplitFrac ("electrons", "class0", eleFracsEBgt15[0]+eleFracsEBgt15[1]+eleFracsEBgt15[2]) ;
00159     _EBgt15lh->setSplitFrac ("electrons", "class1", eleFracsEBgt15[3]) ;
00160   }
00161   else {
00162     throw cms::Exception("BadConfig") << "Only class (non-showering / showering)"
00163                                       << " and fullclass (golden / bigbrem / narrow / showering)" 
00164                                       << " splitting is implemented right now";
00165   }
00166 
00167   if (m_eleIDSwitches.m_useDeltaPhiIn)   _EBgt15lh->addPdf ("electrons", "dPhiVtx",     splitSignalPdfs) ;  
00168   if (m_eleIDSwitches.m_useDeltaEtaCalo) _EBgt15lh->addPdf ("electrons", "dEtaCalo",    splitSignalPdfs) ;
00169   if (m_eleIDSwitches.m_useEoverPOut)    _EBgt15lh->addPdf ("electrons", "EoPout",      splitSignalPdfs) ;
00170   if (m_eleIDSwitches.m_useHoverE)       _EBgt15lh->addPdf ("electrons", "HoE",         splitSignalPdfs) ;
00171   if (m_eleIDSwitches.m_useShapeFisher)  _EBgt15lh->addPdf ("electrons", "shapeFisher", splitSignalPdfs) ;
00172   if (m_eleIDSwitches.m_useSigmaEtaEta)  _EBgt15lh->addPdf ("electrons", "sigmaEtaEta", splitSignalPdfs) ;
00173   if (m_eleIDSwitches.m_useE9overE25)    _EBgt15lh->addPdf ("electrons", "s9s25",       splitSignalPdfs) ;
00174 
00175 
00176   if(backgroundWeightSplitting.compare("fullclass")==0) {
00177     _EBgt15lh->setSplitFrac ("hadrons", "fullclass0", piFracsEBgt15[0]) ;
00178     _EBgt15lh->setSplitFrac ("hadrons", "fullclass1", piFracsEBgt15[1]) ;
00179     _EBgt15lh->setSplitFrac ("hadrons", "fullclass2", piFracsEBgt15[2]) ;
00180     _EBgt15lh->setSplitFrac ("hadrons", "fullclass3", piFracsEBgt15[3]) ;
00181   }
00182   else if(backgroundWeightSplitting.compare("class")==0) {
00183     _EBgt15lh->setSplitFrac ("hadrons", "class0", piFracsEBgt15[0]+piFracsEBgt15[1]+ piFracsEBgt15[2]) ;
00184     _EBgt15lh->setSplitFrac ("hadrons", "class1", piFracsEBgt15[3]) ;
00185   }
00186   else {
00187     throw cms::Exception("BadConfig") << "Only class (non-showering / showering)"
00188                                       << " and fullclass (golden / bigbrem / narrow / showering)" 
00189                                       << " splitting is implemented right now";
00190   }
00191 
00192   if (m_eleIDSwitches.m_useDeltaPhiIn)   _EBgt15lh->addPdf ("hadrons", "dPhiVtx",       splitBackgroundPdfs) ; 
00193   if (m_eleIDSwitches.m_useDeltaEtaCalo) _EBgt15lh->addPdf ("hadrons", "dEtaCalo",      splitBackgroundPdfs) ;
00194   if (m_eleIDSwitches.m_useEoverPOut)    _EBgt15lh->addPdf ("hadrons", "EoPout",        splitBackgroundPdfs) ;
00195   if (m_eleIDSwitches.m_useHoverE)       _EBgt15lh->addPdf ("hadrons", "HoE",           splitBackgroundPdfs) ;
00196   if (m_eleIDSwitches.m_useShapeFisher)  _EBgt15lh->addPdf ("hadrons", "shapeFisher",   splitBackgroundPdfs) ;
00197   if (m_eleIDSwitches.m_useSigmaEtaEta)  _EBgt15lh->addPdf ("hadrons", "sigmaEtaEta",   splitBackgroundPdfs) ;
00198   if (m_eleIDSwitches.m_useE9overE25)    _EBgt15lh->addPdf ("hadrons", "s9s25",         splitBackgroundPdfs) ;
00199 
00200   // ECAL ENDCAP LIKELIHOOD - Pt < 15 GeV
00201   _EElt15lh->initFromDB (calibration) ;
00202 
00203   _EElt15lh->addSpecies ("electrons",eleWeight) ;
00204   _EElt15lh->addSpecies ("hadrons",piWeight) ;
00205 
00206   if(signalWeightSplitting.compare("fullclass")==0) {
00207     _EElt15lh->setSplitFrac ("electrons", "fullclass0", eleFracsEElt15[0]) ;
00208     _EElt15lh->setSplitFrac ("electrons", "fullclass1", eleFracsEElt15[1]) ;
00209     _EElt15lh->setSplitFrac ("electrons", "fullclass2", eleFracsEElt15[2]) ;
00210     _EElt15lh->setSplitFrac ("electrons", "fullclass3", eleFracsEElt15[3]) ;
00211   }
00212   else if(signalWeightSplitting.compare("class")==0) {
00213     _EElt15lh->setSplitFrac ("electrons", "class0", eleFracsEElt15[0]+eleFracsEElt15[1]+eleFracsEElt15[2]) ;
00214     _EElt15lh->setSplitFrac ("electrons", "class1", eleFracsEElt15[3]) ;
00215   }
00216   else {
00217     throw cms::Exception("BadConfig") << "Only class (non-showering / showering)"
00218                                       << " and fullclass (golden / bigbrem / narrow / showering)" 
00219                                       << " splitting is implemented right now";
00220   }
00221 
00222   if (m_eleIDSwitches.m_useDeltaPhiIn)   _EElt15lh->addPdf ("electrons", "dPhiVtx",     splitSignalPdfs) ;  
00223   if (m_eleIDSwitches.m_useDeltaEtaCalo) _EElt15lh->addPdf ("electrons", "dEtaCalo",    splitSignalPdfs) ;
00224   if (m_eleIDSwitches.m_useEoverPOut)    _EElt15lh->addPdf ("electrons", "EoPout",      splitSignalPdfs) ;
00225   if (m_eleIDSwitches.m_useHoverE)       _EElt15lh->addPdf ("electrons", "HoE",         splitSignalPdfs) ;
00226   if (m_eleIDSwitches.m_useShapeFisher)  _EElt15lh->addPdf ("electrons", "shapeFisher", splitSignalPdfs) ;
00227   if (m_eleIDSwitches.m_useSigmaEtaEta)  _EElt15lh->addPdf ("electrons", "sigmaEtaEta", splitSignalPdfs) ;
00228   if (m_eleIDSwitches.m_useE9overE25)    _EElt15lh->addPdf ("electrons", "s9s25",       splitSignalPdfs) ;
00229 
00230   if(backgroundWeightSplitting.compare("fullclass")==0) {
00231     _EElt15lh->setSplitFrac ("hadrons", "fullclass0", piFracsEElt15[0]) ;
00232     _EElt15lh->setSplitFrac ("hadrons", "fullclass1", piFracsEElt15[1]) ;
00233     _EElt15lh->setSplitFrac ("hadrons", "fullclass2", piFracsEElt15[2]) ;
00234     _EElt15lh->setSplitFrac ("hadrons", "fullclass3", piFracsEElt15[3]) ;
00235   }
00236   else if(backgroundWeightSplitting.compare("class")==0) {
00237     _EElt15lh->setSplitFrac ("hadrons", "class0", piFracsEElt15[0]+piFracsEElt15[1]+piFracsEElt15[2]) ;
00238     _EElt15lh->setSplitFrac ("hadrons", "class1", piFracsEElt15[3]) ;
00239   }
00240   else {
00241     throw cms::Exception("BadConfig") << "Only class (non-showering / showering)"
00242                                       << " and fullclass (golden / bigbrem / narrow / showering)" 
00243                                       << " splitting is implemented right now";
00244   }
00245 
00246 
00247   if (m_eleIDSwitches.m_useDeltaPhiIn)   _EElt15lh->addPdf ("hadrons", "dPhiVtx",     splitBackgroundPdfs) ;  
00248   if (m_eleIDSwitches.m_useDeltaEtaCalo) _EElt15lh->addPdf ("hadrons", "dEtaCalo",    splitBackgroundPdfs) ;
00249   if (m_eleIDSwitches.m_useEoverPOut)    _EElt15lh->addPdf ("hadrons", "EoPout",      splitBackgroundPdfs) ;
00250   if (m_eleIDSwitches.m_useHoverE)       _EElt15lh->addPdf ("hadrons", "HoE",         splitBackgroundPdfs) ;
00251   if (m_eleIDSwitches.m_useShapeFisher)  _EElt15lh->addPdf ("hadrons", "shapeFisher", splitBackgroundPdfs) ;
00252   if (m_eleIDSwitches.m_useSigmaEtaEta)  _EElt15lh->addPdf ("hadrons", "sigmaEtaEta", splitBackgroundPdfs) ;
00253   if (m_eleIDSwitches.m_useE9overE25)    _EElt15lh->addPdf ("hadrons", "s9s25",       splitBackgroundPdfs) ;
00254 
00255 
00256   // ECAL ENDCAP LIKELIHOOD - Pt >= 15 GeV
00257   _EEgt15lh->initFromDB (calibration) ;
00258 
00259   _EEgt15lh->addSpecies ("electrons",eleWeight) ;
00260   _EEgt15lh->addSpecies ("hadrons",piWeight) ;
00261 
00262   if(signalWeightSplitting.compare("fullclass")==0) {
00263     _EEgt15lh->setSplitFrac ("electrons", "fullclass0", eleFracsEEgt15[0]) ;
00264     _EEgt15lh->setSplitFrac ("electrons", "fullclass1", eleFracsEEgt15[1]) ;
00265     _EEgt15lh->setSplitFrac ("electrons", "fullclass2", eleFracsEEgt15[2]) ;
00266     _EEgt15lh->setSplitFrac ("electrons", "fullclass3", eleFracsEEgt15[3]) ;
00267   }
00268   else if(signalWeightSplitting.compare("class")==0) {
00269     _EEgt15lh->setSplitFrac ("electrons", "class0", eleFracsEEgt15[0]+eleFracsEEgt15[1]+eleFracsEEgt15[2]) ;
00270     _EEgt15lh->setSplitFrac ("electrons", "class1", eleFracsEEgt15[3]) ;
00271   }
00272   else {
00273     throw cms::Exception("BadConfig") << "Only class (non-showering / showering)"
00274                                       << " and fullclass (golden / bigbrem / narrow / showering)" 
00275                                       << " splitting is implemented right now";
00276   }
00277 
00278   if (m_eleIDSwitches.m_useDeltaPhiIn)   _EEgt15lh->addPdf ("electrons", "dPhiVtx",     splitSignalPdfs) ;
00279   if (m_eleIDSwitches.m_useDeltaEtaCalo) _EEgt15lh->addPdf ("electrons", "dEtaCalo",    splitSignalPdfs) ;
00280   if (m_eleIDSwitches.m_useEoverPOut)    _EEgt15lh->addPdf ("electrons", "EoPout",      splitSignalPdfs) ;
00281   if (m_eleIDSwitches.m_useHoverE)       _EEgt15lh->addPdf ("electrons", "HoE",         splitSignalPdfs) ;
00282   if (m_eleIDSwitches.m_useShapeFisher)  _EEgt15lh->addPdf ("electrons", "shapeFisher", splitSignalPdfs) ;
00283   if (m_eleIDSwitches.m_useSigmaEtaEta)  _EEgt15lh->addPdf ("electrons", "sigmaEtaEta", splitSignalPdfs) ;
00284   if (m_eleIDSwitches.m_useE9overE25)    _EEgt15lh->addPdf ("electrons", "s9s25",       splitSignalPdfs) ;
00285 
00286   if(backgroundWeightSplitting.compare("fullclass")==0) {
00287     _EEgt15lh->setSplitFrac ("hadrons", "fullclass0", piFracsEEgt15[0]) ;
00288     _EEgt15lh->setSplitFrac ("hadrons", "fullclass1", piFracsEEgt15[1]) ;
00289     _EEgt15lh->setSplitFrac ("hadrons", "fullclass2", piFracsEEgt15[2]) ;
00290     _EEgt15lh->setSplitFrac ("hadrons", "fullclass3", piFracsEEgt15[3]) ;
00291   }
00292   else if(backgroundWeightSplitting.compare("class")==0) {
00293     _EEgt15lh->setSplitFrac ("hadrons", "class0", piFracsEEgt15[0]+piFracsEEgt15[1]+piFracsEEgt15[2]) ;
00294     _EEgt15lh->setSplitFrac ("hadrons", "class1", piFracsEEgt15[3]) ;
00295   }
00296   else {
00297     throw cms::Exception("BadConfig") << "Only class (non-showering / showering)"
00298                                       << " and fullclass (golden / bigbrem / narrow / showering)" 
00299                                       << " splitting is implemented right now";
00300   }
00301 
00302   // initialise the fisher coefficients from cfi
00303   m_fisherEBLt15 = fisherEBLt15;
00304   m_fisherEBGt15 = fisherEBGt15;
00305   m_fisherEELt15 = fisherEELt15;
00306   m_fisherEEGt15 = fisherEEGt15;
00307 
00308   if (m_eleIDSwitches.m_useDeltaPhiIn)   _EEgt15lh->addPdf ("hadrons", "dPhiVtx",     splitBackgroundPdfs) ; 
00309   if (m_eleIDSwitches.m_useDeltaEtaCalo) _EEgt15lh->addPdf ("hadrons", "dEtaCalo",    splitBackgroundPdfs) ;
00310   if (m_eleIDSwitches.m_useEoverPOut)    _EEgt15lh->addPdf ("hadrons", "EoPout",      splitBackgroundPdfs) ;
00311   if (m_eleIDSwitches.m_useHoverE)       _EEgt15lh->addPdf ("hadrons", "HoE",         splitBackgroundPdfs) ;
00312   if (m_eleIDSwitches.m_useShapeFisher)  _EEgt15lh->addPdf ("hadrons", "shapeFisher", splitBackgroundPdfs) ;
00313   if (m_eleIDSwitches.m_useSigmaEtaEta)  _EEgt15lh->addPdf ("hadrons", "sigmaEtaEta", splitBackgroundPdfs) ;
00314   if (m_eleIDSwitches.m_useE9overE25)    _EEgt15lh->addPdf ("hadrons", "s9s25",       splitBackgroundPdfs) ;
00315 }
00316 
00317 
00318 
00319 // --------------------------------------------------------
00320 
00321 
00322 
00323 void 
00324 ElectronLikelihood::getInputVar (const reco::GsfElectron &electron, 
00325                                  std::vector<float> &measurements, 
00326                                  EcalClusterLazyTools myEcalCluster) const 
00327 {
00328 
00329   // the variables entering the likelihood
00330   if (m_eleIDSwitches.m_useDeltaPhiIn) measurements.push_back ( electron.deltaPhiSuperClusterTrackAtVtx () ) ;
00331   if (m_eleIDSwitches.m_useDeltaEtaCalo) measurements.push_back ( electron.deltaEtaSeedClusterTrackAtCalo () ) ;
00332   if (m_eleIDSwitches.m_useEoverPOut) measurements.push_back ( electron.eSeedClusterOverPout () ) ;
00333   if (m_eleIDSwitches.m_useHoverE) measurements.push_back ( electron.hadronicOverEm () ) ;
00334   if (m_eleIDSwitches.m_useShapeFisher) measurements.push_back ( CalculateFisher(electron, myEcalCluster) ) ;
00335   std::vector<float> vCov = myEcalCluster.covariances(*(electron.superCluster()->seed())) ;
00336   if (m_eleIDSwitches.m_useSigmaEtaEta) measurements.push_back ( sqrt (vCov[0]) );
00337   if (m_eleIDSwitches.m_useE9overE25)   measurements.push_back ( myEcalCluster.e3x3(*(electron.superCluster()->seed()))/
00338                                                                  myEcalCluster.e5x5(*(electron.superCluster()->seed()))) ;
00339 }
00340 
00341 
00342 
00343 // --------------------------------------------------------
00344 
00345 
00346 
00347 double 
00348 ElectronLikelihood::CalculateFisher(const reco::GsfElectron &electron,
00349                                     EcalClusterLazyTools myEcalCluster) const
00350 {
00351 
00352   // the variables entering the shape fisher
00353   double s9s25, sigmaEtaEta, etaLat, a20;
00354 
00355   s9s25=myEcalCluster.e3x3(*(electron.superCluster()->seed()))/
00356         myEcalCluster.e5x5(*(electron.superCluster()->seed())) ;
00357   std::vector<float> vCov = myEcalCluster.covariances(*(electron.superCluster()->seed())) ;
00358   sigmaEtaEta=sqrt (vCov[0]) ;
00359   std::vector<float> vLat = myEcalCluster.lat(*(electron.superCluster()->seed())) ;
00360   etaLat=vLat[0] ;
00361   a20=myEcalCluster.zernike20 (*(electron.superCluster()->seed())) ;
00362 
00363 
00364   // evaluate the Fisher discriminant
00365   double clusterShapeFisher;
00366   std::vector<DetId> vecId=electron.superCluster()->getHitsByDetId () ;
00367   EcalSubdetector subdet = EcalSubdetector (vecId[0].subdetId ()) ;
00368   
00369   if (subdet==EcalBarrel) {
00370     if (electron.pt()<15.) {
00371       clusterShapeFisher = m_fisherEBLt15[0] + 
00372         m_fisherEBLt15[1] * sigmaEtaEta +
00373         m_fisherEBLt15[2] * s9s25 +
00374         m_fisherEBLt15[3] * etaLat +
00375         m_fisherEBLt15[4] * a20;
00376     }
00377     else {
00378       clusterShapeFisher = m_fisherEBGt15[0] + 
00379         m_fisherEBGt15[1] * sigmaEtaEta +
00380         m_fisherEBGt15[2] * s9s25 +
00381         m_fisherEBGt15[3] * etaLat +
00382         m_fisherEBGt15[4] * a20;
00383     }
00384   }
00385   else if (subdet==EcalEndcap) {
00386     if (electron.pt()<15.) {
00387       clusterShapeFisher = m_fisherEELt15[0] + 
00388         m_fisherEELt15[1] * sigmaEtaEta +
00389         m_fisherEELt15[2] * s9s25 +
00390         m_fisherEELt15[3] * etaLat +
00391         m_fisherEELt15[4] * a20;
00392     }
00393     else {
00394       clusterShapeFisher = m_fisherEEGt15[0] + 
00395         m_fisherEEGt15[1] * sigmaEtaEta +
00396         m_fisherEEGt15[2] * s9s25 +
00397         m_fisherEEGt15[3] * etaLat +
00398         m_fisherEEGt15[4] * a20;
00399     }
00400   }
00401   else {
00402     clusterShapeFisher = -999 ;
00403     edm::LogWarning ("ElectronLikelihood") << "Undefined electron, eta = " << electron.eta () << "!" ;
00404   }
00405   return clusterShapeFisher;
00406 }
00407 
00408 
00409 
00410 // --------------------------------------------------------
00411 
00412 
00413 
00414 float 
00415 ElectronLikelihood::result (const reco::GsfElectron &electron, 
00416                             EcalClusterLazyTools myEcalCluster) const 
00417 {
00418 
00419   //=======================================================
00420   // used classification:
00421   // golden                     =>  0
00422   // big brem                   => 10
00423   // narrow                     => 20
00424   // showering nbrem 0          => 30
00425   // showering nbrem 1          => 31
00426   // showering nbrem 2          => 32
00427   // showering nbrem 3          => 33
00428   // showering nbrem 4 ou plus  => 34
00429   // cracks                     => 40
00430   // endcap                     => barrel + 100
00431   //=======================================================
00432 
00433   std::vector<float> measurements ;
00434   getInputVar (electron, measurements, myEcalCluster) ;
00435 
00436   // Split using only the 10^1 bit (golden/big brem/narrow/showering)
00437   int bitVal=-1 ;
00438   int gsfclass=electron.classification () ;
00439   if (gsfclass<99) 
00440     bitVal=int (gsfclass)/10 ; 
00441   else
00442     bitVal=int (int (gsfclass)%100)/10 ;
00443   
00444   // temporary: crack electrons goes into Class3 (showering)
00445   if (bitVal==4) bitVal=3 ;
00446   if (bitVal<0 || bitVal>3)
00447     throw cms::Exception ("ElectronLikelihood") << "ERROR! electron class " << gsfclass << " is not accepted!\n" ;
00448 
00449   char className[20] ;
00450   if(m_signalWeightSplitting.compare("fullclass")==0) {
00451     sprintf (className,"fullclass%d",bitVal) ;
00452   }
00453   else if(m_signalWeightSplitting.compare("class")==0) {
00454     int classVal = (bitVal<3) ? 0 : 1;
00455     sprintf (className,"class%d",classVal);
00456   }
00457   else {
00458     throw cms::Exception("BadConfig") << "Only class (non-showering / showering)"
00459                                       << " and fullclass (golden / bigbrem / narrow / showering)" 
00460                                       << " splitting is implemented right now";
00461   }
00462 
00463   reco::SuperClusterRef sclusRef = electron.superCluster() ;
00464   std::vector<DetId> vecId=sclusRef->getHitsByDetId () ;
00465   EcalSubdetector subdet = EcalSubdetector (vecId[0].subdetId ()) ;
00466   float thisPt =  electron.pt();
00467 
00468   if (subdet==EcalBarrel && thisPt<15.)
00469     return _EBlt15lh->getRatio ("electrons",measurements,std::string (className)) ;
00470   else if (subdet==EcalBarrel && thisPt>=15.)
00471     return _EBgt15lh->getRatio ("electrons",measurements,std::string (className)) ;
00472   else if (subdet==EcalEndcap && thisPt<15.)
00473     return _EElt15lh->getRatio ("electrons",measurements,std::string (className)) ;
00474   else if (subdet==EcalEndcap && thisPt>=15.)
00475     return _EEgt15lh->getRatio ("electrons",measurements,std::string (className)) ;
00476   else return -999. ;
00477 }
00478 

Generated on Tue Jun 9 17:43:28 2009 for CMSSW by  doxygen 1.5.4