15#ifndef dealii_sparse_direct_h
16#define dealii_sparse_direct_h
29#ifdef DEAL_II_WITH_UMFPACK
41#ifdef SuiteSparse_long
42 using suitesparse_index = SuiteSparse_long;
149 template <
class Matrix>
156 template <
class Matrix>
260 solve(
Vector<std::complex<double>> &rhs_and_solution,
283 template <
class Matrix>
285 solve(
const Matrix &matrix,
292 template <
class Matrix>
294 solve(
const Matrix &matrix,
295 Vector<std::complex<double>> &rhs_and_solution,
301 template <
class Matrix>
303 solve(
const Matrix &matrix,
310 template <
class Matrix>
312 solve(
const Matrix &matrix,
313 BlockVector<std::complex<double>> &rhs_and_solution,
329 <<
"UMFPACK routine " << arg1 <<
" returned error status " << arg2 <<
'.'
331 << (
"A complete list of error codes can be found in the file "
332 "<bundled/umfpack/UMFPACK/Include/umfpack.h>."
334 "That said, the two most common errors that can happen are "
335 "that your matrix cannot be factorized because it is "
336 "rank deficient, and that UMFPACK runs out of memory "
337 "because your problem is too large."
339 "The first of these cases most often happens if you "
340 "forget terms in your bilinear form necessary to ensure "
341 "that the matrix has full rank, or if your equation has a "
342 "spatially variable coefficient (or nonlinearity) that is "
343 "supposed to be strictly positive but, for whatever "
344 "reasons, is negative or zero. In either case, you probably "
345 "want to check your assembly procedure. Similarly, a "
346 "matrix can be rank deficient if you forgot to apply the "
347 "appropriate boundary conditions. For example, the "
348 "Laplace equation for a problem where only Neumann boundary "
349 "conditions are posed (or where you forget to apply Dirichlet "
350 "boundary conditions) has exactly one eigenvalue equal to zero "
351 "and its rank is therefore deficient by one. Finally, the matrix "
352 "may be rank deficient because you are using a quadrature "
353 "formula with too few quadrature points."
355 "The other common situation is that you run out of memory. "
356 "On a typical laptop or desktop, it should easily be possible "
357 "to solve problems with 100,000 unknowns in 2d. If you are "
358 "solving problems with many more unknowns than that, in "
359 "particular if you are in 3d, then you may be running out "
360 "of memory and you will need to consider iterative "
361 "solvers instead of the direct solver employed by "
396 template <
typename number>
400 template <
typename number>
404 template <
typename number>
419 std::vector<types::suitesparse_index>
Ap;
420 std::vector<types::suitesparse_index>
Ai;
421 std::vector<double>
Ax;
422 std::vector<double>
Az;
void * numeric_decomposition
~SparseDirectUMFPACK() override
void initialize(const SparsityPattern &sparsity_pattern)
void Tvmult(Vector< double > &dst, const Vector< double > &src) const
void * symbolic_decomposition
void solve(Vector< double > &rhs_and_solution, const bool transpose=false) const
void sort_arrays(const SparseMatrixEZ< number > &)
void factorize(const Matrix &matrix)
std::vector< double > control
void vmult(Vector< double > &dst, const Vector< double > &src) const
std::vector< types::suitesparse_index > Ap
std::vector< types::suitesparse_index > Ai
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcUMFPACKError(std::string arg1, int arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
unsigned int global_dof_index
long int suitesparse_index