Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:26:44

0001 
0002 
0003 // system include files
0004 #include <memory>
0005 #include <iostream>
0006 
0007 // user include files
0008 #include "FWCore/Framework/interface/Frameworkfwd.h"
0009 #include "FWCore/Framework/interface/global/EDFilter.h"
0010 
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/MakerMacros.h"
0013 
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 
0016 #include "DataFormats/TrackReco/interface/Track.h"
0017 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0018 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
0019 #include "DataFormats/METReco/interface/PFMET.h"
0020 #include "DataFormats/METReco/interface/PFMETFwd.h"
0021 #include "DataFormats/METReco/interface/PFMETCollection.h"
0022 //
0023 // class declaration
0024 //
0025 
0026 class ChargedHadronTrackResolutionFilter : public edm::global::EDFilter<> {
0027 public:
0028   explicit ChargedHadronTrackResolutionFilter(const edm::ParameterSet&);
0029   ~ChargedHadronTrackResolutionFilter() override;
0030 
0031 private:
0032   bool filter(edm::StreamID iID, edm::Event&, const edm::EventSetup&) const override;
0033 
0034   // ----------member data ---------------------------
0035 
0036   edm::EDGetTokenT<reco::PFCandidateCollection> tokenPFCandidates_;
0037   edm::EDGetTokenT<reco::PFMETCollection> tokenPFMET_;
0038   const bool taggingMode_;
0039   const double ptMin_;
0040   const double metSignifMin_;
0041   const double p1_;
0042   const double p2_;
0043   const double p3_;
0044   const bool debug_;
0045 };
0046 
0047 //
0048 // constants, enums and typedefs
0049 //
0050 
0051 //
0052 // static data member definitions
0053 //
0054 
0055 //
0056 // constructors and destructor
0057 //
0058 ChargedHadronTrackResolutionFilter::ChargedHadronTrackResolutionFilter(const edm::ParameterSet& iConfig)
0059     : tokenPFCandidates_(consumes<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("PFCandidates"))),
0060       tokenPFMET_(consumes<reco::PFMETCollection>(iConfig.getParameter<edm::InputTag>("PFMET"))),
0061       taggingMode_(iConfig.getParameter<bool>("taggingMode")),
0062       ptMin_(iConfig.getParameter<double>("ptMin")),
0063       metSignifMin_(iConfig.getParameter<double>("MetSignifMin")),
0064       p1_(iConfig.getParameter<double>("p1")),
0065       p2_(iConfig.getParameter<double>("p2")),
0066       p3_(iConfig.getParameter<double>("p3")),
0067       debug_(iConfig.getParameter<bool>("debug")) {
0068   produces<bool>();
0069 }
0070 
0071 ChargedHadronTrackResolutionFilter::~ChargedHadronTrackResolutionFilter() {}
0072 
0073 //
0074 // member functions
0075 //
0076 
0077 // ------------ method called on each new Event  ------------
0078 bool ChargedHadronTrackResolutionFilter::filter(edm::StreamID iID,
0079                                                 edm::Event& iEvent,
0080                                                 const edm::EventSetup& iSetup) const {
0081   using namespace std;
0082   using namespace edm;
0083 
0084   Handle<reco::PFCandidateCollection> pfCandidates;
0085   iEvent.getByToken(tokenPFCandidates_, pfCandidates);
0086   Handle<reco::PFMETCollection> pfMET;
0087   iEvent.getByToken(tokenPFMET_, pfMET);
0088 
0089   bool foundBadTrack = false;
0090   if (debug_)
0091     cout << "starting loop over pfCandidates" << endl;
0092 
0093   for (unsigned i = 0; i < pfCandidates->size(); ++i) {
0094     const reco::PFCandidate& cand = (*pfCandidates)[i];
0095 
0096     if (std::abs(cand.pdgId()) == 211) {
0097       if (cand.trackRef().isNull())
0098         continue;
0099 
0100       const reco::TrackRef trackref = cand.trackRef();
0101       const double Pt = trackref->pt();
0102       const double DPt = trackref->ptError();
0103       if (Pt < ptMin_)
0104         continue;
0105       if (debug_)
0106         cout << "charged hadron track pT > " << Pt << " GeV - "
0107              << " dPt: " << DPt << " GeV - algorithm: " << trackref->algo() << std::endl;
0108 
0109       const double P = trackref->p();
0110 
0111       const unsigned int LostHits = trackref->numberOfLostHits();
0112 
0113       if ((DPt / Pt) > (p1_ * sqrt(p2_ * p2_ / P + p3_ * p3_) / (1. + LostHits))) {
0114         const double MET_px = pfMET->begin()->px();
0115         const double MET_py = pfMET->begin()->py();
0116         const double MET_et = pfMET->begin()->et();
0117         const double MET_sumEt = pfMET->begin()->sumEt();
0118         const double hadron_px = cand.px();
0119         const double hadron_py = cand.py();
0120         if (MET_sumEt == 0)
0121           continue;
0122         const double MET_signif = MET_et / MET_sumEt;
0123         const double MET_et_corr =
0124             sqrt((MET_px + hadron_px) * (MET_px + hadron_px) + (MET_py + hadron_py) * (MET_py + hadron_py));
0125         const double MET_signif_corr = MET_et_corr / MET_sumEt;
0126         if (debug_)
0127           std::cout << "MET signif before: " << MET_signif << " - after: " << MET_signif_corr
0128                     << " - reduction factor: " << MET_signif / MET_signif_corr << endl;
0129 
0130         if ((MET_signif / MET_signif_corr) > metSignifMin_) {
0131           foundBadTrack = true;
0132 
0133           if (debug_) {
0134             cout << cand << endl;
0135             cout << "charged hadron \t"
0136                  << "track pT = " << Pt << " +/- " << DPt;
0137             cout << endl;
0138             cout << "MET: " << MET_et << " MET phi: " << pfMET->begin()->phi() << " MET sumet: " << MET_sumEt
0139                  << " MET significance: " << MET_et / MET_sumEt << endl;
0140             cout << "MET_px: " << MET_px << " MET_py: " << MET_py << " hadron_px: " << hadron_px
0141                  << " hadron_py: " << hadron_py << endl;
0142             cout << "corrected: " << sqrt(pow((MET_px + hadron_px), 2) + pow((MET_py + hadron_py), 2)) << endl;
0143           }
0144         }
0145       }
0146     }
0147   }  // end loop over PF candidates
0148 
0149   bool pass = !foundBadTrack;
0150 
0151   iEvent.put(std::make_unique<bool>(pass));
0152 
0153   return taggingMode_ || pass;
0154 }
0155 
0156 //define this as a plug-in
0157 DEFINE_FWK_MODULE(ChargedHadronTrackResolutionFilter);