CMS 3D CMS Logo

FWCaloClusterProxyBuilder.cc
Go to the documentation of this file.
7 
8 #include "TEveBoxSet.h"
9 #include "TEveStraightLineSet.h"
10 
11 class FWCaloClusterProxyBuilder : public FWHeatmapProxyBuilderTemplate<reco::CaloCluster> {
12 public:
14  ~FWCaloClusterProxyBuilder(void) override {}
15 
17 
18  FWCaloClusterProxyBuilder(const FWCaloClusterProxyBuilder &) = delete; // stop default
19  const FWCaloClusterProxyBuilder &operator=(const FWCaloClusterProxyBuilder &) = delete; // stop default
20 
21 private:
24  long layer;
26  bool heatmap;
27  bool z_plus;
28  bool z_minus;
30 
31  void setItem(const FWEventItem *iItem) override;
32 
33  void build(const FWEventItem *iItem, TEveElementList *product, const FWViewContext *vc) override;
34  void build(const reco::CaloCluster &iData,
35  unsigned int iIndex,
36  TEveElement &oItemHolder,
37  const FWViewContext *) override;
38 };
39 
42  if (iItem) {
43  iItem->getConfig()->assertParam("Cluster(0)/RecHit(1)", false);
44  iItem->getConfig()->assertParam("EnableTimeFilter", false);
45  iItem->getConfig()->assertParam("TimeLowerBound(ns)", 0.01, 0.0, 75.0);
46  iItem->getConfig()->assertParam("TimeUpperBound(ns)", 0.01, 0.0, 75.0);
47  }
48 }
49 
50 void FWCaloClusterProxyBuilder::build(const FWEventItem *iItem, TEveElementList *product, const FWViewContext *vc) {
51  iItem->getEvent()->getByLabel(edm::InputTag("hgcalLayerClusters", "timeLayerCluster"), TimeValueMapHandle);
53  timeLowerBound = std::min(item()->getConfig()->value<double>("TimeLowerBound(ns)"),
54  item()->getConfig()->value<double>("TimeUpperBound(ns)"));
55  timeUpperBound = std::max(item()->getConfig()->value<double>("TimeLowerBound(ns)"),
56  item()->getConfig()->value<double>("TimeUpperBound(ns)"));
57  } else {
58  std::cerr << "Warning: couldn't locate 'timeLayerCluster' ValueMap in root file." << std::endl;
59  }
60 
61  layer = item()->getConfig()->value<long>("Layer");
62  saturation_energy = item()->getConfig()->value<double>("EnergyCutOff");
63  heatmap = item()->getConfig()->value<bool>("Heatmap");
64  z_plus = item()->getConfig()->value<bool>("Z+");
65  z_minus = item()->getConfig()->value<bool>("Z-");
66  enableTimeFilter = item()->getConfig()->value<bool>("EnableTimeFilter");
67 
68  FWHeatmapProxyBuilderTemplate::build(iItem, product, vc);
69 }
70 
72  unsigned int iIndex,
73  TEveElement &oItemHolder,
74  const FWViewContext *) {
76  const float time = TimeValueMapHandle->get(iIndex).first;
77  if (time < timeLowerBound || time > timeUpperBound)
78  return;
79  }
80 
81  std::vector<std::pair<DetId, float>> clusterDetIds = iData.hitsAndFractions();
82 
83  bool h_hex(false);
84  TEveBoxSet *hex_boxset = new TEveBoxSet();
85  if (!heatmap)
86  hex_boxset->UseSingleColor();
87  hex_boxset->SetPickable(true);
88  hex_boxset->Reset(TEveBoxSet::kBT_Hex, true, 64);
89  hex_boxset->SetAntiFlick(true);
90 
91  bool h_box(false);
92  TEveBoxSet *boxset = new TEveBoxSet();
93  if (!heatmap)
94  boxset->UseSingleColor();
95  boxset->SetPickable(true);
96  boxset->Reset(TEveBoxSet::kBT_FreeBox, true, 64);
97  boxset->SetAntiFlick(true);
98 
99  for (std::vector<std::pair<DetId, float>>::iterator it = clusterDetIds.begin(), itEnd = clusterDetIds.end();
100  it != itEnd;
101  ++it) {
102  const uint8_t type = ((it->first >> 28) & 0xF);
103 
104  const float *corners = item()->getGeom()->getCorners(it->first);
105  if (corners == nullptr)
106  continue;
107 
108  // HGCal
109  if (iData.algo() == 8 || (type >= 8 && type <= 10)) {
110  if (heatmap && hitmap->find(it->first) == hitmap->end())
111  continue;
112 
113  const bool z = (it->first >> 25) & 0x1;
114 
115  // discard everything thats not at the side that we are intersted in
116  if (((z_plus & z_minus) != 1) && (((z_plus | z_minus) == 0) || !(z == z_minus || z == !z_plus)))
117  continue;
118 
119  const float *parameters = item()->getGeom()->getParameters(it->first);
120  const float *shapes = item()->getGeom()->getShapePars(it->first);
121 
122  if (parameters == nullptr || shapes == nullptr)
123  continue;
124 
125  const int total_points = parameters[0];
126  const bool isScintillator = (total_points == 4);
127 
128  uint8_t ll = layer;
129  if (layer > 0) {
130  if (layer > 28) {
131  if (type == 8) {
132  continue;
133  }
134  ll -= 28;
135  } else {
136  if (type != 8) {
137  continue;
138  }
139  }
140 
141  if (ll != ((it->first >> (isScintillator ? 17 : 20)) & 0x1F))
142  continue;
143  }
144 
145  // seed and cluster position
146  if (iData.seed().rawId() == it->first.rawId()) {
147  TEveStraightLineSet *marker = new TEveStraightLineSet;
148  marker->SetLineWidth(1);
149 
150  // center of RecHit
151  const float center[3] = {corners[total_points * 3 + 0],
152  corners[total_points * 3 + 1],
153  corners[total_points * 3 + 2] + shapes[3] * 0.5f};
154 
155  // draw 3D cross
156  const float crossScale = 1.0f + fmin(iData.energy(), 5.0f);
157  marker->AddLine(center[0] - crossScale, center[1], center[2], center[0] + crossScale, center[1], center[2]);
158  marker->AddLine(center[0], center[1] - crossScale, center[2], center[0], center[1] + crossScale, center[2]);
159  marker->AddLine(center[0], center[1], center[2] - crossScale, center[0], center[1], center[2] + crossScale);
160 
161  oItemHolder.AddElement(marker);
162 
163  TEveStraightLineSet *position_marker = new TEveStraightLineSet;
164  position_marker->SetLineWidth(2);
165  position_marker->SetLineColor(kOrange);
166  auto const &pos = iData.position();
167  const float position_crossScale = crossScale * 0.5;
168  position_marker->AddLine(
169  pos.x() - position_crossScale, pos.y(), pos.z(), pos.x() + position_crossScale, pos.y(), pos.z());
170  position_marker->AddLine(
171  pos.x(), pos.y() - position_crossScale, pos.z(), pos.x(), pos.y() + position_crossScale, pos.z());
172 
173  oItemHolder.AddElement(position_marker);
174  }
175 
176  const float energy =
177  fmin((item()->getConfig()->value<bool>("Cluster(0)/RecHit(1)") ? hitmap->at(it->first)->energy()
178  : iData.energy()) /
180  1.0f);
181  const uint8_t colorFactor = gradient_steps * energy;
182 
183  // Scintillator
184  if (isScintillator) {
185  const int total_vertices = 3 * total_points;
186 
187  std::vector<float> pnts(24);
188  for (int i = 0; i < total_points; ++i) {
189  pnts[i * 3 + 0] = corners[i * 3];
190  pnts[i * 3 + 1] = corners[i * 3 + 1];
191  pnts[i * 3 + 2] = corners[i * 3 + 2];
192 
193  pnts[(i * 3 + 0) + total_vertices] = corners[i * 3];
194  pnts[(i * 3 + 1) + total_vertices] = corners[i * 3 + 1];
195  pnts[(i * 3 + 2) + total_vertices] = corners[i * 3 + 2] + shapes[3];
196  }
197  boxset->AddBox(&pnts[0]);
198  if (heatmap) {
199  energy ? boxset->DigitColor(gradient[0][colorFactor], gradient[1][colorFactor], gradient[2][colorFactor])
200  : boxset->DigitColor(64, 64, 64);
201  }
202 
203  h_box = true;
204  }
205  // Silicon
206  else {
207  constexpr int offset = 9;
208 
209  float centerX = (corners[6] + corners[6 + offset]) / 2;
210  float centerY = (corners[7] + corners[7 + offset]) / 2;
211  float radius = fabs(corners[6] - corners[6 + offset]) / 2;
212  hex_boxset->AddHex(TEveVector(centerX, centerY, corners[2]), radius, 90.0, shapes[3]);
213  if (heatmap) {
214  energy ? hex_boxset->DigitColor(gradient[0][colorFactor], gradient[1][colorFactor], gradient[2][colorFactor])
215  : hex_boxset->DigitColor(64, 64, 64);
216  }
217 
218  h_hex = true;
219  }
220  }
221  // Not HGCal
222  else {
223  h_box = true;
224 
225  std::vector<float> pnts(24);
226  fireworks::energyTower3DCorners(corners, (*it).second, pnts);
227  boxset->AddBox(&pnts[0]);
228  }
229  }
230 
231  if (h_hex) {
232  hex_boxset->RefitPlex();
233 
234  hex_boxset->CSCTakeAnyParentAsMaster();
235  if (!heatmap) {
236  hex_boxset->CSCApplyMainColorToMatchingChildren();
237  hex_boxset->CSCApplyMainTransparencyToMatchingChildren();
238  hex_boxset->SetMainColor(item()->modelInfo(iIndex).displayProperties().color());
239  hex_boxset->SetMainTransparency(item()->defaultDisplayProperties().transparency());
240  }
241  oItemHolder.AddElement(hex_boxset);
242  }
243 
244  if (h_box) {
245  boxset->RefitPlex();
246 
247  boxset->CSCTakeAnyParentAsMaster();
248  if (!heatmap) {
249  boxset->CSCApplyMainColorToMatchingChildren();
250  boxset->CSCApplyMainTransparencyToMatchingChildren();
251  boxset->SetMainColor(item()->modelInfo(iIndex).displayProperties().color());
252  boxset->SetMainTransparency(item()->defaultDisplayProperties().transparency());
253  }
254  oItemHolder.AddElement(boxset);
255  }
256 }
257 
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:154
#define REGISTER_PROXYBUILDER_METHODS()
#define REGISTER_FWPROXYBUILDER(_name_, _type_, _purpose_, _view_)
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
AlgoId algo() const
algorithm identifier
Definition: CaloCluster.h:190
void setItem(const FWEventItem *iItem) override
FWProxyBuilderConfiguration * getConfig() const
Definition: FWEventItem.h:150
const_reference_type get(ProductID id, size_t idx) const
Definition: ValueMap.h:144
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:219
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:92
void setItem(const FWEventItem *iItem) override
const float * getCorners(unsigned int id) const
Definition: FWGeometry.cc:461
double f[11][100]
static constexpr uint8_t gradient[3][gradient_steps]
const edm::EventBase * getEvent() const
Definition: FWEventItem.h:131
FWGenericParameter< T > * assertParam(const std::string &name, T def)
double energy() const
cluster energy
Definition: CaloCluster.h:149
void energyTower3DCorners(const float *corners, float scale, std::vector< float > &, bool reflect=false)
edm::Handle< edm::ValueMap< std::pair< float, float > > > TimeValueMapHandle
const FWGeometry * getGeom() const
Definition: FWEventItem.cc:548
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
const FWCaloClusterProxyBuilder & operator=(const FWCaloClusterProxyBuilder &)=delete
const float * getShapePars(unsigned int id) const
Definition: FWGeometry.cc:483
bool isValid() const
Definition: HandleBase.h:70
const float * getParameters(unsigned int id) const
Definition: FWGeometry.cc:472
const FWEventItem * item() const
std::unordered_map< DetId, const HGCRecHit *> * hitmap