1*bf2c3715SXin Li // This file is triangularView of Eigen, a lightweight C++ template library
2*bf2c3715SXin Li // for linear algebra.
3*bf2c3715SXin Li //
4*bf2c3715SXin Li // Copyright (C) 2008-2009 Gael Guennebaud <[email protected]>
5*bf2c3715SXin Li //
6*bf2c3715SXin Li // This Source Code Form is subject to the terms of the Mozilla
7*bf2c3715SXin Li // Public License v. 2.0. If a copy of the MPL was not distributed
8*bf2c3715SXin Li // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9*bf2c3715SXin Li
10*bf2c3715SXin Li #ifdef EIGEN_TEST_PART_100
11*bf2c3715SXin Li # define EIGEN_NO_DEPRECATED_WARNING
12*bf2c3715SXin Li #endif
13*bf2c3715SXin Li
14*bf2c3715SXin Li #include "main.h"
15*bf2c3715SXin Li
16*bf2c3715SXin Li
triangular_deprecated(const MatrixType & m)17*bf2c3715SXin Li template<typename MatrixType> void triangular_deprecated(const MatrixType &m)
18*bf2c3715SXin Li {
19*bf2c3715SXin Li Index rows = m.rows();
20*bf2c3715SXin Li Index cols = m.cols();
21*bf2c3715SXin Li MatrixType m1, m2, m3, m4;
22*bf2c3715SXin Li m1.setRandom(rows,cols);
23*bf2c3715SXin Li m2.setRandom(rows,cols);
24*bf2c3715SXin Li m3 = m1; m4 = m2;
25*bf2c3715SXin Li // deprecated method:
26*bf2c3715SXin Li m1.template triangularView<Eigen::Upper>().swap(m2);
27*bf2c3715SXin Li // use this method instead:
28*bf2c3715SXin Li m3.template triangularView<Eigen::Upper>().swap(m4.template triangularView<Eigen::Upper>());
29*bf2c3715SXin Li VERIFY_IS_APPROX(m1,m3);
30*bf2c3715SXin Li VERIFY_IS_APPROX(m2,m4);
31*bf2c3715SXin Li // deprecated method:
32*bf2c3715SXin Li m1.template triangularView<Eigen::Lower>().swap(m4);
33*bf2c3715SXin Li // use this method instead:
34*bf2c3715SXin Li m3.template triangularView<Eigen::Lower>().swap(m2.template triangularView<Eigen::Lower>());
35*bf2c3715SXin Li VERIFY_IS_APPROX(m1,m3);
36*bf2c3715SXin Li VERIFY_IS_APPROX(m2,m4);
37*bf2c3715SXin Li }
38*bf2c3715SXin Li
39*bf2c3715SXin Li
triangular_square(const MatrixType & m)40*bf2c3715SXin Li template<typename MatrixType> void triangular_square(const MatrixType& m)
41*bf2c3715SXin Li {
42*bf2c3715SXin Li typedef typename MatrixType::Scalar Scalar;
43*bf2c3715SXin Li typedef typename NumTraits<Scalar>::Real RealScalar;
44*bf2c3715SXin Li typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
45*bf2c3715SXin Li
46*bf2c3715SXin Li RealScalar largerEps = 10*test_precision<RealScalar>();
47*bf2c3715SXin Li
48*bf2c3715SXin Li Index rows = m.rows();
49*bf2c3715SXin Li Index cols = m.cols();
50*bf2c3715SXin Li
51*bf2c3715SXin Li MatrixType m1 = MatrixType::Random(rows, cols),
52*bf2c3715SXin Li m2 = MatrixType::Random(rows, cols),
53*bf2c3715SXin Li m3(rows, cols),
54*bf2c3715SXin Li m4(rows, cols),
55*bf2c3715SXin Li r1(rows, cols),
56*bf2c3715SXin Li r2(rows, cols);
57*bf2c3715SXin Li VectorType v2 = VectorType::Random(rows);
58*bf2c3715SXin Li
59*bf2c3715SXin Li MatrixType m1up = m1.template triangularView<Upper>();
60*bf2c3715SXin Li MatrixType m2up = m2.template triangularView<Upper>();
61*bf2c3715SXin Li
62*bf2c3715SXin Li if (rows*cols>1)
63*bf2c3715SXin Li {
64*bf2c3715SXin Li VERIFY(m1up.isUpperTriangular());
65*bf2c3715SXin Li VERIFY(m2up.transpose().isLowerTriangular());
66*bf2c3715SXin Li VERIFY(!m2.isLowerTriangular());
67*bf2c3715SXin Li }
68*bf2c3715SXin Li
69*bf2c3715SXin Li // VERIFY_IS_APPROX(m1up.transpose() * m2, m1.upper().transpose().lower() * m2);
70*bf2c3715SXin Li
71*bf2c3715SXin Li // test overloaded operator+=
72*bf2c3715SXin Li r1.setZero();
73*bf2c3715SXin Li r2.setZero();
74*bf2c3715SXin Li r1.template triangularView<Upper>() += m1;
75*bf2c3715SXin Li r2 += m1up;
76*bf2c3715SXin Li VERIFY_IS_APPROX(r1,r2);
77*bf2c3715SXin Li
78*bf2c3715SXin Li // test overloaded operator=
79*bf2c3715SXin Li m1.setZero();
80*bf2c3715SXin Li m1.template triangularView<Upper>() = m2.transpose() + m2;
81*bf2c3715SXin Li m3 = m2.transpose() + m2;
82*bf2c3715SXin Li VERIFY_IS_APPROX(m3.template triangularView<Lower>().transpose().toDenseMatrix(), m1);
83*bf2c3715SXin Li
84*bf2c3715SXin Li // test overloaded operator=
85*bf2c3715SXin Li m1.setZero();
86*bf2c3715SXin Li m1.template triangularView<Lower>() = m2.transpose() + m2;
87*bf2c3715SXin Li VERIFY_IS_APPROX(m3.template triangularView<Lower>().toDenseMatrix(), m1);
88*bf2c3715SXin Li
89*bf2c3715SXin Li VERIFY_IS_APPROX(m3.template triangularView<Lower>().conjugate().toDenseMatrix(),
90*bf2c3715SXin Li m3.conjugate().template triangularView<Lower>().toDenseMatrix());
91*bf2c3715SXin Li
92*bf2c3715SXin Li m1 = MatrixType::Random(rows, cols);
93*bf2c3715SXin Li for (int i=0; i<rows; ++i)
94*bf2c3715SXin Li while (numext::abs2(m1(i,i))<RealScalar(1e-1)) m1(i,i) = internal::random<Scalar>();
95*bf2c3715SXin Li
96*bf2c3715SXin Li Transpose<MatrixType> trm4(m4);
97*bf2c3715SXin Li // test back and forward substitution with a vector as the rhs
98*bf2c3715SXin Li m3 = m1.template triangularView<Upper>();
99*bf2c3715SXin Li VERIFY(v2.isApprox(m3.adjoint() * (m1.adjoint().template triangularView<Lower>().solve(v2)), largerEps));
100*bf2c3715SXin Li m3 = m1.template triangularView<Lower>();
101*bf2c3715SXin Li VERIFY(v2.isApprox(m3.transpose() * (m1.transpose().template triangularView<Upper>().solve(v2)), largerEps));
102*bf2c3715SXin Li m3 = m1.template triangularView<Upper>();
103*bf2c3715SXin Li VERIFY(v2.isApprox(m3 * (m1.template triangularView<Upper>().solve(v2)), largerEps));
104*bf2c3715SXin Li m3 = m1.template triangularView<Lower>();
105*bf2c3715SXin Li VERIFY(v2.isApprox(m3.conjugate() * (m1.conjugate().template triangularView<Lower>().solve(v2)), largerEps));
106*bf2c3715SXin Li
107*bf2c3715SXin Li // test back and forward substitution with a matrix as the rhs
108*bf2c3715SXin Li m3 = m1.template triangularView<Upper>();
109*bf2c3715SXin Li VERIFY(m2.isApprox(m3.adjoint() * (m1.adjoint().template triangularView<Lower>().solve(m2)), largerEps));
110*bf2c3715SXin Li m3 = m1.template triangularView<Lower>();
111*bf2c3715SXin Li VERIFY(m2.isApprox(m3.transpose() * (m1.transpose().template triangularView<Upper>().solve(m2)), largerEps));
112*bf2c3715SXin Li m3 = m1.template triangularView<Upper>();
113*bf2c3715SXin Li VERIFY(m2.isApprox(m3 * (m1.template triangularView<Upper>().solve(m2)), largerEps));
114*bf2c3715SXin Li m3 = m1.template triangularView<Lower>();
115*bf2c3715SXin Li VERIFY(m2.isApprox(m3.conjugate() * (m1.conjugate().template triangularView<Lower>().solve(m2)), largerEps));
116*bf2c3715SXin Li
117*bf2c3715SXin Li // check M * inv(L) using in place API
118*bf2c3715SXin Li m4 = m3;
119*bf2c3715SXin Li m1.transpose().template triangularView<Eigen::Upper>().solveInPlace(trm4);
120*bf2c3715SXin Li VERIFY_IS_APPROX(m4 * m1.template triangularView<Eigen::Lower>(), m3);
121*bf2c3715SXin Li
122*bf2c3715SXin Li // check M * inv(U) using in place API
123*bf2c3715SXin Li m3 = m1.template triangularView<Upper>();
124*bf2c3715SXin Li m4 = m3;
125*bf2c3715SXin Li m3.transpose().template triangularView<Eigen::Lower>().solveInPlace(trm4);
126*bf2c3715SXin Li VERIFY_IS_APPROX(m4 * m1.template triangularView<Eigen::Upper>(), m3);
127*bf2c3715SXin Li
128*bf2c3715SXin Li // check solve with unit diagonal
129*bf2c3715SXin Li m3 = m1.template triangularView<UnitUpper>();
130*bf2c3715SXin Li VERIFY(m2.isApprox(m3 * (m1.template triangularView<UnitUpper>().solve(m2)), largerEps));
131*bf2c3715SXin Li
132*bf2c3715SXin Li // VERIFY(( m1.template triangularView<Upper>()
133*bf2c3715SXin Li // * m2.template triangularView<Upper>()).isUpperTriangular());
134*bf2c3715SXin Li
135*bf2c3715SXin Li // test swap
136*bf2c3715SXin Li m1.setOnes();
137*bf2c3715SXin Li m2.setZero();
138*bf2c3715SXin Li m2.template triangularView<Upper>().swap(m1.template triangularView<Eigen::Upper>());
139*bf2c3715SXin Li m3.setZero();
140*bf2c3715SXin Li m3.template triangularView<Upper>().setOnes();
141*bf2c3715SXin Li VERIFY_IS_APPROX(m2,m3);
142*bf2c3715SXin Li VERIFY_RAISES_STATIC_ASSERT(m1.template triangularView<Eigen::Lower>().swap(m2.template triangularView<Eigen::Upper>()));
143*bf2c3715SXin Li
144*bf2c3715SXin Li m1.setRandom();
145*bf2c3715SXin Li m3 = m1.template triangularView<Upper>();
146*bf2c3715SXin Li Matrix<Scalar, MatrixType::ColsAtCompileTime, Dynamic> m5(cols, internal::random<int>(1,20)); m5.setRandom();
147*bf2c3715SXin Li Matrix<Scalar, Dynamic, MatrixType::RowsAtCompileTime> m6(internal::random<int>(1,20), rows); m6.setRandom();
148*bf2c3715SXin Li VERIFY_IS_APPROX(m1.template triangularView<Upper>() * m5, m3*m5);
149*bf2c3715SXin Li VERIFY_IS_APPROX(m6*m1.template triangularView<Upper>(), m6*m3);
150*bf2c3715SXin Li
151*bf2c3715SXin Li m1up = m1.template triangularView<Upper>();
152*bf2c3715SXin Li VERIFY_IS_APPROX(m1.template selfadjointView<Upper>().template triangularView<Upper>().toDenseMatrix(), m1up);
153*bf2c3715SXin Li VERIFY_IS_APPROX(m1up.template selfadjointView<Upper>().template triangularView<Upper>().toDenseMatrix(), m1up);
154*bf2c3715SXin Li VERIFY_IS_APPROX(m1.template selfadjointView<Upper>().template triangularView<Lower>().toDenseMatrix(), m1up.adjoint());
155*bf2c3715SXin Li VERIFY_IS_APPROX(m1up.template selfadjointView<Upper>().template triangularView<Lower>().toDenseMatrix(), m1up.adjoint());
156*bf2c3715SXin Li
157*bf2c3715SXin Li VERIFY_IS_APPROX(m1.template selfadjointView<Upper>().diagonal(), m1.diagonal());
158*bf2c3715SXin Li
159*bf2c3715SXin Li m3.setRandom();
160*bf2c3715SXin Li const MatrixType& m3c(m3);
161*bf2c3715SXin Li VERIFY( is_same_type(m3c.template triangularView<Lower>(),m3.template triangularView<Lower>().template conjugateIf<false>()) );
162*bf2c3715SXin Li VERIFY( is_same_type(m3c.template triangularView<Lower>().conjugate(),m3.template triangularView<Lower>().template conjugateIf<true>()) );
163*bf2c3715SXin Li VERIFY_IS_APPROX(m3.template triangularView<Lower>().template conjugateIf<true>().toDenseMatrix(),
164*bf2c3715SXin Li m3.conjugate().template triangularView<Lower>().toDenseMatrix());
165*bf2c3715SXin Li VERIFY_IS_APPROX(m3.template triangularView<Lower>().template conjugateIf<false>().toDenseMatrix(),
166*bf2c3715SXin Li m3.template triangularView<Lower>().toDenseMatrix());
167*bf2c3715SXin Li
168*bf2c3715SXin Li VERIFY( is_same_type(m3c.template selfadjointView<Lower>(),m3.template selfadjointView<Lower>().template conjugateIf<false>()) );
169*bf2c3715SXin Li VERIFY( is_same_type(m3c.template selfadjointView<Lower>().conjugate(),m3.template selfadjointView<Lower>().template conjugateIf<true>()) );
170*bf2c3715SXin Li VERIFY_IS_APPROX(m3.template selfadjointView<Lower>().template conjugateIf<true>().toDenseMatrix(),
171*bf2c3715SXin Li m3.conjugate().template selfadjointView<Lower>().toDenseMatrix());
172*bf2c3715SXin Li VERIFY_IS_APPROX(m3.template selfadjointView<Lower>().template conjugateIf<false>().toDenseMatrix(),
173*bf2c3715SXin Li m3.template selfadjointView<Lower>().toDenseMatrix());
174*bf2c3715SXin Li
175*bf2c3715SXin Li }
176*bf2c3715SXin Li
177*bf2c3715SXin Li
triangular_rect(const MatrixType & m)178*bf2c3715SXin Li template<typename MatrixType> void triangular_rect(const MatrixType& m)
179*bf2c3715SXin Li {
180*bf2c3715SXin Li typedef typename MatrixType::Scalar Scalar;
181*bf2c3715SXin Li typedef typename NumTraits<Scalar>::Real RealScalar;
182*bf2c3715SXin Li enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime };
183*bf2c3715SXin Li
184*bf2c3715SXin Li Index rows = m.rows();
185*bf2c3715SXin Li Index cols = m.cols();
186*bf2c3715SXin Li
187*bf2c3715SXin Li MatrixType m1 = MatrixType::Random(rows, cols),
188*bf2c3715SXin Li m2 = MatrixType::Random(rows, cols),
189*bf2c3715SXin Li m3(rows, cols),
190*bf2c3715SXin Li m4(rows, cols),
191*bf2c3715SXin Li r1(rows, cols),
192*bf2c3715SXin Li r2(rows, cols);
193*bf2c3715SXin Li
194*bf2c3715SXin Li MatrixType m1up = m1.template triangularView<Upper>();
195*bf2c3715SXin Li MatrixType m2up = m2.template triangularView<Upper>();
196*bf2c3715SXin Li
197*bf2c3715SXin Li if (rows>1 && cols>1)
198*bf2c3715SXin Li {
199*bf2c3715SXin Li VERIFY(m1up.isUpperTriangular());
200*bf2c3715SXin Li VERIFY(m2up.transpose().isLowerTriangular());
201*bf2c3715SXin Li VERIFY(!m2.isLowerTriangular());
202*bf2c3715SXin Li }
203*bf2c3715SXin Li
204*bf2c3715SXin Li // test overloaded operator+=
205*bf2c3715SXin Li r1.setZero();
206*bf2c3715SXin Li r2.setZero();
207*bf2c3715SXin Li r1.template triangularView<Upper>() += m1;
208*bf2c3715SXin Li r2 += m1up;
209*bf2c3715SXin Li VERIFY_IS_APPROX(r1,r2);
210*bf2c3715SXin Li
211*bf2c3715SXin Li // test overloaded operator=
212*bf2c3715SXin Li m1.setZero();
213*bf2c3715SXin Li m1.template triangularView<Upper>() = 3 * m2;
214*bf2c3715SXin Li m3 = 3 * m2;
215*bf2c3715SXin Li VERIFY_IS_APPROX(m3.template triangularView<Upper>().toDenseMatrix(), m1);
216*bf2c3715SXin Li
217*bf2c3715SXin Li
218*bf2c3715SXin Li m1.setZero();
219*bf2c3715SXin Li m1.template triangularView<Lower>() = 3 * m2;
220*bf2c3715SXin Li VERIFY_IS_APPROX(m3.template triangularView<Lower>().toDenseMatrix(), m1);
221*bf2c3715SXin Li
222*bf2c3715SXin Li m1.setZero();
223*bf2c3715SXin Li m1.template triangularView<StrictlyUpper>() = 3 * m2;
224*bf2c3715SXin Li VERIFY_IS_APPROX(m3.template triangularView<StrictlyUpper>().toDenseMatrix(), m1);
225*bf2c3715SXin Li
226*bf2c3715SXin Li
227*bf2c3715SXin Li m1.setZero();
228*bf2c3715SXin Li m1.template triangularView<StrictlyLower>() = 3 * m2;
229*bf2c3715SXin Li VERIFY_IS_APPROX(m3.template triangularView<StrictlyLower>().toDenseMatrix(), m1);
230*bf2c3715SXin Li m1.setRandom();
231*bf2c3715SXin Li m2 = m1.template triangularView<Upper>();
232*bf2c3715SXin Li VERIFY(m2.isUpperTriangular());
233*bf2c3715SXin Li VERIFY(!m2.isLowerTriangular());
234*bf2c3715SXin Li m2 = m1.template triangularView<StrictlyUpper>();
235*bf2c3715SXin Li VERIFY(m2.isUpperTriangular());
236*bf2c3715SXin Li VERIFY(m2.diagonal().isMuchSmallerThan(RealScalar(1)));
237*bf2c3715SXin Li m2 = m1.template triangularView<UnitUpper>();
238*bf2c3715SXin Li VERIFY(m2.isUpperTriangular());
239*bf2c3715SXin Li m2.diagonal().array() -= Scalar(1);
240*bf2c3715SXin Li VERIFY(m2.diagonal().isMuchSmallerThan(RealScalar(1)));
241*bf2c3715SXin Li m2 = m1.template triangularView<Lower>();
242*bf2c3715SXin Li VERIFY(m2.isLowerTriangular());
243*bf2c3715SXin Li VERIFY(!m2.isUpperTriangular());
244*bf2c3715SXin Li m2 = m1.template triangularView<StrictlyLower>();
245*bf2c3715SXin Li VERIFY(m2.isLowerTriangular());
246*bf2c3715SXin Li VERIFY(m2.diagonal().isMuchSmallerThan(RealScalar(1)));
247*bf2c3715SXin Li m2 = m1.template triangularView<UnitLower>();
248*bf2c3715SXin Li VERIFY(m2.isLowerTriangular());
249*bf2c3715SXin Li m2.diagonal().array() -= Scalar(1);
250*bf2c3715SXin Li VERIFY(m2.diagonal().isMuchSmallerThan(RealScalar(1)));
251*bf2c3715SXin Li // test swap
252*bf2c3715SXin Li m1.setOnes();
253*bf2c3715SXin Li m2.setZero();
254*bf2c3715SXin Li m2.template triangularView<Upper>().swap(m1.template triangularView<Eigen::Upper>());
255*bf2c3715SXin Li m3.setZero();
256*bf2c3715SXin Li m3.template triangularView<Upper>().setOnes();
257*bf2c3715SXin Li VERIFY_IS_APPROX(m2,m3);
258*bf2c3715SXin Li }
259*bf2c3715SXin Li
bug_159()260*bf2c3715SXin Li void bug_159()
261*bf2c3715SXin Li {
262*bf2c3715SXin Li Matrix3d m = Matrix3d::Random().triangularView<Lower>();
263*bf2c3715SXin Li EIGEN_UNUSED_VARIABLE(m)
264*bf2c3715SXin Li }
265*bf2c3715SXin Li
EIGEN_DECLARE_TEST(triangular)266*bf2c3715SXin Li EIGEN_DECLARE_TEST(triangular)
267*bf2c3715SXin Li {
268*bf2c3715SXin Li int maxsize = (std::min)(EIGEN_TEST_MAX_SIZE,20);
269*bf2c3715SXin Li for(int i = 0; i < g_repeat ; i++)
270*bf2c3715SXin Li {
271*bf2c3715SXin Li int r = internal::random<int>(2,maxsize); TEST_SET_BUT_UNUSED_VARIABLE(r)
272*bf2c3715SXin Li int c = internal::random<int>(2,maxsize); TEST_SET_BUT_UNUSED_VARIABLE(c)
273*bf2c3715SXin Li
274*bf2c3715SXin Li CALL_SUBTEST_1( triangular_square(Matrix<float, 1, 1>()) );
275*bf2c3715SXin Li CALL_SUBTEST_2( triangular_square(Matrix<float, 2, 2>()) );
276*bf2c3715SXin Li CALL_SUBTEST_3( triangular_square(Matrix3d()) );
277*bf2c3715SXin Li CALL_SUBTEST_4( triangular_square(Matrix<std::complex<float>,8, 8>()) );
278*bf2c3715SXin Li CALL_SUBTEST_5( triangular_square(MatrixXcd(r,r)) );
279*bf2c3715SXin Li CALL_SUBTEST_6( triangular_square(Matrix<float,Dynamic,Dynamic,RowMajor>(r, r)) );
280*bf2c3715SXin Li
281*bf2c3715SXin Li CALL_SUBTEST_7( triangular_rect(Matrix<float, 4, 5>()) );
282*bf2c3715SXin Li CALL_SUBTEST_8( triangular_rect(Matrix<double, 6, 2>()) );
283*bf2c3715SXin Li CALL_SUBTEST_9( triangular_rect(MatrixXcf(r, c)) );
284*bf2c3715SXin Li CALL_SUBTEST_5( triangular_rect(MatrixXcd(r, c)) );
285*bf2c3715SXin Li CALL_SUBTEST_6( triangular_rect(Matrix<float,Dynamic,Dynamic,RowMajor>(r, c)) );
286*bf2c3715SXin Li
287*bf2c3715SXin Li CALL_SUBTEST_100( triangular_deprecated(Matrix<float, 5, 7>()) );
288*bf2c3715SXin Li CALL_SUBTEST_100( triangular_deprecated(MatrixXd(r,c)) );
289*bf2c3715SXin Li }
290*bf2c3715SXin Li
291*bf2c3715SXin Li CALL_SUBTEST_1( bug_159() );
292*bf2c3715SXin Li }
293