16 #ifndef dealii_sparse_direct_h 17 #define dealii_sparse_direct_h 20 #include <deal.II/base/config.h> 22 #include <deal.II/base/exceptions.h> 23 #include <deal.II/base/subscriptor.h> 24 #include <deal.II/base/thread_management.h> 26 #include <deal.II/lac/block_sparse_matrix.h> 27 #include <deal.II/lac/sparse_matrix.h> 28 #include <deal.II/lac/sparse_matrix_ez.h> 29 #include <deal.II/lac/vector.h> 31 #ifdef DEAL_II_WITH_UMFPACK 35 DEAL_II_NAMESPACE_OPEN
43 #ifdef SuiteSparse_long 152 template <
class Matrix>
159 template <
class Matrix>
262 template <
class Matrix>
264 solve(
const Matrix & matrix,
271 template <
class Matrix>
273 solve(
const Matrix & matrix,
290 <<
"UMFPACK routine " << arg1 <<
" returned error status " << arg2 <<
"." 292 << (
"A complete list of error codes can be found in the file " 293 "<bundled/umfpack/UMFPACK/Include/umfpack.h>." 295 "That said, the two most common errors that can happen are " 296 "that your matrix cannot be factorized because it is " 297 "rank deficient, and that UMFPACK runs out of memory " 298 "because your problem is too large." 300 "The first of these cases most often happens if you " 301 "forget terms in your bilinear form necessary to ensure " 302 "that the matrix has full rank, or if your equation has a " 303 "spatially variable coefficient (or nonlinearity) that is " 304 "supposed to be strictly positive but, for whatever " 305 "reasons, is negative or zero. In either case, you probably " 306 "want to check your assembly procedure. Similarly, a " 307 "matrix can be rank deficient if you forgot to apply the " 308 "appropriate boundary conditions. For example, the " 309 "Laplace equation without boundary conditions has a " 310 "single zero eigenvalue and its rank is therefore " 313 "The other common situation is that you run out of memory. " 314 "On a typical laptop or desktop, it should easily be possible " 315 "to solve problems with 100,000 unknowns in 2d. If you are " 316 "solving problems with many more unknowns than that, in " 317 "particular if you are in 3d, then you may be running out " 318 "of memory and you will need to consider iterative " 319 "solvers instead of the direct solver employed by " 339 void *numeric_decomposition;
353 template <
typename number>
357 template <
typename number>
361 template <
typename number>
368 std::vector<types::suitesparse_index>
Ap;
369 std::vector<types::suitesparse_index> Ai;
370 std::vector<double> Ax;
378 DEAL_II_NAMESPACE_CLOSE
380 #endif // dealii_sparse_direct_h static ::ExceptionBase & ExcUMFPACKError(std::string arg1, int arg2)
void * symbolic_decomposition
#define DeclException2(Exception2, type1, type2, outsequence)
long int suitesparse_index
void factorize(const Matrix &matrix)
void vmult(Vector< double > &dst, const Vector< double > &src) const
void solve(Vector< double > &rhs_and_solution, const bool transpose=false) const
types::global_dof_index size_type
SymmetricTensor< rank_, dim, Number > transpose(const SymmetricTensor< rank_, dim, Number > &t)
std::vector< double > control
unsigned int global_dof_index
void initialize(const SparsityPattern &sparsity_pattern)
std::vector< types::suitesparse_index > Ap
~SparseDirectUMFPACK() override
void Tvmult(Vector< double > &dst, const Vector< double > &src) const
void sort_arrays(const SparseMatrixEZ< number > &)