Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ > Class Template Reference

A versatile sparse matrix representation where each element is a block. More...

+ Inheritance diagram for Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >:

Public Types

enum  {
  Options ,
  Flags ,
  BlockSize ,
  RowsAtCompileTime ,
  ColsAtCompileTime ,
  MaxRowsAtCompileTime ,
  MaxColsAtCompileTime ,
  IsVectorAtCompileTime ,
  IsColMajor
}
 
typedef Matrix< RealScalar, _BlockAtCompileTime, _BlockAtCompileTime, IsColMajor ? ColMajor :RowMajorBlockRealScalar
 
typedef Matrix< Scalar, _BlockAtCompileTime, _BlockAtCompileTime, IsColMajor ? ColMajor :RowMajorBlockScalar
 
typedef std::conditional_t< _BlockAtCompileTime==Dynamic, Scalar, BlockScalarBlockScalarReturnType
 
typedef internal::ref_selector< BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ > >::type Nested
 
typedef BlockSparseMatrix< Scalar, BlockSize, IsColMajor ? ColMajor :RowMajor, StorageIndexPlainObject
 
typedef NumTraits< Scalar >::Real RealScalar
 
typedef Scalar_ Scalar
 
typedef StorageIndex_ StorageIndex
 
- Public Types inherited from Eigen::SparseMatrixBase< BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ > >
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, Eigen::Transpose< const Derived > >, Transpose< const Derived > > AdjointReturnType
 
typedef Transpose< const Derived > ConstTransposeReturnType
 
typedef Matrix< StorageIndex, Dynamic, 1 > IndexVector
 
typedef internal::add_const_on_value_type_if_arithmetic< typename internal::packet_traits< Scalar >::type >::type PacketReturnType
 
typedef internal::packet_traits< Scalar >::type PacketScalar
 
typedef SparseMatrix< Scalar, Flags &RowMajorBit ? RowMajor :ColMajor, StorageIndexPlainObject
 
typedef internal::traits< Derived >::Scalar Scalar
 
typedef Matrix< Scalar, Dynamic, 1 > ScalarVector
 
typedef SparseMatrixBase StorageBaseType
 
typedef internal::traits< Derived >::StorageIndex StorageIndex
 
typedef internal::traits< Derived >::StorageKind StorageKind
 
typedef Transpose< Derived > TransposeReturnType
 
typedef Scalar value_type
 
- Public Types inherited from Eigen::EigenBase< class >
typedef Eigen::Index Index
 
typedef internal::traits< Derived >::StorageKind StorageKind
 

Public Member Functions

Index blockCols () const
 
Index blockColsIndex (Index bj) const
 
Index blockInnerIndex (Index bi) const
 
Index blockInnerSize (Index bi) const
 
Index blockOuterIndex (Index bj) const
 
Index blockOuterSize (Index bj) const
 
Index blockPtr (Index id) const
 
Index blockRows () const
 
Index blockRowsIndex (Index bi) const
 
 BlockSparseMatrix ()
 
 BlockSparseMatrix (const BlockSparseMatrix &other)
 Copy-constructor. More...
 
template<typename MatrixType >
 BlockSparseMatrix (const MatrixType &spmat)
 Constructor from a sparse matrix. More...
 
 BlockSparseMatrix (Index brow, Index bcol)
 Construct and resize. More...
 
Map< const BlockScalarcoeff (Index brow, Index bcol) const
 
Ref< BlockScalarcoeffRef (Index brow, Index bcol)
 
Index cols () const
 
Index innerBlocks () const
 
StorageIndexinnerIndexPtr ()
 
const StorageIndexinnerIndexPtr () const
 
Index innerSize () const
 
Index innerToBlock (Index inner) const
 
bool isCompressed () const
 for compatibility purposes with the SparseMatrix class More...
 
Index nonZeros () const
 
Index nonZerosBlocks () const
 
template<typename VecType >
BlockSparseTimeDenseProduct< BlockSparseMatrix, VecType > operator* (const VecType &lhs) const
 
BlockSparseMatrixoperator= (BlockSparseMatrix other)
 
template<typename MatrixType >
BlockSparseMatrixoperator= (const MatrixType &spmat)
 Assignment from a sparse matrix with the same storage order. More...
 
Index outerBlocks () const
 
StorageIndexouterIndexPtr ()
 
const StorageIndexouterIndexPtr () const
 
Index outerSize () const
 
Index outerToBlock (Index outer) const
 
void reserve (const Index nonzerosblocks)
 Allocate the internal array of pointers to blocks and their inner indices. More...
 
void resize (Index brow, Index bcol)
 Set the number of rows and columns blocks. More...
 
Index rows () const
 
void setBlockLayout (const VectorXi &rowBlocks, const VectorXi &colBlocks)
 Set the row and column block layouts,. More...
 
void setBlockSize (Index blockSize)
 set the block size at runtime for fixed-size block layout More...
 
template<typename MatrixType >
void setBlockStructure (const MatrixType &blockPattern)
 Set the nonzero block pattern of the matrix. More...
 
template<typename InputIterator >
void setFromTriplets (const InputIterator &begin, const InputIterator &end)
 Fill values in a matrix from a triplet list. More...
 
BlockScalarReturnTypevaluePtr ()
 
 ~BlockSparseMatrix ()
 
- Public Member Functions inherited from Eigen::SparseMatrixBase< BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ > >
const AdjointReturnType adjoint () const
 
RealScalar blueNorm () const
 
Index cols () const
 
const SparseMatrixBase< Derived >::template CwiseProductDenseReturnType< OtherDerived >::Type cwiseProduct (const MatrixBase< OtherDerived > &other) const
 
const CwiseProductDenseReturnType< OtherDerived >::Type cwiseProduct (const MatrixBase< OtherDerived > &other) const
 
internal::traits< Derived >::Scalar dot (const MatrixBase< OtherDerived > &other) const
 
Scalar dot (const MatrixBase< OtherDerived > &other) const
 
internal::traits< Derived >::Scalar dot (const SparseMatrixBase< OtherDerived > &other) const
 
Scalar dot (const SparseMatrixBase< OtherDerived > &other) const
 
const internal::eval< Derived >::type eval () const
 
Index innerSize () const
 
bool isApprox (const MatrixBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
 
bool isApprox (const SparseMatrixBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
 
bool isRValue () const
 
bool isVector () const
 
Derived & markAsRValue ()
 
RealScalar norm () const
 
const Product< Derived, OtherDerived > operator* (const DiagonalBase< OtherDerived > &other) const
 
const Product< Derived, OtherDerived > operator* (const MatrixBase< OtherDerived > &other) const
 
const Product< Derived, OtherDerived, AliasFreeProductoperator* (const SparseMatrixBase< OtherDerived > &other) const
 
Derived & operator*= (const Scalar &other)
 
Derived & operator*= (const SparseMatrixBase< OtherDerived > &other)
 
Derived & operator+= (const DiagonalBase< OtherDerived > &other)
 
Derived & operator+= (const EigenBase< OtherDerived > &other)
 
Derived & operator+= (const SparseMatrixBase< OtherDerived > &other)
 
Derived & operator-= (const DiagonalBase< OtherDerived > &other)
 
Derived & operator-= (const EigenBase< OtherDerived > &other)
 
Derived & operator-= (const SparseMatrixBase< OtherDerived > &other)
 
Derived & operator/= (const Scalar &other)
 
Derived & operator= (const Derived &other)
 
Derived & operator= (const EigenBase< OtherDerived > &other)
 
Derived & operator= (const ReturnByValue< OtherDerived > &other)
 
Derived & operator= (const SparseMatrixBase< OtherDerived > &other)
 
Index outerSize () const
 
const SparseView< Derived > pruned (const Scalar &reference=Scalar(0), const RealScalar &epsilon=NumTraits< Scalar >::dummy_precision()) const
 
Index rows () const
 
SelfAdjointViewReturnType< UpLo >::Type selfadjointView ()
 
SparseMatrixBase< Derived >::template SelfAdjointViewReturnType< UpLo >::Type selfadjointView ()
 
SparseMatrixBase< Derived >::template ConstSelfAdjointViewReturnType< UpLo >::Type selfadjointView () const
 
ConstSelfAdjointViewReturnType< UpLo >::Type selfadjointView () const
 
Index size () const
 
 SparseMatrixBase ()
 
RealScalar squaredNorm () const
 
Scalar sum () const
 
DenseMatrixType toDense () const
 
TransposeReturnType transpose ()
 
const ConstTransposeReturnType transpose () const
 
const TriangularView< const Derived, Mode > triangularView () const
 
SparseSymmetricPermutationProduct< Derived, Upper|LowertwistedBy (const PermutationMatrix< Dynamic, Dynamic, StorageIndex > &perm) const
 
- Public Member Functions inherited from Eigen::EigenBase< class >
void addTo (Dest &dst) const
 
void applyThisOnTheLeft (Dest &dst) const
 
void applyThisOnTheRight (Dest &dst) const
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
Derived & const_cast_derived () const
 
const Derived & const_derived () const
 
Derived & derived ()
 
const Derived & derived () const
 
void evalTo (Dest &dst) const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 
void subTo (Dest &dst) const
 

Protected Member Functions

Map< BlockScalarinsert (Index brow, Index bcol)
 
- Protected Member Functions inherited from Eigen::SparseMatrixBase< BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ > >
Derived & assign (const OtherDerived &other)
 
void assignGeneric (const OtherDerived &other)
 

Protected Attributes

StorageIndexm_blockPtr
 
Index m_blockSize
 
StorageIndexm_indices
 
Index m_innerBSize
 
StorageIndexm_innerOffset
 
Index m_nonzeros
 
Index m_nonzerosblocks
 
Index m_outerBSize
 
StorageIndexm_outerIndex
 
StorageIndexm_outerOffset
 
Scalarm_values
 
- Protected Attributes inherited from Eigen::SparseMatrixBase< BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ > >
bool m_isRValue
 

Additional Inherited Members

- Public Attributes inherited from Eigen::SparseMatrixBase< BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ > >
 ColsAtCompileTime
 
 Flags
 
 InnerSizeAtCompileTime
 
 IsRowMajor
 
 IsVectorAtCompileTime
 
 MaxColsAtCompileTime
 
 MaxRowsAtCompileTime
 
 MaxSizeAtCompileTime
 
 NumDimensions
 
 RowsAtCompileTime
 
 SizeAtCompileTime
 
- Static Protected Member Functions inherited from Eigen::SparseMatrixBase< BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ > >
static StorageIndex convert_index (const Index idx)
 

Detailed Description

template<typename Scalar_, int _BlockAtCompileTime, int Options_, typename StorageIndex_>
class Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >

A versatile sparse matrix representation where each element is a block.

This class provides routines to manipulate block sparse matrices stored in a BSR-like representation. There are two main types :

  1. All blocks have the same number of rows and columns, called block size in the following. In this case, if this block size is known at compile time, it can be given as a template parameter like
    BlockSparseMatrix<Scalar, 3, ColMajor> bmat(b_rows, b_cols);
    Here, bmat is a b_rows x b_cols block sparse matrix where each coefficient is a 3x3 dense matrix. If the block size is fixed but will be given at runtime,
    BlockSparseMatrix<Scalar, Dynamic, ColMajor> bmat(b_rows, b_cols);
    bmat.setBlockSize(block_size);
  2. The second case is for variable-block sparse matrices. Here each block has its own dimensions. The only restriction is that all the blocks in a row (resp. a column) should have the same number of rows (resp. of columns). It is thus required in this case to describe the layout of the matrix by calling setBlockLayout(rowBlocks, colBlocks).

In any of the previous case, the matrix can be filled by calling setFromTriplets(). A regular sparse matrix can be converted to a block sparse matrix and vice versa. It is obviously required to describe the block layout beforehand by calling either setBlockSize() for fixed-size blocks or setBlockLayout for variable-size blocks.

Template Parameters
Scalar_The Scalar type
_BlockAtCompileTimeThe block layout option. It takes the following values Dynamic : block size known at runtime a numeric number : fixed-size block known at compile time

Definition at line 286 of file BlockSparseMatrix.h.

Member Typedef Documentation

◆ BlockRealScalar

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
typedef Matrix<RealScalar, _BlockAtCompileTime, _BlockAtCompileTime,IsColMajor ? ColMajor : RowMajor> Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::BlockRealScalar

Definition at line 306 of file BlockSparseMatrix.h.

◆ BlockScalar

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
typedef Matrix<Scalar, _BlockAtCompileTime, _BlockAtCompileTime,IsColMajor ? ColMajor : RowMajor> Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::BlockScalar

Definition at line 305 of file BlockSparseMatrix.h.

◆ BlockScalarReturnType

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
typedef std::conditional_t<_BlockAtCompileTime==Dynamic, Scalar, BlockScalar> Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::BlockScalarReturnType

Definition at line 307 of file BlockSparseMatrix.h.

◆ Nested

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
typedef internal::ref_selector<BlockSparseMatrix<Scalar_, _BlockAtCompileTime, Options_, StorageIndex_> >::type Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::Nested

Definition at line 292 of file BlockSparseMatrix.h.

◆ PlainObject

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
typedef BlockSparseMatrix<Scalar, BlockSize, IsColMajor ? ColMajor : RowMajor, StorageIndex> Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::PlainObject

Definition at line 308 of file BlockSparseMatrix.h.

◆ RealScalar

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
typedef NumTraits<Scalar>::Real Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::RealScalar

Definition at line 290 of file BlockSparseMatrix.h.

◆ Scalar

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
typedef Scalar_ Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::Scalar

Definition at line 289 of file BlockSparseMatrix.h.

◆ StorageIndex

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
typedef StorageIndex_ Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::StorageIndex

Definition at line 291 of file BlockSparseMatrix.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
anonymous enum
Enumerator
Options 
Flags 
BlockSize 
RowsAtCompileTime 
ColsAtCompileTime 
MaxRowsAtCompileTime 
MaxColsAtCompileTime 
IsVectorAtCompileTime 
IsColMajor 

Definition at line 294 of file BlockSparseMatrix.h.

Constructor & Destructor Documentation

◆ BlockSparseMatrix() [1/4]

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::BlockSparseMatrix ( )
inline

◆ BlockSparseMatrix() [2/4]

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::BlockSparseMatrix ( Index  brow,
Index  bcol 
)
inline

Construct and resize.

Definition at line 322 of file BlockSparseMatrix.h.

323  : m_innerBSize(IsColMajor ? brow : bcol),
324  m_outerBSize(IsColMajor ? bcol : brow),
326  m_values(0),m_blockPtr(0),m_indices(0),
328  { }

◆ BlockSparseMatrix() [3/4]

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::BlockSparseMatrix ( const BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ > &  other)
inline

Copy-constructor.

Definition at line 333 of file BlockSparseMatrix.h.

334  : m_innerBSize(other.m_innerBSize),m_outerBSize(other.m_outerBSize),
335  m_nonzerosblocks(other.m_nonzerosblocks),m_nonzeros(other.m_nonzeros),
336  m_blockPtr(0),m_blockSize(other.m_blockSize)
337  {
338  // should we allow copying between variable-size blocks and fixed-size blocks ??
339  eigen_assert(m_blockSize == BlockSize && " CAN NOT COPY BETWEEN FIXED-SIZE AND VARIABLE-SIZE BLOCKS");
340 
341  std::copy(other.m_innerOffset, other.m_innerOffset+m_innerBSize+1, m_innerOffset);
342  std::copy(other.m_outerOffset, other.m_outerOffset+m_outerBSize+1, m_outerOffset);
343  std::copy(other.m_values, other.m_values+m_nonzeros, m_values);
344 
345  if(m_blockSize != Dynamic)
346  std::copy(other.m_blockPtr, other.m_blockPtr+m_nonzerosblocks, m_blockPtr);
347 
348  std::copy(other.m_indices, other.m_indices+m_nonzerosblocks, m_indices);
349  std::copy(other.m_outerIndex, other.m_outerIndex+m_outerBSize, m_outerIndex);
350  }
#define eigen_assert(x)

◆ ~BlockSparseMatrix()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::~BlockSparseMatrix ( )
inline

Definition at line 375 of file BlockSparseMatrix.h.

376  {
377  delete[] m_outerIndex;
378  delete[] m_innerOffset;
379  delete[] m_outerOffset;
380  delete[] m_indices;
381  delete[] m_blockPtr;
382  delete[] m_values;
383  }

◆ BlockSparseMatrix() [4/4]

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
template<typename MatrixType >
Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::BlockSparseMatrix ( const MatrixType spmat)
inline

Constructor from a sparse matrix.

Definition at line 391 of file BlockSparseMatrix.h.

392  {
393  EIGEN_STATIC_ASSERT((m_blockSize != Dynamic), THIS_METHOD_IS_ONLY_FOR_FIXED_SIZE);
394 
395  *this = spmat;
396  }
#define EIGEN_STATIC_ASSERT(X, MSG)

Member Function Documentation

◆ blockCols()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::blockCols ( ) const
inline
Returns
the number of columns grouped by blocks

Definition at line 772 of file BlockSparseMatrix.h.

773  {
774  return (IsColMajor ? m_outerBSize : m_innerBSize);
775  }

◆ blockColsIndex()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::blockColsIndex ( Index  bj) const
inline
Returns
the starting index of the bj col block

Definition at line 886 of file BlockSparseMatrix.h.

887  {
888  return IsColMajor ? blockOuterIndex(bj) : blockInnerIndex(bj);
889  }
Index blockInnerIndex(Index bi) const
Index blockOuterIndex(Index bj) const

◆ blockInnerIndex()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::blockInnerIndex ( Index  bi) const
inline

Definition at line 895 of file BlockSparseMatrix.h.

896  {
897  return (m_blockSize == Dynamic) ? m_innerOffset[bi] : (bi * m_blockSize);
898  }

◆ blockInnerSize()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::blockInnerSize ( Index  bi) const
inline

Definition at line 901 of file BlockSparseMatrix.h.

902  {
903  return (m_blockSize == Dynamic) ? (m_innerOffset[bi+1] - m_innerOffset[bi]) : m_blockSize;
904  }

◆ blockOuterIndex()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::blockOuterIndex ( Index  bj) const
inline

Definition at line 891 of file BlockSparseMatrix.h.

892  {
893  return (m_blockSize == Dynamic) ? m_outerOffset[bj] : (bj * m_blockSize);
894  }

◆ blockOuterSize()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::blockOuterSize ( Index  bj) const
inline

Definition at line 905 of file BlockSparseMatrix.h.

906  {
907  return (m_blockSize == Dynamic) ? (m_outerOffset[bj+1]- m_outerOffset[bj]) : m_blockSize;
908  }

◆ blockPtr()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::blockPtr ( Index  id) const
inline
Returns
the starting position of the block id in the array of values

Definition at line 938 of file BlockSparseMatrix.h.

939  {
940  if(m_blockSize == Dynamic) return m_blockPtr[id];
941  else return id * m_blockSize * m_blockSize;
942  //return blockDynIdx(id, std::conditional_t<(BlockSize==Dynamic), internal::true_type, internal::false_type>());
943  }

◆ blockRows()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::blockRows ( ) const
inline
Returns
the number of rows grouped by blocks

Definition at line 767 of file BlockSparseMatrix.h.

768  {
769  return (IsColMajor ? m_innerBSize : m_outerBSize);
770  }

◆ blockRowsIndex()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::blockRowsIndex ( Index  bi) const
inline
Returns
the starting index of the bi row block

Definition at line 878 of file BlockSparseMatrix.h.

879  {
880  return IsColMajor ? blockInnerIndex(bi) : blockOuterIndex(bi);
881  }

◆ coeff()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Map<const BlockScalar> Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::coeff ( Index  brow,
Index  bcol 
) const
inline
Returns
the value of the (i,j) block as an Eigen Dense Matrix

Definition at line 834 of file BlockSparseMatrix.h.

835  {
836  eigen_assert(brow < blockRows() && "BLOCK ROW INDEX OUT OF BOUNDS");
837  eigen_assert(bcol < blockCols() && "BLOCK COLUMN OUT OF BOUNDS");
838 
839  StorageIndex rsize = IsColMajor ? blockInnerSize(brow): blockOuterSize(bcol);
840  StorageIndex csize = IsColMajor ? blockOuterSize(bcol) : blockInnerSize(brow);
841  StorageIndex inner = IsColMajor ? brow : bcol;
842  StorageIndex outer = IsColMajor ? bcol : brow;
843  StorageIndex offset = m_outerIndex[outer];
844  while(offset < m_outerIndex[outer+1] && m_indices[offset] != inner) offset++;
845  if(m_indices[offset] == inner)
846  {
847  return Map<const BlockScalar> (&(m_values[blockPtr(offset)]), rsize, csize);
848  }
849  else
850 // return BlockScalar::Zero(rsize, csize);
851  eigen_assert("NOT YET SUPPORTED");
852  }
Index blockPtr(Index id) const
Index blockInnerSize(Index bi) const
Index blockOuterSize(Index bj) const

◆ coeffRef()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Ref<BlockScalar> Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::coeffRef ( Index  brow,
Index  bcol 
)
inline
Returns
a reference to the (i,j) block as an Eigen Dense Matrix

Definition at line 808 of file BlockSparseMatrix.h.

809  {
810  eigen_assert(brow < blockRows() && "BLOCK ROW INDEX OUT OF BOUNDS");
811  eigen_assert(bcol < blockCols() && "BLOCK nzblocksFlagCOLUMN OUT OF BOUNDS");
812 
813  StorageIndex rsize = IsColMajor ? blockInnerSize(brow): blockOuterSize(bcol);
814  StorageIndex csize = IsColMajor ? blockOuterSize(bcol) : blockInnerSize(brow);
815  StorageIndex inner = IsColMajor ? brow : bcol;
816  StorageIndex outer = IsColMajor ? bcol : brow;
817  StorageIndex offset = m_outerIndex[outer];
818  while(offset < m_outerIndex[outer+1] && m_indices[offset] != inner)
819  offset++;
820  if(m_indices[offset] == inner)
821  {
822  return Map<BlockScalar>(&(m_values[blockPtr(offset)]), rsize, csize);
823  }
824  else
825  {
826  //FIXME the block does not exist, Insert it !!!!!!!!!
827  eigen_assert("DYNAMIC INSERTION IS NOT YET SUPPORTED");
828  }
829  }

◆ cols()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::cols ( void  ) const
inline
Returns
the number of cols

Definition at line 749 of file BlockSparseMatrix.h.

750  {
751 // return blockCols();
752  return (IsColMajor ? outerSize() : innerSize());
753  }

◆ innerBlocks()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::innerBlocks ( ) const
inline

Definition at line 778 of file BlockSparseMatrix.h.

778 { return m_innerBSize; }

◆ innerIndexPtr() [1/2]

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
StorageIndex* Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::innerIndexPtr ( )
inline

Definition at line 868 of file BlockSparseMatrix.h.

868 {return m_indices; }

◆ innerIndexPtr() [2/2]

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
const StorageIndex* Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::innerIndexPtr ( ) const
inline

Definition at line 869 of file BlockSparseMatrix.h.

869 {return m_indices; }

◆ innerSize()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::innerSize ( ) const
inline

Definition at line 755 of file BlockSparseMatrix.h.

756  {
758  else return (m_innerBSize * m_blockSize) ;
759  }

◆ innerToBlock()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::innerToBlock ( Index  inner) const
inline
Returns
the block index where inner belongs to

Definition at line 793 of file BlockSparseMatrix.h.

794  {
795  eigen_assert(inner < innerSize() && "OUTER INDEX OUT OF BOUNDS");
796 
797  if(m_blockSize != Dynamic)
798  return (inner / m_blockSize); // Integer division
799 
800  StorageIndex b_inner = 0;
801  while(m_innerOffset[b_inner] <= inner) ++b_inner;
802  return b_inner - 1;
803  }

◆ insert()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Map<BlockScalar> Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::insert ( Index  brow,
Index  bcol 
)
protected

◆ isCompressed()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
bool Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::isCompressed ( ) const
inline

for compatibility purposes with the SparseMatrix class

Definition at line 874 of file BlockSparseMatrix.h.

874 {return true;}

◆ nonZeros()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::nonZeros ( ) const
inline
Returns
the total number of nonzero elements, including eventual explicit zeros in blocks

Definition at line 864 of file BlockSparseMatrix.h.

864 { return m_nonzeros; }

◆ nonZerosBlocks()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::nonZerosBlocks ( ) const
inline
Returns
the number of nonzero blocks

Definition at line 862 of file BlockSparseMatrix.h.

862 { return m_nonzerosblocks; }

◆ operator*()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
template<typename VecType >
BlockSparseTimeDenseProduct<BlockSparseMatrix, VecType> Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::operator* ( const VecType &  lhs) const
inline

Definition at line 856 of file BlockSparseMatrix.h.

857  {
858  return BlockSparseTimeDenseProduct<BlockSparseMatrix, VecType>(*this, lhs);
859  }

◆ operator=() [1/2]

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
BlockSparseMatrix& Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::operator= ( BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >  other)
inline

Definition at line 367 of file BlockSparseMatrix.h.

368  {
369  //Copy-and-swap paradigm ... avoid leaked data if thrown
370  swap(*this, other);
371  return *this;
372  }
friend void swap(BlockSparseMatrix &first, BlockSparseMatrix &second)

◆ operator=() [2/2]

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
template<typename MatrixType >
BlockSparseMatrix& Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::operator= ( const MatrixType spmat)
inline

Assignment from a sparse matrix with the same storage order.

Convert from a sparse matrix to block sparse matrix.

Warning
Before calling this function, tt is necessary to call either setBlockLayout() (matrices with variable-size blocks) or setBlockSize() (for fixed-size blocks).

Definition at line 407 of file BlockSparseMatrix.h.

408  {
409  eigen_assert((m_innerBSize != 0 && m_outerBSize != 0)
410  && "Trying to assign to a zero-size matrix, call resize() first");
411  eigen_assert(((MatrixType::Options&RowMajorBit) != IsColMajor) && "Wrong storage order");
412  typedef SparseMatrix<bool,MatrixType::Options,typename MatrixType::Index> MatrixPatternType;
413  MatrixPatternType blockPattern(blockRows(), blockCols());
414  m_nonzeros = 0;
415 
416  // First, compute the number of nonzero blocks and their locations
417  for(StorageIndex bj = 0; bj < m_outerBSize; ++bj)
418  {
419  // Browse each outer block and compute the structure
420  std::vector<bool> nzblocksFlag(m_innerBSize,false); // Record the existing blocks
421  blockPattern.startVec(bj);
422  for(StorageIndex j = blockOuterIndex(bj); j < blockOuterIndex(bj+1); ++j)
423  {
424  typename MatrixType::InnerIterator it_spmat(spmat, j);
425  for(; it_spmat; ++it_spmat)
426  {
427  StorageIndex bi = innerToBlock(it_spmat.index()); // Index of the current nonzero block
428  if(!nzblocksFlag[bi])
429  {
430  // Save the index of this nonzero block
431  nzblocksFlag[bi] = true;
432  blockPattern.insertBackByOuterInnerUnordered(bj, bi) = true;
433  // Compute the total number of nonzeros (including explicit zeros in blocks)
435  }
436  }
437  } // end current outer block
438  }
439  blockPattern.finalize();
440 
441  // Allocate the internal arrays
442  setBlockStructure(blockPattern);
443 
444  for(StorageIndex nz = 0; nz < m_nonzeros; ++nz) m_values[nz] = Scalar(0);
445  for(StorageIndex bj = 0; bj < m_outerBSize; ++bj)
446  {
447  // Now copy the values
448  for(StorageIndex j = blockOuterIndex(bj); j < blockOuterIndex(bj+1); ++j)
449  {
450  // Browse the outer block column by column (for column-major matrices)
451  typename MatrixType::InnerIterator it_spmat(spmat, j);
452  for(; it_spmat; ++it_spmat)
453  {
454  StorageIndex idx = 0; // Position of this block in the column block
455  StorageIndex bi = innerToBlock(it_spmat.index()); // Index of the current nonzero block
456  // Go to the inner block where this element belongs to
457  while(bi > m_indices[m_outerIndex[bj]+idx]) ++idx; // Not expensive for ordered blocks
458  StorageIndex idxVal;// Get the right position in the array of values for this element
459  if(m_blockSize == Dynamic)
460  {
461  // Offset from all blocks before ...
462  idxVal = m_blockPtr[m_outerIndex[bj]+idx];
463  // ... and offset inside the block
464  idxVal += (j - blockOuterIndex(bj)) * blockOuterSize(bj) + it_spmat.index() - m_innerOffset[bi];
465  }
466  else
467  {
468  // All blocks before
469  idxVal = (m_outerIndex[bj] + idx) * m_blockSize * m_blockSize;
470  // inside the block
471  idxVal += (j - blockOuterIndex(bj)) * m_blockSize + (it_spmat.index()%m_blockSize);
472  }
473  // Insert the value
474  m_values[idxVal] = it_spmat.value();
475  } // end of this column
476  } // end of this block
477  } // end of this outer block
478 
479  return *this;
480  }
Index innerToBlock(Index inner) const
void setBlockStructure(const MatrixType &blockPattern)
Set the nonzero block pattern of the matrix.
std::ptrdiff_t j

◆ outerBlocks()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::outerBlocks ( ) const
inline

Definition at line 777 of file BlockSparseMatrix.h.

777 { return m_outerBSize; }

◆ outerIndexPtr() [1/2]

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
StorageIndex* Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::outerIndexPtr ( )
inline

Definition at line 870 of file BlockSparseMatrix.h.

870 {return m_outerIndex; }

◆ outerIndexPtr() [2/2]

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
const StorageIndex* Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::outerIndexPtr ( ) const
inline

Definition at line 871 of file BlockSparseMatrix.h.

871 {return m_outerIndex; }

◆ outerSize()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::outerSize ( ) const
inline

Definition at line 761 of file BlockSparseMatrix.h.

762  {
764  else return (m_outerBSize * m_blockSize) ;
765  }

◆ outerToBlock()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::outerToBlock ( Index  outer) const
inline
Returns
the block index where outer belongs to

Definition at line 781 of file BlockSparseMatrix.h.

782  {
783  eigen_assert(outer < outerSize() && "OUTER INDEX OUT OF BOUNDS");
784 
785  if(m_blockSize != Dynamic)
786  return (outer / m_blockSize); // Integer division
787 
788  StorageIndex b_outer = 0;
789  while(m_outerOffset[b_outer] <= outer) ++b_outer;
790  return b_outer - 1;
791  }

◆ reserve()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
void Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::reserve ( const Index  nonzerosblocks)
inline

Allocate the internal array of pointers to blocks and their inner indices.

Note
For fixed-size blocks, call setBlockSize() to set the block. And For variable-size blocks, call setBlockLayout() before using this function
Parameters
nonzerosblocksNumber of nonzero blocks. The total number of nonzeros is is computed in setBlockLayout() for variable-size blocks
See also
setBlockSize()

Definition at line 598 of file BlockSparseMatrix.h.

599  {
600  eigen_assert((m_innerBSize != 0 && m_outerBSize != 0) &&
601  "TRYING TO RESERVE ZERO-SIZE MATRICES, CALL resize() first");
602 
603  //FIXME Should free if already allocated
605 
606  m_nonzerosblocks = nonzerosblocks;
607  if(m_blockSize != Dynamic)
608  {
609  m_nonzeros = nonzerosblocks * (m_blockSize * m_blockSize);
610  m_blockPtr = 0;
611  }
612  else
613  {
614  // m_nonzeros is already computed in setBlockLayout()
616  }
618  m_values = new Scalar[m_nonzeros];
619  }

◆ resize()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
void Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::resize ( Index  brow,
Index  bcol 
)
inline

Set the number of rows and columns blocks.

Definition at line 540 of file BlockSparseMatrix.h.

541  {
542  m_innerBSize = IsColMajor ? brow : bcol;
543  m_outerBSize = IsColMajor ? bcol : brow;
544  }

◆ rows()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::rows ( void  ) const
inline
Returns
the number of rows

Definition at line 740 of file BlockSparseMatrix.h.

741  {
742 // return blockRows();
743  return (IsColMajor ? innerSize() : outerSize());
744  }

◆ setBlockLayout()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
void Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::setBlockLayout ( const VectorXi rowBlocks,
const VectorXi colBlocks 
)
inline

Set the row and column block layouts,.

This function set the size of each row and column block. So this function should be used only for blocks with variable size.

Parameters
rowBlocks: Number of rows per row block
colBlocks: Number of columns per column block
See also
resize(), setBlockSize()

Definition at line 565 of file BlockSparseMatrix.h.

566  {
567  const VectorXi& innerBlocks = IsColMajor ? rowBlocks : colBlocks;
568  const VectorXi& outerBlocks = IsColMajor ? colBlocks : rowBlocks;
569  eigen_assert(m_innerBSize == innerBlocks.size() && "CHECK THE NUMBER OF ROW OR COLUMN BLOCKS");
570  eigen_assert(m_outerBSize == outerBlocks.size() && "CHECK THE NUMBER OF ROW OR COLUMN BLOCKS");
571  m_outerBSize = outerBlocks.size();
572  // starting index of blocks... cumulative sums
575  m_innerOffset[0] = 0;
576  m_outerOffset[0] = 0;
577  std::partial_sum(&innerBlocks[0], &innerBlocks[m_innerBSize-1]+1, &m_innerOffset[1]);
578  std::partial_sum(&outerBlocks[0], &outerBlocks[m_outerBSize-1]+1, &m_outerOffset[1]);
579 
580  // Compute the total number of nonzeros
581  m_nonzeros = 0;
582  for(StorageIndex bj = 0; bj < m_outerBSize; ++bj)
583  for(StorageIndex bi = 0; bi < m_innerBSize; ++bi)
584  m_nonzeros += outerBlocks[bj] * innerBlocks[bi];
585 
586  }
Matrix< int, Dynamic, 1 > VectorXi

◆ setBlockSize()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
void Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::setBlockSize ( Index  blockSize)
inline

set the block size at runtime for fixed-size block layout

Call this only for fixed-size blocks

Definition at line 551 of file BlockSparseMatrix.h.

552  {
553  m_blockSize = blockSize;
554  }

◆ setBlockStructure()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
template<typename MatrixType >
void Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::setBlockStructure ( const MatrixType blockPattern)
inline

Set the nonzero block pattern of the matrix.

Given a sparse matrix describing the nonzero block pattern, this function prepares the internal pointers for values. After calling this function, any nonzero block (bi, bj) can be set with a simple call to coeffRef(bi,bj).

Warning
Before calling this function, tt is necessary to call either setBlockLayout() (matrices with variable-size blocks) or setBlockSize() (for fixed-size blocks).
Parameters
blockPatternSparse matrix of boolean elements describing the block structure
See also
setBlockLayout()
setBlockSize()

Definition at line 500 of file BlockSparseMatrix.h.

501  {
502  resize(blockPattern.rows(), blockPattern.cols());
503  reserve(blockPattern.nonZeros());
504 
505  // Browse the block pattern and set up the various pointers
506  m_outerIndex[0] = 0;
507  if(m_blockSize == Dynamic) m_blockPtr[0] = 0;
508  for(StorageIndex nz = 0; nz < m_nonzeros; ++nz) m_values[nz] = Scalar(0);
509  for(StorageIndex bj = 0; bj < m_outerBSize; ++bj)
510  {
511  //Browse each outer block
512 
513  //First, copy and save the indices of nonzero blocks
514  //FIXME : find a way to avoid this ...
515  std::vector<int> nzBlockIdx;
516  typename MatrixType::InnerIterator it(blockPattern, bj);
517  for(; it; ++it)
518  {
519  nzBlockIdx.push_back(it.index());
520  }
521  std::sort(nzBlockIdx.begin(), nzBlockIdx.end());
522 
523  // Now, fill block indices and (eventually) pointers to blocks
524  for(StorageIndex idx = 0; idx < nzBlockIdx.size(); ++idx)
525  {
526  StorageIndex offset = m_outerIndex[bj]+idx; // offset in m_indices
527  m_indices[offset] = nzBlockIdx[idx];
528  if(m_blockSize == Dynamic)
529  m_blockPtr[offset] = m_blockPtr[offset-1] + blockInnerSize(nzBlockIdx[idx]) * blockOuterSize(bj);
530  // There is no blockPtr for fixed-size blocks... not needed !???
531  }
532  // Save the pointer to the next outer block
533  m_outerIndex[bj+1] = m_outerIndex[bj] + nzBlockIdx.size();
534  }
535  }
void reserve(const Index nonzerosblocks)
Allocate the internal array of pointers to blocks and their inner indices.
void resize(Index brow, Index bcol)
Set the number of rows and columns blocks.

◆ setFromTriplets()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
template<typename InputIterator >
void Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::setFromTriplets ( const InputIterator &  begin,
const InputIterator &  end 
)
inline

Fill values in a matrix from a triplet list.

Each triplet item has a block stored in an Eigen dense matrix. The InputIterator class should provide the functions row(), col() and value()

Note
For fixed-size blocks, call setBlockSize() before this function.

FIXME Do not accept duplicates

Definition at line 633 of file BlockSparseMatrix.h.

634  {
635  eigen_assert((m_innerBSize!=0 && m_outerBSize !=0) && "ZERO BLOCKS, PLEASE CALL resize() before");
636 
637  /* First, sort the triplet list
638  * FIXME This can be unnecessarily expensive since only the inner indices have to be sorted
639  * The best approach is like in SparseMatrix::setFromTriplets()
640  */
641  internal::TripletComp<InputIterator, IsColMajor> tripletcomp;
642  std::sort(begin, end, tripletcomp);
643 
644  /* Count the number of rows and column blocks,
645  * and the number of nonzero blocks per outer dimension
646  */
647  VectorXi rowBlocks(m_innerBSize); // Size of each block row
648  VectorXi colBlocks(m_outerBSize); // Size of each block column
649  rowBlocks.setZero(); colBlocks.setZero();
650  VectorXi nzblock_outer(m_outerBSize); // Number of nz blocks per outer vector
651  VectorXi nz_outer(m_outerBSize); // Number of nz per outer vector...for variable-size blocks
652  nzblock_outer.setZero();
653  nz_outer.setZero();
654  for(InputIterator it(begin); it !=end; ++it)
655  {
656  eigen_assert(it->row() >= 0 && it->row() < this->blockRows() && it->col() >= 0 && it->col() < this->blockCols());
657  eigen_assert((it->value().rows() == it->value().cols() && (it->value().rows() == m_blockSize))
658  || (m_blockSize == Dynamic));
659 
660  if(m_blockSize == Dynamic)
661  {
662  eigen_assert((rowBlocks[it->row()] == 0 || rowBlocks[it->row()] == it->value().rows()) &&
663  "NON CORRESPONDING SIZES FOR ROW BLOCKS");
664  eigen_assert((colBlocks[it->col()] == 0 || colBlocks[it->col()] == it->value().cols()) &&
665  "NON CORRESPONDING SIZES FOR COLUMN BLOCKS");
666  rowBlocks[it->row()] =it->value().rows();
667  colBlocks[it->col()] = it->value().cols();
668  }
669  nz_outer(IsColMajor ? it->col() : it->row()) += it->value().rows() * it->value().cols();
670  nzblock_outer(IsColMajor ? it->col() : it->row())++;
671  }
672  // Allocate member arrays
673  if(m_blockSize == Dynamic) setBlockLayout(rowBlocks, colBlocks);
674  StorageIndex nzblocks = nzblock_outer.sum();
675  reserve(nzblocks);
676 
677  // Temporary markers
678  VectorXi block_id(m_outerBSize); // To be used as a block marker during insertion
679 
680  // Setup outer index pointers and markers
681  m_outerIndex[0] = 0;
682  if (m_blockSize == Dynamic) m_blockPtr[0] = 0;
683  for(StorageIndex bj = 0; bj < m_outerBSize; ++bj)
684  {
685  m_outerIndex[bj+1] = m_outerIndex[bj] + nzblock_outer(bj);
686  block_id(bj) = m_outerIndex[bj];
687  if(m_blockSize==Dynamic)
688  {
689  m_blockPtr[m_outerIndex[bj+1]] = m_blockPtr[m_outerIndex[bj]] + nz_outer(bj);
690  }
691  }
692 
693  // Fill the matrix
694  for(InputIterator it(begin); it!=end; ++it)
695  {
696  StorageIndex outer = IsColMajor ? it->col() : it->row();
697  StorageIndex inner = IsColMajor ? it->row() : it->col();
698  m_indices[block_id(outer)] = inner;
699  StorageIndex block_size = it->value().rows()*it->value().cols();
700  StorageIndex nz_marker = blockPtr(block_id[outer]);
701  memcpy(&(m_values[nz_marker]), it->value().data(), block_size * sizeof(Scalar));
702  if(m_blockSize == Dynamic)
703  {
704  m_blockPtr[block_id(outer)+1] = m_blockPtr[block_id(outer)] + block_size;
705  }
706  block_id(outer)++;
707  }
708 
709  // An alternative when the outer indices are sorted...no need to use an array of markers
710 // for(Index bcol = 0; bcol < m_outerBSize; ++bcol)
711 // {
712 // Index id = 0, id_nz = 0, id_nzblock = 0;
713 // for(InputIterator it(begin); it!=end; ++it)
714 // {
715 // while (id<bcol) // one pass should do the job unless there are empty columns
716 // {
717 // id++;
718 // m_outerIndex[id+1]=m_outerIndex[id];
719 // }
720 // m_outerIndex[id+1] += 1;
721 // m_indices[id_nzblock]=brow;
722 // Index block_size = it->value().rows()*it->value().cols();
723 // m_blockPtr[id_nzblock+1] = m_blockPtr[id_nzblock] + block_size;
724 // id_nzblock++;
725 // memcpy(&(m_values[id_nz]),it->value().data(), block_size*sizeof(Scalar));
726 // id_nz += block_size;
727 // }
728 // while(id < m_outerBSize-1) // Empty columns at the end
729 // {
730 // id++;
731 // m_outerIndex[id+1]=m_outerIndex[id];
732 // }
733 // }
734  }
void setBlockLayout(const VectorXi &rowBlocks, const VectorXi &colBlocks)
Set the row and column block layouts,.
static const lastp1_t end

◆ valuePtr()

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
BlockScalarReturnType* Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::valuePtr ( )
inline

Definition at line 866 of file BlockSparseMatrix.h.

866 {return static_cast<BlockScalarReturnType *>(m_values);}
std::conditional_t< _BlockAtCompileTime==Dynamic, Scalar, BlockScalar > BlockScalarReturnType

Member Data Documentation

◆ m_blockPtr

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
StorageIndex* Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::m_blockPtr
protected

Definition at line 967 of file BlockSparseMatrix.h.

◆ m_blockSize

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::m_blockSize
protected

Definition at line 970 of file BlockSparseMatrix.h.

◆ m_indices

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
StorageIndex* Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::m_indices
protected

Definition at line 968 of file BlockSparseMatrix.h.

◆ m_innerBSize

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::m_innerBSize
protected

Definition at line 960 of file BlockSparseMatrix.h.

◆ m_innerOffset

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
StorageIndex* Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::m_innerOffset
protected

Definition at line 962 of file BlockSparseMatrix.h.

◆ m_nonzeros

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::m_nonzeros
protected

Definition at line 965 of file BlockSparseMatrix.h.

◆ m_nonzerosblocks

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::m_nonzerosblocks
protected

Definition at line 964 of file BlockSparseMatrix.h.

◆ m_outerBSize

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Index Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::m_outerBSize
protected

Definition at line 961 of file BlockSparseMatrix.h.

◆ m_outerIndex

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
StorageIndex* Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::m_outerIndex
protected

Definition at line 969 of file BlockSparseMatrix.h.

◆ m_outerOffset

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
StorageIndex* Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::m_outerOffset
protected

Definition at line 963 of file BlockSparseMatrix.h.

◆ m_values

template<typename Scalar_ , int _BlockAtCompileTime, int Options_, typename StorageIndex_ >
Scalar* Eigen::BlockSparseMatrix< Scalar_, _BlockAtCompileTime, Options_, StorageIndex_ >::m_values
protected

Definition at line 966 of file BlockSparseMatrix.h.


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