CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
L2TauJetsProvider.cc
Go to the documentation of this file.
2 #include "Math/GenVector/VectorUtil.h"
7 //
8 // class decleration
9 //
10 using namespace reco;
11 using namespace std;
12 using namespace edm;
13 using namespace l1extra;
14 
16 {
17  jetSrc = iConfig.getParameter<vtag>("JetSrc");
18  for( vtag::const_iterator s = jetSrc.begin(); s != jetSrc.end(); ++ s ) {
19  edm::EDGetTokenT<CaloJetCollection> aToken = consumes<CaloJetCollection>( *s );
20  jetSrcToken.push_back(aToken);
21  }
22  l1ParticlesTau = consumes<L1JetParticleCollection>(iConfig.getParameter<InputTag>("L1ParticlesTau"));
23  l1ParticlesJet = consumes<L1JetParticleCollection>(iConfig.getParameter<InputTag>("L1ParticlesJet"));
24  tauTrigger = consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<InputTag>("L1TauTrigger"));
25  mEt_Min = iConfig.getParameter<double>("EtMin");
26 
27  produces<CaloJetCollection>();
28 }
29 
31 
33 {
34 
35  using namespace edm;
36  using namespace std;
37  using namespace reco;
38  using namespace trigger;
39  using namespace l1extra;
40 
41 
42  //Getting all the L1Seeds
43 
44 
45  //Getting the Collections of L2ReconstructedJets from L1Seeds
46  //and removing the collinear jets
47  myL2L1JetsMap.clear();
48  int iL1Jet = 0;
49  typedef std::vector<edm::EDGetTokenT<reco::CaloJetCollection> > vtag_token;
50  for( vtag_token::const_iterator s = jetSrcToken.begin(); s != jetSrcToken.end(); ++ s ) {
52  iEvent.getByToken( * s, tauJets );
53  CaloJetCollection::const_iterator iTau = tauJets->begin();
54  if(iTau != tauJets->end()){
55  //Create a Map to associate to every Jet its L1SeedId, i.e. 0,1,2 or 3
56  myL2L1JetsMap.insert(std::pair<int, const CaloJet>(iL1Jet, *(iTau)));
57  }
58  iL1Jet++;
59  }
60  std::auto_ptr<CaloJetCollection> tauL2jets(new CaloJetCollection);
61 
62  //Loop over the jetSrc to select the proper jets
63 
64  double deltaR = 0.1;
65  double matchingR = 0.01;
66  //Loop over the Map to find which jets has fired the trigger
67  //myL1Tau is the Collection of L1TauCandidates (from 0 to max 4 elements)
68  //get the list of trigger candidates from the HLTL1SeedGT filter
71 
72  iEvent.getByToken(l1ParticlesTau, tauColl);
73  iEvent.getByToken(l1ParticlesJet, jetColl);
74 
75  const L1JetParticleCollection myL1Tau = *(tauColl.product());
76 
77  const L1JetParticleCollection myL1Jet = *(jetColl.product());
78 
80  myL1Obj.reserve(8);
81 
82  for(unsigned int i=0;i<myL1Tau.size();i++)
83  {
84  myL1Obj.push_back(myL1Tau[i]);
85  }
86  for(unsigned int j=0;j<myL1Jet.size();j++)
87  {
88  myL1Obj.push_back(myL1Jet[j]);
89  }
90 
91 
93  if(iEvent.getByToken(tauTrigger,l1TriggeredTaus)){
94 
95 
96  tauCandRefVec.clear();
97  jetCandRefVec.clear();
98 
99  l1TriggeredTaus->getObjects( trigger::TriggerL1TauJet,tauCandRefVec);
100  l1TriggeredTaus->getObjects( trigger::TriggerL1CenJet,jetCandRefVec);
101 
102  for( unsigned int iL1Tau=0; iL1Tau <tauCandRefVec.size();iL1Tau++)
103  {
104  for(unsigned int iJet=0;iJet<myL1Obj.size();iJet++)
105  {
106  //Find the relative L2TauJets, to see if it has been reconstructed
107  std::map<int, const reco::CaloJet>::const_iterator myL2itr = myL2L1JetsMap.find(iJet);
108  if(myL2itr!=myL2L1JetsMap.end()){
109  //Calculate the DeltaR between L1TauCandidate and L1Tau which fired the trigger
110  if(&tauCandRefVec[iL1Tau])
111  deltaR = ROOT::Math::VectorUtil::DeltaR(myL1Obj[iJet].p4().Vect(), (tauCandRefVec[iL1Tau]->p4()).Vect());
112  if(deltaR < matchingR ) {
113  // Getting back from the map the L2TauJet
114  const CaloJet myL2TauJet = myL2itr->second;
115  if(myL2TauJet.pt() > mEt_Min) tauL2jets->push_back(myL2TauJet);
116  myL2L1JetsMap.erase(myL2itr->first);
117  break;
118 
119  }
120  }
121 
122  }
123  }
124 
125  for(unsigned int iL1Tau=0; iL1Tau <jetCandRefVec.size();iL1Tau++)
126  {
127  for(unsigned int iJet=0;iJet<myL1Obj.size();iJet++)
128  {
129  //Find the relative L2TauJets, to see if it has been reconstructed
130  std::map<int, const reco::CaloJet>::const_iterator myL2itr = myL2L1JetsMap.find(iJet);
131  if(myL2itr!=myL2L1JetsMap.end()){
132  //Calculate the DeltaR between L1TauCandidate and L1Tau which fired the trigger
133  if(&jetCandRefVec[iL1Tau])
134  deltaR = ROOT::Math::VectorUtil::DeltaR(myL1Obj[iJet].p4().Vect(), (jetCandRefVec[iL1Tau]->p4()).Vect());
135  if(deltaR < matchingR) {
136  // Getting back from the map the L2TauJet
137  const CaloJet myL2TauJet = myL2itr->second;
138 
139  if(myL2TauJet.pt() > mEt_Min) tauL2jets->push_back(myL2TauJet);
140  myL2L1JetsMap.erase(myL2itr->first);
141  break;
142 
143  }
144  }
145 
146  }
147  }
148 
149  }
150  // std::cout <<"Size of L2 jets "<<tauL2jets->size()<<std::endl;
151 
152  iEvent.put(tauL2jets);
153 
154 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
Jets made from CaloTowers.
Definition: CaloJet.h:29
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
double deltaR(const T1 &t1, const T2 &t2)
Definition: deltaR.h:48
virtual double pt() const
transverse momentum
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
double p4[4]
Definition: TauolaWrapper.h:92
int j
Definition: DBlmapReader.cc:9
L2TauJetsProvider(const edm::ParameterSet &)
virtual void produce(edm::Event &, const edm::EventSetup &) override
std::vector< edm::InputTag > vtag
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects