preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
EigenHelperFunctions.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Eigen/Core>
4#include <iterator>
5#include <vector>
6#include "utils/algorithm.hpp"
7#include "utils/assertion.hpp"
8
9namespace precice::utils {
10
11void shiftSetFirst(Eigen::MatrixXd &A, const Eigen::VectorXd &v);
12
13void appendFront(Eigen::MatrixXd &A, Eigen::VectorXd &v);
14
15void removeColumnFromMatrix(Eigen::MatrixXd &A, int col);
16
17void append(Eigen::VectorXd &v, double value);
18
19template <typename Derived1>
20void append(
21 Eigen::MatrixXd &A,
22 const Eigen::PlainObjectBase<Derived1> &B)
23{
24 int n = A.rows(), m = A.cols();
25 if (n <= 0 && m <= 0) {
26 A = B;
27 } else {
28 PRECICE_ASSERT(B.rows() == n, B.rows(), A.rows());
29 A.conservativeResize(n, m + B.cols());
30 for (auto i = m; i < m + B.cols(); i++)
31 A.col(i) = B.col(i - m);
32 }
33}
34
35template <typename Derived1>
36void append(
37 Eigen::VectorXd &v,
38 const Eigen::PlainObjectBase<Derived1> &app)
39{
40 int n = v.size();
41 if (n <= 0) {
42 v = app;
43 } else {
44 PRECICE_ASSERT(app.cols() == 1, app.cols());
45 v.conservativeResize(n + app.size());
46 for (int i = 0; i < app.size(); i++)
47 v(n + i) = app(i);
48 }
49}
50
56template <typename Derived>
57auto firstN(const Eigen::PlainObjectBase<Derived> &val, unsigned n) -> const Eigen::Map<const Eigen::Matrix<typename Derived::Scalar, 1, Eigen::Dynamic>>
58{
59 return {val.data(), std::min<Eigen::Index>(n, val.size())};
60}
61
62template <typename Derived, typename Iter = const typename Derived::Scalar *, typename Size = typename std::iterator_traits<Iter>::difference_type>
63const RangePreview<Iter> previewRange(Size n, const Eigen::PlainObjectBase<Derived> &eigen)
64{
65 auto begin = eigen.data();
66 auto end = begin;
67 std::advance(end, eigen.size());
68 return {n, begin, end};
69}
70
71template <typename DerivedLHS, typename DerivedRHS>
72bool componentWiseLess(const Eigen::PlainObjectBase<DerivedLHS> &lhs, const Eigen::PlainObjectBase<DerivedRHS> &rhs)
73{
74 if (lhs.size() != rhs.size())
75 return false;
76
77 for (int i = 0; i < rhs.size(); ++i) {
78 if (lhs[i] != rhs[i]) {
79 return lhs[i] < rhs[i];
80 }
81 }
82 return false;
83}
84
86 template <typename DerivedLHS, typename DerivedRHS>
87 bool operator()(const Eigen::PlainObjectBase<DerivedLHS> &lhs, const Eigen::PlainObjectBase<DerivedRHS> &rhs) const
88 {
89 return componentWiseLess(lhs, rhs);
90 }
91};
92
93} // namespace precice::utils
T advance(T... args)
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
T make_unique(T... args)
contains precice-related utilities.
void removeColumnFromMatrix(Eigen::MatrixXd &A, int col)
void appendFront(Eigen::MatrixXd &A, Eigen::VectorXd &v)
void append(Eigen::VectorXd &v, double value)
const RangePreview< Iter > previewRange(Size n, const Range &range)
auto firstN(const Eigen::PlainObjectBase< Derived > &val, unsigned n) -> const Eigen::Map< const Eigen::Matrix< typename Derived::Scalar, 1, Eigen::Dynamic > >
void shiftSetFirst(Eigen::MatrixXd &A, const Eigen::VectorXd &v)
bool componentWiseLess(const Eigen::PlainObjectBase< DerivedLHS > &lhs, const Eigen::PlainObjectBase< DerivedRHS > &rhs)
bool operator()(const Eigen::PlainObjectBase< DerivedLHS > &lhs, const Eigen::PlainObjectBase< DerivedRHS > &rhs) const
The RangePreview object used as a lazy proxy struct for proviewing the content of a Range.