CMS 3D CMS Logo

EcalIntercalibConstantsMC_PayloadInspector.cc
Go to the documentation of this file.
8 // the data format of the condition to be inspected
10 
11 #include "TH2F.h"
12 #include "TCanvas.h"
13 #include "TStyle.h"
14 #include "TLine.h"
15 #include "TLatex.h"
16 
17 #include <memory>
18 #include <sstream>
19 
20 namespace {
21  enum {kEBChannels = 61200, kEEChannels = 14648};
22  enum {MIN_IETA = 1, MIN_IPHI = 1, MAX_IETA = 85, MAX_IPHI = 360, EBhistEtaMax = 171}; // barrel lower and upper bounds on eta and phi
23  enum {IX_MIN = 1, IY_MIN = 1, IX_MAX = 100, IY_MAX = 100, EEhistXMax = 220}; // endcaps lower and upper bounds on x and y
24 
25  /*******************************************************
26 
27  2d histogram of ECAL barrel Intercalib Constants MC of 1 IOV
28 
29  *******************************************************/
30 
31  // inherit from one of the predefined plot class: Histogram2D
32  class EcalIntercalibConstantsMCEBMap : public cond::payloadInspector::Histogram2D<EcalIntercalibConstantsMC> {
33 
34  public:
35  EcalIntercalibConstantsMCEBMap() : cond::payloadInspector::Histogram2D<EcalIntercalibConstantsMC>("ECAL Barrel Intercalib Constants MC- map ",
36  "iphi", MAX_IPHI, MIN_IPHI, MAX_IPHI + 1,
37  "ieta", EBhistEtaMax, -MAX_IETA, MAX_IETA + 1) {
38  Base::setSingleIov( true );
39  }
40 
41  // Histogram2D::fill (virtual) needs be overridden - the implementation should use fillWithValue
42  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
43 
44  for (auto const & iov: iovs) {
45  std::shared_ptr<EcalIntercalibConstantsMC> payload = Base::fetchPayload( std::get<1>(iov) );
46  if( payload.get() ){
47  // looping over the EB channels, via the dense-index, mapped into EBDetId's
48  if (payload->barrelItems().empty()) return false;
49 
50  // set to -1 for ieta 0 (no crystal)
51  for(int iphi = MIN_IPHI; iphi < MAX_IPHI+1; iphi++) fillWithValue(iphi, 0, -1);
52 
53  for(int cellid = EBDetId::MIN_HASH; cellid < EBDetId::kSizeForDenseIndexing; ++cellid) {
54  uint32_t rawid = EBDetId::unhashIndex(cellid);
55 
56  // check the existence of ECAL Intercalib Constants, for a given ECAL barrel channel
57  EcalFloatCondObjectContainer::const_iterator value_ptr = payload->find(rawid);
58  if (value_ptr == payload->end())
59  continue; // cell absent from payload
60 
61  float weight = (float)(*value_ptr);
62 
63  // fill the Histogram2D here
64  fillWithValue( (EBDetId(rawid)).iphi() , (EBDetId(rawid)).ieta(), weight);
65  }// loop over cellid
66 
67 
68 
69  }// if payload.get()
70  }// loop over IOV's (1 in this case)
71 
72  return true;
73 
74  }//fill method
75  };
76 
77  /*******************************************************
78 
79  2d histogram of ECAL EndCaps Intercalib Constants MC of 1 IOV
80 
81  *******************************************************/
82 
83  class EcalIntercalibConstantsMCEEMap : public cond::payloadInspector::Histogram2D<EcalIntercalibConstantsMC> {
84 
85  private:
86  int EEhistSplit = 20;
87 
88  public:
89  EcalIntercalibConstantsMCEEMap() : cond::payloadInspector::Histogram2D<EcalIntercalibConstantsMC>( "ECAL Endcap Intercalib Constants MC- map ",
90  "ix", EEhistXMax, IX_MIN, EEhistXMax + 1,
91  "iy", IY_MAX, IY_MIN, IY_MAX + 1) {
92  Base::setSingleIov( true );
93  }
94 
95  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
96 
97  for (auto const & iov: iovs) {
98  std::shared_ptr<EcalIntercalibConstantsMC> payload = Base::fetchPayload( std::get<1>(iov) );
99  if( payload.get() ){
100  if (payload->endcapItems().empty()) return false;
101 
102  // set to -1 everywhwere
103  for(int ix = IX_MIN; ix < EEhistXMax + 1; ix++)
104  for(int iy = IY_MAX; iy < IY_MAX + 1; iy++)
105  fillWithValue(ix, iy, -1);
106 
107  for (int cellid = 0; cellid < EEDetId::kSizeForDenseIndexing; ++cellid){ // loop on EE cells
108  if (EEDetId::validHashIndex(cellid)){
109  uint32_t rawid = EEDetId::unhashIndex(cellid);
110  EcalFloatCondObjectContainer::const_iterator value_ptr = payload->find(rawid);
111  if (value_ptr == payload->end())
112  continue; // cell absent from payload
113 
114  float weight = (float)(*value_ptr);
115  EEDetId myEEId(rawid);
116  if(myEEId.zside() == -1)
117  fillWithValue(myEEId.ix(), myEEId.iy(), weight);
118  else
119  fillWithValue(myEEId.ix() + IX_MAX + EEhistSplit, myEEId.iy(), weight);
120  } // validDetId
121  } // loop over cellid
122 
123  } // payload
124  } // loop over IOV's (1 in this case)
125  return true;
126  }// fill method
127 
128  };
129 
130 
131  /*************************************************
132  2d plot of ECAL IntercalibConstantsMC of 1 IOV
133  *************************************************/
134  class EcalIntercalibConstantsMCPlot : public cond::payloadInspector::PlotImage<EcalIntercalibConstantsMC> {
135 
136  public:
137  EcalIntercalibConstantsMCPlot() : cond::payloadInspector::PlotImage<EcalIntercalibConstantsMC>("ECAL Intercalib Constants MC - map ") {
138  setSingleIov( true );
139  }
140 
141  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
142  TH2F* barrel = new TH2F("EB", "mean EB", MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
143  TH2F* endc_p = new TH2F("EE+", "mean EE+", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
144  TH2F* endc_m = new TH2F("EE-", "mean EE-", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
145 
146  auto iov = iovs.front();
147  std::shared_ptr<EcalIntercalibConstantsMC> payload = fetchPayload( std::get<1>(iov) );
148  unsigned int run = std::get<0>(iov);
149 
150  if( payload.get() ){
151 
152  if (payload->barrelItems().empty())
153  return false;
154 
155  fillEBMap_SingleIOV<EcalIntercalibConstantsMC>(payload, barrel);
156 
157 
158  if (payload->endcapItems().empty())
159  return false;
160 
161  fillEEMap_SingleIOV<EcalIntercalibConstantsMC>(payload, endc_m, endc_p);
162 
163  } // payload
164 
165  gStyle->SetPalette(1);
166  gStyle->SetOptStat(0);
167  TCanvas canvas("CC map","CC map", 1600, 450);
168  TLatex t1;
169  t1.SetNDC();
170  t1.SetTextAlign(26);
171  t1.SetTextSize(0.04);
172  t1.DrawLatex(0.5, 0.96, Form("Ecal IntercalibConstants MC, IOV %i", run));
173 
174  float xmi[3] = {0.0 , 0.24, 0.76};
175  float xma[3] = {0.24, 0.76, 1.00};
176  TPad** pad = new TPad*;
177  for (int obj = 0; obj < 3; obj++) {
178  pad[obj] = new TPad(Form("p_%i", obj),Form("p_%i", obj), xmi[obj], 0.0, xma[obj], 0.94);
179  pad[obj]->Draw();
180  }
181  // EcalDrawMaps ICMap;
182  pad[0]->cd();
183  // ICMap.DrawEE(endc_m, 0., 2.);
184  DrawEE(endc_m, 0., 2.5);
185  pad[1]->cd();
186  // ICMap.DrawEB(barrel, 0., 2.);
187  DrawEB(barrel, 0., 2.5);
188  pad[2]->cd();
189  // ICMap.DrawEE(endc_p, 0., 2.);
190  DrawEE(endc_p, 0., 2.5);
191 
192  std::string ImageName(m_imageFileName);
193  canvas.SaveAs(ImageName.c_str());
194  return true;
195  }// fill method
196  };
197 
198  /*****************************************************************
199  2d plot of ECAL IntercalibConstantsMC difference between 2 IOVs
200  *****************************************************************/
201  class EcalIntercalibConstantsMCDiff : public cond::payloadInspector::PlotImage<EcalIntercalibConstantsMC> {
202 
203  public:
204  EcalIntercalibConstantsMCDiff() : cond::payloadInspector::PlotImage<EcalIntercalibConstantsMC>("ECAL Intercalib Constants MC difference ") {
205  setSingleIov(false);
206  }
207 
208  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
209  TH2F* barrel = new TH2F("EB", "mean EB", MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
210  TH2F* endc_p = new TH2F("EE+", "mean EE+", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
211  TH2F* endc_m = new TH2F("EE-", "mean EE-", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
212  float pEBmin, pEEmin, pEBmax, pEEmax;
213  pEBmin = 10.;
214  pEEmin = 10.;
215  pEBmax = -10.;
216  pEEmax = -10.;
217 
218  unsigned int run[2], irun = 0;
219  float pEB[kEBChannels], pEE[kEEChannels];
220  for ( auto const & iov: iovs) {
221 
222  std::shared_ptr<EcalIntercalibConstantsMC> payload = fetchPayload( std::get<1>(iov) );
223  run[irun] = std::get<0>(iov);
224 
225  if( payload.get() ){
226 
227  if (payload->barrelItems().empty())
228  return false;
229 
230  fillEBMap_DiffIOV<EcalIntercalibConstantsMC>(payload, barrel, irun, pEB, pEBmin, pEBmax);
231 
232  if (payload->endcapItems().empty())
233  return false;
234 
235  fillEEMap_DiffIOV<EcalIntercalibConstantsMC>(payload, endc_m, endc_p, irun, pEE, pEEmin, pEEmax);
236 
237  }// payload
238  irun++;
239 
240  }// loop over IOVs
241 
242  gStyle->SetPalette(1);
243  gStyle->SetOptStat(0);
244  TCanvas canvas("CC map","CC map", 1600, 450);
245  TLatex t1;
246  t1.SetNDC();
247  t1.SetTextAlign(26);
248  t1.SetTextSize(0.04);
249  t1.DrawLatex(0.5, 0.96, Form("Ecal IntercalibConstants MC Diff, IOV %i - %i", run[1], run[0]));
250 
251  float xmi[3] = {0.0 , 0.24, 0.76};
252  float xma[3] = {0.24, 0.76, 1.00};
253  TPad** pad = new TPad*;
254 
255  for (int obj = 0; obj < 3; obj++) {
256  pad[obj] = new TPad(Form("p_%i", obj),Form("p_%i", obj), xmi[obj], 0.0, xma[obj], 0.94);
257  pad[obj]->Draw();
258  }
259 
260  pad[0]->cd();
261  DrawEE(endc_m, pEEmin, pEEmax);
262  pad[1]->cd();
263  DrawEB(barrel, pEBmin, pEBmax);
264  pad[2]->cd();
265  DrawEE(endc_p, pEEmin, pEEmax);
266 
267  std::string ImageName(m_imageFileName);
268  canvas.SaveAs(ImageName.c_str());
269  return true;
270  }// fill method
271  };
272 
273 /*******************************************************
274  2d plot of Ecal Intercalib Constants MC Summary of 1 IOV
275  *******************************************************/
276 class EcalIntercalibConstantsMCSummaryPlot: public cond::payloadInspector::PlotImage<EcalIntercalibConstantsMC>{
277  public:
278  EcalIntercalibConstantsMCSummaryPlot():
279  cond::payloadInspector::PlotImage<EcalIntercalibConstantsMC>("Ecal Intercalib Constants MC Summary - map "){
280  setSingleIov(true);
281  }
282 
283  bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs)override {
284  auto iov=iovs.front();
285  std::shared_ptr <EcalIntercalibConstantsMC> payload = fetchPayload(std::get<1> (iov));
286  unsigned int run=std::get<0> (iov);
287  TH2F* align;
288  int NbRows;
289 
290  if(payload.get()){
291  NbRows=2;
292  align=new TH2F("", "", 0, 0, 0, 0, 0, 0);
293 
294 
295  float mean_x_EB=0.0f;
296  float mean_x_EE=0.0f;
297 
298  float rms_EB=0.0f;
299  float rms_EE=0.0f;
300 
301  int num_x_EB=0;
302  int num_x_EE=0;
303 
304 
305  payload->summary(mean_x_EB, rms_EB, num_x_EB, mean_x_EE, rms_EE, num_x_EE);
306  fillTableWithSummary(align, "Ecal Intercalib Constants MC", mean_x_EB, rms_EB, num_x_EB, mean_x_EE, rms_EE, num_x_EE);
307 
308  }else
309  return false;
310 
311  gStyle->SetPalette(1);
312  gStyle->SetOptStat(0);
313  TCanvas canvas("CC map", "CC map", 1000, 1000);
314  TLatex t1;
315  t1.SetNDC();
316  t1.SetTextAlign(26);
317  t1.SetTextSize(0.04);
318  t1.SetTextColor(2);
319  t1.DrawLatex(0.5, 0.96,Form("Ecal Intercalib Constants MC Summary, IOV %i", run));
320 
321 
322  TPad* pad = new TPad("pad", "pad", 0.0, 0.0, 1.0, 0.94);
323  pad->Draw();
324  pad->cd();
325  align->Draw("TEXT");
326 
327  drawTable(NbRows,4);
328 
329  std::string ImageName(m_imageFileName);
330  canvas.SaveAs(ImageName.c_str());
331 
332  return true;
333  }
334 };
335 
336 
337 
338 }
339 
340 
341 // Register the classes as boost python plugin
343  PAYLOAD_INSPECTOR_CLASS( EcalIntercalibConstantsMCEBMap);
344  PAYLOAD_INSPECTOR_CLASS( EcalIntercalibConstantsMCEEMap);
345  PAYLOAD_INSPECTOR_CLASS( EcalIntercalibConstantsMCPlot );
346  PAYLOAD_INSPECTOR_CLASS( EcalIntercalibConstantsMCDiff );
347  PAYLOAD_INSPECTOR_CLASS( EcalIntercalibConstantsMCSummaryPlot );
348 }
void fillTableWithSummary(TH2F *&align, std::string title, const float &mean_x_EB, const float &rms_EB, const int &num_x_EB, const float &mean_x_EE, const float &rms_EE, const int &num_x_EE)
Definition: weight.py:1
void fillWithValue(float xvalue, float yvalue, float weight=1)
static EEDetId unhashIndex(int hi)
Definition: EEDetId.cc:65
static bool validHashIndex(int i)
Definition: EEDetId.h:239
void DrawEE(TH2F *endc, float min, float max)
Definition: EcalDrawUtils.h:29
void DrawEB(TH2F *ebmap, float min, float max)
Definition: EcalDrawUtils.h:4
bool fill(const std::vector< std::tuple< cond::Time_t, cond::Hash > > &iovs) override
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
static const int MIN_HASH
Definition: EBDetId.h:149
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
std::vector< Item >::const_iterator const_iterator
static EBDetId unhashIndex(int hi)
get a DetId from a compact index for arrays
Definition: EBDetId.h:110
Definition: plugin.cc:24
def canvas(sub, attr)
Definition: svgfig.py:482
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)
void drawTable(int nbRows, int nbColumns)