HyperHDG
|
A namespace containing the different wrapper functions. More...
Classes | |
struct | LAPACKexception |
Exception to be thrown if LAPACK's solve fails. More... | |
Typedefs | |
template<unsigned int hyEdge_dim, unsigned int space_dim, typename TPCC::boundaries bndT = TPCC::boundaries::both, typename index_t = unsigned int> | |
using | tpcc_t = TPCC::Lexicographic< space_dim, hyEdge_dim, bndT, index_t, unsigned int, unsigned int > |
Type of a tensor product chain complex. More... | |
template<unsigned int hyEdge_dim, unsigned int space_dim> | |
using | tpcc_elem_t = TPCC::Element< space_dim, hyEdge_dim, unsigned int, unsigned int > |
Type of an element of a tensor product chain complex. More... | |
Functions | |
template<unsigned int system_size, unsigned int n_rhs_cols = 1, typename lapack_float_t > | |
std::array< lapack_float_t, system_size *n_rhs_cols > | lapack_solve (std::array< lapack_float_t, system_size *system_size > &dense_mat, std::array< lapack_float_t, system_size *n_rhs_cols > &rhs) |
Solve local system of equations. More... | |
template<unsigned int n_rows, unsigned int n_cols, typename lapack_float_t > | |
lapack_float_t | lapack_det (std::array< lapack_float_t, n_rows *n_cols > &dense_mat) |
Determinant of a rectangular system. More... | |
template<unsigned int n_rows, unsigned int n_cols, typename lapack_float_t > | |
std::array< lapack_float_t, n_rows *n_rows > | lapack_qr_decomp_q (std::array< lapack_float_t, n_rows *n_cols > &dense_mat) |
Matrix Q of QR decomposition. More... | |
template<unsigned int n_rows, unsigned int n_cols, typename lapack_float_t > | |
std::array< lapack_float_t, n_rows *n_cols > & | lapack_qr_decomp_r (std::array< lapack_float_t, n_rows *n_cols > &dense_mat) |
Matrix R of QR decomposition. More... | |
template<unsigned int n_rows, unsigned int n_cols, typename lapack_float_t > | |
void | lapack_qr_decomp (std::array< lapack_float_t, n_rows *n_cols > &dense_mat, std::array< lapack_float_t, n_rows *n_rows > &mat_q) |
Matrices Q and R of QR decomposition. More... | |
template<unsigned int n_rows, unsigned int n_cols, typename lapack_float_t > | |
void | lapack_qr_decomp (std::array< lapack_float_t, n_rows *n_cols > &dense_mat, std::array< lapack_float_t, n_rows *n_rows > &mat_q, std::array< lapack_float_t, n_cols *n_cols > &mat_r) |
Matrices Q and R of QR decomposition. More... | |
void | lapack_solve (int system_size, int n_rhs_cols, double *mat_a, double *rhs_b) |
Solve local system of equations with double floating point numbers — DO NOT USE. More... | |
void | lapack_solve (int system_size, int n_rhs_cols, float *mat_a, float *rhs_b) |
Solve local system of equations with float floating point numbers — DO NOT USE. More... | |
void | lapack_qr (int n_rows, int n_cols, double *mat_a, double *tau) |
QR decomposition in double floating point arithmetic — DO NOT USE. More... | |
void | lapack_qr (int n_rows, int n_cols, float *mat_a, float *tau) |
QR decomposition in float floating point arithmetic — DO NOT USE. More... | |
template<unsigned int n_rows, unsigned int n_cols, unsigned int rank, typename lapack_float_t > | |
std::array< lapack_float_t, n_rows *n_rows > | get_q_from_lapack_qr_result (const std::array< lapack_float_t, n_rows *n_cols > &dense_mat, const std::array< lapack_float_t, rank > &tau) |
Matrix Q of QR decomposition — DO NOT USE. More... | |
template<unsigned int n_rows, unsigned int n_cols, unsigned int rank, typename lapack_float_t > | |
void | get_q_from_lapack_qr_result (const std::array< lapack_float_t, n_rows *n_cols > &dense_mat, const std::array< lapack_float_t, rank > &tau, std::array< lapack_float_t, n_rows *n_rows > &mat_q) |
Matrix Q of QR decomposition — DO NOT USE. More... | |
template<unsigned int n_rows, unsigned int n_cols, typename lapack_float_t > | |
std::array< lapack_float_t, n_rows *n_cols > & | get_r_from_lapack_qr_result (std::array< lapack_float_t, n_rows *n_cols > &dense_mat) |
Matrix R of QR decomposition — DO NOT USE. More... | |
template<unsigned int n_rows, unsigned int n_cols, typename lapack_float_t > | |
void | get_r_from_lapack_qr_result (const std::array< lapack_float_t, n_rows *n_cols > &lapack_mat, std::array< lapack_float_t, n_cols *n_cols > &mat_r) |
Matrix R of QR decomposition — DO NOT USE. More... | |
unsigned int | binomial (unsigned int n, unsigned int k) |
Binomial coefficient of unsigned integers. More... | |
template<unsigned int hyEdge_dim, unsigned int space_dim, typename TPCC::boundaries bndT = TPCC::boundaries::both, typename index_t = unsigned int> | |
tpcc_t< hyEdge_dim, space_dim, bndT, index_t > | create_tpcc (const SmallVec< space_dim, index_t > &vec) |
Create a tensor product chain complex. More... | |
template<typename TPCC::boundaries bndT> | |
auto | tpcc_faces (const auto &elements) |
Create a tensor product chain complex associated to the facets. More... | |
template<typename index_t = unsigned int> | |
index_t | n_elements (const auto &tpcc) |
Return the element of given index the TPCC. More... | |
auto | get_element (const auto &tpcc, const auto index) |
Return the element of given index the TPCC. More... | |
template<typename index_t = unsigned int> | |
index_t | get_index (const auto &tpcc, const auto &elem) |
Return index of given element within TPCC. More... | |
template<typename index_t = unsigned int> | |
index_t | get_index_in_slice (const auto &tpcc, const auto &elem) |
Return index of given element within TPCC. More... | |
auto | get_face (const auto &elem, const unsigned int index) |
Return i-th element facet. More... | |
unsigned int | elem_orientation (const auto &elem) |
Determine the orientation of an element. More... | |
unsigned int | exterior_direction (const auto &elem, const unsigned int index) |
Return the index-th orthonormal direction to element. More... | |
unsigned int | exterior_coordinate (const auto &elem, const unsigned int index) |
Return coordinate value with respect to index-th orthonormal direction on element. More... | |
unsigned int | interior_coordinate (const auto &elem, const unsigned int index) |
Return coordinate value with respect to index-th orthonormal direction of element. More... | |
A namespace containing the different wrapper functions.
using Wrapper::tpcc_elem_t = typedef TPCC::Element<space_dim, hyEdge_dim, unsigned int, unsigned int> |
Type of an element of a tensor product chain complex.
using Wrapper::tpcc_t = typedef TPCC::Lexicographic<space_dim, hyEdge_dim, bndT, index_t, unsigned int, unsigned int> |
Type of a tensor product chain complex.
unsigned int Wrapper::binomial | ( | unsigned int | n, |
unsigned int | k | ||
) |
Binomial coefficient of unsigned integers.
tpcc_t<hyEdge_dim, space_dim, bndT, index_t> Wrapper::create_tpcc | ( | const SmallVec< space_dim, index_t > & | vec | ) |
Create a tensor product chain complex.
unsigned int Wrapper::elem_orientation | ( | const auto & | elem | ) |
Determine the orientation of an element.
unsigned int Wrapper::exterior_coordinate | ( | const auto & | elem, |
const unsigned int | index | ||
) |
Return coordinate value with respect to index-th orthonormal direction on element.
unsigned int Wrapper::exterior_direction | ( | const auto & | elem, |
const unsigned int | index | ||
) |
Return the index-th orthonormal direction to element.
auto Wrapper::get_element | ( | const auto & | tpcc, |
const auto | index | ||
) |
Return the element of given index the TPCC.
auto Wrapper::get_face | ( | const auto & | elem, |
const unsigned int | index | ||
) |
Return i-th element facet.
index_t Wrapper::get_index | ( | const auto & | tpcc, |
const auto & | elem | ||
) |
Return index of given element within TPCC.
index_t Wrapper::get_index_in_slice | ( | const auto & | tpcc, |
const auto & | elem | ||
) |
Return index of given element within TPCC.
|
inline |
Matrix Q of QR decomposition — DO NOT USE.
Return matrix Q of the Householder QR decomposition of the matrix. The matrix is supposed to contain the Householder QR decomposition which is encoded in the style of LAPACK. The function returns matrix Q in standard encoding.
n_rows | Number of rows of the matrix whose determinant should be calculated. |
n_cols | Number of columns of the matrix whose determinant should be calculated. |
float_t | Floating type which this function should be executed with. Only float and double are supported. |
dense_mat | Array comprising the QR decomposed matrix. |
tau | Auxiliary data from LAPACK to encode matrix Q. |
mat_q | Matrix Q of Householder QR decomposition. |
|
inline |
Matrix Q of QR decomposition — DO NOT USE.
Return matrix Q of the Householder QR decomposition of the matrix. The matrix is supposed to contain the Householder QR decomposition which is encoded in the style of LAPACK. The function returns matrix Q in standard encoding.
n_rows | Number of rows of the matrix whose determinant should be calculated. |
n_cols | Number of columns of the matrix whose determinant should be calculated. |
float_t | Floating type which this function should be executed with. Only float and double are supported. |
dense_mat | Array comprising the QR decomposed matrix. |
tau | Auxiliary data from LAPACK to encode matrix Q. |
mat_q | Matrix to be filled with entries of the return value. |
mat_q | Matrix Q of Householder QR decomposition. |
|
inline |
Matrix R of QR decomposition — DO NOT USE.
Return matrix R of the Householder QR decomposition of the matrix. The matrix is supposed to contain the Householder QR decomposition which is encoded in the style of LAPACK. The function returns matrix R in standard encoding, but with dimensions n_cols
times n_cols
.
n_rows | Number of rows of the matrix whose determinant should be calculated. |
n_cols | Number of columns of the matrix whose determinant should be calculated. |
float_t | Floating type which this function should be executed with. Only float and double are supported. |
lapack_mat | Array comprising the QR decomposed matrix. |
mat_r | Matrix to be filled with the entries of the return value. |
mat_r | Matrix R of Householder QR decomposition. |
|
inline |
Matrix R of QR decomposition — DO NOT USE.
Return matrix R of the Householder QR decomposition of the matrix. The matrix is supposed to contain the Householder QR decomposition which is encoded in the style of LAPACK. The function returns matrix R in standard encoding.
The function also manipulates the input matrix to be R.
n_rows | Number of rows of the matrix whose determinant should be calculated. |
n_cols | Number of columns of the matrix whose determinant should be calculated. |
float_t | Floating type which this function should be executed with. Only float and double are supported. |
dense_mat | Array comprising the QR decomposed matrix. |
mat_r | Matrix R of Householder QR decomposition. |
unsigned int Wrapper::interior_coordinate | ( | const auto & | elem, |
const unsigned int | index | ||
) |
Return coordinate value with respect to index-th orthonormal direction of element.
lapack_float_t Wrapper::lapack_det | ( | std::array< lapack_float_t, n_rows *n_cols > & | dense_mat | ) |
Determinant of a rectangular system.
Calculate the generalized determinant of a rectangular matrix. If the matrix is square, this is the standard determinant. The determinant is determined by doing a QR decomposition based on Householder transformations. Thus det(Q) = +1, if the number of rows if even and det(Q) = -1, if the number of rows is odd. This number is multiplied by the diagonal entries of R.
The matrix is destroyed using this function. Its entries might be used to generate matrices Q and R of the QR descomposition.
n_rows | Number of rows of the matrix whose determinant should be calculated. |
n_cols | Number of columns of the matrix whose determinant should be calculated. |
float_t | Floating type which this function should be executed with. Only float and double are supported. |
dense_mat | Array comprising the matrix describing the linear system of equations. |
determinant | Generalized determinant of the matrix. |
|
inline |
QR decomposition in double
floating point arithmetic — DO NOT USE.
Perform QR decomposition of a given matrix using Householder transformations. The matrices Q and R are stored according to the LAPACK function in an encoded way in matrix mat_a and vector tau.
n_rows | Number of rows of matrix. |
n_cols | Number of columns of matrix. |
mat_a | Matrix that is to be decomposed. |
tau | Space for return value, i.e. vector with auxiliary darta. |
mat_a | Encoded QR decomposition of the matrix. |
tau | Auxiliary parameters needed to reconstruct matrix Q. |
|
inline |
QR decomposition in float
floating point arithmetic — DO NOT USE.
Perform QR decomposition of a given matrix using Householder transformations. The matrices Q and R are stored according to the LAPACK function in an encoded way in matrix mat_a and vector tau.
n_rows | Number of rows of matrix. |
n_cols | Number of columns of matrix. |
mat_a | Matrix that is to be decomposed. |
tau | Space for return value, i.e. vector with auxiliary darta. |
mat_a | Encoded QR decomposition of the matrix. |
tau | Auxiliary parameters needed to reconstruct matrix Q. |
void Wrapper::lapack_qr_decomp | ( | std::array< lapack_float_t, n_rows *n_cols > & | dense_mat, |
std::array< lapack_float_t, n_rows *n_rows > & | mat_q | ||
) |
Matrices Q and R of QR decomposition.
Return matriices Q and R of the Householder QR decomposition of the matrix.
n_rows | Number of rows of the matrix whose determinant should be calculated. |
n_cols | Number of columns of the matrix whose determinant should be calculated. |
float_t | Floating type which this function should be executed with. Only float and double are supported. |
dense_mat | Array comprising the matrix describing the linear system of equations. |
mat_q | Reference to empy matrix to be filled. |
dense_mat | Matrix R of Householder QR decomposition. |
mat_q | Matrix Q of Householder QR decomposition. |
void Wrapper::lapack_qr_decomp | ( | std::array< lapack_float_t, n_rows *n_cols > & | dense_mat, |
std::array< lapack_float_t, n_rows *n_rows > & | mat_q, | ||
std::array< lapack_float_t, n_cols *n_cols > & | mat_r | ||
) |
Matrices Q and R of QR decomposition.
Return matriices Q and R of the Householder QR decomposition of the matrix.
n_rows | Number of rows of the matrix whose determinant should be calculated. |
n_cols | Number of columns of the matrix whose determinant should be calculated. |
float_t | Floating type which this function should be executed with. Only float and double are supported. |
dense_mat | Array comprising the matrix describing the linear system of equations. |
mat_q | Reference to empy matrix to be filled. |
mat_r | Reference to empy matrix to be filled. |
dense_mat | Matrix R of Householder QR decomposition. |
mat_q | Matrix Q of Householder QR decomposition. |
mat_r | Square system of size n_cols that contains respective part of R. |
std::array< lapack_float_t, n_rows *n_rows > Wrapper::lapack_qr_decomp_q | ( | std::array< lapack_float_t, n_rows *n_cols > & | dense_mat | ) |
Matrix Q of QR decomposition.
Return matrix Q of the Householder QR decomposition of the matrix. The matrix is destroyed using this function. Its entries might be used to generate matrices Q and R of the QR descomposition.
n_rows | Number of rows of the matrix whose determinant should be calculated. |
n_cols | Number of columns of the matrix whose determinant should be calculated. |
float_t | Floating type which this function should be executed with. Only float and double are supported. |
dense_mat | Array comprising the matrix describing the linear system of equations. |
mat_q | Matrix Q of Householder QR decomposition. |
std::array< lapack_float_t, n_rows *n_cols > & Wrapper::lapack_qr_decomp_r | ( | std::array< lapack_float_t, n_rows *n_cols > & | dense_mat | ) |
Matrix R of QR decomposition.
Return matrix R of the Householder QR decomposition of the matrix. The matrix is destroyed using this function. Its entries are the same as the ones of R of the QR descomposition.
n_rows | Number of rows of the matrix whose determinant should be calculated. |
n_cols | Number of columns of the matrix whose determinant should be calculated. |
float_t | Floating type which this function should be executed with. Only float and double are supported. |
dense_mat | Array comprising the matrix describing the linear system of equations. |
dense_mat | Matrix R of Householder QR decomposition. |
mat_r | Matrix R of Householder QR decomposition. |
|
inline |
Solve local system of equations with double
floating point numbers — DO NOT USE.
Solve linear (dense) system of equations \(Ax=b\), where $A$ is an \(n \times n\) square matrix, which enters as a pointer to double
, \(n\) is provided via system_size
, and the double
pointer rhs_b
is both, input (i.e., \(b\)) and output (i.e., \(x\)) of the function.
Independent of const
expressions of functions using lapack_solve
one should not use mat_a
after calling this function. The input that has been within rhs_b
will have been replaced by the solution of the system of equations.
system_size | Size of the system of equations. |
n_rhs_cols | Number of right-hand sides the system is solved for. |
mat_a | Pointer to the matrix describing the linear system of equations. |
rhs_b | Pointer to the right-hand side of the system. |
rhs_b | Pointer to the solution of the system of equations. |
|
inline |
Solve local system of equations with float
floating point numbers — DO NOT USE.
Solve linear (dense) system of equations \(Ax=b\), where $A$ is an \(n \times n\) square matrix, which enters as a pointer to float
, \(n\) is provided via system_size
, and the float
pointer rhs_b
is both, input (i.e., \(b\)) and output (i.e., \(x\)) of the function.
Independent of const
expressions of functions using lapack_solve
one should not use mat_a
after calling this function. The input that has been within rhs_b
will have been replaced by the solution of the system of equations.
system_size | Size of the system of equations. |
n_rhs_cols | Number of right-hand sides the system is solved for. |
mat_a | Pointer to the matrix describing the linear system of equations. |
rhs_b | Pointer to the right-hand side of the system. |
rhs_b | Pointer to the solution of the system of equations. |
std::array< lapack_float_t, system_size *n_rhs_cols > Wrapper::lapack_solve | ( | std::array< lapack_float_t, system_size *system_size > & | dense_mat, |
std::array< lapack_float_t, system_size *n_rhs_cols > & | rhs | ||
) |
Solve local system of equations.
Solve linear (dense) system of equations \(Ax=b\), where $A$ is an \(n \times n\) square matrix, which enters as a std::array
of double
, \(n\) is provided via system_size
, and the std::array
of double
rhs_b
is both, input (i.e., \(b\)) and output (i.e., \(x\)) of the function.
Independent of const
expressions of functions using lapack_solve
one should not use mat_a
after calling this function. The input that has been within rhs_b
will have been replaced by the solution of the system of equations.
system_size | Size of the system of equations. |
n_rhs_cols | Number of columns of the right hand side matrix. Defaults to 1. |
float_t | Floating type which this function should be executed with. Only float and double are supported. |
dense_mat | Array comprising the matrix describing the linear system of equations. |
rhs | Array comprising the right-hand side of the system. |
rhs_b | Array comprising the solution of the system of equations. |
index_t Wrapper::n_elements | ( | const auto & | tpcc | ) |
Return the element of given index the TPCC.
auto Wrapper::tpcc_faces | ( | const auto & | elements | ) |
Create a tensor product chain complex associated to the facets.