src
JetMETCorrections
FFTJetObjects
interface
L2AbsScaleCalculator.h
Go to the documentation of this file.
1
#ifndef JetMETCorrections_FFTJetObjects_L2AbsScaleCalculator_h
2
#define JetMETCorrections_FFTJetObjects_L2AbsScaleCalculator_h
3
4
#include <cmath>
5
#include <cassert>
6
7
#include "
JetMETCorrections/FFTJetObjects/interface/AbsFFTSpecificScaleCalculator.h
"
8
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
9
#include "
FWCore/Utilities/interface/Exception.h
"
10
11
class
L2AbsScaleCalculator
:
public
AbsFFTSpecificScaleCalculator
{
12
public
:
13
inline
explicit
L2AbsScaleCalculator
(
const
edm::ParameterSet
& ps)
14
:
m_radiusFactor
(ps.getParameter<double>(
"radiusFactor"
)),
15
m_zeroPtLog
(ps.getParameter<double>(
"zeroPtLog"
)),
16
m_takePtLog
(ps.getParameter<
bool
>(
"takePtLog"
)) {}
17
18
inline
~L2AbsScaleCalculator
()
override
{}
19
20
inline
void
mapFFTJet
(
const
reco::Jet
&
/* jet */
,
21
const
reco::FFTJet<float>
& fftJet,
22
const
math::XYZTLorentzVector
& current,
23
double
*
buf
,
24
const
unsigned
dim)
const override
{
25
if
(dim != 2)
26
throw
cms::Exception
(
"FFTJetBadConfig"
) <<
"In L2AbsScaleCalculator::mapFFTJet: "
27
<<
"invalid table dimensionality: "
<< dim << std::endl;
28
assert
(
buf
);
29
const
double
radius
= fftJet.
f_recoScale
();
30
const
double
pt
= current.pt();
31
buf
[0] =
radius
*
m_radiusFactor
;
32
if
(
m_takePtLog
) {
33
if
(
pt
> 0.0)
34
buf
[1] =
log
(
pt
);
35
else
36
buf
[1] =
m_zeroPtLog
;
37
}
else
38
buf
[1] =
pt
;
39
}
40
41
private
:
42
double
m_radiusFactor
;
43
double
m_zeroPtLog
;
44
bool
m_takePtLog
;
45
};
46
47
#endif // JetMETCorrections_FFTJetObjects_L2AbsScaleCalculator_h
Exception
Definition:
hltDiff.cc:245
L2AbsScaleCalculator::m_radiusFactor
double m_radiusFactor
Definition:
L2AbsScaleCalculator.h:42
reco::Jet
Base class for all types of Jets.
Definition:
Jet.h:20
L2AbsScaleCalculator::~L2AbsScaleCalculator
~L2AbsScaleCalculator() override
Definition:
L2AbsScaleCalculator.h:18
cms::cuda::assert
assert(be >=bs)
DiDispStaMuonMonitor_cfi.pt
pt
Definition:
DiDispStaMuonMonitor_cfi.py:39
ParameterSet.h
math::XYZTLorentzVector
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition:
LorentzVector.h:29
visDQMUpload.buf
buf
Definition:
visDQMUpload.py:153
CosmicsPD_Skims.radius
radius
Definition:
CosmicsPD_Skims.py:135
L2AbsScaleCalculator::m_zeroPtLog
double m_zeroPtLog
Definition:
L2AbsScaleCalculator.h:43
L2AbsScaleCalculator
Definition:
L2AbsScaleCalculator.h:11
L2AbsScaleCalculator::mapFFTJet
void mapFFTJet(const reco::Jet &, const reco::FFTJet< float > &fftJet, const math::XYZTLorentzVector ¤t, double *buf, const unsigned dim) const override
Definition:
L2AbsScaleCalculator.h:20
electrons_cff.bool
bool
Definition:
electrons_cff.py:319
reco::FFTJet< float >
L2AbsScaleCalculator::m_takePtLog
bool m_takePtLog
Definition:
L2AbsScaleCalculator.h:44
Exception.h
AbsFFTSpecificScaleCalculator
Definition:
AbsFFTSpecificScaleCalculator.h:13
edm::ParameterSet
Definition:
ParameterSet.h:47
L2AbsScaleCalculator::L2AbsScaleCalculator
L2AbsScaleCalculator(const edm::ParameterSet &ps)
Definition:
L2AbsScaleCalculator.h:13
dqm-mbProfile.log
log
Definition:
dqm-mbProfile.py:17
AbsFFTSpecificScaleCalculator.h
reco::FFTJet::f_recoScale
Real f_recoScale() const
Definition:
FFTJet.h:85
Generated for CMSSW Reference Manual by
1.8.14