Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:09:30

0001 #ifndef HLTCATopTagFilter_h
0002 #define HLTCATopTagFilter_h
0003 
0004 // system include files
0005 #include <memory>
0006 #include <vector>
0007 #include <sstream>
0008 
0009 // user include files
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/MakerMacros.h"
0012 #include "FWCore/PluginManager/interface/ModuleDef.h"
0013 #include "FWCore/Framework/interface/Frameworkfwd.h"
0014 #include "DataFormats/BTauReco/interface/CATopJetTagInfo.h"
0015 #include "FWCore/Utilities/interface/InputTag.h"
0016 #include "DataFormats/Math/interface/deltaR.h"
0017 #include "DataFormats/JetReco/interface/BasicJet.h"
0018 #include "DataFormats/Candidate/interface/CompositeCandidate.h"
0019 #include "CommonTools/CandUtils/interface/AddFourMomenta.h"
0020 #include "DataFormats/Candidate/interface/CandMatchMap.h"
0021 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0022 #include "HLTrigger/HLTcore/interface/HLTFilter.h"
0023 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0024 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0025 #include <Math/VectorUtil.h>
0026 
0027 class CATopJetHelperUser {
0028 public:
0029   CATopJetHelperUser(double TopMass) : TopMass_(TopMass) {}
0030 
0031   reco::CATopJetProperties operator()(reco::Jet const& ihardJet) const;
0032 
0033 protected:
0034   double TopMass_;
0035 };
0036 
0037 struct GreaterByPtCandPtrUser {
0038   bool operator()(const edm::Ptr<reco::Candidate>& t1, const edm::Ptr<reco::Candidate>& t2) const {
0039     return t1->pt() > t2->pt();
0040   }
0041 };
0042 
0043 //
0044 // class declaration
0045 //
0046 
0047 class HLTCATopTagFilter : public HLTFilter {
0048 public:
0049   explicit HLTCATopTagFilter(const edm::ParameterSet&);
0050   ~HLTCATopTagFilter() override;
0051   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0052   bool hltFilter(edm::Event&,
0053                  const edm::EventSetup&,
0054                  trigger::TriggerFilterObjectWithRefs& filterobject) const override;
0055 
0056 private:
0057   // ----------member data ---------------------------
0058 
0059   edm::InputTag src_;
0060   edm::InputTag pfsrc_;
0061   const edm::EDGetTokenT<reco::BasicJetCollection> inputToken_;
0062   const edm::EDGetTokenT<reco::PFJetCollection> inputPFToken_;
0063   double TopMass_;
0064   double minTopMass_;
0065   double maxTopMass_;
0066   double minMinMass_;
0067 };
0068 
0069 inline reco::CATopJetProperties CATopJetHelperUser::operator()(reco::Jet const& ihardJet) const {
0070   reco::CATopJetProperties properties;
0071   // Get subjets
0072   reco::Jet::Constituents subjets = ihardJet.getJetConstituents();
0073   properties.nSubJets = subjets.size();  // number of subjets
0074   properties.topMass = ihardJet.mass();  // jet mass
0075   properties.minMass = 999999.;          // minimum mass pairing
0076 
0077   // Require at least three subjets in all cases, if not, untagged
0078   if (properties.nSubJets >= 3) {
0079     // Take the highest 3 pt subjets for cuts
0080     sort(subjets.begin(), subjets.end(), GreaterByPtCandPtrUser());
0081 
0082     // Now look at the subjets that were formed
0083     for (int isub = 0; isub < 2; ++isub) {
0084       // Get this subjet
0085       reco::Jet::Constituent icandJet = subjets[isub];
0086 
0087       // Now look at the "other" subjets than this one, form the minimum invariant mass
0088       // pairing, as well as the "closest" combination to the W mass
0089       for (int jsub = isub + 1; jsub < 3; ++jsub) {
0090         // Get the second subjet
0091         reco::Jet::Constituent jcandJet = subjets[jsub];
0092 
0093         reco::Candidate::LorentzVector wCand = icandJet->p4() + jcandJet->p4();
0094 
0095         // Get the candidate mass
0096         double imw = wCand.mass();
0097 
0098         // Find the minimum mass pairing.
0099         if (fabs(imw) < properties.minMass) {
0100           properties.minMass = imw;
0101         }
0102       }  // end second loop over subjets
0103     }    // end first loop over subjets
0104   }      // endif 3 subjets
0105 
0106   if (properties.minMass == 999999)
0107     properties.minMass = -1;
0108 
0109   return properties;
0110 }
0111 
0112 #endif