CMS 3D CMS Logo

PFJetsMaxInvMassModule.cc
Go to the documentation of this file.
2 #include "Math/GenVector/VectorUtil.h"
6 
7 //
8 // class declaration
9 //
11  pfJetSrc_ ( consumes<reco::PFJetCollection>(iConfig.getParameter<edm::InputTag>("PFJetSrc") ) ),
12  maxInvMassPairOnly_ ( iConfig.getParameter<bool>("maxInvMassPairOnly") ),
13  removeMaxInvMassPair_ ( iConfig.getParameter<bool>("removeMaxInvMassPair") )
14 {
15  produces<reco::PFJetCollection>();
16 }
18 
20 {
21  std::unique_ptr<reco::PFJetCollection> addPFJets(new reco::PFJetCollection);
22 
24  iEvent.getByToken(pfJetSrc_, jets);
25 
26  unsigned iCan = 0;
27  unsigned jCan = 0;
28  double mjj_max = 0;
29 
30  if(jets->size() > 1){
31  for(unsigned i = 0; i < jets->size(); i++){
32  for(unsigned j = i+1; j < jets->size(); j++){
33  double test = ((*jets)[i].p4()+(*jets)[j].p4()).M();
34  if (test > mjj_max) {
35  mjj_max = test;
36  iCan = i;
37  jCan = j;
38  }
39  }
40  }
41 
42  if (maxInvMassPairOnly_) {
43  addPFJets->push_back((*jets)[iCan]);
44  addPFJets->push_back((*jets)[jCan]);
45  }
46  else if (removeMaxInvMassPair_) {
47  for (unsigned i = 0; i < jets->size(); i++) {
48  if (i != iCan && i != jCan) addPFJets->push_back((*jets)[i]);
49  }
50  }
51  }
52 
53  iEvent.put(std::move(addPFJets));
54 }
55 
57 {
59  desc.add<edm::InputTag>("PFJetSrc", edm::InputTag(""))->setComment("Input collection of PFJets that pass a filter" );
60  desc.add<bool>("maxInvMassPairOnly", true)->setComment("Add only max mjj pair");
61  desc.add<bool>("removeMaxInvMassPair", false)->setComment("Remove max mjj pair and keep all other jets");
62  descriptions.setComment("This module produces a collection of PFJets that are cross-cleaned with respect to PFTaus passing a HLT filter.");
63  descriptions.add ("PFJetsMaxInvMassModule",desc);
64 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
int iEvent
Definition: GenABIO.cc:230
vector< PseudoJet > jets
ParameterDescriptionBase * add(U const &iLabel, T const &value)
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.
def move(src, dest)
Definition: eostools.py:510