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_coarse.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2002 - 2020 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_coarse_h
17#define dealii_mg_coarse_h
18
19
20#include <deal.II/base/config.h>
21
25
27
29
32
37template <class VectorType = Vector<double>>
39{
40public:
45
50
54 void
56
60 void
62
66 void
67 operator()(const unsigned int level,
68 VectorType & dst,
69 const VectorType & src) const override;
70
71private:
78};
79
80
81
88template <class VectorType,
89 class SolverType,
90 class MatrixType,
91 class PreconditionerType>
93{
94public:
99
105 const MatrixType & matrix,
106 const PreconditionerType &precondition);
107
112 void
114 const MatrixType & matrix,
115 const PreconditionerType &precondition);
116
120 void
122
127 virtual void
128 operator()(const unsigned int level,
129 VectorType & dst,
130 const VectorType & src) const override;
131
132private:
139 MatrixType,
140 PreconditionerType>>
142
146 SmartPointer<const MatrixType,
149 MatrixType,
150 PreconditionerType>>
152
156 SmartPointer<const PreconditionerType,
159 MatrixType,
160 PreconditionerType>>
162};
163
164
165
174template <typename number = double, class VectorType = Vector<number>>
176{
177public:
182
186 void
188
189 void
190 operator()(const unsigned int level,
191 VectorType & dst,
192 const VectorType & src) const override;
193
194private:
199};
200
207template <typename number = double, class VectorType = Vector<number>>
208class MGCoarseGridSVD : public MGCoarseGridBase<VectorType>
209{
210public:
214 MGCoarseGridSVD() = default;
215
219 void
220 initialize(const FullMatrix<number> &A, const double threshold = 0);
221
222 void
223 operator()(const unsigned int level,
224 VectorType & dst,
225 const VectorType & src) const;
226
230 void
231 log() const;
232
233private:
238};
239
242#ifndef DOXYGEN
243/* ------------------ Functions for MGCoarseGridApplySmoother -----------*/
244template <class VectorType>
246 : coarse_smooth(nullptr)
247{}
248
249template <class VectorType>
251 const MGSmootherBase<VectorType> &coarse_smooth)
252 : coarse_smooth(nullptr)
253{
254 initialize(coarse_smooth);
255}
256
257
258template <class VectorType>
259void
261 const MGSmootherBase<VectorType> &coarse_smooth_)
262{
263 coarse_smooth =
266 typeid(*this).name());
267}
268
269
270template <class VectorType>
271void
273{
274 coarse_smooth = nullptr;
275}
276
277
278template <class VectorType>
279void
281 VectorType & dst,
282 const VectorType & src) const
283{
284 coarse_smooth->smooth(level, dst, src);
285}
286
287/* ------------------ Functions for MGCoarseGridIterativeSolver ------------ */
288
289template <class VectorType,
290 class SolverType,
291 class MatrixType,
292 class PreconditionerType>
295 MatrixType,
296 PreconditionerType>::MGCoarseGridIterativeSolver()
297 : solver(0, typeid(*this).name())
298 , matrix(0, typeid(*this).name())
299 , preconditioner(0, typeid(*this).name())
300{}
301
302
303
304template <class VectorType,
305 class SolverType,
306 class MatrixType,
307 class PreconditionerType>
310 MatrixType,
311 PreconditionerType>::
312 MGCoarseGridIterativeSolver(SolverType & solver,
313 const MatrixType & matrix,
314 const PreconditionerType &preconditioner)
315 : solver(&solver, typeid(*this).name())
316 , matrix(&matrix, typeid(*this).name())
317 , preconditioner(&preconditioner, typeid(*this).name())
318{}
319
320
321
322template <class VectorType,
323 class SolverType,
324 class MatrixType,
325 class PreconditionerType>
326void
328 VectorType,
330 MatrixType,
331 PreconditionerType>::initialize(SolverType & solver_,
332 const MatrixType & matrix_,
333 const PreconditionerType &preconditioner_)
334{
335 solver = &solver_;
336 matrix = &matrix_;
337 preconditioner = &preconditioner_;
338}
339
340
341
342template <class VectorType,
343 class SolverType,
344 class MatrixType,
345 class PreconditionerType>
346void
349 MatrixType,
350 PreconditionerType>::clear()
351{
352 solver = 0;
353 matrix = 0;
354 preconditioner = 0;
355}
356
357
358
359template <class VectorType,
360 class SolverType,
361 class MatrixType,
362 class PreconditionerType>
363void
366 MatrixType,
367 PreconditionerType>::
368operator()(const unsigned int /*level*/,
369 VectorType & dst,
370 const VectorType &src) const
371{
372 Assert(solver != nullptr, ExcNotInitialized());
373 Assert(matrix != nullptr, ExcNotInitialized());
374 Assert(preconditioner != nullptr, ExcNotInitialized());
375 solver->solve(*matrix, dst, src, *preconditioner);
376}
377
378
379
380/* ------------------ Functions for MGCoarseGridHouseholder ------------ */
381
382template <typename number, class VectorType>
384 const FullMatrix<number> *A)
385{
386 if (A != nullptr)
387 householder.initialize(*A);
388}
389
390
391
392template <typename number, class VectorType>
393void
395 const FullMatrix<number> &A)
396{
397 householder.initialize(A);
398}
399
400
401
402template <typename number, class VectorType>
403void
405operator()(const unsigned int /*level*/,
406 VectorType & dst,
407 const VectorType &src) const
408{
409 householder.least_squares(dst, src);
410}
411
412//---------------------------------------------------------------------------
413
414
415
416template <typename number, class VectorType>
417void
419 double threshold)
420{
421 matrix.reinit(A.n_rows(), A.n_cols());
422 matrix = A;
423 matrix.compute_inverse_svd(threshold);
424}
425
426
427template <typename number, class VectorType>
428void
429MGCoarseGridSVD<number, VectorType>::operator()(const unsigned int /*level*/,
430 VectorType & dst,
431 const VectorType &src) const
432{
433 matrix.vmult(dst, src);
434}
435
436
437template <typename number, class VectorType>
438void
440{
441 const unsigned int n = std::min(matrix.n_rows(), matrix.n_cols());
442
443 for (unsigned int i = 0; i < n; ++i)
444 deallog << ' ' << matrix.singular_value(i);
445 deallog << std::endl;
446}
447
448
449#endif // DOXYGEN
450
452
453#endif
void initialize(const MGSmootherBase< VectorType > &coarse_smooth)
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const override
SmartPointer< const MGSmootherBase< VectorType >, MGCoarseGridApplySmoother< VectorType > > coarse_smooth
Definition: mg_coarse.h:77
MGCoarseGridApplySmoother(const MGSmootherBase< VectorType > &coarse_smooth)
void initialize(const FullMatrix< number > &A)
Householder< number > householder
Definition: mg_coarse.h:198
MGCoarseGridHouseholder(const FullMatrix< number > *A=nullptr)
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const override
SmartPointer< const MatrixType, MGCoarseGridIterativeSolver< VectorType, SolverType, MatrixType, PreconditionerType > > matrix
Definition: mg_coarse.h:151
SmartPointer< SolverType, MGCoarseGridIterativeSolver< VectorType, SolverType, MatrixType, PreconditionerType > > solver
Definition: mg_coarse.h:141
SmartPointer< const PreconditionerType, MGCoarseGridIterativeSolver< VectorType, SolverType, MatrixType, PreconditionerType > > preconditioner
Definition: mg_coarse.h:161
void initialize(SolverType &solver, const MatrixType &matrix, const PreconditionerType &precondition)
virtual void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const override
MGCoarseGridIterativeSolver(SolverType &solver, const MatrixType &matrix, const PreconditionerType &precondition)
void log() const
MGCoarseGridSVD()=default
LAPACKFullMatrix< number > matrix
Definition: mg_coarse.h:237
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const
void initialize(const FullMatrix< number > &A, const double threshold=0)
#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()
LogStream deallog
Definition: logstream.cc:37
@ matrix
Contents is actually a matrix.
static const char A
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)