CMS 3D CMS Logo

L1TDTTPGClient.cc
Go to the documentation of this file.
2 
11 #include "TRandom.h"
12 
13 #include <TF1.h>
14 #include <cstdio>
15 #include <sstream>
16 #include <cmath>
17 #include <TProfile.h>
18 #include <TProfile2D.h>
19 
20 using namespace edm;
21 using namespace std;
22 
24 {
25  parameters_=ps;
26  initialize();
27 }
28 
30  LogInfo("TriggerDQM")<<"[TriggerDQM]: ending... ";
31 }
32 
33 //--------------------------------------------------------
35 
36  counterLS_=0;
37  counterEvt_=0;
38 
39  // get back-end interface
40  //dbe_ = Service<DQMStore>().operator->();
41 
42  // base folder for the contents of this job
43  monitorName_ = parameters_.getUntrackedParameter<string>("monitorName","");
44 // cout << "Monitor name = " << monitorName_ << endl;
45  prescaleLS_ = parameters_.getUntrackedParameter<int>("prescaleLS", -1);
46 // cout << "DQM lumi section prescale = " << prescaleLS_ << " lumi section(s)"<< endl;
47  prescaleEvt_ = parameters_.getUntrackedParameter<int>("prescaleEvt", -1);
48 // cout << "DQM event prescale = " << prescaleEvt_ << " events(s)"<< endl;
49  output_dir_ = parameters_.getUntrackedParameter<string>("output_dir","");
50 // cout << "DQM output dir = " << output_dir_ << endl;
51  input_dir_ = parameters_.getUntrackedParameter<string>("input_dir","");
52 // cout << "DQM input dir = " << input_dir_ << endl;
53 
54  LogInfo( "TriggerDQM");
55 }
56 
57 //--------------------------------------------------------
59 
60  ibooker.setCurrentFolder(output_dir_);
61 
62  // booking
63 
64  dttpgphmapcorrf = ibooker.book2D("DT_TPG_phi_map_corr_frac",
65  "Fraction of correlated best triggers per station",20,1,21,12,0,12);
66  dttpgphmap2ndf = ibooker.book2D("DT_TPG_phi_map_2nd_frac",
67  "Fraction of second tracks per station",20,1,21,12,0,12);
68  dttpgphmapbxf[0] = ibooker.book2D("DT_TPG_phi_map_bx-1_frac",
69  "Fraction of triggers per station (BX=-1)",20,1,21,12,0,12);
70  dttpgphmapbxf[1] = ibooker.book2D("DT_TPG_phi_map_bx0_frac",
71  "Fraction of triggers per station (BX=0)",20,1,21,12,0,12);
72  dttpgphmapbxf[2] = ibooker.book2D("DT_TPG_phi_map_bx+1_frac",
73  "Fraction of triggers per station (BX=1)",20,1,21,12,0,12);
74  setMapPhLabel(dttpgphmapcorrf);
75  setMapPhLabel(dttpgphmap2ndf);
76  setMapPhLabel(dttpgphmapbxf[0]);
77  setMapPhLabel(dttpgphmapbxf[1]);
78  setMapPhLabel(dttpgphmapbxf[2]);
79 
80  dttpgthmaphf = ibooker.book2D("DT_TPG_theta_map_corr_frac",
81  "Fraction of H quality best triggers per station",15,1,16,12,0,12);
82  dttpgthmapbxf[0] = ibooker.book2D("DT_TPG_theta_map_bx-1_frac",
83  "Fraction of triggers per station (BX=-1)",15,1,16,12,0,12);
84  dttpgthmapbxf[1] = ibooker.book2D("DT_TPG_theta_map_bx0_frac",
85  "Fraction of triggers per station (BX=0)",15,1,16,12,0,12);
86  dttpgthmapbxf[2] = ibooker.book2D("DT_TPG_theta_map_bx+1_frac",
87  "Fraction of triggers per station (BX=1)",15,1,16,12,0,12);
88  setMapThLabel(dttpgthmaphf);
89  setMapThLabel(dttpgthmapbxf[0]);
90  setMapThLabel(dttpgthmapbxf[1]);
91  setMapThLabel(dttpgthmapbxf[2]);
92 
93 // cout << "L1TDTTPGClient::analyze" << endl;
94  counterEvt_++;
95  if (prescaleEvt_<1) return;
96  if (prescaleEvt_>0 && counterEvt_%prescaleEvt_ != 0) return;
97 
98  string nName = "DT_TPG_phi_best_map_corr";
99  string dName = "DT_TPG_phi_best_map";
100  makeRatioHisto(igetter,dttpgphmapcorrf,nName,dName);
101  dName = "DT_TPG_phi_map";
102  nName = "DT_TPG_phi_map_2nd";
103  makeRatioHisto(igetter,dttpgphmap2ndf,nName,dName);
104  nName = "DT_TPG_phi_map_bx-1";
105  makeRatioHisto(igetter,dttpgphmapbxf[0],nName,dName);
106  nName = "DT_TPG_phi_map_bx0";
107  makeRatioHisto(igetter,dttpgphmapbxf[1],nName,dName);
108  nName = "DT_TPG_phi_map_bx+1";
109  makeRatioHisto(igetter,dttpgphmapbxf[2],nName,dName);
110 
111  nName = "DT_TPG_theta_best_map_h";
112  dName = "DT_TPG_theta_best_map";
113  makeRatioHisto(igetter,dttpgthmaphf,nName,dName);
114  dName = "DT_TPG_theta_map";
115  nName = "DT_TPG_theta_map_bx-1";
116  makeRatioHisto(igetter,dttpgthmapbxf[0],nName,dName);
117  nName = "DT_TPG_theta_map_bx0";
118  makeRatioHisto(igetter,dttpgthmapbxf[1],nName,dName);
119  nName = "DT_TPG_theta_map_bx+1";
120  makeRatioHisto(igetter,dttpgthmapbxf[2],nName,dName);
121 }
122 
123 void L1TDTTPGClient::makeRatioHisto(DQMStore::IGetter &igetter, MonitorElement *ratioME, string &nName, string &dName)
124 {
125 
126  TH2F *numerator;
127  TH2F *denominator;
128 
129  denominator = this->get2DHisto(input_dir_+"/"+dName,igetter);
130  numerator = this->get2DHisto(input_dir_+"/"+nName,igetter);
131 
132  if (numerator && denominator) {
133 
134  TH2F * ratio = ratioME->getTH2F();
135  if (ratio) {
136  ratio->Divide(numerator,denominator);
137  }
138  else {
139  LogInfo("TriggerDQM") << "[TriggerDQM]: ratio histo named \"" << ratioME->getName() << "\" not found!" << endl;
140  }
141  }
142  else {
143  if (!numerator)
144  LogInfo("TriggerDQM") << "[TriggerDQM]: numerator histo \"" << nName << "\" not found!" << endl;
145  if (!denominator)
146  LogInfo("TriggerDQM") << "[TriggerDQM]: denominator histo \"" << dName << "\" not found!" << endl;
147  }
148 
149 }
150 
151 TH1F * L1TDTTPGClient::get1DHisto(string meName, DQMStore::IGetter &igetter)
152 {
153 
154  MonitorElement * me_ = igetter.get(meName);
155 
156  if (!me_) {
157  LogInfo("TriggerDQM") << "ME NOT FOUND.";
158  return nullptr;
159  }
160 
161  return me_->getTH1F();
162 }
163 
164 TH2F * L1TDTTPGClient::get2DHisto(string meName, DQMStore::IGetter &igetter)
165 {
166 
167 
168  MonitorElement * me_ = igetter.get(meName);
169 
170  if (!me_) {
171  LogInfo("TriggerDQM") << "ME NOT FOUND.";
172  return nullptr;
173  }
174 
175  return me_->getTH2F();
176 }
177 
178 
179 
180 TProfile2D * L1TDTTPGClient::get2DProfile(string meName, DQMStore::IGetter &igetter)
181 {
182 
183 
184  MonitorElement * me_ = igetter.get(meName);
185 
186  if (!me_) {
187  LogInfo("TriggerDQM") << "ME NOT FOUND.";
188  return nullptr;
189  }
190 
191  return me_->getTProfile2D();
192 }
193 
194 
195 TProfile * L1TDTTPGClient::get1DProfile(string meName, DQMStore::IGetter &igetter)
196 {
197 
198 
199  MonitorElement * me_ = igetter.get(meName);
200 
201  if (!me_) {
202  LogInfo("TriggerDQM") << "ME NOT FOUND.";
203  return nullptr;
204  }
205 
206  return me_->getTProfile();
207 }
208 
209 
211 {
212 
213  me->setAxisTitle("DTTF Sector",2);
214  for(int i=0;i<5;i++){
215  ostringstream wheel;
216  wheel << i-2;
217  me->setBinLabel(1+i*4,"Wheel "+ wheel.str(),1);
218  }
219 
220 }
221 
222 
224 {
225 
226  me->setAxisTitle("DTTF Sector",2);
227  for(int i=0;i<5;i++){
228  ostringstream wheel;
229  wheel << i-2;
230  me->setBinLabel(1+i*3,"Wheel "+ wheel.str(),1);
231  }
232 
233 }
static AlgebraicMatrix initialize()
const std::string & getName(void) const
get name of ME
void setMapThLabel(MonitorElement *me)
TProfile2D * getTProfile2D(void) const
TProfile * get1DProfile(std::string meName, DQMStore::IGetter &igetter)
MonitorElement * get(const std::string &path)
Definition: DQMStore.cc:305
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
TH1F * get1DHisto(std::string meName, DQMStore::IGetter &igetter)
TProfile2D * get2DProfile(std::string meName, DQMStore::IGetter &igetter)
TH2F * get2DHisto(std::string meName, DQMStore::IGetter &igetter)
void dqmEndJob(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) override
~L1TDTTPGClient() override
Destructor.
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
L1TDTTPGClient(const edm::ParameterSet &ps)
Constructor.
TH1F * getTH1F(void) const
void makeRatioHisto(DQMStore::IGetter &igetter, MonitorElement *ratioME, std::string &nName, std::string &dName)
HLT enums.
TProfile * getTProfile(void) const
void setMapPhLabel(MonitorElement *me)
TH2F * getTH2F(void) const
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)