CMS 3D CMS Logo

FWCandidateHGCalLegoProxyBuilder.cc
Go to the documentation of this file.
1 #ifndef _FWCANDIDATEHGCALLEGOPROXYBUILDER_H_
2 #define _FWCANDIDATEHGCALLEGOPROXYBUILDER_H_
3 
4 // -*- C++ -*-
5 //
6 // Package: Candidates
7 // Class : FWCandidateHGCalLegoProxyBuilder
8 //
9 //
10 
11 // User include files
16 
18 
19 #include <math.h>
20 
21 //-----------------------------------------------------------------------------
22 // FWCandidateHGCalLegoProxyBuilder
23 //-----------------------------------------------------------------------------
24 class FWCandidateHGCalLegoProxyBuilder : public FWSimpleProxyBuilderTemplate<reco::HGCalMultiCluster>
25 {
26  public:
27  // ---------------- Constructor(s)/Destructor ----------------------
30 
31  // --------------------- Member Functions --------------------------
32  virtual bool havePerViewProduct( FWViewType::EType ) const override { return true; }
33  virtual void scaleProduct( TEveElementList*, FWViewType::EType, const FWViewContext* ) override;
34  virtual void localModelChanges( const FWModelId&, TEveElement*, FWViewType::EType,
35  const FWViewContext* ) override;
36 
38 
39  private:
40  // ----------------------- Data Members ----------------------------
43 
44  // --------------------- Member Functions --------------------------
46  void build( const reco::HGCalMultiCluster&, unsigned int, TEveElement&, const FWViewContext* ) override;
47 };
48 #endif
49 //=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_
50 
51 
52 //______________________________________________________________________________
53 void
55  const FWViewContext *vc )
56 {
57  for( TEveElement::List_i i = parent->BeginChildren(); i != parent->EndChildren(); ++i )
58  {
59  if( (*i)->HasChildren() )
60  {
61  TEveElement *el = (*i)->FirstChild(); // There is only one child
62  FWLegoCandidate *candidate = dynamic_cast<FWLegoCandidate*> (el);
63  candidate->updateScale( vc, context() );
64  }
65  }
66 }
67 
68 //______________________________________________________________________________
69 void
72 {
73  // Line set marker is nto the same color as line, have to fix it here
74  if( (parent)->HasChildren() )
75  {
76  TEveElement *el = (parent)->FirstChild(); // There is only one child
77  FWLegoCandidate *candidate = dynamic_cast<FWLegoCandidate*> (el);
78  const FWDisplayProperties &dp = item()->modelInfo( iId.index() ).displayProperties();
79  candidate->SetMarkerColor( dp.color() );
80  candidate->ElementChanged();
81  }
82 }
83 
84 //______________________________________________________________________________
85 void
87  TEveElement &oItemHolder, const FWViewContext *vc )
88 {
89  const auto & clusters = iData.clusters();
90 
91  for (const auto & c : clusters)
92  {
93  auto pt = c->energy()/std::cosh(c->eta());
94  FWLegoCandidate *candidate = new FWLegoCandidate( vc, context(), c->energy(),
95  pt, pt,
96  c->eta(), c->phi() );
97 
98  candidate->SetMarkerColor( item()->defaultDisplayProperties().color() );
99  context().voteMaxEtAndEnergy( pt, c->energy() );
100  setupAddElement( candidate, &oItemHolder );
101  }
102 }
103 
104 //______________________________________________________________________________
type
Definition: HCALResponse.h:21
const fireworks::Context & context() const
virtual bool havePerViewProduct(FWViewType::EType) const override
#define REGISTER_PROXYBUILDER_METHODS()
#define REGISTER_FWPROXYBUILDER(_name_, _type_, _purpose_, _view_)
void voteMaxEtAndEnergy(float Et, float energy) const
Definition: Context.cc:183
void setupAddElement(TEveElement *el, TEveElement *parent, bool set_color=true) const
const FWEventItem * item() const
Color_t color() const
int index() const
Definition: FWModelId.h:49
const edm::PtrVector< reco::BasicCluster > & clusters() const
virtual void localModelChanges(const FWModelId &, TEveElement *, FWViewType::EType, const FWViewContext *) override
auto dp
Definition: deltaR.h:22
void updateScale(const FWViewContext *vc, const fireworks::Context &)
const FWCandidateHGCalLegoProxyBuilder & operator=(const FWCandidateHGCalLegoProxyBuilder &)
virtual void scaleProduct(TEveElementList *, FWViewType::EType, const FWViewContext *) override
ModelInfo modelInfo(int iIndex) const
Definition: FWEventItem.cc:537