eigen_ref.cpp
Go to the documentation of this file.
1 /*
2  * Copyright 2014-2019, CNRS
3  * Copyright 2018-2022, INRIA
4  */
5 
6 #include <iostream>
7 
8 #include "eigenpy/eigenpy.hpp"
10 
11 using namespace Eigen;
12 using namespace eigenpy;
13 
14 template <typename MatType>
15 void printMatrix(const Eigen::Ref<const MatType> mat) {
16  if (MatType::IsVectorAtCompileTime) std::cout << "isVector" << std::endl;
17  std::cout << "input size: cols " << mat.cols() << " rows " << mat.rows()
18  << std::endl;
19  std::cout << mat << std::endl;
20 }
21 
22 template <typename VecType>
23 void printVector(const Eigen::Ref<const VecType>& vec) {
24  EIGEN_STATIC_ASSERT_VECTOR_ONLY(VecType);
26 }
27 
28 template <typename MatType>
29 void setOnes(Eigen::Ref<MatType> mat) {
30  mat.setOnes();
31 }
32 
33 template <typename VecType>
34 VecType copyVectorFromConstRef(const Eigen::Ref<const VecType> vec) {
35  std::cout << "copyVectorFromConstRef::vec: " << vec.transpose() << std::endl;
36  return VecType(vec);
37 }
38 
39 template <typename MatType>
40 Eigen::Ref<MatType> getBlock(Eigen::Ref<MatType> mat, Eigen::DenseIndex i,
41  Eigen::DenseIndex j, Eigen::DenseIndex n,
42  Eigen::DenseIndex m) {
43  return mat.block(i, j, n, m);
44 }
45 
46 template <typename MatType>
47 Eigen::Ref<MatType> editBlock(Eigen::Ref<MatType> mat, Eigen::DenseIndex i,
48  Eigen::DenseIndex j, Eigen::DenseIndex n,
49  Eigen::DenseIndex m) {
50  typename Eigen::Ref<MatType>::BlockXpr B = mat.block(i, j, n, m);
51  int k = 0;
52  for (int i = 0; i < B.rows(); ++i) {
53  for (int j = 0; j < B.cols(); ++j) {
54  B(i, j) = k++;
55  }
56  }
57  std::cout << "B:\n" << B << std::endl;
58  return mat;
59 }
60 
61 template <typename MatType>
62 void fill(Eigen::Ref<MatType> mat, const typename MatType::Scalar& value) {
63  mat.fill(value);
64 }
65 
67 template <typename MatType>
68 Eigen::Ref<MatType> getRefToStatic(const int rows, const int cols) {
69  static MatType mat(rows, cols);
70  std::cout << "create ref to matrix of size (" << rows << "," << cols << ")\n";
71  return mat;
72 }
73 
74 template <typename MatType>
75 Eigen::Ref<MatType> asRef(Eigen::Ref<MatType> mat) {
76  std::cout << "create Ref to input mutable Ref\n";
77  return Eigen::Ref<MatType>(mat);
78 }
79 
80 template <typename MatType>
81 const Eigen::Ref<const MatType> asConstRef(Eigen::Ref<const MatType> mat) {
82  return Eigen::Ref<const MatType>(mat);
83 }
84 
85 struct modify_block {
86  MatrixXd J;
87  modify_block() : J(10, 10) { J.setZero(); }
88  void modify(int n, int m) { call(J.topLeftCorner(n, m)); }
89  virtual void call(Eigen::Ref<MatrixXd> mat) = 0;
90 };
91 
92 struct modify_block_wrap : modify_block, bp::wrapper<modify_block> {
93  modify_block_wrap() : modify_block() {}
94  void call(Eigen::Ref<MatrixXd> mat) { this->get_override("call")(mat); }
95 };
96 
98  MatrixXd J;
99  Eigen::Ref<MatrixXd> Jref;
100  has_ref_member() : J(4, 4), Jref(J.topRightCorner(3, 3)) { J.setZero(); }
101 };
102 
104  namespace bp = boost::python;
106 
107  bp::def("printMatrix", printMatrix<Vector3d>);
108  bp::def("printMatrix", printMatrix<VectorXd>);
109  bp::def("printMatrix", printMatrix<MatrixXd>);
110 
111  bp::def("printVector", printVector<VectorXd>);
112  bp::def("printRowVector", printVector<RowVectorXd>);
113 
114  bp::def("setOnes", setOnes<Vector3d>);
115  bp::def("setOnes", setOnes<VectorXd>);
116  bp::def("setOnes", setOnes<MatrixXd>);
117 
118  bp::def("fillVec3", fill<Vector3d>);
119  bp::def("fillVec", fill<VectorXd>);
120  bp::def("fill", fill<MatrixXd>);
121 
122  bp::def("getRefToStatic", getRefToStatic<MatrixXd>);
123  bp::def("asRef", asRef<MatrixXd>);
124  bp::def("asConstRef", asConstRef<MatrixXd>);
125 
126  bp::def("getBlock", &getBlock<MatrixXd>);
127  bp::def("editBlock", &editBlock<MatrixXd>);
128 
129  bp::def("copyVectorFromConstRef", &copyVectorFromConstRef<VectorXd>);
130  bp::def("copyRowVectorFromConstRef", &copyVectorFromConstRef<RowVectorXd>);
131 
132  bp::class_<modify_block_wrap, boost::noncopyable>("modify_block",
133  bp::init<>())
134  .def_readonly("J", &modify_block::J)
135  .def("modify", &modify_block::modify)
136  .def("call", bp::pure_virtual(&modify_block_wrap::call));
137 
138  bp::class_<has_ref_member, boost::noncopyable>("has_ref_member", bp::init<>())
139  .def_readonly("J", &has_ref_member::J)
140  .add_property(
141  "Jref",
142  bp::make_getter(&has_ref_member::Jref,
143  bp::return_value_policy<bp::return_by_value>()));
144  // can't return Eigen::Ref by reference but by value
145  // (def_readonly creates a by-reference getter)
146 }
boost::python
Definition: alignment.hpp:49
Eigen
Definition: complex.cpp:7
fill
void fill(Eigen::Ref< MatType > mat, const typename MatType::Scalar &value)
Definition: eigen_ref.cpp:62
editBlock
Eigen::Ref< MatType > editBlock(Eigen::Ref< MatType > mat, Eigen::DenseIndex i, Eigen::DenseIndex j, Eigen::DenseIndex n, Eigen::DenseIndex m)
Definition: eigen_ref.cpp:47
modify_block::modify_block
modify_block()
Definition: eigen_ref.cpp:87
modify_block_wrap::call
void call(Eigen::Ref< MatrixXd > mat)
Definition: eigen_ref.cpp:94
eigenpy::enableEigenPy
void EIGENPY_DLLAPI enableEigenPy()
Definition: eigenpy.cpp:43
test_CholmodSimplicialLDLT.B
B
Definition: test_CholmodSimplicialLDLT.py:18
modify_block::modify
void modify(int n, int m)
Definition: eigen_ref.cpp:88
test_complex.rows
int rows
Definition: test_complex.py:4
eigen-from-python.hpp
test_matrix.vec
vec
Definition: test_matrix.py:180
printMatrix
void printMatrix(const Eigen::Ref< const MatType > mat)
Definition: eigen_ref.cpp:15
BOOST_PYTHON_MODULE
BOOST_PYTHON_MODULE(eigen_ref)
Definition: eigen_ref.cpp:103
has_ref_member::J
MatrixXd J
Definition: eigen_ref.cpp:98
eigenpy
Definition: alignment.hpp:14
modify_block_wrap
Definition: eigen_ref.cpp:92
has_ref_member::Jref
Eigen::Ref< MatrixXd > Jref
Definition: eigen_ref.cpp:99
getRefToStatic
Eigen::Ref< MatType > getRefToStatic(const int rows, const int cols)
Get ref to a static matrix of size ( rows, cols )
Definition: eigen_ref.cpp:68
getBlock
Eigen::Ref< MatType > getBlock(Eigen::Ref< MatType > mat, Eigen::DenseIndex i, Eigen::DenseIndex j, Eigen::DenseIndex n, Eigen::DenseIndex m)
Definition: eigen_ref.cpp:40
asConstRef
const Eigen::Ref< const MatType > asConstRef(Eigen::Ref< const MatType > mat)
Definition: eigen_ref.cpp:81
test_eigen_ref.mat
mat
Definition: test_eigen_ref.py:137
test_matrix.value
float value
Definition: test_matrix.py:161
has_ref_member::has_ref_member
has_ref_member()
Definition: eigen_ref.cpp:100
eigenpy.hpp
asRef
Eigen::Ref< MatType > asRef(Eigen::Ref< MatType > mat)
Definition: eigen_ref.cpp:75
printVector
void printVector(const Eigen::Ref< const VecType > &vec)
Definition: eigen_ref.cpp:23
test_sparse_matrix.m
m
Definition: test_sparse_matrix.py:5
copyVectorFromConstRef
VecType copyVectorFromConstRef(const Eigen::Ref< const VecType > vec)
Definition: eigen_ref.cpp:34
modify_block::J
MatrixXd J
Definition: eigen_ref.cpp:86
test_complex.cols
int cols
Definition: test_complex.py:5
setOnes
void setOnes(Eigen::Ref< MatType > mat)
Definition: eigen_ref.cpp:29
modify_block_wrap::modify_block_wrap
modify_block_wrap()
Definition: eigen_ref.cpp:93
eigenpy::call
Allows a template specialization.
Definition: expose.hpp:15


eigenpy
Author(s): Justin Carpentier, Nicolas Mansard
autogenerated on Sat Nov 2 2024 02:14:45