MyraMath
zLDLHSolver3


Source: tests/dense/ZLDLHSolver.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>
16 
17 // Algorithms.
19 #include <myramath/dense/gemm.h>
20 #include <myramath/dense/herk.h>
25 #include <myramath/dense/inverse.h>
26 #include <myramath/dense/gesvd.h>
27 
28 // Solver class under test.
30 
31 // Serialization testing.
33 
34 // Reporting.
35 #include <tests/myratest.h>
37 
38 using namespace myra;
39 
40 namespace {
41 
42 template<class Precision> void test(int I, int J, Precision tolerance)
43  {
44  typedef std::complex<Precision> Number;
45  myra::out() << typestring<Number>() << std::endl;
46  // Construct random hermitian matrix, A = G+G'
47  auto G = Matrix<Number>::random(I,I);
48  auto A = G + hermitian(G);
49 
50  // Restrict A to lower triangle, generate LDL' solver for it.
52  typedef ZLDLHSolver<Precision> Solver;
53  Solver solver(L);
54  // Total of 8 solve() variations to validiate, side['L',R'] x op['N','T','H','C']
55  // Solve A*X=B
56  {
57  auto X = Matrix<Number>::random(I,J);
58  auto B = A*X;
59  solver.solve(B,'L','N');
60  Precision error = frobenius(B-X)/frobenius(X);
61  myra::out() << " |inv(A)*B-X| = " << error << std::endl;
62  REQUIRE(error < tolerance);
63  }
64  // Solve transpose(A)*X=B
65  {
66  auto X = Matrix<Number>::random(I,J);
67  auto B = transpose(A)*X;
68  solver.solve(B,'L','T');
69  Precision error = frobenius(B-X)/frobenius(X);
70  myra::out() << " |inv(A^T)*B-X| = " << error << std::endl;
71  REQUIRE(error < tolerance);
72  }
73  // Solve hermitian(A)*X=B
74  {
75  auto X = Matrix<Number>::random(I,J);
76  auto B = hermitian(A)*X;
77  solver.solve(B,'L','H');
78  Precision error = frobenius(B-X)/frobenius(X);
79  myra::out() << " |inv(A^H)*B-X| = " << error << std::endl;
80  REQUIRE(error < tolerance);
81  }
82  // Solve conjugate(A)*X=B
83  {
84  auto X = Matrix<Number>::random(I,J);
85  auto B = conjugate(A)*X;
86  solver.solve(B,'L','C');
87  Precision error = frobenius(B-X)/frobenius(X);
88  myra::out() << " |inv(A^C)*B-X| = " << error << std::endl;
89  REQUIRE(error < tolerance);
90  }
91 
92  // Solve X*A=B
93  {
94  auto X = Matrix<Number>::random(J,I);
95  auto B = X*A;
96  solver.solve(B,'R','N');
97  Precision error = frobenius(B-X)/frobenius(X);
98  myra::out() << " |B*inv(A)-X| = " << error << std::endl;
99  REQUIRE(error < tolerance);
100  }
101  // Solve X*transpose(A)=B
102  {
103  auto X = Matrix<Number>::random(J,I);
104  auto B = X*transpose(A);
105  solver.solve(B,'R','T');
106  Precision error = frobenius(B-X)/frobenius(X);
107  myra::out() << " |B*inv(A^T)-X| = " << error << std::endl;
108  REQUIRE(error < tolerance);
109  }
110  // Solve X*hermitian(A)=B
111  {
112  auto X = Matrix<Number>::random(J,I);
113  auto B = X*hermitian(A);
114  solver.solve(B,'R','H');
115  Precision error = frobenius(B-X)/frobenius(X);
116  myra::out() << " |B*inv(A^H)-X| = " << error << std::endl;
117  REQUIRE(error < tolerance);
118  }
119  // Solve X*conjugate(A)=B
120  {
121  auto X = Matrix<Number>::random(J,I);
122  auto B = X*conjugate(A);
123  solver.solve(B,'R','C');
124  Precision error = frobenius(B-X)/frobenius(X);
125  myra::out() << " |B*inv(A^C)-X| = " << error << std::endl;
126  REQUIRE(error < tolerance);
127  }
128  // Generate symmetric schur complement.
129  {
130  // ZLDLHSolver is well suited for computing symmetric schur complements,
131  // S = Bh*inv(A)*B. Try it out and compare to a reference calculation.
132  // (In the case of sparse A/B, this approach is dramatically more efficient)
133  // Reference answer.
134  auto B = Matrix<Number>::random(I,J);
135  LowerMatrix<Number> S1 (hermitian(B)*inverse(A)*B);
136  // Structure exploiting answer, approximately half the work and inplace on B.
137  Precision one(1);
138  int n_plus = solver.inertia().first;
139  int n_minus = solver.inertia().second;
140  LowerMatrix<Number> S2(J);
141  solver.solveL(B,'L','N');
142  herk_inplace(S2,B.top(n_plus),'H',one,one);
143  herk_inplace(S2,B.bottom(n_minus),'H',-one,one);
144  // Compare/report.
145  Precision error = frobenius(S1-S2) / frobenius(S1);
146  myra::out() << " |B'*inv(A)*B - (L\\B)'*I*(L\\B)| = " << error << std::endl;
147  REQUIRE(error < tolerance);
148  }
149  // Generate unsymmetric schur complement.
150  {
151  // ZLTDHSolver can also do unsymmetric schur complements, S = Ch*inv(A)*B.
152  // Not especially efficient, but still tricky so it's demonstrated here.
153  // (In the case of sparse A/B/C, this approach is dramatically more efficient)
154  // Reference answer.
155  auto B = Matrix<Number>::random(I,J);
156  auto C = Matrix<Number>::random(I,J+3);
157  Matrix<Number> S1 = (hermitian(C)*inverse(A)*B);
158  // Structure exploiting answer.
159  Number one(1);
160  int n_plus = solver.inertia().first;
161  int n_minus = solver.inertia().second;
162  Matrix<Number> S2(J+3,J);
163  solver.solveL(B,'L','N');
164  solver.solveL(C,'L','N');
165  gemm_inplace(S2,C.top(n_plus),'H',B.top(n_plus),'N',one,one);
166  gemm_inplace(S2,C.bottom(n_minus),'H',B.bottom(n_minus),'N',-one,one);
167  // Compare/report.
168  Precision error = frobenius(S1-S2) / frobenius(S1);
169  myra::out() << " |C'*inv(A)*B - (L\\C)'*I*(L\\B)| = " << error << std::endl;
170  REQUIRE(error < tolerance);
171  }
172  // Save solver to disk, then solve A*X=B again.
173  {
174  VectorStream inout;
175  solver.write(inout);
176  Solver saved(inout);
177  auto X = Matrix<Number>::random(I,J);
178  auto B = A*X;
179  saved.solve(B,'L','N');
180  Precision error = frobenius(B-X)/frobenius(X);
181  myra::out() << " |inv(A)*B-X| (saved) = " << error << std::endl;
182  REQUIRE(error < tolerance);
183  }
184  }
185 
186 template<class Precision> void test3(Precision tolerance)
187  {
188  typedef std::complex<Precision> Number;
189  myra::out() << typestring<Number>() << std::endl;
190  // Construct A = [0 1 0; 1 0 0; 0 0 1], a problematic example.
191  auto A = Matrix<Number>::zeros(3,3);
192  Precision one(1);
193  Precision two(2);
194  Precision pi = std::acos(-one);
195  Precision phase = two*pi*random<Precision>();
196  Number phasor(std::cos(phase),std::sin(phase));
197  A(1,0) = phasor;
198  A(0,1) = conjugate(phasor);
199  A(2,2) = one;
200  // Restrict A to lower triangle, generate LDL' solver for it.
201  LowerMatrix<Number> L(A);
202  typedef ZLDLHSolver<Precision> Solver;
203  Solver solver(L);
204  // Solve A*X=B
205  auto X = Matrix<Number>::random(3,2);
206  auto B = A*X;
207  solver.solve(B,'L','N');
208  Precision error = frobenius(B-X)/frobenius(X);
209  myra::out() << " |inv(A)*B-X| = " << error << std::endl;
210  REQUIRE(error < tolerance);
211  }
212 
213 } // namespace
214 
215 ADD_TEST("cLDLHSolver","[dense][solver]")
216  { test<float>(43,14,1.0e-3f); }
217 
218 ADD_TEST("zLDLHSolver","[dense][solver]")
219  { test<double>(82,45,1.0e-10); }
220 
221 ADD_TEST("cLDLHSolver3","[dense][solver]")
222  { test3<float>(1.0e-3f); }
223 
224 ADD_TEST("zLDLHSolver3","[dense][solver]")
225  { test3<double>(1.0e-10); }
Returns a conjugated copy of a Matrix or Vector. Or, conjugate one inplace.
Tabulates an IxJ matrix. Allows random access, has column major layout to be compatible with BLAS/LAP...
Definition: bdsqr.h:20
static Matrix< Number > zeros(int I, int J)
Generates a zeros Matrix of specified size.
Definition: Matrix.cpp:357
Routines for hermitian rank-k updates, a specialized form of Matrix*Matrix multiplication.
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
Wraps a std::vector<char>, presents it as both an InputStream and OutputStream. Useful for hygienic u...
Definition: VectorStream.h:22
Range construct for a lower triangular matrix stored in rectangular packed format.
Definition: syntax.dox:1
Specialized container for a lower triangular matrix, O(N^2/2) storage. Used by symmetry exploiting ma...
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.
Factors A = L*D*Lh, where D is a "sign matrix" of the form [I 0; 0 -I]. Presents solve methods...
Definition: ZLDLHSolver.h:30
A stream that serialize/deserializes to std::vector<char> buffer.
General purpose dense matrix container, O(i*j) storage.
Routines for computing singular value decompositions.
Overwrites a LowerMatrix, DiagonalMatrix, or square Matrix with its own inverse. Or, returns it as a copy.
A solver suitable for a complex hermitian Matrix (indefinite is fine). Encapsulates hetrf() ...
Returns a hermitian copy of a Matrix. The inplace version only works on a square operand.
Simplistic random number functions.
Streams that serialize/deserialize to/from files.
Stores a lower triangular matrix in rectangular packed format.
Definition: conjugate.h:22
Variety of routines all for dense Matrix*Matrix multiplies. Delegates to the BLAS.
CMatrixRange< Number > top(int i) const
Returns a MatrixRange over the i topmost rows, this(0:i,:)
Definition: Matrix.cpp:207


Results: [PASS]

std::complex<double>
|inv(A)*B-X| = 3.92676e-17


Go back to Summary of /test programs.