CMS 3D CMS Logo

FWTracksterHitsProxyBuilder.cc
Go to the documentation of this file.
12 
13 #include "TEveBoxSet.h"
14 #include "TEveStraightLineSet.h"
15 
17 public:
19  ~FWTracksterHitsProxyBuilder(void) override {}
20 
22 
23  FWTracksterHitsProxyBuilder(const FWTracksterHitsProxyBuilder &) = delete; // stop default
24  const FWTracksterHitsProxyBuilder &operator=(const FWTracksterHitsProxyBuilder &) = delete; // stop default
25 
26 private:
30  long layer_;
32  bool heatmap_;
33  bool z_plus_;
34  bool z_minus_;
39 
40  void setItem(const FWEventItem *iItem) override;
41 
42  void build(const FWEventItem *iItem, TEveElementList *product, const FWViewContext *vc) override;
43  void build(const ticl::Trackster &iData,
44  unsigned int iIndex,
45  TEveElement &oItemHolder,
46  const FWViewContext *) override;
47 };
48 
51  if (iItem) {
52  iItem->getConfig()->assertParam("Cluster(0)/RecHit(1)", false);
53  iItem->getConfig()->assertParam("EnableSeedLines", false);
54  iItem->getConfig()->assertParam("EnablePositionLines", false);
55  iItem->getConfig()->assertParam("EnableEdges", false);
56  iItem->getConfig()->assertParam("EnableTimeFilter", false);
57  iItem->getConfig()->assertParam("TimeLowerBound(ns)", 0.01, 0.0, 75.0);
58  iItem->getConfig()->assertParam("TimeUpperBound(ns)", 0.01, 0.0, 75.0);
59  }
60 }
61 
62 void FWTracksterHitsProxyBuilder::build(const FWEventItem *iItem, TEveElementList *product, const FWViewContext *vc) {
63  iItem->getEvent()->getByLabel(edm::InputTag("hgcalLayerClusters", "timeLayerCluster"), TimeValueMapHandle_);
64  iItem->getEvent()->getByLabel(edm::InputTag("hgcalLayerClusters"), layerClustersHandle_);
66  timeLowerBound_ = item()->getConfig()->value<double>("TimeLowerBound(ns)");
67  timeUpperBound_ = item()->getConfig()->value<double>("TimeUpperBound(ns)");
69  edm::LogWarning("InvalidParameters")
70  << "lower time bound is larger than upper time bound. Maybe opposite is desired?";
71  }
72  } else {
73  edm::LogWarning("DataNotFound|InvalidData") << "couldn't locate 'timeLayerCluster' ValueMap in root file.";
74  }
75 
77  edm::LogWarning("DataNotFound|InvalidData") << "couldn't locate 'timeLayerCluster' ValueMap in root file.";
78  }
79 
80  layer_ = item()->getConfig()->value<long>("Layer");
81  saturation_energy_ = item()->getConfig()->value<double>("EnergyCutOff");
82  heatmap_ = item()->getConfig()->value<bool>("Heatmap");
83  z_plus_ = item()->getConfig()->value<bool>("Z+");
84  z_minus_ = item()->getConfig()->value<bool>("Z-");
85  enableTimeFilter_ = item()->getConfig()->value<bool>("EnableTimeFilter");
86  enableSeedLines_ = item()->getConfig()->value<bool>("EnableSeedLines");
87  enablePositionLines_ = item()->getConfig()->value<bool>("EnablePositionLines");
88  enableEdges_ = item()->getConfig()->value<bool>("EnableEdges");
89 
90  FWHeatmapProxyBuilderTemplate::build(iItem, product, vc);
91 }
92 
94  unsigned int iIndex,
95  TEveElement &oItemHolder,
96  const FWViewContext *) {
98  const float time = TimeValueMapHandle_->get(iIndex).first;
99  if (time < timeLowerBound_ || time > timeUpperBound_)
100  return;
101  }
102 
103  const ticl::Trackster &trackster = iData;
104  const size_t N = trackster.vertices().size();
105  const std::vector<reco::CaloCluster> &layerClusters = *layerClustersHandle_;
106 
107  TEveBoxSet *hex_boxset = new TEveBoxSet();
108  if (!heatmap_)
109  hex_boxset->UseSingleColor();
110  hex_boxset->SetPickable(true);
111  hex_boxset->Reset(TEveBoxSet::kBT_Hex, true, 64);
112 
113  TEveBoxSet *boxset = new TEveBoxSet();
114  if (!heatmap_)
115  boxset->UseSingleColor();
116  boxset->SetPickable(true);
117  boxset->Reset(TEveBoxSet::kBT_FreeBox, true, 64);
118 
119  TEveStraightLineSet *seed_marker = nullptr;
120  if (enableSeedLines_) {
121  seed_marker = new TEveStraightLineSet("seeds");
122  seed_marker->SetLineWidth(2);
123  seed_marker->SetLineColor(kOrange + 10);
124  seed_marker->GetLinePlex().Reset(sizeof(TEveStraightLineSet::Line_t), 2 * N);
125  }
126 
127  TEveStraightLineSet *position_marker = nullptr;
128  if (enablePositionLines_) {
129  position_marker = new TEveStraightLineSet("positions");
130  position_marker->SetLineWidth(2);
131  position_marker->SetLineColor(kOrange);
132  position_marker->GetLinePlex().Reset(sizeof(TEveStraightLineSet::Line_t), 2 * N);
133  }
134 
135  for (size_t i = 0; i < N; ++i) {
136  const reco::CaloCluster layerCluster = layerClusters[trackster.vertices(i)];
137  std::vector<std::pair<DetId, float>> clusterDetIds = layerCluster.hitsAndFractions();
138 
139  for (std::vector<std::pair<DetId, float>>::iterator it = clusterDetIds.begin(), itEnd = clusterDetIds.end();
140  it != itEnd;
141  ++it) {
142  const float *corners = item()->getGeom()->getCorners(it->first);
143  if (corners == nullptr)
144  continue;
145 
146  if (heatmap_ && hitmap->find(it->first) == hitmap->end())
147  continue;
148 
149  const float *parameters = item()->getGeom()->getParameters(it->first);
150  const float *shapes = item()->getGeom()->getShapePars(it->first);
151 
152  if (parameters == nullptr || shapes == nullptr)
153  continue;
154 
155  const int total_points = parameters[0];
156  const int layer = parameters[1];
157  const int zside = parameters[2];
158  const bool isSilicon = parameters[3];
159 
160  // discard everything that's not at the side that we are intersted in
161  auto const z_selection_is_on = z_plus_ ^ z_minus_;
162  auto const z_plus_selection_ok = z_plus_ && (zside == 1);
163  auto const z_minus_selection_ok = z_minus_ && (zside == -1);
164  if (!z_minus_ && !z_plus_)
165  break;
166  if (z_selection_is_on && !(z_plus_selection_ok || z_minus_selection_ok))
167  break;
168 
169  if (layer_ > 0 && layer != layer_)
170  break;
171 
172  // seed and cluster position
173  if (layerCluster.seed().rawId() == it->first.rawId()) {
174  const float crossScale = 0.2f + fmin(layerCluster.energy(), 5.0f);
175  if (enableSeedLines_) {
176  // center of RecHit
177  const float center[3] = {corners[total_points * 3 + 0],
178  corners[total_points * 3 + 1],
179  corners[total_points * 3 + 2] + shapes[3] * 0.5f};
180 
181  // draw 3D cross
182  seed_marker->AddLine(
183  center[0] - crossScale, center[1], center[2], center[0] + crossScale, center[1], center[2]);
184  seed_marker->AddLine(
185  center[0], center[1] - crossScale, center[2], center[0], center[1] + crossScale, center[2]);
186  }
187 
188  if (enablePositionLines_) {
189  auto const &pos = layerCluster.position();
190  const float position_crossScale = crossScale * 0.5;
191  position_marker->AddLine(
192  pos.x() - position_crossScale, pos.y(), pos.z(), pos.x() + position_crossScale, pos.y(), pos.z());
193  position_marker->AddLine(
194  pos.x(), pos.y() - position_crossScale, pos.z(), pos.x(), pos.y() + position_crossScale, pos.z());
195  }
196  }
197 
198  const float energy =
199  fmin((item()->getConfig()->value<bool>("Cluster(0)/RecHit(1)") ? hitmap->at(it->first)->energy()
200  : layerCluster.energy()) /
202  1.0f);
203  const uint8_t colorFactor = gradient_steps * energy;
204  auto transparency = item()->modelInfo(iIndex).displayProperties().transparency();
205  UChar_t alpha = (255 * (100 - transparency)) / 100;
206 
207  // Scintillator
208  if (!isSilicon) {
209  const int total_vertices = 3 * total_points;
210 
211  std::vector<float> pnts(24);
212  for (int i = 0; i < total_points; ++i) {
213  pnts[i * 3 + 0] = corners[i * 3];
214  pnts[i * 3 + 1] = corners[i * 3 + 1];
215  pnts[i * 3 + 2] = corners[i * 3 + 2];
216 
217  pnts[(i * 3 + 0) + total_vertices] = corners[i * 3];
218  pnts[(i * 3 + 1) + total_vertices] = corners[i * 3 + 1];
219  pnts[(i * 3 + 2) + total_vertices] = corners[i * 3 + 2] + shapes[3];
220  }
221  boxset->AddBox(&pnts[0]);
222  if (heatmap_) {
223  energy
224  ? boxset->DigitColor(gradient[0][colorFactor], gradient[1][colorFactor], gradient[2][colorFactor], alpha)
225  : boxset->DigitColor(64, 64, 64, alpha);
226  }
227  }
228  // Silicon
229  else {
230  constexpr int offset = 9;
231 
232  float centerX = (corners[6] + corners[6 + offset]) / 2;
233  float centerY = (corners[7] + corners[7 + offset]) / 2;
234  float radius = fabs(corners[6] - corners[6 + offset]) / 2;
235  hex_boxset->AddHex(TEveVector(centerX, centerY, corners[2]), radius, shapes[2], shapes[3]);
236  if (heatmap_) {
237  energy ? hex_boxset->DigitColor(
238  gradient[0][colorFactor], gradient[1][colorFactor], gradient[2][colorFactor], alpha)
239  : hex_boxset->DigitColor(64, 64, 64, alpha);
240  } else {
241  hex_boxset->CSCApplyMainColorToMatchingChildren();
242  hex_boxset->CSCApplyMainTransparencyToMatchingChildren();
243  hex_boxset->SetMainColor(item()->modelInfo(iIndex).displayProperties().color());
244  hex_boxset->SetMainTransparency(item()->defaultDisplayProperties().transparency());
245  }
246  }
247  } // End of loop over rechits of a single layercluster
248  } // End loop over the layerclusters of the trackster
249 
250  hex_boxset->RefitPlex();
251  boxset->RefitPlex();
252  setupAddElement(hex_boxset, &oItemHolder);
253  setupAddElement(boxset, &oItemHolder);
254 
255  if (enableSeedLines_)
256  oItemHolder.AddElement(seed_marker);
257 
259  oItemHolder.AddElement(position_marker);
260 
261  if (enableEdges_) {
262  auto &edges = trackster.edges();
263 
264  TEveStraightLineSet *adjacent_marker = new TEveStraightLineSet("adj_edges");
265  adjacent_marker->SetLineWidth(2);
266  adjacent_marker->SetLineColor(kYellow);
267  adjacent_marker->GetLinePlex().Reset(sizeof(TEveStraightLineSet::Line_t), edges.size());
268 
269  TEveStraightLineSet *non_adjacent_marker = new TEveStraightLineSet("non_adj_edges");
270  non_adjacent_marker->SetLineWidth(2);
271  non_adjacent_marker->SetLineColor(kRed);
272  non_adjacent_marker->GetLinePlex().Reset(sizeof(TEveStraightLineSet::Line_t), edges.size());
273 
274  for (auto edge : edges) {
275  auto doublet = std::make_pair(layerClusters[edge[0]], layerClusters[edge[1]]);
276 
277  int layerIn = item()->getGeom()->getParameters(doublet.first.seed())[1];
278  int layerOut = item()->getGeom()->getParameters(doublet.second.seed())[1];
279 
280  const bool isAdjacent = std::abs(layerOut - layerIn) == 1;
281 
282  // draw 3D cross
283  if (layer_ == 0 || std::abs(layerIn - layer_) == 0 || std::abs(layerOut - layer_) == 0) {
284  if (isAdjacent)
285  adjacent_marker->AddLine(doublet.first.x(),
286  doublet.first.y(),
287  doublet.first.z(),
288  doublet.second.x(),
289  doublet.second.y(),
290  doublet.second.z());
291  else
292  non_adjacent_marker->AddLine(doublet.first.x(),
293  doublet.first.y(),
294  doublet.first.z(),
295  doublet.second.x(),
296  doublet.second.y(),
297  doublet.second.z());
298  }
299  } // End of loop over all edges of the trackster
300  oItemHolder.AddElement(adjacent_marker);
301  oItemHolder.AddElement(non_adjacent_marker);
302  }
303 }
304 
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:154
float alpha
Definition: AMPTWrapper.h:105
#define REGISTER_PROXYBUILDER_METHODS()
#define REGISTER_FWPROXYBUILDER(_name_, _type_, _purpose_, _view_)
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
void setItem(const FWEventItem *iItem) override
void setupAddElement(TEveElement *el, TEveElement *parent, bool set_color=true) const
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
Char_t transparency() const
edm::Handle< edm::ValueMap< std::pair< float, float > > > TimeValueMapHandle_
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:97
int zside(DetId const &)
std::vector< std::array< unsigned int, 2 > > & edges()
Definition: Trackster.h:59
edm::Handle< std::vector< reco::CaloCluster > > layerClustersHandle_
const FWDisplayProperties & displayProperties() const
Definition: FWEventItem.h:64
const float * getCorners(unsigned int id) const
Definition: FWGeometry.cc:439
void setItem(const FWEventItem *iItem) override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
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
std::vector< unsigned int > & vertices()
Definition: Trackster.h:57
const FWGeometry * getGeom() const
Definition: FWEventItem.cc:548
#define N
Definition: blowfish.cc:9
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
const float * getShapePars(unsigned int id) const
Definition: FWGeometry.cc:461
bool isValid() const
Definition: HandleBase.h:70
const FWTracksterHitsProxyBuilder & operator=(const FWTracksterHitsProxyBuilder &)=delete
const float * getParameters(unsigned int id) const
Definition: FWGeometry.cc:450
Log< level::Warning, false > LogWarning
ModelInfo modelInfo(int iIndex) const
Definition: FWEventItem.cc:446
const FWEventItem * item() const
std::unordered_map< DetId, const HGCRecHit *> * hitmap