Main Page
Namespaces
Classes
Package Documentation
JetMETCorrections
Modules
plugins
QGLikelihoodSystematicsDBReader.cc
Go to the documentation of this file.
1
#include <memory>
2
#include "
FWCore/Framework/interface/Frameworkfwd.h
"
3
#include "
FWCore/Framework/interface/EDAnalyzer.h
"
4
#include "
FWCore/Framework/interface/Event.h
"
5
#include "
FWCore/Framework/interface/MakerMacros.h
"
6
#include "
FWCore/Framework/interface/EventSetup.h
"
7
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
8
#include "
FWCore/Framework/interface/ESHandle.h
"
9
#include "
CondFormats/JetMETObjects/interface/QGLikelihoodObject.h
"
10
#include "
JetMETCorrections/Objects/interface/JetCorrectionsRecord.h
"
11
#include "
CondFormats/DataRecord/interface/QGLikelihoodSystematicsRcd.h
"
12
13
class
QGLikelihoodSystematicsDBReader
:
public
edm::EDAnalyzer
{
14
public
:
15
explicit
QGLikelihoodSystematicsDBReader
(
const
edm::ParameterSet
&);
16
~QGLikelihoodSystematicsDBReader
()
override
{};
17
18
private
:
19
void
beginJob
()
override
{};
20
void
analyze
(
const
edm::Event
&,
const
edm::EventSetup
&)
override
;
21
void
endJob
()
override
{};
22
23
std::string
mPayloadName
;
24
bool
mCreateTextFile
,
mPrintScreen
;
25
};
26
27
28
QGLikelihoodSystematicsDBReader::QGLikelihoodSystematicsDBReader
(
const
edm::ParameterSet
& iConfig){
29
mPayloadName
= iConfig.
getUntrackedParameter
<
std::string
>(
"payloadName"
);
30
mPrintScreen
= iConfig.
getUntrackedParameter
<
bool
>(
"printScreen"
);
31
mCreateTextFile
= iConfig.
getUntrackedParameter
<
bool
>(
"createTextFile"
);
32
}
33
34
35
void
QGLikelihoodSystematicsDBReader::analyze
(
const
edm::Event
&
iEvent
,
const
edm::EventSetup
& iSetup){
36
edm::LogInfo
(
"UserOutput"
) <<
"Inspecting QGLikelihood payload with label:"
<<
mPayloadName
<< std::endl;
37
edm::ESHandle<QGLikelihoodSystematicsObject>
QGLSysPar;
38
QGLikelihoodSystematicsRcd
const
& rcdhandle = iSetup.
get
<
QGLikelihoodSystematicsRcd
>();
39
rcdhandle.
get
(
mPayloadName
, QGLSysPar);
40
41
std::vector<QGLikelihoodSystematicsObject::Entry>
const
&
data
= QGLSysPar->
data
;
42
edm::LogInfo
(
"UserOutput"
) <<
"There are "
<< data.size() <<
" entries (categories with parameters for smearing):"
<< std::endl;
43
for
(
auto
idata = data.begin(); idata != data.end(); ++idata){
44
int
qgBin = idata->systCategory.QGIndex;
45
double
etaMin
= idata->systCategory.EtaMin;
46
double
etaMax
= idata->systCategory.EtaMax;
47
double
rhoMin
= idata->systCategory.RhoMin;
48
double
rhoMax
= idata->systCategory.RhoMax;
49
double
ptMin
= idata->systCategory.PtMin;
50
double
ptMax
= idata->systCategory.PtMax;
51
double
a
= idata->a;
52
double
b
= idata->b;
53
double
lmin = idata->lmin;
54
double
lmax = idata->lmax;
55
56
char
buff[1000];
57
sprintf(buff,
"qg=%1d, ptMin=%8.2f, ptMax=%8.2f, etaMin=%3.1f, etaMax=%3.1f, rhoMin=%6.2f, rhoMax=%6.2f, a=%7.3f, b=%7.3f, lmin=%6.2f, lmax=%6.2f"
, qgBin, ptMin, ptMax, etaMin, etaMax, rhoMin, rhoMax, a, b, lmin, lmax);
58
edm::LogVerbatim
(
"UserOutput"
) << buff << std::endl;
59
}
60
}
61
62
DEFINE_FWK_MODULE
(
QGLikelihoodSystematicsDBReader
);
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
QGLikelihoodSystematicsDBReader::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition:
QGLikelihoodSystematicsDBReader.cc:35
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition:
MakerMacros.h:17
Event.h
egammaForCoreTracking_cff.rhoMax
rhoMax
Definition:
egammaForCoreTracking_cff.py:32
QGLikelihoodSystematicsDBReader::endJob
void endJob() override
Definition:
QGLikelihoodSystematicsDBReader.cc:21
MakerMacros.h
EventSetup.h
QGLikelihoodSystematicsDBReader::beginJob
void beginJob() override
Definition:
QGLikelihoodSystematicsDBReader.cc:19
QGLikelihoodObject.h
QGLikelihoodSystematicsDBReader::mCreateTextFile
bool mCreateTextFile
Definition:
QGLikelihoodSystematicsDBReader.cc:24
ALCARECOTkAlBeamHalo_cff.etaMin
etaMin
GeV.
Definition:
ALCARECOTkAlBeamHalo_cff.py:32
QGLikelihoodSystematicsDBReader::~QGLikelihoodSystematicsDBReader
~QGLikelihoodSystematicsDBReader() override
Definition:
QGLikelihoodSystematicsDBReader.cc:16
Frameworkfwd.h
QGLikelihoodSystematicsRcd.h
ParameterSet.h
QGLikelihoodSystematicsObject::data
std::vector< Entry > data
Definition:
QGLikelihoodObject.h:44
iEvent
int iEvent
Definition:
GenABIO.cc:230
edm::ESHandle
Definition:
DTSurvey.h:22
ESHandle.h
edm::EventSetup
Definition:
EventSetup.h:54
edm::LogVerbatim
Definition:
MessageLogger.h:271
ALCARECOTkAlBeamHalo_cff.etaMax
etaMax
Definition:
ALCARECOTkAlBeamHalo_cff.py:33
ptMin
float ptMin
Definition:
PhotonIDValueMapProducer.cc:139
QGLikelihoodSystematicsDBReader::mPrintScreen
bool mPrintScreen
Definition:
QGLikelihoodSystematicsDBReader.cc:24
edm::EDAnalyzer
Definition:
EDAnalyzer.h:28
EDAnalyzer.h
edm::LogInfo
Definition:
MessageLogger.h:238
QGLikelihoodSystematicsRcd
Definition:
QGLikelihoodSystematicsRcd.h:23
AlignmentTrackSelector_cfi.ptMax
ptMax
Definition:
AlignmentTrackSelector_cfi.py:12
b
double b
Definition:
hdecay.h:120
QGLikelihoodSystematicsDBReader
Definition:
QGLikelihoodSystematicsDBReader.cc:13
data
char data[epos_bytes_allocation]
Definition:
EPOS_Wrapper.h:82
a
double a
Definition:
hdecay.h:121
edm::EventSetup::get
T get() const
Definition:
EventSetup.h:68
QGLikelihoodSystematicsDBReader::QGLikelihoodSystematicsDBReader
QGLikelihoodSystematicsDBReader(const edm::ParameterSet &)
Definition:
QGLikelihoodSystematicsDBReader.cc:28
edm::ParameterSet
Definition:
ParameterSet.h:36
edm::Event
Definition:
Event.h:70
electronConversionRejectionValidator.rhoMin
rhoMin
Definition:
electronConversionRejectionValidator.py:59
edm::eventsetup::EventSetupRecord::get
bool get(HolderT &iHolder) const
Definition:
EventSetupRecord.h:109
QGLikelihoodSystematicsDBReader::mPayloadName
std::string mPayloadName
Definition:
QGLikelihoodSystematicsDBReader.cc:21
JetCorrectionsRecord.h
Generated for CMSSW Reference Manual by
1.8.11