MyraMath
ssytrdU


Source: tests/dense/sytrd.cpp

1 // ========================================================================= //
2 // This file is part of MyraMath, copyright (c) 2014-2019 by Ryan A Chilton //
3 // and distributed by MyraCore, LLC. See LICENSE.txt for license terms. //
4 // ========================================================================= //
5 
11 // Containers.
13 #include <myramath/dense/Matrix.h>
15 #include <myramath/dense/Vector.h>
18 
19 // Algorithms.
20 #include <myramath/dense/diag.h>
21 #include <myramath/dense/gemm.h>
22 #include <myramath/dense/ormqr.h>
23 #include <myramath/dense/ormlq.h>
24 #include <myramath/dense/sytrd.h>
25 #include <myramath/dense/syev.h>
26 #include <myramath/dense/stev.h>
30 
31 // Reporting.
32 #include <tests/myratest.h>
33 
34 using namespace myra;
35 
36 namespace {
37 
38 // Tests sytrd('L'ower)
39 template<class Precision> void testL(int N, Precision tolerance)
40  {
41  myra::out() << typestring<Precision>() << std::endl;
42  // Make random real symmetric matrix A, a few copies.
43  auto A = Matrix<Precision>::random(N,N);
44  A += transpose(A);
45  auto A1 = A;
46  auto A2 = A;
47  // Tridiagonalize A = Q*T*Q'.
48  auto tau = sytrd_inplace('L',A);
49  // Extract T factor from diagonal/subdiagonal of A.
50  auto T = Matrix<Precision>::zeros(N,N);
51  auto T0 = Vector<Precision>::zeros(N);
52  auto T1 = Vector<Precision>::zeros(N-1);
53  for (int n = 0; n < N; ++n)
54  T0(n) = T(n,n) = A(n,n);
55  for (int n = 0; n < N-1; ++n)
56  T1(n) = T(n+1,n) = T(n,n+1) = A(n+1,n);
57  // Verify A and T have same eigenvalues, using syev() and stev()
58  auto lambda_A = syev_inplace(A1);
59  auto lambda_T = stev(T0,T1).second;
60  Precision error1 = frobenius(lambda_A-lambda_T) / frobenius(lambda_A);
61  myra::out() << "|eig(A)-eig(T)| = " << error1 << std::endl;
62  REQUIRE(error1 < tolerance);
63  // Verify Q'*A*Q = T (ormqr)
64  auto Q = A.bottom(N-1).left(N-1);
65  ormqr_inplace(Q,tau,A2.bottom(N-1),'L','T');
66  ormqr_inplace(Q,tau,A2.right(N-1),'R','N');
67  Precision error2 = frobenius(A2-T) / frobenius(T);
68  myra::out() << "|Q'*A*Q - T| = " << error2 << std::endl;
69  REQUIRE(error2 < tolerance);
70  }
71 
72 // Tests sytrd('U'pper)
73 template<class Precision> void testU(int N, Precision tolerance)
74  {
75  myra::out() << typestring<Precision>() << std::endl;
76  // Make random real symmetric matrix A, a few copies.
77  auto A = Matrix<Precision>::random(N,N);
78  A += transpose(A);
79  auto A1 = A;
80  auto A2 = A;
81  // Tridiagonalize A = Q'*T*Q.
82  auto tau = sytrd_inplace('U',A);
83  // Extract T factor from diagonal/superdiagonal of A.
84  auto T = Matrix<Precision>::zeros(N,N);
85  auto T0 = Vector<Precision>::zeros(N);
86  auto T1 = Vector<Precision>::zeros(N-1);
87  for (int n = 0; n < N; ++n)
88  T0(n) = T(n,n) = A(n,n);
89  for (int n = 0; n < N-1; ++n)
90  T1(n) = T(n+1,n) = T(n,n+1) = A(n,n+1);
91  // Verify A and T have same eigenvalues, using syev() and stev()
92  auto lambda_A = syev_inplace(A1);
93  auto lambda_T = stev(T0,T1).second;
94  Precision error1 = frobenius(lambda_A-lambda_T) / frobenius(lambda_A);
95  myra::out() << "|eig(A)-eig(T)| = " << error1 << std::endl;
96  REQUIRE(error1 < tolerance);
97  // Verify Q*A*Q' = T (ormlq)
98  auto Q = A.right(N-1).top(N-1);
99  ormlq_inplace(Q,tau,A2.bottom(N-1),'L','N');
100  ormlq_inplace(Q,tau,A2.right(N-1),'R','T');
101  Precision error2 = frobenius(A2-T) / frobenius(T);
102  myra::out() << "|Q*A*Q' - T| = " << error2 << std::endl;
103  // Note this check will fail with reference LAPACK/BLAS, since
104  // they use a QL decomposition instead of an LQ decomposition.
105  // REQUIRE(error2 < tolerance);
106  }
107 
108 } // namespace
109 
110 ADD_TEST("ssytrdL","[dense][lapack]")
111  { testL<float>(157,1.0e-4f); }
112 
113 ADD_TEST("dsytrdL","[dense][lapack]")
114  { testL<double>(157,1.0e-10); }
115 
116 ADD_TEST("ssytrdU","[dense][lapack]")
117  { testU<float>(157,1.0e-4f); }
118 
119 ADD_TEST("dsytrdU","[dense][lapack]")
120  { testU<double>(157,1.0e-10); }
121 
CMatrixRange< Number > right(int j) const
Returns a MatrixRange over the j rightmost columns, this(:,J-j:J)
Definition: Matrix.cpp:281
Interface class for representing subranges of DiagonalMatrix&#39;s.
Computes all eigenpairs of real symmetric Matrix.
Interface class for representing subranges of dense Matrix&#39;s.
Routines to apply a Householder Q, as previously factored via gelqf_inplace()
static Matrix< Number > zeros(int I, int J)
Generates a zeros Matrix of specified size.
Definition: Matrix.cpp:357
Routines for computing Frobenius norms of various algebraic containers.
static Matrix< Number > random(int I, int J)
Generates a random Matrix of specified size.
Definition: Matrix.cpp:353
CMatrixRange< Number > bottom(int i) const
Returns a MatrixRange over the i bottommost rows, this(I-i:I,:)
Definition: Matrix.cpp:225
Householder tridiagonalization of a real symmetric Matrix.
Replaces small values with explicit zeros.
Definition: syntax.dox:1
Returns a transposed copy of a Matrix. The inplace version only works on a square operand...
Various utility functions/classes related to scalar Number types.
static Vector< Number > zeros(int N)
Generates a zeros Vector of specified size.
Definition: Vector.cpp:280
General purpose dense matrix container, O(i*j) storage.
Container for either a column vector or row vector (depends upon the usage context) ...
Returns the diagonal of a dense Matrix or LowerMatrix.
Routines to apply a Householder Q, as previously factored via geqrf_inplace()
Computes all eigenpairs of real symmetric tridiagonal matrix.
Container for a diagonal matrix, O(n) storage. Used by SVD, row/column scaling, etc.
Variety of routines all for dense Matrix*Matrix multiplies. Delegates to the BLAS.


Results: [PASS]

float
|eig(A)-eig(T)| = 9.00067e-07
|Q*A*Q' - T| = 9.041e-07


Go back to Summary of /test programs.