aboutsummaryrefslogtreecommitdiffhomepage
path: root/test/bandmatrix.cpp
blob: 2f228cdd239407714bf59a74463b1e4ae01472ec (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
// This file is triangularView 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"

template<typename MatrixType> void bandmatrix(const MatrixType& _m)
{
  typedef typename MatrixType::Index Index;
  typedef typename MatrixType::Scalar Scalar;
  typedef typename NumTraits<Scalar>::Real RealScalar;
  typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrixType;

  Index rows = _m.rows();
  Index cols = _m.cols();
  Index supers = _m.supers();
  Index subs = _m.subs();

  MatrixType m(rows,cols,supers,subs);

  DenseMatrixType dm1(rows,cols);
  dm1.setZero();

  m.diagonal().setConstant(123);
  dm1.diagonal().setConstant(123);
  for (int i=1; i<=m.supers();++i)
  {
    m.diagonal(i).setConstant(static_cast<RealScalar>(i));
    dm1.diagonal(i).setConstant(static_cast<RealScalar>(i));
  }
  for (int i=1; i<=m.subs();++i)
  {
    m.diagonal(-i).setConstant(-static_cast<RealScalar>(i));
    dm1.diagonal(-i).setConstant(-static_cast<RealScalar>(i));
  }
  //std::cerr << m.m_data << "\n\n" << m.toDense() << "\n\n" << dm1 << "\n\n\n\n";
  VERIFY_IS_APPROX(dm1,m.toDenseMatrix());

  for (int i=0; i<cols; ++i)
  {
    m.col(i).setConstant(static_cast<RealScalar>(i+1));
    dm1.col(i).setConstant(static_cast<RealScalar>(i+1));
  }
  Index d = std::min(rows,cols);
  Index a = std::max<Index>(0,cols-d-supers);
  Index b = std::max<Index>(0,rows-d-subs);
  if(a>0) dm1.block(0,d+supers,rows,a).setZero();
  dm1.block(0,supers+1,cols-supers-1-a,cols-supers-1-a).template triangularView<Upper>().setZero();
  dm1.block(subs+1,0,rows-subs-1-b,rows-subs-1-b).template triangularView<Lower>().setZero();
  if(b>0) dm1.block(d+subs,0,b,cols).setZero();
  //std::cerr << m.m_data << "\n\n" << m.toDense() << "\n\n" << dm1 << "\n\n";
  VERIFY_IS_APPROX(dm1,m.toDenseMatrix());

}

void test_bandmatrix()
{
  typedef BandMatrix<float>::Index Index;

  for(int i = 0; i < 10*g_repeat ; i++) {
    Index rows = ei_random<Index>(1,10);
    Index cols = ei_random<Index>(1,10);
    Index sups = ei_random<Index>(0,cols-1);
    Index subs = ei_random<Index>(0,rows-1);
    CALL_SUBTEST(bandmatrix(BandMatrix<float>(rows,cols,sups,subs)) );
  }
}