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() == reco::CaloCluster::hgcal_em || iData.algo() == reco::CaloCluster::hgcal_had ||
110  (type >= 8 && type <= 10)) {
111  if (heatmap && hitmap->find(it->first) == hitmap->end())
112  continue;
113 
114  const bool z = (it->first >> 25) & 0x1;
115 
116  // discard everything thats not at the side that we are intersted in
117  if (((z_plus & z_minus) != 1) && (((z_plus | z_minus) == 0) || !(z == z_minus || z == !z_plus)))
118  continue;
119 
120  const float *parameters = item()->getGeom()->getParameters(it->first);
121  const float *shapes = item()->getGeom()->getShapePars(it->first);
122 
123  if (parameters == nullptr || shapes == nullptr)
124  continue;
125 
126  const int total_points = parameters[0];
127  const bool isScintillator = (total_points == 4);
128 
129  uint8_t ll = layer;
130  if (layer > 0) {
131  if (layer > 28) {
132  if (type == 8) {
133  continue;
134  }
135  ll -= 28;
136  } else {
137  if (type != 8) {
138  continue;
139  }
140  }
141 
142  if (ll != ((it->first >> (isScintillator ? 17 : 20)) & 0x1F))
143  continue;
144  }
145 
146  // seed and cluster position
147  if (iData.seed().rawId() == it->first.rawId()) {
148  TEveStraightLineSet *marker = new TEveStraightLineSet;
149  marker->SetLineWidth(1);
150 
151  // center of RecHit
152  const float center[3] = {corners[total_points * 3 + 0],
153  corners[total_points * 3 + 1],
154  corners[total_points * 3 + 2] + shapes[3] * 0.5f};
155 
156  // draw 3D cross
157  const float crossScale = 1.0f + fmin(iData.energy(), 5.0f);
158  marker->AddLine(center[0] - crossScale, center[1], center[2], center[0] + crossScale, center[1], center[2]);
159  marker->AddLine(center[0], center[1] - crossScale, center[2], center[0], center[1] + crossScale, center[2]);
160  marker->AddLine(center[0], center[1], center[2] - crossScale, center[0], center[1], center[2] + crossScale);
161 
162  oItemHolder.AddElement(marker);
163 
164  TEveStraightLineSet *position_marker = new TEveStraightLineSet;
165  position_marker->SetLineWidth(2);
166  position_marker->SetLineColor(kOrange);
167  auto const &pos = iData.position();
168  const float position_crossScale = crossScale * 0.5;
169  position_marker->AddLine(
170  pos.x() - position_crossScale, pos.y(), pos.z(), pos.x() + position_crossScale, pos.y(), pos.z());
171  position_marker->AddLine(
172  pos.x(), pos.y() - position_crossScale, pos.z(), pos.x(), pos.y() + position_crossScale, pos.z());
173 
174  oItemHolder.AddElement(position_marker);
175  }
176 
177  const float energy =
178  fmin((item()->getConfig()->value<bool>("Cluster(0)/RecHit(1)") ? hitmap->at(it->first)->energy()
179  : iData.energy()) /
181  1.0f);
182  const uint8_t colorFactor = gradient_steps * energy;
183 
184  // Scintillator
185  if (isScintillator) {
186  const int total_vertices = 3 * total_points;
187 
188  std::vector<float> pnts(24);
189  for (int i = 0; i < total_points; ++i) {
190  pnts[i * 3 + 0] = corners[i * 3];
191  pnts[i * 3 + 1] = corners[i * 3 + 1];
192  pnts[i * 3 + 2] = corners[i * 3 + 2];
193 
194  pnts[(i * 3 + 0) + total_vertices] = corners[i * 3];
195  pnts[(i * 3 + 1) + total_vertices] = corners[i * 3 + 1];
196  pnts[(i * 3 + 2) + total_vertices] = corners[i * 3 + 2] + shapes[3];
197  }
198  boxset->AddBox(&pnts[0]);
199  if (heatmap) {
200  energy ? boxset->DigitColor(gradient[0][colorFactor], gradient[1][colorFactor], gradient[2][colorFactor])
201  : boxset->DigitColor(64, 64, 64);
202  }
203 
204  h_box = true;
205  }
206  // Silicon
207  else {
208  constexpr int offset = 9;
209 
210  float centerX = (corners[6] + corners[6 + offset]) / 2;
211  float centerY = (corners[7] + corners[7 + offset]) / 2;
212  float radius = fabs(corners[6] - corners[6 + offset]) / 2;
213  hex_boxset->AddHex(TEveVector(centerX, centerY, corners[2]), radius, shapes[2], shapes[3]);
214  if (heatmap) {
215  energy ? hex_boxset->DigitColor(gradient[0][colorFactor], gradient[1][colorFactor], gradient[2][colorFactor])
216  : hex_boxset->DigitColor(64, 64, 64);
217  }
218 
219  h_hex = true;
220  }
221  }
222  // Not HGCal
223  else {
224  h_box = true;
225 
226  std::vector<float> pnts(24);
227  fireworks::energyTower3DCorners(corners, (*it).second, pnts);
228  boxset->AddBox(&pnts[0]);
229  }
230  }
231 
232  if (h_hex) {
233  hex_boxset->RefitPlex();
234 
235  hex_boxset->CSCTakeAnyParentAsMaster();
236  if (!heatmap) {
237  hex_boxset->CSCApplyMainColorToMatchingChildren();
238  hex_boxset->CSCApplyMainTransparencyToMatchingChildren();
239  hex_boxset->SetMainColor(item()->modelInfo(iIndex).displayProperties().color());
240  hex_boxset->SetMainTransparency(item()->defaultDisplayProperties().transparency());
241  }
242  oItemHolder.AddElement(hex_boxset);
243  }
244 
245  if (h_box) {
246  boxset->RefitPlex();
247 
248  boxset->CSCTakeAnyParentAsMaster();
249  if (!heatmap) {
250  boxset->CSCApplyMainColorToMatchingChildren();
251  boxset->CSCApplyMainTransparencyToMatchingChildren();
252  boxset->SetMainColor(item()->modelInfo(iIndex).displayProperties().color());
253  boxset->SetMainTransparency(item()->defaultDisplayProperties().transparency());
254  }
255  oItemHolder.AddElement(boxset);
256  }
257 }
258 
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:153
#define REGISTER_PROXYBUILDER_METHODS()
#define REGISTER_FWPROXYBUILDER(_name_, _type_, _purpose_, _view_)
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:209
AlgoId algo() const
algorithm identifier
Definition: CaloCluster.h:189
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:218
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:97
void setItem(const FWEventItem *iItem) override
const float * getCorners(unsigned int id) const
Definition: FWGeometry.cc:439
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:148
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:461
bool isValid() const
Definition: HandleBase.h:70
const float * getParameters(unsigned int id) const
Definition: FWGeometry.cc:450
const FWEventItem * item() const
std::unordered_map< DetId, const HGCRecHit *> * hitmap