Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
mg_matrix.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2003 - 2021 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_matrix_h
17#define dealii_mg_matrix_h
18
19#include <deal.II/base/config.h>
20
22
25#include <deal.II/lac/vector.h>
26
28
29#include <memory>
30
32
35
36namespace mg
37{
43 template <typename VectorType = Vector<double>>
44 class Matrix : public MGMatrixBase<VectorType>
45 {
46 public:
50 Matrix() = default;
51
56 template <typename MatrixType>
58
63 template <typename MatrixType>
64 void
66
70 void
71 reset();
72
76 const LinearOperator<VectorType> &operator[](unsigned int level) const;
77
78 virtual void
79 vmult(const unsigned int level,
80 VectorType & dst,
81 const VectorType & src) const override;
82 virtual void
83 vmult_add(const unsigned int level,
84 VectorType & dst,
85 const VectorType & src) const override;
86 virtual void
87 Tvmult(const unsigned int level,
88 VectorType & dst,
89 const VectorType & src) const override;
90 virtual void
91 Tvmult_add(const unsigned int level,
92 VectorType & dst,
93 const VectorType & src) const override;
94 virtual unsigned int
95 get_minlevel() const override;
96 virtual unsigned int
97 get_maxlevel() const override;
98
102 std::size_t
103 memory_consumption() const;
104
105 private:
107 };
108
109} // namespace mg
110
111
120template <typename MatrixType, typename number>
121class MGMatrixSelect : public MGMatrixBase<Vector<number>>
122{
123public:
128 MGMatrixSelect(const unsigned int row = 0,
129 const unsigned int col = 0,
131
136 void
138
142 void
143 select_block(const unsigned int row, const unsigned int col);
144
148 virtual void
149 vmult(const unsigned int level,
150 Vector<number> & dst,
151 const Vector<number> &src) const;
152
156 virtual void
157 vmult_add(const unsigned int level,
158 Vector<number> & dst,
159 const Vector<number> &src) const;
160
164 virtual void
165 Tvmult(const unsigned int level,
166 Vector<number> & dst,
167 const Vector<number> &src) const;
168
172 virtual void
173 Tvmult_add(const unsigned int level,
174 Vector<number> & dst,
175 const Vector<number> &src) const;
176
177private:
186 unsigned int row;
190 unsigned int col;
191};
192
195/*----------------------------------------------------------------------*/
196
197namespace mg
198{
199 template <typename VectorType>
200 template <typename MatrixType>
201 inline void
203 {
204 matrices.resize(p.min_level(), p.max_level());
205 for (unsigned int level = p.min_level(); level <= p.max_level(); ++level)
206 {
207 // Workaround: Unfortunately, not every "p[level]" object has a
208 // rich enough interface to populate reinit_(domain|range)_vector.
209 // Thus, apply an empty LinearOperator exemplar.
210 matrices[level] =
211 linear_operator<VectorType>(LinearOperator<VectorType>(),
213 p[level]));
214 }
215 }
216
217
218
219 template <typename VectorType>
220 inline void
222 {
223 matrices.resize(0, 0);
224 }
225
226
227
228 template <typename VectorType>
229 template <typename MatrixType>
231 {
232 initialize(p);
233 }
234
235
236
237 template <typename VectorType>
239 operator[](unsigned int level) const
240 {
241 return matrices[level];
242 }
243
244
245
246 template <typename VectorType>
247 void
249 VectorType & dst,
250 const VectorType & src) const
251 {
252 matrices[level].vmult(dst, src);
253 }
254
255
256
257 template <typename VectorType>
258 void
260 VectorType & dst,
261 const VectorType & src) const
262 {
263 matrices[level].vmult_add(dst, src);
264 }
265
266
267
268 template <typename VectorType>
269 void
271 VectorType & dst,
272 const VectorType & src) const
273 {
274 matrices[level].Tvmult(dst, src);
275 }
276
277
278
279 template <typename VectorType>
280 void
282 VectorType & dst,
283 const VectorType & src) const
284 {
285 matrices[level].Tvmult_add(dst, src);
286 }
287
288
289
290 template <typename VectorType>
291 unsigned int
293 {
294 return matrices.min_level();
295 }
296
297
298
299 template <typename VectorType>
300 unsigned int
302 {
303 return matrices.max_level();
304 }
305
306
307
308 template <typename VectorType>
309 inline std::size_t
311 {
312 return sizeof(*this) + matrices->memory_consumption();
313 }
314} // namespace mg
315
316
317/*----------------------------------------------------------------------*/
318
319template <typename MatrixType, typename number>
321 const unsigned int col,
323 : matrix(p, typeid(*this).name())
324 , row(row)
325 , col(col)
326{}
327
328
329
330template <typename MatrixType, typename number>
331void
333{
334 matrix = p;
335}
336
337
338
339template <typename MatrixType, typename number>
340void
342 const unsigned int bcol)
343{
344 row = brow;
345 col = bcol;
346}
347
348
349
350template <typename MatrixType, typename number>
351void
353 Vector<number> & dst,
354 const Vector<number> &src) const
355{
357
359 m[level].block(row, col).vmult(dst, src);
360}
361
362
363
364template <typename MatrixType, typename number>
365void
367 Vector<number> & dst,
368 const Vector<number> &src) const
369{
371
373 m[level].block(row, col).vmult_add(dst, src);
374}
375
376
377
378template <typename MatrixType, typename number>
379void
381 Vector<number> & dst,
382 const Vector<number> &src) const
383{
385
387 m[level].block(row, col).Tvmult(dst, src);
388}
389
390
391
392template <typename MatrixType, typename number>
393void
395 Vector<number> & dst,
396 const Vector<number> &src) const
397{
399
401 m[level].block(row, col).Tvmult_add(dst, src);
402}
403
405
406#endif
unsigned int max_level() const
unsigned int min_level() const
void set_matrix(MGLevelObject< MatrixType > *M)
Definition: mg_matrix.h:332
virtual void vmult(const unsigned int level, Vector< number > &dst, const Vector< number > &src) const
Definition: mg_matrix.h:352
virtual void vmult_add(const unsigned int level, Vector< number > &dst, const Vector< number > &src) const
Definition: mg_matrix.h:366
MGMatrixSelect(const unsigned int row=0, const unsigned int col=0, MGLevelObject< MatrixType > *matrix=0)
Definition: mg_matrix.h:320
void select_block(const unsigned int row, const unsigned int col)
Definition: mg_matrix.h:341
virtual void Tvmult(const unsigned int level, Vector< number > &dst, const Vector< number > &src) const
Definition: mg_matrix.h:380
unsigned int row
Definition: mg_matrix.h:186
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSelect< MatrixType, number > > matrix
Definition: mg_matrix.h:182
unsigned int col
Definition: mg_matrix.h:190
virtual void Tvmult_add(const unsigned int level, Vector< number > &dst, const Vector< number > &src) const
Definition: mg_matrix.h:394
virtual unsigned int get_minlevel() const override
Definition: mg_matrix.h:292
virtual void vmult(const unsigned int level, VectorType &dst, const VectorType &src) const override
Definition: mg_matrix.h:248
virtual void vmult_add(const unsigned int level, VectorType &dst, const VectorType &src) const override
Definition: mg_matrix.h:259
std::size_t memory_consumption() const
Definition: mg_matrix.h:310
virtual void Tvmult(const unsigned int level, VectorType &dst, const VectorType &src) const override
Definition: mg_matrix.h:270
void initialize(const MGLevelObject< MatrixType > &M)
Definition: mg_matrix.h:202
void reset()
Definition: mg_matrix.h:221
virtual unsigned int get_maxlevel() const override
Definition: mg_matrix.h:301
MGLevelObject< LinearOperator< VectorType > > matrices
Definition: mg_matrix.h:106
virtual void Tvmult_add(const unsigned int level, VectorType &dst, const VectorType &src) const override
Definition: mg_matrix.h:281
Matrix()=default
const LinearOperator< VectorType > & operator[](unsigned int level) const
Definition: mg_matrix.h:239
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
unsigned int level
Definition: grid_out.cc:4590
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcNotInitialized()
@ matrix
Contents is actually a matrix.
T & get_underlying_value(T &p)
Definition: utilities.h:1422
Definition: mg.h:82