MFEM  v4.6.0
Finite element discretization library
multigrid.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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  : cycleType(CycleType::VCYCLE), preSmoothingSteps(1), postSmoothingSteps(1),
19  nrhs(0)
20 {}
21 
23  const Array<Solver*>& smoothers_,
24  const Array<bool>& ownedOperators_,
25  const Array<bool>& ownedSmoothers_)
26  : Solver(operators_.Last()->Height(), operators_.Last()->Width()),
27  cycleType(CycleType::VCYCLE), preSmoothingSteps(1), postSmoothingSteps(1),
28  nrhs(0)
29 {
30  operators_.Copy(operators);
31  smoothers_.Copy(smoothers);
32  ownedOperators_.Copy(ownedOperators);
33  ownedSmoothers_.Copy(ownedSmoothers);
34 }
35 
37 {
38  for (int i = 0; i < operators.Size(); ++i)
39  {
40  if (ownedOperators[i])
41  {
42  delete operators[i];
43  }
44  if (ownedSmoothers[i])
45  {
46  delete smoothers[i];
47  }
48  }
49  EraseVectors();
50 }
51 
52 void MultigridBase::InitVectors() const
53 {
54  if (X.NumRows() > 0 && X.NumCols() > 0) { EraseVectors(); }
55  const int M = NumLevels();
56  X.SetSize(M, nrhs);
57  Y.SetSize(M, nrhs);
58  R.SetSize(M, nrhs);
59  Z.SetSize(M, nrhs);
60  for (int i = 0; i < X.NumRows(); ++i)
61  {
62  const int n = operators[i]->Height();
63  for (int j = 0; j < X.NumCols(); ++j)
64  {
65  X(i, j) = new Vector(n);
66  Y(i, j) = new Vector(n);
67  R(i, j) = new Vector(n);
68  Z(i, j) = new Vector(n);
69  }
70  }
71 }
72 
73 void MultigridBase::EraseVectors() const
74 {
75  for (int i = 0; i < X.NumRows(); ++i)
76  {
77  for (int j = 0; j < X.NumCols(); ++j)
78  {
79  delete X(i, j);
80  delete Y(i, j);
81  delete R(i, j);
82  delete Z(i, j);
83  }
84  }
85 }
86 
88  bool ownOperator, bool ownSmoother)
89 {
90  height = op->Height();
91  width = op->Width();
92  operators.Append(op);
93  smoothers.Append(smoother);
94  ownedOperators.Append(ownOperator);
95  ownedSmoothers.Append(ownSmoother);
96 }
97 
98 void MultigridBase::SetCycleType(CycleType cycleType_, int preSmoothingSteps_,
99  int postSmoothingSteps_)
100 {
101  cycleType = cycleType_;
102  preSmoothingSteps = preSmoothingSteps_;
103  postSmoothingSteps = postSmoothingSteps_;
104 }
105 
106 void MultigridBase::Mult(const Vector& x, Vector& y) const
107 {
108  Array<const Vector*> X_(1);
109  Array<Vector*> Y_(1);
110  X_[0] = &x;
111  Y_[0] = &y;
112  ArrayMult(X_, Y_);
113 }
114 
116  Array<Vector*>& Y_) const
117 {
118  MFEM_ASSERT(operators.Size() > 0,
119  "Multigrid solver does not have operators set!");
120  MFEM_ASSERT(X_.Size() == Y_.Size(),
121  "Number of columns mismatch in MultigridBase::Mult!");
122  if (iterative_mode)
123  {
124  MFEM_WARNING("Multigrid solver does not use iterative_mode and ignores "
125  "the initial guess!");
126  }
127 
128  // Add capacity as necessary
129  nrhs = X_.Size();
130  if (X.NumCols() < nrhs) { InitVectors(); }
131 
132  // Perform a single cycle
133  const int M = NumLevels();
134  for (int j = 0; j < nrhs; ++j)
135  {
136  MFEM_ASSERT(X_[j] && Y_[j], "Missing Vector in MultigridBase::Mult!");
137  *X(M - 1, j) = *X_[j];
138  *Y(M - 1, j) = 0.0;
139  }
140  Cycle(M - 1);
141  for (int j = 0; j < nrhs; ++j)
142  {
143  *Y_[j] = *Y(M - 1, j);
144  }
145 }
146 
147 void MultigridBase::SmoothingStep(int level, bool zero, bool transpose) const
148 {
149  // y = y + S (x - A y) or y = y + S^T (x - A y)
150  if (zero)
151  {
152  Array<Vector *> X_(X[level], nrhs), Y_(Y[level], nrhs);
153  GetSmootherAtLevel(level)->ArrayMult(X_, Y_);
154  }
155  else
156  {
157  Array<Vector *> Y_(Y[level], nrhs), R_(R[level], nrhs),
158  Z_(Z[level], nrhs);
159  for (int j = 0; j < nrhs; ++j)
160  {
161  *R_[j] = *X(level, j);
162  }
163  GetOperatorAtLevel(level)->ArrayAddMult(Y_, R_, -1.0);
164  if (transpose)
165  {
166  GetSmootherAtLevel(level)->ArrayMultTranspose(R_, Z_);
167  }
168  else
169  {
170  GetSmootherAtLevel(level)->ArrayMult(R_, Z_);
171  }
172  for (int j = 0; j < nrhs; ++j)
173  {
174  *Y_[j] += *Z_[j];
175  }
176  }
177 }
178 
179 void MultigridBase::Cycle(int level) const
180 {
181  // Coarse solve
182  if (level == 0)
183  {
184  SmoothingStep(0, true, false);
185  return;
186  }
187 
188  // Pre-smooth
189  for (int i = 0; i < preSmoothingSteps; ++i)
190  {
191  SmoothingStep(level, (cycleType == CycleType::VCYCLE && i == 0), false);
192  }
193 
194  // Compute residual and restrict
195  {
196  Array<Vector *> Y_(Y[level], nrhs), R_(R[level], nrhs),
197  X_(X[level - 1], nrhs);
198  for (int j = 0; j < nrhs; ++j)
199  {
200  *R_[j] = *X(level, j);
201  }
202  GetOperatorAtLevel(level)->ArrayAddMult(Y_, R_, -1.0);
203  GetProlongationAtLevel(level - 1)->ArrayMultTranspose(R_, X_);
204  for (int j = 0; j < nrhs; ++j)
205  {
206  *Y(level - 1, j) = 0.0;
207  }
208  }
209 
210  // Corrections
211  Cycle(level - 1);
213  {
214  Cycle(level - 1);
215  }
216 
217  // Prolongate and add
218  {
219  Array<Vector *> Y_(Y[level - 1], nrhs), Z_(Z[level], nrhs);
220  GetProlongationAtLevel(level - 1)->ArrayMult(Y_, Z_);
221  for (int j = 0; j < nrhs; ++j)
222  {
223  *Y(level, j) += *Z_[j];
224  }
225  }
226 
227  // Post-smooth
228  for (int i = 0; i < postSmoothingSteps; ++i)
229  {
230  SmoothingStep(level, false, true);
231  }
232 }
233 
235  : MultigridBase()
236 {}
237 
239  const Array<Solver*>& smoothers_,
240  const Array<Operator*>& prolongations_,
241  const Array<bool>& ownedOperators_,
242  const Array<bool>& ownedSmoothers_,
243  const Array<bool>& ownedProlongations_)
244  : MultigridBase(operators_, smoothers_, ownedOperators_, ownedSmoothers_)
245 {
246  prolongations_.Copy(prolongations);
247  ownedProlongations_.Copy(ownedProlongations);
248 }
249 
251 {
252  for (int i = 0; i < prolongations.Size(); ++i)
253  {
254  if (ownedProlongations[i])
255  {
256  delete prolongations[i];
257  }
258  }
259 }
260 
263  : MultigridBase(), fespaces(fespaces_)
264 {}
265 
267 {
268  for (int i = 0; i < bfs.Size(); ++i)
269  {
270  delete bfs[i];
271  }
272  for (int i = 0; i < essentialTrueDofs.Size(); ++i)
273  {
274  delete essentialTrueDofs[i];
275  }
276 }
277 
279  OperatorHandle& A,
280  Vector& X, Vector& B)
281 {
282  bfs.Last()->FormLinearSystem(*essentialTrueDofs.Last(), x, b, A, X, B);
283 }
284 
286  const Vector& b, Vector& x)
287 {
288  bfs.Last()->RecoverFEMSolution(X, b, x);
289 }
290 
291 } // namespace mfem
Array< bool > ownedOperators
Definition: multigrid.hpp:37
Array2D< Vector * > Z
Definition: multigrid.hpp:44
virtual ~GeometricMultigrid()
Destructor.
Definition: multigrid.cpp:266
Array2D< Vector * > X
Definition: multigrid.hpp:44
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
Array< BilinearForm * > bfs
Definition: multigrid.hpp:170
Array2D< Vector * > Y
Definition: multigrid.hpp:44
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:686
Array< Array< int > * > essentialTrueDofs
Definition: multigrid.hpp:169
Array< Operator * > prolongations
Definition: multigrid.hpp:138
Array< Operator * > operators
Definition: multigrid.hpp:35
Array< Solver * > smoothers
Definition: multigrid.hpp:36
Array< bool > ownedSmoothers
Definition: multigrid.hpp:38
Array2D< Vector * > R
Definition: multigrid.hpp:44
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:759
virtual void ArrayMult(const Array< const Vector *> &X, Array< Vector *> &Y) const
Operator application on a matrix: Y=A(X).
Definition: operator.cpp:66
double b
Definition: lissajous.cpp:42
const Solver * GetSmootherAtLevel(int level) const
Returns smoother at given level.
Definition: multigrid.hpp:95
const Operator * GetOperatorAtLevel(int level) const
Returns operator at given level.
Definition: multigrid.hpp:75
virtual void Mult(const Vector &x, Vector &y) const override
Application of the multigrid as a preconditioner.
Definition: multigrid.cpp:106
CycleType cycleType
Definition: multigrid.hpp:40
virtual void ArrayMult(const Array< const Vector *> &X_, Array< Vector *> &Y_) const override
Operator application on a matrix: Y=A(X).
Definition: multigrid.cpp:115
void SetCycleType(CycleType cycleType_, int preSmoothingSteps_, int postSmoothingSteps_)
Set cycle type and number of pre- and post-smoothing steps used by Mult.
Definition: multigrid.cpp:98
virtual ~MultigridBase()
Destructor.
Definition: multigrid.cpp:36
Array< bool > ownedProlongations
Definition: multigrid.hpp:139
virtual ~Multigrid()
Destructor.
Definition: multigrid.cpp:250
MultigridBase()
Constructs an empty multigrid hierarchy.
Definition: multigrid.cpp:17
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
virtual void ArrayMultTranspose(const Array< const Vector *> &X, Array< Vector *> &Y) const
Action of the transpose operator on a matrix: Y=A^t(X).
Definition: operator.cpp:78
void RecoverFineFEMSolution(const Vector &X, const Vector &b, Vector &x)
Recover the solution of a linear system formed with FormFineLinearSystem()
Definition: multigrid.cpp:285
int NumLevels() const
Returns the number of levels.
Definition: multigrid.hpp:69
void FormFineLinearSystem(Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)
Definition: multigrid.cpp:278
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition: array.hpp:864
virtual void ArrayAddMult(const Array< const Vector *> &X, Array< Vector *> &Y, const double a=1.0) const
Operator application on a matrix: Y+=A(X) (default) or Y+=a*A(X).
Definition: operator.cpp:90
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
Vector data type.
Definition: vector.hpp:58
void AddLevel(Operator *op, Solver *smoother, bool ownOperator, bool ownSmoother)
Adds a level to the multigrid operator hierarchy.
Definition: multigrid.cpp:87
Multigrid()
Constructs an empty multigrid hierarchy.
Definition: multigrid.cpp:234
Base class for solvers.
Definition: operator.hpp:682
GeometricMultigrid(const FiniteElementSpaceHierarchy &fespaces_)
Definition: multigrid.cpp:262
Abstract operator.
Definition: operator.hpp:24
Abstract base class for Multigrid solvers.
Definition: multigrid.hpp:25
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28