Mantid
Loading...
Searching...
No Matches
EigenMatrixView.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 +
8
9namespace Mantid::CurveFitting {
10// EigenMatrix_View Constructors
11// default constructor
12EigenMatrix_View::EigenMatrix_View() : m_view({}, 0, 0, dynamic_stride(0, 0)), m_isConst(false) {}
13
14// constructor: array->matrix view
22EigenMatrix_View::EigenMatrix_View(double *base, const size_t nTotalRows, size_t nTotalCols, size_t nElements_1,
23 size_t nElements_2, const size_t startElement_1, const size_t startElement_2)
24 : m_view(base, nTotalRows, nTotalCols, dynamic_stride(nTotalRows, 1)) {
25 if (nElements_1 == SIZE_T_NULL && nElements_2 == SIZE_T_NULL)
26 // if both are default, exit as m_view is initialised as such
27 return;
28 initialiseMatrix(nTotalRows, nTotalCols, nElements_1, nElements_2);
29 new (&m_view) map_type(base + (startElement_2 * nTotalRows) + startElement_1, nElements_1, nElements_2,
30 dynamic_stride(nTotalRows, 1));
31}
32
33// constructor: matrix->matrix view
39EigenMatrix_View::EigenMatrix_View(Eigen::MatrixXd &matrix, size_t nElements_1, size_t nElements_2,
40 const size_t startElement_1, const size_t startElement_2)
41 : m_view(matrix.data(), matrix.rows(), matrix.cols(), dynamic_stride(matrix.outerStride(), matrix.innerStride())) {
42 if (nElements_1 == SIZE_T_NULL && nElements_2 == SIZE_T_NULL)
43 // if both are default, exit as m_view is initialised as such
44 return;
45 initialiseMatrix(matrix.rows(), matrix.cols(), nElements_1, nElements_2);
46 new (&m_view) map_type(matrix.data() + (startElement_2 * matrix.rows()) + startElement_1, nElements_1, nElements_2,
47 dynamic_stride(matrix.outerStride(), matrix.innerStride()));
48}
49
50// constructor: map->matrix view
56EigenMatrix_View::EigenMatrix_View(map_type &matrix, size_t nElements_1, size_t nElements_2,
57 const size_t startElement_1, const size_t startElement_2)
58 : m_view(matrix.data(), matrix.rows(), matrix.cols(), dynamic_stride(matrix.outerStride(), matrix.innerStride())) {
59 if (nElements_1 == SIZE_T_NULL && nElements_2 == SIZE_T_NULL)
60 // if both are default, exit as m_view is initialised as such
61 return;
62 initialiseMatrix(matrix.rows(), matrix.cols(), nElements_1, nElements_2);
63 new (&m_view) map_type(matrix.data() + (startElement_2 * matrix.rows()) + startElement_1, nElements_1, nElements_2,
64 dynamic_stride(matrix.outerStride(), matrix.innerStride()));
65}
66
75EigenMatrix_View::EigenMatrix_View(const double *base, const size_t nTotalRows, size_t nTotalCols, size_t nElements_1,
76 size_t nElements_2, const size_t startElement_1, const size_t startElement_2)
77 : m_view({}, 0, 0, dynamic_stride(0, 0)), m_isConst(true) {
78 initialiseMatrix(nTotalRows, nTotalCols, nElements_1, nElements_2);
79 new (&m_view) const_map_type(base + (startElement_2 * nTotalRows) + startElement_1, nElements_1, nElements_2,
80 dynamic_stride(nTotalRows, 1));
81}
82
89EigenMatrix_View::EigenMatrix_View(const Eigen::MatrixXd &matrix, size_t nElements_1, size_t nElements_2,
90 const size_t startElement_1, const size_t startElement_2)
91 : m_view({}, 0, 0, dynamic_stride(0, 0)), m_isConst(true) {
92 initialiseMatrix(matrix.rows(), matrix.cols(), nElements_1, nElements_2);
93 new (&m_view) const_map_type(matrix.data() + (startElement_2 * matrix.rows()) + startElement_1, nElements_1,
94 nElements_2, dynamic_stride(matrix.rows(), 1));
95}
96
103EigenMatrix_View::EigenMatrix_View(const map_type &matrix, size_t nElements_1, size_t nElements_2,
104 const size_t startElement_1, const size_t startElement_2)
105 : m_view({}, 0, 0, dynamic_stride(0, 0)), m_isConst(true) {
106 initialiseMatrix(matrix.rows(), matrix.cols(), nElements_1, nElements_2);
107 new (&m_view) const_map_type(matrix.data() + (startElement_2 * matrix.rows()) + startElement_1, nElements_1,
108 nElements_2, dynamic_stride(matrix.outerStride(), matrix.innerStride()));
109}
110
114 : m_view({}, 0, 0, dynamic_stride(0, 0)), m_isConst(true) {
115 new (&m_view)
116 const_map_type(v.matrix_inspector().data(), v.rows(), v.cols(), dynamic_stride(v.outerStride(), v.innerStride()));
117}
118
122 : m_view(v.matrix_mutator().data(), v.rows(), v.cols(), dynamic_stride(v.outerStride(), v.innerStride())),
123 m_isConst(false) {}
124
127 if (!m_isConst) {
128 return m_view;
129 } else {
130 throw std::runtime_error("Matrix is const matrix, a const matrix cannot be mutated.");
131 }
132}
133
136 new (&m_view) map_type(V.m_view.data(), V.m_view.rows(), V.m_view.cols(),
137 dynamic_stride(V.m_view.outerStride(), V.m_view.innerStride()));
138 return *this;
139}
140
142 m_isConst = V.m_isConst;
143 new (&m_view) map_type(V.m_view.data(), V.m_view.rows(), V.m_view.cols(),
144 dynamic_stride(V.m_view.outerStride(), V.m_view.innerStride()));
145 return *this;
146}
147} // namespace Mantid::CurveFitting
#define SIZE_T_NULL
EigenMatrix_View & operator=(EigenMatrix_View &V)
void initialiseMatrix(const size_t nTotalRows, const size_t nTotalCols, size_t &nElements_1, size_t &nElements_2)
Eigen::Map< Eigen::MatrixXd, 0, dynamic_stride > map_type
Eigen::Stride< Eigen::Dynamic, Eigen::Dynamic > dynamic_stride
Eigen::Map< const Eigen::MatrixXd, 0, dynamic_stride > const_map_type