Reference documentation for deal.II version 8.4.2
mg_block_smoother.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2005 - 2015 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__mg_block_smoother_h
17 #define dealii__mg_block_smoother_h
18 
19 
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_base.h>
26 #include <deal.II/base/mg_level_object.h>
27 #include <vector>
28 
29 DEAL_II_NAMESPACE_OPEN
30 
31 /*
32  * MGSmootherBase is defined in mg_base.h
33  */
34 
37 
46 template <typename MatrixType, class RelaxationType, typename number>
48  : public MGSmootherBase<BlockVector<number> >
49 {
50 public:
55  const unsigned int steps = 1,
56  const bool variable = false,
57  const bool symmetric = false,
58  const bool transpose = false,
59  const bool reverse = false);
60 
74  template <class MGMatrixType, class MGRelaxationType>
75  void initialize (const MGMatrixType &matrices,
76  const MGRelaxationType &smoothers);
77 
81  void clear ();
82 
86  void set_steps (const unsigned int);
87 
91  void set_variable (const bool);
92 
96  void set_symmetric (const bool);
97 
101  void set_transpose (const bool);
102 
106  void set_reverse (const bool);
107 
113  virtual void smooth (const unsigned int level,
115  const BlockVector<number> &rhs) const;
116 private:
121 
126 
130  unsigned int steps;
131 
135  bool variable;
136 
140  bool symmetric;
141 
142  /*
143  * Transposed?
144  */
145  bool transpose;
146 
150  bool reverse;
151 
156 
157 };
158 
161 //---------------------------------------------------------------------------
162 
163 #ifndef DOXYGEN
164 
165 template <typename MatrixType, class RelaxationType, typename number>
166 inline
169  const unsigned int steps,
170  const bool variable,
171  const bool symmetric,
172  const bool transpose,
173  const bool reverse)
174  :
175  steps(steps),
176  variable(variable),
177  symmetric(symmetric),
178  transpose(transpose),
179  reverse(reverse),
180  mem(mem)
181 {}
182 
183 
184 template <typename MatrixType, class RelaxationType, typename number>
185 inline void
187 {
188  unsigned int i=matrices.min_level(),
189  max_level=matrices.max_level();
190  for (; i<=max_level; ++i)
191  {
192  smoothers[i] = 0;
193  matrices[i] = 0;
194  }
195 }
196 
197 
198 template <typename MatrixType, class RelaxationType, typename number>
199 template <class MGMatrixType, class MGRelaxationType>
200 inline void
202  const MGRelaxationType &s)
203 {
204  const unsigned int min = m.min_level();
205  const unsigned int max = m.max_level();
206 
207  matrices.resize(min, max);
208  smoothers.resize(min, max);
209 
210  for (unsigned int i=min; i<=max; ++i)
211  {
212  matrices[i] = &m[i];
213  smoothers[i] = &s[i];
214  }
215 }
216 
217 template <typename MatrixType, class RelaxationType, typename number>
218 inline void
220 set_steps (const unsigned int s)
221 {
222  steps = s;
223 }
224 
225 
226 template <typename MatrixType, class RelaxationType, typename number>
227 inline void
229 set_variable (const bool flag)
230 {
231  variable = flag;
232 }
233 
234 
235 template <typename MatrixType, class RelaxationType, typename number>
236 inline void
238 set_symmetric (const bool flag)
239 {
240  symmetric = flag;
241 }
242 
243 
244 template <typename MatrixType, class RelaxationType, typename number>
245 inline void
247 set_transpose (const bool flag)
248 {
249  transpose = flag;
250 }
251 
252 
253 template <typename MatrixType, class RelaxationType, typename number>
254 inline void
256 set_reverse (const bool flag)
257 {
258  reverse = flag;
259 }
260 
261 
262 template <typename MatrixType, class RelaxationType, typename number>
263 inline void
266  const BlockVector<number> &rhs) const
267 {
268  deallog.push("Smooth");
269 
270  unsigned int maxlevel = matrices.max_level();
271  unsigned int steps2 = steps;
272 
273  if (variable)
274  steps2 *= (1<<(maxlevel-level));
275 
278  r->reinit(u);
279  d->reinit(u);
280 
281  bool T = transpose;
282  if (symmetric && (steps2 % 2 == 0))
283  T = false;
284 
285  for (unsigned int i=0; i<steps2; ++i)
286  {
287  if (T)
288  {
289  matrices[level].vmult(*r,u);
290  r->sadd(-1.,1.,rhs);
291  smoothers[level].Tvmult(*d, *r);
292  }
293  else
294  {
295  matrices[level].vmult(*r,u);
296  r->sadd(-1.,1.,rhs);
297  smoothers[level].vmult(*d, *r);
298  }
299  u += *d;
300  if (symmetric)
301  T = !T;
302  }
303 
304  mem.free(r);
305  mem.free(d);
306  deallog.pop();
307 }
308 
309 #endif // DOXYGEN
310 
311 DEAL_II_NAMESPACE_CLOSE
312 
313 #endif
unsigned int max_level() const
void pop()
Definition: logstream.cc:300
virtual VectorType * alloc()=0
MGLevelObject< PointerMatrix< MatrixType, BlockVector< number > > > matrices
void set_symmetric(const bool)
void initialize(const MGMatrixType &matrices, const MGRelaxationType &smoothers)
void set_reverse(const bool)
void set_variable(const bool)
void set_transpose(const bool)
unsigned int min_level() const
virtual void free(const VectorType *const)=0
virtual void smooth(const unsigned int level, BlockVector< number > &u, const BlockVector< number > &rhs) const
void set_steps(const unsigned int)
void reinit(const unsigned int n_blocks, const size_type block_size=0, const bool omit_zeroing_entries=false)
void sadd(const value_type s, const BlockVectorBase &V)
VectorMemory< BlockVector< number > > & mem
void push(const std::string &text)
Definition: logstream.cc:288
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)
void resize(const unsigned int new_minlevel, const unsigned int new_maxlevel)
MGLevelObject< PointerMatrix< RelaxationType, BlockVector< number > > > smoothers
unsigned int steps