CMS 3D CMS Logo

L2TauJetsMerger.cc
Go to the documentation of this file.
2 #include "Math/GenVector/VectorUtil.h"
5 //
6 // class decleration
7 //
8 using namespace reco;
9 using namespace std;
10 using namespace edm;
11 
13  : jetSrc(iConfig.getParameter<vtag>("JetSrc")), mEt_Min(iConfig.getParameter<double>("EtMin")) {
14  for (vtag::const_iterator it = jetSrc.begin(); it != jetSrc.end(); ++it) {
15  edm::EDGetTokenT<CaloJetCollection> aToken = consumes<CaloJetCollection>(*it);
16  jetSrc_token.push_back(aToken);
17  }
18 
19  produces<CaloJetCollection>();
20 }
21 
23 
25  using namespace edm;
26  using namespace std;
27  using namespace reco;
28 
29  //Getting the Collections of L2ReconstructedJets from L1Seeds
30  //and removing the collinear jets
31  CaloJetCollection myTmpJets;
32 
33  for (vtoken_cjets::const_iterator s = jetSrc_token.begin(); s != jetSrc_token.end(); ++s) {
35  iEvent.getByToken(*s, tauJets);
36  for (CaloJetCollection::const_iterator iTau = tauJets->begin(); iTau != tauJets->end(); ++iTau) {
37  if (iTau->et() > mEt_Min) {
38  //Add the Pdg Id here
39  CaloJet myJet = *iTau;
40  myJet.setPdgId(15);
41  myTmpJets.push_back(myJet);
42  }
43  }
44  }
45 
46  std::unique_ptr<CaloJetCollection> tauL2jets(new CaloJetCollection);
47 
48  //Removing the collinear jets correctly!
49 
50  //First sort the jets you have merged
52  std::sort(myTmpJets.begin(), myTmpJets.end(), sorter);
53 
54  //Remove Collinear Jets by prefering the highest ones!
55  while (!myTmpJets.empty()) {
56  tauL2jets->push_back(myTmpJets[0]);
58  for (unsigned int i = 1; i < myTmpJets.size(); ++i) {
59  double DR2 = ROOT::Math::VectorUtil::DeltaR2(myTmpJets[0].p4(), myTmpJets[i].p4());
60  if (DR2 > 0.1 * 0.1)
61  tmp.push_back(myTmpJets[i]);
62  }
63  myTmpJets.swap(tmp);
64  }
65 
66  iEvent.put(std::move(tauL2jets));
67 }
68 
71  std::vector<edm::InputTag> inputTags;
72  inputTags.push_back(edm::InputTag("hltAkIsoTau1Regional"));
73  inputTags.push_back(edm::InputTag("hltAkIsoTau2Regional"));
74  inputTags.push_back(edm::InputTag("hltAkIsoTau3Regional"));
75  inputTags.push_back(edm::InputTag("hltAkIsoTau4Regional"));
76  desc.add<std::vector<edm::InputTag> >("JetSrc", inputTags)->setComment("CaloJet collections to merge");
77  desc.add<double>("EtMin", 20.0)->setComment("Minimal ET of jet to merge");
78  descriptions.setComment("Merges CaloJet collections removing duplicates");
79  descriptions.add("L2TauJetsMerger", desc);
80 }
81 
L2TauJetsMerger(const edm::ParameterSet &)
Jets made from CaloTowers.
Definition: CaloJet.h:27
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
~L2TauJetsMerger() override
int iEvent
Definition: GenABIO.cc:224
vtoken_cjets jetSrc_token
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const vtag jetSrc
void setComment(std::string const &value)
inputTags
All input tags are specified in this pset for convenience.
const double mEt_Min
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
std::vector< edm::InputTag > vtag
tmp
align.sh
Definition: createJobs.py:716
void setPdgId(int pdgId) final
def move(src, dest)
Definition: eostools.py:511
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects