|
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\}}\)
Go to the documentation of this file.
16 #ifndef dealii_solver_richardson_h
17 #define dealii_solver_richardson_h
62 template <
class VectorType = Vector<
double>>
110 template <
typename MatrixType,
typename PreconditionerType>
112 solve(
const MatrixType &
A,
115 const PreconditionerType &preconditioner);
120 template <
typename MatrixType,
typename PreconditionerType>
125 const PreconditionerType &preconditioner);
151 virtual typename VectorType::value_type
165 template <
class VectorType>
168 const bool use_preconditioned_residual)
170 , use_preconditioned_residual(use_preconditioned_residual)
174 template <
class VectorType>
177 const AdditionalData & data)
179 , additional_data(data)
184 template <
class VectorType>
186 const AdditionalData &data)
188 , additional_data(data)
193 template <
class VectorType>
194 template <
typename MatrixType,
typename PreconditionerType>
199 const PreconditionerType &preconditioner)
205 unsigned int iter = 0;
227 preconditioner.vmult(
d, r);
231 last_criterion = criterion(r,
d);
232 conv = this->iteration_status(iter, last_criterion, x);
236 x.add(additional_data.omega,
d);
237 print_vectors(iter, x, r,
d);
250 template <
class VectorType>
251 template <
typename MatrixType,
typename PreconditionerType>
256 const PreconditionerType &preconditioner)
261 unsigned int iter = 0;
283 preconditioner.Tvmult(
d, r);
285 last_criterion = criterion(r,
d);
286 conv = this->iteration_status(iter, last_criterion, x);
290 x.add(additional_data.omega,
d);
291 print_vectors(iter, x, r,
d);
304 template <
class VectorType>
314 template <
class VectorType>
315 inline typename VectorType::value_type
319 if (!additional_data.use_preconditioned_residual)
326 template <
class VectorType>
330 additional_data.omega = om;
bool use_preconditioned_residual
virtual VectorType::value_type criterion(const VectorType &r, const VectorType &d) const
virtual ~SolverRichardson() override=default
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
@ iterate
Continue iteration.
#define DEAL_II_NAMESPACE_OPEN
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
AdditionalData additional_data
void set_omega(const double om=1.)
@ success
Stop iteration, goal reached.
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
AdditionalData(const double omega=1, const bool use_preconditioned_residual=false)
#define DEAL_II_NAMESPACE_CLOSE
void Tsolve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
#define AssertThrow(cond, exc)
SolverRichardson(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
T max(const T &t, const MPI_Comm &mpi_communicator)