CMS 3D CMS Logo

EcalIntercalibConstants_PayloadInspector.cc
Go to the documentation of this file.
7 
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 of 1 IOV
28 
29  *******************************************************/
30 
31  // inherit from one of the predefined plot class: Histogram2D
32  class EcalIntercalibConstantsEBMap : public cond::payloadInspector::Histogram2D<EcalIntercalibConstants> {
33 
34  public:
35  EcalIntercalibConstantsEBMap() : cond::payloadInspector::Histogram2D<EcalIntercalibConstants>("ECAL Barrel Intercalib Constants - 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<EcalIntercalibConstants> 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  // set to -1 for ieta 0 (no crystal)
50  for(int iphi = MIN_IPHI; iphi < MAX_IPHI+1; iphi++) fillWithValue(iphi, 0, -1);
51 
52  for(int cellid = EBDetId::MIN_HASH; cellid < EBDetId::kSizeForDenseIndexing; ++cellid) {
53  uint32_t rawid = EBDetId::unhashIndex(cellid);
54 
55  // check the existence of ECAL Intercalib Constants, for a given ECAL barrel channel
56  EcalFloatCondObjectContainer::const_iterator value_ptr = payload->find(rawid);
57  if (value_ptr == payload->end()) continue; // cell absent from payload
58  float weight = (float)(*value_ptr);
59 
60  // fill the Histogram2D here
61  fillWithValue( (EBDetId(rawid)).iphi() , (EBDetId(rawid)).ieta(), weight);
62  }// loop over cellid
63  }// if payload.get()
64  }// loop over IOV's (1 in this case)
65 
66  return true;
67  }// fill method
68  };
69 
70  class EcalIntercalibConstantsEEMap : public cond::payloadInspector::Histogram2D<EcalIntercalibConstants> {
71 
72  private:
73  int EEhistSplit = 20;
74 
75  public:
76  EcalIntercalibConstantsEEMap() : cond::payloadInspector::Histogram2D<EcalIntercalibConstants>( "ECAL Endcap Intercalib Constants - map ",
77  "ix", EEhistXMax, IX_MIN, EEhistXMax + 1,
78  "iy", IY_MAX, IY_MIN, IY_MAX + 1) {
79  Base::setSingleIov( true );
80  }
81 
82  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
83 
84  for (auto const & iov: iovs) {
85  std::shared_ptr<EcalIntercalibConstants> payload = Base::fetchPayload( std::get<1>(iov) );
86  if( payload.get() ){
87  if (payload->endcapItems().empty()) return false;
88 
89  // set to -1 everywhwere
90  for(int ix = IX_MIN; ix < EEhistXMax + 1; ix++)
91  for(int iy = IY_MAX; iy < IY_MAX + 1; iy++)
92  fillWithValue(ix, iy, -1);
93 
94  for (int cellid = 0; cellid < EEDetId::kSizeForDenseIndexing; ++cellid){ // loop on EE cells
95  if (EEDetId::validHashIndex(cellid)){
96  uint32_t rawid = EEDetId::unhashIndex(cellid);
97  EcalFloatCondObjectContainer::const_iterator value_ptr = payload->find(rawid);
98  if (value_ptr == payload->end()) continue; // cell absent from payload
99  float weight = (float)(*value_ptr);
100  EEDetId myEEId(rawid);
101  if(myEEId.zside() == -1)
102  fillWithValue(myEEId.ix(), myEEId.iy(), weight);
103  else
104  fillWithValue(myEEId.ix() + IX_MAX + EEhistSplit, myEEId.iy(), weight);
105  } // validDetId
106  } // loop over cellid
107  } // payload
108  } // loop over IOV's (1 in this case)
109  return true;
110  }// fill method
111  };
112 
113  /*************************************************
114  2d plot of ECAL IntercalibConstants of 1 IOV
115  *************************************************/
116  class EcalIntercalibConstantsPlot : public cond::payloadInspector::PlotImage<EcalIntercalibConstants> {
117 
118  public:
119  EcalIntercalibConstantsPlot() : cond::payloadInspector::PlotImage<EcalIntercalibConstants>("ECAL Intercalib Constants - map ") {
120  setSingleIov( true );
121  }
122 
123  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
124  TH2F* barrel = new TH2F("EB", "mean EB", MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
125  TH2F* endc_p = new TH2F("EE+", "mean EE+", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
126  TH2F* endc_m = new TH2F("EE-", "mean EE-", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
127 
128  auto iov = iovs.front();
129  std::shared_ptr<EcalIntercalibConstants> payload = fetchPayload( std::get<1>(iov) );
130  unsigned int run = std::get<0>(iov);
131  if( payload.get() ){
132  if (payload->barrelItems().empty()) return false;
133  for(int cellid = EBDetId::MIN_HASH; cellid < EBDetId::kSizeForDenseIndexing; ++cellid) {
134  uint32_t rawid = EBDetId::unhashIndex(cellid);
135  EcalFloatCondObjectContainer::const_iterator value_ptr = payload->find(rawid);
136  if (value_ptr == payload->end()) continue; // cell absent from payload
137  float weight = (float)(*value_ptr);
138  Double_t phi = (Double_t)(EBDetId(rawid)).iphi() - 0.5;
139  Double_t eta = (Double_t)(EBDetId(rawid)).ieta();
140  if(eta > 0.) eta = eta - 0.5; // 0.5 to 84.5
141  else eta = eta + 0.5; // -84.5 to -0.5
142  barrel->Fill(phi, eta, weight);
143  }// loop over cellid
144 
145  if (payload->endcapItems().empty()) return false;
146  // looping over the EE channels
147  for(int iz = -1; iz < 2; iz = iz + 2) // -1 or +1
148  for(int iy = IY_MIN; iy < IY_MAX+IY_MIN; iy++)
149  for(int ix = IX_MIN; ix < IX_MAX+IX_MIN; ix++)
150  if(EEDetId::validDetId(ix, iy, iz)) {
151  EEDetId myEEId = EEDetId(ix, iy, iz, EEDetId::XYMODE);
152  uint32_t rawid = myEEId.rawId();
153  EcalFloatCondObjectContainer::const_iterator value_ptr = payload->find(rawid);
154  if (value_ptr == payload->end()) continue; // cell absent from payload
155  float weight = (float)(*value_ptr);
156  if(iz == 1)
157  endc_p->Fill(ix, iy, weight);
158  else
159  endc_m->Fill(ix, iy, weight);
160  } // validDetId
161  } // payload
162 
163  gStyle->SetPalette(1);
164  gStyle->SetOptStat(0);
165  TCanvas canvas("CC map","CC map", 1600, 450);
166  TLatex t1;
167  t1.SetNDC();
168  t1.SetTextAlign(26);
169  t1.SetTextSize(0.05);
170  t1.DrawLatex(0.5, 0.96, Form("Ecal IntercalibConstants, IOV %i", run));
171 
172  float xmi[3] = {0.0 , 0.24, 0.76};
173  float xma[3] = {0.24, 0.76, 1.00};
174  TPad** pad = new TPad*;
175  for (int obj = 0; obj < 3; obj++) {
176  pad[obj] = new TPad(Form("p_%i", obj),Form("p_%i", obj), xmi[obj], 0.0, xma[obj], 0.94);
177  pad[obj]->Draw();
178  }
179  // EcalDrawMaps ICMap;
180  pad[0]->cd();
181  // ICMap.DrawEE(endc_m, 0., 2.);
182  DrawEE(endc_m, 0., 2.5);
183  pad[1]->cd();
184  // ICMap.DrawEB(barrel, 0., 2.);
185  DrawEB(barrel, 0., 2.5);
186  pad[2]->cd();
187  // ICMap.DrawEE(endc_p, 0., 2.);
188  DrawEE(endc_p, 0., 2.5);
189 
190  std::string ImageName(m_imageFileName);
191  canvas.SaveAs(ImageName.c_str());
192  return true;
193  }// fill method
194  };
195 
196  /*****************************************************************
197  2d plot of ECAL IntercalibConstants difference between 2 IOVs
198  *****************************************************************/
199  class EcalIntercalibConstantsDiff : public cond::payloadInspector::PlotImage<EcalIntercalibConstants> {
200 
201  public:
202  EcalIntercalibConstantsDiff() : cond::payloadInspector::PlotImage<EcalIntercalibConstants>("ECAL Intercalib Constants difference ") {
203  setSingleIov(false);
204  }
205 
206  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
207  TH2F* barrel = new TH2F("EB", "mean EB", MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
208  TH2F* endc_p = new TH2F("EE+", "mean EE+", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
209  TH2F* endc_m = new TH2F("EE-", "mean EE-", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
210  float pEBmin, pEEmin, pEBmax, pEEmax;
211  pEBmin = 10.;
212  pEEmin = 10.;
213  pEBmax = -10.;
214  pEEmax = -10.;
215 
216  unsigned int run[2], irun = 0;
217  float pEB[kEBChannels], pEE[kEEChannels];
218  for ( auto const & iov: iovs) {
219  std::shared_ptr<EcalIntercalibConstants> payload = fetchPayload( std::get<1>(iov) );
220  run[irun] = std::get<0>(iov);
221  if( payload.get() ){
222  if (payload->barrelItems().empty()) return false;
223  for(int cellid = EBDetId::MIN_HASH; cellid < EBDetId::kSizeForDenseIndexing; ++cellid) {
224  uint32_t rawid = EBDetId::unhashIndex(cellid);
225  EcalFloatCondObjectContainer::const_iterator value_ptr = payload->find(rawid);
226  if (value_ptr == payload->end()) continue; // cell absent from payload
227  float weight = (float)(*value_ptr);
228  if(irun == 0)
229  pEB[cellid] = weight;
230  else {
231  Double_t phi = (Double_t)(EBDetId(rawid)).iphi() - 0.5;
232  Double_t eta = (Double_t)(EBDetId(rawid)).ieta();
233  if(eta > 0.) eta = eta - 0.5; // 0.5 to 84.5
234  else eta = eta + 0.5; // -84.5 to -0.5
235  double diff = weight - pEB[cellid];
236  if(diff < pEBmin) pEBmin = diff;
237  if(diff > pEBmax) pEBmax = diff;
238  barrel->Fill(phi, eta, diff);
239  }
240  }// loop over cellid
241 
242  if (payload->endcapItems().empty()) return false;
243  // looping over the EE channels
244  for(int iz = -1; iz < 2; iz = iz + 2) // -1 or +1
245  for(int iy = IY_MIN; iy < IY_MAX+IY_MIN; iy++)
246  for(int ix = IX_MIN; ix < IX_MAX+IX_MIN; ix++)
247  if(EEDetId::validDetId(ix, iy, iz)) {
248  EEDetId myEEId = EEDetId(ix, iy, iz, EEDetId::XYMODE);
249  uint32_t cellid = myEEId.hashedIndex();
250  uint32_t rawid = myEEId.rawId();
251  EcalFloatCondObjectContainer::const_iterator value_ptr = payload->find(rawid);
252  if (value_ptr == payload->end()) continue; // cell absent from payload
253  float weight = (float)(*value_ptr);
254  if(irun == 0)
255  pEE[cellid] = weight;
256  else {
257  double diff = weight - pEE[cellid];
258  if(diff < pEEmin) pEEmin = diff;
259  if(diff > pEEmax) pEEmax = diff;
260  if(iz == 1)
261  endc_p->Fill(ix, iy, diff);
262  else
263  endc_m->Fill(ix, iy, diff);
264  }
265  } // validDetId
266  } // payload
267  irun++;
268  } // loop over IOVs
269 
270  gStyle->SetPalette(1);
271  gStyle->SetOptStat(0);
272  TCanvas canvas("CC map","CC map", 1600, 450);
273  TLatex t1;
274  t1.SetNDC();
275  t1.SetTextAlign(26);
276  t1.SetTextSize(0.05);
277  t1.DrawLatex(0.5, 0.96, Form("Ecal IntercalibConstants, IOV %i - %i", run[1], run[0]));
278 
279  float xmi[3] = {0.0 , 0.24, 0.76};
280  float xma[3] = {0.24, 0.76, 1.00};
281  TPad** pad = new TPad*;
282  for (int obj = 0; obj < 3; obj++) {
283  pad[obj] = new TPad(Form("p_%i", obj),Form("p_%i", obj), xmi[obj], 0.0, xma[obj], 0.94);
284  pad[obj]->Draw();
285  }
286  pad[0]->cd();
287  DrawEE(endc_m, pEEmin, pEEmax);
288  pad[1]->cd();
289  DrawEB(barrel, pEBmin, pEBmax);
290  pad[2]->cd();
291  DrawEE(endc_p, pEEmin, pEEmax);
292 
293  std::string ImageName(m_imageFileName);
294  canvas.SaveAs(ImageName.c_str());
295  return true;
296  }// fill method
297  };
298 
299 } // close namespace
300 
301 // Register the classes as boost python plugin
303  PAYLOAD_INSPECTOR_CLASS( EcalIntercalibConstantsEBMap);
304  PAYLOAD_INSPECTOR_CLASS( EcalIntercalibConstantsEEMap);
305  PAYLOAD_INSPECTOR_CLASS( EcalIntercalibConstantsPlot);
306  PAYLOAD_INSPECTOR_CLASS( EcalIntercalibConstantsDiff);
307 }
static const int XYMODE
Definition: EEDetId.h:339
Definition: weight.py:1
void fillWithValue(float xvalue, float yvalue, float weight=1)
static EEDetId unhashIndex(int hi)
Definition: EEDetId.cc:99
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
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
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:156
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
int hashedIndex() const
Definition: EEDetId.h:182
static bool validDetId(int crystal_ix, int crystal_iy, int iz)
Definition: EEDetId.h:248
std::vector< Item >::const_iterator const_iterator
static EBDetId unhashIndex(int hi)
get a DetId from a compact index for arrays
Definition: EBDetId.h:114
Definition: plugin.cc:24
def canvas(sub, attr)
Definition: svgfig.py:481
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)