aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/SparseLU/SparseLU_Matrix.h
blob: 9f2dcaa5656f369481c3571b767c30e3f219a627 (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
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
// Copyright (C) 2012 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/>.

#ifndef EIGEN_SPARSELU_MATRIX_H
#define EIGEN_SPARSELU_MATRIX_H

/** \ingroup SparseLU_Module
 * \brief a class to manipulate the L supernodal factor from the SparseLU factorization
 * 
 * This class  contain the data to easily store 
 * and manipulate the supernodes during the factorization and solution phase of Sparse LU. 
 * Only the lower triangular matrix has supernodes.
 * 
 * NOTE : This class corresponds to the SCformat structure in SuperLU
 * 
 */
/* TO DO
 * InnerIterator as for sparsematrix 
 * SuperInnerIterator to iterate through all supernodes 
 * Function for triangular solve
 */
template <typename _Scalar, typename _Index>
class SuperNodalMatrix
{
  public:
    typedef _Scalar Scalar; 
    typedef _Index Index;
    typedef Matrix<Index,Dynamic,1> IndexVector; 
    typedef Matrix<Scalar,Dynamic,1> ScalarVector;
  public:
    SuperNodalMatrix()
    {
      
    }
    SuperNodalMatrix(int m, int n,  ScalarVector& nzval, IndexVector& nzval_colptr, IndexVector& rowind, 
             IndexVector& rowind_colptr, IndexVector& col_to_sup, IndexVector& sup_to_col )
    {
      setInfos(m, n, nzval, nzval_colptr, rowind, rowind_colptr, col_to_sup, sup_to_col);
    }
    
    ~SuperNodalMatrix()
    {
      
    }
    /**
     * Set appropriate pointers for the lower triangular supernodal matrix
     * These infos are available at the end of the numerical factorization
     * FIXME This class will be modified such that it can be use in the course 
     * of the factorization.
     */
    void setInfos(int m, int n, ScalarVector& nzval, IndexVector& nzval_colptr, IndexVector& rowind, 
             IndexVector& rowind_colptr, IndexVector& col_to_sup, IndexVector& sup_to_col )
    {
      m_row = m;
      m_col = n; 
      m_nzval = nzval.data(); 
      m_nzval_colptr = nzval_colptr.data(); 
      m_rowind = rowind.data(); 
      m_rowind_colptr = rowind_colptr.data(); 
      m_nsuper = col_to_sup(n); 
      m_col_to_sup = col_to_sup.data(); 
      m_sup_to_col = sup_to_col.data(); 
      
    }
    
    /**
     * Number of rows
     */
    int rows()
    {
      return m_row;
    }
    
    /**
     * Number of columns
     */
    int cols()
    {
      return m_col;
    }
    
    /**
     * Return the array of nonzero values packed by column
     * 
     * The size is nnz
     */
    Scalar* valuePtr()
    {
      return m_nzval; 
    }
    
    const Scalar* valuePtr() const 
    {
      return m_nzval; 
    }
    /**
     * Return the pointers to the beginning of each column in \ref valuePtr()
     */
    Index* colIndexPtr()
    {
      return m_nzval_colptr; 
    }
    
    const Index* colIndexPtr() const
    {
      return m_nzval_colptr; 
    }
    
    /**
     * Return the array of compressed row indices of all supernodes
     */
    Index* rowIndex()
    {
      return m_rowind; 
    }
    
    const Index* rowIndex() const
    {
      return m_rowind; 
    }
    
    /**
     * Return the location in \em rowvaluePtr() which starts each column
     */
    Index* rowIndexPtr()
    {
      return m_rowind_colptr; 
    }
    
    const Index* rowIndexPtr() const 
    {
      return m_rowind_colptr; 
    }
    
    /** 
     * Return the array of column-to-supernode mapping 
     */
    Index* colToSup()
    {
      return m_col_to_sup;       
    }
    
    const Index* colToSup() const
    {
      return m_col_to_sup;       
    }
    /**
     * Return the array of supernode-to-column mapping
     */
    Index* supToCol()
    {
      return m_sup_to_col;
    }
    
    const Index* supToCol() const 
    {
      return m_sup_to_col;
    }
    
    /**
     * Return the number of supernodes
     */
    int nsuper() const 
    {
      return m_nsuper; 
    }
    class InnerIterator; 
    class SuperNodeIterator;
    
  protected:
    Index m_row; // Number of rows
    Index m_col; // Number of columns 
    Index m_nsuper; // Number of supernodes 
    Scalar* m_nzval; //array of nonzero values packed by column
    Index* m_nzval_colptr; //nzval_colptr[j] Stores the location in nzval[] which starts column j 
    Index* m_rowind; // Array of compressed row indices of rectangular supernodes
    Index* m_rowind_colptr; //rowind_colptr[j] stores the location in rowind[] which starts column j
    Index* m_col_to_sup; // col_to_sup[j] is the supernode number to which column j belongs
    Index* m_sup_to_col; //sup_to_col[s] points to the starting column of the s-th supernode
    
  private :
};

/**
  * \brief InnerIterator class to iterate over nonzero values in the triangular supernodal matrix
  * 
  */
template<typename Scalar, typename Index>
class SuperNodalMatrix<Scalar,Index>::InnerIterator
{
  public:
     InnerIterator(const SuperNodalMatrix& mat, Index outer)
      : m_matrix(mat),
        m_outer(outer), 
        m_idval(mat.colIndexPtr()[outer]),
        m_startval(m_idval),
        m_endval(mat.colIndexPtr()[outer+1]),
        m_idrow(mat.rowIndexPtr()[outer]),
        m_startidrow(m_idrow),
        m_endidrow(mat.rowIndexPtr()[outer+1])
    {}
    inline InnerIterator& operator++()
    { 
      m_idval++; 
      m_idrow++ ;
      return *this; 
    }
    inline Scalar value() const { return m_matrix.valuePtr()[m_idval]; }
    
    inline Scalar& valueRef() { return const_cast<Scalar&>(m_matrix.valuePtr()[m_idval]); }
    
    inline Index index() const { return m_matrix.rowIndex()[m_idrow]; }
    inline Index row() const { return index(); }
    inline Index col() const { return m_outer; }
    
    inline Index supIndex() const { return m_matrix.colToSup()[m_outer]; }
    
    inline operator bool() const 
    { 
      return ( (m_idval < m_endval) && (m_idval > m_startval) && 
                (m_idrow < m_endidrow) && (m_idrow > m_startidrow) ); 
    }
    
  protected:
        const SuperNodalMatrix& m_matrix; // Supernodal lower triangular matrix 
        const Index m_outer; // Current column 
        Index m_idval; //Index to browse the values in the current column
        const Index m_startval; // Start of the column value 
        const Index m_endval; // End of the column value 
        Index m_idrow;  //Index to browse the row indices 
        const Index m_startidrow; // Start of the row indices of the current column value
        const Index m_endidrow; // End of the row indices of the current column value
};

/**
 * \brief Iterator class to iterate over Supernodes in the triangular supernodal matrix
 * 
 * The final goal is to use this class when dealing with supernodes during numerical factorization
 */
template<typename Scalar, typename Index>
class SuperNodalMatrix<Scalar,Index>::SuperNodeIterator
{
  public: 
    SuperNodeIterator(const SuperNodalMatrix& mat)
    {
      
    }
    SuperNodeIterator(const SuperNodalMatrix& mat, Index supno)
    {
      
    }
    
    /*
     * Available Methods : 
     * Browse all supernodes (operator ++ )
     * Number of supernodes 
     * Columns of the current supernode
     * triangular matrix of the current supernode 
     * rectangular part of the current supernode
     */
  protected:
    const SuperNodalMatrix& m_matrix; // Supernodal lower triangular matrix 
    Index m_idsup;  // Index to browse all supernodes
    const Index m_nsuper; // Number of all supernodes
    Index m_startidsup; 
    Index m_endidsup; 
  
}; 
#endif