CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/RecoTauTag/RecoTau/interface/HPSPFRecoTauAlgorithm.h

Go to the documentation of this file.
00001 /*
00002 Hadrons + Strips Tau Identification Algorithm
00003 ---------------------------------------------
00004 Michail Bachtis
00005 University of Wisconsin-Madison 
00006 bachtis@cern.ch
00007 */
00008 
00009 #ifndef RecoTauTag_RecoTau_HPSPFTauAlgorithm
00010 #define RecoTauTag_RecoTau_HPSPFTauAlgorithm
00011 
00012 #include "RecoTauTag/TauTagTools/interface/PFCandidateStripMerger.h"
00013 #include "RecoTauTag/RecoTau/interface/PFRecoTauAlgorithmBase.h"
00014 #include "DataFormats/MuonReco/interface/Muon.h"
00015 #include "DataFormats/MuonReco/interface/MuonFwd.h"
00016 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00017 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
00018 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00019 #include "FWCore/Utilities/interface/Exception.h"
00020 #include "TrackingTools/IPTools/interface/IPTools.h"
00021 #include "RecoTauTag/TauTagTools/interface/TauTagTools.h"
00022 #include "RecoTauTag/TauTagTools/interface/PFTauElementsOperators.h"
00023 
00024 
00025 class HPSPFRecoTauAlgorithm : public PFRecoTauAlgorithmBase
00026 {
00027  public:
00028   //  Constructors are following the PFRecoTauProducer scheme. 
00029   //  The input is a ParameterSet of the algorithm configuration
00030   HPSPFRecoTauAlgorithm();
00031   HPSPFRecoTauAlgorithm(const edm::ParameterSet&);
00032   ~HPSPFRecoTauAlgorithm();
00033 
00034   //Basic Method that creates the taus 
00035   reco::PFTau buildPFTau(const reco::PFTauTagInfoRef&,const reco::Vertex&); 
00036 
00037  private:   
00038   //*Private Members *//
00039   PFCandidateMergerBase * candidateMerger_;
00040 
00041 
00042   //* Helper Methods *//
00043 
00044   //Creators of the Decay Modes
00045   void buildOneProng(const reco::PFTauTagInfoRef&,const reco::PFCandidateRefVector& );
00046   void buildOneProngStrip(const reco::PFTauTagInfoRef&,const std::vector<reco::PFCandidateRefVector>&,const reco::PFCandidateRefVector&);
00047   void buildOneProngTwoStrips(const reco::PFTauTagInfoRef&,const std::vector<reco::PFCandidateRefVector>&,const reco::PFCandidateRefVector&);
00048   void buildThreeProngs(const reco::PFTauTagInfoRef&,const reco::PFCandidateRefVector&);
00049 
00050   //Narrowness selection
00051   bool isNarrowTau(const reco::PFTau&,double);
00052 
00053   //Associate Isolation Candidates
00054   void associateIsolationCandidates(reco::PFTau&,double);
00055 
00056   reco::PFTau getBestTauCandidate(reco::PFTauCollection&); //get the best tau if we have overlap 
00057   void applyMuonRejection(reco::PFTau&);
00058 
00059   //Apply electron Rejection
00060   void applyElectronRejection(reco::PFTau&,double);
00061 
00062   //Method to create a candidate from the merged EM Candidates vector;
00063   math::XYZTLorentzVector createMergedLorentzVector(const reco::PFCandidateRefVector&);
00064 
00065   void removeCandidateFromRefVector(const reco::PFCandidateRef&,reco::PFCandidateRefVector&); 
00066   void applyMassConstraint(math::XYZTLorentzVector&,double );
00067 
00068   bool refitThreeProng(reco::PFTau&); 
00069 
00070 
00071 
00072   //Configure the algorithm!
00073   void configure(const edm::ParameterSet&);
00074 
00075   //* Variables for configuration*//
00076   
00077   //Which merging algorithm to use
00078   std::string emMerger_; 
00079 
00080   //Overlap Removal Criterion for overlaping decay modes
00081   std::string overlapCriterion_;
00082 
00083   //Decay Mode activation
00084   bool doOneProngs_;
00085   bool doOneProngStrips_;
00086   bool doOneProngTwoStrips_;
00087   bool doThreeProngs_;
00088 
00089 
00090   //Minimum Tau Pt for the reconstruction
00091   double tauThreshold_;
00092   //Lead Pion Threshold
00093   double leadPionThreshold_;
00094 
00095   //strip threshold
00096   double stripPtThreshold_;
00097   //Isolation Cone sizes for different particles
00098   double chargeIsolationCone_;
00099   double gammaIsolationCone_;
00100   double neutrHadrIsolationCone_;
00101 
00102   //Use a solid cone or an annulus?
00103   bool   useIsolationAnnulus_;
00104 
00105   //Mass Windows 
00106   std::vector<double> oneProngStripMassWindow_;
00107   std::vector<double> oneProngTwoStripsMassWindow_;
00108   std::vector<double> oneProngTwoStripsPi0MassWindow_;
00109   std::vector<double> threeProngMassWindow_;
00110 
00111   //Matching Cone between PFTau , PFJet 
00112   double matchingCone_;
00113 
00114 
00115   //Cone Narrowness. Unfortunately we have to remove TFormula! 
00116   //The code supports A(DeltaR or Angle) < coneParameter_/(ET or ENERGY)
00117   std::string coneSizeFormula_;
00118   std::string coneMetric_;
00119   double minSignalCone_;
00120   double maxSignalCone_;
00121 
00122 
00123   TFormula coneSizeFormula;
00124 
00125   reco::PFTauCollection pfTaus_;
00126 
00127 
00128 
00129   class HPSTauPtSorter {
00130   public:
00131 
00132     HPSTauPtSorter()
00133       {
00134       }
00135 
00136 
00137     ~HPSTauPtSorter()
00138       {}
00139 
00140     bool operator()(reco::PFTau a , reco::PFTau b) {
00141       return (a.pt() > b.pt());
00142     }
00143   };
00144 
00145 
00146   class HPSTauIsolationSorter {
00147   public:
00148 
00149     HPSTauIsolationSorter()
00150       {
00151       }
00152 
00153 
00154     ~HPSTauIsolationSorter()
00155       {}
00156 
00157     bool operator()(reco::PFTau a , reco::PFTau b) {
00158       return (a.isolationPFGammaCandsEtSum()+a.isolationPFChargedHadrCandsPtSum())<
00159         (b.isolationPFGammaCandsEtSum()+b.isolationPFChargedHadrCandsPtSum());
00160     }
00161   };
00162 
00163   
00164 };
00165 #endif
00166 
00167 
00168