src
Fireworks
Candidates
src
FWLegoCandidate.cc
Go to the documentation of this file.
1
#include "
Fireworks/Candidates/interface/FWLegoCandidate.h
"
2
3
#include "
Fireworks/Core/interface/Context.h
"
4
#include "
Fireworks/Core/interface/FWViewContext.h
"
5
#include "
Fireworks/Core/interface/FWViewEnergyScale.h
"
6
7
//______________________________________________________________________________
8
FWLegoCandidate::FWLegoCandidate
(
9
const
FWViewContext
*vc,
const
fireworks::Context
&
context
,
float
energy
,
float
et
,
float
pt
,
float
eta
,
float
phi
)
10
: m_energy(
energy
), m_et(
et
), m_pt(
pt
), m_eta(
eta
), m_phi(
phi
) {
11
float
base
= 0.001;
// Floor offset 1%
12
13
// First vertical line
14
FWViewEnergyScale
*caloScale = vc->
getEnergyScale
();
15
float
val
= caloScale->
getPlotEt
() ?
m_pt
:
m_energy
;
// Use pt instead of et
16
17
AddLine(
m_eta
,
m_phi
,
base
,
m_eta
,
m_phi
,
base
+
val
* caloScale->
getScaleFactorLego
());
18
19
AddMarker(0, 1.
f
);
20
SetMarkerStyle(3);
21
SetMarkerSize(0.01);
22
SetDepthTest(
false
);
23
24
// Circle pt
25
const
unsigned
int
nLineSegments = 20;
26
const
double
radius
=
log
(1 +
m_pt
) /
log
(10) / 30.f;
27
//const double radius = m_pt / 100.f;
28
const
double
twoPi
= 2 *
TMath::Pi
();
29
30
for
(
unsigned
int
iPhi = 0; iPhi < nLineSegments; ++iPhi) {
31
AddLine(
m_eta
+
radius
*
cos
(
twoPi
/ nLineSegments * iPhi),
32
m_phi
+
radius
*
sin
(
twoPi
/ nLineSegments * iPhi),
33
base
,
34
m_eta
+
radius
*
cos
(
twoPi
/ nLineSegments * (iPhi + 1)),
35
m_phi
+
radius
*
sin
(
twoPi
/ nLineSegments * (iPhi + 1)),
36
base
);
37
}
38
}
39
40
//______________________________________________________________________________
41
void
FWLegoCandidate::updateScale
(
const
FWViewContext
*vc,
const
fireworks::Context
&
context
) {
42
FWViewEnergyScale
*caloScale = vc->
getEnergyScale
();
43
float
val
= caloScale->
getPlotEt
() ?
m_pt
:
m_energy
;
// Use pt instead of et
44
float
scaleFac = caloScale->
getScaleFactorLego
();
45
46
// Resize first line
47
TEveChunkManager::iterator li(GetLinePlex());
48
li.next();
49
TEveStraightLineSet::Line_t &
l
= *(TEveStraightLineSet::Line_t *)li();
50
l
.fV2[2] =
l
.fV1[2] +
val
* scaleFac;
51
52
// Move end point (marker)
53
TEveChunkManager::iterator mi(GetMarkerPlex());
54
mi.next();
55
TEveStraightLineSet::Marker_t &
m
= *(TEveStraightLineSet::Marker_t *)mi();
56
m
.fV[2] =
l
.fV2[2];
// Set to new top of line
57
}
Pi
const double Pi
Definition:
CosmicMuonParameters.h:18
FWLegoCandidate::FWLegoCandidate
FWLegoCandidate()
Definition:
FWLegoCandidate.h:38
FWLegoCandidate::m_energy
float m_energy
Definition:
FWLegoCandidate.h:49
PVValHelper::phi
Definition:
PVValidationHelpers.h:69
FWLegoCandidate.h
edmMakeDummyCfis.base
base
Definition:
edmMakeDummyCfis.py:22
funct::sin
Sin< T >::type sin(const T &t)
Definition:
Sin.h:22
FWLegoCandidate::m_phi
float m_phi
Definition:
FWLegoCandidate.h:53
PVValHelper::eta
Definition:
PVValidationHelpers.h:70
FWViewContext
Definition:
FWViewContext.h:32
Context.h
hcalRecHitTable_cff.energy
energy
Definition:
hcalRecHitTable_cff.py:13
FWLegoCandidate::m_eta
float m_eta
Definition:
FWLegoCandidate.h:52
DiDispStaMuonMonitor_cfi.pt
pt
Definition:
DiDispStaMuonMonitor_cfi.py:39
FWViewEnergyScale::getPlotEt
bool getPlotEt() const
Definition:
FWViewEnergyScale.h:47
FWViewEnergyScale::getScaleFactorLego
float getScaleFactorLego() const
Definition:
FWViewEnergyScale.h:45
FWViewContext::getEnergyScale
FWViewEnergyScale * getEnergyScale() const
Definition:
FWViewContext.cc:25
CosmicsPD_Skims.radius
radius
Definition:
CosmicsPD_Skims.py:135
FWLegoCandidate::m_pt
float m_pt
Definition:
FWLegoCandidate.h:51
funct::cos
Cos< T >::type cos(const T &t)
Definition:
Cos.h:22
visualization-live-secondInstance_cfg.m
m
Definition:
visualization-live-secondInstance_cfg.py:84
fireworks::Context
Definition:
Context.h:41
f
double f[11][100]
Definition:
MuScleFitUtils.cc:78
visDQMUpload.context
context
Definition:
visDQMUpload.py:30
FWLegoCandidate::updateScale
void updateScale(const FWViewContext *vc, const fireworks::Context &)
Definition:
FWLegoCandidate.cc:41
FWViewEnergyScale
Definition:
FWViewEnergyScale.h:34
FWViewEnergyScale.h
FWViewContext.h
dqm-mbProfile.log
log
Definition:
dqm-mbProfile.py:17
Geom::twoPi
constexpr double twoPi()
Definition:
Pi.h:32
heppy_batch.val
val
Definition:
heppy_batch.py:351
l1tnanotables_cff.et
et
Definition:
l1tnanotables_cff.py:32
MainPageGenerator.l
l
Definition:
MainPageGenerator.py:429
Generated for CMSSW Reference Manual by
1.8.14