RecoEcal
EgammaCoreTools
src
SCDynamicDPhiParametersHelper.cc
Go to the documentation of this file.
1
// Implementation of the SC dynamic dPhi parameters interface
2
3
#include "
RecoEcal/EgammaCoreTools/interface/SCDynamicDPhiParametersHelper.h
"
4
5
#include <algorithm>
6
#include <utility>
7
8
using namespace
reco
;
9
10
SCDynamicDPhiParametersHelper::SCDynamicDPhiParametersHelper
(
EcalSCDynamicDPhiParameters
&
params
,
11
const
edm::ParameterSet
&iConfig)
12
: parameters_(
params
) {
13
// dynamic dPhi parameters
14
// clear the vector in case the EcalMustacheSCParameters had been initialised before
15
if
(!
parameters_
.
dynamicDPhiParametersCollection_
.empty()) {
16
parameters_
.
dynamicDPhiParametersCollection_
.clear();
17
}
18
const
auto
dynamicDPhiPSets = iConfig.
getParameter
<std::vector<edm::ParameterSet>>(
"dynamicDPhiParameterSets"
);
19
for
(
const
auto
&pSet : dynamicDPhiPSets) {
20
EcalSCDynamicDPhiParameters::DynamicDPhiParameters
dynDPhiParams({pSet.getParameter<
double
>(
"eMin"
),
21
pSet.getParameter<
double
>(
"etaMin"
),
22
pSet.getParameter<
double
>(
"yoffset"
),
23
pSet.getParameter<
double
>(
"scale"
),
24
pSet.getParameter<
double
>(
"xoffset"
),
25
pSet.getParameter<
double
>(
"width"
),
26
pSet.getParameter<
double
>(
"saturation"
),
27
pSet.getParameter<
double
>(
"cutoff"
)});
28
addDynamicDPhiParameters
(dynDPhiParams);
29
sortDynamicDPhiParametersCollection
();
30
}
31
}
32
33
void
SCDynamicDPhiParametersHelper::addDynamicDPhiParameters
(
34
const
EcalSCDynamicDPhiParameters::DynamicDPhiParameters
&dynDPhiParams) {
35
parameters_
.
dynamicDPhiParametersCollection_
.emplace_back(dynDPhiParams);
36
}
37
38
void
SCDynamicDPhiParametersHelper::sortDynamicDPhiParametersCollection
() {
39
std::sort
(
parameters_
.
dynamicDPhiParametersCollection_
.begin(),
40
parameters_
.
dynamicDPhiParametersCollection_
.end(),
41
[](
const
EcalSCDynamicDPhiParameters::DynamicDPhiParameters
&
p1
,
42
const
EcalSCDynamicDPhiParameters::DynamicDPhiParameters
&
p2
) {
43
const
auto
p1Mins = std::make_pair(
p1
.eMin,
p1
.etaMin);
44
const
auto
p2Mins = std::make_pair(
p2
.eMin,
p2
.etaMin);
45
return
p1Mins < p2Mins;
46
});
47
}
CalibrationSummaryClient_cfi.params
params
Definition:
CalibrationSummaryClient_cfi.py:14
SCDynamicDPhiParametersHelper.h
reco::SCDynamicDPhiParametersHelper::SCDynamicDPhiParametersHelper
SCDynamicDPhiParametersHelper(EcalSCDynamicDPhiParameters ¶ms, const edm::ParameterSet &iConfig)
Definition:
SCDynamicDPhiParametersHelper.cc:10
reco::SCDynamicDPhiParametersHelper::sortDynamicDPhiParametersCollection
void sortDynamicDPhiParametersCollection()
Definition:
SCDynamicDPhiParametersHelper.cc:38
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:45
EcalSCDynamicDPhiParameters::DynamicDPhiParameters
Definition:
EcalSCDynamicDPhiParameters.h:18
EcalSCDynamicDPhiParameters
Definition:
EcalSCDynamicDPhiParameters.h:13
p2
double p2[4]
Definition:
TauolaWrapper.h:90
edm::ParameterSet
Definition:
ParameterSet.h:47
jetUpdater_cfi.sort
sort
Definition:
jetUpdater_cfi.py:29
reco::SCDynamicDPhiParametersHelper::addDynamicDPhiParameters
void addDynamicDPhiParameters(const EcalSCDynamicDPhiParameters::DynamicDPhiParameters &dynDPhiParams)
Definition:
SCDynamicDPhiParametersHelper.cc:33
p1
double p1[4]
Definition:
TauolaWrapper.h:89
EcalSCDynamicDPhiParameters::dynamicDPhiParametersCollection_
std::vector< DynamicDPhiParameters > dynamicDPhiParametersCollection_
Definition:
EcalSCDynamicDPhiParameters.h:45
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition:
ParameterSet.h:303
reco::SCDynamicDPhiParametersHelper::parameters_
EcalSCDynamicDPhiParameters & parameters_
Definition:
SCDynamicDPhiParametersHelper.h:18
Generated for CMSSW Reference Manual by
1.8.16