MyraMath
chemmLower


Source: tests/dense/hemm3.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/hemm.h>
23 
24 // Reporting.
25 #include <tests/myratest.h>
26 
27 using namespace myra;
28 
29 namespace {
30 
31 template<class Precision> void test(int I, int J, Precision tolerance)
32  {
33  typedef std::complex<Precision> Number;
34  myra::out() << typestring<Number>() << std::endl;
35  // Make random hermitian matrix A.
36  auto A = Matrix<Number>::random(I,I);
37  A = A + hermitian(A);
38  Number alpha = random<Number>();
39  Number beta = random<Number>();
40  // Check hemm('L')
41  {
42  auto B = Matrix<Number>::random(I,J);
43  auto C1 = Matrix<Number>::random(I,J);
44  Matrix<Number> C2 = C1;
46  gemm_inplace(C1, A, 'N', B, 'N', alpha, beta);
47  hemm_inplace('L', C2, L, B, alpha, beta);
48  Precision error = frobenius(C1-C2);
49  myra::out() << " |hemm('L')-gemm()| = " << error << std::endl;
50  REQUIRE(error < tolerance);
51  }
52  // Check hemm('R')
53  {
54  auto B = Matrix<Number>::random(J,I);
55  auto C1 = Matrix<Number>::random(J,I);
56  Matrix<Number> C2 = C1;
58  gemm_inplace(C1, B, 'N', A, 'N', alpha, beta);
59  hemm_inplace('R', C2, L, B, alpha, beta);
60  Precision error = frobenius(C1-C2);
61  myra::out() << " |hemm('R')-gemm()| = " << error << std::endl;
62  REQUIRE(error < tolerance);
63  }
64  }
65 
66 } // namespace
67 
68 ADD_TEST("chemmLower","[dense][blas]")
69  { test<float> (53,28, 1.0e-4f); }
70 
71 ADD_TEST("zhemmLower","[dense][blas]")
72  { test<double> (53,28, 1.0e-8); }
73 
Tabulates an IxJ matrix. Allows random access, has column major layout to be compatible with BLAS/LAP...
Definition: bdsqr.h:20
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
Routines for hermitian Matrix * dense Matrix multiplication.
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...
Various utility functions/classes related to scalar Number types.
General purpose dense matrix container, O(i*j) storage.
Returns a hermitian copy of a Matrix. The inplace version only works on a square operand.
Simplistic random number functions.
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.


Results: [PASS]

std::complex<float>
|hemm('L')-gemm()| = 4.53846e-05
|hemm('R')-gemm()| = 4.78509e-05


Go back to Summary of /test programs.