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") ),
14  mEt_Min( iConfig.getParameter<double>("EtMin") )
15 {
16  for(vtag::const_iterator it = jetSrc.begin(); it != jetSrc.end(); ++it) {
17  edm::EDGetTokenT<CaloJetCollection> aToken = consumes<CaloJetCollection>(*it);
18  jetSrc_token.push_back(aToken);
19  }
20 
21  produces<CaloJetCollection>();
22 }
23 
25 
27 {
28 
29  using namespace edm;
30  using namespace std;
31  using namespace reco;
32 
33  //Getting the Collections of L2ReconstructedJets from L1Seeds
34  //and removing the collinear jets
35  CaloJetCollection myTmpJets;
36 
37  int iL1Jet = 0;
38  for( vtoken_cjets::const_iterator s = jetSrc_token.begin(); s != jetSrc_token.end(); ++ s ) {
40  iEvent.getByToken( * s, tauJets );
41  for(CaloJetCollection::const_iterator iTau = tauJets->begin(); iTau !=tauJets->end(); ++iTau)
42  {
43  if(iTau->et() > mEt_Min) {
44 
45  //Add the Pdg Id here
46  CaloJet myJet = *iTau;
47  myJet.setPdgId(15);
48  myTmpJets.push_back(myJet);
49  }
50  }
51  iL1Jet++;
52  }
53 
54  std::unique_ptr<CaloJetCollection> tauL2jets(new CaloJetCollection);
55 
56  //Removing the collinear jets correctly!
57 
58  //First sort the jets you have merged
60  std::sort(myTmpJets.begin(),myTmpJets.end(),sorter);
61 
62  //Remove Collinear Jets by prefering the highest ones!
63  while(!myTmpJets.empty()) {
64  tauL2jets->push_back(myTmpJets[0]);
66  for(unsigned int i=1 ;i<myTmpJets.size();++i) {
67  double DR = ROOT::Math::VectorUtil::DeltaR(myTmpJets[0].p4(),myTmpJets[i].p4());
68  if(DR>0.1)
69  tmp.push_back(myTmpJets[i]);
70  }
71  myTmpJets.swap(tmp);
72  }
73 
74 
75  iEvent.put(std::move(tauL2jets));
76 
77 }
78 
80 {
81 
83  std::vector<edm::InputTag> inputTags;
84  inputTags.push_back( edm::InputTag("hltAkIsoTau1Regional") );
85  inputTags.push_back( edm::InputTag("hltAkIsoTau2Regional") );
86  inputTags.push_back( edm::InputTag("hltAkIsoTau3Regional") );
87  inputTags.push_back( edm::InputTag("hltAkIsoTau4Regional") );
88  desc.add< std::vector<edm::InputTag> >("JetSrc",inputTags)->setComment("CaloJet collections to merge");
89  desc.add<double>("EtMin",20.0)->setComment("Minimal ET of jet to merge");
90  descriptions.setComment("Merges CaloJet collections removing duplicates");
91  descriptions.add("L2TauJetsMerger",desc);
92 
93 }
void setComment(std::string const &value)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:127
L2TauJetsMerger(const edm::ParameterSet &)
Jets made from CaloTowers.
Definition: CaloJet.h:29
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
~L2TauJetsMerger() override
int iEvent
Definition: GenABIO.cc:230
double p4[4]
Definition: TauolaWrapper.h:92
vtoken_cjets jetSrc_token
const vtag jetSrc
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void setComment(std::string const &value)
const double mEt_Min
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
fixed size matrix
HLT enums.
std::vector< edm::InputTag > vtag
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
void setPdgId(int pdgId) final
def move(src, dest)
Definition: eostools.py:510
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects