Reference documentation for deal.II version 8.4.2
multigrid.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2016 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__multigrid_h
17 #define dealii__multigrid_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/subscriptor.h>
22 #include <deal.II/base/smartpointer.h>
23 #include <deal.II/dofs/dof_handler.h>
24 #include <deal.II/lac/sparse_matrix.h>
25 #include <deal.II/lac/vector.h>
26 #include <deal.II/multigrid/mg_base.h>
27 #include <deal.II/base/mg_level_object.h>
28 
29 #include <vector>
30 
31 DEAL_II_NAMESPACE_OPEN
32 
35 
63 template <typename VectorType>
64 class Multigrid : public Subscriptor
65 {
66 public:
70  enum Cycle
71  {
78  };
79 
80  typedef VectorType vector_type;
81  typedef const VectorType const_vector_type;
82 
92  template <int dim>
93  Multigrid(const DoFHandler<dim> &mg_dof_handler,
99  Cycle cycle = v_cycle);
100 
106  Multigrid(const unsigned int minlevel,
107  const unsigned int maxlevel,
113  Cycle cycle = v_cycle);
114 
118  void reinit (const unsigned int minlevel,
119  const unsigned int maxlevel);
120 
125  void cycle ();
126 
137  void vcycle ();
138 
153 
168 
172  unsigned int get_maxlevel() const;
173 
177  unsigned int get_minlevel() const;
178 
185  void set_maxlevel (const unsigned int);
186 
202  void set_minlevel (const unsigned int level,
203  bool relative = false);
204 
208  void set_cycle(Cycle);
209 
214  void set_debug (const unsigned int);
215 
216 private:
217 
224  void level_v_step (const unsigned int level);
225 
233  void level_step (const unsigned int level, Cycle cycle);
234 
239 
243  unsigned int minlevel;
244 
248  unsigned int maxlevel;
249 
250 public:
256 
261 
262 private:
267 
272 
273 
278 
283 
288 
293 
298 
305 
313 
320 
327 
331  unsigned int debug;
332 
333  template<int dim, class OtherVectorType, class TRANSFER> friend class PreconditionMG;
334 };
335 
336 
349 template<int dim, typename VectorType, class TRANSFER>
351 {
352 public:
357  PreconditionMG(const DoFHandler<dim> &dof_handler,
359  const TRANSFER &transfer);
360 
364  bool empty () const;
365 
372  template<class OtherVectorType>
373  void vmult (OtherVectorType &dst,
374  const OtherVectorType &src) const;
375 
380  template<class OtherVectorType>
381  void vmult_add (OtherVectorType &dst,
382  const OtherVectorType &src) const;
383 
389  template<class OtherVectorType>
390  void Tvmult (OtherVectorType &dst,
391  const OtherVectorType &src) const;
392 
398  template<class OtherVectorType>
399  void Tvmult_add (OtherVectorType &dst,
400  const OtherVectorType &src) const;
401 
402 private:
407 
412 
417 };
418 
421 #ifndef DOXYGEN
422 /* --------------------------- inline functions --------------------- */
423 
424 
425 template <typename VectorType>
426 template <int dim>
433  Cycle cycle)
434  :
435  cycle_type(cycle),
436  minlevel(0),
437  maxlevel(mg_dof_handler.get_triangulation().n_global_levels()-1),
442  matrix(&matrix, typeid(*this).name()),
443  coarse(&coarse, typeid(*this).name()),
444  transfer(&transfer, typeid(*this).name()),
445  pre_smooth(&pre_smooth, typeid(*this).name()),
446  post_smooth(&post_smooth, typeid(*this).name()),
447  edge_down(0, typeid(*this).name()),
448  edge_up(0, typeid(*this).name()),
449  debug(0)
450 {}
451 
452 
453 
454 template <typename VectorType>
455 inline
456 unsigned int
458 {
459  return maxlevel;
460 }
461 
462 
463 
464 template <typename VectorType>
465 inline
466 unsigned int
468 {
469  return minlevel;
470 }
471 
472 
473 /* --------------------------- inline functions --------------------- */
474 
475 
476 template<int dim, typename VectorType, class TRANSFER>
480  const TRANSFER &transfer)
481  :
482  dof_handler(&dof_handler),
483  multigrid(&mg),
485 {}
486 
487 template<int dim, typename VectorType, class TRANSFER>
488 inline bool
490 {
491  return false;
492 }
493 
494 template<int dim, typename VectorType, class TRANSFER>
495 template<class OtherVectorType>
496 void
498 (OtherVectorType &dst,
499  const OtherVectorType &src) const
500 {
501  transfer->copy_to_mg(*dof_handler,
502  multigrid->defect,
503  src);
504  multigrid->cycle();
505 
506  transfer->copy_from_mg(*dof_handler,
507  dst,
508  multigrid->solution);
509 }
510 
511 
512 template<int dim, typename VectorType, class TRANSFER>
513 template<class OtherVectorType>
514 void
516 (OtherVectorType &dst,
517  const OtherVectorType &src) const
518 {
519  transfer->copy_to_mg(*dof_handler,
520  multigrid->defect,
521  src);
522  multigrid->cycle();
523  transfer->copy_from_mg_add(*dof_handler,
524  dst,
525  multigrid->solution);
526 }
527 
528 
529 template<int dim, typename VectorType, class TRANSFER>
530 template<class OtherVectorType>
531 void
533 (OtherVectorType &,
534  const OtherVectorType &) const
535 {
536  Assert(false, ExcNotImplemented());
537 }
538 
539 
540 template<int dim, typename VectorType, class TRANSFER>
541 template<class OtherVectorType>
542 void
544 (OtherVectorType &,
545  const OtherVectorType &) const
546 {
547  Assert(false, ExcNotImplemented());
548 }
549 
550 #endif // DOXYGEN
551 
552 DEAL_II_NAMESPACE_CLOSE
553 
554 #endif
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > matrix
Definition: multigrid.h:277
void vcycle()
void cycle()
SmartPointer< const TRANSFER, PreconditionMG< dim, VectorType, TRANSFER > > transfer
Definition: multigrid.h:416
SmartPointer< const MGMatrixBase< VectorType > > edge_in
Definition: multigrid.h:312
void vmult(OtherVectorType &dst, const OtherVectorType &src) const
void set_edge_matrices(const MGMatrixBase< VectorType > &edge_out, const MGMatrixBase< VectorType > &edge_in)
Definition: mg.h:81
unsigned int get_maxlevel() const
void set_maxlevel(const unsigned int)
MGLevelObject< VectorType > t
Definition: multigrid.h:266
bool empty() const
SmartPointer< const MGCoarseGridBase< VectorType >, Multigrid< VectorType > > coarse
Definition: multigrid.h:282
unsigned int debug
Definition: multigrid.h:331
Multigrid(const DoFHandler< dim > &mg_dof_handler, const MGMatrixBase< VectorType > &matrix, const MGCoarseGridBase< VectorType > &coarse, const MGTransferBase< VectorType > &transfer, const MGSmootherBase< VectorType > &pre_smooth, const MGSmootherBase< VectorType > &post_smooth, Cycle cycle=v_cycle)
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > edge_up
Definition: multigrid.h:326
void set_debug(const unsigned int)
void level_v_step(const unsigned int level)
SmartPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > post_smooth
Definition: multigrid.h:297
unsigned int maxlevel
Definition: multigrid.h:248
#define Assert(cond, exc)
Definition: exceptions.h:294
The F-cycle.
Definition: multigrid.h:77
void Tvmult(OtherVectorType &dst, const OtherVectorType &src) const
SmartPointer< Multigrid< VectorType >, PreconditionMG< dim, VectorType, TRANSFER > > multigrid
Definition: multigrid.h:411
unsigned int get_minlevel() const
SmartPointer< const MGMatrixBase< VectorType > > edge_out
Definition: multigrid.h:304
void vmult_add(OtherVectorType &dst, const OtherVectorType &src) const
The W-cycle.
Definition: multigrid.h:75
SmartPointer< const MGTransferBase< VectorType >, Multigrid< VectorType > > transfer
Definition: multigrid.h:287
void set_minlevel(const unsigned int level, bool relative=false)
The V-cycle.
Definition: multigrid.h:73
void set_edge_flux_matrices(const MGMatrixBase< VectorType > &edge_down, const MGMatrixBase< VectorType > &edge_up)
SmartPointer< const DoFHandler< dim >, PreconditionMG< dim, VectorType, TRANSFER > > dof_handler
Definition: multigrid.h:406
void Tvmult_add(OtherVectorType &dst, const OtherVectorType &src) const
unsigned int minlevel
Definition: multigrid.h:243
MGLevelObject< VectorType > solution
Definition: multigrid.h:260
void set_cycle(Cycle)
void level_step(const unsigned int level, Cycle cycle)
SmartPointer< const MGSmootherBase< VectorType >, Multigrid< VectorType > > pre_smooth
Definition: multigrid.h:292
SmartPointer< const MGMatrixBase< VectorType >, Multigrid< VectorType > > edge_down
Definition: multigrid.h:319
MGLevelObject< VectorType > defect2
Definition: multigrid.h:271
MGLevelObject< VectorType > defect
Definition: multigrid.h:255
PreconditionMG(const DoFHandler< dim > &dof_handler, Multigrid< VectorType > &mg, const TRANSFER &transfer)
Cycle cycle_type
Definition: multigrid.h:238
void reinit(const unsigned int minlevel, const unsigned int maxlevel)