CMS 3D CMS Logo

FWCandidateLegoProxyBuilder.cc
Go to the documentation of this file.
1 #ifndef _FWCANDIDATELEGOPROXYBUILDER_H_
2 #define _FWCANDIDATELEGOPROXYBUILDER_H_
3 
4 // -*- C++ -*-
5 //
6 // Package: Candidates
7 // Class : FWCandidateLegoProxyBuilder
8 //
9 // Implementation:
10 // <Notes on implementation>
11 //
12 // Original Author: Simon Harris
13 // Created: 24/06/2011
14 //
15 
16 // User include files
21 
23 
24 //-----------------------------------------------------------------------------
25 // FWCandidateLegoProxyBuilder
26 //-----------------------------------------------------------------------------
28 public:
29  // ---------------- Constructor(s)/Destructor ----------------------
32 
33  // --------------------- Member Functions --------------------------
34  bool havePerViewProduct(FWViewType::EType) const override { return true; }
35  void scaleProduct(TEveElementList *, FWViewType::EType, const FWViewContext *) override;
36  void localModelChanges(const FWModelId &, TEveElement *, FWViewType::EType, const FWViewContext *) override;
37 
39 
40  // ----------------------- Data Members ----------------------------
43 
44 private:
45  // --------------------- Member Functions --------------------------
47  void build(const reco::Candidate &, unsigned int, TEveElement &, const FWViewContext *) override;
48 };
49 #endif
50 //=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_
51 
52 //______________________________________________________________________________
55  const FWViewContext *vc) {
56  for (TEveElement::List_i i = parent->BeginChildren(); i != parent->EndChildren(); ++i) {
57  if ((*i)->HasChildren()) {
58  TEveElement *el = (*i)->FirstChild(); // There is only one child
59  FWLegoCandidate *candidate = dynamic_cast<FWLegoCandidate *>(el);
60  candidate->updateScale(vc, context());
61  }
62  }
63 }
64 
65 //______________________________________________________________________________
67  TEveElement *parent,
69  const FWViewContext *vc) {
70  // Line set marker is nto the same color as line, have to fix it here
71  if ((parent)->HasChildren()) {
72  TEveElement *el = (parent)->FirstChild(); // There is only one child
73  FWLegoCandidate *candidate = dynamic_cast<FWLegoCandidate *>(el);
74  candidate->SetMarkerColor(item()->modelInfo(iId.index()).displayProperties().color());
75  candidate->ElementChanged();
76  }
77 }
78 
79 //______________________________________________________________________________
81  unsigned int iIndex,
82  TEveElement &oItemHolder,
83  const FWViewContext *vc) {
84  FWLegoCandidate *candidate =
85  new FWLegoCandidate(vc, context(), iData.energy(), iData.et(), iData.pt(), iData.eta(), iData.phi());
86 
87  candidate->SetMarkerColor(item()->defaultDisplayProperties().color());
88  context().voteMaxEtAndEnergy(iData.et(), iData.energy());
89 
90  setupAddElement(candidate, &oItemHolder);
91 }
92 
93 //______________________________________________________________________________
96  "Candidates",
virtual double energy() const =0
energy
#define REGISTER_PROXYBUILDER_METHODS()
#define REGISTER_FWPROXYBUILDER(_name_, _type_, _purpose_, _view_)
void setupAddElement(TEveElement *el, TEveElement *parent, bool set_color=true) const
virtual double et() const =0
transverse energy
virtual double pt() const =0
transverse momentum
void scaleProduct(TEveElementList *, FWViewType::EType, const FWViewContext *) override
int index() const
Definition: FWModelId.h:41
const fireworks::Context & context() const
bool havePerViewProduct(FWViewType::EType) const override
void voteMaxEtAndEnergy(float Et, float energy) const
Definition: Context.cc:162
void updateScale(const FWViewContext *vc, const fireworks::Context &)
const FWCandidateLegoProxyBuilder & operator=(const FWCandidateLegoProxyBuilder &)=delete
void localModelChanges(const FWModelId &, TEveElement *, FWViewType::EType, const FWViewContext *) override
const FWEventItem * item() const
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity