aboutsummaryrefslogtreecommitdiffhomepage
path: root/test/sparse.h
blob: 6cd07fc0aa6fd942aa64623ff66303ad644727c2 (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
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
//
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.

#ifndef EIGEN_TESTSPARSE_H
#define EIGEN_TESTSPARSE_H

#define EIGEN_YES_I_KNOW_SPARSE_MODULE_IS_NOT_STABLE_YET

#include "main.h"

#if EIGEN_HAS_CXX11

#ifdef min
#undef min
#endif

#ifdef max
#undef max
#endif

#include <unordered_map>
#define EIGEN_UNORDERED_MAP_SUPPORT

#endif

#include <Eigen/Cholesky>
#include <Eigen/LU>
#include <Eigen/Sparse>

enum {
  ForceNonZeroDiag = 1,
  MakeLowerTriangular = 2,
  MakeUpperTriangular = 4,
  ForceRealDiag = 8
};

/* Initializes both a sparse and dense matrix with same random values,
 * and a ratio of \a density non zero entries.
 * \param flags is a union of ForceNonZeroDiag, MakeLowerTriangular and MakeUpperTriangular
 *        allowing to control the shape of the matrix.
 * \param zeroCoords and nonzeroCoords allows to get the coordinate lists of the non zero,
 *        and zero coefficients respectively.
 */
template<typename Scalar,int Opt1,int Opt2,typename StorageIndex> void
initSparse(double density,
           Matrix<Scalar,Dynamic,Dynamic,Opt1>& refMat,
           SparseMatrix<Scalar,Opt2,StorageIndex>& sparseMat,
           int flags = 0,
           std::vector<Matrix<StorageIndex,2,1> >* zeroCoords = 0,
           std::vector<Matrix<StorageIndex,2,1> >* nonzeroCoords = 0)
{
  enum { IsRowMajor = SparseMatrix<Scalar,Opt2,StorageIndex>::IsRowMajor };
  sparseMat.setZero();
  //sparseMat.reserve(int(refMat.rows()*refMat.cols()*density));
  sparseMat.reserve(VectorXi::Constant(IsRowMajor ? refMat.rows() : refMat.cols(), int((1.5*density)*(IsRowMajor?refMat.cols():refMat.rows()))));
  
  for(Index j=0; j<sparseMat.outerSize(); j++)
  {
    //sparseMat.startVec(j);
    for(Index i=0; i<sparseMat.innerSize(); i++)
    {
      Index ai(i), aj(j);
      if(IsRowMajor)
        std::swap(ai,aj);
      Scalar v = (internal::random<double>(0,1) < density) ? internal::random<Scalar>() : Scalar(0);
      if ((flags&ForceNonZeroDiag) && (i==j))
      {
        // FIXME: the following is too conservative
        v = internal::random<Scalar>()*Scalar(3.);
        v = v*v;
        if(numext::real(v)>0) v += Scalar(5);
        else                  v -= Scalar(5);
      }
      if ((flags & MakeLowerTriangular) && aj>ai)
        v = Scalar(0);
      else if ((flags & MakeUpperTriangular) && aj<ai)
        v = Scalar(0);

      if ((flags&ForceRealDiag) && (i==j))
        v = numext::real(v);

      if (v!=Scalar(0))
      {
        //sparseMat.insertBackByOuterInner(j,i) = v;
        sparseMat.insertByOuterInner(j,i) = v;
        if (nonzeroCoords)
          nonzeroCoords->push_back(Matrix<StorageIndex,2,1> (ai,aj));
      }
      else if (zeroCoords)
      {
        zeroCoords->push_back(Matrix<StorageIndex,2,1> (ai,aj));
      }
      refMat(ai,aj) = v;
    }
  }
  //sparseMat.finalize();
}

template<typename Scalar,int Opt1,int Opt2,typename Index> void
initSparse(double density,
           Matrix<Scalar,Dynamic,Dynamic, Opt1>& refMat,
           DynamicSparseMatrix<Scalar, Opt2, Index>& sparseMat,
           int flags = 0,
           std::vector<Matrix<Index,2,1> >* zeroCoords = 0,
           std::vector<Matrix<Index,2,1> >* nonzeroCoords = 0)
{
  enum { IsRowMajor = DynamicSparseMatrix<Scalar,Opt2,Index>::IsRowMajor };
  sparseMat.setZero();
  sparseMat.reserve(int(refMat.rows()*refMat.cols()*density));
  for(int j=0; j<sparseMat.outerSize(); j++)
  {
    sparseMat.startVec(j); // not needed for DynamicSparseMatrix
    for(int i=0; i<sparseMat.innerSize(); i++)
    {
      int ai(i), aj(j);
      if(IsRowMajor)
        std::swap(ai,aj);
      Scalar v = (internal::random<double>(0,1) < density) ? internal::random<Scalar>() : Scalar(0);
      if ((flags&ForceNonZeroDiag) && (i==j))
      {
        v = internal::random<Scalar>()*Scalar(3.);
        v = v*v + Scalar(5.);
      }
      if ((flags & MakeLowerTriangular) && aj>ai)
        v = Scalar(0);
      else if ((flags & MakeUpperTriangular) && aj<ai)
        v = Scalar(0);

      if ((flags&ForceRealDiag) && (i==j))
        v = numext::real(v);

      if (v!=Scalar(0))
      {
        sparseMat.insertBackByOuterInner(j,i) = v;
        if (nonzeroCoords)
          nonzeroCoords->push_back(Matrix<Index,2,1> (ai,aj));
      }
      else if (zeroCoords)
      {
        zeroCoords->push_back(Matrix<Index,2,1> (ai,aj));
      }
      refMat(ai,aj) = v;
    }
  }
  sparseMat.finalize();
}

template<typename Scalar,int Options,typename Index> void
initSparse(double density,
           Matrix<Scalar,Dynamic,1>& refVec,
           SparseVector<Scalar,Options,Index>& sparseVec,
           std::vector<int>* zeroCoords = 0,
           std::vector<int>* nonzeroCoords = 0)
{
  sparseVec.reserve(int(refVec.size()*density));
  sparseVec.setZero();
  for(int i=0; i<refVec.size(); i++)
  {
    Scalar v = (internal::random<double>(0,1) < density) ? internal::random<Scalar>() : Scalar(0);
    if (v!=Scalar(0))
    {
      sparseVec.insertBack(i) = v;
      if (nonzeroCoords)
        nonzeroCoords->push_back(i);
    }
    else if (zeroCoords)
        zeroCoords->push_back(i);
    refVec[i] = v;
  }
}

template<typename Scalar,int Options,typename Index> void
initSparse(double density,
           Matrix<Scalar,1,Dynamic>& refVec,
           SparseVector<Scalar,Options,Index>& sparseVec,
           std::vector<int>* zeroCoords = 0,
           std::vector<int>* nonzeroCoords = 0)
{
  sparseVec.reserve(int(refVec.size()*density));
  sparseVec.setZero();
  for(int i=0; i<refVec.size(); i++)
  {
    Scalar v = (internal::random<double>(0,1) < density) ? internal::random<Scalar>() : Scalar(0);
    if (v!=Scalar(0))
    {
      sparseVec.insertBack(i) = v;
      if (nonzeroCoords)
        nonzeroCoords->push_back(i);
    }
    else if (zeroCoords)
        zeroCoords->push_back(i);
    refVec[i] = v;
  }
}


#include <unsupported/Eigen/SparseExtra>
#endif // EIGEN_TESTSPARSE_H