Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:43

0001 #include "RecoTauTag/HLTProducers/interface/L1HLTJetsMatching.h"
0002 #include "Math/GenVector/VectorUtil.h"
0003 #include "DataFormats/Common/interface/Handle.h"
0004 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0005 #include "DataFormats/L1Trigger/interface/L1JetParticle.h"
0006 #include "DataFormats/L1Trigger/interface/L1JetParticleFwd.h"
0007 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0008 #include "FWCore/Utilities/interface/EDMException.h"
0009 #include "FWCore/Utilities/interface/InputTag.h"
0010 //
0011 // class decleration
0012 //
0013 using namespace reco;
0014 using namespace std;
0015 using namespace edm;
0016 using namespace l1extra;
0017 
0018 L1HLTJetsMatching::L1HLTJetsMatching(const edm::ParameterSet& iConfig) {
0019   jetSrc = consumes<edm::View<reco::Candidate> >(iConfig.getParameter<InputTag>("JetSrc"));
0020   tauTrigger = consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<InputTag>("L1TauTrigger"));
0021   mEt_Min = iConfig.getParameter<double>("EtMin");
0022 
0023   produces<CaloJetCollection>();
0024 }
0025 
0026 void L1HLTJetsMatching::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iES) const {
0027   using namespace edm;
0028   using namespace std;
0029   using namespace reco;
0030   using namespace trigger;
0031   using namespace l1extra;
0032 
0033   typedef std::vector<LeafCandidate> LeafCandidateCollection;
0034 
0035   unique_ptr<CaloJetCollection> tauL2jets(new CaloJetCollection);
0036 
0037   constexpr double matchingR2 = 0.5 * 0.5;
0038 
0039   //Getting HLT jets to be matched
0040   edm::Handle<edm::View<Candidate> > tauJets;
0041   iEvent.getByToken(jetSrc, tauJets);
0042 
0043   //        std::cout <<"Size of input jet collection "<<tauJets->size()<<std::endl;
0044 
0045   edm::Handle<trigger::TriggerFilterObjectWithRefs> l1TriggeredTaus;
0046   iEvent.getByToken(tauTrigger, l1TriggeredTaus);
0047 
0048   std::vector<l1extra::L1JetParticleRef> tauCandRefVec;
0049   std::vector<l1extra::L1JetParticleRef> jetCandRefVec;
0050 
0051   l1TriggeredTaus->getObjects(trigger::TriggerL1TauJet, tauCandRefVec);
0052   l1TriggeredTaus->getObjects(trigger::TriggerL1CenJet, jetCandRefVec);
0053   math::XYZPoint a(0., 0., 0.);
0054   CaloJet::Specific f;
0055 
0056   for (unsigned int iL1Tau = 0; iL1Tau < tauCandRefVec.size(); iL1Tau++) {
0057     for (unsigned int iJet = 0; iJet < tauJets->size(); iJet++) {
0058       //Find the relative L2TauJets, to see if it has been reconstructed
0059       const Candidate& myJet = (*tauJets)[iJet];
0060       double deltaR2 = ROOT::Math::VectorUtil::DeltaR2(myJet.p4().Vect(), (tauCandRefVec[iL1Tau]->p4()).Vect());
0061       if (deltaR2 < matchingR2) {
0062         //       LeafCandidate myLC(myJet);
0063         CaloJet myCaloJet(myJet.p4(), a, f);
0064         if (myJet.pt() > mEt_Min) {
0065           //          tauL2LC->push_back(myLC);
0066           tauL2jets->push_back(myCaloJet);
0067         }
0068         break;
0069       }
0070     }
0071   }
0072 
0073   for (unsigned int iL1Tau = 0; iL1Tau < jetCandRefVec.size(); iL1Tau++) {
0074     for (unsigned int iJet = 0; iJet < tauJets->size(); iJet++) {
0075       const Candidate& myJet = (*tauJets)[iJet];
0076       //Find the relative L2TauJets, to see if it has been reconstructed
0077       double deltaR2 = ROOT::Math::VectorUtil::DeltaR2(myJet.p4().Vect(), (jetCandRefVec[iL1Tau]->p4()).Vect());
0078       if (deltaR2 < matchingR2) {
0079         //       LeafCandidate myLC(myJet);
0080         CaloJet myCaloJet(myJet.p4(), a, f);
0081         if (myJet.pt() > mEt_Min) {
0082           //tauL2LC->push_back(myLC);
0083           tauL2jets->push_back(myCaloJet);
0084         }
0085         break;
0086       }
0087     }
0088   }
0089 
0090   //std::cout <<"Size of L1HLT matched jets "<<tauL2jets->size()<<std::endl;
0091 
0092   iEvent.put(std::move(tauL2jets));
0093   // iEvent.put(std::move(tauL2LC));
0094 }
0095 
0096 #include "FWCore/Framework/interface/MakerMacros.h"
0097 DEFINE_FWK_MODULE(L1HLTJetsMatching);