xref: /aosp_15_r20/external/eigen/test/sparse.h (revision bf2c37156dfe67e5dfebd6d394bad8b2ab5804d4)
1*bf2c3715SXin Li // This file is part of Eigen, a lightweight C++ template library
2*bf2c3715SXin Li // for linear algebra.
3*bf2c3715SXin Li //
4*bf2c3715SXin Li // Copyright (C) 2008-2011 Gael Guennebaud <[email protected]>
5*bf2c3715SXin Li //
6*bf2c3715SXin Li // This Source Code Form is subject to the terms of the Mozilla
7*bf2c3715SXin Li // Public License v. 2.0. If a copy of the MPL was not distributed
8*bf2c3715SXin Li // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9*bf2c3715SXin Li 
10*bf2c3715SXin Li #ifndef EIGEN_TESTSPARSE_H
11*bf2c3715SXin Li #define EIGEN_TESTSPARSE_H
12*bf2c3715SXin Li 
13*bf2c3715SXin Li #define EIGEN_YES_I_KNOW_SPARSE_MODULE_IS_NOT_STABLE_YET
14*bf2c3715SXin Li 
15*bf2c3715SXin Li #include "main.h"
16*bf2c3715SXin Li 
17*bf2c3715SXin Li #if EIGEN_HAS_CXX11
18*bf2c3715SXin Li 
19*bf2c3715SXin Li #ifdef min
20*bf2c3715SXin Li #undef min
21*bf2c3715SXin Li #endif
22*bf2c3715SXin Li 
23*bf2c3715SXin Li #ifdef max
24*bf2c3715SXin Li #undef max
25*bf2c3715SXin Li #endif
26*bf2c3715SXin Li 
27*bf2c3715SXin Li #include <unordered_map>
28*bf2c3715SXin Li #define EIGEN_UNORDERED_MAP_SUPPORT
29*bf2c3715SXin Li 
30*bf2c3715SXin Li #endif
31*bf2c3715SXin Li 
32*bf2c3715SXin Li #include <Eigen/Cholesky>
33*bf2c3715SXin Li #include <Eigen/LU>
34*bf2c3715SXin Li #include <Eigen/Sparse>
35*bf2c3715SXin Li 
36*bf2c3715SXin Li enum {
37*bf2c3715SXin Li   ForceNonZeroDiag = 1,
38*bf2c3715SXin Li   MakeLowerTriangular = 2,
39*bf2c3715SXin Li   MakeUpperTriangular = 4,
40*bf2c3715SXin Li   ForceRealDiag = 8
41*bf2c3715SXin Li };
42*bf2c3715SXin Li 
43*bf2c3715SXin Li /* Initializes both a sparse and dense matrix with same random values,
44*bf2c3715SXin Li  * and a ratio of \a density non zero entries.
45*bf2c3715SXin Li  * \param flags is a union of ForceNonZeroDiag, MakeLowerTriangular and MakeUpperTriangular
46*bf2c3715SXin Li  *        allowing to control the shape of the matrix.
47*bf2c3715SXin Li  * \param zeroCoords and nonzeroCoords allows to get the coordinate lists of the non zero,
48*bf2c3715SXin Li  *        and zero coefficients respectively.
49*bf2c3715SXin Li  */
50*bf2c3715SXin Li template<typename Scalar,int Opt1,int Opt2,typename StorageIndex> void
51*bf2c3715SXin Li initSparse(double density,
52*bf2c3715SXin Li            Matrix<Scalar,Dynamic,Dynamic,Opt1>& refMat,
53*bf2c3715SXin Li            SparseMatrix<Scalar,Opt2,StorageIndex>& sparseMat,
54*bf2c3715SXin Li            int flags = 0,
55*bf2c3715SXin Li            std::vector<Matrix<StorageIndex,2,1> >* zeroCoords = 0,
56*bf2c3715SXin Li            std::vector<Matrix<StorageIndex,2,1> >* nonzeroCoords = 0)
57*bf2c3715SXin Li {
58*bf2c3715SXin Li   enum { IsRowMajor = SparseMatrix<Scalar,Opt2,StorageIndex>::IsRowMajor };
59*bf2c3715SXin Li   sparseMat.setZero();
60*bf2c3715SXin Li   //sparseMat.reserve(int(refMat.rows()*refMat.cols()*density));
61*bf2c3715SXin Li   sparseMat.reserve(VectorXi::Constant(IsRowMajor ? refMat.rows() : refMat.cols(), int((1.5*density)*(IsRowMajor?refMat.cols():refMat.rows()))));
62*bf2c3715SXin Li 
63*bf2c3715SXin Li   for(Index j=0; j<sparseMat.outerSize(); j++)
64*bf2c3715SXin Li   {
65*bf2c3715SXin Li     //sparseMat.startVec(j);
66*bf2c3715SXin Li     for(Index i=0; i<sparseMat.innerSize(); i++)
67*bf2c3715SXin Li     {
68*bf2c3715SXin Li       Index ai(i), aj(j);
69*bf2c3715SXin Li       if(IsRowMajor)
70*bf2c3715SXin Li         std::swap(ai,aj);
71*bf2c3715SXin Li       Scalar v = (internal::random<double>(0,1) < density) ? internal::random<Scalar>() : Scalar(0);
72*bf2c3715SXin Li       if ((flags&ForceNonZeroDiag) && (i==j))
73*bf2c3715SXin Li       {
74*bf2c3715SXin Li         // FIXME: the following is too conservative
75*bf2c3715SXin Li         v = internal::random<Scalar>()*Scalar(3.);
76*bf2c3715SXin Li         v = v*v;
77*bf2c3715SXin Li         if(numext::real(v)>0) v += Scalar(5);
78*bf2c3715SXin Li         else                  v -= Scalar(5);
79*bf2c3715SXin Li       }
80*bf2c3715SXin Li       if ((flags & MakeLowerTriangular) && aj>ai)
81*bf2c3715SXin Li         v = Scalar(0);
82*bf2c3715SXin Li       else if ((flags & MakeUpperTriangular) && aj<ai)
83*bf2c3715SXin Li         v = Scalar(0);
84*bf2c3715SXin Li 
85*bf2c3715SXin Li       if ((flags&ForceRealDiag) && (i==j))
86*bf2c3715SXin Li         v = numext::real(v);
87*bf2c3715SXin Li 
88*bf2c3715SXin Li       if (v!=Scalar(0))
89*bf2c3715SXin Li       {
90*bf2c3715SXin Li         //sparseMat.insertBackByOuterInner(j,i) = v;
91*bf2c3715SXin Li         sparseMat.insertByOuterInner(j,i) = v;
92*bf2c3715SXin Li         if (nonzeroCoords)
93*bf2c3715SXin Li           nonzeroCoords->push_back(Matrix<StorageIndex,2,1> (ai,aj));
94*bf2c3715SXin Li       }
95*bf2c3715SXin Li       else if (zeroCoords)
96*bf2c3715SXin Li       {
97*bf2c3715SXin Li         zeroCoords->push_back(Matrix<StorageIndex,2,1> (ai,aj));
98*bf2c3715SXin Li       }
99*bf2c3715SXin Li       refMat(ai,aj) = v;
100*bf2c3715SXin Li     }
101*bf2c3715SXin Li   }
102*bf2c3715SXin Li   //sparseMat.finalize();
103*bf2c3715SXin Li }
104*bf2c3715SXin Li 
105*bf2c3715SXin Li template<typename Scalar,int Opt1,int Opt2,typename Index> void
106*bf2c3715SXin Li initSparse(double density,
107*bf2c3715SXin Li            Matrix<Scalar,Dynamic,Dynamic, Opt1>& refMat,
108*bf2c3715SXin Li            DynamicSparseMatrix<Scalar, Opt2, Index>& sparseMat,
109*bf2c3715SXin Li            int flags = 0,
110*bf2c3715SXin Li            std::vector<Matrix<Index,2,1> >* zeroCoords = 0,
111*bf2c3715SXin Li            std::vector<Matrix<Index,2,1> >* nonzeroCoords = 0)
112*bf2c3715SXin Li {
113*bf2c3715SXin Li   enum { IsRowMajor = DynamicSparseMatrix<Scalar,Opt2,Index>::IsRowMajor };
114*bf2c3715SXin Li   sparseMat.setZero();
115*bf2c3715SXin Li   sparseMat.reserve(int(refMat.rows()*refMat.cols()*density));
116*bf2c3715SXin Li   for(int j=0; j<sparseMat.outerSize(); j++)
117*bf2c3715SXin Li   {
118*bf2c3715SXin Li     sparseMat.startVec(j); // not needed for DynamicSparseMatrix
119*bf2c3715SXin Li     for(int i=0; i<sparseMat.innerSize(); i++)
120*bf2c3715SXin Li     {
121*bf2c3715SXin Li       int ai(i), aj(j);
122*bf2c3715SXin Li       if(IsRowMajor)
123*bf2c3715SXin Li         std::swap(ai,aj);
124*bf2c3715SXin Li       Scalar v = (internal::random<double>(0,1) < density) ? internal::random<Scalar>() : Scalar(0);
125*bf2c3715SXin Li       if ((flags&ForceNonZeroDiag) && (i==j))
126*bf2c3715SXin Li       {
127*bf2c3715SXin Li         v = internal::random<Scalar>()*Scalar(3.);
128*bf2c3715SXin Li         v = v*v + Scalar(5.);
129*bf2c3715SXin Li       }
130*bf2c3715SXin Li       if ((flags & MakeLowerTriangular) && aj>ai)
131*bf2c3715SXin Li         v = Scalar(0);
132*bf2c3715SXin Li       else if ((flags & MakeUpperTriangular) && aj<ai)
133*bf2c3715SXin Li         v = Scalar(0);
134*bf2c3715SXin Li 
135*bf2c3715SXin Li       if ((flags&ForceRealDiag) && (i==j))
136*bf2c3715SXin Li         v = numext::real(v);
137*bf2c3715SXin Li 
138*bf2c3715SXin Li       if (v!=Scalar(0))
139*bf2c3715SXin Li       {
140*bf2c3715SXin Li         sparseMat.insertBackByOuterInner(j,i) = v;
141*bf2c3715SXin Li         if (nonzeroCoords)
142*bf2c3715SXin Li           nonzeroCoords->push_back(Matrix<Index,2,1> (ai,aj));
143*bf2c3715SXin Li       }
144*bf2c3715SXin Li       else if (zeroCoords)
145*bf2c3715SXin Li       {
146*bf2c3715SXin Li         zeroCoords->push_back(Matrix<Index,2,1> (ai,aj));
147*bf2c3715SXin Li       }
148*bf2c3715SXin Li       refMat(ai,aj) = v;
149*bf2c3715SXin Li     }
150*bf2c3715SXin Li   }
151*bf2c3715SXin Li   sparseMat.finalize();
152*bf2c3715SXin Li }
153*bf2c3715SXin Li 
154*bf2c3715SXin Li template<typename Scalar,int Options,typename Index> void
155*bf2c3715SXin Li initSparse(double density,
156*bf2c3715SXin Li            Matrix<Scalar,Dynamic,1>& refVec,
157*bf2c3715SXin Li            SparseVector<Scalar,Options,Index>& sparseVec,
158*bf2c3715SXin Li            std::vector<int>* zeroCoords = 0,
159*bf2c3715SXin Li            std::vector<int>* nonzeroCoords = 0)
160*bf2c3715SXin Li {
161*bf2c3715SXin Li   sparseVec.reserve(int(refVec.size()*density));
162*bf2c3715SXin Li   sparseVec.setZero();
163*bf2c3715SXin Li   for(int i=0; i<refVec.size(); i++)
164*bf2c3715SXin Li   {
165*bf2c3715SXin Li     Scalar v = (internal::random<double>(0,1) < density) ? internal::random<Scalar>() : Scalar(0);
166*bf2c3715SXin Li     if (v!=Scalar(0))
167*bf2c3715SXin Li     {
168*bf2c3715SXin Li       sparseVec.insertBack(i) = v;
169*bf2c3715SXin Li       if (nonzeroCoords)
170*bf2c3715SXin Li         nonzeroCoords->push_back(i);
171*bf2c3715SXin Li     }
172*bf2c3715SXin Li     else if (zeroCoords)
173*bf2c3715SXin Li         zeroCoords->push_back(i);
174*bf2c3715SXin Li     refVec[i] = v;
175*bf2c3715SXin Li   }
176*bf2c3715SXin Li }
177*bf2c3715SXin Li 
178*bf2c3715SXin Li template<typename Scalar,int Options,typename Index> void
179*bf2c3715SXin Li initSparse(double density,
180*bf2c3715SXin Li            Matrix<Scalar,1,Dynamic>& refVec,
181*bf2c3715SXin Li            SparseVector<Scalar,Options,Index>& sparseVec,
182*bf2c3715SXin Li            std::vector<int>* zeroCoords = 0,
183*bf2c3715SXin Li            std::vector<int>* nonzeroCoords = 0)
184*bf2c3715SXin Li {
185*bf2c3715SXin Li   sparseVec.reserve(int(refVec.size()*density));
186*bf2c3715SXin Li   sparseVec.setZero();
187*bf2c3715SXin Li   for(int i=0; i<refVec.size(); i++)
188*bf2c3715SXin Li   {
189*bf2c3715SXin Li     Scalar v = (internal::random<double>(0,1) < density) ? internal::random<Scalar>() : Scalar(0);
190*bf2c3715SXin Li     if (v!=Scalar(0))
191*bf2c3715SXin Li     {
192*bf2c3715SXin Li       sparseVec.insertBack(i) = v;
193*bf2c3715SXin Li       if (nonzeroCoords)
194*bf2c3715SXin Li         nonzeroCoords->push_back(i);
195*bf2c3715SXin Li     }
196*bf2c3715SXin Li     else if (zeroCoords)
197*bf2c3715SXin Li         zeroCoords->push_back(i);
198*bf2c3715SXin Li     refVec[i] = v;
199*bf2c3715SXin Li   }
200*bf2c3715SXin Li }
201*bf2c3715SXin Li 
202*bf2c3715SXin Li 
203*bf2c3715SXin Li #include <unsupported/Eigen/SparseExtra>
204*bf2c3715SXin Li #endif // EIGEN_TESTSPARSE_H
205