CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch1/src/RecoTauTag/RecoTau/interface/AntiElectronIDMVA.h

Go to the documentation of this file.
00001 //--------------------------------------------------------------------------------------------------
00002 // $Id $
00003 //
00004 // AntiElectronIDMVA
00005 //
00006 // Helper Class for applying MVA anti-electron discrimination
00007 //
00008 // Authors: L.Bianchini
00009 //--------------------------------------------------------------------------------------------------
00010 
00011 /*
00012   proposed WP:    epsilonB ~ 18%  epsilonS ~ 91% wrt signal taus passing discr. ag. electrons Medium. 
00013   bool pass = 
00014   (abs(TauEta)<1.5 && TauSignalPFGammaCands==0 && MVAValue(...)>0.054) ||
00015   (abs(TauEta)<1.5 && TauSignalPFGammaCands>0  && TauHasGsf>0.5 && MVAValue(...)>0.060) ||
00016   (abs(TauEta)<1.5 && TauSignalPFGammaCands>0  && TauHasGsf<0.5 && MVAValue(...)>0.054) ||
00017   (abs(TauEta)>1.5 && TauSignalPFGammaCands==0 && MVAValue(...)>0.060) ||
00018   (abs(TauEta)>1.5 && TauSignalPFGammaCands>0  && TauHasGsf>0.5 && MVAValue(...)>0.053) ||
00019   (abs(TauEta)>1.5 && TauSignalPFGammaCands>0  && TauHasGsf<0.5 && MVAValue(...)>0.049);
00020 */
00021 
00022 
00023 #ifndef RECOTAUTAG_RECOTAU_AntiElectronIDMVA_H
00024 #define RECOTAUTAG_RECOTAU_AntiElectronIDMVA_H
00025 
00026 #include "DataFormats/TauReco/interface/PFTau.h"
00027 
00028 #include "TMVA/Tools.h"
00029 #include "TMVA/Reader.h"
00030 
00031 #include <vector>
00032 
00033 class AntiElectronIDMVA {
00034   public:
00035 
00036     AntiElectronIDMVA();
00037     ~AntiElectronIDMVA(); 
00038 
00039     void   Initialize(std::string methodName,
00040                       std::string oneProng0Pi0_BL,
00041                       std::string oneProng1pi0wGSF_BL,
00042                       std::string oneProng1pi0woGSF_BL,
00043                       std::string oneProng0Pi0_EC,
00044                       std::string oneProng1pi0wGSF_EC,
00045                       std::string oneProng1pi0woGSF_EC
00046                       );
00047 
00048     // RECOMMENDED:
00049     double MVAValue(Float_t TauEta, Float_t TauPt,
00050                     Float_t TauSignalPFChargedCands, Float_t TauSignalPFGammaCands, 
00051                     Float_t TauLeadPFChargedHadrMva, 
00052                     Float_t TauLeadPFChargedHadrHoP, Float_t TauLeadPFChargedHadrEoP, 
00053                     Float_t TauHasGsf, Float_t TauVisMass,  Float_t TauEmFraction,
00054                     std::vector<Float_t>* GammasdEta, std::vector<Float_t>* GammasdPhi, std::vector<Float_t>* GammasPt
00055                     );
00056 
00057     /* 
00058     where:
00059 
00060     TauEta                  = myTau->eta();
00061     TauPt                   = myTau->pt();
00062     TauSignalPFChargedCands = myTau->signalPFChargedHadrCands().size();
00063     TauSignalPFGammaCands   = myTau->signalPFGammaCands().size();
00064     TauLeadPFChargedHadrMva = myTau->electronPreIDOutput();
00065     TauLeadPFChargedHadrHoP = myTau->leadPFChargedHadrCand()->hcalEnergy()/myTau->leadPFChargedHadrCand()->p();
00066     TauLeadPFChargedHadrEoP = myTau->leadPFChargedHadrCand()->ecalEnergy()/myTau->leadPFChargedHadrCand()->p();
00067     TauHasGsf               = (myTau->leadPFChargedHadrCand()->gsfTrackRef()).isNonnull();
00068     TauVisMass              = myTau->mass();
00069     TauEmFraction           = myTau->emFraction();
00070   
00071     GammasdEta     = new std::vector< float >();
00072     GammasdPhi     = new std::vector< float >();
00073     GammasPt       = new std::vector< float >();
00074     
00075     for(unsigned int k = 0 ; k < (myTau->signalPFGammaCands()).size() ; k++){
00076     reco::PFCandidateRef gamma = (myTau->signalPFGammaCands()).at(k);
00077     if( (myTau->leadPFChargedHadrCand()).isNonnull() ){
00078         GammasdEta->push_back( gamma->eta() - myTau->leadPFChargedHadrCand()->eta() );
00079         GammasdPhi->push_back( gamma->phi() - myTau->leadPFChargedHadrCand()->phi() );
00080     }
00081     else{
00082         GammasdEta->push_back( gamma->eta() - myTau->eta() );
00083         GammasdPhi->push_back( gamma->phi() - myTau->phi() );
00084     }
00085      GammasPt->push_back(  gamma->pt() );
00086     }
00087     */
00088 
00089     double MVAValue(Float_t TauEta,  Float_t TauPt,
00090                     Float_t TauSignalPFChargedCands, Float_t TauSignalPFGammaCands, 
00091                     Float_t TauLeadPFChargedHadrMva, 
00092                     Float_t TauLeadPFChargedHadrHoP , Float_t TauLeadPFChargedHadrEoP, 
00093                     Float_t TauHasGsf, Float_t TauVisMass,  Float_t TauEmFraction,
00094                     Float_t GammaEtaMom, Float_t GammaPhiMom, Float_t GammaEnFrac
00095                     );
00096     /*
00097       see AntiElectronIDMVA.cc for GammaEtaMom,GammaPhiMom,GammaEnFrac
00098     */
00099 
00100 
00101     double MVAValue(const reco::PFTauRef& thePFTauRef);
00102 
00103  private:
00104 
00105     Bool_t isInitialized_;
00106     std::string methodName_;
00107     TMVA::Reader* fTMVAReader_[6];
00108     Float_t TauSignalPFGammaCands_; 
00109     Float_t TauVisMass_; 
00110     Float_t GammadEta_; 
00111     Float_t GammadPhi_; 
00112     Float_t GammadPt_;
00113     Float_t TauLeadPFChargedHadrMva_;
00114     Float_t TauLeadPFChargedHadrHoP_;
00115     Float_t TauLeadPFChargedHadrEoP_;
00116     Float_t TauEmFraction_;
00117     
00118 };
00119 
00120 #endif