CMS 3D CMS Logo

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 
36 public:
38  FWSimpleProxyBuilder(typeid(T)) {
39  }
40 
41  //virtual ~FWHeatmapProxyBuilderTemplate();
42 
43  // ---------- const member functions ---------------------
44 
45  // ---------- static member functions --------------------
46 
47  // ---------- member functions ---------------------------
48 
49 protected:
50  std::map<DetId, const HGCRecHit*> hitmap;
51 
52  static constexpr uint8_t gradient_steps = 9;
53  static constexpr uint8_t gradient[3][gradient_steps] = {
54  {static_cast<uint8_t>(0.2082*255), static_cast<uint8_t>(0.0592*255), static_cast<uint8_t>(0.0780*255),
55  static_cast<uint8_t>(0.0232*255), static_cast<uint8_t>(0.1802*255), static_cast<uint8_t>(0.5301*255),
56  static_cast<uint8_t>(0.8186*255), static_cast<uint8_t>(0.9956*255), static_cast<uint8_t>(0.9764*255)},
57 
58  {static_cast<uint8_t>(0.1664*255), static_cast<uint8_t>(0.3599*255), static_cast<uint8_t>(0.5041*255),
59  static_cast<uint8_t>(0.6419*255), static_cast<uint8_t>(0.7178*255), static_cast<uint8_t>(0.7492*255),
60  static_cast<uint8_t>(0.7328*255), static_cast<uint8_t>(0.7862*255), static_cast<uint8_t>(0.9832*255)},
61 
62  {static_cast<uint8_t>(0.5293*255), static_cast<uint8_t>(0.8684*255), static_cast<uint8_t>(0.8385*255),
63  static_cast<uint8_t>(0.7914*255), static_cast<uint8_t>(0.6425*255), static_cast<uint8_t>(0.4662*255),
64  static_cast<uint8_t>(0.3499*255), static_cast<uint8_t>(0.1968*255), static_cast<uint8_t>(0.0539*255)}
65  };
66 
67  const T& modelData(int index) { return *reinterpret_cast<const T*>(m_helper.offsetObject(item()->modelData(index))); }
68 
69  void setItem(const FWEventItem *iItem) override
70  {
72  if (iItem)
73  {
74  iItem->getConfig()->keepEntries(true);
75  iItem->getConfig()->assertParam("Layer", 0L, 0L, 52L);
76  iItem->getConfig()->assertParam("EnergyCutOff", 0.5, 0.2, 5.0);
77  iItem->getConfig()->assertParam("Heatmap", true);
78  iItem->getConfig()->assertParam("Z+", true);
79  iItem->getConfig()->assertParam("Z-", true);
80  }
81  }
82 
83  void build(const FWEventItem *iItem, TEveElementList *product, const FWViewContext *vc) override
84  {
85  if (item()->getConfig()->template value<bool>("Heatmap"))
86  {
87  hitmap.clear();
88 
89  const edm::EventBase *event = iItem->getEvent();
90 
91  edm::Handle<HGCRecHitCollection> recHitHandleEE;
92  edm::Handle<HGCRecHitCollection> recHitHandleFH;
93  edm::Handle<HGCRecHitCollection> recHitHandleBH;
94 
95  event->getByLabel( edm::InputTag( "HGCalRecHit", "HGCEERecHits" ), recHitHandleEE );
96  event->getByLabel( edm::InputTag( "HGCalRecHit", "HGCHEFRecHits" ), recHitHandleFH );
97  event->getByLabel( edm::InputTag( "HGCalRecHit", "HGCHEBRecHits" ), recHitHandleBH );
98 
99  const auto& rechitsEE = *recHitHandleEE;
100  const auto& rechitsFH = *recHitHandleFH;
101  const auto& rechitsBH = *recHitHandleBH;
102 
103  for (unsigned int i = 0; i < rechitsEE.size(); ++i) {
104  hitmap[rechitsEE[i].detid().rawId()] = &rechitsEE[i];
105  }
106  for (unsigned int i = 0; i < rechitsFH.size(); ++i) {
107  hitmap[rechitsFH[i].detid().rawId()] = &rechitsFH[i];
108  }
109  for (unsigned int i = 0; i < rechitsBH.size(); ++i) {
110  hitmap[rechitsBH[i].detid().rawId()] = &rechitsBH[i];
111  }
112  }
113 
114  FWSimpleProxyBuilder::build(iItem, product, vc);
115  }
116 
118  void build(const void *iData, unsigned int iIndex, TEveElement& oItemHolder, const FWViewContext* context) override
119  {
120  if(nullptr!=iData) {
121  build(*reinterpret_cast<const T*> (iData), iIndex, oItemHolder, context);
122  }
123  }
124 
126  void buildViewType(const void *iData, unsigned int iIndex, TEveElement& oItemHolder, FWViewType::EType viewType, const FWViewContext* context) override
127  {
128  if(nullptr!=iData) {
129  buildViewType(*reinterpret_cast<const T*> (iData), iIndex, oItemHolder, viewType, context);
130  }
131  }
135  virtual void build(const T& iData, unsigned int iIndex,TEveElement& oItemHolder, const FWViewContext*)
136  {
137  throw std::runtime_error("virtual build(const T&, unsigned int, TEveElement&, const FWViewContext*) not implemented by inherited class.");
138  }
139 
140  virtual void buildViewType(const T& iData, unsigned int iIndex,TEveElement& oItemHolder, FWViewType::EType viewType, const FWViewContext*)
141  {
142  throw std::runtime_error("virtual buildViewType(const T&, unsigned int, TEveElement&, FWViewType::EType, const FWViewContext*) not implemented by inherited class");
143  };
144 private:
145  FWHeatmapProxyBuilderTemplate(const FWHeatmapProxyBuilderTemplate&) = delete; // stop default
146 
147  const FWHeatmapProxyBuilderTemplate& operator=(const FWHeatmapProxyBuilderTemplate&) = delete; // stop default
148 
149 
150 
151  // ---------- member data --------------------------------
152 
153 };
154 
155 
156 #endif
const fireworks::Context & context() const
FWProxyBuilderConfiguration * getConfig() const
Definition: FWEventItem.h:169
std::map< DetId, const HGCRecHit * > hitmap
void setItem(const FWEventItem *iItem) override
static uint8_t gradient[3][gradient_steps]
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)
FWGenericParameter< T > * assertParam(const std::string &name, T def)
const edm::EventBase * getEvent() const
Definition: FWEventItem.h:148
FWSimpleProxyHelper m_helper
void buildViewType(const FWEventItem *iItem, TEveElementList *product, FWViewType::EType viewType, const FWViewContext *) override
virtual void build(const T &iData, unsigned int iIndex, TEveElement &oItemHolder, const FWViewContext *)
long double T
#define constexpr