aboutsummaryrefslogtreecommitdiffhomepage
path: root/test/product_trsolve.cpp
blob: 352f5475123bd698c48e6d163373caf7c796437d (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
//
// 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"

#define VERIFY_TRSM(TRI,XB) { \
    (XB).setRandom(); ref = (XB); \
    (TRI).solveInPlace(XB); \
    VERIFY_IS_APPROX((TRI).toDenseMatrix() * (XB), ref); \
  }

#define VERIFY_TRSM_ONTHERIGHT(TRI,XB) { \
    (XB).setRandom(); ref = (XB); \
    (TRI).transpose().template solveInPlace<OnTheRight>(XB.transpose()); \
    VERIFY_IS_APPROX((XB).transpose() * (TRI).transpose().toDenseMatrix(), ref.transpose()); \
  }

template<typename Scalar,int Size, int Cols> void trsolve(int size=Size,int cols=Cols)
{
  typedef typename NumTraits<Scalar>::Real RealScalar;

  Matrix<Scalar,Size,Size,ColMajor> cmLhs(size,size);
  Matrix<Scalar,Size,Size,RowMajor> rmLhs(size,size);

  enum {  colmajor = Size==1 ? RowMajor : ColMajor,
          rowmajor = Cols==1 ? ColMajor : RowMajor };
  Matrix<Scalar,Size,Cols,colmajor> cmRhs(size,cols), ref(size,cols);
  Matrix<Scalar,Size,Cols,rowmajor> rmRhs(size,cols);

  cmLhs.setRandom(); cmLhs *= static_cast<RealScalar>(0.1); cmLhs.diagonal().array() += static_cast<RealScalar>(1);
  rmLhs.setRandom(); rmLhs *= static_cast<RealScalar>(0.1); rmLhs.diagonal().array() += static_cast<RealScalar>(1);

  VERIFY_TRSM(cmLhs.conjugate().template triangularView<Lower>(), cmRhs);
  VERIFY_TRSM(cmLhs            .template triangularView<Upper>(), cmRhs);
  VERIFY_TRSM(cmLhs            .template triangularView<Lower>(), rmRhs);
  VERIFY_TRSM(cmLhs.conjugate().template triangularView<Upper>(), rmRhs);

  VERIFY_TRSM(cmLhs.conjugate().template triangularView<UnitLower>(), cmRhs);
  VERIFY_TRSM(cmLhs            .template triangularView<UnitUpper>(), rmRhs);

  VERIFY_TRSM(rmLhs            .template triangularView<Lower>(), cmRhs);
  VERIFY_TRSM(rmLhs.conjugate().template triangularView<UnitUpper>(), rmRhs);


  VERIFY_TRSM_ONTHERIGHT(cmLhs.conjugate().template triangularView<Lower>(), cmRhs);
  VERIFY_TRSM_ONTHERIGHT(cmLhs            .template triangularView<Upper>(), cmRhs);
  VERIFY_TRSM_ONTHERIGHT(cmLhs            .template triangularView<Lower>(), rmRhs);
  VERIFY_TRSM_ONTHERIGHT(cmLhs.conjugate().template triangularView<Upper>(), rmRhs);

  VERIFY_TRSM_ONTHERIGHT(cmLhs.conjugate().template triangularView<UnitLower>(), cmRhs);
  VERIFY_TRSM_ONTHERIGHT(cmLhs            .template triangularView<UnitUpper>(), rmRhs);

  VERIFY_TRSM_ONTHERIGHT(rmLhs            .template triangularView<Lower>(), cmRhs);
  VERIFY_TRSM_ONTHERIGHT(rmLhs.conjugate().template triangularView<UnitUpper>(), rmRhs);
}

void test_product_trsolve()
{
  for(int i = 0; i < g_repeat ; i++)
  {
    // matrices
    CALL_SUBTEST_1((trsolve<float,Dynamic,Dynamic>(ei_random<int>(1,320),ei_random<int>(1,320))));
    CALL_SUBTEST_2((trsolve<double,Dynamic,Dynamic>(ei_random<int>(1,320),ei_random<int>(1,320))));
    CALL_SUBTEST_3((trsolve<std::complex<float>,Dynamic,Dynamic>(ei_random<int>(1,200),ei_random<int>(1,200))));
    CALL_SUBTEST_4((trsolve<std::complex<double>,Dynamic,Dynamic>(ei_random<int>(1,200),ei_random<int>(1,200))));

    // vectors
    CALL_SUBTEST_5((trsolve<std::complex<double>,Dynamic,1>(ei_random<int>(1,320))));
    CALL_SUBTEST_6((trsolve<float,1,1>()));
    CALL_SUBTEST_7((trsolve<float,1,2>()));
    CALL_SUBTEST_8((trsolve<std::complex<float>,4,1>()));
  }
}