SparseLU_SupernodalMatrix.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
5 // Copyright (C) 2012 Gael Guennebaud <gael.guennebaud@inria.fr>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_SPARSELU_SUPERNODAL_MATRIX_H
12 #define EIGEN_SPARSELU_SUPERNODAL_MATRIX_H
13 
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 namespace internal {
18 
29 /* TODO
30  * InnerIterator as for sparsematrix
31  * SuperInnerIterator to iterate through all supernodes
32  * Function for triangular solve
33  */
34 template <typename Scalar_, typename StorageIndex_>
35 class MappedSuperNodalMatrix
36 {
37  public:
38  typedef Scalar_ Scalar;
39  typedef StorageIndex_ StorageIndex;
40  typedef Matrix<StorageIndex,Dynamic,1> IndexVector;
41  typedef Matrix<Scalar,Dynamic,1> ScalarVector;
42  public:
43  MappedSuperNodalMatrix()
44  {
45 
46  }
47  MappedSuperNodalMatrix(Index m, Index n, ScalarVector& nzval, IndexVector& nzval_colptr, IndexVector& rowind,
48  IndexVector& rowind_colptr, IndexVector& col_to_sup, IndexVector& sup_to_col )
49  {
50  setInfos(m, n, nzval, nzval_colptr, rowind, rowind_colptr, col_to_sup, sup_to_col);
51  }
52 
53  ~MappedSuperNodalMatrix()
54  {
55 
56  }
63  void setInfos(Index m, Index n, ScalarVector& nzval, IndexVector& nzval_colptr, IndexVector& rowind,
64  IndexVector& rowind_colptr, IndexVector& col_to_sup, IndexVector& sup_to_col )
65  {
66  m_row = m;
67  m_col = n;
68  m_nzval = nzval.data();
69  m_nzval_colptr = nzval_colptr.data();
70  m_rowind = rowind.data();
71  m_rowind_colptr = rowind_colptr.data();
72  m_nsuper = col_to_sup(n);
73  m_col_to_sup = col_to_sup.data();
74  m_sup_to_col = sup_to_col.data();
75  }
76 
80  Index rows() const { return m_row; }
81 
85  Index cols() const { return m_col; }
86 
92  Scalar* valuePtr() { return m_nzval; }
93 
94  const Scalar* valuePtr() const
95  {
96  return m_nzval;
97  }
101  StorageIndex* colIndexPtr()
102  {
103  return m_nzval_colptr;
104  }
105 
106  const StorageIndex* colIndexPtr() const
107  {
108  return m_nzval_colptr;
109  }
110 
114  StorageIndex* rowIndex() { return m_rowind; }
115 
116  const StorageIndex* rowIndex() const
117  {
118  return m_rowind;
119  }
120 
124  StorageIndex* rowIndexPtr() { return m_rowind_colptr; }
125 
126  const StorageIndex* rowIndexPtr() const
127  {
128  return m_rowind_colptr;
129  }
130 
134  StorageIndex* colToSup() { return m_col_to_sup; }
135 
136  const StorageIndex* colToSup() const
137  {
138  return m_col_to_sup;
139  }
143  StorageIndex* supToCol() { return m_sup_to_col; }
144 
145  const StorageIndex* supToCol() const
146  {
147  return m_sup_to_col;
148  }
149 
153  Index nsuper() const
154  {
155  return m_nsuper;
156  }
157 
158  class InnerIterator;
159  template<typename Dest>
160  void solveInPlace( MatrixBase<Dest>&X) const;
161  template<bool Conjugate, typename Dest>
162  void solveTransposedInPlace( MatrixBase<Dest>&X) const;
163 
164 
165 
166 
167 
168  protected:
169  Index m_row; // Number of rows
170  Index m_col; // Number of columns
171  Index m_nsuper; // Number of supernodes
172  Scalar* m_nzval; //array of nonzero values packed by column
173  StorageIndex* m_nzval_colptr; //nzval_colptr[j] Stores the location in nzval[] which starts column j
174  StorageIndex* m_rowind; // Array of compressed row indices of rectangular supernodes
175  StorageIndex* m_rowind_colptr; //rowind_colptr[j] stores the location in rowind[] which starts column j
176  StorageIndex* m_col_to_sup; // col_to_sup[j] is the supernode number to which column j belongs
177  StorageIndex* m_sup_to_col; //sup_to_col[s] points to the starting column of the s-th supernode
178 
179  private :
180 };
181 
186 template<typename Scalar, typename StorageIndex>
187 class MappedSuperNodalMatrix<Scalar,StorageIndex>::InnerIterator
188 {
189  public:
190  InnerIterator(const MappedSuperNodalMatrix& mat, Index outer)
191  : m_matrix(mat),
192  m_outer(outer),
193  m_supno(mat.colToSup()[outer]),
194  m_idval(mat.colIndexPtr()[outer]),
195  m_startidval(m_idval),
196  m_endidval(mat.colIndexPtr()[outer+1]),
197  m_idrow(mat.rowIndexPtr()[mat.supToCol()[mat.colToSup()[outer]]]),
198  m_endidrow(mat.rowIndexPtr()[mat.supToCol()[mat.colToSup()[outer]]+1])
199  {}
200  inline InnerIterator& operator++()
201  {
202  m_idval++;
203  m_idrow++;
204  return *this;
205  }
206  inline Scalar value() const { return m_matrix.valuePtr()[m_idval]; }
207 
208  inline Scalar& valueRef() { return const_cast<Scalar&>(m_matrix.valuePtr()[m_idval]); }
209 
210  inline Index index() const { return m_matrix.rowIndex()[m_idrow]; }
211  inline Index row() const { return index(); }
212  inline Index col() const { return m_outer; }
213 
214  inline Index supIndex() const { return m_supno; }
215 
216  inline operator bool() const
217  {
218  return ( (m_idval < m_endidval) && (m_idval >= m_startidval)
219  && (m_idrow < m_endidrow) );
220  }
221 
222  protected:
223  const MappedSuperNodalMatrix& m_matrix; // Supernodal lower triangular matrix
224  const Index m_outer; // Current column
225  const Index m_supno; // Current SuperNode number
226  Index m_idval; // Index to browse the values in the current column
227  const Index m_startidval; // Start of the column value
228  const Index m_endidval; // End of the column value
229  Index m_idrow; // Index to browse the row indices
230  Index m_endidrow; // End index of row indices of the current column
231 };
232 
237 template<typename Scalar, typename Index_>
238 template<typename Dest>
239 void MappedSuperNodalMatrix<Scalar,Index_>::solveInPlace( MatrixBase<Dest>&X) const
240 {
241  /* Explicit type conversion as the Index type of MatrixBase<Dest> may be wider than Index */
242 // eigen_assert(X.rows() <= NumTraits<Index>::highest());
243 // eigen_assert(X.cols() <= NumTraits<Index>::highest());
244  Index n = int(X.rows());
245  Index nrhs = Index(X.cols());
246  const Scalar * Lval = valuePtr(); // Nonzero values
247  Matrix<Scalar,Dynamic,Dest::ColsAtCompileTime, ColMajor> work(n, nrhs); // working vector
248  work.setZero();
249  for (Index k = 0; k <= nsuper(); k ++)
250  {
251  Index fsupc = supToCol()[k]; // First column of the current supernode
252  Index istart = rowIndexPtr()[fsupc]; // Pointer index to the subscript of the current column
253  Index nsupr = rowIndexPtr()[fsupc+1] - istart; // Number of rows in the current supernode
254  Index nsupc = supToCol()[k+1] - fsupc; // Number of columns in the current supernode
255  Index nrow = nsupr - nsupc; // Number of rows in the non-diagonal part of the supernode
256  Index irow; //Current index row
257 
258  if (nsupc == 1 )
259  {
260  for (Index j = 0; j < nrhs; j++)
261  {
262  InnerIterator it(*this, fsupc);
263  ++it; // Skip the diagonal element
264  for (; it; ++it)
265  {
266  irow = it.row();
267  X(irow, j) -= X(fsupc, j) * it.value();
268  }
269  }
270  }
271  else
272  {
273  // The supernode has more than one column
274  Index luptr = colIndexPtr()[fsupc];
275  Index lda = colIndexPtr()[fsupc+1] - luptr;
276 
277  // Triangular solve
278  Map<const Matrix<Scalar,Dynamic,Dynamic, ColMajor>, 0, OuterStride<> > A( &(Lval[luptr]), nsupc, nsupc, OuterStride<>(lda) );
279  typename Dest::RowsBlockXpr U = X.derived().middleRows(fsupc, nsupc);
280  U = A.template triangularView<UnitLower>().solve(U);
281  // Matrix-vector product
282  new (&A) Map<const Matrix<Scalar,Dynamic,Dynamic, ColMajor>, 0, OuterStride<> > ( &(Lval[luptr+nsupc]), nrow, nsupc, OuterStride<>(lda) );
283  work.topRows(nrow).noalias() = A * U;
284 
285  //Begin Scatter
286  for (Index j = 0; j < nrhs; j++)
287  {
288  Index iptr = istart + nsupc;
289  for (Index i = 0; i < nrow; i++)
290  {
291  irow = rowIndex()[iptr];
292  X(irow, j) -= work(i, j); // Scatter operation
293  work(i, j) = Scalar(0);
294  iptr++;
295  }
296  }
297  }
298  }
299 }
300 
301 template<typename Scalar, typename Index_>
302 template<bool Conjugate, typename Dest>
303 void MappedSuperNodalMatrix<Scalar,Index_>::solveTransposedInPlace( MatrixBase<Dest>&X) const
304 {
305  using numext::conj;
306  Index n = int(X.rows());
307  Index nrhs = Index(X.cols());
308  const Scalar * Lval = valuePtr(); // Nonzero values
309  Matrix<Scalar,Dynamic,Dest::ColsAtCompileTime, ColMajor> work(n, nrhs); // working vector
310  work.setZero();
311  for (Index k = nsuper(); k >= 0; k--)
312  {
313  Index fsupc = supToCol()[k]; // First column of the current supernode
314  Index istart = rowIndexPtr()[fsupc]; // Pointer index to the subscript of the current column
315  Index nsupr = rowIndexPtr()[fsupc+1] - istart; // Number of rows in the current supernode
316  Index nsupc = supToCol()[k+1] - fsupc; // Number of columns in the current supernode
317  Index nrow = nsupr - nsupc; // Number of rows in the non-diagonal part of the supernode
318  Index irow; //Current index row
319 
320  if (nsupc == 1 )
321  {
322  for (Index j = 0; j < nrhs; j++)
323  {
324  InnerIterator it(*this, fsupc);
325  ++it; // Skip the diagonal element
326  for (; it; ++it)
327  {
328  irow = it.row();
329  X(fsupc,j) -= X(irow, j) * (Conjugate?conj(it.value()):it.value());
330  }
331  }
332  }
333  else
334  {
335  // The supernode has more than one column
336  Index luptr = colIndexPtr()[fsupc];
337  Index lda = colIndexPtr()[fsupc+1] - luptr;
338 
339  //Begin Gather
340  for (Index j = 0; j < nrhs; j++)
341  {
342  Index iptr = istart + nsupc;
343  for (Index i = 0; i < nrow; i++)
344  {
345  irow = rowIndex()[iptr];
346  work.topRows(nrow)(i,j)= X(irow,j); // Gather operation
347  iptr++;
348  }
349  }
350 
351  // Matrix-vector product with transposed submatrix
352  Map<const Matrix<Scalar,Dynamic,Dynamic, ColMajor>, 0, OuterStride<> > A( &(Lval[luptr+nsupc]), nrow, nsupc, OuterStride<>(lda) );
353  typename Dest::RowsBlockXpr U = X.derived().middleRows(fsupc, nsupc);
354  if(Conjugate)
355  U = U - A.adjoint() * work.topRows(nrow);
356  else
357  U = U - A.transpose() * work.topRows(nrow);
358 
359  // Triangular solve (of transposed diagonal block)
360  new (&A) Map<const Matrix<Scalar,Dynamic,Dynamic, ColMajor>, 0, OuterStride<> > ( &(Lval[luptr]), nsupc, nsupc, OuterStride<>(lda) );
361  if(Conjugate)
362  U = A.adjoint().template triangularView<UnitUpper>().solve(U);
363  else
364  U = A.transpose().template triangularView<UnitUpper>().solve(U);
365 
366  }
367 
368  }
369 }
370 
371 
372 } // end namespace internal
373 
374 } // end namespace Eigen
375 
376 #endif // EIGEN_SPARSELU_MATRIX_H
Matrix3f m
int n
RowXpr row(Index i)
This is the const version of row(). *‍/.
ColXpr col(Index i)
This is the const version of col().
MatrixXcf A
bfloat16 operator++(bfloat16 &a)
Definition: BFloat16.h:298
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
std::ptrdiff_t j