Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ > Class Template Reference

a class to manipulate the L supernodal factor from the SparseLU factorization More...

#include <SparseLU_SupernodalMatrix.h>

Classes

class  InnerIterator
 InnerIterator class to iterate over nonzero values of the current column in the supernodal matrix L. More...
 

Public Types

typedef Scalar_ Scalar
 
typedef StorageIndex_ StorageIndex
 
typedef Matrix< StorageIndex, Dynamic, 1 > IndexVector
 
typedef Matrix< Scalar, Dynamic, 1 > ScalarVector
 

Public Member Functions

 MappedSuperNodalMatrix ()
 
 MappedSuperNodalMatrix (Index m, Index n, ScalarVector &nzval, IndexVector &nzval_colptr, IndexVector &rowind, IndexVector &rowind_colptr, IndexVector &col_to_sup, IndexVector &sup_to_col)
 
 ~MappedSuperNodalMatrix ()
 
void setInfos (Index m, Index n, ScalarVector &nzval, IndexVector &nzval_colptr, IndexVector &rowind, IndexVector &rowind_colptr, IndexVector &col_to_sup, IndexVector &sup_to_col)
 
Index rows () const
 
Index cols () const
 
ScalarvaluePtr ()
 
const ScalarvaluePtr () const
 
StorageIndexcolIndexPtr ()
 
const StorageIndexcolIndexPtr () const
 
StorageIndexrowIndex ()
 
const StorageIndexrowIndex () const
 
StorageIndexrowIndexPtr ()
 
const StorageIndexrowIndexPtr () const
 
StorageIndexcolToSup ()
 
const StorageIndexcolToSup () const
 
StorageIndexsupToCol ()
 
const StorageIndexsupToCol () const
 
Index nsuper () const
 
template<typename Dest >
void solveInPlace (MatrixBase< Dest > &X) const
 Solve with the supernode triangular matrix. More...
 
template<bool Conjugate, typename Dest >
void solveTransposedInPlace (MatrixBase< Dest > &X) const
 

Protected Attributes

Index m_row
 
Index m_col
 
Index m_nsuper
 
Scalarm_nzval
 
StorageIndexm_nzval_colptr
 
StorageIndexm_rowind
 
StorageIndexm_rowind_colptr
 
StorageIndexm_col_to_sup
 
StorageIndexm_sup_to_col
 

Detailed Description

template<typename Scalar_, typename StorageIndex_>
class Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >

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

Member Typedef Documentation

◆ IndexVector

template<typename Scalar_ , typename StorageIndex_ >
typedef Matrix<StorageIndex, Dynamic, 1> Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::IndexVector

◆ Scalar

template<typename Scalar_ , typename StorageIndex_ >
typedef Scalar_ Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::Scalar

◆ ScalarVector

template<typename Scalar_ , typename StorageIndex_ >
typedef Matrix<Scalar, Dynamic, 1> Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::ScalarVector

◆ StorageIndex

template<typename Scalar_ , typename StorageIndex_ >
typedef StorageIndex_ Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::StorageIndex

Constructor & Destructor Documentation

◆ MappedSuperNodalMatrix() [1/2]

template<typename Scalar_ , typename StorageIndex_ >
Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::MappedSuperNodalMatrix ( )
inline
44 {}

◆ MappedSuperNodalMatrix() [2/2]

template<typename Scalar_ , typename StorageIndex_ >
Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::MappedSuperNodalMatrix ( Index  m,
Index  n,
ScalarVector nzval,
IndexVector nzval_colptr,
IndexVector rowind,
IndexVector rowind_colptr,
IndexVector col_to_sup,
IndexVector sup_to_col 
)
inline
46  {
47  setInfos(m, n, nzval, nzval_colptr, rowind, rowind_colptr, col_to_sup, sup_to_col);
48  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
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
int * m
Definition: level2_cplx_impl.h:294

References m, n, and Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::setInfos().

◆ ~MappedSuperNodalMatrix()

template<typename Scalar_ , typename StorageIndex_ >
Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::~MappedSuperNodalMatrix ( )
inline
50 {}

Member Function Documentation

◆ colIndexPtr() [1/2]

template<typename Scalar_ , typename StorageIndex_ >
StorageIndex* Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::colIndexPtr ( )
inline

Return the pointers to the beginning of each column in valuePtr()

91 { return m_nzval_colptr; }
StorageIndex * m_nzval_colptr
Definition: SparseLU_SupernodalMatrix.h:138

References Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::m_nzval_colptr.

◆ colIndexPtr() [2/2]

template<typename Scalar_ , typename StorageIndex_ >
const StorageIndex* Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::colIndexPtr ( ) const
inline

◆ cols()

◆ colToSup() [1/2]

template<typename Scalar_ , typename StorageIndex_ >
StorageIndex* Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::colToSup ( )
inline

Return the array of column-to-supernode mapping

112 { return m_col_to_sup; }
StorageIndex * m_col_to_sup
Definition: SparseLU_SupernodalMatrix.h:141

References Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::m_col_to_sup.

◆ colToSup() [2/2]

template<typename Scalar_ , typename StorageIndex_ >
const StorageIndex* Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::colToSup ( ) const
inline

◆ nsuper()

template<typename Scalar_ , typename StorageIndex_ >
Index Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::nsuper ( ) const
inline

Return the number of supernodes

125 { return m_nsuper; }
Index m_nsuper
Definition: SparseLU_SupernodalMatrix.h:136

References Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::m_nsuper.

◆ rowIndex() [1/2]

template<typename Scalar_ , typename StorageIndex_ >
StorageIndex* Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::rowIndex ( )
inline

Return the array of compressed row indices of all supernodes

98 { return m_rowind; }
StorageIndex * m_rowind
Definition: SparseLU_SupernodalMatrix.h:139

References Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::m_rowind.

◆ rowIndex() [2/2]

template<typename Scalar_ , typename StorageIndex_ >
const StorageIndex* Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::rowIndex ( ) const
inline

◆ rowIndexPtr() [1/2]

template<typename Scalar_ , typename StorageIndex_ >
StorageIndex* Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::rowIndexPtr ( )
inline

Return the location in rowvaluePtr() which starts each column

105 { return m_rowind_colptr; }
StorageIndex * m_rowind_colptr
Definition: SparseLU_SupernodalMatrix.h:140

References Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::m_rowind_colptr.

◆ rowIndexPtr() [2/2]

template<typename Scalar_ , typename StorageIndex_ >
const StorageIndex* Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::rowIndexPtr ( ) const
inline

◆ rows()

◆ setInfos()

template<typename Scalar_ , typename StorageIndex_ >
void Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::setInfos ( Index  m,
Index  n,
ScalarVector nzval,
IndexVector nzval_colptr,
IndexVector rowind,
IndexVector rowind_colptr,
IndexVector col_to_sup,
IndexVector sup_to_col 
)
inline

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.

58  {
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  }
Scalar * m_nzval
Definition: SparseLU_SupernodalMatrix.h:137
StorageIndex * m_sup_to_col
Definition: SparseLU_SupernodalMatrix.h:142

References Eigen::PlainObjectBase< Derived >::data(), m, Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::m_col, Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::m_col_to_sup, Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::m_nsuper, Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::m_nzval, Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::m_nzval_colptr, Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::m_row, Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::m_rowind, Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::m_rowind_colptr, Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::m_sup_to_col, and n.

Referenced by Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::MappedSuperNodalMatrix().

◆ solveInPlace()

template<typename Scalar , typename Index_ >
template<typename Dest >
void Eigen::internal::MappedSuperNodalMatrix< Scalar, Index_ >::solveInPlace ( MatrixBase< Dest > &  X) const

Solve with the supernode triangular matrix.

199  {
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
206  Matrix<Scalar, Dynamic, Dest::ColsAtCompileTime, ColMajor> work(n, nrhs); // working vector
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 }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
SCALAR Scalar
Definition: bench_gemm.cpp:45
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
StorageIndex * colIndexPtr()
Definition: SparseLU_SupernodalMatrix.h:91
StorageIndex * rowIndex()
Definition: SparseLU_SupernodalMatrix.h:98
StorageIndex * supToCol()
Definition: SparseLU_SupernodalMatrix.h:118
Scalar_ Scalar
Definition: SparseLU_SupernodalMatrix.h:38
Index nsuper() const
Definition: SparseLU_SupernodalMatrix.h:125
Scalar * valuePtr()
Definition: SparseLU_SupernodalMatrix.h:85
StorageIndex * rowIndexPtr()
Definition: SparseLU_SupernodalMatrix.h:105
#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
char char char int int * k
Definition: level2_impl.h:374
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
double U
Swimming speed.
Definition: two_d_variable_diff_adapt.cc:53
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References i, int(), j, k, lda, n, Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::InnerIterator::row(), Eigen::PlainObjectBase< Derived >::setZero(), RachelsAdvectionDiffusion::U, Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::InnerIterator::value(), and X.

◆ solveTransposedInPlace()

template<typename Scalar , typename Index_ >
template<bool Conjugate, typename Dest >
void Eigen::internal::MappedSuperNodalMatrix< Scalar, Index_ >::solveTransposedInPlace ( MatrixBase< Dest > &  X) const
256  {
257  using numext::conj;
258  Index n = int(X.rows());
259  Index nrhs = Index(X.cols());
260  const Scalar* Lval = valuePtr(); // Nonzero values
261  Matrix<Scalar, Dynamic, Dest::ColsAtCompileTime, ColMajor> work(n, nrhs); // working vector
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 }
AnnoyingScalar conj(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:133
squared absolute value
Definition: GlobalFunctions.h:87
const AutoDiffScalar< DerType > & conj(const AutoDiffScalar< DerType > &x)
Definition: AutoDiffScalar.h:482

References conj(), Eigen::conj(), i, int(), j, k, lda, n, Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::InnerIterator::row(), Eigen::PlainObjectBase< Derived >::setZero(), RachelsAdvectionDiffusion::U, Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::InnerIterator::value(), and X.

◆ supToCol() [1/2]

template<typename Scalar_ , typename StorageIndex_ >
StorageIndex* Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::supToCol ( )
inline

Return the array of supernode-to-column mapping

118 { return m_sup_to_col; }

References Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::m_sup_to_col.

◆ supToCol() [2/2]

template<typename Scalar_ , typename StorageIndex_ >
const StorageIndex* Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::supToCol ( ) const
inline

◆ valuePtr() [1/2]

template<typename Scalar_ , typename StorageIndex_ >
Scalar* Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::valuePtr ( )
inline

Return the array of nonzero values packed by column

The size is nnz

85 { return m_nzval; }

References Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::m_nzval.

◆ valuePtr() [2/2]

template<typename Scalar_ , typename StorageIndex_ >
const Scalar* Eigen::internal::MappedSuperNodalMatrix< Scalar_, StorageIndex_ >::valuePtr ( ) const
inline

Member Data Documentation

◆ m_col

◆ m_col_to_sup

◆ m_nsuper

◆ m_nzval

◆ m_nzval_colptr

◆ m_row

◆ m_rowind

◆ m_rowind_colptr

◆ m_sup_to_col


The documentation for this class was generated from the following file: