CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AlCaDiJetsProducer.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 #include <string>
4 // user include files
13 
31 
34 #include<iostream>
35 
36 //
37 // class declaration
38 //
39 
41  public:
42  explicit AlCaDiJetsProducer(const edm::ParameterSet&);
44  virtual void beginJob() ;
45  virtual void produce(edm::Event &, const edm::EventSetup&);
46  virtual void endJob();
47  private:
49 
50  // ----------member data ---------------------------
51 
53  double minPtJet_;
55 
60  //edm::EDGetTokenT<edm::TriggerResults> tok_TrigRes_;
63 };
64 
65 AlCaDiJetsProducer::AlCaDiJetsProducer(const edm::ParameterSet& iConfig) : nAll_(0), nSelect_(0) {
66  // Take input
67  labelPFJet_ = iConfig.getParameter<edm::InputTag>("PFjetInput");
68  labelHBHE_ = iConfig.getParameter<edm::InputTag>("HBHEInput");
69  labelHF_ = iConfig.getParameter<edm::InputTag>("HFInput");
70  labelHO_ = iConfig.getParameter<edm::InputTag>("HOInput");
71  //labelTrigger_ = iConfig.getParameter<edm::InputTag>("TriggerResults");
72  labelPFCandidate_= iConfig.getParameter<edm::InputTag>("particleFlowInput");
73  labelVertex_ = iConfig.getParameter<edm::InputTag>("VertexInput");
74  minPtJet_ = iConfig.getParameter<double>("MinPtJet");
75 
76  tok_PFJet_ = consumes<reco::PFJetCollection>(labelPFJet_);
77  tok_HBHE_ = consumes<edm::SortedCollection<HBHERecHit,edm::StrictWeakOrdering<HBHERecHit>>>(labelHBHE_);
78  tok_HF_ = consumes<edm::SortedCollection<HFRecHit,edm::StrictWeakOrdering<HFRecHit>>>(labelHF_);
79  tok_HO_ = consumes<edm::SortedCollection<HORecHit,edm::StrictWeakOrdering<HORecHit>>>(labelHO_);
80  //tok_TrigRes_= consumes<edm::TriggerResults>(labelTrigger_);
81  tok_PFCand_ = consumes<reco::PFCandidateCollection>(labelPFCandidate_);
82  tok_Vertex_ = consumes<reco::VertexCollection>(labelVertex_);
83 
84  // register your products
85  produces<reco::PFJetCollection>(labelPFJet_.encode());
86  produces<edm::SortedCollection<HBHERecHit,edm::StrictWeakOrdering<HBHERecHit>>>(labelHBHE_.encode());
87  produces<edm::SortedCollection<HFRecHit,edm::StrictWeakOrdering<HFRecHit>>>(labelHF_.encode());
88  produces<edm::SortedCollection<HORecHit,edm::StrictWeakOrdering<HORecHit>>>(labelHO_.encode());
89  //produces<edm::TriggerResults>(labelTrigger_.encode());
90  produces<reco::PFCandidateCollection>(labelPFCandidate_.encode());
91  produces<reco::VertexCollection>(labelVertex_.encode());
92 }
93 
95 
97 
99  edm::LogInfo("AlcaDiJets") << "Accepts " << nSelect_ << " events from a total of " << nAll_ << " events";
100 }
101 
103  if (jt.size()<2) return false;
104  if (((jt.at(0)).pt())<minPtJet_) return false;
105  return true;
106 }
107 // ------------ method called to produce the data ------------
109  nAll_++;
110 
111  // Access the collections from iEvent
113  iEvent.getByToken(tok_PFJet_,pfjet);
114  if (!pfjet.isValid()) {
115  edm::LogWarning("AlCaDiJets") << "AlCaDiJetsProducer: Error! can't get product " << labelPFJet_;
116  return ;
117  }
118  const reco::PFJetCollection pfjets = *(pfjet.product());
119 
121  iEvent.getByToken(tok_PFCand_,pfc);
122  if (!pfc.isValid()) {
123  edm::LogWarning("AlCaDiJets") << "AlCaDiJetsProducer: Error! can't get product " << labelPFCandidate_;
124  return ;
125  }
126  const reco::PFCandidateCollection pfcand = *(pfc.product());
127 
129  iEvent.getByToken(tok_Vertex_,vt);
130  if (!vt.isValid()) {
131  edm::LogWarning("AlCaDiJets") << "AlCaDiJetsProducer: Error! can't get product " << labelVertex_;
132  return ;
133  }
134  const reco::VertexCollection vtx = *(vt.product());
135 
137  iEvent.getByToken(tok_HBHE_,hbhe);
138  if (!hbhe.isValid()) {
139  edm::LogWarning("AlCaDiJets") << "AlCaDiJetsProducer: Error! can't get product " << labelHBHE_;
140  return ;
141  }
143 
145  iEvent.getByToken(tok_HO_,ho);
146  if(!ho.isValid()) {
147  edm::LogWarning("AlCaDiJets") << "AlCaDiJetsProducer: Error! can't get product " << labelHO_;
148  return ;
149  }
151 
153  iEvent.getByToken(tok_HF_,hf);
154  if(!hf.isValid()) {
155  edm::LogWarning("AlCaDiJets") << "AlCaDiJetsProducer: Error! can't get product " << labelHF_;
156  return ;
157  }
159 
160  // See if this event is useful
161  bool accept = select(pfjets);
162  if (accept) {
163  nSelect_++;
164 
165  //Copy from standard place
166  std::auto_ptr<reco::PFJetCollection> miniPFjetCollection(new reco::PFJetCollection);
167  for(reco::PFJetCollection::const_iterator pfjetItr=pfjets.begin();
168  pfjetItr!=pfjets.end(); pfjetItr++) {
169  miniPFjetCollection->push_back(*pfjetItr);
170  }
171 
172  std::auto_ptr<reco::PFCandidateCollection> miniPFCandCollection(new reco::PFCandidateCollection);
173  for(reco::PFCandidateCollection::const_iterator pfcItr=pfcand.begin();
174  pfcItr!=pfcand.end(); pfcItr++) {
175  miniPFCandCollection->push_back(*pfcItr);
176  }
177 
178  std::auto_ptr<reco::VertexCollection> miniVtxCollection(new reco::VertexCollection);
179  for(reco::VertexCollection::const_iterator vtxItr=vtx.begin();
180  vtxItr!=vtx.end(); vtxItr++) {
181  miniVtxCollection->push_back(*vtxItr);
182  }
183 
184  std::auto_ptr<edm::SortedCollection<HBHERecHit,edm::StrictWeakOrdering<HBHERecHit>>> miniHBHECollection(new edm::SortedCollection<HBHERecHit,edm::StrictWeakOrdering<HBHERecHit>>);
185  for(edm::SortedCollection<HBHERecHit,edm::StrictWeakOrdering<HBHERecHit> >::const_iterator hbheItr=Hithbhe.begin();
186  hbheItr!=Hithbhe.end(); hbheItr++) {
187  miniHBHECollection->push_back(*hbheItr);
188  }
189 
190  std::auto_ptr<edm::SortedCollection<HORecHit,edm::StrictWeakOrdering<HORecHit>>> miniHOCollection(new edm::SortedCollection<HORecHit,edm::StrictWeakOrdering<HORecHit>>);
191  for(edm::SortedCollection<HORecHit,edm::StrictWeakOrdering<HORecHit> >::const_iterator hoItr=Hitho.begin();
192  hoItr!=Hitho.end(); hoItr++) {
193  miniHOCollection->push_back(*hoItr);
194  }
195 
196  std::auto_ptr<edm::SortedCollection<HFRecHit,edm::StrictWeakOrdering<HFRecHit>>> miniHFCollection(new edm::SortedCollection<HFRecHit,edm::StrictWeakOrdering<HFRecHit>>);
197  for(edm::SortedCollection<HFRecHit,edm::StrictWeakOrdering<HFRecHit> >::const_iterator hfItr=Hithf.begin();
198  hfItr!=Hithf.end(); hfItr++) {
199  miniHFCollection->push_back(*hfItr);
200  }
201 
202  //Put them in the event
203  iEvent.put( miniPFjetCollection, labelPFJet_.encode());
204  iEvent.put( miniHBHECollection, labelHBHE_.encode());
205  iEvent.put( miniHFCollection, labelHF_.encode());
206  iEvent.put( miniHOCollection, labelHO_.encode());
207  //iEvent.put( miniTriggerCollection, labelTrigger_.encode());
208  iEvent.put( miniPFCandCollection, labelPFCandidate_.encode());
209  iEvent.put( miniVtxCollection, labelVertex_.encode());
210  }
211  return;
212 
213 }
214 
T getParameter(std::string const &) const
edm::EDGetTokenT< edm::SortedCollection< HBHERecHit, edm::StrictWeakOrdering< HBHERecHit > > > tok_HBHE_
bool select(reco::PFJetCollection)
edm::InputTag labelPFCandidate_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< edm::SortedCollection< HORecHit, edm::StrictWeakOrdering< HORecHit > > > tok_HO_
edm::InputTag labelHBHE_
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< reco::PFJetCollection > tok_PFJet_
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:25
edm::EDGetTokenT< edm::SortedCollection< HFRecHit, edm::StrictWeakOrdering< HFRecHit > > > tok_HF_
std::string encode() const
Definition: InputTag.cc:164
return((rh^lh)&mask)
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
edm::InputTag labelPFJet_
edm::EDGetTokenT< reco::VertexCollection > tok_Vertex_
bool isValid() const
Definition: HandleBase.h:75
edm::EDGetTokenT< reco::PFCandidateCollection > tok_PFCand_
const_iterator end() const
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
T const * product() const
Definition: Handle.h:81
std::vector< PFJet > PFJetCollection
collection of PFJet objects
AlCaDiJetsProducer(const edm::ParameterSet &)
virtual void beginJob()
edm::InputTag labelVertex_
const_iterator begin() const
virtual void produce(edm::Event &, const edm::EventSetup &)