Mantid
Loading...
Searching...
No Matches
Framework
CurveFitting
src
EigenVectorView.cpp
Go to the documentation of this file.
1
// Mantid Repository : https://github.com/mantidproject/mantid
2
//
3
// Copyright © 2022 ISIS Rutherford Appleton Laboratory UKRI,
4
// NScD Oak Ridge National Laboratory, European Spallation Source,
5
// Institut Laue - Langevin & CSNS, Institute of High Energy Physics, CAS
6
// SPDX - License - Identifier: GPL - 3.0 +
7
#include "
MantidCurveFitting/EigenVectorView.h
"
8
9
namespace
Mantid::CurveFitting
{
10
// EigenVector_View Constructors
11
// default constructor
12
EigenVector_View::EigenVector_View
() : m_view({}, 0,
dynamic_stride
(0, 0)) {}
13
14
// constructor: map->vector view
18
EigenVector_View::EigenVector_View
(
vec_map_type
&vector,
const
int
nElements,
const
size_t
startElement)
19
: m_view(vector.data(), vector.size(),
dynamic_stride
(0, 1)) {
20
if
(nElements == -1) {
21
// if nElements is default, do nothing as m_view is initialised as such
22
}
else
{
23
new
(&
m_view
)
vec_map_type
(vector.data() + startElement, nElements,
dynamic_stride
(0, 1));
24
}
25
}
26
27
// constructor: vector->vector view
31
EigenVector_View::EigenVector_View
(Eigen::VectorXd &vector,
const
int
nElements,
const
size_t
startElement)
32
: m_view(vector.data(), vector.size(),
dynamic_stride
(0, 1)) {
33
if
(nElements == -1) {
34
// if nElements is default, do nothing as m_view is initialised as such
35
}
else
{
36
new
(&
m_view
)
vec_map_type
(vector.data() + startElement, nElements,
dynamic_stride
(0, 1));
37
}
38
}
39
40
// constructor: array->vector view
44
EigenVector_View::EigenVector_View
(
double
*base,
const
size_t
nElements,
const
size_t
startElement)
45
: m_view(base + startElement, nElements,
dynamic_stride
(0, 1)) {}
46
47
// CONST constructor: map->vector view
51
EigenVector_View::EigenVector_View
(
const
vec_map_type
&vector,
const
size_t
nElements,
const
size_t
startElement)
52
: m_view({}, 0,
dynamic_stride
(0, 0)) {
53
new
(&m_view)
vec_const_map_type
(vector.data() + startElement, nElements,
dynamic_stride
(0, 1));
54
}
55
56
// CONST constructor: vector->vector view
60
EigenVector_View::EigenVector_View
(
const
Eigen::VectorXd &vector,
const
size_t
nElements,
const
size_t
startElement)
61
: m_view({}, 0,
dynamic_stride
(0, 0)), m_isConst(
true
) {
62
new
(&
m_view
)
vec_const_map_type
(vector.data() + startElement, nElements,
dynamic_stride
(0, 1));
63
}
64
65
// CONST constructor: array->vector view
69
EigenVector_View::EigenVector_View
(
const
double
*base,
const
size_t
nElements,
const
size_t
startElement)
70
: m_view({}, 0,
dynamic_stride
(0, 0)), m_isConst(
true
) {
71
new
(&
m_view
)
vec_const_map_type
(base + startElement, nElements,
dynamic_stride
(0, 1));
72
}
73
74
// CONST copy constructor
76
EigenVector_View::EigenVector_View
(
const
EigenVector_View
&v) : m_view({}, 0,
dynamic_stride
(0, 0)), m_isConst(
true
) {
77
new
(&
m_view
)
vec_const_map_type
(v.vector_inspector().data(), v.size(),
dynamic_stride
(0, 1));
78
}
79
80
// copy constructor
82
EigenVector_View::EigenVector_View
(
EigenVector_View
&v)
83
: m_view(v.vector_mutator().data(), v.size(),
dynamic_stride
(0, 1)), m_isConst(false) {}
84
86
vec_map_type
&
EigenVector_View::vector_mutator
() {
87
if
(!
m_isConst
) {
88
return
m_view
;
89
}
else
{
90
throw
std::runtime_error(
"Vector is const vector, cannot mutate const vector."
);
91
}
92
}
93
94
EigenVector_View
&
EigenVector_View::operator=
(
EigenVector_View
&v) {
95
new
(&
m_view
)
96
vec_map_type
(v.
m_view
.data(), v.
m_view
.size(),
dynamic_stride
(v.
m_view
.outerStride(), v.
m_view
.innerStride()));
97
98
return
*
this
;
99
}
100
101
EigenVector_View
&
EigenVector_View::operator=
(
EigenVector_View
&&v) {
102
new
(&
m_view
)
103
vec_map_type
(v.m_view.data(), v.m_view.size(),
dynamic_stride
(v.m_view.outerStride(), v.m_view.innerStride()));
104
return
*
this
;
105
}
106
}
// namespace Mantid::CurveFitting
EigenVectorView.h
Mantid::CurveFitting::EigenVector_View
Definition
EigenVectorView.h:17
Mantid::CurveFitting::EigenVector_View::EigenVector_View
EigenVector_View()
Definition
EigenVectorView.cpp:12
Mantid::CurveFitting::EigenVector_View::m_view
vec_map_type m_view
Definition
EigenVectorView.h:56
Mantid::CurveFitting::EigenVector_View::m_isConst
bool m_isConst
Definition
EigenVectorView.h:57
Mantid::CurveFitting::EigenVector_View::operator=
EigenVector_View & operator=(EigenVector_View &V)
Definition
EigenVectorView.cpp:94
Mantid::CurveFitting::EigenVector_View::vector_mutator
vec_map_type & vector_mutator()
Definition
EigenVectorView.cpp:86
Mantid::CurveFitting
Definition
IFunction1D.h:19
Mantid::CurveFitting::vec_map_type
Eigen::Map< Eigen::VectorXd, 0, dynamic_stride > vec_map_type
Definition
EigenVectorView.h:14
Mantid::CurveFitting::vec_const_map_type
Eigen::Map< const Eigen::VectorXd, 0, dynamic_stride > vec_const_map_type
Definition
EigenVectorView.h:15
Mantid::CurveFitting::dynamic_stride
Eigen::Stride< Eigen::Dynamic, Eigen::Dynamic > dynamic_stride
Definition
EigenComplexVector.h:22
Generated by
1.9.8