CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalPulseCovariancesPyWrapper.cc
Go to the documentation of this file.
5 #include "TH2F.h"
6 #include "TCanvas.h"
7 #include "TStyle.h"
8 #include "TLine.h"
11 
14 
15 #include <string>
16 #include <sstream>
17 #include <algorithm>
18 #include <numeric>
19 #include <iterator>
20 #include <boost/ref.hpp>
21 #include <boost/bind.hpp>
22 #include <boost/function.hpp>
23 #include <boost/iterator/transform_iterator.hpp>
24 #include <fstream>
25 
27 
28 #define THERAW 0
29 
30 namespace {
31  struct Printer {
32  void doit(EcalPulseCovariance const & item) {
33  for(int s = 0; s<EcalPulseShape::TEMPLATESAMPLES; ++s)
34  ss << item.val(THERAW,s);
35  ss << " ";
36  }
37  std::stringstream ss;
38  };
39 }
40 
41 namespace cond {
42 
43 
44  namespace ecalpulsecovariance {
47 
48  float average(EcalPulseCovariances const & pulseshapes, Quantity q) {
49  return std::accumulate(
50  boost::make_transform_iterator(pulseshapes.barrelItems().begin(),bind(&EcalPulseCovariance::val,_1,THERAW,q-1)),
51  boost::make_transform_iterator(pulseshapes.barrelItems().end(),bind(&EcalPulseCovariance::val,_1,THERAW,q-1)),
52  0.)/float(pulseshapes.barrelItems().size());
53  }
54 
55  void extractAverage(EcalPulseCovariances const & pulseshapes, Quantity q, std::vector<int> const &, std::vector<float> & result) {
56  result.resize(1);
57  result[0] = average(pulseshapes,q);
58  }
59 
60  void extractSuperModules(EcalPulseCovariances const & pulseshapes, Quantity q, std::vector<int> const & which, std::vector<float> & result) {
61  // bho...
62  }
63 
64  void extractSingleChannel(EcalPulseCovariances const & pulseshapes, Quantity q, std::vector<int> const & which, std::vector<float> & result) {
65  for (unsigned int i=0; i<which.size();i++) {
66  // absolutely arbitraty
67  if ((unsigned int) (which[i])< pulseshapes.barrelItems().size())
68  result.push_back( pulseshapes.barrelItems()[which[i]].val(THERAW,q-1));
69  }
70  }
71 
72  typedef boost::function<void(EcalPulseCovariances const & pulseshapes, Quantity q, std::vector<int> const & which, std::vector<float> & result)> PulseCovarianceExtractor;
73  }
74 
75  template<>
77 
80  std::vector<int> m_which;
81 
82  ecalpulsecovariance::Quantity const & quantity() const { return m_quantity;}
83  ecalpulsecovariance::How const & how() const { return m_how;}
84  std::vector<int> const & which() const { return m_which;}
85 
86 
89  void set_which(std::vector<int> & i) { m_which.swap(i);}
90  };
91 
92 
93  template<>
94  class ValueExtractor<EcalPulseCovariances>: public BaseValueExtractor<EcalPulseCovariances> {
95  public:
96 
102  };
103  return fun[how];
104  }
105 
108  static What what() { return What();}
109 
112  : m_what(what)
113  {
114  // here one can make stuff really complicated... (select mean rms, 12,6,1)
115  // ask to make average on selected channels...
116  }
117 
118  void compute(Class const & it) override{
119  std::vector<float> res;
120  extractor(m_what.how())(it,m_what.quantity(),m_what.which(),res);
121  swap(res);
122  }
123 
124  private:
126 
127  };
128 
129 
130  template<>
133  std::stringstream ss;
134  EcalCondHeader header;
135  ss<<EcalPulseCovariancesXMLTranslator::dumpXML(header,object());
136  return ss.str();
137  } // dump
138 
139 
140  class EcalPulseCovariancesHelper: public EcalPyWrapperHelper<EcalPulseCovariance>{
141  public:
143  protected:
145  type_vValues getValues( const std::vector<EcalPulseCovariance> & vItems) override
146  {
147  type_vValues vValues(total_values);
148 
149  for(int s = 0; s<EcalPulseShape::TEMPLATESAMPLES; ++s) vValues[s].first = Form("sample_0.vs.sample_%d",s);
150 
151  //get info:
152  for(std::vector<EcalPulseCovariance>::const_iterator iItems = vItems.begin(); iItems != vItems.end(); ++iItems) {
153  for(int s = 0; s<EcalPulseShape::TEMPLATESAMPLES; ++s) vValues[s].second = iItems->val(THERAW,s);
154  }
155  return vValues;
156  }
157  };
158 
159  template<>
161  std::stringstream ss;
163  ss << helper.printBarrelsEndcaps(object().barrelItems(), object().endcapItems());
164  return ss.str();
165  } // summary
166 
167 
168  // return the real name of the file including extension...
169  template<>
171  std::string const &,
172  std::vector<int> const&,
173  std::vector<float> const& ) const {
174  gStyle->SetPalette(1);
175  // TCanvas canvas("CC map","CC map",840,600);
176  TCanvas canvas("CC map","CC map",800,1200);
177 
178  float xmi[3] = {0.0 , 0.22, 0.78};
179  float xma[3] = {0.22, 0.78, 1.00};
180  TPad*** pad = new TPad**[6];
181  for (int s = 0; s < 6; s++) {
182  pad[s] = new TPad*[3];
183  for (int obj = 0; obj < 3; obj++) {
184  float yma = 1.- (0.17 * s);
185  float ymi = yma - 0.15;
186  pad[s][obj] = new TPad(Form("p_%i_%i", obj, s),Form("p_%i_%i", obj, s),
187  xmi[obj], ymi, xma[obj], yma);
188  pad[s][obj]->Draw();
189  }
190  }
191 
192  const int kSides = 2;
193  const int kBarlRings = EBDetId::MAX_IETA;
194  const int kBarlWedges = EBDetId::MAX_IPHI;
195  const int kEndcWedgesX = EEDetId::IX_MAX;
196  const int kEndcWedgesY = EEDetId::IY_MAX;
197 
198  TH2F** barrel_s = new TH2F*[EcalPulseShape::TEMPLATESAMPLES];
199  TH2F** endc_p_s = new TH2F*[EcalPulseShape::TEMPLATESAMPLES];
200  TH2F** endc_m_s = new TH2F*[EcalPulseShape::TEMPLATESAMPLES];
201  for(int s = 0; s<EcalPulseShape::TEMPLATESAMPLES; ++s) {
202  barrel_s[s] = new TH2F(Form("EBs%i",s),Form("sample %i EB",s),360,0,360, 170, -85,85);
203  endc_p_s[s] = new TH2F(Form("EE+s%i",s),Form("sample %i EE+",s),100,1,101,100,1,101);
204  endc_m_s[s] = new TH2F(Form("EE-s%i",s),Form("sample %i EE-",s),100,1,101,100,1,101);
205  }
206 
207  for (int sign=0; sign < kSides; sign++) {
208  int thesign = sign==1 ? 1:-1;
209 
210  for (int ieta=0; ieta<kBarlRings; ieta++) {
211  for (int iphi=0; iphi<kBarlWedges; iphi++) {
212  EBDetId id((ieta+1)*thesign, iphi+1);
213  float y = -1 - ieta;
214  if(sign == 1) y = ieta;
215  for(int s = 0; s<EcalPulseShape::TEMPLATESAMPLES; ++s) {
216  barrel_s[s]->Fill(iphi, y, fabs(object()[id.rawId()].covval[THERAW][s]));
217  }
218  } // iphi
219  } // ieta
220 
221  for (int ix=0; ix<kEndcWedgesX; ix++) {
222  for (int iy=0; iy<kEndcWedgesY; iy++) {
223  if (! EEDetId::validDetId(ix+1,iy+1,thesign)) continue;
224  EEDetId id(ix+1,iy+1,thesign);
225  if (thesign==1) {
226  for(int s = 0; s<EcalPulseShape::TEMPLATESAMPLES; ++s)
227  endc_p_s[s]->Fill(ix+1,iy+1,fabs(object()[id.rawId()].covval[THERAW][s]));
228  }
229  else {
230  for(int s = 0; s<EcalPulseShape::TEMPLATESAMPLES; ++s)
231  endc_m_s[s]->Fill(ix+1,iy+1,fabs(object()[id.rawId()].covval[THERAW][s]));
232  }
233  } // iy
234  } // ix
235  } // side
236 
237  //canvas.cd(1);
238  TLine* l = new TLine(0., 0., 0., 0.);
239  l->SetLineWidth(1);
240  int ixSectorsEE[202] = {
241  62, 62, 61, 61, 60, 60, 59, 59, 58, 58, 56, 56, 46, 46, 44, 44, 43, 43, 42, 42,
242  41, 41, 40, 40, 41, 41, 42, 42, 43, 43, 44, 44, 46, 46, 56, 56, 58, 58, 59, 59,
243  60, 60, 61, 61, 62, 62, 0,101,101, 98, 98, 96, 96, 93, 93, 88, 88, 86, 86, 81,
244  81, 76, 76, 66, 66, 61, 61, 41, 41, 36, 36, 26, 26, 21, 21, 16, 16, 14, 14, 9,
245  9, 6, 6, 4, 4, 1, 1, 4, 4, 6, 6, 9, 9, 14, 14, 16, 16, 21, 21, 26,
246  26, 36, 36, 41, 41, 61, 61, 66, 66, 76, 76, 81, 81, 86, 86, 88, 88, 93, 93, 96,
247  96, 98, 98,101,101, 0, 62, 66, 66, 71, 71, 81, 81, 91, 91, 93, 0, 62, 66, 66,
248  91, 91, 98, 0, 58, 61, 61, 66, 66, 71, 71, 76, 76, 81, 81, 0, 51, 51, 0, 44,
249  41, 41, 36, 36, 31, 31, 26, 26, 21, 21, 0, 40, 36, 36, 11, 11, 4, 0, 40, 36,
250  36, 31, 31, 21, 21, 11, 11, 9, 0, 46, 46, 41, 41, 36, 36, 0, 56, 56, 61, 61, 66, 66};
251 
252  int iySectorsEE[202] = {
253  51, 56, 56, 58, 58, 59, 59, 60, 60, 61, 61, 62, 62, 61, 61, 60, 60, 59, 59, 58,
254  58, 56, 56, 46, 46, 44, 44, 43, 43, 42, 42, 41, 41, 40, 40, 41, 41, 42, 42, 43,
255  43, 44, 44, 46, 46, 51, 0, 51, 61, 61, 66, 66, 76, 76, 81, 81, 86, 86, 88, 88,
256  93, 93, 96, 96, 98, 98,101,101, 98, 98, 96, 96, 93, 93, 88, 88, 86, 86, 81, 81,
257  76, 76, 66, 66, 61, 61, 41, 41, 36, 36, 26, 26, 21, 21, 16, 16, 14, 14, 9, 9,
258  6, 6, 4, 4, 1, 1, 4, 4, 6, 6, 9, 9, 14, 14, 16, 16, 21, 21, 26, 26,
259  36, 36, 41, 41, 51, 0, 46, 46, 41, 41, 36, 36, 31, 31, 26, 26, 0, 51, 51, 56,
260  56, 61, 61, 0, 61, 61, 66, 66, 71, 71, 76, 76, 86, 86, 88, 0, 62,101, 0, 61,
261  61, 66, 66, 71, 71, 76, 76, 86, 86, 88, 0, 51, 51, 56, 56, 61, 61, 0, 46, 46,
262  41, 41, 36, 36, 31, 31, 26, 26, 0, 40, 31, 31, 16, 16, 6, 0, 40, 31, 31, 16, 16, 6};
263 
264  int ipad=0;
265  // plot only the measured ones, not the extrapolated
266  for(int s = 0; s<7; ++s) {
267  // plot only the extrapolated ones
268  // for(int s = 7; s<12; ++s) {
269  // do not plot the maximum sample, which is 1 by default
270  if(s==2) continue;
271  pad[ipad][0]->cd();
272  pad[ipad][0]->SetLogz(1);
273  endc_m_s[s]->SetStats(0);
274  endc_m_s[s]->SetMaximum(1.0);
275  endc_m_s[s]->SetMinimum(1.0e-7);
276  endc_m_s[s]->Draw("colz");
277  for ( int i=0; i<201; i=i+1) {
278  if ( (ixSectorsEE[i]!=0 || iySectorsEE[i]!=0) &&
279  (ixSectorsEE[i+1]!=0 || iySectorsEE[i+1]!=0) ) {
280  l->DrawLine(ixSectorsEE[i], iySectorsEE[i],
281  ixSectorsEE[i+1], iySectorsEE[i+1]);
282  l->SetLineWidth(0.2);
283  }
284  }
285  //canvas.cd(2);
286  pad[ipad][1]->cd();
287  pad[ipad][1]->SetLogz(1);
288  barrel_s[s]->SetStats(0);
289  barrel_s[s]->SetMaximum(1.0);
290  barrel_s[s]->SetMinimum(1.0e-7);
291  barrel_s[s]->Draw("colz");
292  for(int i = 0; i <17; i++) {
293  Double_t x = 20.+ (i *20);
294  l = new TLine(x,-85.,x,86.);
295  l->Draw();
296  l->SetLineWidth(0.2);
297  }
298  l = new TLine(0.,0.,360.,0.);
299  l->Draw();
300  //canvas.cd(3);
301  pad[ipad][2]->cd();
302  pad[ipad][2]->SetLogz(1);
303  endc_p_s[s]->SetStats(0);
304  endc_p_s[s]->SetMaximum(1.0);
305  endc_p_s[s]->SetMinimum(1.0e-7);
306  endc_p_s[s]->Draw("colz");
307  for ( int i=0; i<201; i=i+1) {
308  if ( (ixSectorsEE[i]!=0 || iySectorsEE[i]!=0) &&
309  (ixSectorsEE[i+1]!=0 || iySectorsEE[i+1]!=0) ) {
310  l->DrawLine(ixSectorsEE[i], iySectorsEE[i],
311  ixSectorsEE[i+1], iySectorsEE[i+1]);
312  l->SetLineWidth(0.2);
313  }
314  }
315  ipad++;
316  }
317 
318  canvas.SaveAs(filename.c_str());
319  return filename;
320  } // plot
321 }
322 
323 namespace condPython {
324  template<>
326  using namespace boost::python;
327  enum_<cond::ecalpulsecovariance::Quantity>("Quantity")
329  .value("sample_1",cond::ecalpulsecovariance::sample_1)
330  .value("sample_2",cond::ecalpulsecovariance::sample_2)
331  .value("sample_3",cond::ecalpulsecovariance::sample_3)
332  .value("sample_4",cond::ecalpulsecovariance::sample_4)
333  .value("sample_5",cond::ecalpulsecovariance::sample_5)
334  .value("sample_6",cond::ecalpulsecovariance::sample_6)
335  .value("sample_7",cond::ecalpulsecovariance::sample_7)
336  .value("sample_8",cond::ecalpulsecovariance::sample_8)
337  .value("sample_9",cond::ecalpulsecovariance::sample_9)
338  ;
339  enum_<cond::ecalpulsecovariance::How>("How")
341  .value("bySuperModule",cond::ecalpulsecovariance::bySuperModule)
342  .value("all",cond::ecalpulsecovariance::all)
343  ;
344 
346  class_<What>("What",init<>())
347  .def("set_quantity",&What::set_quantity)
348  .def("set_how",&What::set_how)
349  .def("set_which",&What::set_which)
350  .def("quantity",&What::quantity, return_value_policy<copy_const_reference>())
351  .def("how",&What::how, return_value_policy<copy_const_reference>())
352  .def("which",&What::which, return_value_policy<copy_const_reference>())
353  ;
354  }
355 }
356 
357 
358 
int i
Definition: DBlmapReader.cc:9
static std::string dumpXML(const EcalCondHeader &header, const EcalPulseCovariances &record)
std::string plot(std::string const &, std::string const &, std::vector< int > const &, std::vector< float > const &) const
ExtractWhat< Class > What
static const int kBarlRings
const Items & barrelItems() const
std::vector< std::pair< std::string, float > > type_vValues
double sign(double x)
static const int TEMPLATESAMPLES
#define PYTHON_WRAPPER(_class, _name)
static const int kSides
float val(int i, int j) const
ecalpulsecovariance::How const & how() const
boost::function< void(EcalPulseCovariances const &pulseshapes, Quantity q, std::vector< int > const &which, std::vector< float > &result)> PulseCovarianceExtractor
def canvas
Definition: svgfig.py:481
static const int kBarlWedges
U second(std::pair< T, U > const &p)
T x() const
Cartesian x coordinate.
std::string summary() const
static const int kEndcWedgesX
float average(EcalPulseCovariances const &pulseshapes, Quantity q)
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
void set_quantity(ecalpulsecovariance::Quantity i)
tuple result
Definition: query.py:137
void swap(std::vector< float > &v)
void extractAverage(EcalPulseCovariances const &pulseshapes, Quantity q, std::vector< int > const &, std::vector< float > &result)
static const int IX_MAX
Definition: EEDetId.h:302
static const int MAX_IPHI
Definition: EBDetId.h:144
ecalpulsecovariance::Quantity const & quantity() const
static bool validDetId(int crystal_ix, int crystal_iy, int iz)
Definition: EEDetId.h:248
static const int MAX_IETA
Definition: EBDetId.h:143
static ecalpulsecovariance::PulseCovarianceExtractor & extractor(ecalpulsecovariance::How how)
std::string printBarrelsEndcaps(const std::vector< T > &barrelItems, const std::vector< T > &endcapItems)
tuple filename
Definition: lut2db_cfg.py:20
static const int IY_MAX
Definition: EEDetId.h:306
void extractSingleChannel(EcalPulseCovariances const &pulseshapes, Quantity q, std::vector< int > const &which, std::vector< float > &result)
std::string dump() const
type_vValues getValues(const std::vector< EcalPulseCovariance > &vItems) override
void extractSuperModules(EcalPulseCovariances const &pulseshapes, Quantity q, std::vector< int > const &which, std::vector< float > &result)
static const int kEndcWedgesY