CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DQMHcalDiJetsAlCaReco.cc
Go to the documentation of this file.
1 /*
2  * \file DQMHcalPhiSymAlCaReco.cc
3  *
4  * \author Olga Kodolova
5  *
6  * $Date: 2010/04/23 16:35:17 $
7  * $Revision: 1.7 $
8  *
9  *
10  * Description: Monitoring of Phi Symmetry Calibration Stream
11 */
12 
17 
18 // DQM include files
19 
21 
22 // work on collections
23 
24 
34 
35 
36 // #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
37 // #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
38 // #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
39 
42 
43 using namespace std;
44 using namespace edm;
45 using namespace reco;
46 
47 // ******************************************
48 // constructors
49 // *****************************************
50 
52 eventCounter_(0)
53 {
55 //
56 // Input from configurator file
57 //
58  folderName_ = iConfig.getUntrackedParameter<string>("FolderName","ALCAStreamHcalDiJets");
59 
60 
61  jets_= iConfig.getParameter<edm::InputTag>("jetsInput");
62  ec_= iConfig.getParameter<edm::InputTag>("ecInput");
63  hbhe_= iConfig.getParameter<edm::InputTag>("hbheInput");
64  ho_= iConfig.getParameter<edm::InputTag>("hoInput");
65  hf_= iConfig.getParameter<edm::InputTag>("hfInput");
66 
67  saveToFile_= iConfig.getUntrackedParameter<bool>("SaveToFile",false);
68  fileName_= iConfig.getUntrackedParameter<string>("FileName","MonitorAlCaHcalDiJets.root");
69 
70 }
71 
73 {}
74 
75 //--------------------------------------------------------
77 
78 
79  // create and cd into new folder
81 
82  // book some histograms 1D
83  hiDistrRecHitEnergyEBEE_ = dbe_->book1D("RecHitEnergyEBEE", "the number of hits inside jets", 100, 0, 800);
85  hiDistrRecHitEnergyEBEE_->setAxisTitle("# rechits", 2);
86 
87  hiDistrRecHitEnergyHBHE_ = dbe_->book1D("RecHitEnergyHBHE", "the number of hits inside jets", 100,0, 800);
89  hiDistrRecHitEnergyHBHE_->setAxisTitle("# rechits", 2);
90 
91  hiDistrRecHitEnergyHF_ = dbe_->book1D("RecHitEnergyHF", "the number of hits inside jets", 150,0, 1500);
93  hiDistrRecHitEnergyHF_->setAxisTitle("# rechits", 2);
94 
95  hiDistrRecHitEnergyHO_ = dbe_->book1D("RecHitEnergyHO", "the number of hits inside jets", 100,0, 100);
97  hiDistrRecHitEnergyHO_->setAxisTitle("# rechits", 2);
98 
99  hiDistrProbeJetEnergy_ = dbe_->book1D("ProbeJetEnergy", "the energy of probe jets", 250,0, 2500);
100  hiDistrProbeJetEnergy_->setAxisTitle("E, GeV", 1);
101  hiDistrProbeJetEnergy_->setAxisTitle("# jets", 2);
102 
103  hiDistrProbeJetEta_ = dbe_->book1D("ProbeJetEta", "the number of probe jets", 100, -5., 5.);
104  hiDistrProbeJetEta_->setAxisTitle("#eta", 1);
105  hiDistrProbeJetEta_->setAxisTitle("# jets", 2);
106 
107  hiDistrProbeJetPhi_ = dbe_->book1D("ProbeJetPhi", "the number of probe jets", 50, -3.14, 3.14);
108  hiDistrProbeJetPhi_->setAxisTitle("#phi", 1);
109  hiDistrProbeJetPhi_->setAxisTitle("# jets", 2);
110 
111  hiDistrTagJetEnergy_ = dbe_->book1D("TagJetEnergy", "the energy of tsg jets", 250,0, 2500);
112  hiDistrTagJetEnergy_->setAxisTitle("E, GeV", 1);
113  hiDistrTagJetEnergy_->setAxisTitle("# jets", 2);
114 
115  hiDistrTagJetEta_ = dbe_->book1D("TagJetEta", "the number of tag jets", 100, -5., 5.);
116  hiDistrTagJetEta_->setAxisTitle("#eta", 1);
117  hiDistrTagJetEta_->setAxisTitle("# jets", 2);
118 
119  hiDistrTagJetPhi_ = dbe_->book1D("TagJetPhi", "the number of tag jets", 50, -3.14, 3.14);
120  hiDistrTagJetPhi_->setAxisTitle("#phi", 1);
121  hiDistrTagJetPhi_->setAxisTitle("# jets", 2);
122 
123  hiDistrEtThirdJet_ = dbe_->book1D("EtThirdJet", "Et of the third jet", 90, 0, 90);
124 
125 
126 //==================================================================================
127 
128 
129 }
130 
131 //--------------------------------------------------------
132 void DQMHcalDiJetsAlCaReco::beginRun(const edm::Run& r, const EventSetup& context) {
133 
134 }
135 
136 //--------------------------------------------------------
138  const EventSetup& context) {
139 
140 }
141 
142 //-------------------------------------------------------------
143 
145  const EventSetup& iSetup ){
146 
147 
148  eventCounter_++;
149 
150  CaloJet jet1, jet2, jet3;
151  Float_t etVetoJet;
152 
154  iEvent.getByLabel(jets_,jets);
155 
156  if(!jets.isValid()){
157  LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't getjet product!" << std::endl;
158  return ;
159  }
160 
161 
162 
163  if(jets->size()>1){
164  jet1 = (*jets)[0];
165  jet2 = (*jets)[1];
166  if(fabs(jet1.eta())>fabs(jet2.eta())){
167  CaloJet jet = jet1;
168  jet1 = jet2;
169  jet2 = jet;
170  }
171  // if(fabs(jet1.eta())>eta_1 || (fabs(jet2.eta())-jet_R) < eta_2){ return;}
172  } else {return;}
173 
175  hiDistrTagJetEta_->Fill(jet1.eta());
176  hiDistrTagJetPhi_->Fill(jet1.phi());
177 
179  hiDistrProbeJetEta_->Fill(jet2.eta());
180  hiDistrProbeJetPhi_->Fill(jet2.phi());
181 
182  if(jets->size()>2){
183  jet3 = (*jets)[2];
184  etVetoJet = jet3.et();
185  hiDistrEtThirdJet_->Fill(etVetoJet);
186  } else { etVetoJet = 0.; hiDistrEtThirdJet_->Fill(etVetoJet); }
187 
188 
189 
191  iEvent.getByLabel(ec_,ec);
192 
193  if(!ec.isValid()){
194  LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get ec product!" << std::endl;
195  return ;
196  }
197 
198 
199  for(EcalRecHitCollection::const_iterator ecItr = (*ec).begin();
200  ecItr != (*ec).end(); ++ecItr)
201  {
202  hiDistrRecHitEnergyEBEE_->Fill(ecItr->energy());
203  }
204 
205 
206 
207 
209  iEvent.getByLabel(hbhe_, hbhe);
210 
211  if(!hbhe.isValid()){
212  LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get hbhe product!" << std::endl;
213  return ;
214  }
215 
216 
217  for(HBHERecHitCollection::const_iterator hbheItr=hbhe->begin();
218  hbheItr!=hbhe->end(); hbheItr++)
219  {
220  hiDistrRecHitEnergyHBHE_->Fill(hbheItr->energy());
221  }
222 
223 
225  iEvent.getByLabel(ho_, ho);
226 
227  if(!ho.isValid()){
228  LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get ho product!" << std::endl;
229  return ;
230  }
231 
232 
233  for(HORecHitCollection::const_iterator hoItr=ho->begin();
234  hoItr!=ho->end(); hoItr++)
235  {
236  hiDistrRecHitEnergyHO_->Fill(hoItr->energy());
237 
238  }
239 
240 
241 
242 
244  iEvent.getByLabel(hf_, hf);
245 
246  if(!hf.isValid()){
247  LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get hf product!" << std::endl;
248  return ;
249  }
250 
251 
252  for(HFRecHitCollection::const_iterator hfItr=hf->begin();
253  hfItr!=hf->end(); hfItr++)
254  {
255  hiDistrRecHitEnergyHF_->Fill(hfItr->energy());
256  }
257 
258 } //analyze
259 
260 
261 
262 
263 //--------------------------------------------------------
265  const EventSetup& context) {
266 }
267 //--------------------------------------------------------
268 void DQMHcalDiJetsAlCaReco::endRun(const Run& r, const EventSetup& context){
269 
270 }
271 //--------------------------------------------------------
273 
274  if (saveToFile_) {
275  dbe_->save(fileName_);
276  }
277 
278 }
279 
280 
#define LogDebug(id)
MonitorElement * hiDistrTagJetEnergy_
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
void beginRun(const edm::Run &r, const edm::EventSetup &c)
MonitorElement * hiDistrRecHitEnergyHF_
Jets made from CaloTowers.
Definition: CaloJet.h:30
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:717
virtual double et() const
transverse energy
std::string fileName_
Output file name if required.
MonitorElement * hiDistrTagJetPhi_
MonitorElement * hiDistrTagJetEta_
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE")
Definition: DQMStore.cc:2113
std::vector< EcalRecHit >::const_iterator const_iterator
MonitorElement * hiDistrRecHitEnergyEBEE_
virtual double eta() const
momentum pseudorapidity
void Fill(long long x)
virtual double energy() const
energy
int iEvent
Definition: GenABIO.cc:243
void beginLuminosityBlock(const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &context)
vector< PseudoJet > jets
MonitorElement * hiDistrProbeJetEta_
void endLuminosityBlock(const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &c)
bool isValid() const
Definition: HandleBase.h:76
DQMHcalDiJetsAlCaReco(const edm::ParameterSet &)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
edm::InputTag jets_
object to monitor
MonitorElement * hiDistrProbeJetPhi_
bool saveToFile_
Write to file.
MonitorElement * hiDistrRecHitEnergyHBHE_
MonitorElement * hiDistrProbeJetEnergy_
MonitorElement * hiDistrRecHitEnergyHO_
std::string folderName_
DQM folder name.
virtual double phi() const
momentum azimuthal angle
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
MonitorElement * hiDistrEtThirdJet_
void analyze(const edm::Event &e, const edm::EventSetup &c)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:429
Definition: Run.h:33
void endRun(const edm::Run &r, const edm::EventSetup &c)