nomalloc.cpp
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
5 // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 // discard stack allocation as that too bypasses malloc
12 #define EIGEN_STACK_ALLOCATION_LIMIT 0
13 // heap allocation will raise an assert if enabled at runtime
14 #define EIGEN_RUNTIME_NO_MALLOC
15 
16 #include "main.h"
17 #include <Eigen/Cholesky>
18 #include <Eigen/Eigenvalues>
19 #include <Eigen/LU>
20 #include <Eigen/QR>
21 #include <Eigen/SVD>
22 
23 template<typename MatrixType> void nomalloc(const MatrixType& m)
24 {
25  /* this test check no dynamic memory allocation are issued with fixed-size matrices
26  */
27  typedef typename MatrixType::Scalar Scalar;
28 
29  Index rows = m.rows();
30  Index cols = m.cols();
31 
32  MatrixType m1 = MatrixType::Random(rows, cols),
33  m2 = MatrixType::Random(rows, cols),
34  m3(rows, cols);
35 
36  Scalar s1 = internal::random<Scalar>();
37 
38  Index r = internal::random<Index>(0, rows-1),
39  c = internal::random<Index>(0, cols-1);
40 
41  VERIFY_IS_APPROX((m1+m2)*s1, s1*m1+s1*m2);
42  VERIFY_IS_APPROX((m1+m2)(r,c), (m1(r,c))+(m2(r,c)));
43  VERIFY_IS_APPROX(m1.cwiseProduct(m1.block(0,0,rows,cols)), (m1.array()*m1.array()).matrix());
44  VERIFY_IS_APPROX((m1*m1.transpose())*m2, m1*(m1.transpose()*m2));
45 
46  m2.col(0).noalias() = m1 * m1.col(0);
47  m2.col(0).noalias() -= m1.adjoint() * m1.col(0);
48  m2.col(0).noalias() -= m1 * m1.row(0).adjoint();
49  m2.col(0).noalias() -= m1.adjoint() * m1.row(0).adjoint();
50 
51  m2.row(0).noalias() = m1.row(0) * m1;
52  m2.row(0).noalias() -= m1.row(0) * m1.adjoint();
53  m2.row(0).noalias() -= m1.col(0).adjoint() * m1;
54  m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint();
55  VERIFY_IS_APPROX(m2,m2);
56 
57  m2.col(0).noalias() = m1.template triangularView<Upper>() * m1.col(0);
58  m2.col(0).noalias() -= m1.adjoint().template triangularView<Upper>() * m1.col(0);
59  m2.col(0).noalias() -= m1.template triangularView<Upper>() * m1.row(0).adjoint();
60  m2.col(0).noalias() -= m1.adjoint().template triangularView<Upper>() * m1.row(0).adjoint();
61 
62  m2.row(0).noalias() = m1.row(0) * m1.template triangularView<Upper>();
63  m2.row(0).noalias() -= m1.row(0) * m1.adjoint().template triangularView<Upper>();
64  m2.row(0).noalias() -= m1.col(0).adjoint() * m1.template triangularView<Upper>();
65  m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint().template triangularView<Upper>();
66  VERIFY_IS_APPROX(m2,m2);
67 
68  m2.col(0).noalias() = m1.template selfadjointView<Upper>() * m1.col(0);
69  m2.col(0).noalias() -= m1.adjoint().template selfadjointView<Upper>() * m1.col(0);
70  m2.col(0).noalias() -= m1.template selfadjointView<Upper>() * m1.row(0).adjoint();
71  m2.col(0).noalias() -= m1.adjoint().template selfadjointView<Upper>() * m1.row(0).adjoint();
72 
73  m2.row(0).noalias() = m1.row(0) * m1.template selfadjointView<Upper>();
74  m2.row(0).noalias() -= m1.row(0) * m1.adjoint().template selfadjointView<Upper>();
75  m2.row(0).noalias() -= m1.col(0).adjoint() * m1.template selfadjointView<Upper>();
76  m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint().template selfadjointView<Upper>();
77  VERIFY_IS_APPROX(m2,m2);
78 
79  m2.template selfadjointView<Lower>().rankUpdate(m1.col(0),-1);
80  m2.template selfadjointView<Upper>().rankUpdate(m1.row(0),-1);
81  m2.template selfadjointView<Lower>().rankUpdate(m1.col(0), m1.col(0)); // rank-2
82 
83  // The following fancy matrix-matrix products are not safe yet regarding static allocation
84  m2.template selfadjointView<Lower>().rankUpdate(m1);
85  m2 += m2.template triangularView<Upper>() * m1;
86  m2.template triangularView<Upper>() = m2 * m2;
87  m1 += m1.template selfadjointView<Lower>() * m2;
88  VERIFY_IS_APPROX(m2,m2);
89 }
90 
91 template<typename Scalar>
93 {
94  const int maxSize = 16;
95  const int size = 12;
96 
97  typedef Eigen::Matrix<Scalar,
99  0,
100  maxSize, maxSize> Matrix;
101 
102  typedef Eigen::Matrix<Scalar,
103  Eigen::Dynamic, 1,
104  0,
105  maxSize, 1> Vector;
106 
109  0,
110  maxSize, maxSize> ComplexMatrix;
111 
112  const Matrix A(Matrix::Random(size, size)), B(Matrix::Random(size, size));
113  Matrix X(size,size);
114  const ComplexMatrix complexA(ComplexMatrix::Random(size, size));
115  const Matrix saA = A.adjoint() * A;
116  const Vector b(Vector::Random(size));
117  Vector x(size);
118 
119  // Cholesky module
121  X = LLT.solve(B);
122  x = LLT.solve(b);
124  X = LDLT.solve(B);
125  x = LDLT.solve(b);
126 
127  // Eigenvalues module
128  Eigen::HessenbergDecomposition<ComplexMatrix> hessDecomp; hessDecomp.compute(complexA);
129  Eigen::ComplexSchur<ComplexMatrix> cSchur(size); cSchur.compute(complexA);
130  Eigen::ComplexEigenSolver<ComplexMatrix> cEigSolver; cEigSolver.compute(complexA);
131  Eigen::EigenSolver<Matrix> eigSolver; eigSolver.compute(A);
132  Eigen::SelfAdjointEigenSolver<Matrix> saEigSolver(size); saEigSolver.compute(saA);
133  Eigen::Tridiagonalization<Matrix> tridiag; tridiag.compute(saA);
134 
135  // LU module
137  X = ppLU.solve(B);
138  x = ppLU.solve(b);
139  Eigen::FullPivLU<Matrix> fpLU; fpLU.compute(A);
140  X = fpLU.solve(B);
141  x = fpLU.solve(b);
142 
143  // QR module
145  X = hQR.solve(B);
146  x = hQR.solve(b);
148  X = cpQR.solve(B);
149  x = cpQR.solve(b);
151  // FIXME X = fpQR.solve(B);
152  x = fpQR.solve(b);
153 
154  // SVD module
156 }
157 
159  // default constructors:
160  Eigen::MatrixXd A;
161  Eigen::VectorXd v;
162  // explicit zero-sized:
163  Eigen::ArrayXXd A0(0,0);
164  Eigen::ArrayXd v0(0);
165 
166  // assigning empty objects to each other:
167  A=A0;
168  v=v0;
169 }
170 
171 template<typename MatrixType> void test_reference(const MatrixType& m) {
172  typedef typename MatrixType::Scalar Scalar;
173  enum { Flag = MatrixType::IsRowMajor ? Eigen::RowMajor : Eigen::ColMajor};
174  enum { TransposeFlag = !MatrixType::IsRowMajor ? Eigen::RowMajor : Eigen::ColMajor};
175  Index rows = m.rows(), cols=m.cols();
178  // Dynamic reference:
180  typedef Eigen::Ref<const MatrixXT > RefT;
181 
182  Ref r1(m);
183  Ref r2(m.block(rows/3, cols/4, rows/2, cols/2));
184  RefT r3(m.transpose());
185  RefT r4(m.topLeftCorner(rows/2, cols/2).transpose());
186 
187  VERIFY_RAISES_ASSERT(RefT r5(m));
188  VERIFY_RAISES_ASSERT(Ref r6(m.transpose()));
189  VERIFY_RAISES_ASSERT(Ref r7(Scalar(2) * m));
190 
191  // Copy constructors shall also never malloc
192  Ref r8 = r1;
193  RefT r9 = r3;
194 
195  // Initializing from a compatible Ref shall also never malloc
197 
198  // Initializing from an incompatible Ref will malloc:
199  typedef Eigen::Ref<const MatrixX, Aligned> RefAligned;
200  VERIFY_RAISES_ASSERT(RefAligned r12=r10);
201  VERIFY_RAISES_ASSERT(Ref r13=r10); // r10 has more dynamic strides
202 
203 }
204 
206 {
207  // create some dynamic objects
208  Eigen::MatrixXd M1 = MatrixXd::Random(3,3);
209  Ref<const MatrixXd> R1 = 2.0*M1; // Ref requires temporary
210 
211  // from here on prohibit malloc:
212  Eigen::internal::set_is_malloc_allowed(false);
213 
214  // check that our operator new is indeed called:
215  VERIFY_RAISES_ASSERT(MatrixXd dummy(MatrixXd::Random(3,3)));
217  CALL_SUBTEST_2(nomalloc(Matrix4d()) );
219 
220  // Check decomposition modules with dynamic matrices that have a known compile-time max size (ctms)
221  CALL_SUBTEST_4(ctms_decompositions<float>());
222 
224 
227  CALL_SUBTEST_8(Ref<MatrixXd> R2 = M1.topRows<2>(); test_reference(R2));
228 }
Matrix< SCALARB, Dynamic, Dynamic, opt_B > B
Definition: bench_gemm.cpp:49
HouseholderQR & compute(const EigenBase< InputType > &matrix)
Matrix3f m
Robust Cholesky decomposition of a matrix with pivoting.
Definition: LDLT.h:59
const Solve< LLT< _MatrixType, _UpLo >, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: SolverBase.h:106
SCALAR Scalar
Definition: bench_gemm.cpp:46
#define VERIFY_RAISES_ASSERT(a)
Definition: main.h:340
#define CALL_SUBTEST_6(FUNC)
#define CALL_SUBTEST_4(FUNC)
EigenSolver & compute(const EigenBase< InputType > &matrix, bool computeEigenvectors=true)
Computes eigendecomposition of given matrix.
PartialPivLU & compute(const EigenBase< InputType > &matrix)
Definition: PartialPivLU.h:131
Householder rank-revealing QR decomposition of a matrix with full pivoting.
Scalar * b
Definition: benchVecAdd.cpp:17
LDLT & compute(const EigenBase< InputType > &matrix)
#define CALL_SUBTEST_3(FUNC)
ComplexEigenSolver & compute(const EigenBase< InputType > &matrix, bool computeEigenvectors=true)
Computes eigendecomposition of given matrix.
MatrixType m2(n_dims)
EIGEN_DEVICE_FUNC SelfAdjointEigenSolver & compute(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix.
#define CALL_SUBTEST_7(FUNC)
FullPivLU & compute(const EigenBase< InputType > &matrix)
Definition: FullPivLU.h:120
Scalar Scalar * c
Definition: benchVecAdd.cpp:17
Computes eigenvalues and eigenvectors of selfadjoint matrices.
LU decomposition of a matrix with partial pivoting, and related features.
MatrixXf M1
MatrixXf MatrixType
Tridiagonalization & compute(const EigenBase< InputType > &matrix)
Computes tridiagonal decomposition of given matrix.
void nomalloc(const MatrixType &m)
Definition: nomalloc.cpp:23
Tridiagonal decomposition of a selfadjoint matrix.
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:48
ComplexSchur & compute(const EigenBase< InputType > &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.
EIGEN_DECLARE_TEST(nomalloc)
Definition: nomalloc.cpp:205
FullPivHouseholderQR & compute(const EigenBase< InputType > &matrix)
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
void test_reference(const MatrixType &m)
Definition: nomalloc.cpp:171
#define VERIFY_IS_APPROX(a, b)
ColPivHouseholderQR & compute(const EigenBase< InputType > &matrix)
HessenbergDecomposition & compute(const EigenBase< InputType > &matrix)
Computes Hessenberg decomposition of given matrix.
#define CALL_SUBTEST_1(FUNC)
Matrix3d m1
Definition: IOFormat.cpp:2
void test_zerosized()
Definition: nomalloc.cpp:158
Householder rank-revealing QR decomposition of a matrix with column-pivoting.
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition: LLT.h:66
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
Array< int, Dynamic, 1 > v
#define CALL_SUBTEST_8(FUNC)
static const double r2
static const double r3
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:281
#define CALL_SUBTEST_5(FUNC)
static const double r1
Reduces a square matrix to Hessenberg form by an orthogonal similarity transformation.
LU decomposition of a matrix with complete pivoting, and related features.
static const double v0
Householder QR decomposition of a matrix.
Two-sided Jacobi SVD decomposition of a rectangular matrix.
JacobiSVD & compute(const MatrixType &matrix, unsigned int computationOptions)
Method performing the decomposition of given matrix using custom options.
Definition: JacobiSVD.h:666
LLT & compute(const EigenBase< InputType > &matrix)
#define CALL_SUBTEST_2(FUNC)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > Matrix
const int Dynamic
Definition: Constants.h:22
Computes eigenvalues and eigenvectors of general matrices.
Definition: EigenSolver.h:64
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
The matrix class, also used for vectors and row-vectors.
Performs a complex Schur decomposition of a real or complex square matrix.
Definition: ComplexSchur.h:51
#define X
Definition: icosphere.cpp:20
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Computes eigenvalues and eigenvectors of general complex matrices.
Eigen::Matrix< double, Eigen::Dynamic, 1 > Vector
void ctms_decompositions()
Definition: nomalloc.cpp:92


gtsam
Author(s):
autogenerated on Tue Jul 4 2023 02:34:57