CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTJetCollectionsForLeptonPlusJets.cc
Go to the documentation of this file.
2 
4 
9 
11 
13 
17 
18 
19 
20 
21 template <typename jetType>
23  hltLeptonTag(iConfig.getParameter< edm::InputTag > ("HltLeptonTag")),
24  sourceJetTag(iConfig.getParameter< edm::InputTag > ("SourceJetTag")),
25  minDeltaR_(iConfig.getParameter< double > ("minDeltaR"))
26 {
27  using namespace edm;
28  using namespace std;
29  typedef vector<RefVector<vector<jetType>,jetType,refhelper::FindUsingAdvance<vector<jetType>,jetType> > > JetCollectionVector;
30  m_theLeptonToken = consumes<trigger::TriggerFilterObjectWithRefs>(hltLeptonTag);
31  m_theJetToken = consumes<std::vector<jetType>>(sourceJetTag);
32  produces<JetCollectionVector> ();
33 }
34 
35 template <typename jetType>
37 {
38  // do anything here that needs to be done at desctruction time
39  // (e.g. close files, deallocate resources etc.)
40 
41 }
42 
43 template <typename jetType>
44 void
47  desc.add<edm::InputTag> ("HltLeptonTag", edm::InputTag("triggerFilterObjectWithRefs"));
48  desc.add<edm::InputTag> ("SourceJetTag", edm::InputTag("caloJetCollection"));
49  desc.add<double> ("minDeltaR", 0.5);
51 }
52 
53 //
54 // member functions
55 //
56 
57 
58 // ------------ method called to produce the data ------------
59 // template <typename T>
60 template <typename jetType>
61 void
63 {
64  using namespace edm;
65  using namespace std;
66 
67  typedef vector<RefVector<vector<jetType>,jetType,refhelper::FindUsingAdvance<vector<jetType>,jetType> > > JetCollectionVector;
68  typedef vector<jetType> JetCollection;
71 
73  iEvent.getByToken(m_theLeptonToken,PrevFilterOutput);
74 
75  //its easier on the if statement flow if I try everything at once, shouldnt add to timing
76  // Electrons can be stored as objects of types TriggerCluster, TriggerElectron, or TriggerPhoton
77  vector<Ref<reco::RecoEcalCandidateCollection> > clusCands;
78  PrevFilterOutput->getObjects(trigger::TriggerCluster,clusCands);
79 
80  vector<Ref<reco::ElectronCollection> > eleCands;
81  PrevFilterOutput->getObjects(trigger::TriggerElectron,eleCands);
82 
83  trigger::VRphoton photonCands;
84  PrevFilterOutput->getObjects(trigger::TriggerPhoton, photonCands);
85 
86  vector<reco::RecoChargedCandidateRef> muonCands;
87  PrevFilterOutput->getObjects(trigger::TriggerMuon,muonCands);
88 
89  Handle<JetCollection> theJetCollectionHandle;
90  iEvent.getByToken(m_theJetToken, theJetCollectionHandle);
91 
92  const JetCollection & theJetCollection = *theJetCollectionHandle;
93 
94  auto_ptr < JetCollectionVector > allSelections(new JetCollectionVector());
95 
96  if(!clusCands.empty()){ // try trigger clusters
97  for(size_t candNr=0;candNr<clusCands.size();candNr++){
98  JetRefVector refVector;
99  for (unsigned int j = 0; j < theJetCollection.size(); j++) {
100  if (deltaR(clusCands[candNr]->superCluster()->position(),theJetCollection[j]) > minDeltaR_) refVector.push_back(JetRef(theJetCollectionHandle, j));
101  }
102  allSelections->push_back(refVector);
103  }
104  }
105 
106  if(!eleCands.empty()){ // try electrons
107  for(size_t candNr=0;candNr<eleCands.size();candNr++){
108  JetRefVector refVector;
109  for (unsigned int j = 0; j < theJetCollection.size(); j++) {
110  if (deltaR(eleCands[candNr]->superCluster()->position(),theJetCollection[j]) > minDeltaR_) refVector.push_back(JetRef(theJetCollectionHandle, j));
111  }
112  allSelections->push_back(refVector);
113  }
114  }
115 
116  if(!photonCands.empty()){ // try photons
117  for(size_t candNr=0;candNr<photonCands.size();candNr++){
118  JetRefVector refVector;
119  for (unsigned int j = 0; j < theJetCollection.size(); j++) {
120  if (deltaR(photonCands[candNr]->superCluster()->position(),theJetCollection[j]) > minDeltaR_) refVector.push_back(JetRef(theJetCollectionHandle, j));
121  }
122  allSelections->push_back(refVector);
123  }
124  }
125 
126  if(!muonCands.empty()){ // muons
127  for(size_t candNr=0;candNr<muonCands.size();candNr++){
128  JetRefVector refVector;
129  for (unsigned int j = 0; j < theJetCollection.size(); j++) {
130  if (deltaR(muonCands[candNr]->p4(),theJetCollection[j]) > minDeltaR_) refVector.push_back(JetRef(theJetCollectionHandle, j));
131  }
132  allSelections->push_back(refVector);
133  }
134  }
135 
136 
137 
138 
139  iEvent.put(allSelections);
140 
141  return;
142 
143 }
144 
std::string defaultModuleLabel()
std::vector< Jet > JetCollection
Definition: Jet.h:52
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > m_theLeptonToken
edm::Ref< JetCollection > JetRef
Definition: Jet.h:54
edm::RefVector< JetCollection > JetRefVector
Definition: Jet.h:55
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
double p4[4]
Definition: TauolaWrapper.h:92
int j
Definition: DBlmapReader.cc:9
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
static int position[264][3]
Definition: ReadPGInfo.cc:509
edm::EDGetTokenT< std::vector< jetType > > m_theJetToken
virtual void produce(edm::Event &, const edm::EventSetup &)
std::vector< reco::RecoEcalCandidateRef > VRphoton