Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:30

0001 #ifndef HLTJetsCleanedFromLeadingLeptons_h
0002 #define HLTJetsCleanedFromLeadingLeptons_h
0003 
0004 #include <vector>
0005 
0006 #include "FWCore/Framework/interface/Event.h"
0007 #include "FWCore/Framework/interface/stream/EDProducer.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0010 
0011 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0012 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0013 
0014 /**
0015  * \class HLTJetsCleanedFromLeadingLeptons
0016  * \author Andrey Popov, inspired by code by Lukasz Kreczko
0017  * \brief Produces a collection of jets cleaned against leading leptons
0018  * 
0019  * Leptons (muons or electrons) are read from results of a previous HLT filter. They are ordered in
0020  * energy, and the user can configure how many leading leptons are used to clean jets. If the
0021  * requested number is larger than the total number of leptons, all of them are used.
0022  * 
0023  * The plugin loops over the given collection of jets and exclude ones that are close to one of the
0024  * leading leptons. References to surviving jets are stored in the same format as expected by the
0025  * HLTJetCollectionsFilter plugin.
0026  */
0027 template <typename JetType>
0028 class HLTJetsCleanedFromLeadingLeptons : public edm::stream::EDProducer<> {
0029 public:
0030   typedef std::vector<JetType> JetCollection;
0031   typedef edm::Ref<JetCollection> JetRef;
0032   typedef edm::RefVector<JetCollection> JetRefVector;
0033 
0034   typedef std::vector<edm::RefVector<JetCollection, JetType, edm::refhelper::FindUsingAdvance<JetCollection, JetType>>>
0035       JetCollectionVector;
0036   //^ This is the type expected by HLTJetCollectionsFilter
0037 
0038 private:
0039   /**
0040      * \class EtaPhiE
0041      * \brief An auxiliary class to store momentum parametrised in eta, phi, and energy
0042      * 
0043      * It is useful to deal with muons and ECAL superclusters on a common basis.
0044      */
0045   class EtaPhiE {
0046   public:
0047     /// Constructor
0048     EtaPhiE(double eta, double phi, double e);
0049 
0050   public:
0051     /// Returns pseudorapidity
0052     double eta() const;
0053 
0054     /// Returns azimuthal angle
0055     double phi() const;
0056 
0057     /// Returns energy
0058     double e() const;
0059 
0060     /// A comparison operator to sort a collection of objects of this type
0061     bool operator<(EtaPhiE const &rhs) const;
0062 
0063   private:
0064     /// Pseudorapidity and azimuthal angle
0065     double etaValue, phiValue;
0066 
0067     /// Energy
0068     double eValue;
0069   };
0070 
0071 public:
0072   /// Constructor
0073   HLTJetsCleanedFromLeadingLeptons(edm::ParameterSet const &iConfig);
0074 
0075 public:
0076   /// Describes configuration of the plugin
0077   static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0078 
0079   /// Produces jets cleaned against leptons
0080   void produce(edm::Event &iEvent, edm::EventSetup const &iSetup) override;
0081 
0082 private:
0083   /// Token to identify a collection of leptons that pass an HLT filter
0084   edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs> leptonToken;
0085 
0086   /// Token to access a collection of jets
0087   edm::EDGetTokenT<std::vector<JetType>> jetToken;
0088 
0089   /// A square of the minimal allowed angular separation between a lepton and a jet
0090   double minDeltaR2;
0091 
0092   /**
0093      * \brief Number of leading leptons against which the jets are cleaned
0094      * 
0095      * If the number is larger than the total number of leptons, the jets are cleaned against all
0096      * leptons.
0097      */
0098   unsigned numLeptons;
0099 };
0100 
0101 // Implementation
0102 
0103 #include "DataFormats/Common/interface/Handle.h"
0104 
0105 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0106 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0107 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0108 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
0109 
0110 #include "DataFormats/Math/interface/deltaR.h"
0111 
0112 template <typename JetType>
0113 HLTJetsCleanedFromLeadingLeptons<JetType>::EtaPhiE::EtaPhiE(double eta, double phi, double e)
0114     : etaValue(eta), phiValue(phi), eValue(e) {}
0115 
0116 template <typename JetType>
0117 double HLTJetsCleanedFromLeadingLeptons<JetType>::EtaPhiE::eta() const {
0118   return etaValue;
0119 }
0120 
0121 template <typename JetType>
0122 double HLTJetsCleanedFromLeadingLeptons<JetType>::EtaPhiE::phi() const {
0123   return phiValue;
0124 }
0125 
0126 template <typename JetType>
0127 double HLTJetsCleanedFromLeadingLeptons<JetType>::EtaPhiE::e() const {
0128   return eValue;
0129 }
0130 
0131 template <typename JetType>
0132 bool HLTJetsCleanedFromLeadingLeptons<JetType>::EtaPhiE::operator<(EtaPhiE const &rhs) const {
0133   return (eValue < rhs.eValue);
0134 }
0135 
0136 template <typename JetType>
0137 HLTJetsCleanedFromLeadingLeptons<JetType>::HLTJetsCleanedFromLeadingLeptons(edm::ParameterSet const &iConfig)
0138     : minDeltaR2(std::pow(iConfig.getParameter<double>("minDeltaR"), 2)),
0139       numLeptons(iConfig.getParameter<unsigned>("numLeptons")) {
0140   leptonToken = consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<edm::InputTag>("leptons"));
0141   jetToken = consumes<std::vector<JetType>>(iConfig.getParameter<edm::InputTag>("jets"));
0142 
0143   produces<JetCollectionVector>();
0144 }
0145 
0146 template <typename JetType>
0147 void HLTJetsCleanedFromLeadingLeptons<JetType>::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0148   edm::ParameterSetDescription desc;
0149 
0150   desc.add<edm::InputTag>("leptons", edm::InputTag("triggerFilterObjectWithRefs"))
0151       ->setComment("A collection of leptons that pass an HLT filter");
0152   desc.add<edm::InputTag>("jets", edm::InputTag("jetCollection"))->setComment("A collection of jets");
0153   desc.add<double>("minDeltaR", 0.3)->setComment("Minimal allowed angular separation between a jet and a lepton");
0154   desc.add<unsigned>("numLeptons", 1)->setComment("Number of leading leptons against which the jets are cleaned");
0155 
0156   descriptions.add(defaultModuleLabel<HLTJetsCleanedFromLeadingLeptons<JetType>>(), desc);
0157 }
0158 
0159 template <typename JetType>
0160 void HLTJetsCleanedFromLeadingLeptons<JetType>::produce(edm::Event &iEvent, edm::EventSetup const &iSetup) {
0161   // Read results of the lepton filter
0162   edm::Handle<trigger::TriggerFilterObjectWithRefs> filterOutput;
0163   iEvent.getByToken(leptonToken, filterOutput);
0164 
0165   // Momenta of the leptons that passed the filter will be pushed into a vector
0166   std::vector<EtaPhiE> leptonMomenta;
0167 
0168   // First, assume these are muons and try to store their momenta
0169   trigger::VRmuon muons;
0170   filterOutput->getObjects(trigger::TriggerMuon, muons);
0171 
0172   for (auto const &muRef : muons)  // the collection might be empty
0173     leptonMomenta.emplace_back(muRef->eta(), muRef->phi(), muRef->energy());
0174 
0175   // Then get the momenta as if these are electrons. Electrons are tricky because they can be
0176   //stored with three types of trigger objects: TriggerElectron, TriggerPhoton, and
0177   //TriggerCluster. Try them one by one.
0178   trigger::VRelectron electrons;
0179   filterOutput->getObjects(trigger::TriggerElectron, electrons);
0180 
0181   for (auto const &eRef : electrons)  // the collection might be empty
0182   {
0183     auto const &sc = eRef->superCluster();
0184     leptonMomenta.emplace_back(sc->eta(), sc->phi(), sc->energy());
0185   }
0186 
0187   trigger::VRphoton photons;
0188   filterOutput->getObjects(trigger::TriggerPhoton, photons);
0189 
0190   for (auto const &eRef : photons)  // the collection might be empty
0191   {
0192     auto const &sc = eRef->superCluster();
0193     leptonMomenta.emplace_back(sc->eta(), sc->phi(), sc->energy());
0194   }
0195 
0196   trigger::VRphoton clusters;
0197   filterOutput->getObjects(trigger::TriggerCluster, clusters);
0198 
0199   for (auto const &eRef : clusters)  // the collection might be empty
0200   {
0201     auto const &sc = eRef->superCluster();
0202     leptonMomenta.emplace_back(sc->eta(), sc->phi(), sc->energy());
0203   }
0204 
0205   // Make sure the momenta are sorted
0206   std::sort(leptonMomenta.rbegin(), leptonMomenta.rend());
0207 
0208   // Read the source collection of jets
0209   edm::Handle<JetCollection> jetHandle;
0210   iEvent.getByToken(jetToken, jetHandle);
0211   JetCollection const &jets = *jetHandle;
0212 
0213   // Put references to jets that are not matched to leptons into a dedicated collection
0214   JetRefVector cleanedJetRefs;
0215   unsigned const numLeptonsToLoop = std::min<unsigned>(leptonMomenta.size(), numLeptons);
0216 
0217   for (unsigned iJet = 0; iJet < jets.size(); ++iJet) {
0218     bool overlap = false;
0219 
0220     for (unsigned iLepton = 0; iLepton < numLeptonsToLoop; ++iLepton)
0221       if (reco::deltaR2(leptonMomenta.at(iLepton), jets.at(iJet)) < minDeltaR2) {
0222         overlap = true;
0223         break;
0224       }
0225 
0226     if (not overlap)
0227       cleanedJetRefs.push_back(JetRef(jetHandle, iJet));
0228   }
0229 
0230   // Store the collection in the event
0231   std::unique_ptr<JetCollectionVector> product(new JetCollectionVector);
0232   //^ Have to use the depricated unique_ptr here because this is what edm::Event::put expects
0233   product->emplace_back(cleanedJetRefs);
0234   iEvent.put(std::move(product));
0235 }
0236 
0237 #endif  // HLTJetsCleanedFromLeadingLeptons_h