1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2018 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_mg_base_h
17 #define dealii_mg_base_h
18 
19 /*
20  * This file contains some abstract base classes
21  * used by Multigrid.
22  */
23 
24 #include <deal.II/base/config.h>
25 
26 #include <deal.II/base/smartpointer.h>
27 #include <deal.II/base/subscriptor.h>
28 
29 #include <deal.II/lac/vector.h>
30 
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 /*!@addtogroup mg */
35 /*@{*/
36 
37 
38 /**
39  * Multilevel matrix base. This class sets up the interface needed by
40  * multilevel algorithms. It has no relation to the actual matrix type and
41  * takes the vector class as only template argument.
42  *
43  * Usually, the derived class mg::Matrix, which operates on an MGLevelObject
44  * of matrices, will be sufficient for applications.
45  */
46 template <typename VectorType>
47 class MGMatrixBase : public Subscriptor
48 {
49 public:
50   /*
51    * Virtual destructor.
52    */
53   virtual ~MGMatrixBase() override = default;
54 
55   /**
56    * Matrix-vector-multiplication on a certain level.
57    */
58   virtual void
59   vmult(const unsigned int level,
60         VectorType &       dst,
61         const VectorType & src) const = 0;
62 
63   /**
64    * Adding matrix-vector-multiplication on a certain level.
65    */
66   virtual void
67   vmult_add(const unsigned int level,
68             VectorType &       dst,
69             const VectorType & src) const = 0;
70 
71   /**
72    * Transpose matrix-vector-multiplication on a certain level.
73    */
74   virtual void
75   Tvmult(const unsigned int level,
76          VectorType &       dst,
77          const VectorType & src) const = 0;
78 
79   /**
80    * Adding transpose matrix-vector-multiplication on a certain level.
81    */
82   virtual void
83   Tvmult_add(const unsigned int level,
84              VectorType &       dst,
85              const VectorType & src) const = 0;
86 
87   /**
88    * Return the minimal level for which matrices are stored.
89    */
90   virtual unsigned int
91   get_minlevel() const = 0;
92 
93   /**
94    * Return the minimal level for which matrices are stored.
95    */
96   virtual unsigned int
97   get_maxlevel() const = 0;
98 };
99 
100 
101 /**
102  * Base class for coarse grid solvers.  This defines the virtual parenthesis
103  * operator, being the interface used by multigrid methods. Any implementation
104  * will be done by derived classes.
105  */
106 template <typename VectorType>
107 class MGCoarseGridBase : public Subscriptor
108 {
109 public:
110   /**
111    * Virtual destructor.
112    */
113   virtual ~MGCoarseGridBase() override = default;
114 
115   /**
116    * Solution operator.
117    */
118   virtual void
119   operator()(const unsigned int level,
120              VectorType &       dst,
121              const VectorType & src) const = 0;
122 };
123 
124 
125 /**
126  * Base class used to declare the operations needed by a concrete class
127  * implementing prolongation and restriction of vectors in the multigrid
128  * context. This class is abstract and has no implementation of these
129  * operations.
130  *
131  * There are several derived classes, reflecting the fact that vector types
132  * and numbering of the fine-grid discretization and of the multi-level
133  * implementation are independent.
134  *
135  * If you use multigrid for a single PDE or for your complete system of
136  * equations, you will use MGTransferPrebuilt together with Multigrid. The
137  * vector types used on the fine grid as well as for the multilevel operations
138  * may be Vector or BlockVector. In both cases, MGTransferPrebuilt will
139  * operate on all components of the solution.
140  *
141  * @note For the following, it is important to realize the difference between
142  * a solution
143  * @ref GlossComponent "component"
144  * and a solution
145  * @ref GlossBlock "block".
146  * The distinction only applies if vector valued elements are used, but is
147  * quite important then. This is reflected in the fact that it is not possible
148  * right now to use transfer classes based on MGTransferComponentBase for
149  * genuine vector valued elements, but descendants of MGTransferBlockBase
150  * would have to be applied. In the following text, we will use the term
151  * <em>block</em>, but remark that it might refer to components as well.
152  *
153  * @todo update the following documentation, since it does not reflect the
154  * latest changes in structure.
155  *
156  * For mixed systems, it may be required to do multigrid only for a single
157  * component or for some components. The classes MGTransferSelect and
158  * MGTransferBlock handle these cases.
159  *
160  * MGTransferSelect is used if you use multigrid (on Vector objects) for a
161  * single component, possibly grouped using <tt>mg_target_component</tt>.
162  *
163  * The class MGTransferBlock handles the case where your multigrid method
164  * operates on BlockVector objects. These can contain all or a consecutive set
165  * of the blocks of the complete system. Since most smoothers cannot operate
166  * on block structures, it is not clear whether this case is really useful.
167  * Therefore, a tested implementation of this case will be supplied when
168  * needed.
169  */
170 template <typename VectorType>
171 class MGTransferBase : public Subscriptor
172 {
173 public:
174   /**
175    * Destructor. Does nothing here, but needs to be declared virtual anyway.
176    */
177   virtual ~MGTransferBase() override = default;
178 
179   /**
180    * Prolongate a vector from level <tt>to_level-1</tt> to level
181    * <tt>to_level</tt>. The previous content of <tt>dst</tt> is overwritten.
182    *
183    * @arg src is a vector with as many elements as there are degrees of
184    * freedom on the coarser level involved.
185    *
186    * @arg dst has as many elements as there are degrees of freedom on the
187    * finer level.
188    */
189   virtual void
190   prolongate(const unsigned int to_level,
191              VectorType &       dst,
192              const VectorType & src) const = 0;
193 
194   /**
195    * Restrict a vector from level <tt>from_level</tt> to level
196    * <tt>from_level-1</tt> and add this restriction to <tt>dst</tt>. If the
197    * region covered by cells on level <tt>from_level</tt> is smaller than that
198    * of level <tt>from_level-1</tt> (local refinement), then some degrees of
199    * freedom in <tt>dst</tt> are active and will not be altered. For the other
200    * degrees of freedom, the result of the restriction is added.
201    *
202    * @arg src is a vector with as many elements as there are degrees of
203    * freedom on the finer level
204    *
205    * @arg dst has as many elements as there are degrees of freedom on the
206    * coarser level.
207    */
208   virtual void
209   restrict_and_add(const unsigned int from_level,
210                    VectorType &       dst,
211                    const VectorType & src) const = 0;
212 };
213 
214 
215 
216 /**
217  * Base class for multigrid smoothers. Does nothing but defining the interface
218  * used by multigrid methods.
219  *
220  * The smoother interface provides two methods, a smooth() method and an
221  * apply() method. The multigrid preconditioner interfaces distinguish between
222  * the two for efficiency reasons: Upon entry to the preconditioner operations,
223  * the vector @p u needs to be set to zero and smoothing starts by a simple
224  * application of the smoother on the @p rhs vector. This method is provided by
225  * the apply() method of this class. It is the same as first setting @p u to
226  * zero and then calling smooth(), but for many classes the separate apply()
227  * interface is more efficient because it can skip one matrix-vector product.
228  *
229  * In the multigrid preconditioner interfaces, the apply() method is used for
230  * the pre-smoothing operation because the previous content in the solution
231  * vector needs to be overwritten for a new incoming residual. On the other
232  * hand, all subsequent operations need to smooth the content already present
233  * in the vector @p u given the right hand side, which is done by smooth().
234  */
235 template <typename VectorType>
236 class MGSmootherBase : public Subscriptor
237 {
238 public:
239   /**
240    * Virtual destructor.
241    */
242   virtual ~MGSmootherBase() override = default;
243 
244   /**
245    * Release matrices.
246    */
247   virtual void
248   clear() = 0;
249 
250   /**
251    * Smoothing function that smooths the content in @p u given the right hand
252    * side vector @p rhs. This is the function used in multigrid methods.
253    */
254   virtual void
255   smooth(const unsigned int level,
256          VectorType &       u,
257          const VectorType & rhs) const = 0;
258 
259   /**
260    * As opposed to the smooth() function, this function applies the action of
261    * the smoothing, overwriting the previous content in the vector u. This
262    * function must be equivalent to the following code
263    * @code
264    * u = 0;
265    * smooth(level, u, rhs);
266    * @endcode
267    * but can usually be implemented more efficiently than the former. If a
268    * particular smoother does not override the apply() method, the default
269    * implementation as described here is used.
270    *
271    * In the multigrid preconditioner interfaces, the apply() method is used for
272    * the pre-smoothing operation because the previous content in the solution
273    * vector needs to be overwritten for a new incoming residual. On the other
274    * hand, all subsequent operations need to smooth the content already present
275    * in the vector @p u given the right hand side, which is done by smooth().
276    */
277   virtual void
278   apply(const unsigned int level, VectorType &u, const VectorType &rhs) const;
279 };
280 
281 /*@}*/
282 
283 DEAL_II_NAMESPACE_CLOSE
284 
285 #endif
286