Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "RecoTauTag/HLTProducers/interface/L2TauJetsMerger.h"
0002 #include "Math/GenVector/VectorUtil.h"
0003 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0004 #include "FWCore/Utilities/interface/EDMException.h"
0005 //
0006 // class decleration
0007 //
0008 using namespace reco;
0009 using namespace std;
0010 using namespace edm;
0011 
0012 L2TauJetsMerger::L2TauJetsMerger(const edm::ParameterSet& iConfig)
0013     : jetSrc(iConfig.getParameter<vtag>("JetSrc")), mEt_Min(iConfig.getParameter<double>("EtMin")) {
0014   for (vtag::const_iterator it = jetSrc.begin(); it != jetSrc.end(); ++it) {
0015     edm::EDGetTokenT<CaloJetCollection> aToken = consumes<CaloJetCollection>(*it);
0016     jetSrc_token.push_back(aToken);
0017   }
0018 
0019   produces<CaloJetCollection>();
0020 }
0021 
0022 L2TauJetsMerger::~L2TauJetsMerger() {}
0023 
0024 void L2TauJetsMerger::produce(edm::StreamID iSId, edm::Event& iEvent, const edm::EventSetup& iES) const {
0025   using namespace edm;
0026   using namespace std;
0027   using namespace reco;
0028 
0029   //Getting the Collections of L2ReconstructedJets from L1Seeds
0030   //and removing the collinear jets
0031   CaloJetCollection myTmpJets;
0032 
0033   for (vtoken_cjets::const_iterator s = jetSrc_token.begin(); s != jetSrc_token.end(); ++s) {
0034     edm::Handle<CaloJetCollection> tauJets;
0035     iEvent.getByToken(*s, tauJets);
0036     for (CaloJetCollection::const_iterator iTau = tauJets->begin(); iTau != tauJets->end(); ++iTau) {
0037       if (iTau->et() > mEt_Min) {
0038         //Add the Pdg Id here
0039         CaloJet myJet = *iTau;
0040         myJet.setPdgId(15);
0041         myTmpJets.push_back(myJet);
0042       }
0043     }
0044   }
0045 
0046   std::unique_ptr<CaloJetCollection> tauL2jets(new CaloJetCollection);
0047 
0048   //Removing the collinear jets correctly!
0049 
0050   //First sort the jets you have merged
0051   SorterByPt sorter;
0052   std::sort(myTmpJets.begin(), myTmpJets.end(), sorter);
0053 
0054   //Remove Collinear Jets by prefering the highest ones!
0055   while (!myTmpJets.empty()) {
0056     tauL2jets->push_back(myTmpJets[0]);
0057     CaloJetCollection tmp;
0058     for (unsigned int i = 1; i < myTmpJets.size(); ++i) {
0059       double DR2 = ROOT::Math::VectorUtil::DeltaR2(myTmpJets[0].p4(), myTmpJets[i].p4());
0060       if (DR2 > 0.1 * 0.1)
0061         tmp.push_back(myTmpJets[i]);
0062     }
0063     myTmpJets.swap(tmp);
0064   }
0065 
0066   iEvent.put(std::move(tauL2jets));
0067 }
0068 
0069 void L2TauJetsMerger::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0070   edm::ParameterSetDescription desc;
0071   std::vector<edm::InputTag> inputTags;
0072   inputTags.push_back(edm::InputTag("hltAkIsoTau1Regional"));
0073   inputTags.push_back(edm::InputTag("hltAkIsoTau2Regional"));
0074   inputTags.push_back(edm::InputTag("hltAkIsoTau3Regional"));
0075   inputTags.push_back(edm::InputTag("hltAkIsoTau4Regional"));
0076   desc.add<std::vector<edm::InputTag> >("JetSrc", inputTags)->setComment("CaloJet collections to merge");
0077   desc.add<double>("EtMin", 20.0)->setComment("Minimal ET of jet to merge");
0078   descriptions.setComment("Merges CaloJet collections removing duplicates");
0079   descriptions.add("L2TauJetsMerger", desc);
0080 }
0081 
0082 #include "FWCore/Framework/interface/MakerMacros.h"
0083 DEFINE_FWK_MODULE(L2TauJetsMerger);