Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
sparse_ilu.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2018 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #ifndef dealii_sparse_ilu_h
18 #define dealii_sparse_ilu_h
19 
20 
21 #include <deal.II/base/config.h>
22 
23 #include <deal.II/lac/exceptions.h>
26 
28 
60 template <typename number>
61 class SparseILU : public SparseLUDecomposition<number>
62 {
63 public:
68 
75  SparseILU() = default;
76 
82 
100  template <typename somenumber>
101  void
103  const AdditionalData & parameters = AdditionalData());
104 
111  template <typename somenumber>
112  void
113  vmult(Vector<somenumber> &dst, const Vector<somenumber> &src) const;
114 
115 
122  template <typename somenumber>
123  void
124  Tvmult(Vector<somenumber> &dst, const Vector<somenumber> &src) const;
125 
126 
131  std::size_t
132  memory_consumption() const override;
133 
143  double,
144  << "The strengthening parameter " << arg1
145  << " is not greater or equal than zero!");
150  size_type,
151  << "While computing the ILU decomposition, the algorithm "
152  "found a zero pivot on the diagonal of row "
153  << arg1
154  << ". This must stop the ILU algorithm because it means "
155  "that the matrix for which you try to compute a "
156  "decomposition is singular.");
158 };
159 
161 //---------------------------------------------------------------------------
162 
163 
165 
166 #endif // dealii_sparse_ilu_h
sparse_matrix.h
SparseILU::Tvmult
void Tvmult(Vector< somenumber > &dst, const Vector< somenumber > &src) const
SparseMatrix
Definition: sparse_matrix.h:497
SparseLUDecomposition::AdditionalData
Definition: sparse_decomposition.h:144
SparseILU::SparseILU
SparseILU()=default
exceptions.h
sparse_decomposition.h
SparseILU
Definition: sparse_ilu.h:61
SparseLUDecomposition
Definition: sparse_decomposition.h:110
SparseLUDecomposition::size_type
typename SparseMatrix< number >::size_type size_type
Definition: sparse_decomposition.h:126
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
LAPACKSupport::matrix
@ matrix
Contents is actually a matrix.
Definition: lapack_support.h:60
SparseILU::memory_consumption
std::size_t memory_consumption() const override
SparseILU::ExcInvalidStrengthening
static ::ExceptionBase & ExcInvalidStrengthening(double arg1)
SparseILU::initialize
void initialize(const SparseMatrix< somenumber > &matrix, const AdditionalData &parameters=AdditionalData())
DeclException1
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:518
SparseILU::AdditionalData
typename SparseLUDecomposition< number >::AdditionalData AdditionalData
Definition: sparse_ilu.h:81
SparseILU::vmult
void vmult(Vector< somenumber > &dst, const Vector< somenumber > &src) const
config.h
SparseILU::ExcZeroPivot
static ::ExceptionBase & ExcZeroPivot(size_type arg1)
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
Vector
Definition: mapping_q1_eulerian.h:32