CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTJetCollectionsForBoostedLeptonPlusJets.cc
Go to the documentation of this file.
1 #include <string>
2 #include <vector>
3 
4 #include "TVector3.h"
5 
20 
21 
22 
23 
24 template <typename jetType>
26  hltLeptonTag(iConfig.getParameter< edm::InputTag > ("HltLeptonTag")),
27  sourceJetTag(iConfig.getParameter< edm::InputTag > ("SourceJetTag")),
28  minDeltaR_(iConfig.getParameter< double > ("minDeltaR"))
29 {
30  using namespace edm;
31  using namespace std;
32  typedef vector<RefVector<vector<jetType>,jetType,refhelper::FindUsingAdvance<vector<jetType>,jetType> > > JetCollectionVector;
33  m_theLeptonToken = consumes<trigger::TriggerFilterObjectWithRefs>(hltLeptonTag);
34  m_theJetToken = consumes<std::vector<jetType>>(sourceJetTag);
35  produces<JetCollectionVector> ();
36 }
37 
38 template <typename jetType>
40 {
41  // do anything here that needs to be done at desctruction time
42  // (e.g. close files, deallocate resources etc.)
43 
44 }
45 
46 template <typename jetType>
47 void
50  desc.add<edm::InputTag> ("HltLeptonTag", edm::InputTag("triggerFilterObjectWithRefs"));
51  //(2)
52  desc.add<edm::InputTag> ("SourceJetTag", edm::InputTag("jetCollection"));
53  //(2)
54  desc.add<double> ("minDeltaR", 0.5);
56 }
57 
58 //
59 // member functions
60 //
61 
62 
63 // ------------ method called to produce the data ------------
64 // template <typename T>
65 template <typename jetType>
66 void
68 {
69  using namespace edm;
70  using namespace std;
71  //(3)
72  using namespace reco;
73  //(3)
74 
75  typedef vector<RefVector<vector<jetType>,jetType,refhelper::FindUsingAdvance<vector<jetType>,jetType> > > JetCollectionVector;
76  typedef vector<jetType> JetCollection;
79  //(4)
81  //(4)
82 
84  iEvent.getByToken(m_theLeptonToken,PrevFilterOutput);
85 
86  //its easier on the if statement flow if I try everything at once, shouldnt add to timing
87  // Electrons can be stored as objects of types TriggerCluster, TriggerElectron, or TriggerPhoton
88  vector<Ref<reco::RecoEcalCandidateCollection> > clusCands;
89  PrevFilterOutput->getObjects(trigger::TriggerCluster,clusCands);
90 
91  vector<Ref<reco::ElectronCollection> > eleCands;
92  PrevFilterOutput->getObjects(trigger::TriggerElectron,eleCands);
93 
94  trigger::VRphoton photonCands;
95  PrevFilterOutput->getObjects(trigger::TriggerPhoton, photonCands);
96 
97  vector<reco::RecoChargedCandidateRef> muonCands;
98  PrevFilterOutput->getObjects(trigger::TriggerMuon,muonCands);
99 
100  Handle<JetCollection> theJetCollectionHandle;
101  iEvent.getByToken(m_theJetToken, theJetCollectionHandle);
102 
103  const JetCollection & theJetCollection = *theJetCollectionHandle;
104 
105  auto_ptr < JetCollectionVector > allSelections(new JetCollectionVector());
106 
107  if(!clusCands.empty()){ // try trigger clusters
108  for(size_t candNr=0;candNr<clusCands.size();candNr++){
109  JetRefVector refVector;
110  for (unsigned int j = 0; j < theJetCollection.size(); j++) {
111  if (deltaR(clusCands[candNr]->superCluster()->position(),theJetCollection[j]) > minDeltaR_) refVector.push_back(JetRef(theJetCollectionHandle, j));
112  else{
113  unsigned int w =0 ;
114  std::vector<reco::PFCandidatePtr> pfConstituents = theJetCollection[j].getPFConstituents();
115  for(std::vector<reco::PFCandidatePtr>::const_iterator i_candidate = pfConstituents.begin(); i_candidate != pfConstituents.end(); ++i_candidate){
116  TVector3 ClusP(clusCands[candNr]->p4().Px(),clusCands[candNr]->p4().Py(), clusCands[candNr]->p4().Pz());
117  TVector3 PFJetConstP((*i_candidate)->px(),(*i_candidate)->py(),(*i_candidate)->pz());
118  double deltaRPFConste = ClusP.DeltaR(PFJetConstP);
119  if(deltaRPFConste < 0.001 && w==0){
120  const_cast<LorentzVector&>(theJetCollection[j].p4()) = theJetCollection[j].p4() - clusCands[candNr]->p4();
121  w ++;
122  } //if
123  }//for constituents
124  refVector.push_back(JetRef(theJetCollectionHandle, j));
125  }//else
126  }
127  allSelections->push_back(refVector);
128  }
129  }
130 
131  if(!eleCands.empty()){ // try electrons
132  for(size_t candNr=0;candNr<eleCands.size();candNr++){
133  JetRefVector refVector;
134  for (unsigned int j = 0; j < theJetCollection.size(); j++) {
135  if (deltaR(eleCands[candNr]->superCluster()->position(),theJetCollection[j]) > minDeltaR_) refVector.push_back(JetRef(theJetCollectionHandle, j));
136  else{
137  unsigned int w =0 ;
138  std::vector<reco::PFCandidatePtr> pfConstituents = theJetCollection[j].getPFConstituents();
139  for(std::vector<reco::PFCandidatePtr>::const_iterator i_candidate = pfConstituents.begin(); i_candidate != pfConstituents.end(); ++i_candidate){
140  TVector3 ElectronP(eleCands[candNr]->p4().Px(),eleCands[candNr]->p4().Py(), eleCands[candNr]->p4().Pz());
141  TVector3 PFJetConstP((*i_candidate)->px(),(*i_candidate)->py(),(*i_candidate)->pz());
142  double deltaRPFConste = ElectronP.DeltaR(PFJetConstP);
143  if(deltaRPFConste < 0.001 && w==0){
144  const_cast<LorentzVector&>(theJetCollection[j].p4()) = theJetCollection[j].p4() - eleCands[candNr]->p4();
145  w ++;
146  } //if
147  }//for constituents
148  refVector.push_back(JetRef(theJetCollectionHandle, j));
149  }//else
150  }//for jet collection
151 
152  allSelections->push_back(refVector);
153  }
154  }
155 
156  if(!photonCands.empty()){ // try photons
157  for(size_t candNr=0;candNr<photonCands.size();candNr++){
158  JetRefVector refVector;
159  for (unsigned int j = 0; j < theJetCollection.size(); j++) {
160  if (deltaR(photonCands[candNr]->superCluster()->position(),theJetCollection[j]) > minDeltaR_) refVector.push_back(JetRef(theJetCollectionHandle, j));
161  else{
162  unsigned int w =0 ;
163  std::vector<reco::PFCandidatePtr> pfConstituents = theJetCollection[j].getPFConstituents();
164  for(std::vector<reco::PFCandidatePtr>::const_iterator i_candidate = pfConstituents.begin(); i_candidate != pfConstituents.end(); ++i_candidate){
165  TVector3 PhotonP(photonCands[candNr]->p4().Px(),photonCands[candNr]->p4().Py(), photonCands[candNr]->p4().Pz());
166  TVector3 PFJetConstP((*i_candidate)->px(),(*i_candidate)->py(),(*i_candidate)->pz());
167  double deltaRPFConste = PhotonP.DeltaR(PFJetConstP);
168  if(deltaRPFConste < 0.001 && w==0){
169  const_cast<LorentzVector&>(theJetCollection[j].p4()) = theJetCollection[j].p4() - photonCands[candNr]->p4();
170  w ++;
171  } //if
172  }//for constituents
173  refVector.push_back(JetRef(theJetCollectionHandle, j));
174  }//else
175  }//for jet collection
176  allSelections->push_back(refVector);
177  }
178  }
179 
180  if(!muonCands.empty()){ // muons
181  for(size_t candNr=0;candNr<muonCands.size();candNr++){
182  JetRefVector refVector;
183  for (unsigned int j = 0; j < theJetCollection.size(); j++) {
184  if (deltaR(muonCands[candNr]->p4(),theJetCollection[j]) > minDeltaR_) refVector.push_back(JetRef(theJetCollectionHandle, j));
185  else{
186  unsigned int w =0 ;
187  std::vector<reco::PFCandidatePtr> pfConstituents = theJetCollection[j].getPFConstituents();
188  for(std::vector<reco::PFCandidatePtr>::const_iterator i_candidate = pfConstituents.begin(); i_candidate != pfConstituents.end(); ++i_candidate){
189  TVector3 MuP(muonCands[candNr]->p4().Px(),muonCands[candNr]->p4().Py(), muonCands[candNr]->p4().Pz());
190  TVector3 PFJetConstP((*i_candidate)->px(),(*i_candidate)->py(),(*i_candidate)->pz());
191  double deltaRPFConste = MuP.DeltaR(PFJetConstP);
192  if(deltaRPFConste < 0.001 && w==0){
193  const_cast<LorentzVector&>(theJetCollection[j].p4()) = theJetCollection[j].p4() - muonCands[candNr]->p4();
194  w ++;
195  } //if
196  }//for constituents
197  refVector.push_back(JetRef(theJetCollectionHandle, j));
198  }//else
199  }
200  allSelections->push_back(refVector);
201  }
202  }
203 
204 
205 
206 
207  iEvent.put(allSelections);
208 
209  return;
210 
211 }
212 
std::string defaultModuleLabel()
const double w
Definition: UKUtility.cc:23
std::vector< Jet > JetCollection
Definition: Jet.h:52
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
virtual void produce(edm::Event &, const edm::EventSetup &)
edm::EDGetTokenT< std::vector< jetType > > m_theJetToken
edm::Ref< JetCollection > JetRef
Definition: Jet.h:54
edm::RefVector< JetCollection > JetRefVector
Definition: Jet.h:55
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
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
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > m_theLeptonToken
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > LorentzVector
Definition: analysisEnums.h:9
static int position[264][3]
Definition: ReadPGInfo.cc:509
std::vector< reco::RecoEcalCandidateRef > VRphoton
math::PtEtaPhiELorentzVectorF LorentzVector