LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
lf::assemble::COOMatrix< SCALAR > Class Template Reference

A temporary data structure for matrix in COO format. More...

#include <lf/assemble/coomatrix.h>

Public Types

using Scalar = SCALAR
 
using Triplet = Eigen::Triplet<SCALAR>
 
using TripletVec = std::vector<Triplet>
 
using Index = Eigen::Index
 

Public Member Functions

 COOMatrix (size_type num_rows, size_type num_cols)
 
 COOMatrix (const COOMatrix &)=default
 
 COOMatrix (COOMatrix &&) noexcept=default
 
COOMatrixoperator= (const COOMatrix &)=default
 
COOMatrixoperator= (COOMatrix &&) noexcept=default
 
 ~COOMatrix ()=default
 
Index rows () const
 return number of rows
 
Index cols () const
 return number of column
 
void AddToEntry (gdof_idx_t i, gdof_idx_t j, SCALAR increment)
 Add a value to the specified entry.
 
void setZero ()
 Erase all entries of the matrix.
 
template<typename PREDICATE>
void setZero (PREDICATE &&pred)
 Erase specific entries of the COO matrix, that is, set them to zero.
 
TripletVectriplets ()
 Gives access to the vector of triplets.
 
const TripletVectriplets () const
 
template<typename VECTOR>
Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 > MatVecMult (SCALAR alpha, const VECTOR &vec) const
 Computes the product of a (scaled) vector with the matrix in COO format.
 
template<typename VECTOR, typename RESULTVECTOR>
void MatVecMult (SCALAR alpha, const VECTOR &vec, RESULTVECTOR &resvec) const
 In-situ computes the product of a vector with the matrix in COO format.
 
Eigen::SparseMatrix< ScalarmakeSparse () const
 Create an Eigen::SparseMatrix out of the COO format.
 
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > makeDense () const
 Create an Eigen::MatrixX from the COO format.
 
void PrintInfo (std::ostream &o) const
 Ouput of triplet list describing COO matrix.
 

Private Attributes

size_type rows_
 
size_type cols_
 
TripletVec triplets_ {}
 

Friends

template<typename SCALARTYPE>
std::ostream & operator<< (std::ostream &o, const COOMatrix< SCALARTYPE > &mat)
 Output operator for COO matrix.
 

Detailed Description

template<typename SCALAR>
class lf::assemble::COOMatrix< SCALAR >

A temporary data structure for matrix in COO format.

Template Parameters
SCALARbasic scalar type for the matrix

This class provides a container for a matrix in triplet format, also known as COO format see Wikipedia. It essentially manages a vector of Eigen triplets_.

This matrix format allows efficient incremental matrix construction in the context of local assembly.

type requirements for template arguments

SCALAR must be a type that can serve as a scalar type for Eigen::Matrix.

Warning
Objects of this type are no matrix objects in the sense of Eigen. This means you cannot use them instead of an Eigen matrix type, but you have to convert them into a sparse Eigen matrix before by calling COOMatrix::makeSparse().

sample usage:

// Obtain a shared pointed to a mesh
std::shared_ptr<lf::mesh::Mesh> mesh_p{
// Initialization of index mapping for linear finite elements
lf::assemble::UniformFEDofHandler loc_glob_map(
mesh_p, {{lf::base::RefEl::kPoint(), 1}});
// Initialize ENTITYMATRIXPROVIDER object for local computations
lf::uscalfe::LinearFELaplaceElementMatrix loc_mat_laplace{};
// Dimension of finite element space
const lf::assemble::size_type N_dofs(loc_glob_map.NumDofs());
// Matrix in triplet format holding Galerkin matrix
lf::assemble::COOMatrix<double> mat(N_dofs, N_dofs);
// Building the Galerkin matrix (trial space = test space)
0, loc_glob_map, loc_glob_map, loc_mat_laplace, mat);
// the `mat` object now contains the Galerkin matrix in triplet/COO format
// Initialize an Eigen sparse matrix object from `mat`
Eigen::SparseMatrix<double> stiffness_matrix = {mat.makeSparse()};

Definition at line 52 of file coomatrix.h.

Member Typedef Documentation

◆ Index

template<typename SCALAR>
using lf::assemble::COOMatrix< SCALAR >::Index = Eigen::Index

Definition at line 57 of file coomatrix.h.

◆ Scalar

template<typename SCALAR>
using lf::assemble::COOMatrix< SCALAR >::Scalar = SCALAR

Definition at line 54 of file coomatrix.h.

◆ Triplet

template<typename SCALAR>
using lf::assemble::COOMatrix< SCALAR >::Triplet = Eigen::Triplet<SCALAR>

Definition at line 55 of file coomatrix.h.

◆ TripletVec

template<typename SCALAR>
using lf::assemble::COOMatrix< SCALAR >::TripletVec = std::vector<Triplet>

Definition at line 56 of file coomatrix.h.

Constructor & Destructor Documentation

◆ COOMatrix() [1/3]

template<typename SCALAR>
lf::assemble::COOMatrix< SCALAR >::COOMatrix ( size_type num_rows,
size_type num_cols )
inline

Set up zero matrix of a given size

Definition at line 60 of file coomatrix.h.

References cols_, and rows_.

Referenced by COOMatrix(), COOMatrix(), operator<<, operator=(), operator=(), and ~COOMatrix().

◆ COOMatrix() [2/3]

template<typename SCALAR>
lf::assemble::COOMatrix< SCALAR >::COOMatrix ( const COOMatrix< SCALAR > & )
default

References COOMatrix().

◆ COOMatrix() [3/3]

template<typename SCALAR>
lf::assemble::COOMatrix< SCALAR >::COOMatrix ( COOMatrix< SCALAR > && )
defaultnoexcept

References COOMatrix().

◆ ~COOMatrix()

template<typename SCALAR>
lf::assemble::COOMatrix< SCALAR >::~COOMatrix ( )
default

References COOMatrix().

Member Function Documentation

◆ AddToEntry()

template<typename SCALAR>
void lf::assemble::COOMatrix< SCALAR >::AddToEntry ( gdof_idx_t i,
gdof_idx_t j,
SCALAR increment )
inline

Add a value to the specified entry.

Parameters
irow index
jcolumn index
increment

Adds another index-value triplet.

Note
the size information for the matrix is adjusted to the passed indices. This ensures that the internally stored numbers of columns and rows are always bigger than the indices of any triplet. When manipulating the triplet information directly, one has to make sure that the size information remains valid.

Definition at line 87 of file coomatrix.h.

References cols_, rows_, and triplets_.

Referenced by lf::assemble::FixFlaggedSolutionCompAlt(), and lf::assemble::FixFlaggedSolutionComponents().

◆ cols()

template<typename SCALAR>
Index lf::assemble::COOMatrix< SCALAR >::cols ( ) const
inlinenodiscard

◆ makeDense()

template<typename SCALAR>
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > lf::assemble::COOMatrix< SCALAR >::makeDense ( ) const
inlinenodiscard

Create an Eigen::MatrixX from the COO format.

Returns
A dense matrix representing the COO matrix

This method is mainly meant for debugging

Definition at line 188 of file coomatrix.h.

References cols_, rows_, and triplets_.

◆ makeSparse()

template<typename SCALAR>
Eigen::SparseMatrix< Scalar > lf::assemble::COOMatrix< SCALAR >::makeSparse ( ) const
inlinenodiscard

Create an Eigen::SparseMatrix out of the COO format.

Returns
The created sparse matrix
Note
This method can be called multiple times and it does not modify the data stored in this COOMatrix.

Definition at line 172 of file coomatrix.h.

References cols_, rows_, and triplets_.

◆ MatVecMult() [1/2]

template<typename SCALAR>
template<typename VECTOR>
Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 > lf::assemble::COOMatrix< SCALAR >::MatVecMult ( SCALAR alpha,
const VECTOR & vec ) const
nodiscard

Computes the product of a (scaled) vector with the matrix in COO format.

Template Parameters
VECTORa basic vector type for the argument vector
Parameters
alphascalar with which to multiply the argument vector before the matrix x vector multiplication.
vecconstant reference to a vector of type VECTOR
Returns
result vector, a dense vector of Eigen

Requirements for type VECTOR

An object of type VECTOR must provide a method Size() telling the length of the vector and operator [] for read access to vector entries.

Definition at line 236 of file coomatrix.h.

References cols_, rows_, and triplets_.

Referenced by lf::assemble::FixFlaggedSolutionComponents().

◆ MatVecMult() [2/2]

template<typename SCALAR>
template<typename VECTOR, typename RESULTVECTOR>
void lf::assemble::COOMatrix< SCALAR >::MatVecMult ( SCALAR alpha,
const VECTOR & vec,
RESULTVECTOR & resvec ) const

In-situ computes the product of a vector with the matrix in COO format.

Template Parameters
VECTORa basic vector type for the argument vector
RESULTVECTORanother vector type for returning the result
Parameters
alphascalar with which to multiply the argument vector before the matrix x vector multiplication.
vecconstant reference to a vector of type VECTOR
resvecreference to an object of type RESULTVECTOR.
Note
the result of the matrix-vector product will be added to the entries of result!

Requirements for types VECTOR and RESULTVECTOR

An object of type VECTOR or RESULTVECTOR must provide a method Size() telling the length of the vector and operator [] for read/write access to vector entries.

Definition at line 250 of file coomatrix.h.

References cols_, rows_, and triplets_.

◆ operator=() [1/2]

template<typename SCALAR>
COOMatrix & lf::assemble::COOMatrix< SCALAR >::operator= ( const COOMatrix< SCALAR > & )
default

References COOMatrix().

◆ operator=() [2/2]

template<typename SCALAR>
COOMatrix & lf::assemble::COOMatrix< SCALAR >::operator= ( COOMatrix< SCALAR > && )
defaultnoexcept

References COOMatrix().

◆ PrintInfo()

template<typename SCALAR>
void lf::assemble::COOMatrix< SCALAR >::PrintInfo ( std::ostream & o) const
inline

Ouput of triplet list describing COO matrix.

Parameters
ooutput stream

Definition at line 205 of file coomatrix.h.

References cols_, rows_, and triplets_.

Referenced by operator<<.

◆ rows()

template<typename SCALAR>
Index lf::assemble::COOMatrix< SCALAR >::rows ( ) const
inlinenodiscard

◆ setZero() [1/2]

template<typename SCALAR>
void lf::assemble::COOMatrix< SCALAR >::setZero ( )
inline

Erase all entries of the matrix.

This method clears the vector of triplets, effectively setting the matrix to zero. It does not affect the size information about the matrix.

Definition at line 98 of file coomatrix.h.

References triplets_.

Referenced by lf::assemble::FixFlaggedSolutionCompAlt(), and lf::assemble::FixFlaggedSolutionComponents().

◆ setZero() [2/2]

template<typename SCALAR>
template<typename PREDICATE>
void lf::assemble::COOMatrix< SCALAR >::setZero ( PREDICATE && pred)
inline

Erase specific entries of the COO matrix, that is, set them to zero.

Template Parameters
PREDICATEa predicate type compliant with std::function<bool(gdof_idx_t,gdof_idx_t)>
Parameters
predselector predicate

Removes all triplets trp for which pred(trp.row(),trp.col()) is true. This amounts to setting the corresponding matrix entries to zero.

Definition at line 109 of file coomatrix.h.

References triplets_.

◆ triplets() [1/2]

template<typename SCALAR>
TripletVec & lf::assemble::COOMatrix< SCALAR >::triplets ( )
inlinenodiscard

Gives access to the vector of triplets.

Returns
reference to internal vector of triplets

Use of this method is deprecated. Use setZero(pred) and AddToEntry() instead.

Definition at line 124 of file coomatrix.h.

References triplets_.

◆ triplets() [2/2]

template<typename SCALAR>
const TripletVec & lf::assemble::COOMatrix< SCALAR >::triplets ( ) const
inlinenodiscard

Definition at line 125 of file coomatrix.h.

References triplets_.

Friends And Related Symbol Documentation

◆ operator<<

template<typename SCALAR>
template<typename SCALARTYPE>
std::ostream & operator<< ( std::ostream & o,
const COOMatrix< SCALARTYPE > & mat )
friend

Output operator for COO matrix.

This function prints matrix size and the list of triplets

Definition at line 229 of file coomatrix.h.

References COOMatrix(), and PrintInfo().

Member Data Documentation

◆ cols_

template<typename SCALAR>
size_type lf::assemble::COOMatrix< SCALAR >::cols_
private

dimensions of matrix

Definition at line 223 of file coomatrix.h.

Referenced by AddToEntry(), cols(), COOMatrix(), makeDense(), makeSparse(), MatVecMult(), MatVecMult(), and PrintInfo().

◆ rows_

template<typename SCALAR>
size_type lf::assemble::COOMatrix< SCALAR >::rows_
private

◆ triplets_

template<typename SCALAR>
TripletVec lf::assemble::COOMatrix< SCALAR >::triplets_ {}
private

COO format data

Definition at line 224 of file coomatrix.h.

Referenced by AddToEntry(), makeDense(), makeSparse(), MatVecMult(), MatVecMult(), PrintInfo(), setZero(), setZero(), triplets(), and triplets().


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