Matrix Formats
We provide different frontends for high-performance linear algebra packages like Petsc and Slepc and Eigen. Respective functions are available if the libraries are included.
For dense operations (e.g., determining the full spectrum), a dense matrix with Eigen is returned by create_EigenDenseMatrix(…). Eigen can be combined using different LAPACK libraries such as Intel’s MKL. For sparse operations (e.g., Krylov space methods), we recommend using Petsc and Slepc. The libraries exploit the sparsity of the matrices and utilize the full power of a computing cluster by distributing tasks over several nodes. create_PetscSparseMatrix(…) returns the sparse matrix that can be used for different routines in Slepc. Equivalently, create_PetscShellMatrix(…) returns a matrix shell allowing matrix-free multiplications. Typically, the memory consumption is reduced by the system size but with increased run time.
Furthermore, create_DenseMatrix() and create_SparseMatrix() return the matrix in a dense format using std::vector
and a sparse matrix using the SparseMatrix class.
The functions are available without Petsc or Eigen.
All return functions are parallelized – if possible. Functions using native cpp container and Eigen work with OpenMP. The Petsc frontends use MPI.
-
ShellMatrix<Operator<T, ScalarType, OT>> danceq::internal::Operator::create_ShellMatrix(const bool show_info = false, const double overhead_in_GB_per_core = 2.) const
Returns a ShellMatrix for matrix-free multiplication.
Creates a ShellMatrix that allows a matrix-free matrix-vector multiplication. The ShellMatrix works together with the Vector class. An instance of the Vector class with the correct memory layout can be retrieved using the function create_Vector(). Importantly, the ShellMatrix (and Vector class) is compatible with MPI and openMP. The syntax is equivalent in all cases. Additional memory for the communication between processes is requiered whe MPI is enabled. This can be limited by overhead_in_GB_per_core.
auto matrix = H.create_ShellMatrix();
- Parameters:
show_info – Prints information like communication pattern used by MPI if
true
(Default:true
)overhead_in_GB_per_core – Maximal memory per core in GB used during the multiplication (Default: 2.)
- Returns:
matrix
-
Mat danceq::internal::Operator::create_PetscSparseMatrix(const bool show_info = false)
Returns a Petsc Mat as a MATMPIAIJ type.
Creates and fills a sparse matrix over several MPI processes using in a MATMPIAIJ format from Petsc. This functions is a wrapper for create_PetscSparseMatrix_reference(…).
auto matrix = H.create_PetscSparseMatrix();
- Parameters:
show_info – Prints information like memory usage if
true
- Returns:
matrix
-
Mat danceq::internal::Operator::create_PetscShellMatrix(const bool show_info = false, const double overhead_in_GB_per_core = 2.) const
Returns a Petsc Mat as a MATSHELL type for matrix-free multiplication.
Creates a shell matrix over several MPI processes using in a MATSHELL format from Petsc. It is possible to control the maximal memory overhead per process during the multiplication using overhead_in_GB_per_core. More detailed information can be found here. This functions is a wrapper for create_PetscShellMatrix_reference(…).
auto matrix = H.create_PetscShellMatrix();
- Parameters:
show_info – Prints information like communication pattern used by MPI if
true
(Default:false
)overhead_in_GB_per_core – Maximal memory per core in GB used during the multiplication (Default: 2.)
- Returns:
matrix
-
EigenMatrixType danceq::internal::Operator::create_EigenDenseMatrix(void) const
Returns a dense matix with Eigen.
Creates and fills a dense matrix with Eigen:
auto matrix = H.create_EigenDenseMatrix();
- Returns:
matrix
-
std::vector<std::vector<ScalarType>> danceq::internal::Operator::create_DenseMatrix(void) const
Returns the operator in a dense format.
The operator is returned as a dense matrix in the form of a two-dimensional
std::vector
of size dim times dim. The function is compatible with openMP.auto matrix = H.create_DenseMatrix();
- Returns:
error_code
-
SparseMatrix<ScalarType> danceq::internal::Operator::create_SparseMatrix(void)
Returns the matrix in a sparse format.
The operator is returned as a SparseMatrix. The function is compatible with openMP. We do not provide a MPI version as this should be done using Petsc via create_PetscSparseMatrix(…).
auto matrix = H.create_SparseMatrix();
- Returns:
matrix