MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
multigrid.cpp
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 #include "multigrid.hpp"
13 
14 namespace mfem
15 {
16 
18  : fespaces(fespaces_), cycleType(CycleType::VCYCLE), preSmoothingSteps(1),
19  postSmoothingSteps(1)
20 {}
21 
23 {
24  for (int i = 0; i < operators.Size(); ++i)
25  {
26  if (ownedOperators[i])
27  {
28  delete operators[i];
29  }
30  if (ownedSmoothers[i])
31  {
32  delete smoothers[i];
33  }
34  delete X[i];
35  delete Y[i];
36  delete R[i];
37  delete Z[i];
38  }
39 
40  operators.DeleteAll();
41  smoothers.DeleteAll();
42  X.DeleteAll();
43  Y.DeleteAll();
44  R.DeleteAll();
45  Z.DeleteAll();
46 
47  for (int i = 0; i < bfs.Size(); ++i)
48  {
49  delete bfs[i];
50  }
51 
52  bfs.DeleteAll();
53 
54  for (int i = 0; i < essentialTrueDofs.Size(); ++i)
55  {
56  delete essentialTrueDofs[i];
57  }
58 
59  essentialTrueDofs.DeleteAll();
60 }
61 
62 void Multigrid::AddLevel(Operator* opr, Solver* smoother, bool ownOperator,
63  bool ownSmoother)
64 {
65  operators.Append(opr);
66  smoothers.Append(smoother);
67  ownedOperators.Append(ownOperator);
68  ownedSmoothers.Append(ownSmoother);
69  width = opr->Width();
70  height = opr->Height();
71 
72  X.Append(new Vector(height));
73  *X.Last() = 0.0;
74  Y.Append(new Vector(height));
75  *Y.Last() = 0.0;
76  R.Append(new Vector(height));
77  *R.Last() = 0.0;
78  Z.Append(new Vector(height));
79  *Z.Last() = 0.0;
80 }
81 
82 int Multigrid::NumLevels() const { return operators.Size(); }
83 
84 int Multigrid::GetFinestLevelIndex() const { return NumLevels() - 1; }
85 
86 const Operator* Multigrid::GetOperatorAtLevel(int level) const
87 {
88  return operators[level];
89 }
90 
92 {
93  return operators[level];
94 }
95 
97 {
98  return GetOperatorAtLevel(operators.Size() - 1);
99 }
100 
102 {
103  return GetOperatorAtLevel(operators.Size() - 1);
104 }
105 
107 {
108  return smoothers[level];
109 }
110 
112 {
113  return smoothers[level];
114 }
115 
116 void Multigrid::SetCycleType(CycleType cycleType_, int preSmoothingSteps_,
117  int postSmoothingSteps_)
118 {
119  cycleType = cycleType_;
120  preSmoothingSteps = preSmoothingSteps_;
121  postSmoothingSteps = postSmoothingSteps_;
122 }
123 
124 void Multigrid::Mult(const Vector& x, Vector& y) const
125 {
126  MFEM_ASSERT(NumLevels() > 0, "");
127  *X.Last() = x;
128  *Y.Last() = 0.0;
129  Cycle(GetFinestLevelIndex());
130  y = *Y.Last();
131 }
132 
134 {
135  MFEM_ABORT("SetOperator not supported in Multigrid");
136 }
137 
138 void Multigrid::SmoothingStep(int level, bool transpose) const
139 {
140  GetOperatorAtLevel(level)->Mult(*Y[level], *R[level]); // r = A x
141  subtract(*X[level], *R[level], *R[level]); // r = b - A x
142  if (transpose)
143  {
144  GetSmootherAtLevel(level)->MultTranspose(*R[level], *Z[level]); // z = S r
145  }
146  else
147  {
148  GetSmootherAtLevel(level)->Mult(*R[level], *Z[level]); // z = S r
149  }
150  add(*Y[level], 1.0, *Z[level], *Y[level]); // x = x + S (b - A x)
151 }
152 
153 void Multigrid::Cycle(int level) const
154 {
155  if (level == 0)
156  {
157  GetSmootherAtLevel(level)->Mult(*X[level], *Y[level]);
158  return;
159  }
160 
161  for (int i = 0; i < preSmoothingSteps; i++)
162  {
163  SmoothingStep(level, false);
164  }
165 
166  // Compute residual
167  GetOperatorAtLevel(level)->Mult(*Y[level], *R[level]);
168  subtract(*X[level], *R[level], *R[level]);
169 
170  // Restrict residual
171  fespaces.GetProlongationAtLevel(level - 1)->MultTranspose(*R[level],
172  *X[level - 1]);
173 
174  // Init zeros
175  *Y[level - 1] = 0.0;
176 
177  // Corrections
178  int corrections = 1;
179  if (cycleType == CycleType::WCYCLE)
180  {
181  corrections = 2;
182  }
183  for (int correction = 0; correction < corrections; ++correction)
184  {
185  Cycle(level - 1);
186  }
187 
188  // Prolongate
189  fespaces.GetProlongationAtLevel(level - 1)->Mult(*Y[level - 1], *R[level]);
190 
191  // Add update
192  *Y[level] += *R[level];
193 
194  // Post-smooth
195  for (int i = 0; i < postSmoothingSteps; i++)
196  {
197  SmoothingStep(level, true);
198  }
199 }
200 
202  Vector& X, Vector& B)
203 {
204  bfs.Last()->FormLinearSystem(*essentialTrueDofs.Last(), x, b, A, X, B);
205 }
206 
208  Vector& x)
209 {
210  bfs.Last()->RecoverFEMSolution(X, b, x);
211 }
212 
213 } // namespace mfem
Array< BilinearForm * > bfs
Definition: multigrid.hpp:37
Multigrid(const FiniteElementSpaceHierarchy &fespaces_)
Constructs an empty multigrid for the given FiniteElementSpaceHierarchy.
Definition: multigrid.cpp:17
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:71
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:261
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: operator.hpp:92
const FiniteElementSpaceHierarchy & fespaces
Definition: multigrid.hpp:35
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:65
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:726
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
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition: vector.cpp:408
virtual ~Multigrid()
Destructor.
Definition: multigrid.cpp:22
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
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
Operator * GetProlongationAtLevel(int level) const
Returns the prolongation operator from the finite element space at level to the finite element space ...
Abstract operator.
Definition: operator.hpp:24
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
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