CMS 3D CMS Logo

DQMHcalDiJetsAlCaReco.cc
Go to the documentation of this file.
1 /*
2  * \file DQMHcalPhiSymAlCaReco.cc
3  *
4  * \author Olga Kodolova
5  *
6  *
7  *
8  * Description: Monitoring of Phi Symmetry Calibration Stream
9 */
10 
15 
16 // DQM include files
17 
19 
20 // work on collections
21 
26 
27 
28 // #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
29 // #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
30 // #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
31 
33 
34 using namespace std;
35 using namespace edm;
36 using namespace reco;
37 
38 // ******************************************
39 // constructors
40 // *****************************************
41 
43 eventCounter_(0)
44 {
45 //
46 // Input from configurator file
47 //
48  folderName_ = iConfig.getUntrackedParameter<string>("FolderName","ALCAStreamHcalDiJets");
49 
50 
51  jets_= consumes<CaloJetCollection>(iConfig.getParameter<edm::InputTag>("jetsInput"));
52  ec_= consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("ecInput"));
53  hbhe_= consumes<HBHERecHitCollection>(iConfig.getParameter<edm::InputTag>("hbheInput"));
54  ho_= consumes<HORecHitCollection>(iConfig.getParameter<edm::InputTag>("hoInput"));
55  hf_= consumes<HFRecHitCollection>(iConfig.getParameter<edm::InputTag>("hfInput"));
56 
57  saveToFile_= iConfig.getUntrackedParameter<bool>("SaveToFile",false);
58  fileName_= iConfig.getUntrackedParameter<string>("FileName","MonitorAlCaHcalDiJets.root");
59 
60 }
61 
63 {}
64 
65 //--------------------------------------------------------
66 void DQMHcalDiJetsAlCaReco::bookHistograms(DQMStore::IBooker & ibooker, edm::Run const & /* iRun*/, edm::EventSetup const & /* iSetup */)
67 {
68 
69  // create and cd into new folder
71 
72  // book some histograms 1D
73  hiDistrRecHitEnergyEBEE_ = ibooker.book1D("RecHitEnergyEBEE", "the number of hits inside jets", 100, 0, 800);
75  hiDistrRecHitEnergyEBEE_->setAxisTitle("# rechits", 2);
76 
77  hiDistrRecHitEnergyHBHE_ = ibooker.book1D("RecHitEnergyHBHE", "the number of hits inside jets", 100,0, 800);
79  hiDistrRecHitEnergyHBHE_->setAxisTitle("# rechits", 2);
80 
81  hiDistrRecHitEnergyHF_ = ibooker.book1D("RecHitEnergyHF", "the number of hits inside jets", 150,0, 1500);
83  hiDistrRecHitEnergyHF_->setAxisTitle("# rechits", 2);
84 
85  hiDistrRecHitEnergyHO_ = ibooker.book1D("RecHitEnergyHO", "the number of hits inside jets", 100,0, 100);
87  hiDistrRecHitEnergyHO_->setAxisTitle("# rechits", 2);
88 
89  hiDistrProbeJetEnergy_ = ibooker.book1D("ProbeJetEnergy", "the energy of probe jets", 250,0, 2500);
92 
93  hiDistrProbeJetEta_ = ibooker.book1D("ProbeJetEta", "the number of probe jets", 100, -5., 5.);
95  hiDistrProbeJetEta_->setAxisTitle("# jets", 2);
96 
97  hiDistrProbeJetPhi_ = ibooker.book1D("ProbeJetPhi", "the number of probe jets", 50, -3.14, 3.14);
99  hiDistrProbeJetPhi_->setAxisTitle("# jets", 2);
100 
101  hiDistrTagJetEnergy_ = ibooker.book1D("TagJetEnergy", "the energy of tsg jets", 250,0, 2500);
102  hiDistrTagJetEnergy_->setAxisTitle("E, GeV", 1);
103  hiDistrTagJetEnergy_->setAxisTitle("# jets", 2);
104 
105  hiDistrTagJetEta_ = ibooker.book1D("TagJetEta", "the number of tag jets", 100, -5., 5.);
106  hiDistrTagJetEta_->setAxisTitle("#eta", 1);
107  hiDistrTagJetEta_->setAxisTitle("# jets", 2);
108 
109  hiDistrTagJetPhi_ = ibooker.book1D("TagJetPhi", "the number of tag jets", 50, -3.14, 3.14);
110  hiDistrTagJetPhi_->setAxisTitle("#phi", 1);
111  hiDistrTagJetPhi_->setAxisTitle("# jets", 2);
112 
113  hiDistrEtThirdJet_ = ibooker.book1D("EtThirdJet", "Et of the third jet", 90, 0, 90);
114 
115 
116 //==================================================================================
117 
118 
119 }
120 
121 //-------------------------------------------------------------
122 
124  const EventSetup& iSetup ){
125 
126 
127  eventCounter_++;
128 
129  CaloJet jet1, jet2, jet3;
130  Float_t etVetoJet;
131 
133  iEvent.getByToken(jets_,jets);
134 
135  if(!jets.isValid()){
136  LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't getjet product!" << std::endl;
137  return ;
138  }
139 
140 
141 
142  if(jets->size()>1){
143  jet1 = (*jets)[0];
144  jet2 = (*jets)[1];
145  if(fabs(jet1.eta())>fabs(jet2.eta())){
146  CaloJet jet = jet1;
147  jet1 = jet2;
148  jet2 = jet;
149  }
150  // if(fabs(jet1.eta())>eta_1 || (fabs(jet2.eta())-jet_R) < eta_2){ return;}
151  } else {return;}
152 
154  hiDistrTagJetEta_->Fill(jet1.eta());
155  hiDistrTagJetPhi_->Fill(jet1.phi());
156 
158  hiDistrProbeJetEta_->Fill(jet2.eta());
159  hiDistrProbeJetPhi_->Fill(jet2.phi());
160 
161  if(jets->size()>2){
162  jet3 = (*jets)[2];
163  etVetoJet = jet3.et();
164  hiDistrEtThirdJet_->Fill(etVetoJet);
165  } else { etVetoJet = 0.; hiDistrEtThirdJet_->Fill(etVetoJet); }
166 
167 
168 
170  iEvent.getByToken(ec_,ec);
171 
172  if(!ec.isValid()){
173  LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get ec product!" << std::endl;
174  return ;
175  }
176 
177 
178  for(EcalRecHitCollection::const_iterator ecItr = (*ec).begin();
179  ecItr != (*ec).end(); ++ecItr)
180  {
181  hiDistrRecHitEnergyEBEE_->Fill(ecItr->energy());
182  }
183 
184 
185 
186 
188  iEvent.getByToken(hbhe_, hbhe);
189 
190  if(!hbhe.isValid()){
191  LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get hbhe product!" << std::endl;
192  return ;
193  }
194 
195 
196  for(HBHERecHitCollection::const_iterator hbheItr=hbhe->begin();
197  hbheItr!=hbhe->end(); hbheItr++)
198  {
199  hiDistrRecHitEnergyHBHE_->Fill(hbheItr->energy());
200  }
201 
202 
204  iEvent.getByToken(ho_, ho);
205 
206  if(!ho.isValid()){
207  LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get ho product!" << std::endl;
208  return ;
209  }
210 
211 
213  hoItr!=ho->end(); hoItr++)
214  {
215  hiDistrRecHitEnergyHO_->Fill(hoItr->energy());
216 
217  }
218 
219 
220 
221 
223  iEvent.getByToken(hf_, hf);
224 
225  if(!hf.isValid()){
226  LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get hf product!" << std::endl;
227  return ;
228  }
229 
230 
232  hfItr!=hf->end(); hfItr++)
233  {
234  hiDistrRecHitEnergyHF_->Fill(hfItr->energy());
235  }
236 
237 } //analyze
238 
239 
240 
#define LogDebug(id)
edm::EDGetTokenT< reco::CaloJetCollection > jets_
object to monitor
MonitorElement * hiDistrTagJetEnergy_
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
void analyze(const edm::Event &e, const edm::EventSetup &c) override
MonitorElement * hiDistrRecHitEnergyHF_
double eta() const final
momentum pseudorapidity
Jets made from CaloTowers.
Definition: CaloJet.h:29
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
std::string fileName_
Output file name if required.
MonitorElement * hiDistrTagJetPhi_
MonitorElement * hiDistrTagJetEta_
std::vector< EcalRecHit >::const_iterator const_iterator
MonitorElement * hiDistrRecHitEnergyEBEE_
return((rh^lh)&mask)
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:230
double et() const final
transverse energy
edm::EDGetTokenT< HORecHitCollection > ho_
vector< PseudoJet > jets
edm::EDGetTokenT< HFRecHitCollection > hf_
double energy() const final
energy
MonitorElement * hiDistrProbeJetEta_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
bool isValid() const
Definition: HandleBase.h:74
DQMHcalDiJetsAlCaReco(const edm::ParameterSet &)
const_iterator end() const
MonitorElement * hiDistrProbeJetPhi_
edm::EDGetTokenT< HBHERecHitCollection > hbhe_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
bool saveToFile_
Write to file.
fixed size matrix
HLT enums.
edm::EDGetTokenT< EcalRecHitCollection > ec_
MonitorElement * hiDistrRecHitEnergyHBHE_
MonitorElement * hiDistrProbeJetEnergy_
MonitorElement * hiDistrRecHitEnergyHO_
std::string folderName_
DQM folder name.
double phi() const final
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_
const_iterator begin() const
Definition: Run.h:44
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override