Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include <cmath>
0002 #include <string>
0003 #include <vector>
0004 
0005 #include "CommonTools/Utils/interface/PtComparator.h"
0006 #include "DataFormats/Common/interface/Handle.h"
0007 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0008 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0009 #include "DataFormats/JetReco/interface/PFJet.h"
0010 #include "DataFormats/Math/interface/deltaR.h"
0011 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0012 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
0013 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
0014 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
0015 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0016 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0018 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0019 #include "HLTrigger/JetMET/interface/HLTJetCollectionsForBoostedLeptonPlusJets.h"
0020 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0021 
0022 template <typename jetType>
0023 HLTJetCollectionsForBoostedLeptonPlusJets<jetType>::HLTJetCollectionsForBoostedLeptonPlusJets(
0024     const edm::ParameterSet& iConfig)
0025     : hltLeptonTag(iConfig.getParameter<edm::InputTag>("HltLeptonTag")),
0026       sourceJetTag(iConfig.getParameter<edm::InputTag>("SourceJetTag")),
0027       // minimum delta-R^2 threshold with sign
0028       minDeltaR2_(iConfig.getParameter<double>("minDeltaR") * std::abs(iConfig.getParameter<double>("minDeltaR"))) {
0029   using namespace edm;
0030   using namespace std;
0031 
0032   typedef vector<RefVector<vector<jetType>, jetType, refhelper::FindUsingAdvance<vector<jetType>, jetType>>>
0033       JetCollectionVector;
0034   typedef vector<jetType> JetCollection;
0035 
0036   m_theLeptonToken = consumes<trigger::TriggerFilterObjectWithRefs>(hltLeptonTag);
0037   m_theJetToken = consumes<std::vector<jetType>>(sourceJetTag);
0038   produces<JetCollectionVector>();
0039   produces<JetCollection>();
0040 }
0041 
0042 template <typename jetType>
0043 HLTJetCollectionsForBoostedLeptonPlusJets<jetType>::~HLTJetCollectionsForBoostedLeptonPlusJets() {
0044   // do anything here that needs to be done at desctruction time
0045   // (e.g. close files, deallocate resources etc.)
0046 }
0047 
0048 template <typename jetType>
0049 void HLTJetCollectionsForBoostedLeptonPlusJets<jetType>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0050   edm::ParameterSetDescription desc;
0051   desc.add<edm::InputTag>("HltLeptonTag", edm::InputTag("triggerFilterObjectWithRefs"));
0052   //(2)
0053   desc.add<edm::InputTag>("SourceJetTag", edm::InputTag("jetCollection"));
0054   //(2)
0055   desc.add<double>("minDeltaR", 0.5);
0056   descriptions.add(defaultModuleLabel<HLTJetCollectionsForBoostedLeptonPlusJets<jetType>>(), desc);
0057 }
0058 
0059 //
0060 // member functions
0061 //
0062 
0063 // ------------ method called to produce the data  ------------
0064 // template <typename T>
0065 template <typename jetType>
0066 void HLTJetCollectionsForBoostedLeptonPlusJets<jetType>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0067   using namespace edm;
0068   using namespace std;
0069   //(3)
0070   using namespace reco;
0071   //(3)
0072   typedef vector<RefVector<vector<jetType>, jetType, refhelper::FindUsingAdvance<vector<jetType>, jetType>>>
0073       JetCollectionVector;
0074   typedef vector<jetType> JetCollection;
0075   typedef edm::Ref<JetCollection> JetRef;
0076   typedef edm::RefVector<JetCollection> JetRefVector;
0077 
0078   Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
0079   iEvent.getByToken(m_theLeptonToken, PrevFilterOutput);
0080 
0081   //its easier on the if statement flow if I try everything at once, shouldnt add to timing
0082   // Electrons can be stored as objects of types TriggerCluster, TriggerElectron, or TriggerPhoton
0083   vector<reco::RecoChargedCandidateRef> muonCands;
0084   PrevFilterOutput->getObjects(trigger::TriggerMuon, muonCands);
0085 
0086   vector<Ref<reco::ElectronCollection>> eleCands;
0087   PrevFilterOutput->getObjects(trigger::TriggerElectron, eleCands);
0088 
0089   trigger::VRphoton photonCands;
0090   PrevFilterOutput->getObjects(trigger::TriggerPhoton, photonCands);
0091 
0092   vector<Ref<reco::RecoEcalCandidateCollection>> clusCands;
0093   PrevFilterOutput->getObjects(trigger::TriggerCluster, clusCands);
0094 
0095   Handle<JetCollection> theJetCollectionHandle;
0096   iEvent.getByToken(m_theJetToken, theJetCollectionHandle);
0097 
0098   typename JetCollection::const_iterator jet;
0099 
0100   unique_ptr<JetCollection> allSelections(new JetCollection);
0101   unique_ptr<JetCollectionVector> product(new JetCollectionVector);
0102 
0103   std::vector<size_t> usedCands;
0104 
0105   if (!muonCands.empty()) {  // muons
0106     for (jet = theJetCollectionHandle->begin(); jet != theJetCollectionHandle->end(); jet++) {
0107       //const jetType* referenceJet = &*jet;
0108       jetType cleanedJet = *jet;  //copy original jet
0109       for (size_t candNr = 0; candNr < muonCands.size(); candNr++) {
0110         if (std::find(usedCands.begin(), usedCands.end(), candNr) != usedCands.end())
0111           continue;
0112         if (reco::deltaR2((*muonCands[candNr]), cleanedJet) <= minDeltaR2_) {
0113           std::vector<edm::Ptr<reco::PFCandidate>> pfConstituents = cleanedJet.getPFConstituents();
0114           for (std::vector<edm::Ptr<reco::PFCandidate>>::const_iterator i_candidate = pfConstituents.begin();
0115                i_candidate != pfConstituents.end();
0116                ++i_candidate) {
0117             if (reco::deltaR2((*muonCands[candNr]), (**i_candidate)) < 1e-6) {
0118               cleanedJet.setP4(cleanedJet.p4() - muonCands[candNr]->p4());
0119               usedCands.push_back(candNr);
0120               break;
0121             }  //if constituent matched
0122           }    //for constituents
0123         }      //if dR<min
0124       }        //for cands
0125       allSelections->push_back(cleanedJet);
0126     }  //for jets
0127   }    //if cands
0128 
0129   if (!eleCands.empty()) {  // electrons
0130     for (jet = theJetCollectionHandle->begin(); jet != theJetCollectionHandle->end(); jet++) {
0131       //const jetType* referenceJet = &*jet;
0132       jetType cleanedJet = *jet;  //copy original jet
0133       for (size_t candNr = 0; candNr < eleCands.size(); candNr++) {
0134         if (std::find(usedCands.begin(), usedCands.end(), candNr) != usedCands.end())
0135           continue;
0136         if (reco::deltaR2((*eleCands[candNr]), cleanedJet) <= minDeltaR2_) {
0137           std::vector<edm::Ptr<reco::PFCandidate>> pfConstituents = cleanedJet.getPFConstituents();
0138           for (std::vector<edm::Ptr<reco::PFCandidate>>::const_iterator i_candidate = pfConstituents.begin();
0139                i_candidate != pfConstituents.end();
0140                ++i_candidate) {
0141             if (reco::deltaR2((*eleCands[candNr]), (**i_candidate)) < 1e-6) {
0142               cleanedJet.setP4(cleanedJet.p4() - eleCands[candNr]->p4());
0143               usedCands.push_back(candNr);
0144               break;
0145             }  //if constituent matched
0146           }    //for constituents
0147         }      //if dR<min
0148       }        //for cands
0149       allSelections->push_back(cleanedJet);
0150     }  //for jets
0151   }    //if cands
0152 
0153   if (!photonCands.empty()) {  // photons
0154     for (jet = theJetCollectionHandle->begin(); jet != theJetCollectionHandle->end(); jet++) {
0155       //const jetType* referenceJet = &*jet;
0156       jetType cleanedJet = *jet;  //copy original jet
0157       for (size_t candNr = 0; candNr < photonCands.size(); candNr++) {
0158         if (std::find(usedCands.begin(), usedCands.end(), candNr) != usedCands.end())
0159           continue;
0160         if (reco::deltaR2((*photonCands[candNr]), cleanedJet) <= minDeltaR2_) {
0161           std::vector<edm::Ptr<reco::PFCandidate>> pfConstituents = cleanedJet.getPFConstituents();
0162           for (std::vector<edm::Ptr<reco::PFCandidate>>::const_iterator i_candidate = pfConstituents.begin();
0163                i_candidate != pfConstituents.end();
0164                ++i_candidate) {
0165             if (reco::deltaR2((*photonCands[candNr]), (**i_candidate)) < 1e-6) {
0166               cleanedJet.setP4(cleanedJet.p4() - photonCands[candNr]->p4());
0167               usedCands.push_back(candNr);
0168               break;
0169             }  //if constituent matched
0170           }    //for constituents
0171         }      //if dR<min
0172       }        //for cands
0173       allSelections->push_back(cleanedJet);
0174     }  //for jets
0175   }    //if cands
0176 
0177   if (!clusCands.empty()) {  // trigger clusters
0178     for (jet = theJetCollectionHandle->begin(); jet != theJetCollectionHandle->end(); jet++) {
0179       //const jetType* referenceJet = &*jet;
0180       jetType cleanedJet = *jet;  //copy original jet
0181       for (size_t candNr = 0; candNr < clusCands.size(); candNr++) {
0182         if (std::find(usedCands.begin(), usedCands.end(), candNr) != usedCands.end())
0183           continue;
0184         if (reco::deltaR2((*clusCands[candNr]), cleanedJet) <= minDeltaR2_) {
0185           std::vector<edm::Ptr<reco::PFCandidate>> pfConstituents = cleanedJet.getPFConstituents();
0186           for (std::vector<edm::Ptr<reco::PFCandidate>>::const_iterator i_candidate = pfConstituents.begin();
0187                i_candidate != pfConstituents.end();
0188                ++i_candidate) {
0189             if (reco::deltaR2((*clusCands[candNr]), (**i_candidate)) < 1e-6) {
0190               cleanedJet.setP4(cleanedJet.p4() - clusCands[candNr]->p4());
0191               usedCands.push_back(candNr);
0192               break;
0193             }  //if constituent matched
0194           }    //for constituents
0195         }      //if dR<min
0196       }        //for cands
0197       allSelections->push_back(cleanedJet);
0198     }  //for jets
0199   }    //if cands
0200 
0201   NumericSafeGreaterByPt<jetType> compJets;
0202   // reorder cleaned jets
0203   std::sort(allSelections->begin(), allSelections->end(), compJets);
0204   edm::OrphanHandle<JetCollection> cleanedJetHandle = iEvent.put(std::move(allSelections));
0205 
0206   JetCollection const& jets = *cleanedJetHandle;
0207 
0208   JetRefVector cleanedJetRefs;
0209   cleanedJetRefs.reserve(jets.size());
0210   for (unsigned iJet = 0; iJet < jets.size(); ++iJet) {
0211     cleanedJetRefs.push_back(JetRef(cleanedJetHandle, iJet));
0212   }
0213 
0214   product->emplace_back(cleanedJetRefs);
0215   iEvent.put(std::move(product));
0216 
0217   return;
0218 }