CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FWLegoEvePFCandidate.cc
Go to the documentation of this file.
1 #include "TEveCaloData.h"
2 
9 
10 
12  m_energy(0.f),
13  m_et(0.f)
14 {
15  m_et = iData.et();
16  m_energy = iData.energy();
17 
18  float base = 0.001; // flour offset 1%
19 
20  // first vertical line , which is et/energy
21  FWViewEnergyScale* caloScale = vc->getEnergyScale();
22  float val = caloScale->getPlotEt() ? m_et : m_energy;
23  AddLine(iData.eta(),iData.phi(), base,
24  iData.eta(),iData.phi(), base + val*caloScale->getScaleFactorLego());
25 
26 
27  AddMarker(0, 1.f);
28  SetMarkerStyle(3);
29  SetMarkerSize(0.01);
30  SetDepthTest(false);
31 
32  // circle pt
33  const unsigned int nLineSegments = 20;
34  float circleScalingFactor = 50;
35  const double jetRadius = iData.pt()/circleScalingFactor;
36 
37  for ( unsigned int iphi = 0; iphi < nLineSegments; ++iphi ) {
38  AddLine(iData.eta()+jetRadius*cos(2*M_PI/nLineSegments*iphi),
39  iData.phi()+jetRadius*sin(2*M_PI/nLineSegments*iphi),
40  base,
41  iData.eta()+jetRadius*cos(2*M_PI/nLineSegments*(iphi+1)),
42  iData.phi()+jetRadius*sin(2*M_PI/nLineSegments*(iphi+1)),
43  base);
44  }
45 }
46 
47 void
49 {
50  FWViewEnergyScale* caloScale = vc->getEnergyScale();
51  float val = caloScale->getPlotEt() ? m_et : m_energy;
52 
53  // printf("update scale %f \n", getScale(vc, context)); fflush(stdout);
54  float scaleFac = caloScale->getScaleFactorLego();
55  // resize first line
56  TEveChunkManager::iterator li(GetLinePlex());
57  li.next();
58  TEveStraightLineSet::Line_t& l = * (TEveStraightLineSet::Line_t*) li();
59  l.fV2[2] = l.fV1[2] + val*scaleFac;
60 
61  // move end point
62  TEveChunkManager::iterator mi(GetMarkerPlex());
63  mi.next();
64  TEveStraightLineSet::Marker_t& m = * (TEveStraightLineSet::Marker_t*) mi();
65  m.fV[2] = l.fV2[2];
66 }
67 
tuple base
Main Program
Definition: newFWLiteAna.py:91
virtual double energy() const final
energy
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
FWViewEnergyScale * getEnergyScale() const
virtual double phi() const final
momentum azimuthal angle
void updateScale(const FWViewContext *, const fireworks::Context &)
float getScaleFactorLego() const
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double f[11][100]
#define M_PI
FWLegoEvePFCandidate(const reco::PFCandidate &pfc, const FWViewContext *, const fireworks::Context &)
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
virtual double et() const final
transverse energy
virtual double eta() const final
momentum pseudorapidity
bool getPlotEt() const
virtual double pt() const final
transverse momentum