Reference documentation for deal.II version 8.5.1
sparse_direct.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 2016 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__sparse_direct_h
17 #define dealii__sparse_direct_h
18 
19 
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>
28 
29 #ifdef DEAL_II_WITH_UMFPACK
30 # include <umfpack.h>
31 #endif
32 #ifndef SuiteSparse_long
33 # define SuiteSparse_long long int
34 #endif
35 
36 DEAL_II_NAMESPACE_OPEN
37 
83 {
84 public:
89 
95  {};
96 
97 
103 
108 
120  void initialize (const SparsityPattern &sparsity_pattern);
121 
139  template <class Matrix>
140  void factorize (const Matrix &matrix);
141 
145  template <class Matrix>
146  void initialize(const Matrix &matrix,
147  const AdditionalData additional_data = AdditionalData());
148 
169  void vmult (Vector<double> &dst,
170  const Vector<double> &src) const;
171 
175  void vmult (BlockVector<double> &dst,
176  const BlockVector<double> &src) const;
177 
182  void Tvmult (Vector<double> &dst,
183  const Vector<double> &src) const;
184 
188  void Tvmult (BlockVector<double> &dst,
189  const BlockVector<double> &src) const;
190 
195  size_type m () const;
196 
201  size_type n () const;
202 
229  void solve (Vector<double> &rhs_and_solution,
230  const bool transpose = false) const;
231 
235  void solve (BlockVector<double> &rhs_and_solution,
236  const bool transpose = false) const;
237 
244  template <class Matrix>
245  void solve (const Matrix &matrix,
246  Vector<double> &rhs_and_solution,
247  const bool transpose = false);
248 
252  template <class Matrix>
253  void solve (const Matrix &matrix,
254  BlockVector<double> &rhs_and_solution,
255  const bool transpose = false);
256 
266  DeclException2 (ExcUMFPACKError, char *, int,
267  << "UMFPACK routine " << arg1
268  << " returned error status " << arg2 << "."
269  << "\n\n"
270  << ("A complete list of error codes can be found in the file "
271  "<bundled/umfpack/UMFPACK/Include/umfpack.h>."
272  "\n\n"
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."
277  "\n\n"
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 "
289  "deficient by one."
290  "\n\n"
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 "
298  "UMFPACK."));
299 
300 private:
305 
310 
317  void *numeric_decomposition;
318 
322  void clear ();
323 
330  template <typename number>
331  void sort_arrays (const SparseMatrixEZ<number> &);
332 
333  template <typename number>
334  void sort_arrays (const SparseMatrix<number> &);
335 
336  template <typename number>
337  void sort_arrays (const BlockSparseMatrix<number> &);
338 
344  std::vector<SuiteSparse_long> Ap;
345  std::vector<SuiteSparse_long> Ai;
346  std::vector<double> Ax;
347 
351  std::vector<double> control;
352 };
353 
354 DEAL_II_NAMESPACE_CLOSE
355 
356 #endif // dealii__sparse_direct_h
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
types::global_dof_index size_type
Definition: sparse_direct.h:88
size_type n() const
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
Definition: types.h:88
std::vector< SuiteSparse_long > Ap
static ::ExceptionBase & ExcUMFPACKError(char *arg1, int arg2)
std::vector< double > control
void initialize(const SparsityPattern &sparsity_pattern)
size_type m() const
void Tvmult(Vector< double > &dst, const Vector< double > &src) const
void sort_arrays(const SparseMatrixEZ< number > &)