Discontinuous Galerkin Library
#include "dg/algorithm.h"
|
Classes | |
class | dg::Operator< T > |
A square nxn matrix. More... | |
Functions | |
template<class ContainerType > | |
auto | dg::asDenseMatrix (const std::vector< const ContainerType * > &in) |
Lightweight DenseMatrix for dg::blas2::gemv . More... | |
template<class ContainerType > | |
auto | dg::asDenseMatrix (const std::vector< const ContainerType * > &in, unsigned size) |
Lightweight DenseMatrix for dg::blas2::gemv . More... | |
template<class ContainerType > | |
std::vector< const ContainerType * > | dg::asPointers (const std::vector< ContainerType > &in) |
Convert a vector of vectors to a vector of pointers. More... | |
template<class T > | |
T | dg::create::lu_pivot (dg::Operator< T > &m, std::vector< unsigned > &p) |
LU Decomposition with partial pivoting. More... | |
template<class T > | |
void | dg::create::lu_solve (const dg::Operator< T > &lu, const std::vector< unsigned > &p, std::vector< T > &b) |
Solve the linear system with the LU decomposition. More... | |
template<class T > | |
dg::Operator< T > | dg::create::inverse (const dg::Operator< T > &in) |
Compute the inverse of a square matrix. More... | |
template<class T > | |
dg::Operator< T > | dg::invert (const dg::Operator< T > &in) |
Alias for dg::create::inverse . Compute inverse of square matrix. More... | |
auto dg::asDenseMatrix | ( | const std::vector< const ContainerType * > & | in | ) |
Lightweight DenseMatrix for dg::blas2::gemv
.
The philosophy is that a column matrix is represented by a std::vector
of pointers and can be multiplied with a coefficient vector
\[ \vec y = V \vec c = \sum_i c_i \vec v_{i} \]
where \( v_i\) are the columns of \( V\)
in | a collection of pointers that form the columns of the dense matrix |
std::vector
is to be interpreted as a dense matrix and call the correct implementation.ContainerType | Any class for which a specialization of TensorTraits exists and which fulfills the requirements of the there defined data and execution policies derived from AnyVectorTag and AnyPolicyTag . Among others
ContainerTypes in the argument list, then TensorTraits must exist for all of them |
auto dg::asDenseMatrix | ( | const std::vector< const ContainerType * > & | in, |
unsigned | size | ||
) |
Lightweight DenseMatrix for dg::blas2::gemv
.
The philosophy is that a column matrix is represented by a std::vector
of pointers and can be multiplied with a coefficient vector
\[ \vec y = V \vec c = \sum_i c_i \vec v_{i} \]
where \( v_i\) are the columns of \( V\)
in | a collection of pointers that form the columns of the dense matrix |
std::vector
is to be interpreted as a dense matrix and call the correct implementation.ContainerType | Any class for which a specialization of TensorTraits exists and which fulfills the requirements of the there defined data and execution policies derived from AnyVectorTag and AnyPolicyTag . Among others
ContainerTypes in the argument list, then TensorTraits must exist for all of them |
size | only the first size pointers are used in the matrix (i.e. the number of columns is size ) |
std::vector< const ContainerType * > dg::asPointers | ( | const std::vector< ContainerType > & | in | ) |
Convert a vector of vectors to a vector of pointers.
A convenience function that can be used in combination with asDenseMatrix
in | a collection of vectors that form the columns of the dense matrix |
ContainerType | Any class for which a specialization of TensorTraits exists and which fulfills the requirements of the there defined data and execution policies derived from AnyVectorTag and AnyPolicyTag . Among others
ContainerTypes in the argument list, then TensorTraits must exist for all of them |
dg::Operator< T > dg::create::inverse | ( | const dg::Operator< T > & | in | ) |
Compute the inverse of a square matrix.
using lu decomposition in combination with our accurate scalar products
T | value type |
in | input matrix |
std::runtime_error | if in is singular |
dg::Operator< T > dg::invert | ( | const dg::Operator< T > & | in | ) |
Alias for dg::create::inverse
. Compute inverse of square matrix.
using lu decomposition in combination with our accurate scalar products
T | value type |
in | input matrix |
std::runtime_error | if in is singular |
T dg::create::lu_pivot | ( | dg::Operator< T > & | m, |
std::vector< unsigned > & | p | ||
) |
LU Decomposition with partial pivoting.
T | value type |
std::runtime_error | if the matrix is singular |
m | contains lu decomposition of input on output (inplace transformation) |
p | contains the pivot elements on output |
m
dg::exblas
which makes it quite robust against almost singular matricesdg::create::lu_solve
void dg::create::lu_solve | ( | const dg::Operator< T > & | lu, |
const std::vector< unsigned > & | p, | ||
std::vector< T > & | b | ||
) |
Solve the linear system with the LU decomposition.
T | value type |
lu | result of dg::create::lu_pivot |
p | pivot vector from dg::create::lu_pivot |
b | right hand side (contains solution on output) |