Reference documentation for deal.II version 9.3.3
\(\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\}}\)
solver_relaxation.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2010 - 2020 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_solver_relaxation_h
17#define dealii_solver_relaxation_h
18
19
20#include <deal.II/base/config.h>
21
24
25#include <deal.II/lac/solver.h>
27
29
56template <typename VectorType = Vector<double>>
57class SolverRelaxation : public SolverBase<VectorType>
58{
59public:
65 {};
66
71 const AdditionalData &data = AdditionalData());
72
78 template <typename MatrixType, class RelaxationType>
79 void
80 solve(const MatrixType & A,
81 VectorType & x,
82 const VectorType & b,
83 const RelaxationType &R);
84};
85
86//----------------------------------------------------------------------//
87
88template <class VectorType>
90 const AdditionalData &)
91 : SolverBase<VectorType>(cn)
92{}
93
94
95
96template <class VectorType>
97template <typename MatrixType, class RelaxationType>
98void
100 VectorType & x,
101 const VectorType & b,
102 const RelaxationType &R)
103{
106
107 // Memory allocation
108 typename VectorMemory<VectorType>::Pointer Vr(mem);
109 VectorType & r = *Vr;
110 r.reinit(x);
111 typename VectorMemory<VectorType>::Pointer Vd(mem);
112 VectorType & d = *Vd;
113 d.reinit(x);
114
115 LogStream::Prefix prefix("Relaxation");
116
117 int iter = 0;
118 // Main loop
119 for (; conv == SolverControl::iterate; iter++)
120 {
121 // Compute residual
122 A.vmult(r, x);
123 r.sadd(-1., 1., b);
124
125 // The required norm of the
126 // (preconditioned)
127 // residual is computed in
128 // criterion() and stored
129 // in res.
130 conv = this->iteration_status(iter, r.l2_norm(), x);
131 if (conv != SolverControl::iterate)
132 break;
133 R.step(x, b);
134 }
135
136 // in case of failure: throw exception
138 SolverControl::NoConvergence(iter, r.l2_norm()));
139 // otherwise exit as normal
140}
141
142
144
145#endif
@ iterate
Continue iteration.
@ success
Stop iteration, goal reached.
SolverRelaxation(SolverControl &cn, const AdditionalData &data=AdditionalData())
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const RelaxationType &R)
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define AssertThrow(cond, exc)
Definition: exceptions.h:1575
static const char A
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)