CMS 3D CMS Logo

PFJetsMaxInvMassModule.cc
Go to the documentation of this file.
2 #include "Math/GenVector/VectorUtil.h"
8 
9 //
10 // class declaration
11 //
13  : pfJetSrc_(consumes<reco::PFJetCollection>(iConfig.getParameter<edm::InputTag>("PFJetSrc"))),
14  maxInvMassPairOnly_(iConfig.getParameter<bool>("maxInvMassPairOnly")),
15  removeMaxInvMassPair_(iConfig.getParameter<bool>("removeMaxInvMassPair")) {
16  produces<reco::PFJetCollection>();
17 }
18 
20  std::unique_ptr<reco::PFJetCollection> addPFJets(new reco::PFJetCollection);
21 
23  iEvent.getByToken(pfJetSrc_, jets);
24 
25  unsigned iCan = 0;
26  unsigned jCan = 0;
27  double m2jj_max = 0;
28 
29  if (jets->size() > 1) {
30  for (unsigned i = 0; i < jets->size(); i++) {
31  for (unsigned j = i + 1; j < jets->size(); j++) {
32  double test = ((*jets)[i].p4() + (*jets)[j].p4()).M2();
33  if (test > m2jj_max) {
34  m2jj_max = test;
35  iCan = i;
36  jCan = j;
37  }
38  }
39  }
40 
41  if (maxInvMassPairOnly_) {
42  addPFJets->push_back((*jets)[iCan]);
43  addPFJets->push_back((*jets)[jCan]);
44  } else if (removeMaxInvMassPair_) {
45  for (unsigned i = 0; i < jets->size(); i++) {
46  if (i != iCan && i != jCan)
47  addPFJets->push_back((*jets)[i]);
48  }
49  }
50  }
51 
52  iEvent.put(std::move(addPFJets));
53 }
54 
57  desc.add<edm::InputTag>("PFJetSrc", edm::InputTag(""))->setComment("Input collection of PFJets that pass a filter");
58  desc.add<bool>("maxInvMassPairOnly", true)->setComment("Add only max mjj pair");
59  desc.add<bool>("removeMaxInvMassPair", false)->setComment("Remove max mjj pair and keep all other jets");
60  descriptions.setComment(
61  "This module produces a collection of PFJets that are cross-cleaned with respect to PFTaus passing a HLT "
62  "filter.");
63  descriptions.add("PFJetsMaxInvMassModule", desc);
64 }
65 
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void setComment(std::string const &value)
PFJetsMaxInvMassModule(const edm::ParameterSet &)
const edm::EDGetTokenT< reco::PFJetCollection > pfJetSrc_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< PFJet > PFJetCollection
collection of PFJet objects
fixed size matrix
HLT enums.
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
def move(src, dest)
Definition: eostools.py:511