OpenBEM
Open-source framework for electromagnetic simulation with the boundary element method.
Loading...
Searching...
No Matches
bem::MatrixBase< T > Class Template Referenceabstract

Base class for for matrix algebra containers. More...

#include <base.hpp>

+ Inheritance diagram for bem::MatrixBase< T >:

Detailed Description

template<typename T>
class bem::MatrixBase< T >

Base class for for matrix algebra containers.

Template Parameters
T- Data type to be stored in the matrix (e.g., float, double, std::complex).

Definition at line 45 of file base.hpp.

Member Function Documentation

◆ num_rows()

◆ num_cols()

◆ resize()

template<typename T >
virtual void bem::MatrixBase< T >::resize ( Index  rows,
Index  cols 
)
pure virtual

Resizes the matrix to a new number of rows and columns.

Parameters
[in]rows- New number of rows.
[in]cols- New number of columns.

Implemented in bem::EigenMatrixBase< T, MatrixType >, bem::EigenMatrixBase< T, Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor > >, and bem::EigenMatrixBase< T, Eigen::SparseMatrix< T, Eigen::ColMajor > >.

◆ clear()

◆ value()

template<typename T >
virtual T bem::MatrixBase< T >::value ( Index  row,
Index  col 
) const
pure virtual

Returns the matrix value at the specified row and column.

Parameters
[in]row- Row index.
[in]col- Column index.
Returns
Value at the specified position in the matrix.

Implemented in bem::EigenDenseMatrix< T, StorageOrder >, and bem::EigenSparseMatrix< T, StorageOrder >.

◆ set_value()

template<typename T >
virtual void bem::MatrixBase< T >::set_value ( Index  row,
Index  col,
const T a 
)
pure virtual

Sets the matrix value at the specified row and column.

Parameters
[in]row- Row index.
[in]col- Column index.
[in]a- Value to set.

Implemented in bem::EigenDenseMatrix< T, StorageOrder >, and bem::EigenSparseMatrix< T, StorageOrder >.

◆ add_ax()

template<typename T >
virtual void bem::MatrixBase< T >::add_ax ( const MatrixBase< T > &  x,
const T a = T(1) 
)
pure virtual

Computes \( \mathbf{M} = \mathbf{M} + a\mathbf{X} \) where \( \mathbf{M} \) is this matrix, \( a \) is a scalar, and \( \mathbf{X} \) is a matrix.

Parameters
[in]x- Matrix to scale and add, must have the same dimensions as this matrix.
[in]a- Scalar which which to scale x.

Implemented in bem::EigenMatrixBase< T, MatrixType >, bem::EigenMatrixBase< T, Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor > >, and bem::EigenMatrixBase< T, Eigen::SparseMatrix< T, Eigen::ColMajor > >.

◆ scale()

template<typename T >
virtual void bem::MatrixBase< T >::scale ( const T a)
pure virtual

Scales all matrix entries by a given value.

Parameters
[in]a- Value by which to scale.

Implemented in bem::EigenDenseMatrix< T, StorageOrder >, and bem::EigenSparseMatrix< T, StorageOrder >.

◆ set_zero()

◆ data() [1/2]

template<typename T >
virtual const T * bem::MatrixBase< T >::data ( ) const
pure virtual

Returns a read-only pointer to the underlying raw data.

Returns
Read-only pointer to the raw data.

Implemented in bem::EigenDenseMatrix< T, StorageOrder >, and bem::EigenSparseMatrix< T, StorageOrder >.

◆ data() [2/2]

template<typename T >
virtual T * bem::MatrixBase< T >::data ( )
pure virtual

Returns a writable pointer to the underlying raw data - use with care.

Returns
Writable pointer to the raw data.

Implemented in bem::EigenDenseMatrix< T, StorageOrder >, and bem::EigenSparseMatrix< T, StorageOrder >.

◆ mat_solve()

template<typename T >
virtual void bem::MatrixBase< T >::mat_solve ( MatrixBase< T > &  x,
const MatrixBase< T > &  b 
)
pure virtual

Solves \( \mathbf{M}\mathbf{X} = \mathbf{B} \) for matrix \( \mathbf{X} \) with a direct solver, where \( \mathbf{M} \) is this matrix, and \( \mathbf{B} \) is a given right-hand side matrix.

Parameters
[out]x- Solution.
[in]b- Right-hand side matrix, must have the same number of rows as this matrix.

Implemented in bem::EigenDenseMatrix< T, StorageOrder >, and bem::EigenSparseMatrix< T, StorageOrder >.

◆ size()

template<typename T >
virtual Index bem::MatrixBase< T >::size ( ) const
inlinevirtual

◆ add_value()

template<typename T >
virtual void bem::MatrixBase< T >::add_value ( Index  row,
Index  col,
const T a 
)
inlinevirtual

Adds to the matrix value at the specified row and column.

Parameters
[in]row- Row index.
[in]col- Column index.
[in]a- Value to add.

Reimplemented in bem::EigenDenseMatrix< T, StorageOrder >, and bem::EigenSparseMatrix< T, StorageOrder >.

Definition at line 154 of file base.hpp.

◆ preallocate() [1/2]

template<typename T >
virtual void bem::MatrixBase< T >::preallocate ( const std::vector< Index > &  nnz)
inlinevirtual

Preallocates memory for a given number of non-zero values per row.

Parameters
[in]nnz- Vector containing the number of non-zero values for each row.

Reimplemented in bem::EigenSparseMatrix< T, StorageOrder >.

Definition at line 165 of file base.hpp.

◆ preallocate() [2/2]

template<typename T >
virtual void bem::MatrixBase< T >::preallocate ( Index  num_entries)
inlinevirtual

Preallocates memory for a given total number of non-zero values.

Parameters
[in]num_entries- Total number of non-zero values.

Reimplemented in bem::EigenSparseMatrix< T, StorageOrder >.

Definition at line 172 of file base.hpp.

◆ assemble()

template<typename T >
virtual void bem::MatrixBase< T >::assemble ( )
inlinevirtual

Assembles cached data into the matrix.

Reimplemented in bem::EigenSparseMatrix< T, StorageOrder >.

Definition at line 178 of file base.hpp.

◆ set_identity()

template<typename T >
virtual void bem::MatrixBase< T >::set_identity ( )
inlinevirtual

◆ write_binary()

template<typename T >
virtual void bem::MatrixBase< T >::write_binary ( const std::string  name = "matrix") const
inlinevirtual

Writes the matrix to a binary file.

Definition at line 196 of file base.hpp.

◆ read_binary()

template<typename T >
virtual void bem::MatrixBase< T >::read_binary ( const std::string  name = "matrix")
inlinevirtual

Reads the matrix from a binary file.

Definition at line 211 of file base.hpp.

◆ print()

template<typename T >
virtual void bem::MatrixBase< T >::print ( const std::string  name = "matrix") const
inlinevirtual

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