mirror of
https://gitlab.com/libeigen/eigen.git
synced 2024-12-21 07:19:46 +08:00
139 lines
5.7 KiB
C++
139 lines
5.7 KiB
C++
// This file is triangularView of Eigen, a lightweight C++ template library
|
|
// for linear algebra.
|
|
//
|
|
// Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@gmail.com>
|
|
//
|
|
// Eigen is free software; you can redistribute it and/or
|
|
// modify it under the terms of the GNU Lesser General Public
|
|
// License as published by the Free Software Foundation; either
|
|
// version 3 of the License, or (at your option) any later version.
|
|
//
|
|
// Alternatively, you can redistribute it and/or
|
|
// modify it under the terms of the GNU General Public License as
|
|
// published by the Free Software Foundation; either version 2 of
|
|
// the License, or (at your option) any later version.
|
|
//
|
|
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
|
|
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
|
|
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
|
|
// GNU General Public License for more details.
|
|
//
|
|
// You should have received a copy of the GNU Lesser General Public
|
|
// License and a copy of the GNU General Public License along with
|
|
// Eigen. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
#include "main.h"
|
|
|
|
template<typename MatrixType> void triangular(const MatrixType& m)
|
|
{
|
|
typedef typename MatrixType::Scalar Scalar;
|
|
typedef typename NumTraits<Scalar>::Real RealScalar;
|
|
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
|
|
|
|
RealScalar largerEps = 10*test_precision<RealScalar>();
|
|
|
|
int rows = m.rows();
|
|
int cols = m.cols();
|
|
|
|
MatrixType m1 = MatrixType::Random(rows, cols),
|
|
m2 = MatrixType::Random(rows, cols),
|
|
m3(rows, cols),
|
|
m4(rows, cols),
|
|
r1(rows, cols),
|
|
r2(rows, cols),
|
|
mzero = MatrixType::Zero(rows, cols),
|
|
mones = MatrixType::Ones(rows, cols),
|
|
identity = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
|
::Identity(rows, rows),
|
|
square = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
|
::Random(rows, rows);
|
|
VectorType v1 = VectorType::Random(rows),
|
|
v2 = VectorType::Random(rows),
|
|
vzero = VectorType::Zero(rows);
|
|
|
|
MatrixType m1up = m1.template triangularView<Eigen::UpperTriangular>();
|
|
MatrixType m2up = m2.template triangularView<Eigen::UpperTriangular>();
|
|
|
|
if (rows*cols>1)
|
|
{
|
|
VERIFY(m1up.isUpperTriangular());
|
|
VERIFY(m2up.transpose().isLowerTriangular());
|
|
VERIFY(!m2.isLowerTriangular());
|
|
}
|
|
|
|
// VERIFY_IS_APPROX(m1up.transpose() * m2, m1.upper().transpose().lower() * m2);
|
|
|
|
// test overloaded operator+=
|
|
r1.setZero();
|
|
r2.setZero();
|
|
r1.template triangularView<Eigen::UpperTriangular>() += m1;
|
|
r2 += m1up;
|
|
VERIFY_IS_APPROX(r1,r2);
|
|
|
|
// test overloaded operator=
|
|
m1.setZero();
|
|
m1.template triangularView<Eigen::UpperTriangular>() = (m2.transpose() * m2).lazy();
|
|
m3 = m2.transpose() * m2;
|
|
VERIFY_IS_APPROX(m3.template triangularView<Eigen::LowerTriangular>().transpose().toDense(), m1);
|
|
|
|
// test overloaded operator=
|
|
m1.setZero();
|
|
m1.template triangularView<Eigen::LowerTriangular>() = (m2.transpose() * m2).lazy();
|
|
VERIFY_IS_APPROX(m3.template triangularView<Eigen::LowerTriangular>().toDense(), m1);
|
|
|
|
m1 = MatrixType::Random(rows, cols);
|
|
for (int i=0; i<rows; ++i)
|
|
while (ei_abs2(m1(i,i))<1e-3) m1(i,i) = ei_random<Scalar>();
|
|
|
|
Transpose<MatrixType> trm4(m4);
|
|
// test back and forward subsitution
|
|
m3 = m1.template triangularView<Eigen::UpperTriangular>();
|
|
VERIFY(m2.isApprox(m3.adjoint() * (m1.adjoint().template triangularView<Eigen::LowerTriangular>().solve(m2)), largerEps));
|
|
m3 = m1.template triangularView<Eigen::LowerTriangular>();
|
|
VERIFY(m2.isApprox(m3.transpose() * (m1.transpose().template triangularView<Eigen::UpperTriangular>().solve(m2)), largerEps));
|
|
m3 = m1.template triangularView<Eigen::UpperTriangular>();
|
|
VERIFY(m2.isApprox(m3 * (m1.template triangularView<Eigen::UpperTriangular>().solve(m2)), largerEps));
|
|
m3 = m1.template triangularView<Eigen::LowerTriangular>();
|
|
VERIFY(m2.isApprox(m3.conjugate() * (m1.conjugate().template triangularView<Eigen::LowerTriangular>().solve(m2)), largerEps));
|
|
|
|
// check M * inv(L) using in place API
|
|
m4 = m3;
|
|
m3.transpose().template triangularView<Eigen::UpperTriangular>().solveInPlace(trm4);
|
|
VERIFY(m4.cwise().abs().isIdentity(test_precision<RealScalar>()));
|
|
|
|
// check M * inv(U) using in place API
|
|
m3 = m1.template triangularView<Eigen::UpperTriangular>();
|
|
m4 = m3;
|
|
m3.transpose().template triangularView<Eigen::LowerTriangular>().solveInPlace(trm4);
|
|
VERIFY(m4.cwise().abs().isIdentity(test_precision<RealScalar>()));
|
|
|
|
// check solve with unit diagonal
|
|
m3 = m1.template triangularView<Eigen::UnitUpperTriangular>();
|
|
VERIFY(m2.isApprox(m3 * (m1.template triangularView<Eigen::UnitUpperTriangular>().solve(m2)), largerEps));
|
|
|
|
// VERIFY(( m1.template triangularView<Eigen::UpperTriangular>()
|
|
// * m2.template triangularView<Eigen::UpperTriangular>()).isUpperTriangular());
|
|
|
|
// test swap
|
|
m1.setOnes();
|
|
m2.setZero();
|
|
m2.template triangularView<Eigen::UpperTriangular>().swap(m1);
|
|
m3.setZero();
|
|
m3.template triangularView<Eigen::UpperTriangular>().setOnes();
|
|
VERIFY_IS_APPROX(m2,m3);
|
|
|
|
}
|
|
|
|
void test_triangular()
|
|
{
|
|
for(int i = 0; i < g_repeat ; i++) {
|
|
CALL_SUBTEST( triangular(Matrix<float, 1, 1>()) );
|
|
CALL_SUBTEST( triangular(Matrix<float, 2, 2>()) );
|
|
CALL_SUBTEST( triangular(Matrix3d()) );
|
|
CALL_SUBTEST( triangular(MatrixXcf(4, 4)) );
|
|
CALL_SUBTEST( triangular(Matrix<std::complex<float>,8, 8>()) );
|
|
CALL_SUBTEST( triangular(MatrixXd(17,17)) );
|
|
CALL_SUBTEST( triangular(Matrix<float,Dynamic,Dynamic,RowMajor>(5, 5)) );
|
|
}
|
|
}
|