MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
multigrid.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_MULTIGRID
13 #define MFEM_MULTIGRID
14 
15 #include "fespacehierarchy.hpp"
16 #include "bilinearform.hpp"
17 
18 #include "../linalg/operator.hpp"
19 #include "../linalg/handle.hpp"
20 
21 namespace mfem
22 {
23 
24 /// Multigrid solver class
25 class Multigrid : public Solver
26 {
27 public:
28  enum class CycleType
29  {
30  VCYCLE,
31  WCYCLE
32  };
33 
34 protected:
38 
39 private:
40  Array<Operator*> operators;
41  Array<Solver*> smoothers;
42 
43  Array<bool> ownedOperators;
44  Array<bool> ownedSmoothers;
45 
46  CycleType cycleType;
47  int preSmoothingSteps;
48  int postSmoothingSteps;
49 
50  mutable Array<Vector*> X;
51  mutable Array<Vector*> Y;
52  mutable Array<Vector*> R;
53  mutable Array<Vector*> Z;
54 
55 public:
56  /// Constructs an empty multigrid for the given FiniteElementSpaceHierarchy
57  Multigrid(const FiniteElementSpaceHierarchy& fespaces_);
58 
59  /// Destructor
60  virtual ~Multigrid();
61 
62  /// Adds a level to the multigrid operator hierarchy.
63  /** The ownership of the operators and solvers/smoothers may be transferred
64  to the Multigrid by setting the according boolean variables. */
65  void AddLevel(Operator* opr, Solver* smoother, bool ownOperator,
66  bool ownSmoother);
67 
68  /// Returns the number of levels
69  int NumLevels() const;
70 
71  /// Returns the index of the finest level
72  int GetFinestLevelIndex() const;
73 
74  /// Returns operator at given level
75  const Operator* GetOperatorAtLevel(int level) const;
76 
77  /// Returns operator at given level
78  Operator* GetOperatorAtLevel(int level);
79 
80  /// Returns operator at finest level
81  const Operator* GetOperatorAtFinestLevel() const;
82 
83  /// Returns operator at finest level
85 
86  /// Returns smoother at given level
87  Solver* GetSmootherAtLevel(int level) const;
88 
89  /// Returns smoother at given level
90  Solver* GetSmootherAtLevel(int level);
91 
92  /// Set the cycle type and number of pre- and post-smoothing steps used by Mult
93  void SetCycleType(CycleType cycleType_, int preSmoothingSteps_,
94  int postSmoothingSteps_);
95 
96  /// Application of the multigrid as a preconditioner
97  virtual void Mult(const Vector& x, Vector& y) const override;
98 
99  /// Not supported for multigrid
100  virtual void SetOperator(const Operator& op) override;
101 
102  /// Form the linear system A X = B, corresponding to the operator on the finest level
104  Vector& B);
105 
106  /// Recover the solution of a linear system formed with FormFineLinearSystem()
107  void RecoverFineFEMSolution(const Vector& X, const Vector& b, Vector& x);
108 
109 private:
110  /// Application of a smoothing step at particular level
111  void SmoothingStep(int level, bool transpose) const;
112 
113  /// Application of a cycle at particular level
114  void Cycle(int level) const;
115 };
116 
117 } // namespace mfem
118 
119 #endif
Multigrid solver class.
Definition: multigrid.hpp:25
Array< BilinearForm * > bfs
Definition: multigrid.hpp:37
Multigrid(const FiniteElementSpaceHierarchy &fespaces_)
Constructs an empty multigrid for the given FiniteElementSpaceHierarchy.
Definition: multigrid.cpp:17
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
const FiniteElementSpaceHierarchy & fespaces
Definition: multigrid.hpp:35
void FormFineLinearSystem(Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)
Form the linear system A X = B, corresponding to the operator on the finest level.
Definition: multigrid.cpp:201
double b
Definition: lissajous.cpp:42
virtual void SetOperator(const Operator &op) override
Not supported for multigrid.
Definition: multigrid.cpp:133
virtual ~Multigrid()
Destructor.
Definition: multigrid.cpp:22
void AddLevel(Operator *opr, Solver *smoother, bool ownOperator, bool ownSmoother)
Adds a level to the multigrid operator hierarchy.
Definition: multigrid.cpp:62
Solver * GetSmootherAtLevel(int level) const
Returns smoother at given level.
Definition: multigrid.cpp:106
const Operator * GetOperatorAtFinestLevel() const
Returns operator at finest level.
Definition: multigrid.cpp:96
virtual void Mult(const Vector &x, Vector &y) const override
Application of the multigrid as a preconditioner.
Definition: multigrid.cpp:124
int GetFinestLevelIndex() const
Returns the index of the finest level.
Definition: multigrid.cpp:84
Array< Array< int > * > essentialTrueDofs
Definition: multigrid.hpp:36
int NumLevels() const
Returns the number of levels.
Definition: multigrid.cpp:82
const Operator * GetOperatorAtLevel(int level) const
Returns operator at given level.
Definition: multigrid.cpp:86
Vector data type.
Definition: vector.hpp:51
void RecoverFineFEMSolution(const Vector &X, const Vector &b, Vector &x)
Recover the solution of a linear system formed with FormFineLinearSystem()
Definition: multigrid.cpp:207
Base class for solvers.
Definition: operator.hpp:634
Abstract operator.
Definition: operator.hpp:24
void SetCycleType(CycleType cycleType_, int preSmoothingSteps_, int postSmoothingSteps_)
Set the cycle type and number of pre- and post-smoothing steps used by Mult.
Definition: multigrid.cpp:116