16#ifndef dealii_solver_richardson_h
17#define dealii_solver_richardson_h
59template <
class VectorType = Vector<
double>>
107 template <
typename MatrixType,
typename PreconditionerType>
111 const VectorType & b,
112 const PreconditionerType &preconditioner);
117 template <
typename MatrixType,
typename PreconditionerType>
121 const VectorType & b,
122 const PreconditionerType &preconditioner);
137 const VectorType & x,
138 const VectorType & r,
139 const VectorType & d)
const;
148 virtual typename VectorType::value_type
149 criterion(
const VectorType &r,
const VectorType &d)
const;
162template <
class VectorType>
165 const bool use_preconditioned_residual)
167 , use_preconditioned_residual(use_preconditioned_residual)
171template <
class VectorType>
174 const AdditionalData & data)
176 , additional_data(data)
181template <
class VectorType>
183 const AdditionalData &data)
185 , additional_data(data)
190template <
class VectorType>
191template <
typename MatrixType,
typename PreconditionerType>
195 const VectorType & b,
196 const PreconditionerType &preconditioner)
200 double last_criterion = std::numeric_limits<double>::lowest();
202 unsigned int iter = 0;
224 preconditioner.vmult(d, r);
228 last_criterion = criterion(r, d);
229 conv = this->iteration_status(iter, last_criterion, x);
233 x.add(additional_data.omega, d);
234 print_vectors(iter, x, r, d);
247template <
class VectorType>
248template <
typename MatrixType,
typename PreconditionerType>
252 const VectorType & b,
253 const PreconditionerType &preconditioner)
256 double last_criterion = std::numeric_limits<double>::lowest();
258 unsigned int iter = 0;
280 preconditioner.Tvmult(d, r);
282 last_criterion = criterion(r, d);
283 conv = this->iteration_status(iter, last_criterion, x);
287 x.add(additional_data.omega, d);
288 print_vectors(iter, x, r, d);
301template <
class VectorType>
306 const VectorType &)
const
311template <
class VectorType>
312inline typename VectorType::value_type
314 const VectorType &d)
const
316 if (!additional_data.use_preconditioned_residual)
323template <
class VectorType>
327 additional_data.omega = om;
@ iterate
Continue iteration.
@ success
Stop iteration, goal reached.
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
void set_omega(const double om=1.)
SolverRichardson(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
void Tsolve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
virtual ~SolverRichardson() override=default
virtual VectorType::value_type criterion(const VectorType &r, const VectorType &d) const
AdditionalData additional_data
SolverRichardson(SolverControl &cn, const AdditionalData &data=AdditionalData())
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define AssertThrow(cond, exc)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
AdditionalData(const double omega=1, const bool use_preconditioned_residual=false)
bool use_preconditioned_residual