MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
multigrid.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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
14namespace 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
52void 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
73void 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
98void MultigridBase::SetCycleType(CycleType cycleType_, int preSmoothingSteps_,
99 int postSmoothingSteps_)
100{
101 cycleType = cycleType_;
102 preSmoothingSteps = preSmoothingSteps_;
103 postSmoothingSteps = postSmoothingSteps_;
104}
105
106void MultigridBase::Mult(const Vector& x, Vector& y) const
107{
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
147void 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 {
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
179void 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 const Array<Solver*>& smoothers_,
236 const Array<Operator*>& prolongations_,
237 const Array<bool>& ownedOperators_,
238 const Array<bool>& ownedSmoothers_,
239 const Array<bool>& ownedProlongations_)
240 : MultigridBase(operators_, smoothers_, ownedOperators_, ownedSmoothers_)
241{
242 prolongations_.Copy(prolongations);
243 ownedProlongations_.Copy(ownedProlongations);
244}
245
247{
248 for (int i = 0; i < prolongations.Size(); ++i)
249 {
250 if (ownedProlongations[i])
251 {
252 delete prolongations[i];
253 }
254 }
255}
256
258 const FiniteElementSpaceHierarchy& fespaces_)
259 : fespaces(fespaces_)
260{
261 const int nlevels = fespaces.GetNumLevels();
262 ownedProlongations.SetSize(nlevels - 1);
263 ownedProlongations = false;
264
265 prolongations.SetSize(nlevels - 1);
266 for (int level = 0; level < nlevels - 1; ++level)
267 {
269 }
270}
271
273 const FiniteElementSpaceHierarchy& fespaces_,
274 const Array<int> &ess_bdr)
275 : fespaces(fespaces_)
276{
277 bool have_ess_bdr = false;
278 for (int i = 0; i < ess_bdr.Size(); i++)
279 {
280 if (ess_bdr[i]) { have_ess_bdr = true; break; }
281 }
282
283 const int nlevels = fespaces.GetNumLevels();
284 ownedProlongations.SetSize(nlevels - 1);
285 ownedProlongations = have_ess_bdr;
286
287 if (have_ess_bdr)
288 {
289 essentialTrueDofs.SetSize(nlevels);
290 for (int level = 0; level < nlevels; ++level)
291 {
292 essentialTrueDofs[level] = new Array<int>;
294 ess_bdr, *essentialTrueDofs[level]);
295 }
296 }
297
298 prolongations.SetSize(nlevels - 1);
299 for (int level = 0; level < nlevels - 1; ++level)
300 {
301 if (have_ess_bdr)
302 {
305 *essentialTrueDofs[level],
306 *essentialTrueDofs[level + 1]
307 );
308 }
309 else
310 {
312 }
313 }
314}
315
317{
318 for (int i = 0; i < bfs.Size(); ++i)
319 {
320 delete bfs[i];
321 }
322 for (int i = 0; i < essentialTrueDofs.Size(); ++i)
323 {
324 delete essentialTrueDofs[i];
325 }
326}
327
330 Vector& X, Vector& B)
331{
332 bfs.Last()->FormLinearSystem(*essentialTrueDofs.Last(), x, b, A, X, B);
333}
334
336 const Vector& b, Vector& x)
337{
338 bfs.Last()->RecoverFEMSolution(X, b, x);
339}
340
341} // namespace mfem
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
int Size() const
Return the logical size of the array.
Definition array.hpp:144
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:769
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition array.hpp:874
int GetNumLevels() const
Returns the number of levels in the hierarchy.
Operator * GetProlongationAtLevel(int level) const
Returns the prolongation operator from the finite element space at level to the finite element space ...
virtual const FiniteElementSpace & GetFESpaceAtLevel(int level) const
Returns the finite element space at the given level.
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Definition fespace.cpp:588
MFEM_DEPRECATED GeometricMultigrid(const FiniteElementSpaceHierarchy &fespaces_)
Deprecated.
void RecoverFineFEMSolution(const Vector &X, const Vector &b, Vector &x)
Recover the solution of a linear system formed with FormFineLinearSystem()
Array< Array< int > * > essentialTrueDofs
Array< BilinearForm * > bfs
const FiniteElementSpaceHierarchy & fespaces
void FormFineLinearSystem(Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)
virtual ~GeometricMultigrid()
Destructor.
Abstract base class for Multigrid solvers.
Definition multigrid.hpp:26
virtual void ArrayMult(const Array< const Vector * > &X_, Array< Vector * > &Y_) const override
Operator application on a matrix: Y=A(X).
Array2D< Vector * > R
Definition multigrid.hpp:44
Array2D< Vector * > Y
Definition multigrid.hpp:44
int NumLevels() const
Returns the number of levels.
Definition multigrid.hpp:69
virtual ~MultigridBase()
Destructor.
Definition multigrid.cpp:36
Array< Operator * > operators
Definition multigrid.hpp:35
Array2D< Vector * > Z
Definition multigrid.hpp:44
Array< bool > ownedSmoothers
Definition multigrid.hpp:38
const Operator * GetOperatorAtLevel(int level) const
Returns operator at given level.
Definition multigrid.hpp:75
Array< bool > ownedOperators
Definition multigrid.hpp:37
Array2D< Vector * > X
Definition multigrid.hpp:44
void AddLevel(Operator *op, Solver *smoother, bool ownOperator, bool ownSmoother)
Adds a level to the multigrid operator hierarchy.
Definition multigrid.cpp:87
Array< Solver * > smoothers
Definition multigrid.hpp:36
MultigridBase()
Constructs an empty multigrid hierarchy.
Definition multigrid.cpp:17
virtual void Mult(const Vector &x, Vector &y) const override
Application of the multigrid as a preconditioner.
const Solver * GetSmootherAtLevel(int level) const
Returns smoother at given level.
Definition multigrid.hpp:95
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
Multigrid()=default
Constructs an empty multigrid hierarchy.
Array< bool > ownedProlongations
Array< Operator * > prolongations
virtual ~Multigrid()
Destructor.
Pointer to an Operator of a specified type.
Definition handle.hpp:34
Abstract operator.
Definition operator.hpp:25
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0)
Form a constrained linear system using a matrix-free approach.
Definition operator.cpp:114
int width
Dimension of the input / number of columns in the matrix.
Definition operator.hpp:28
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
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 ArrayAddMult(const Array< const Vector * > &X, Array< Vector * > &Y, const real_t a=1.0) const
Operator application on a matrix: Y+=A(X) (default) or Y+=a*A(X).
Definition operator.cpp:90
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition operator.hpp:72
virtual void ArrayMult(const Array< const Vector * > &X, Array< Vector * > &Y) const
Operator application on a matrix: Y=A(X).
Definition operator.cpp:66
Rectangular Operator for imposing essential boundary conditions on the input space using only the act...
Definition operator.hpp:970
Base class for solvers.
Definition operator.hpp:683
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition operator.hpp:686
Vector data type.
Definition vector.hpp:80
real_t b
Definition lissajous.cpp:42