CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalTPGCrystalStatusPyWrapper.cc
Go to the documentation of this file.
4 #include "TROOT.h"
5 #include "TH2F.h"
6 #include "TCanvas.h"
7 #include "TStyle.h"
8 #include "TColor.h"
9 #include "TLine.h"
12 
15 
16 #include <string>
17 #include <sstream>
18 #include <algorithm>
19 #include <numeric>
20 #include <iterator>
21 #include <boost/ref.hpp>
22 #include <boost/bind.hpp>
23 #include <boost/function.hpp>
24 #include <boost/iterator/transform_iterator.hpp>
25 
26 #include <fstream>
27 
29 
30 namespace cond {
31 
32  namespace ecalcond {
33 
36 
38 
39  int bad(Items const & cont) {
40  return std::count_if(cont.begin(),cont.end(),
41  boost::bind(std::greater<int>(),
42  boost::bind(&value_type::getStatusCode,_1),0)
43  );
44  }
45 
46  void extractBarrel(EcalTPGCrystalStatus const & cont, std::vector<int> const &, std::vector<float> & result) {
47  result.resize(1);
48  result[0] = bad(cont.barrelItems());
49  }
50 
51  void extractEndcap(EcalTPGCrystalStatus const & cont, std::vector<int> const &, std::vector<float> & result) {
52  result.resize(1);
53  result[0] = bad(cont.endcapItems());
54  }
55  void extractAll(EcalTPGCrystalStatus const & cont, std::vector<int> const &, std::vector<float> & result) {
56  result.resize(1);
57  result[0] = bad(cont.barrelItems())+bad(cont.endcapItems());
58  }
59 
60  void extractSuperModules(EcalTPGCrystalStatus const & cont, std::vector<int> const & which, std::vector<float> & result) {
61  // bho...
62  }
63 
64  void extractSingleChannel(EcalTPGCrystalStatus const & cont, std::vector<int> const & which, std::vector<float> & result) {
65  result.reserve(which.size());
66  for (unsigned int i=0; i<which.size();i++) {
67  result.push_back(cont[which[i]].getStatusCode());
68  }
69  }
70 
71  typedef boost::function<void(EcalTPGCrystalStatus const & cont, std::vector<int> const & which, std::vector<float> & result)> CondExtractor;
72  } // namespace ecalcond
73 
74  template<>
76 
78  std::vector<int> m_which;
79 
80  ecalcond::How const & how() const { return m_how;}
81  std::vector<int> const & which() const { return m_which;}
82 
83  void set_how(ecalcond::How i) {m_how=i;}
84  void set_which(std::vector<int> & i) { m_which.swap(i);}
85  };
86 
87 
88  template<>
89  class ValueExtractor<EcalTPGCrystalStatus>: public BaseValueExtractor<EcalTPGCrystalStatus> {
90  public:
91 
93  static ecalcond::CondExtractor fun[3] = {
97  };
98  return fun[how];
99  }
100 
103  static What what() { return What();}
104 
107  : m_what(what)
108  {
109  // here one can make stuff really complicated...
110  }
111 
112  void compute(Class const & it) override{
113  std::vector<float> res;
114  extractor(m_what.how())(it,m_what.which(),res);
115  swap(res);
116  }
117 
118  private:
120  };
121 
122 
123  template<>
125  std::stringstream ss;
128  return ss.str();
129  }
130 
131  class EcalTPGCrystalStatusHelper: public EcalPyWrapperHelper<EcalTPGCrystalStatusCode>{
132  public:
133  //change me
135  protected:
136 
137  //change me
139 
140  type_vValues getValues( const std::vector<EcalObject> & vItems) override
141  {
142  //change me
143  //unsigned int totalValues = 2;
144 
145  type_vValues vValues(total_values);
146 
147  //change us
148  vValues[0].first = "[0]StatusCode";
149 
150 
151  vValues[0].second = .0;
152 
153  //get info:
154  unsigned int shift = 0, mask = 1;
155  unsigned int statusCode;
156  for(std::vector<EcalObject>::const_iterator iItems = vItems.begin(); iItems != vItems.end(); ++iItems){
157  //change us
158  statusCode = iItems->getStatusCode();
159  for (shift = 0; shift < total_values; ++shift){
160  mask = 1 << (shift);
161  //std::cout << "; statuscode: " << statusCode;
162  if (statusCode & mask){
163  vValues[shift].second += 1;
164  }
165  }
166  }
167  return vValues;
168  }
169  };
170 
171  template<>
173  std::stringstream ss;
174  // EcalTPGCrystalStatusHelper helper;
175  // ss << helper.printBarrelsEndcaps(object().barrelItems(), object().endcapItems());
176  const int kSides = 2;
177  const int kBarlRings = EBDetId::MAX_IETA;
178  const int kBarlWedges = EBDetId::MAX_IPHI;
179  const int kEndcWedgesX = EEDetId::IX_MAX;
180  const int kEndcWedgesY = EEDetId::IY_MAX;
181 
182  int EB[2] = {0, 0}, EE[2] = {0, 0};
183  for (int sign=0; sign < kSides; sign++) {
184  int thesign = sign==1 ? 1:-1;
185 
186  for (int ieta=0; ieta<kBarlRings; ieta++) {
187  for (int iphi=0; iphi<kBarlWedges; iphi++) {
188  EBDetId id((ieta+1)*thesign, iphi+1);
189  if(object()[id.rawId()].getStatusCode() > 0) EB[sign]++;
190  } // iphi
191  } // ieta
192 
193  for (int ix=0; ix<kEndcWedgesX; ix++) {
194  for (int iy=0; iy<kEndcWedgesY; iy++) {
195  if (! EEDetId::validDetId(ix+1,iy+1,thesign)) continue;
196  EEDetId id(ix+1,iy+1,thesign);
197  if(object()[id.rawId()].getStatusCode() > 0) EE[sign]++;
198  } // iy
199  } // ix
200  } // side
201  ss << " number of masked Crystals EB- " << EB[0] << " EB+ " << EB[1]
202  << " EE- " << EE[0] << " EE+ " << EE[1] << std::endl;
203  return ss.str();
204  }
205 
206  template<>
208  std::string const &,
209  std::vector<int> const&,
210  std::vector<float> const& ) const {
211  gStyle->SetPalette(1);
212  const int TOTAL_PADS = 3;
213 
214  // TCanvas canvas("CC map","CC map",700,800);
215  Double_t w = 600;
216  Double_t h = 600;
217  TCanvas canvas("c", "c", w, h);
218  canvas.SetWindowSize(w + (w - canvas.GetWw()), h + (h - canvas.GetWh()));
219  // const float IMG_SIZE = 1.5;
220  // TCanvas canvas("CC map","CC map",800*IMG_SIZE, 200 * IMG_SIZE);//800, 1200
221 
222  float xmi[3] = {0.0, 0.0, 0.5};
223  float xma[3] = {1.0, 0.5, 1.0};
224  float ymi[3] = {0.5, 0.0, 0.0};
225  float yma[3] = {1.0, 0.5, 0.5};
226 
227 
228  TPad** pad = new TPad*[TOTAL_PADS];
229  for (int obj = 0; obj < TOTAL_PADS; obj++) {
230  pad[obj] = new TPad(Form("p_%i", obj),Form("p_%i", obj),
231  xmi[obj], ymi[obj], xma[obj], yma[obj]);
232  pad[obj]->Draw();
233  }
234 
235  const int kSides = 2;
236  const int kBarlRings = EBDetId::MAX_IETA;
237  const int kBarlWedges = EBDetId::MAX_IPHI;
238  const int kEndcWedgesX = EEDetId::IX_MAX;
239  const int kEndcWedgesY = EEDetId::IY_MAX;
240 
241  TH2F* barrel = new TH2F("EB","EB TPG Crystal Status", 360,0,360, 170, -85,85);
242  TH2F* endc_m = new TH2F("EEm","EE- TPG Crystal Status",100,1,101,100,1,101);
243  TH2F* endc_p = new TH2F("EEp","EE+ TPG Crystal Status",100,1,101,100,1,101);
244 
245  for (int sign=0; sign < kSides; sign++) {
246  int thesign = sign==1 ? 1:-1;
247 
248  for (int ieta=0; ieta<kBarlRings; ieta++) {
249  for (int iphi=0; iphi<kBarlWedges; iphi++) {
250  EBDetId id((ieta+1)*thesign, iphi+1);
251  float y = -1 - ieta;
252  if(sign == 1) y = ieta;
253  barrel->Fill(iphi, y, object()[id.rawId()].getStatusCode());
254  } // iphi
255  } // ieta
256 
257  for (int ix=0; ix<kEndcWedgesX; ix++) {
258  for (int iy=0; iy<kEndcWedgesY; iy++) {
259  if (! EEDetId::validDetId(ix+1,iy+1,thesign)) continue;
260  EEDetId id(ix+1,iy+1,thesign);
261  if (thesign==1) {
262  endc_p->Fill(ix+1,iy+1,object()[id.rawId()].getStatusCode());
263  }
264  else{
265  endc_m->Fill(ix+1,iy+1,object()[id.rawId()].getStatusCode());
266  }
267  } // iy
268  } // ix
269  } // side
270 
271  TLine* l = new TLine(0., 0., 0., 0.);
272  l->SetLineWidth(1);
273  int ixSectorsEE[202] = {
274  62, 62, 61, 61, 60, 60, 59, 59, 58, 58, 56, 56, 46, 46, 44, 44, 43, 43, 42, 42,
275  41, 41, 40, 40, 41, 41, 42, 42, 43, 43, 44, 44, 46, 46, 56, 56, 58, 58, 59, 59,
276  60, 60, 61, 61, 62, 62, 0,101,101, 98, 98, 96, 96, 93, 93, 88, 88, 86, 86, 81,
277  81, 76, 76, 66, 66, 61, 61, 41, 41, 36, 36, 26, 26, 21, 21, 16, 16, 14, 14, 9,
278  9, 6, 6, 4, 4, 1, 1, 4, 4, 6, 6, 9, 9, 14, 14, 16, 16, 21, 21, 26,
279  26, 36, 36, 41, 41, 61, 61, 66, 66, 76, 76, 81, 81, 86, 86, 88, 88, 93, 93, 96,
280  96, 98, 98,101,101, 0, 62, 66, 66, 71, 71, 81, 81, 91, 91, 93, 0, 62, 66, 66,
281  91, 91, 98, 0, 58, 61, 61, 66, 66, 71, 71, 76, 76, 81, 81, 0, 51, 51, 0, 44,
282  41, 41, 36, 36, 31, 31, 26, 26, 21, 21, 0, 40, 36, 36, 11, 11, 4, 0, 40, 36,
283  36, 31, 31, 21, 21, 11, 11, 9, 0, 46, 46, 41, 41, 36, 36, 0, 56, 56, 61, 61, 66, 66};
284 
285  int iySectorsEE[202] = {
286  51, 56, 56, 58, 58, 59, 59, 60, 60, 61, 61, 62, 62, 61, 61, 60, 60, 59, 59, 58,
287  58, 56, 56, 46, 46, 44, 44, 43, 43, 42, 42, 41, 41, 40, 40, 41, 41, 42, 42, 43,
288  43, 44, 44, 46, 46, 51, 0, 51, 61, 61, 66, 66, 76, 76, 81, 81, 86, 86, 88, 88,
289  93, 93, 96, 96, 98, 98,101,101, 98, 98, 96, 96, 93, 93, 88, 88, 86, 86, 81, 81,
290  76, 76, 66, 66, 61, 61, 41, 41, 36, 36, 26, 26, 21, 21, 16, 16, 14, 14, 9, 9,
291  6, 6, 4, 4, 1, 1, 4, 4, 6, 6, 9, 9, 14, 14, 16, 16, 21, 21, 26, 26,
292  36, 36, 41, 41, 51, 0, 46, 46, 41, 41, 36, 36, 31, 31, 26, 26, 0, 51, 51, 56,
293  56, 61, 61, 0, 61, 61, 66, 66, 71, 71, 76, 76, 86, 86, 88, 0, 62,101, 0, 61,
294  61, 66, 66, 71, 71, 76, 76, 86, 86, 88, 0, 51, 51, 56, 56, 61, 61, 0, 46, 46,
295  41, 41, 36, 36, 31, 31, 26, 26, 0, 40, 31, 31, 16, 16, 6, 0, 40, 31, 31, 16, 16, 6};
296 
297  pad[0]->cd();
298  barrel->SetStats(0);
299  barrel->Draw("colz");
300  for(int i = 0; i <17; i++) {
301  Double_t x = 20.+ (i *20);
302  l = new TLine(x,-85.,x,86.);
303  l->Draw();
304  }
305  l = new TLine(0.,0.,360.,0.);
306  l->Draw();
307 
308  pad[1]->cd();
309  endc_m->SetStats(0);
310  endc_m->Draw("col");
311  for ( int i=0; i<201; i=i+1) {
312  if ( (ixSectorsEE[i]!=0 || iySectorsEE[i]!=0) &&
313  (ixSectorsEE[i+1]!=0 || iySectorsEE[i+1]!=0) ) {
314  l->DrawLine(ixSectorsEE[i], iySectorsEE[i],
315  ixSectorsEE[i+1], iySectorsEE[i+1]);
316  l->SetLineWidth(0.2);
317  }
318  }
319 
320  pad[2]->cd();
321  endc_p->SetStats(0);
322  endc_p->Draw("col");
323  for ( int i=0; i<201; i=i+1) {
324  if ( (ixSectorsEE[i]!=0 || iySectorsEE[i]!=0) &&
325  (ixSectorsEE[i+1]!=0 || iySectorsEE[i+1]!=0) ) {
326  l->DrawLine(ixSectorsEE[i], iySectorsEE[i],
327  ixSectorsEE[i+1], iySectorsEE[i+1]);
328  }
329  }
330 
331  canvas.SaveAs(filename.c_str());
332  return filename;
333  } // plot
334 }
335 
void extractSuperModules(Container const &cont, std::vector< int > const &which, std::vector< float > &result)
int i
Definition: DBlmapReader.cc:9
def which
Definition: eostools.py:335
std::string plot(std::string const &, std::string const &, std::vector< int > const &, std::vector< float > const &) const
tuple cont
load Luminosity info ##
Definition: generateEDF.py:622
const double w
Definition: UKUtility.cc:23
ExtractWhat< Class > What
static std::string dumpXML(const EcalCondHeader &header, const EcalTPGCrystalStatus &record)
static const int kBarlRings
const Items & barrelItems() const
std::vector< std::pair< std::string, float > > type_vValues
double sign(double x)
Code getStatusCode() const
return decoded status
boost::function< void(Container const &cont, std::vector< int > const &which, std::vector< float > &result)> CondExtractor
#define PYTHON_WRAPPER(_class, _name)
static const int kSides
def canvas
Definition: svgfig.py:481
int bad(Items const &cont)
static const int kBarlWedges
T x() const
Cartesian x coordinate.
std::string summary() const
static const int kEndcWedgesX
tuple result
Definition: query.py:137
type_vValues getValues(const std::vector< EcalObject > &vItems) override
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
Container::value_type value_type
void swap(std::vector< float > &v)
static const int IX_MAX
Definition: EEDetId.h:302
static const int MAX_IPHI
Definition: EBDetId.h:144
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 ecalcond::CondExtractor & extractor(ecalcond::How how)
void extractSingleChannel(Container const &cont, std::vector< int > const &which, std::vector< float > &result)
tuple filename
Definition: lut2db_cfg.py:20
static unsigned int const shift
static const int IY_MAX
Definition: EEDetId.h:306
void extractBarrel(Container const &cont, std::vector< int > const &, std::vector< float > &result)
const Items & endcapItems() const
std::string dump() const
void extractEndcap(Container const &cont, std::vector< int > const &, std::vector< float > &result)
void extractAll(Container const &cont, std::vector< int > const &, std::vector< float > &result)
static const int kEndcWedgesY