CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
FWHeatmapProxyBuilderTemplate.h
Go to the documentation of this file.
1 #ifndef Fireworks_Core_FWHeatmapProxyBuilderTemplate_h
2 #define Fireworks_Core_FWHeatmapProxyBuilderTemplate_h
3 // -*- C++ -*-
4 //
5 // Package: Core
6 // Class : FWHeatmapProxyBuilderTemplate
7 //
16 //
17 // Original Author: Alex Mourtziapis
18 // Created: Wed Jan 23 14:50:00 EST 2019
19 //
20 
21 // system include files
22 #include <cmath>
23 
24 // user include files
30 
31 // forward declarations
32 
33 template <typename T>
35 public:
37 
38  //virtual ~FWHeatmapProxyBuilderTemplate();
39 
40  // ---------- const member functions ---------------------
41 
42  // ---------- static member functions --------------------
43 
44  // ---------- member functions ---------------------------
45 
46 protected:
47  std::unordered_map<DetId, const HGCRecHit*>* hitmap = new std::unordered_map<DetId, const HGCRecHit*>;
48 
49  static constexpr uint8_t gradient_steps = 9;
50  static constexpr uint8_t gradient[3][gradient_steps] = {{static_cast<uint8_t>(0.2082 * 255),
51  static_cast<uint8_t>(0.0592 * 255),
52  static_cast<uint8_t>(0.0780 * 255),
53  static_cast<uint8_t>(0.0232 * 255),
54  static_cast<uint8_t>(0.1802 * 255),
55  static_cast<uint8_t>(0.5301 * 255),
56  static_cast<uint8_t>(0.8186 * 255),
57  static_cast<uint8_t>(0.9956 * 255),
58  static_cast<uint8_t>(0.9764 * 255)},
59 
60  {static_cast<uint8_t>(0.1664 * 255),
61  static_cast<uint8_t>(0.3599 * 255),
62  static_cast<uint8_t>(0.5041 * 255),
63  static_cast<uint8_t>(0.6419 * 255),
64  static_cast<uint8_t>(0.7178 * 255),
65  static_cast<uint8_t>(0.7492 * 255),
66  static_cast<uint8_t>(0.7328 * 255),
67  static_cast<uint8_t>(0.7862 * 255),
68  static_cast<uint8_t>(0.9832 * 255)},
69 
70  {static_cast<uint8_t>(0.5293 * 255),
71  static_cast<uint8_t>(0.8684 * 255),
72  static_cast<uint8_t>(0.8385 * 255),
73  static_cast<uint8_t>(0.7914 * 255),
74  static_cast<uint8_t>(0.6425 * 255),
75  static_cast<uint8_t>(0.4662 * 255),
76  static_cast<uint8_t>(0.3499 * 255),
77  static_cast<uint8_t>(0.1968 * 255),
78  static_cast<uint8_t>(0.0539 * 255)}};
79 
80  const T& modelData(int index) { return *reinterpret_cast<const T*>(m_helper.offsetObject(item()->modelData(index))); }
81 
82  void setItem(const FWEventItem* iItem) override {
84  if (iItem) {
85  iItem->getConfig()->keepEntries(true);
86  iItem->getConfig()->assertParam("Layer", 0L, 0L, 52L);
87  iItem->getConfig()->assertParam("EnergyCutOff", 0.5, 0.2, 5.0);
88  iItem->getConfig()->assertParam("Heatmap", true);
89  iItem->getConfig()->assertParam("Z+", true);
90  iItem->getConfig()->assertParam("Z-", true);
91  }
92  }
93 
94  void build(const FWEventItem* iItem, TEveElementList* product, const FWViewContext* vc) override {
95  if (item()->getConfig()->template value<bool>("Heatmap")) {
96  hitmap->clear();
97 
98  edm::Handle<HGCRecHitCollection> recHitHandleEE;
99  edm::Handle<HGCRecHitCollection> recHitHandleFH;
100  edm::Handle<HGCRecHitCollection> recHitHandleBH;
101 
102  const edm::EventBase* event = iItem->getEvent();
103 
104  event->getByLabel(edm::InputTag("HGCalRecHit", "HGCEERecHits"), recHitHandleEE);
105  event->getByLabel(edm::InputTag("HGCalRecHit", "HGCHEFRecHits"), recHitHandleFH);
106  event->getByLabel(edm::InputTag("HGCalRecHit", "HGCHEBRecHits"), recHitHandleBH);
107 
108  if (recHitHandleEE.isValid()) {
109  const auto& rechitsEE = *recHitHandleEE;
110 
111  for (unsigned int i = 0; i < rechitsEE.size(); ++i) {
112  (*hitmap)[rechitsEE[i].detid().rawId()] = &rechitsEE[i];
113  }
114  }
115 
116  if (recHitHandleFH.isValid()) {
117  const auto& rechitsFH = *recHitHandleFH;
118 
119  for (unsigned int i = 0; i < rechitsFH.size(); ++i) {
120  (*hitmap)[rechitsFH[i].detid().rawId()] = &rechitsFH[i];
121  }
122  }
123 
124  if (recHitHandleBH.isValid()) {
125  const auto& rechitsBH = *recHitHandleBH;
126 
127  for (unsigned int i = 0; i < rechitsBH.size(); ++i) {
128  (*hitmap)[rechitsBH[i].detid().rawId()] = &rechitsBH[i];
129  }
130  }
131  }
132 
133  FWSimpleProxyBuilder::build(iItem, product, vc);
134  }
135 
137  void build(const void* iData, unsigned int iIndex, TEveElement& oItemHolder, const FWViewContext* context) override {
138  if (nullptr != iData) {
139  build(*reinterpret_cast<const T*>(iData), iIndex, oItemHolder, context);
140  }
141  }
142 
144  void buildViewType(const void* iData,
145  unsigned int iIndex,
146  TEveElement& oItemHolder,
147  FWViewType::EType viewType,
148  const FWViewContext* context) override {
149  if (nullptr != iData) {
150  buildViewType(*reinterpret_cast<const T*>(iData), iIndex, oItemHolder, viewType, context);
151  }
152  }
156  virtual void build(const T& iData, unsigned int iIndex, TEveElement& oItemHolder, const FWViewContext*) {
157  throw std::runtime_error(
158  "virtual build(const T&, unsigned int, TEveElement&, const FWViewContext*) not implemented by inherited "
159  "class.");
160  }
161 
162  virtual void buildViewType(
163  const T& iData, unsigned int iIndex, TEveElement& oItemHolder, FWViewType::EType viewType, const FWViewContext*) {
164  throw std::runtime_error(
165  "virtual buildViewType(const T&, unsigned int, TEveElement&, FWViewType::EType, const FWViewContext*) not "
166  "implemented by inherited class");
167  };
168 
169 public:
170  FWHeatmapProxyBuilderTemplate(const FWHeatmapProxyBuilderTemplate&) = delete; // stop default
171 
172  const FWHeatmapProxyBuilderTemplate& operator=(const FWHeatmapProxyBuilderTemplate&) = delete; // stop default
173 private:
174  // ---------- member data --------------------------------
175 };
176 
177 #endif
const fireworks::Context & context() const
FWProxyBuilderConfiguration * getConfig() const
Definition: FWEventItem.h:150
void setItem(const FWEventItem *iItem) override
void build(const void *iData, unsigned int iIndex, TEveElement &oItemHolder, const FWViewContext *context) override
void build(const FWEventItem *iItem, TEveElementList *product, const FWViewContext *vc) override
virtual void buildViewType(const T &iData, unsigned int iIndex, TEveElement &oItemHolder, FWViewType::EType viewType, const FWViewContext *)
const void * offsetObject(const void *iObj) const
const FWEventItem * item() const
const FWHeatmapProxyBuilderTemplate & operator=(const FWHeatmapProxyBuilderTemplate &)=delete
void buildViewType(const void *iData, unsigned int iIndex, TEveElement &oItemHolder, FWViewType::EType viewType, const FWViewContext *context) override
virtual void setItem(const FWEventItem *iItem)
static constexpr uint8_t gradient[3][gradient_steps]
bool isValid() const
Definition: HandleBase.h:70
FWGenericParameter< T > * assertParam(const std::string &name, T def)
const edm::EventBase * getEvent() const
Definition: FWEventItem.h:131
FWSimpleProxyHelper m_helper
void buildViewType(const FWEventItem *iItem, TEveElementList *product, FWViewType::EType viewType, const FWViewContext *) override
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:92
virtual void build(const T &iData, unsigned int iIndex, TEveElement &oItemHolder, const FWViewContext *)
long double T
std::unordered_map< DetId, const HGCRecHit * > * hitmap