16 #ifndef dealii__sparse_direct_h 17 #define dealii__sparse_direct_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/exceptions.h> 22 #include <deal.II/base/subscriptor.h> 23 #include <deal.II/base/thread_management.h> 24 #include <deal.II/lac/vector.h> 25 #include <deal.II/lac/sparse_matrix.h> 26 #include <deal.II/lac/sparse_matrix_ez.h> 27 #include <deal.II/lac/block_sparse_matrix.h> 29 #ifdef DEAL_II_WITH_UMFPACK 32 #ifndef SuiteSparse_long 33 # define SuiteSparse_long long int 36 DEAL_II_NAMESPACE_OPEN
139 template <
class Matrix>
145 template <
class Matrix>
230 const bool transpose =
false)
const;
236 const bool transpose =
false)
const;
244 template <
class Matrix>
245 void solve (
const Matrix &matrix,
247 const bool transpose =
false);
252 template <
class Matrix>
253 void solve (
const Matrix &matrix,
255 const bool transpose =
false);
267 <<
"UMFPACK routine " << arg1
268 <<
" returned error status " << arg2 <<
"." 270 << (
"A complete list of error codes can be found in the file " 271 "<bundled/umfpack/UMFPACK/Include/umfpack.h>." 273 "That said, the two most common errors that can happen are " 274 "that your matrix cannot be factorized because it is " 275 "rank deficient, and that UMFPACK runs out of memory " 276 "because your problem is too large." 278 "The first of these cases most often happens if you " 279 "forget terms in your bilinear form necessary to ensure " 280 "that the matrix has full rank, or if your equation has a " 281 "spatially variable coefficient (or nonlinearity) that is " 282 "supposed to be strictly positive but, for whatever " 283 "reasons, is negative or zero. In either case, you probably " 284 "want to check your assembly procedure. Similarly, a " 285 "matrix can be rank deficient if you forgot to apply the " 286 "appropriate boundary conditions. For example, the " 287 "Laplace equation without boundary conditions has a " 288 "single zero eigenvalue and its rank is therefore " 291 "The other common situation is that you run out of memory." 292 "On a typical laptop or desktop, it should easily be possible " 293 "to solve problems with 100,000 unknowns in 2d. If you are " 294 "solving problems with many more unknowns than that, in " 295 "particular if you are in 3d, then you may be running out " 296 "of memory and you will need to consider iterative " 297 "solvers instead of the direct solver employed by " 317 void *numeric_decomposition;
330 template <
typename number>
333 template <
typename number>
336 template <
typename number>
344 std::vector<SuiteSparse_long>
Ap;
345 std::vector<SuiteSparse_long> Ai;
346 std::vector<double> Ax;
354 DEAL_II_NAMESPACE_CLOSE
356 #endif // dealii__sparse_direct_h void * symbolic_decomposition
#define DeclException2(Exception2, type1, type2, outsequence)
types::global_dof_index size_type
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
unsigned int global_dof_index
std::vector< SuiteSparse_long > Ap
static ::ExceptionBase & ExcUMFPACKError(char *arg1, int arg2)
std::vector< double > control
void initialize(const SparsityPattern &sparsity_pattern)
void Tvmult(Vector< double > &dst, const Vector< double > &src) const
void sort_arrays(const SparseMatrixEZ< number > &)