aboutsummaryrefslogtreecommitdiffhomepage
path: root/test
diff options
context:
space:
mode:
authorGravatar Adolfo Rodriguez Tsouroukdissian <adolfo.rodriguez@pal-robotics.com>2010-04-21 17:15:57 +0200
committerGravatar Adolfo Rodriguez Tsouroukdissian <adolfo.rodriguez@pal-robotics.com>2010-04-21 17:15:57 +0200
commit28dde19e40a3d758faa94f0fc228857f7b3192ea (patch)
tree2647d94ffd250e7b215e98baddcecb8fb651543a /test
parentfaf8f7732dffd33eeae8afc8aad165668c8b6b2e (diff)
- Added problem size constructor to decompositions that did not have one. It preallocates member data structures.
- Updated unit tests to check above constructor. - In the compute() method of decompositions: Made temporary matrices/vectors class members to avoid heap allocations during compute() (when dynamic matrices are used, of course). These changes can speed up decomposition computation time when a solver instance is used to solve multiple same-sized problems. An added benefit is that the compute() method can now be invoked in contexts were heap allocations are forbidden, such as in real-time control loops. CAVEAT: Not all of the decompositions in the Eigenvalues module have a heap-allocation-free compute() method. A future patch may address this issue, but some required API changes need to be incorporated first.
Diffstat (limited to 'test')
-rw-r--r--test/cholesky.cpp4
-rw-r--r--test/eigensolver_complex.cpp3
-rw-r--r--test/eigensolver_generic.cpp3
-rw-r--r--test/eigensolver_selfadjoint.cpp4
-rw-r--r--test/hessenberg.cpp3
-rw-r--r--test/jacobisvd.cpp3
-rw-r--r--test/lu.cpp4
-rw-r--r--test/qr.cpp3
-rw-r--r--test/qr_colpivoting.cpp3
-rw-r--r--test/qr_fullpivoting.cpp3
-rw-r--r--test/schur_complex.cpp3
-rw-r--r--test/schur_real.cpp3
-rw-r--r--test/svd.cpp3
13 files changed, 42 insertions, 0 deletions
diff --git a/test/cholesky.cpp b/test/cholesky.cpp
index a446f5d73..0ae26c7d5 100644
--- a/test/cholesky.cpp
+++ b/test/cholesky.cpp
@@ -163,4 +163,8 @@ void test_cholesky()
CALL_SUBTEST_7( cholesky_verify_assert<Matrix3d>() );
CALL_SUBTEST_8( cholesky_verify_assert<MatrixXf>() );
CALL_SUBTEST_2( cholesky_verify_assert<MatrixXd>() );
+
+ // Test problem size constructors
+ CALL_SUBTEST_9( LLT<MatrixXf>(10) );
+ CALL_SUBTEST_9( LDLT<MatrixXf>(10) );
}
diff --git a/test/eigensolver_complex.cpp b/test/eigensolver_complex.cpp
index 4e973c877..b3d9ac24b 100644
--- a/test/eigensolver_complex.cpp
+++ b/test/eigensolver_complex.cpp
@@ -63,4 +63,7 @@ void test_eigensolver_complex()
CALL_SUBTEST_3( eigensolver(Matrix<std::complex<float>, 1, 1>()) );
CALL_SUBTEST_4( eigensolver(Matrix3f()) );
}
+
+ // Test problem size constructors
+ CALL_SUBTEST_5(ComplexEigenSolver<MatrixXf>(10));
}
diff --git a/test/eigensolver_generic.cpp b/test/eigensolver_generic.cpp
index 5dc18f141..f24c3b4ed 100644
--- a/test/eigensolver_generic.cpp
+++ b/test/eigensolver_generic.cpp
@@ -89,4 +89,7 @@ void test_eigensolver_generic()
CALL_SUBTEST_2( eigensolver_verify_assert<MatrixXd>() );
CALL_SUBTEST_4( eigensolver_verify_assert<Matrix2d>() );
CALL_SUBTEST_5( eigensolver_verify_assert<MatrixXf>() );
+
+ // Test problem size constructors
+ CALL_SUBTEST_6(EigenSolver<MatrixXf>(10));
}
diff --git a/test/eigensolver_selfadjoint.cpp b/test/eigensolver_selfadjoint.cpp
index f2146cb88..70b3e6791 100644
--- a/test/eigensolver_selfadjoint.cpp
+++ b/test/eigensolver_selfadjoint.cpp
@@ -129,5 +129,9 @@ void test_eigensolver_selfadjoint()
CALL_SUBTEST_6( selfadjointeigensolver(Matrix<double,1,1>()) );
CALL_SUBTEST_7( selfadjointeigensolver(Matrix<double,2,2>()) );
}
+
+ // Test problem size constructors
+ CALL_SUBTEST_8(SelfAdjointEigenSolver<MatrixXf>(10));
+ CALL_SUBTEST_8(Tridiagonalization<MatrixXf>(10));
}
diff --git a/test/hessenberg.cpp b/test/hessenberg.cpp
index cba9c8fda..ec1148bfc 100644
--- a/test/hessenberg.cpp
+++ b/test/hessenberg.cpp
@@ -61,4 +61,7 @@ void test_hessenberg()
CALL_SUBTEST_3(( hessenberg<std::complex<float>,4>() ));
CALL_SUBTEST_4(( hessenberg<float,Dynamic>(ei_random<int>(1,320)) ));
CALL_SUBTEST_5(( hessenberg<std::complex<double>,Dynamic>(ei_random<int>(1,320)) ));
+
+ // Test problem size constructors
+ CALL_SUBTEST_6(HessenbergDecomposition<MatrixXf>(10));
}
diff --git a/test/jacobisvd.cpp b/test/jacobisvd.cpp
index c8a1495d2..bc9d93754 100644
--- a/test/jacobisvd.cpp
+++ b/test/jacobisvd.cpp
@@ -106,4 +106,7 @@ void test_jacobisvd()
CALL_SUBTEST_3(( svd_verify_assert<Matrix3d>() ));
CALL_SUBTEST_9(( svd_verify_assert<MatrixXf>() ));
CALL_SUBTEST_11(( svd_verify_assert<MatrixXd>() ));
+
+ // Test problem size constructors
+ CALL_SUBTEST_12( JacobiSVD<MatrixXf>(10, 20) );
}
diff --git a/test/lu.cpp b/test/lu.cpp
index 37e2990d2..22dca76d2 100644
--- a/test/lu.cpp
+++ b/test/lu.cpp
@@ -212,5 +212,9 @@ void test_lu()
CALL_SUBTEST_6( lu_verify_assert<MatrixXcd>() );
CALL_SUBTEST_7(( lu_non_invertible<Matrix<float,Dynamic,16> >() ));
+
+ // Test problem size constructors
+ CALL_SUBTEST_9( PartialPivLU<MatrixXf>(10) );
+ CALL_SUBTEST_9( FullPivLU<MatrixXf>(10, 20); );
}
}
diff --git a/test/qr.cpp b/test/qr.cpp
index 1ce1f46e6..90209639f 100644
--- a/test/qr.cpp
+++ b/test/qr.cpp
@@ -133,4 +133,7 @@ void test_qr()
CALL_SUBTEST_6(qr_verify_assert<MatrixXd>());
CALL_SUBTEST_7(qr_verify_assert<MatrixXcf>());
CALL_SUBTEST_8(qr_verify_assert<MatrixXcd>());
+
+ // Test problem size constructors
+ CALL_SUBTEST_12(HouseholderQR<MatrixXf>(10, 20));
}
diff --git a/test/qr_colpivoting.cpp b/test/qr_colpivoting.cpp
index 0ca897b52..a34feedd9 100644
--- a/test/qr_colpivoting.cpp
+++ b/test/qr_colpivoting.cpp
@@ -157,4 +157,7 @@ void test_qr_colpivoting()
CALL_SUBTEST_2(qr_verify_assert<MatrixXd>());
CALL_SUBTEST_6(qr_verify_assert<MatrixXcf>());
CALL_SUBTEST_3(qr_verify_assert<MatrixXcd>());
+
+ // Test problem size constructors
+ CALL_SUBTEST_9(ColPivHouseholderQR<MatrixXf>(10, 20));
}
diff --git a/test/qr_fullpivoting.cpp b/test/qr_fullpivoting.cpp
index 7ad3af1fe..82c42c759 100644
--- a/test/qr_fullpivoting.cpp
+++ b/test/qr_fullpivoting.cpp
@@ -139,4 +139,7 @@ void test_qr_fullpivoting()
CALL_SUBTEST_2(qr_verify_assert<MatrixXd>());
CALL_SUBTEST_4(qr_verify_assert<MatrixXcf>());
CALL_SUBTEST_3(qr_verify_assert<MatrixXcd>());
+
+ // Test problem size constructors
+ CALL_SUBTEST_7(FullPivHouseholderQR<MatrixXf>(10, 20));
}
diff --git a/test/schur_complex.cpp b/test/schur_complex.cpp
index 3659a074c..b33411cf2 100644
--- a/test/schur_complex.cpp
+++ b/test/schur_complex.cpp
@@ -64,4 +64,7 @@ void test_schur_complex()
CALL_SUBTEST_2(( schur<MatrixXcf>(ei_random<int>(1,50)) ));
CALL_SUBTEST_3(( schur<Matrix<std::complex<float>, 1, 1> >() ));
CALL_SUBTEST_4(( schur<Matrix<float, 3, 3, Eigen::RowMajor> >() ));
+
+ // Test problem size constructors
+ CALL_SUBTEST_5(ComplexSchur<MatrixXf>(10));
}
diff --git a/test/schur_real.cpp b/test/schur_real.cpp
index c1747e2f5..d0aca4308 100644
--- a/test/schur_real.cpp
+++ b/test/schur_real.cpp
@@ -82,4 +82,7 @@ void test_schur_real()
CALL_SUBTEST_2(( schur<MatrixXd>(ei_random<int>(1,50)) ));
CALL_SUBTEST_3(( schur<Matrix<float, 1, 1> >() ));
CALL_SUBTEST_4(( schur<Matrix<double, 3, 3, Eigen::RowMajor> >() ));
+
+ // Test problem size constructors
+ CALL_SUBTEST_5(RealSchur<MatrixXf>(10));
}
diff --git a/test/svd.cpp b/test/svd.cpp
index 1b41d7032..9f3072d3b 100644
--- a/test/svd.cpp
+++ b/test/svd.cpp
@@ -113,4 +113,7 @@ void test_svd()
CALL_SUBTEST_2( svd_verify_assert<Matrix4d>() );
CALL_SUBTEST_3( svd_verify_assert<MatrixXf>() );
CALL_SUBTEST_4( svd_verify_assert<MatrixXd>() );
+
+ // Test problem size constructors
+ CALL_SUBTEST_9( SVD<MatrixXf>(10, 20) );
}