RecoTracker
TkMSParametrization
src
MultipleScatteringX0Data.cc
Go to the documentation of this file.
1
#include "
MultipleScatteringX0Data.h
"
2
#include "
FWCore/Utilities/interface/Exception.h
"
3
#include "
FWCore/ParameterSet/interface/FileInPath.h
"
4
5
#include "TH2.h"
6
#include "TFile.h"
7
#include "TKey.h"
8
9
#include <iostream>
10
#include <string>
11
#include <cmath>
12
13
using namespace
std
;
14
15
MultipleScatteringX0Data::MultipleScatteringX0Data
() : theData(nullptr) {
16
string
filename
=
fileName
();
17
TFile
theFile
(
filename
.c_str(),
"READ"
);
18
if
(not
theFile
.IsZombie()) {
19
theData
.reset(dynamic_cast<TH2F*>(
theFile
.GetKey(
"h100"
)->ReadObj()));
20
theData
->SetDirectory(
nullptr
);
21
}
22
if
(!
theData
) {
23
throw
cms::Exception
(
"Data not found"
)
24
<<
" ** MultipleScatteringX0Data ** file: "
<<
filename
<<
" <-- no data found!!!"
;
25
}
26
}
27
28
MultipleScatteringX0Data::~MultipleScatteringX0Data
() {}
29
30
string
MultipleScatteringX0Data::fileName
() {
31
string
defName =
"RecoTracker/TkMSParametrization/data/MultipleScatteringX0Data.root"
;
32
edm::FileInPath
f
(defName);
33
return
f
.fullPath();
34
}
35
36
int
MultipleScatteringX0Data::nBinsEta
()
const
{
37
if
(
theData
)
38
return
theData
->GetNbinsX();
39
else
40
return
0;
41
}
42
43
float
MultipleScatteringX0Data::minEta
()
const
{
44
if
(
theData
)
45
return
theData
->GetXaxis()->GetXmin();
46
else
47
return
0;
48
}
49
50
float
MultipleScatteringX0Data::maxEta
()
const
{
51
if
(
theData
)
52
return
theData
->GetXaxis()->GetXmax();
53
else
54
return
0;
55
}
56
57
float
MultipleScatteringX0Data::sumX0atEta
(
float
eta
,
float
r
)
const
{
58
if
(!
theData
)
59
return
0.;
60
eta
= fabs(
eta
);
61
62
int
ieta
=
theData
->GetXaxis()->FindBin(
eta
);
63
int
irad =
theData
->GetYaxis()->FindBin(
r
);
64
65
if
(irad < theData->GetNbinsY()) {
66
return
theData
->GetBinContent(
ieta
, irad);
67
}
else
{
68
float
sumX0 = 0;
69
for
(
int
ir =
theData
->GetNbinsY(); ir > 0; ir--) {
70
float
val
=
theData
->GetBinContent(
ieta
, ir);
71
if
(
val
> sumX0)
72
sumX0 =
val
;
73
}
74
return
sumX0;
75
}
76
}
f
double f[11][100]
Definition:
MuScleFitUtils.cc:78
MultipleScatteringX0Data::minEta
float minEta() const
Definition:
MultipleScatteringX0Data.cc:43
MultipleScatteringX0Data::maxEta
float maxEta() const
Definition:
MultipleScatteringX0Data.cc:50
FileInPath.h
MultipleScatteringX0Data::fileName
std::string fileName()
Definition:
MultipleScatteringX0Data.cc:30
MultipleScatteringX0Data::nBinsEta
int nBinsEta() const
Definition:
MultipleScatteringX0Data.cc:36
edm::FileInPath
Definition:
FileInPath.h:64
interactiveExample.theFile
theFile
Definition:
interactiveExample.py:17
PVValHelper::eta
Definition:
PVValidationHelpers.h:69
corrVsCorr.filename
filename
Definition:
corrVsCorr.py:123
MultipleScatteringX0Data::sumX0atEta
float sumX0atEta(float eta, float r) const override
Definition:
MultipleScatteringX0Data.cc:57
LEDCalibrationChannels.ieta
ieta
Definition:
LEDCalibrationChannels.py:63
alignCSCRings.r
r
Definition:
alignCSCRings.py:93
MultipleScatteringX0Data::~MultipleScatteringX0Data
~MultipleScatteringX0Data() override
Definition:
MultipleScatteringX0Data.cc:28
heppy_batch.val
val
Definition:
heppy_batch.py:351
std
Definition:
JetResolutionObject.h:76
Exception
Definition:
hltDiff.cc:246
MultipleScatteringX0Data::theData
std::unique_ptr< TH2F > theData
Definition:
MultipleScatteringX0Data.h:31
MultipleScatteringX0Data.h
Exception.h
MultipleScatteringX0Data::MultipleScatteringX0Data
MultipleScatteringX0Data()
Definition:
MultipleScatteringX0Data.cc:15
Generated for CMSSW Reference Manual by
1.8.16