iterative-solver 0.0
Matrix.h
1#ifndef LINEARALGEBRA_SRC_MOLPRO_LINALG_ITSOLV_SUBSPACE_MATRIX_H
2#define LINEARALGEBRA_SRC_MOLPRO_LINALG_ITSOLV_SUBSPACE_MATRIX_H
3#include <algorithm>
4#include <cassert>
5#include <cstddef>
6#include <format>
7#include <iomanip>
8#include <sstream>
9#include <stdexcept>
10#include <string_view>
11#include <utility>
12#include <vector>
13
15
16// FIXME Since we are only supporting a few underlying data types Matrix should not be header only
17// FIXME Usage of Eigen for matrix would not be as catastrophic now that the headers are separated
29template <typename T>
30class Matrix {
31protected:
32 class Slice;
33 class CSlice;
34
35public:
36 using value_type = T;
37 using index_type = size_t;
38 using coord_type = std::pair<size_t, size_t>;
39
40public:
41 explicit Matrix(coord_type dims) : m_rows(dims.first), m_cols(dims.second), m_buffer(m_rows * m_cols) {}
43 explicit Matrix(std::vector<T>&& data, coord_type dims)
44 : m_rows(dims.first), m_cols(dims.second), m_buffer(std::move(data)) {
45 if (m_buffer.size() != size())
46 throw std::runtime_error("data buffer is of the wrong size");
47 }
48 explicit Matrix(const std::vector<T>& data, coord_type dims)
49 : m_rows(dims.first), m_cols(dims.second), m_buffer(data) {
50 if (m_buffer.size() != size())
51 throw std::runtime_error("data buffer is of the wrong size");
52 }
53 Matrix() = default;
54 ~Matrix() = default;
55 Matrix(const Matrix<T>&) = default;
56 Matrix(Matrix<T>&&) noexcept = default;
57 Matrix<T>& operator=(const Matrix<T>&) = default;
58 Matrix<T>& operator=(Matrix<T>&&) noexcept = default;
59
61 T& operator()(index_type i, index_type j) { return m_buffer[i * m_cols + j]; }
62
64 T operator()(index_type i, index_type j) const { return m_buffer[i * m_cols + j]; }
65
67 const std::vector<T>& data() const& { return m_buffer; }
69 std::vector<T>&& data() && {
70 m_rows = 0;
71 m_cols = 0;
72 return std::move(m_buffer);
73 }
74
76 bool empty() const { return size() == 0; }
77
79 void clear() { resize({0, 0}); }
80
82 coord_type to_coord(size_t ind) const {
83 if (ind >= size())
84 throw std::out_of_range("index is larger than size");
85 auto row = ind / m_cols;
86 auto col = ind - row * m_cols;
87 return {row, col};
88 }
89
91 void fill(T value) { std::fill(begin(m_buffer), end(m_buffer), value); }
92
99 Slice slice(coord_type upper_left, coord_type bottom_right) {
100 return Slice(*this, std::move(upper_left), std::move(bottom_right));
101 }
102
104 Slice slice() { return slice({0, 0}, dimensions()); }
105
112 CSlice slice(coord_type upper_left, coord_type bottom_right) const {
113 return CSlice(*this, std::move(upper_left), std::move(bottom_right));
114 }
115
117 CSlice slice() const { return slice({0, 0}, dimensions()); }
118
120 Slice row(size_t i) { return slice({i, 0}, {i + 1, m_cols}); }
121 CSlice row(size_t i) const { return slice({i, 0}, {i + 1, m_cols}); }
122
124 Slice col(size_t j) { return slice({0, j}, {m_rows, j + 1}); }
125 CSlice col(size_t j) const { return slice({0, j}, {m_rows, j + 1}); }
126
128 void resize(const coord_type& dims) {
129 if (dims == dimensions())
130 return;
131 if (dims.second == m_cols) {
132 m_rows = dims.first;
133 m_buffer.resize(size());
134 } else {
135 auto m = Matrix<T>(dims);
136 auto upper_left = coord_type{0, 0};
137 auto bottom_right = coord_type{std::min(rows(), m.rows()), std::min(cols(), m.cols())};
138 m.slice(upper_left, bottom_right) = slice(upper_left, bottom_right);
139 std::swap(*this, m);
140 }
141 }
142
145 if (row >= m_rows)
146 throw std::runtime_error("row is out of range");
147 slice({row, 0}, {m_rows - 1, m_cols}) = slice({row + 1, 0}, dimensions());
148 resize({m_rows - 1, m_cols});
149 }
150
153 if (col >= m_cols)
154 throw std::runtime_error("column is out of range");
155 auto m = Matrix<T>({m_rows, m_cols - 1});
156 m.slice({0, 0}, {m_rows, col}) = slice({0, 0}, {m_rows, col});
157 m.slice({0, col}, {m_rows, m.m_cols}) = slice({0, col + 1}, dimensions());
158 std::swap(*this, m);
159 }
160
165 }
166
167 index_type rows() const { return m_rows; }
168 index_type cols() const { return m_cols; }
170 size_t size() const { return m_rows * m_cols; }
171
172protected:
175 std::vector<T> m_buffer;
176
180 class Slice {
181 public:
182 Slice(Matrix<T>& matrix, coord_type upper_left, coord_type bottom_right)
183 : mat(matrix), upl(std::move(upper_left)), btr(std::move(bottom_right)) {
184 if (upl.first > btr.first || upl.second > btr.second)
185 throw std::runtime_error("specified incorrect corners");
186 if (upl.first < 0 || upl.second < 0)
187 throw std::runtime_error("slice is out of range");
188 if (btr.first > matrix.rows() || btr.second > matrix.cols())
189 throw std::runtime_error("slice is out of range");
190 }
191 Slice() = delete;
192 ~Slice() = default;
193 Slice(const Slice&) = delete;
194 Slice(Slice&&) noexcept = default;
195 Slice& operator=(Slice&&) noexcept = default;
196
197 Slice& operator=(const Slice& right) {
198 if (dimensions() != right.dimensions())
199 throw std::runtime_error("attempting to copy slices of different dimensions");
200 for (size_t i = 0; i < dimensions().first; ++i) {
201 for (size_t j = 0; j < dimensions().second; ++j) {
202 mat(upl.first + i, upl.second + j) = right.mat(right.upl.first + i, right.upl.second + j);
203 }
204 }
205 return *this;
206 }
207
208 T& operator()(size_t i, size_t j) { return mat(upl.first + i, upl.second + j); }
209 T operator()(size_t i, size_t j) const { return mat(upl.first + i, upl.second + j); }
210
211 Slice& axpy(T a, const Slice& x) {
212 if (dimensions() != x.dimensions())
213 throw std::runtime_error("attempting to copy slices of different dimensions");
214 for (size_t i = 0; i < dimensions().first; ++i) {
215 for (size_t j = 0; j < dimensions().second; ++j) {
216 mat(upl.first + i, upl.second + j) += a * x.mat(x.upl.first + i, x.upl.second + j);
217 }
218 }
219 return *this;
220 }
221
222 Slice& axpy(T a, const CSlice& x) {
223 axpy(a, x.m_slice);
224 return *this;
225 }
226
228 Slice& scal(T a) {
229 for (size_t i = 0; i < dimensions().first; ++i)
230 for (size_t j = 0; j < dimensions().second; ++j)
231 mat(upl.first + i, upl.second + j) *= a;
232 return *this;
233 }
234
236 Slice& fill(T a) {
237 for (size_t i = 0; i < dimensions().first; ++i)
238 for (size_t j = 0; j < dimensions().second; ++j)
239 mat(upl.first + i, upl.second + j) = a;
240 return *this;
241 }
242
243 Slice& operator=(const CSlice& right) {
244 *this = right.m_slice;
245 return *this;
246 }
247
248 Slice& operator=(const Matrix<T>& right) {
249 *this = right.slice({0, 0}, right.dimensions());
250 return *this;
251 }
252
253 coord_type dimensions() const { return {btr.first - upl.first, btr.second - upl.second}; }
254 size_t rows() const { return btr.first - upl.first; }
255 size_t cols() const { return btr.second - upl.second; }
256
257 protected:
261 };
262
264 class CSlice {
265 friend Slice;
266
267 public:
268 CSlice(const Matrix<T>& matrix, coord_type upper_left, coord_type bottom_right)
269 : m_slice(const_cast<Matrix<T>&>(matrix), std::move(upper_left), std::move(bottom_right)) {}
270 CSlice() = delete;
271 ~CSlice() = default;
272 CSlice(const CSlice&) = delete;
273 CSlice(CSlice&&) noexcept = default;
274 CSlice& operator=(CSlice&&) = default;
275 T operator()(size_t i, size_t j) const { return const_cast<Slice&>(m_slice)(i, j); }
276
277 protected:
279 };
280};
281
282template <class ML, class MR>
283void transpose_copy(ML&& ml, const MR& mr) {
284 assert(ml.rows() == mr.cols() && ml.cols() == mr.rows());
285 for (size_t i = 0; i < ml.rows(); ++i)
286 for (size_t j = 0; j < ml.cols(); ++j)
287 ml(i, j) = mr(j, i);
288}
289
290template <class Mat>
291std::string as_string(const Mat& m, int precision = 6) {
292 auto s = std::stringstream{};
293 s << std::setprecision(precision);
294 auto dims = m.dimensions();
295 if (dims.first * dims.second != 0) {
296 s << "\n[";
297 for (size_t i = 0; i < dims.first; ++i) {
298 s << "[";
299 for (size_t j = 0; j < dims.second; ++j) {
300 s << m(i, j);
301 if (j != dims.second - 1)
302 s << ", ";
303 }
304 s << "]";
305 if (i != dims.first - 1)
306 s << ",\n ";
307 }
308 s << "]";
309 } else {
310 s << "[[]]";
311 }
312 return s.str();
313}
314
315} // namespace molpro::linalg::itsolv::subspace
316
317template<typename T>
318struct std::formatter<molpro::linalg::itsolv::subspace::Matrix<T>> : std::formatter<std::string_view> {
319 using Base = std::formatter<std::string_view>;
320
321 std::size_t precision = 6;
322
323 constexpr auto parse(std::format_parse_context &ctx) {
324 auto it = ctx.begin();
325
326 if (it == ctx.end()) {
327 return it;
328 }
329
330 if (*it != '.') {
331 return Base::parse(ctx);
332 }
333
334 ++it;
335
336 if (it == ctx.end() || ( *it < '0' && *it > '9' )) {
337 throw std::format_error("Format string uses precision quantifier '.' without giving a (numeric) precision");
338 }
339
340 // Parse specified precision
341 // Unfortunately, the stdlib doesn't have a constexpr function for this until C++23
342 precision = 0;
343 while (it != ctx.end() && ( *it >= '0' && *it <= '9' )) {
344 precision *= 10;
345 precision += *it - '0';
346 ++it;
347 }
348
349 // Interpret remaining format context as specification for std::string_view formatter
350 ctx.advance_to(it);
351 return Base::parse(ctx);
352 }
353
354 auto format(const molpro::linalg::itsolv::subspace::Matrix<T> &mat, std::format_context &ctx) const {
355 std::string tmp = molpro::linalg::itsolv::subspace::as_string(mat, precision);
356
357 return Base::format(tmp, ctx);
358 }
359};
360
361#endif // LINEARALGEBRA_SRC_MOLPRO_LINALG_ITSOLV_SUBSPACE_MATRIX_H
Constant slice that cannot be assigned to.
Definition: Matrix.h:264
CSlice(const Matrix< T > &matrix, coord_type upper_left, coord_type bottom_right)
Definition: Matrix.h:268
Proxy mapping to a rectangular slice of the matrix data. Implements simple assignment.
Definition: Matrix.h:180
coord_type upl
upper left corner
Definition: Matrix.h:259
Slice & operator=(const CSlice &right)
Definition: Matrix.h:243
coord_type dimensions() const
Definition: Matrix.h:253
Matrix< T > & mat
matrix being sliced
Definition: Matrix.h:258
Slice & fill(T a)
Fill all elements of the slice with new values.
Definition: Matrix.h:236
Slice(Matrix< T > &matrix, coord_type upper_left, coord_type bottom_right)
Definition: Matrix.h:182
size_t rows() const
Definition: Matrix.h:254
T operator()(size_t i, size_t j) const
Definition: Matrix.h:209
T & operator()(size_t i, size_t j)
Definition: Matrix.h:208
coord_type btr
bottom right corner
Definition: Matrix.h:260
Slice & operator=(const Matrix< T > &right)
Definition: Matrix.h:248
size_t cols() const
Definition: Matrix.h:255
Slice & axpy(T a, const CSlice &x)
Definition: Matrix.h:222
Slice & scal(T a)
Scale all elements of the slice.
Definition: Matrix.h:228
Slice & axpy(T a, const Slice &x)
Definition: Matrix.h:211
Matrix container that allows simple data access, slicing, copying and resizing without loosing data.
Definition: Matrix.h:30
Matrix(const std::vector< T > &data, coord_type dims)
Definition: Matrix.h:48
index_type m_rows
number of rows
Definition: Matrix.h:173
Slice slice()
Access the whole matrix as a slice.
Definition: Matrix.h:104
T operator()(index_type i, index_type j) const
Access element.
Definition: Matrix.h:64
Slice row(size_t i)
Access row slice.
Definition: Matrix.h:120
T value_type
Definition: Matrix.h:36
Slice col(size_t j)
Access column slice.
Definition: Matrix.h:124
std::pair< size_t, size_t > coord_type
Definition: Matrix.h:38
void remove_row(index_type row)
removes a row from the matrix
Definition: Matrix.h:144
Slice slice(coord_type upper_left, coord_type bottom_right)
Access a rectangular slice of the matrix.
Definition: Matrix.h:99
std::vector< T > m_buffer
data buffer
Definition: Matrix.h:175
void remove_row_col(index_type row, index_type col)
removes row and column
Definition: Matrix.h:162
CSlice row(size_t i) const
Definition: Matrix.h:121
Matrix(const Matrix< T > &)=default
CSlice col(size_t j) const
Definition: Matrix.h:125
CSlice slice(coord_type upper_left, coord_type bottom_right) const
Access a constant rectangular slice of the matrix.
Definition: Matrix.h:112
const std::vector< T > & data() const &
Access the underlying data buffer.
Definition: Matrix.h:67
void clear()
Clears all elements and sets dimensions to 0.
Definition: Matrix.h:79
void remove_col(index_type col)
removes a column from the matrix
Definition: Matrix.h:152
coord_type dimensions() const
Definition: Matrix.h:169
Matrix(coord_type dims)
Definition: Matrix.h:41
size_t size() const
Definition: Matrix.h:170
bool empty() const
Returns true if matrix is empty.
Definition: Matrix.h:76
Matrix(std::vector< T > &&data, coord_type dims)
Construct a matrix by taking ownership of an existing data buffer which must be of correct size.
Definition: Matrix.h:43
void fill(T value)
Sets all elements of matrix to value.
Definition: Matrix.h:91
Matrix(Matrix< T > &&) noexcept=default
void resize(const coord_type &dims)
Resize the matrix. The old data is preserved and any new rows/cols are zeroed.
Definition: Matrix.h:128
CSlice slice() const
Access the whole matrix as a slice.
Definition: Matrix.h:117
index_type rows() const
Definition: Matrix.h:167
coord_type to_coord(size_t ind) const
Converts index of 1D data to matrix coordinate.
Definition: Matrix.h:82
size_t index_type
Definition: Matrix.h:37
index_type cols() const
Definition: Matrix.h:168
std::vector< T > && data() &&
Access the raw data buffer of an r-value matrix. The matrix is left empty.
Definition: Matrix.h:69
index_type m_cols
number of columns
Definition: Matrix.h:174
Definition: PSpace.h:7
std::string as_string(const Mat &m, int precision=6)
Definition: Matrix.h:291
void transpose_copy(ML &&ml, const MR &mr)
Definition: Matrix.h:283
auto format(const molpro::linalg::itsolv::subspace::Matrix< T > &mat, std::format_context &ctx) const
Definition: Matrix.h:354
std::formatter< std::string_view > Base
Definition: Matrix.h:319
constexpr auto parse(std::format_parse_context &ctx)
Definition: Matrix.h:323