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 // IWYU pragma: private
15 #include "./InternalHeaderCheck.h"
16 
17 namespace Eigen {
18 namespace internal {
19 
30 /* TODO
31  * InnerIterator as for sparsematrix
32  * SuperInnerIterator to iterate through all supernodes
33  * Function for triangular solve
34  */
35 template <typename Scalar_, typename StorageIndex_>
37  public:
38  typedef Scalar_ Scalar;
39  typedef StorageIndex_ StorageIndex;
42 
43  public:
46  IndexVector& rowind_colptr, IndexVector& col_to_sup, IndexVector& sup_to_col) {
47  setInfos(m, n, nzval, nzval_colptr, rowind, rowind_colptr, col_to_sup, sup_to_col);
48  }
49 
57  void setInfos(Index m, Index n, ScalarVector& nzval, IndexVector& nzval_colptr, IndexVector& rowind,
58  IndexVector& rowind_colptr, IndexVector& col_to_sup, IndexVector& sup_to_col) {
59  m_row = m;
60  m_col = n;
61  m_nzval = nzval.data();
62  m_nzval_colptr = nzval_colptr.data();
63  m_rowind = rowind.data();
64  m_rowind_colptr = rowind_colptr.data();
65  m_nsuper = col_to_sup(n);
66  m_col_to_sup = col_to_sup.data();
67  m_sup_to_col = sup_to_col.data();
68  }
69 
73  Index rows() const { return m_row; }
74 
78  Index cols() const { return m_col; }
79 
85  Scalar* valuePtr() { return m_nzval; }
86 
87  const Scalar* valuePtr() const { return m_nzval; }
92 
93  const StorageIndex* colIndexPtr() const { return m_nzval_colptr; }
94 
99 
100  const StorageIndex* rowIndex() const { return m_rowind; }
101 
106 
107  const StorageIndex* rowIndexPtr() const { return m_rowind_colptr; }
108 
113 
114  const StorageIndex* colToSup() const { return m_col_to_sup; }
119 
120  const StorageIndex* supToCol() const { return m_sup_to_col; }
121 
125  Index nsuper() const { return m_nsuper; }
126 
127  class InnerIterator;
128  template <typename Dest>
130  template <bool Conjugate, typename Dest>
132 
133  protected:
134  Index m_row; // Number of rows
135  Index m_col; // Number of columns
136  Index m_nsuper; // Number of supernodes
137  Scalar* m_nzval; // array of nonzero values packed by column
138  StorageIndex* m_nzval_colptr; // nzval_colptr[j] Stores the location in nzval[] which starts column j
139  StorageIndex* m_rowind; // Array of compressed row indices of rectangular supernodes
140  StorageIndex* m_rowind_colptr; // rowind_colptr[j] stores the location in rowind[] which starts column j
141  StorageIndex* m_col_to_sup; // col_to_sup[j] is the supernode number to which column j belongs
142  StorageIndex* m_sup_to_col; // sup_to_col[s] points to the starting column of the s-th supernode
143 
144  private:
145 };
146 
151 template <typename Scalar, typename StorageIndex>
153  public:
155  : m_matrix(mat),
156  m_outer(outer),
157  m_supno(mat.colToSup()[outer]),
158  m_idval(mat.colIndexPtr()[outer]),
159  m_startidval(m_idval),
160  m_endidval(mat.colIndexPtr()[outer + 1]),
161  m_idrow(mat.rowIndexPtr()[mat.supToCol()[mat.colToSup()[outer]]]),
162  m_endidrow(mat.rowIndexPtr()[mat.supToCol()[mat.colToSup()[outer]] + 1]) {}
164  m_idval++;
165  m_idrow++;
166  return *this;
167  }
168  inline Scalar value() const { return m_matrix.valuePtr()[m_idval]; }
169 
170  inline Scalar& valueRef() { return const_cast<Scalar&>(m_matrix.valuePtr()[m_idval]); }
171 
172  inline Index index() const { return m_matrix.rowIndex()[m_idrow]; }
173  inline Index row() const { return index(); }
174  inline Index col() const { return m_outer; }
175 
176  inline Index supIndex() const { return m_supno; }
177 
178  inline operator bool() const {
179  return ((m_idval < m_endidval) && (m_idval >= m_startidval) && (m_idrow < m_endidrow));
180  }
181 
182  protected:
183  const MappedSuperNodalMatrix& m_matrix; // Supernodal lower triangular matrix
184  const Index m_outer; // Current column
185  const Index m_supno; // Current SuperNode number
186  Index m_idval; // Index to browse the values in the current column
187  const Index m_startidval; // Start of the column value
188  const Index m_endidval; // End of the column value
189  Index m_idrow; // Index to browse the row indices
190  Index m_endidrow; // End index of row indices of the current column
191 };
192 
197 template <typename Scalar, typename Index_>
198 template <typename Dest>
200  /* Explicit type conversion as the Index type of MatrixBase<Dest> may be wider than Index */
201  // eigen_assert(X.rows() <= NumTraits<Index>::highest());
202  // eigen_assert(X.cols() <= NumTraits<Index>::highest());
203  Index n = int(X.rows());
204  Index nrhs = Index(X.cols());
205  const Scalar* Lval = valuePtr(); // Nonzero values
207  work.setZero();
208  for (Index k = 0; k <= nsuper(); k++) {
209  Index fsupc = supToCol()[k]; // First column of the current supernode
210  Index istart = rowIndexPtr()[fsupc]; // Pointer index to the subscript of the current column
211  Index nsupr = rowIndexPtr()[fsupc + 1] - istart; // Number of rows in the current supernode
212  Index nsupc = supToCol()[k + 1] - fsupc; // Number of columns in the current supernode
213  Index nrow = nsupr - nsupc; // Number of rows in the non-diagonal part of the supernode
214  Index irow; // Current index row
215 
216  if (nsupc == 1) {
217  for (Index j = 0; j < nrhs; j++) {
218  InnerIterator it(*this, fsupc);
219  ++it; // Skip the diagonal element
220  for (; it; ++it) {
221  irow = it.row();
222  X(irow, j) -= X(fsupc, j) * it.value();
223  }
224  }
225  } else {
226  // The supernode has more than one column
227  Index luptr = colIndexPtr()[fsupc];
228  Index lda = colIndexPtr()[fsupc + 1] - luptr;
229 
230  // Triangular solve
231  Map<const Matrix<Scalar, Dynamic, Dynamic, ColMajor>, 0, OuterStride<> > A(&(Lval[luptr]), nsupc, nsupc,
232  OuterStride<>(lda));
233  typename Dest::RowsBlockXpr U = X.derived().middleRows(fsupc, nsupc);
234  U = A.template triangularView<UnitLower>().solve(U);
235  // Matrix-vector product
236  new (&A) Map<const Matrix<Scalar, Dynamic, Dynamic, ColMajor>, 0, OuterStride<> >(&(Lval[luptr + nsupc]), nrow,
237  nsupc, OuterStride<>(lda));
238  work.topRows(nrow).noalias() = A * U;
239 
240  // Begin Scatter
241  for (Index j = 0; j < nrhs; j++) {
242  Index iptr = istart + nsupc;
243  for (Index i = 0; i < nrow; i++) {
244  irow = rowIndex()[iptr];
245  X(irow, j) -= work(i, j); // Scatter operation
246  work(i, j) = Scalar(0);
247  iptr++;
248  }
249  }
250  }
251  }
252 }
253 
254 template <typename Scalar, typename Index_>
255 template <bool Conjugate, typename Dest>
257  using numext::conj;
258  Index n = int(X.rows());
259  Index nrhs = Index(X.cols());
260  const Scalar* Lval = valuePtr(); // Nonzero values
262  work.setZero();
263  for (Index k = nsuper(); k >= 0; k--) {
264  Index fsupc = supToCol()[k]; // First column of the current supernode
265  Index istart = rowIndexPtr()[fsupc]; // Pointer index to the subscript of the current column
266  Index nsupr = rowIndexPtr()[fsupc + 1] - istart; // Number of rows in the current supernode
267  Index nsupc = supToCol()[k + 1] - fsupc; // Number of columns in the current supernode
268  Index nrow = nsupr - nsupc; // Number of rows in the non-diagonal part of the supernode
269  Index irow; // Current index row
270 
271  if (nsupc == 1) {
272  for (Index j = 0; j < nrhs; j++) {
273  InnerIterator it(*this, fsupc);
274  ++it; // Skip the diagonal element
275  for (; it; ++it) {
276  irow = it.row();
277  X(fsupc, j) -= X(irow, j) * (Conjugate ? conj(it.value()) : it.value());
278  }
279  }
280  } else {
281  // The supernode has more than one column
282  Index luptr = colIndexPtr()[fsupc];
283  Index lda = colIndexPtr()[fsupc + 1] - luptr;
284 
285  // Begin Gather
286  for (Index j = 0; j < nrhs; j++) {
287  Index iptr = istart + nsupc;
288  for (Index i = 0; i < nrow; i++) {
289  irow = rowIndex()[iptr];
290  work.topRows(nrow)(i, j) = X(irow, j); // Gather operation
291  iptr++;
292  }
293  }
294 
295  // Matrix-vector product with transposed submatrix
296  Map<const Matrix<Scalar, Dynamic, Dynamic, ColMajor>, 0, OuterStride<> > A(&(Lval[luptr + nsupc]), nrow, nsupc,
297  OuterStride<>(lda));
298  typename Dest::RowsBlockXpr U = X.derived().middleRows(fsupc, nsupc);
299  if (Conjugate)
300  U = U - A.adjoint() * work.topRows(nrow);
301  else
302  U = U - A.transpose() * work.topRows(nrow);
303 
304  // Triangular solve (of transposed diagonal block)
305  new (&A) Map<const Matrix<Scalar, Dynamic, Dynamic, ColMajor>, 0, OuterStride<> >(&(Lval[luptr]), nsupc, nsupc,
306  OuterStride<>(lda));
307  if (Conjugate)
308  U = A.adjoint().template triangularView<UnitUpper>().solve(U);
309  else
310  U = A.transpose().template triangularView<UnitUpper>().solve(U);
311  }
312  }
313 }
314 
315 } // end namespace internal
316 
317 } // end namespace Eigen
318 
319 #endif // EIGEN_SPARSELU_MATRIX_H
AnnoyingScalar conj(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:133
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
SCALAR Scalar
Definition: bench_gemm.cpp:45
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
Definition: ForwardDeclarations.h:102
An InnerIterator allows to loop over the element of any matrix expression.
Definition: CoreIterators.h:37
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:96
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
Convenience specialization of Stride to specify only an outer stride See class Map for some examples.
Definition: Stride.h:104
constexpr EIGEN_DEVICE_FUNC const Scalar * data() const
Definition: PlainObjectBase.h:273
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
Definition: CwiseNullaryOp.h:569
InnerIterator class to iterate over nonzero values of the current column in the supernodal matrix L.
Definition: SparseLU_SupernodalMatrix.h:152
const MappedSuperNodalMatrix & m_matrix
Definition: SparseLU_SupernodalMatrix.h:183
Index m_idrow
Definition: SparseLU_SupernodalMatrix.h:189
Index row() const
Definition: SparseLU_SupernodalMatrix.h:173
Scalar value() const
Definition: SparseLU_SupernodalMatrix.h:168
const Index m_outer
Definition: SparseLU_SupernodalMatrix.h:184
InnerIterator & operator++()
Definition: SparseLU_SupernodalMatrix.h:163
Index index() const
Definition: SparseLU_SupernodalMatrix.h:172
Scalar & valueRef()
Definition: SparseLU_SupernodalMatrix.h:170
Index m_idval
Definition: SparseLU_SupernodalMatrix.h:186
Index supIndex() const
Definition: SparseLU_SupernodalMatrix.h:176
Index m_endidrow
Definition: SparseLU_SupernodalMatrix.h:190
const Index m_startidval
Definition: SparseLU_SupernodalMatrix.h:187
const Index m_supno
Definition: SparseLU_SupernodalMatrix.h:185
Index col() const
Definition: SparseLU_SupernodalMatrix.h:174
const Index m_endidval
Definition: SparseLU_SupernodalMatrix.h:188
InnerIterator(const MappedSuperNodalMatrix &mat, Index outer)
Definition: SparseLU_SupernodalMatrix.h:154
a class to manipulate the L supernodal factor from the SparseLU factorization
Definition: SparseLU_SupernodalMatrix.h:36
Matrix< Scalar, Dynamic, 1 > ScalarVector
Definition: SparseLU_SupernodalMatrix.h:41
StorageIndex * m_col_to_sup
Definition: SparseLU_SupernodalMatrix.h:141
StorageIndex * colIndexPtr()
Definition: SparseLU_SupernodalMatrix.h:91
StorageIndex * rowIndex()
Definition: SparseLU_SupernodalMatrix.h:98
const StorageIndex * colIndexPtr() const
Definition: SparseLU_SupernodalMatrix.h:93
Scalar * m_nzval
Definition: SparseLU_SupernodalMatrix.h:137
StorageIndex * supToCol()
Definition: SparseLU_SupernodalMatrix.h:118
void setInfos(Index m, Index n, ScalarVector &nzval, IndexVector &nzval_colptr, IndexVector &rowind, IndexVector &rowind_colptr, IndexVector &col_to_sup, IndexVector &sup_to_col)
Definition: SparseLU_SupernodalMatrix.h:57
Scalar_ Scalar
Definition: SparseLU_SupernodalMatrix.h:38
const StorageIndex * colToSup() const
Definition: SparseLU_SupernodalMatrix.h:114
StorageIndex * colToSup()
Definition: SparseLU_SupernodalMatrix.h:112
const StorageIndex * rowIndexPtr() const
Definition: SparseLU_SupernodalMatrix.h:107
Index nsuper() const
Definition: SparseLU_SupernodalMatrix.h:125
StorageIndex * m_sup_to_col
Definition: SparseLU_SupernodalMatrix.h:142
MappedSuperNodalMatrix(Index m, Index n, ScalarVector &nzval, IndexVector &nzval_colptr, IndexVector &rowind, IndexVector &rowind_colptr, IndexVector &col_to_sup, IndexVector &sup_to_col)
Definition: SparseLU_SupernodalMatrix.h:45
Index m_row
Definition: SparseLU_SupernodalMatrix.h:134
const StorageIndex * supToCol() const
Definition: SparseLU_SupernodalMatrix.h:120
StorageIndex * m_rowind_colptr
Definition: SparseLU_SupernodalMatrix.h:140
void solveTransposedInPlace(MatrixBase< Dest > &X) const
Definition: SparseLU_SupernodalMatrix.h:256
Index m_col
Definition: SparseLU_SupernodalMatrix.h:135
const Scalar * valuePtr() const
Definition: SparseLU_SupernodalMatrix.h:87
Scalar * valuePtr()
Definition: SparseLU_SupernodalMatrix.h:85
MappedSuperNodalMatrix()
Definition: SparseLU_SupernodalMatrix.h:44
const StorageIndex * rowIndex() const
Definition: SparseLU_SupernodalMatrix.h:100
StorageIndex * rowIndexPtr()
Definition: SparseLU_SupernodalMatrix.h:105
void solveInPlace(MatrixBase< Dest > &X) const
Solve with the supernode triangular matrix.
Definition: SparseLU_SupernodalMatrix.h:199
Matrix< StorageIndex, Dynamic, 1 > IndexVector
Definition: SparseLU_SupernodalMatrix.h:40
StorageIndex * m_nzval_colptr
Definition: SparseLU_SupernodalMatrix.h:138
StorageIndex_ StorageIndex
Definition: SparseLU_SupernodalMatrix.h:39
Index rows() const
Definition: SparseLU_SupernodalMatrix.h:73
StorageIndex * m_rowind
Definition: SparseLU_SupernodalMatrix.h:139
Index m_nsuper
Definition: SparseLU_SupernodalMatrix.h:136
Index cols() const
Definition: SparseLU_SupernodalMatrix.h:78
~MappedSuperNodalMatrix()
Definition: SparseLU_SupernodalMatrix.h:50
#define X
Definition: icosphere.cpp:20
return int(ret)+1
const char const int const RealScalar const RealScalar const int * lda
Definition: level2_cplx_impl.h:20
int * m
Definition: level2_cplx_impl.h:294
char char char int int * k
Definition: level2_impl.h:374
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
const AutoDiffScalar< DerType > & conj(const AutoDiffScalar< DerType > &x)
Definition: AutoDiffScalar.h:482
double U
Swimming speed.
Definition: two_d_variable_diff_adapt.cc:53
Definition: Eigen_Colamd.h:49
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2