File indexing completed on 2024-04-06 12:27:44
0001 #include "RecoTauTag/HLTProducers/interface/PFJetsMaxInvMassModule.h"
0002 #include "Math/GenVector/VectorUtil.h"
0003 #include "CommonTools/Utils/interface/PtComparator.h"
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0006 #include "FWCore/Utilities/interface/InputTag.h"
0007 #include "DataFormats/Common/interface/Handle.h"
0008
0009
0010
0011
0012 PFJetsMaxInvMassModule::PFJetsMaxInvMassModule(const edm::ParameterSet& iConfig)
0013 : pfJetSrc_(consumes<reco::PFJetCollection>(iConfig.getParameter<edm::InputTag>("PFJetSrc"))),
0014 maxInvMassPairOnly_(iConfig.getParameter<bool>("maxInvMassPairOnly")),
0015 removeMaxInvMassPair_(iConfig.getParameter<bool>("removeMaxInvMassPair")) {
0016 produces<reco::PFJetCollection>();
0017 }
0018
0019 void PFJetsMaxInvMassModule::produce(edm::StreamID iSId, edm::Event& iEvent, const edm::EventSetup& iES) const {
0020 std::unique_ptr<reco::PFJetCollection> addPFJets(new reco::PFJetCollection);
0021
0022 edm::Handle<reco::PFJetCollection> jets;
0023 iEvent.getByToken(pfJetSrc_, jets);
0024
0025 unsigned iCan = 0;
0026 unsigned jCan = 0;
0027 double m2jj_max = 0;
0028
0029 if (jets->size() > 1) {
0030 for (unsigned i = 0; i < jets->size(); i++) {
0031 for (unsigned j = i + 1; j < jets->size(); j++) {
0032 double test = ((*jets)[i].p4() + (*jets)[j].p4()).M2();
0033 if (test > m2jj_max) {
0034 m2jj_max = test;
0035 iCan = i;
0036 jCan = j;
0037 }
0038 }
0039 }
0040
0041 if (maxInvMassPairOnly_) {
0042 addPFJets->push_back((*jets)[iCan]);
0043 addPFJets->push_back((*jets)[jCan]);
0044 } else if (removeMaxInvMassPair_) {
0045 for (unsigned i = 0; i < jets->size(); i++) {
0046 if (i != iCan && i != jCan)
0047 addPFJets->push_back((*jets)[i]);
0048 }
0049 }
0050 }
0051
0052 iEvent.put(std::move(addPFJets));
0053 }
0054
0055 void PFJetsMaxInvMassModule::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0056 edm::ParameterSetDescription desc;
0057 desc.add<edm::InputTag>("PFJetSrc", edm::InputTag(""))->setComment("Input collection of PFJets that pass a filter");
0058 desc.add<bool>("maxInvMassPairOnly", true)->setComment("Add only max mjj pair");
0059 desc.add<bool>("removeMaxInvMassPair", false)->setComment("Remove max mjj pair and keep all other jets");
0060 descriptions.setComment(
0061 "This module produces a collection of PFJets that are cross-cleaned with respect to PFTaus passing a HLT "
0062 "filter.");
0063 descriptions.add("PFJetsMaxInvMassModule", desc);
0064 }
0065
0066 #include "FWCore/Framework/interface/MakerMacros.h"
0067 DEFINE_FWK_MODULE(PFJetsMaxInvMassModule);