Reference documentation for deal.II version 9.0.0
mg_coarse.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2002 - 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_mg_coarse_h
17 #define dealii_mg_coarse_h
18 
19 
20 #include <deal.II/lac/full_matrix.h>
21 #include <deal.II/lac/householder.h>
22 #include <deal.II/lac/linear_operator.h>
23 #include <deal.II/lac/matrix_lib.h>
24 #include <deal.II/multigrid/mg_base.h>
25 
26 DEAL_II_NAMESPACE_OPEN
27 
30 
37 template <class VectorType = Vector<double> >
38 class MGCoarseGridApplySmoother : public MGCoarseGridBase<VectorType>
39 {
40 public:
45 
50 
54  void clear ();
55 
60 
64  void operator() (const unsigned int level,
65  VectorType &dst,
66  const VectorType &src) const;
67 
68 private:
73 };
74 
75 
89 template <typename SolverType, class VectorType = Vector<double> >
90 class DEAL_II_DEPRECATED MGCoarseGridLACIteration : public MGCoarseGridBase<VectorType>
91 {
92 public:
97 
102  template <typename MatrixType, typename PreconditionerType>
103  MGCoarseGridLACIteration (SolverType &,
104  const MatrixType &,
105  const PreconditionerType &);
106 
111 
115  template <typename MatrixType, typename PreconditionerType>
116  void initialize (SolverType &,
117  const MatrixType &,
118  const PreconditionerType &);
119 
123  void clear ();
124 
129  void operator() (const unsigned int level,
130  VectorType &dst,
131  const VectorType &src) const;
132 
137  template <typename MatrixType>
138  void set_matrix (const MatrixType &);
139 
140 private:
145 
150 
155 
156 };
157 
158 
159 
166 template <class VectorType,
167  class SolverType,
168  class MatrixType,
169  class PreconditionerType>
171 {
172 public:
177 
182  MGCoarseGridIterativeSolver (SolverType &solver,
183  const MatrixType &matrix,
184  const PreconditionerType &precondition);
185 
190  void initialize (SolverType &solver,
191  const MatrixType &matrix,
192  const PreconditionerType &precondition);
193 
197  void clear ();
198 
203  virtual void operator() (const unsigned int level,
204  VectorType &dst,
205  const VectorType &src) const;
206 
207 private:
211  SmartPointer<SolverType,MGCoarseGridIterativeSolver<VectorType,
212  SolverType,
213  MatrixType,
214  PreconditionerType> > solver;
215 
219  SmartPointer<const MatrixType,MGCoarseGridIterativeSolver<VectorType,
220  SolverType,
221  MatrixType,
222  PreconditionerType> > matrix;
223 
227  SmartPointer<const PreconditionerType,MGCoarseGridIterativeSolver<VectorType,
228  SolverType,
229  MatrixType,
230  PreconditionerType> > preconditioner;
231 
232 };
233 
234 
235 
236 
247 template <typename number = double, class VectorType = Vector<number> >
248 class MGCoarseGridHouseholder : public MGCoarseGridBase<VectorType>
249 {
250 public:
254  MGCoarseGridHouseholder (const FullMatrix<number> *A = nullptr);
255 
259  void initialize (const FullMatrix<number> &A);
260 
261  void operator() (const unsigned int level,
262  VectorType &dst,
263  const VectorType &src) const;
264 
265 private:
270 };
271 
280 template <typename number = double, class VectorType = Vector<number> >
281 class MGCoarseGridSVD : public MGCoarseGridBase<VectorType>
282 {
283 public:
287  MGCoarseGridSVD () = default;
288 
292  void initialize (const FullMatrix<number> &A, const double threshold = 0);
293 
294  void operator() (const unsigned int level,
295  VectorType &dst,
296  const VectorType &src) const;
297 
301  void log () const;
302 
303 private:
304 
309 };
310 
313 #ifndef DOXYGEN
314 /* ------------------ Functions for MGCoarseGridApplySmoother -----------*/
315 template <class VectorType>
317  : coarse_smooth(nullptr)
318 {
319 }
320 
321 template <class VectorType>
323  : coarse_smooth(nullptr)
324 {
325  initialize(coarse_smooth);
326 }
327 
328 
329 template <class VectorType>
330 void
332 {
333  coarse_smooth =
335  (&coarse_smooth_,typeid(*this).name());
336 }
337 
338 
339 template <class VectorType>
340 void
342 {
343  coarse_smooth = nullptr;
344 }
345 
346 
347 template <class VectorType>
348 void
350  VectorType &dst,
351  const VectorType &src) const
352 {
353  coarse_smooth->smooth(level, dst, src);
354 }
355 
356 /* ------------------ Functions for MGCoarseGridLACIteration ------------ */
357 
358 
359 template <typename SolverType, class VectorType>
362  :
363  solver(0, typeid(*this).name()),
364  matrix(0),
365  precondition(0)
366 {}
367 
368 
369 template <typename SolverType, class VectorType>
370 template <typename MatrixType, typename PreconditionerType>
372 ::MGCoarseGridLACIteration (SolverType &s,
373  const MatrixType &m,
374  const PreconditionerType &p)
375  :
376  solver(&s, typeid(*this).name())
377 {
378  // Workaround: Unfortunately, not every "m" object has a rich enough
379  // interface to populate reinit_(domain|range)_vector. Thus, supply an
380  // empty LinearOperator exemplar.
381  matrix = linear_operator<VectorType>(LinearOperator<VectorType>(), m);
382  precondition = linear_operator<VectorType>(matrix, p);
383 }
384 
385 
386 template <typename SolverType, class VectorType>
389 {
390  clear();
391 }
392 
393 
394 template <typename SolverType, class VectorType>
395 template <typename MatrixType, typename PreconditionerType>
396 void
398 ::initialize (SolverType &s,
399  const MatrixType &m,
400  const PreconditionerType &p)
401 {
402  solver = &s;
403  // Workaround: Unfortunately, not every "m" object has a rich enough
404  // interface to populate reinit_(domain|range)_vector. Thus, supply an
405  // empty LinearOperator exemplar.
406  matrix = linear_operator<VectorType>(LinearOperator<VectorType>(), m);
407  precondition = linear_operator<VectorType>(matrix, p);
408 }
409 
410 
411 template <typename SolverType, class VectorType>
412 void
414 ::clear()
415 {
416  solver = nullptr;
418  precondition = LinearOperator<VectorType>();
419 
420 }
421 
422 
423 template <typename SolverType, class VectorType>
424 void
426 ::operator() (const unsigned int /* level */,
427  VectorType &dst,
428  const VectorType &src) const
429 {
430  Assert(solver!=nullptr, ExcNotInitialized());
431  solver->solve(matrix, dst, src, precondition);
432 }
433 
434 
435 template <typename SolverType, class VectorType>
436 template <typename MatrixType>
437 void
439 ::set_matrix(const MatrixType &m)
440 {
441  // Workaround: Unfortunately, not every "m" object has a rich enough
442  // interface to populate reinit_(domain|range)_vector. Thus, supply an
443  // empty LinearOperator exemplar.
444  matrix = linear_operator<VectorType>(LinearOperator<VectorType>(), m);
445 }
446 
447 
448 
449 /* ------------------ Functions for MGCoarseGridIterativeSolver ------------ */
450 
451 template <class VectorType,
452  class SolverType,
453  class MatrixType,
454  class PreconditionerType>
457  :
458  solver(0, typeid(*this).name()),
459  matrix(0, typeid(*this).name()),
460  preconditioner(0, typeid(*this).name())
461 {}
462 
463 
464 
465 template <class VectorType,
466  class SolverType,
467  class MatrixType,
468  class PreconditionerType>
470 ::MGCoarseGridIterativeSolver (SolverType &solver,
471  const MatrixType &matrix,
472  const PreconditionerType &preconditioner)
473  :
474  solver (&solver, typeid(*this).name()),
475  matrix (&matrix, typeid(*this).name()),
476  preconditioner (&preconditioner, typeid(*this).name())
477 {}
478 
479 
480 
481 template <class VectorType,
482  class SolverType,
483  class MatrixType,
484  class PreconditionerType>
485 void
487 ::initialize (SolverType &solver_,
488  const MatrixType &matrix_,
489  const PreconditionerType &preconditioner_)
490 {
491  solver = &solver_;
492  matrix = &matrix_;
493  preconditioner = &preconditioner_;
494 
495 }
496 
497 
498 
499 template <class VectorType,
500  class SolverType,
501  class MatrixType,
502  class PreconditionerType>
503 void
505 ::clear ()
506 {
507  solver = 0;
508  matrix = 0;
509  preconditioner = 0;
510 }
511 
512 
513 
514 template <class VectorType,
515  class SolverType,
516  class MatrixType,
517  class PreconditionerType>
518 void
520 ::operator() (const unsigned int /*level*/,
521  VectorType &dst,
522  const VectorType &src) const
523 {
524  Assert(solver!=nullptr, ExcNotInitialized());
525  Assert(matrix!=nullptr, ExcNotInitialized());
526  Assert(preconditioner!=nullptr, ExcNotInitialized());
527  solver->solve(*matrix, dst, src, *preconditioner);
528 }
529 
530 
531 
532 /* ------------------ Functions for MGCoarseGridHouseholder ------------ */
533 
534 template <typename number, class VectorType>
536 (const FullMatrix<number> *A)
537 {
538  if (A != nullptr) householder.initialize(*A);
539 }
540 
541 
542 
543 template <typename number, class VectorType>
544 void
546 {
547  householder.initialize(A);
548 }
549 
550 
551 
552 template <typename number, class VectorType>
553 void
555  VectorType &dst,
556  const VectorType &src) const
557 {
558  householder.least_squares(dst, src);
559 }
560 
561 //---------------------------------------------------------------------------
562 
563 
564 
565 template <typename number, class VectorType>
566 void
568  double threshold)
569 {
570  matrix.reinit(A.n_rows(), A.n_cols());
571  matrix = A;
572  matrix.compute_inverse_svd(threshold);
573 }
574 
575 
576 template <typename number, class VectorType>
577 void
579  const unsigned int /*level*/,
580  VectorType &dst,
581  const VectorType &src) const
582 {
583  matrix.vmult(dst, src);
584 }
585 
586 
587 template <typename number, class VectorType>
588 void
590 {
591  const unsigned int n = std::min(matrix.n_rows(), matrix.n_cols());
592 
593  for (unsigned int i=0; i<n; ++i)
594  deallog << ' ' << matrix.singular_value(i);
595  deallog << std::endl;
596 }
597 
598 
599 #endif // DOXYGEN
600 
601 DEAL_II_NAMESPACE_CLOSE
602 
603 #endif
virtual void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const =0
void initialize(const FullMatrix< number > &A, const double threshold=0)
Contents is actually a matrix.
LinearOperator< VectorType > precondition
Definition: mg_coarse.h:154
static ::ExceptionBase & ExcNotInitialized()
SmartPointer< const PreconditionerType, MGCoarseGridIterativeSolver< VectorType, SolverType, MatrixType, PreconditionerType > > preconditioner
Definition: mg_coarse.h:230
virtual void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const
void initialize(SolverType &, const MatrixType &, const PreconditionerType &)
LAPACKFullMatrix< number > matrix
Definition: mg_coarse.h:308
SmartPointer< SolverType, MGCoarseGridIterativeSolver< VectorType, SolverType, MatrixType, PreconditionerType > > solver
Definition: mg_coarse.h:214
SmartPointer< SolverType, MGCoarseGridLACIteration< SolverType, VectorType > > solver
Definition: mg_coarse.h:144
#define Assert(cond, exc)
Definition: exceptions.h:1142
MGCoarseGridSVD()=default
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const
void initialize(SolverType &solver, const MatrixType &matrix, const PreconditionerType &precondition)
LinearOperator< VectorType > matrix
Definition: mg_coarse.h:149
void initialize(const MGSmootherBase< VectorType > &coarse_smooth)
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const
MGCoarseGridHouseholder(const FullMatrix< number > *A=nullptr)
SmartPointer< const MatrixType, MGCoarseGridIterativeSolver< VectorType, SolverType, MatrixType, PreconditionerType > > matrix
Definition: mg_coarse.h:222
Householder< number > householder
Definition: mg_coarse.h:269
void operator()(const unsigned int level, VectorType &dst, const VectorType &src) const
SmartPointer< const MGSmootherBase< VectorType >, MGCoarseGridApplySmoother< VectorType > > coarse_smooth
Definition: mg_coarse.h:72
void initialize(const FullMatrix< number > &A)
void set_matrix(const MatrixType &)
void log() const