16 #ifndef dealii_solver_richardson_h 17 #define dealii_solver_richardson_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/logstream.h> 22 #include <deal.II/lac/solver.h> 23 #include <deal.II/lac/solver_control.h> 24 #include <deal.II/base/signaling_nan.h> 26 DEAL_II_NAMESPACE_OPEN
60 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>
122 Tsolve (
const MatrixType &A,
125 const PreconditionerType &preconditioner);
140 const VectorType &d)
const;
149 virtual typename VectorType::value_type
criterion(
const VectorType &r,
150 const VectorType &d)
const;
163 template <
class VectorType>
167 const bool use_preconditioned_residual)
170 use_preconditioned_residual(use_preconditioned_residual)
174 template <
class VectorType>
177 const AdditionalData &data)
179 Solver<VectorType> (cn,mem),
180 additional_data(data)
185 template <
class VectorType>
187 const AdditionalData &data)
190 additional_data(data)
195 template <
class VectorType>
196 template <
typename MatrixType,
typename PreconditionerType>
201 const PreconditionerType &preconditioner)
205 double last_criterion = -std::numeric_limits<double>::max();
207 unsigned int iter = 0;
229 preconditioner.vmult(d,r);
233 last_criterion = criterion(r, d);
234 conv = this->iteration_status (iter, last_criterion, x);
238 x.add(additional_data.omega,d);
239 print_vectors(iter,x,r,d);
253 template <
class VectorType>
254 template <
typename MatrixType,
typename PreconditionerType>
259 const PreconditionerType &preconditioner)
262 double last_criterion = -std::numeric_limits<double>::max();
264 unsigned int iter = 0;
286 preconditioner.Tvmult(d,r);
288 last_criterion = criterion(r, d);
289 conv = this->iteration_status (iter, last_criterion, x);
293 x.add(additional_data.omega,d);
294 print_vectors(iter,x,r,d);
307 template <
class VectorType>
312 const VectorType &)
const 317 template <
class VectorType>
319 typename VectorType::value_type
321 const VectorType &d)
const 323 if (!additional_data.use_preconditioned_residual)
330 template <
class VectorType>
334 additional_data.omega=om;
339 DEAL_II_NAMESPACE_CLOSE
virtual ~SolverRichardson()=default
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.
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)