CondFormats
PPSObjects
src
LHCInterpolatedOpticalFunctionsSet.cc
Go to the documentation of this file.
1
// Original Author: Jan Kašpar
2
3
#include "
CondFormats/PPSObjects/interface/LHCInterpolatedOpticalFunctionsSet.h
"
4
5
//----------------------------------------------------------------------------------------------------
6
7
void
LHCInterpolatedOpticalFunctionsSet::initializeSplines
() {
8
const
unsigned
int
num_xi_vals =
m_xi_values
.size();
9
10
m_splines
.resize(
m_fcn_values
.size());
11
for
(
unsigned
int
i
= 0;
i
<
m_fcn_values
.size(); ++
i
)
12
m_splines
[
i
] = std::make_shared<TSpline3>(
""
,
m_xi_values
.data(),
m_fcn_values
[
i
].data(), num_xi_vals);
13
}
14
15
//----------------------------------------------------------------------------------------------------
16
17
void
LHCInterpolatedOpticalFunctionsSet::transport
(
const
LHCInterpolatedOpticalFunctionsSet::Kinematics
&
input
,
18
LHCInterpolatedOpticalFunctionsSet::Kinematics
&
output
,
19
bool
calculateAngles)
const
{
20
const
double
xi
=
input
.xi;
21
22
output
.x =
m_splines
[
exd
]->Eval(
xi
) +
m_splines
[
evx
]->Eval(
xi
) *
input
.x +
m_splines
[
eLx
]->Eval(
xi
) *
input
.th_x +
23
m_splines
[
e14
]->Eval(
xi
) *
input
.th_y;
24
25
output
.th_x = (!calculateAngles) ? 0.
26
:
m_splines
[
expd
]->Eval(
xi
) +
m_splines
[
evpx
]->Eval(
xi
) *
input
.x +
27
m_splines
[
eLpx
]->Eval(
xi
) *
input
.th_x +
m_splines
[
e24
]->Eval(
xi
) *
input
.th_y;
28
29
output
.y =
m_splines
[
eyd
]->Eval(
xi
) +
m_splines
[
evy
]->Eval(
xi
) *
input
.y +
m_splines
[
eLy
]->Eval(
xi
) *
input
.th_y +
30
m_splines
[
e32
]->Eval(
xi
) *
input
.th_x;
31
32
output
.th_y = (!calculateAngles) ? 0.
33
:
m_splines
[
eypd
]->Eval(
xi
) +
m_splines
[
evpy
]->Eval(
xi
) *
input
.y +
34
m_splines
[
eLpy
]->Eval(
xi
) *
input
.th_y +
m_splines
[
e42
]->Eval(
xi
) *
input
.th_x;
35
36
output
.xi =
input
.xi;
37
}
LHCOpticalFunctionsSet::exd
Definition:
LHCOpticalFunctionsSet.h:16
mps_fire.i
i
Definition:
mps_fire.py:355
input
static const std::string input
Definition:
EdmProvDump.cc:48
LHCOpticalFunctionsSet::expd
Definition:
LHCOpticalFunctionsSet.h:16
LHCOpticalFunctionsSet::e42
Definition:
LHCOpticalFunctionsSet.h:16
LHCInterpolatedOpticalFunctionsSet::Kinematics
proton kinematics description
Definition:
LHCInterpolatedOpticalFunctionsSet.h:28
LHCOpticalFunctionsSet::e24
Definition:
LHCOpticalFunctionsSet.h:16
convertSQLitetoXML_cfg.output
output
Definition:
convertSQLitetoXML_cfg.py:32
LHCOpticalFunctionsSet::m_fcn_values
std::vector< std::vector< double > > m_fcn_values
length unit cm
Definition:
LHCOpticalFunctionsSet.h:36
LHCOpticalFunctionsSet::evx
Definition:
LHCOpticalFunctionsSet.h:16
hybridSuperClusters_cfi.xi
xi
Definition:
hybridSuperClusters_cfi.py:10
LHCInterpolatedOpticalFunctionsSet.h
LHCOpticalFunctionsSet::eLx
Definition:
LHCOpticalFunctionsSet.h:16
LHCOpticalFunctionsSet::eLy
Definition:
LHCOpticalFunctionsSet.h:16
LHCOpticalFunctionsSet::evy
Definition:
LHCOpticalFunctionsSet.h:16
LHCOpticalFunctionsSet::e32
Definition:
LHCOpticalFunctionsSet.h:16
LHCOpticalFunctionsSet::eLpy
Definition:
LHCOpticalFunctionsSet.h:16
LHCOpticalFunctionsSet::evpy
Definition:
LHCOpticalFunctionsSet.h:16
LHCOpticalFunctionsSet::eyd
Definition:
LHCOpticalFunctionsSet.h:16
LHCOpticalFunctionsSet::evpx
Definition:
LHCOpticalFunctionsSet.h:16
LHCOpticalFunctionsSet::eLpx
Definition:
LHCOpticalFunctionsSet.h:16
LHCInterpolatedOpticalFunctionsSet::initializeSplines
void initializeSplines()
builds splines from m_*_values fields
Definition:
LHCInterpolatedOpticalFunctionsSet.cc:7
LHCOpticalFunctionsSet::eypd
Definition:
LHCOpticalFunctionsSet.h:16
LHCInterpolatedOpticalFunctionsSet::m_splines
std::vector< std::shared_ptr< const TSpline3 > > m_splines
Definition:
LHCInterpolatedOpticalFunctionsSet.h:43
LHCInterpolatedOpticalFunctionsSet::transport
void transport(const Kinematics &input, Kinematics &output, bool calculateAngles=false) const
transports proton according to the splines
Definition:
LHCInterpolatedOpticalFunctionsSet.cc:17
LHCOpticalFunctionsSet::m_xi_values
std::vector< double > m_xi_values
Definition:
LHCOpticalFunctionsSet.h:35
LHCOpticalFunctionsSet::e14
Definition:
LHCOpticalFunctionsSet.h:16
Generated for CMSSW Reference Manual by
1.8.16