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/smartpointer.h> 22 #include <deal.II/lac/pointer_matrix.h> 23 #include <deal.II/lac/vector_memory.h> 24 #include <deal.II/lac/block_vector.h> 25 #include <deal.II/multigrid/mg_smoother.h> 26 #include <deal.II/base/mg_level_object.h> 29 DEAL_II_NAMESPACE_OPEN
46 template <
typename MatrixType,
class RelaxationType,
typename number>
58 const unsigned int steps = 1,
62 const bool reverse =
false) DEAL_II_DEPRECATED;
86 template <
class MGMatrixType,
class MGRelaxationType>
105 virtual void smooth (
const unsigned int level,
141 template <
typename MatrixType,
class RelaxationType,
typename number>
145 const unsigned int steps,
147 const bool symmetric,
155 template <
typename MatrixType,
class RelaxationType,
typename number>
158 (
const unsigned int steps,
160 const bool symmetric,
165 mem(&this->vector_memory)
169 template <
typename MatrixType,
class RelaxationType,
typename number>
173 unsigned int i=matrices.min_level(),
174 max_level=matrices.max_level();
175 for (; i<=max_level; ++i)
183 template <
typename MatrixType,
class RelaxationType,
typename number>
184 template <
class MGMatrixType,
class MGRelaxationType>
187 const MGRelaxationType &s)
189 const unsigned int min = m.min_level();
190 const unsigned int max = m.max_level();
192 matrices.resize(
min,
max);
193 smoothers.resize(
min,
max);
195 for (
unsigned int i=
min; i<=
max; ++i)
198 smoothers[i] = &s[i];
203 template <
typename MatrixType,
class RelaxationType,
typename number>
212 template <
typename MatrixType,
class RelaxationType,
typename number>
218 + matrices.memory_consumption()
219 + smoothers.memory_consumption()
220 + this->vector_memory.memory_consumption();
224 template <
typename MatrixType,
class RelaxationType,
typename number>
230 deallog.
push(
"Smooth");
232 unsigned int maxlevel = matrices.max_level();
233 unsigned int steps2 = this->steps;
236 steps2 *= (1<<(maxlevel-level));
244 if (this->symmetric && (steps2 % 2 == 0))
247 for (
unsigned int i=0; i<steps2; ++i)
251 matrices[level].vmult(*r,u);
253 smoothers[level].Tvmult(*d, *r);
257 matrices[level].vmult(*r,u);
259 smoothers[level].vmult(*d, *r);
271 DEAL_II_NAMESPACE_CLOSE
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) 1
MGLevelObject< PointerMatrix< MatrixType, BlockVector< number > > > matrices
void initialize(const MGMatrixType &matrices, const MGRelaxationType &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)
void push(const std::string &text)
std::size_t memory_consumption() const
MGLevelObject< PointerMatrix< RelaxationType, BlockVector< number > > > smoothers
VectorizedArray< Number > max(const ::VectorizedArray< Number > &x, const ::VectorizedArray< Number > &y)