16 #ifndef dealii_solver_richardson_h 17 #define dealii_solver_richardson_h 20 #include <deal.II/base/config.h> 22 #include <deal.II/base/logstream.h> 23 #include <deal.II/base/signaling_nan.h> 25 #include <deal.II/lac/solver.h> 26 #include <deal.II/lac/solver_control.h> 28 DEAL_II_NAMESPACE_OPEN
62 template <
class VectorType = Vector<
double>>
110 template <
typename MatrixType,
typename PreconditionerType>
112 solve(
const MatrixType & A,
114 const VectorType & b,
115 const PreconditionerType &preconditioner);
120 template <
typename MatrixType,
typename PreconditionerType>
122 Tsolve(
const MatrixType & A,
124 const VectorType & b,
125 const PreconditionerType &preconditioner);
140 const VectorType & x,
141 const VectorType & r,
142 const VectorType & d)
const;
151 virtual typename VectorType::value_type
152 criterion(
const VectorType &r,
const VectorType &d)
const;
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>
198 const VectorType & b,
199 const PreconditionerType &preconditioner)
203 double last_criterion = -std::numeric_limits<double>::max();
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>
255 const VectorType & b,
256 const PreconditionerType &preconditioner)
259 double last_criterion = -std::numeric_limits<double>::max();
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>
309 const VectorType &)
const 314 template <
class VectorType>
315 inline typename VectorType::value_type
317 const VectorType &d)
const 319 if (!additional_data.use_preconditioned_residual)
326 template <
class VectorType>
330 additional_data.omega = om;
335 DEAL_II_NAMESPACE_CLOSE
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
bool use_preconditioned_residual
void set_omega(const double om=1.)
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
AdditionalData additional_data
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())
Stop iteration, goal reached.
virtual ~SolverRichardson() override=default
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
virtual VectorType::value_type criterion(const VectorType &r, const VectorType &d) const
AdditionalData(const double omega=1, const bool use_preconditioned_residual=false)