16 #ifndef dealii_mg_block_smoother_h 17 #define dealii_mg_block_smoother_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/mg_level_object.h> 22 #include <deal.II/base/smartpointer.h> 23 #include <deal.II/lac/block_vector.h> 24 #include <deal.II/lac/linear_operator.h> 25 #include <deal.II/lac/vector_memory.h> 26 #include <deal.II/multigrid/mg_smoother.h> 29 DEAL_II_NAMESPACE_OPEN
46 template <
typename MatrixType,
class RelaxationType,
typename number>
59 const unsigned int steps = 1,
87 template <
class MGMatrixType,
class MGRelaxationType>
106 virtual void smooth (
const unsigned int level,
142 template <
typename MatrixType,
class RelaxationType,
typename number>
146 const unsigned int steps,
148 const bool symmetric,
156 template <
typename MatrixType,
class RelaxationType,
typename number>
159 (
const unsigned int steps,
161 const bool symmetric,
166 mem(&this->vector_memory)
170 template <
typename MatrixType,
class RelaxationType,
typename number>
174 unsigned int i=matrices.min_level(),
175 max_level=matrices.max_level();
176 for (; i<=max_level; ++i)
184 template <
typename MatrixType,
class RelaxationType,
typename number>
185 template <
class MGMatrixType,
class MGRelaxationType>
188 const MGRelaxationType &s)
190 const unsigned int min = m.min_level();
191 const unsigned int max = m.max_level();
193 matrices.resize(
min,
max);
194 smoothers.resize(
min,
max);
196 for (
unsigned int i=
min; i<=
max; ++i)
201 matrices[i] = linear_operator<BlockVector<number> >(
203 smoothers[i] = linear_operator<BlockVector<number> >(matrices[i], s[i]);
208 template <
typename MatrixType,
class RelaxationType,
typename number>
217 template <
typename MatrixType,
class RelaxationType,
typename number>
223 + matrices.memory_consumption()
224 + smoothers.memory_consumption()
225 + this->vector_memory.memory_consumption();
229 template <
typename MatrixType,
class RelaxationType,
typename number>
237 unsigned int maxlevel = matrices.max_level();
238 unsigned int steps2 = this->steps;
241 steps2 *= (1<<(maxlevel-level));
249 if (this->symmetric && (steps2 % 2 == 0))
252 for (
unsigned int i=0; i<steps2; ++i)
256 matrices[level].vmult(*r,u);
258 smoothers[level].Tvmult(*d, *r);
262 matrices[level].vmult(*r,u);
264 smoothers[level].vmult(*d, *r);
274 DEAL_II_NAMESPACE_CLOSE
void initialize(const MGMatrixType &matrices, const MGRelaxationType &smoothers)
MGLevelObject< LinearOperator< BlockVector< number > > > smoothers
void set_reverse(const bool)
virtual void smooth(const unsigned int level, BlockVector< number > &u, const BlockVector< number > &rhs) const
VectorizedArray< Number > min(const ::VectorizedArray< Number > &x, const ::VectorizedArray< Number > &y)
SmartPointer< VectorMemory< BlockVector< number > >, MGSmootherBlock< MatrixType, RelaxationType, number > > mem
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
MGSmootherBlock(VectorMemory< BlockVector< number > > &mem, const unsigned int steps=1, const bool variable=false, const bool symmetric=false, const bool transpose=false, const bool reverse=false)
MGLevelObject< LinearOperator< BlockVector< number > > > matrices
std::size_t memory_consumption() const
VectorizedArray< Number > max(const ::VectorizedArray< Number > &x, const ::VectorizedArray< Number > &y)