CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 {
29  public:
30  // ---------------- Constructor(s)/Destructor ----------------------
33 
34  // --------------------- Member Functions --------------------------
35  virtual bool havePerViewProduct( FWViewType::EType ) const { return true; }
36  virtual void scaleProduct( TEveElementList*, FWViewType::EType, const FWViewContext* );
37  virtual void localModelChanges( const FWModelId&, TEveElement*, FWViewType::EType,
38  const FWViewContext* );
39 
41 
42  private:
43  // ----------------------- Data Members ----------------------------
46 
47  // --------------------- Member Functions --------------------------
48  void build( const reco::Candidate&, unsigned int, TEveElement&, const FWViewContext* );
49 };
50 #endif
51 //=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_=_
52 
53 
54 //______________________________________________________________________________
55 void
57  const FWViewContext *vc )
58 {
59  for( TEveElement::List_i i = parent->BeginChildren(); i != parent->EndChildren(); ++i )
60  {
61  if( (*i)->HasChildren() )
62  {
63  TEveElement *el = (*i)->FirstChild(); // There is only one child
64  FWLegoCandidate *candidate = dynamic_cast<FWLegoCandidate*> (el);
65  candidate->updateScale( vc, context() );
66  }
67  }
68 }
69 
70 //______________________________________________________________________________
71 void
74 {
75  // Line set marker is nto the same color as line, have to fix it here
76  if( (parent)->HasChildren() )
77  {
78  TEveElement *el = (parent)->FirstChild(); // There is only one child
79  FWLegoCandidate *candidate = dynamic_cast<FWLegoCandidate*> (el);
80  const FWDisplayProperties &dp = item()->modelInfo( iId.index() ).displayProperties();
81  candidate->SetMarkerColor( dp.color() );
82  candidate->ElementChanged();
83  }
84 }
85 
86 //______________________________________________________________________________
87 void
88 FWCandidateLegoProxyBuilder::build( const reco::Candidate &iData, unsigned int iIndex,
89  TEveElement &oItemHolder, const FWViewContext *vc )
90 {
91  FWLegoCandidate *candidate = new FWLegoCandidate( vc, context(), iData.energy(),
92  iData.et(), iData.pt(),
93  iData.eta(), iData.phi() );
94 
95  candidate->SetMarkerColor( item()->defaultDisplayProperties().color() );
96  context().voteMaxEtAndEnergy( iData.et(), iData.energy() );
97 
98  setupAddElement( candidate, &oItemHolder );
99 }
100 
101 //______________________________________________________________________________
type
Definition: HCALResponse.h:22
const fireworks::Context & context() const
int i
Definition: DBlmapReader.cc:9
virtual double energy() const =0
energy
#define REGISTER_FWPROXYBUILDER(_name_, _type_, _purpose_, _view_)
list parent
Definition: dbtoconf.py:74
virtual double et() const =0
transverse energy
virtual void scaleProduct(TEveElementList *, FWViewType::EType, const FWViewContext *)
virtual double pt() const =0
transverse momentum
void voteMaxEtAndEnergy(float Et, float energy) const
Definition: Context.cc:185
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:50
virtual bool havePerViewProduct(FWViewType::EType) const
virtual void localModelChanges(const FWModelId &, TEveElement *, FWViewType::EType, const FWViewContext *)
void updateScale(const FWViewContext *vc, const fireworks::Context &)
const FWCandidateLegoProxyBuilder & operator=(const FWCandidateLegoProxyBuilder &)
ModelInfo modelInfo(int iIndex) const
Definition: FWEventItem.cc:535
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity