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/subscriptor.h> 25 #include <deal.II/base/signaling_nan.h> 27 DEAL_II_NAMESPACE_OPEN
61 template <
class VectorType = Vector<
double> >
111 template<
typename MatrixType,
typename PreconditionerType>
113 solve (
const MatrixType &A,
116 const PreconditionerType &precondition);
121 template<
typename MatrixType,
typename PreconditionerType>
123 Tsolve (
const MatrixType &A,
126 const PreconditionerType &precondition);
141 const VectorType &d)
const;
147 virtual typename VectorType::value_type
criterion();
172 typename VectorType::value_type
res;
180 template <
class VectorType>
184 const bool use_preconditioned_residual)
187 use_preconditioned_residual(use_preconditioned_residual)
191 template <
class VectorType>
194 const AdditionalData &data)
196 Solver<VectorType> (cn,mem),
199 additional_data(data)
204 template <
class VectorType>
206 const AdditionalData &data)
211 additional_data(data),
217 template <
class VectorType>
222 template <
class VectorType>
223 template <
typename MatrixType,
typename PreconditionerType>
228 const PreconditionerType &precondition)
232 double last_criterion = -std::numeric_limits<double>::max();
234 unsigned int iter = 0;
237 Vr = this->memory.alloc();
240 Vd = this->memory.alloc();
244 deallog.
push(
"Richardson");
255 precondition.vmult(d,r);
262 last_criterion = criterion();
263 conv = this->iteration_status (iter, last_criterion, x);
267 x.add(additional_data.omega,d);
268 print_vectors(iter,x,r,d);
275 this->memory.free(Vr);
276 this->memory.free(Vd);
281 this->memory.free(Vr);
282 this->memory.free(Vd);
293 template <
class VectorType>
294 template <
typename MatrixType,
typename PreconditionerType>
299 const PreconditionerType &precondition)
302 double last_criterion = -std::numeric_limits<double>::max();
304 unsigned int iter = 0;
307 Vr = this->memory.alloc();
310 Vd =this-> memory.alloc();
314 deallog.
push(
"RichardsonT");
325 precondition.Tvmult(d,r);
327 last_criterion = criterion();
328 conv = this->iteration_status (iter, last_criterion, x);
332 x.add(additional_data.omega,d);
333 print_vectors(iter,x,r,d);
340 this->memory.free(Vr);
341 this->memory.free(Vd);
347 this->memory.free(Vr);
348 this->memory.free(Vd);
358 template <
class VectorType>
363 const VectorType &)
const 368 template <
class VectorType>
369 inline typename VectorType::value_type
372 if (!additional_data.use_preconditioned_residual)
380 template <
class VectorType>
384 additional_data.omega=om;
389 DEAL_II_NAMESPACE_CLOSE
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
#define AssertThrow(cond, exc)
SolverRichardson(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
Stop iteration, goal reached.
VectorType::value_type res
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
void push(const std::string &text)
void Tsolve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
virtual ~SolverRichardson()
AdditionalData(const double omega=1, const bool use_preconditioned_residual=false)
virtual VectorType::value_type criterion()