CMS 3D CMS Logo

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