CMS 3D CMS Logo

FWECALCaloDataDetailViewBuilder.cc
Go to the documentation of this file.
1 // FIXME - needed to set fixed eta-phi limits. Without the
2 // visible area may change widely depending on energy
3 // deposition availability
4 
5 #include "TEveCaloData.h"
6 #include "TEveViewer.h"
7 #include "TEveCalo.h"
8 #include "TAxis.h"
9 #include "TMath.h"
10 #include "THLimitsFinder.h"
11 #include "TLatex.h"
12 
18 
22 
23 #include "TGeoMatrix.h"
24 #include "TEveTrans.h"
25 
26 #include <utility>
27 
29  // get the hits from the event
30 
32  const EcalRecHitCollection* hits = nullptr;
33 
34  if (fabs(m_eta) < 1.5) {
35  try {
36  edm::InputTag tag("ecalRecHit", "EcalRecHitsEB");
37  m_event->getByLabel(tag, handle_hits);
38  if (handle_hits.isValid()) {
39  hits = &*handle_hits;
40  }
41  } catch (...) {
42  fwLog(fwlog::kWarning) << "FWECALCaloDataDetailViewBuilder::build():: Failed to access EcalRecHitsEB collection."
43  << std::endl;
44  }
45  if (!handle_hits.isValid()) {
46  try {
47  edm::InputTag tag("reducedEcalRecHitsEB");
48  m_event->getByLabel(tag, handle_hits);
49  if (handle_hits.isValid()) {
50  hits = &*handle_hits;
51  }
52 
53  } catch (...) {
55  << "FWECALCaloDataDetailViewBuilder::build():: Failed to access reducedEcalRecHitsEB collection."
56  << std::endl;
57  }
58  }
59  } else {
60  try {
61  edm::InputTag tag("ecalRecHit", "EcalRecHitsEE");
62  m_event->getByLabel(tag, handle_hits);
63  if (handle_hits.isValid())
64  hits = &*handle_hits;
65  } catch (...) {
66  fwLog(fwlog::kWarning) << "FWECALCaloDataDetailViewBuilder::build():: Failed to access ecalRecHitsEE collection."
67  << std::endl;
68  }
69 
70  if (!handle_hits.isValid()) {
71  try {
72  edm::InputTag tag("reducedEcalRecHitsEE");
73  m_event->getByLabel(tag, handle_hits);
74  if (handle_hits.isValid()) {
75  hits = &*handle_hits;
76  }
77 
78  } catch (...) {
80  << "FWECALCaloDataDetailViewBuilder::build():: Failed to access reducedEcalRecHitsEE collection."
81  << std::endl;
82  }
83  }
84  }
85 
86  // data
87  TEveCaloDataVec* data = new TEveCaloDataVec(1 + m_colors.size());
88  data->RefSliceInfo(0).Setup("hits (not clustered)", 0.0, m_defaultColor);
89  for (size_t i = 0; i < m_colors.size(); ++i) {
90  data->RefSliceInfo(i + 1).Setup("hits (clustered)", 0.0, m_colors[i]);
91  }
92 
93  if (handle_hits.isValid()) {
94  fillData(hits, data, xyEE);
95  }
96 
97  // axis
98  Double_t etaMin(0), etaMax(0), phiMin(0), phiMax(0);
99  if (data->Empty()) {
100  fwLog(fwlog::kWarning) << "FWECALCaloDataDetailViewBuilder::build():: No hits found in "
101  << Form("FWECALCaloDataDetailViewBuilder::build():: No hits found in eta[%f] phi[%f] region",
102  m_eta,
103  m_phi)
104  << ".\n";
105 
106  // add dummy background
107  float x = m_size * TMath::DegToRad();
108  if (fabs(m_eta) < 1.5 || xyEE == false) {
109  etaMin = m_eta - x;
110  etaMax = m_eta + x;
111  phiMin = m_phi - x;
112  phiMax = m_phi + x;
113  data->AddTower(etaMin, etaMax, phiMin, phiMax);
114  } else {
115  float theta = TEveCaloData::EtaToTheta(m_eta);
116  float r = TMath::Tan(theta) * 290;
117  phiMin = r * TMath::Cos(m_phi - x) - 300;
118  phiMax = r * TMath::Cos(m_phi + x) + 300;
119  etaMin = r * TMath::Sin(m_phi - x) - 300;
120  etaMax = r * TMath::Sin(m_phi + x) + 300;
121  data->AddTower(TMath::Min(etaMin, etaMax),
125  }
126  data->FillSlice(0, 0.1);
127  }
128 
129  TAxis* eta_axis = nullptr;
130  TAxis* phi_axis = nullptr;
131  data->GetEtaLimits(etaMin, etaMax);
132  data->GetPhiLimits(phiMin, phiMax);
133  // printf("data rng %f %f %f %f\n",etaMin, etaMax, phiMin, phiMax );
134 
135  if (fabs(m_eta) > 1.5 && xyEE) {
136  eta_axis = new TAxis(10, etaMin, etaMax);
137  phi_axis = new TAxis(10, phiMin, phiMax);
138  eta_axis->SetTitle("X[cm]");
139  phi_axis->SetTitle("Y[cm]");
140  phi_axis->SetTitleSize(0.05);
141  eta_axis->SetTitleSize(0.05);
142  } else {
143  std::vector<double> etaBinsWithinLimits;
144  etaBinsWithinLimits.push_back(etaMin);
145  for (unsigned int i = 0; i < 83; ++i)
147  etaBinsWithinLimits.push_back(fw3dlego::xbins[i]);
148  etaBinsWithinLimits.push_back(etaMax);
149 
150  std::vector<double> phiBinsWithinLimits;
151  phiBinsWithinLimits.push_back(phiMin);
152  for (double phi = -M_PI; phi < M_PI; phi += M_PI / 36)
153  if (phi > phiMin && phi < phiMax)
154  phiBinsWithinLimits.push_back(phi);
155  phiBinsWithinLimits.push_back(phiMax);
156 
157  eta_axis = new TAxis((int)etaBinsWithinLimits.size() - 1, &etaBinsWithinLimits[0]);
158  phi_axis = new TAxis((int)phiBinsWithinLimits.size() - 1, &phiBinsWithinLimits[0]);
159 
160  eta_axis->SetTitleFont(122);
161  eta_axis->SetTitle("h");
162  eta_axis->SetTitleSize(0.07);
163  phi_axis->SetTitleFont(122);
164  phi_axis->SetTitle("f");
165  phi_axis->SetTitleSize(0.07);
166  }
167  eta_axis->SetNdivisions(510);
168  phi_axis->SetNdivisions(510);
169  data->SetEtaBins(eta_axis);
170  data->SetPhiBins(phi_axis);
171  return data;
172 }
173 
174 //_______________________________________________________________
176  TEveCaloData* data = buildCaloData(true);
177 
178  double etaMin, etaMax, phiMin, phiMax;
179  data->GetEtaLimits(etaMin, etaMax);
180  data->GetPhiLimits(phiMin, phiMax);
181 
182  // lego
183  TEveCaloLego* lego = new TEveCaloLego(data);
184  lego->SetDrawNumberCellPixels(100);
185  // scale and translate to real world coordinates
186  lego->SetEta(etaMin, etaMax);
187  lego->SetPhiWithRng((phiMin + phiMax) * 0.5, (phiMax - phiMin) * 0.5); // phi range = 2* phiOffset
188  Double_t legoScale = ((etaMax - etaMin) < (phiMax - phiMin)) ? (etaMax - etaMin) : (phiMax - phiMin);
189  lego->InitMainTrans();
190  lego->RefMainTrans().SetScale(legoScale, legoScale, legoScale * 0.5);
191  lego->RefMainTrans().SetPos((etaMax + etaMin) * 0.5, (phiMax + phiMin) * 0.5, 0);
192  lego->SetAutoRebin(kFALSE);
193  lego->Set2DMode(TEveCaloLego::kValSizeOutline);
194  lego->SetName("ECALDetail Lego");
195  return lego;
196 }
197 
198 void FWECALCaloDataDetailViewBuilder::setColor(Color_t color, const std::vector<DetId>& detIds) {
199  m_colors.push_back(color);
200 
201  // get the slice for this group of detIds
202  // note that the zeroth slice is the default one (all else)
203  int slice = m_colors.size();
204  // take a note of which slice these detids are going to go into
205  for (size_t i = 0; i < detIds.size(); ++i)
206  m_detIdsToColor[detIds[i]] = slice;
207 }
208 
210  std::vector<DetId> clusterDetIds;
211  const std::vector<std::pair<DetId, float> >& hitsAndFractions = cluster.hitsAndFractions();
212  for (size_t j = 0; j < hitsAndFractions.size(); ++j) {
213  clusterDetIds.push_back(hitsAndFractions[j].first);
214  }
215 
216  setColor(color, clusterDetIds);
217 }
218 
219 void FWECALCaloDataDetailViewBuilder::showSuperClusters(Color_t color1, Color_t color2) {
220  // get the superclusters from the event
222 
223  if (fabs(m_eta) < 1.5) {
224  try {
225  m_event->getByLabel(edm::InputTag("correctedHybridSuperClusters"), collection);
226  } catch (...) {
227  fwLog(fwlog::kWarning) << "no barrel superclusters are available" << std::endl;
228  }
229  } else {
230  try {
231  m_event->getByLabel(edm::InputTag("correctedMulti5x5SuperClustersWithPreshower"), collection);
232  } catch (...) {
233  fwLog(fwlog::kWarning) << "no endcap superclusters are available" << std::endl;
234  }
235  }
236  if (collection.isValid()) {
237  unsigned int colorIndex = 0;
238  // sort clusters in eta so neighboring clusters have distinct colors
239  reco::SuperClusterCollection sorted = *collection.product();
240  std::sort(sorted.begin(), sorted.end(), superClusterEtaLess);
241  for (size_t i = 0; i < sorted.size(); ++i) {
242  if (!(fabs(sorted[i].eta() - m_eta) < (m_size * 0.0172) && fabs(sorted[i].phi() - m_phi) < (m_size * 0.0172)))
243  continue;
244 
245  if (colorIndex % 2 == 0)
246  showSuperCluster(sorted[i], color1);
247  else
248  showSuperCluster(sorted[i], color2);
249  ++colorIndex;
250  }
251  }
252 }
253 
254 void FWECALCaloDataDetailViewBuilder::fillData(const EcalRecHitCollection* hits, TEveCaloDataVec* data, bool xyEE) {
255  const float barrelCR = m_size * 0.0172; // barrel cell range
256 
257  // loop on all the detids
258  for (EcalRecHitCollection::const_iterator k = hits->begin(), kEnd = hits->end(); k != kEnd; ++k) {
259  // get reco geometry
260  double centerEta = 0;
261  double centerPhi = 0;
262  const float* points = m_geom->getCorners(k->id().rawId());
263  if (points != nullptr) {
264  TEveVector v;
265  int j = 0;
266  for (int i = 0; i < 8; ++i) {
267  v += TEveVector(points[j], points[j + 1], points[j + 2]);
268  j += 3;
269  }
270  centerEta = v.Eta();
271  centerPhi = v.Phi();
272  } else
273  fwLog(fwlog::kInfo) << "cannot get geometry for DetId: " << k->id().rawId() << ". Ignored.\n";
274 
275  double size = k->energy() / cosh(centerEta);
276 
277  // check what slice to put in
278  int slice = 0;
279  std::map<DetId, int>::const_iterator itr = m_detIdsToColor.find(k->id());
280  if (itr != m_detIdsToColor.end())
281  slice = itr->second;
282 
283  // if in the EB
284  if (k->id().subdetId() == EcalBarrel || xyEE == false) {
285  // do phi wrapping
286  if (centerPhi > m_phi + M_PI)
287  centerPhi -= 2 * M_PI;
288  if (centerPhi < m_phi - M_PI)
289  centerPhi += 2 * M_PI;
290 
291  // check if the hit is in the window to be drawn
292  if (!(fabs(centerEta - m_eta) < barrelCR && fabs(centerPhi - m_phi) < barrelCR))
293  continue;
294 
295  double minEta(10), maxEta(-10), minPhi(4), maxPhi(-4);
296  if (points != nullptr) {
297  // calorimeter crystalls have slightly non-symetrical form in eta-phi projection
298  // so if we simply get the largest eta and phi, cells will overlap
299  // therefore we get a smaller eta-phi range representing the inner square
300  // we also should use only points from the inner face of the crystal, since
301  // non-projecting direction of crystals leads to large shift in eta on outter
302  // face.
303  int j = 0;
304  float eps = 0.005;
305  for (unsigned int i = 0; i < 8; ++i) {
306  TEveVector crystal(points[j], points[j + 1], points[j + 2]);
307  j += 3;
308  double eta = crystal.Eta();
309  double phi = crystal.Phi();
310  if (((k->id().subdetId() == EcalBarrel) && (crystal.Perp() > 135)) ||
311  ((k->id().subdetId() == EcalEndcap) && (crystal.Perp() > 155)))
312  continue;
313  if (minEta - eta > eps)
314  minEta = eta;
315  if (eta - minEta > 0 && eta - minEta < eps)
316  minEta = eta;
317  if (eta - maxEta > eps)
318  maxEta = eta;
319  if (maxEta - eta > 0 && maxEta - eta < eps)
320  maxEta = eta;
321  if (minPhi - phi > eps)
322  minPhi = phi;
323  if (phi - minPhi > 0 && phi - minPhi < eps)
324  minPhi = phi;
325  if (phi - maxPhi > eps)
326  maxPhi = phi;
327  if (maxPhi - phi > 0 && maxPhi - phi < eps)
328  maxPhi = phi;
329  }
330  } else {
331  double delta = 0.0172 * 0.5;
332  minEta = centerEta - delta;
333  maxEta = centerEta + delta;
334  minPhi = centerPhi - delta;
335  maxPhi = centerPhi + delta;
336  }
337  if (minPhi >= (m_phi - barrelCR) && maxPhi <= (m_phi + barrelCR) && minEta >= (m_eta - barrelCR) &&
338  maxEta <= (m_eta + barrelCR)) {
339  // printf("add %f %f %f %f \n",minEta, maxEta, minPhi, maxPhi );
340  data->AddTower(minEta, maxEta, minPhi, maxPhi);
341  data->FillSlice(slice, size);
342  }
343  // otherwise in the EE
344  } else if (k->id().subdetId() == EcalEndcap) {
345  // check if the hit is in the window to be drawn
346  double crystalSize = m_size * 0.0172;
347  if (!(fabs(centerEta - m_eta) < (crystalSize) && fabs(centerPhi - m_phi) < (crystalSize)))
348  continue;
349 
350  if (points != nullptr) {
351  double minX(9999), maxX(-9999), minY(9999), maxY(-9999);
352  int j = 0;
353  for (unsigned int i = 0; i < 8; ++i) {
354  TEveVector crystal(points[j], points[j + 1], points[j + 2]);
355  j += 3;
356  double x = crystal.fX;
357  double y = crystal.fY;
358  if (fabs(crystal.fZ) > 330)
359  continue;
360  if (minX - x > 0.01)
361  minX = x;
362  if (x - maxX > 0.01)
363  maxX = x;
364  if (minY - y > 0.01)
365  minY = y;
366  if (y - maxY > 0.01)
367  maxY = y;
368  }
369  data->AddTower(minX, maxX, minY, maxY);
370  // printf("EE add %f %f %f %f \n",minX, maxX, minY, maxY );
371  }
372  data->FillSlice(slice, size);
373  }
374  } // end loop on hits
375 
376  data->DataChanged();
377 }
378 
380  double x0, double y0, Color_t clustered1, Color_t clustered2, Color_t supercluster) {
381  Double_t fontsize = 0.07;
382  TLatex* latex = new TLatex();
383  Double_t x = x0;
384  Double_t y = y0;
385  Double_t boxH = 0.25 * fontsize;
386  Double_t yStep = 0.04;
387 
388  y -= yStep;
389  latex->DrawLatex(x, y, "Energy types:");
390  y -= yStep;
391 
392  Double_t pos[4];
393  pos[0] = x + 0.05;
394  pos[2] = x + 0.20;
395 
396  pos[1] = y;
397  pos[3] = pos[1] + boxH;
399  latex->DrawLatex(x + 0.25, y, "unclustered");
400  y -= yStep;
401  if (clustered1 < 0)
402  return y;
403 
404  pos[1] = y;
405  pos[3] = pos[1] + boxH;
407  latex->DrawLatex(x + 0.25, y, "clustered");
408  y -= yStep;
409  if (clustered2 < 0)
410  return y;
411 
412  pos[1] = y;
413  pos[3] = pos[1] + boxH;
415  latex->DrawLatex(x + 0.25, y, "clustered");
416  y -= yStep;
417  if (supercluster < 0)
418  return y;
419 
420  pos[1] = y;
421  pos[3] = pos[1] + boxH;
422  FWDetailViewBase::drawCanvasBox(pos, supercluster);
423  latex->DrawLatex(x + 0.25, y, "super-cluster");
424  y -= yStep;
425 
426  return y;
427 }
size
Write out results.
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
const double xbins[]
void showSuperCluster(const reco::SuperCluster &cluster, Color_t color=kYellow)
std::vector< EcalRecHit >::const_iterator const_iterator
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:97
const float * getCorners(unsigned int id) const
Definition: FWGeometry.cc:439
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
static bool superClusterEtaLess(const reco::CaloCluster &lhs, const reco::CaloCluster &rhs)
double makeLegend(double x0=0.02, double y0=0.95, Color_t clustered1=kGreen+1, Color_t clustered2=kTeal, Color_t supercluster=kYellow)
void setColor(Color_t color, const std::vector< DetId > &detIds)
#define M_PI
void fillData(const EcalRecHitCollection *hits, TEveCaloDataVec *data, bool xyEE)
#define fwLog(_level_)
Definition: fwLog.h:45
bool isValid() const
Definition: HandleBase.h:70
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:80
void showSuperClusters(Color_t color1=kGreen+2, Color_t color2=kTeal)
Geom::Theta< T > theta() const
if(threadIdxLocalY==0 &&threadIdxLocalX==0)
static void drawCanvasBox(Double_t *pos, Color_t fillCol, Int_t fillType=0, bool bg=kTRUE)