16#ifndef dealii_mg_block_smoother_h
17#define dealii_mg_block_smoother_h
48template <
typename MatrixType,
class RelaxationType,
typename number>
74 template <
class MGMatrixType,
class MGRelaxationType>
136template <
typename MatrixType,
class RelaxationType,
typename number>
138 const unsigned int steps,
140 const bool symmetric,
145 , mem(&this->vector_memory)
149template <
typename MatrixType,
class RelaxationType,
typename number>
153 unsigned int i = matrices.min_level(), max_level = matrices.max_level();
154 for (; i <= max_level; ++i)
162template <
typename MatrixType,
class RelaxationType,
typename number>
163template <
class MGMatrixType,
class MGRelaxationType>
166 const MGMatrixType & m,
167 const MGRelaxationType &s)
169 const unsigned int min = m.min_level();
170 const unsigned int max = m.max_level();
172 matrices.resize(min, max);
173 smoothers.resize(min, max);
175 for (
unsigned int i = min; i <=
max; ++i)
180 matrices[i] = linear_operator<BlockVector<number>>(
182 smoothers[i] = linear_operator<BlockVector<number>>(matrices[i], s[i]);
187template <
typename MatrixType,
class RelaxationType,
typename number>
196template <
typename MatrixType,
class RelaxationType,
typename number>
200 return sizeof(*this) + matrices.memory_consumption() +
201 smoothers.memory_consumption() +
202 this->vector_memory.memory_consumption();
206template <
typename MatrixType,
class RelaxationType,
typename number>
209 const unsigned int level,
215 unsigned int maxlevel = matrices.max_level();
216 unsigned int steps2 = this->steps;
219 steps2 *= (1 << (maxlevel -
level));
227 if (this->symmetric && (steps2 % 2 == 0))
230 for (
unsigned int i = 0; i < steps2; ++i)
234 matrices[
level].vmult(*r, u);
235 r->sadd(-1., 1., rhs);
236 smoothers[
level].Tvmult(*d, *r);
240 matrices[
level].vmult(*r, u);
241 r->sadd(-1., 1., rhs);
242 smoothers[
level].vmult(*d, *r);
virtual void smooth(const unsigned int level, BlockVector< number > &u, const BlockVector< number > &rhs) const
MGLevelObject< LinearOperator< BlockVector< number > > > smoothers
std::size_t memory_consumption() const
SmartPointer< VectorMemory< BlockVector< number > >, MGSmootherBlock< MatrixType, RelaxationType, number > > mem
void initialize(const MGMatrixType &matrices, const MGRelaxationType &smoothers)
void set_reverse(const bool)
MGLevelObject< LinearOperator< BlockVector< number > > > matrices
MGSmootherBlock(const unsigned int steps=1, const bool variable=false, const bool symmetric=false, const bool transpose=false, const bool reverse=false)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)