MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ex26p.cpp
Go to the documentation of this file.
1// MFEM Example 26 - Parallel Version
2//
3// Compile with: make ex26p
4//
5// Sample runs: mpirun -np 4 ex26p -m ../data/star.mesh
6// mpirun -np 4 ex26p -m ../data/fichera.mesh
7// mpirun -np 4 ex26p -m ../data/beam-hex.mesh
8//
9// Device sample runs:
10// mpirun -np 4 ex26p -d cuda
11// mpirun -np 4 ex26p -d occa-cuda
12// mpirun -np 4 ex26p -d raja-omp
13// mpirun -np 4 ex26p -d ceed-cpu
14// mpirun -np 4 ex26p -d ceed-cuda
15//
16// Description: This example code demonstrates the use of MFEM to define a
17// simple finite element discretization of the Laplace problem
18// -Delta u = 1 with homogeneous Dirichlet boundary conditions
19// as in Example 1.
20//
21// It highlights on the creation of a hierarchy of discretization
22// spaces with partial assembly and the construction of an
23// efficient multigrid preconditioner for the iterative solver.
24//
25// We recommend viewing Example 1 before viewing this example.
26
27#include "mfem.hpp"
28#include <fstream>
29#include <iostream>
30
31using namespace std;
32using namespace mfem;
33
34// Class for constructing a multigrid preconditioner for the diffusion operator.
35// This example multigrid preconditioner class demonstrates the creation of the
36// parallel diffusion bilinear forms and operators using partial assembly for
37// all spaces except the coarsest one in the ParFiniteElementSpaceHierarchy.
38// The multigrid uses a PCG solver preconditioned with AMG on the coarsest level
39// and second order Chebyshev accelerated smoothers on the other levels.
40class DiffusionMultigrid : public GeometricMultigrid
41{
42private:
44 HypreBoomerAMG* amg;
45
46public:
47 // Constructs a diffusion multigrid for the ParFiniteElementSpaceHierarchy
48 // and the array of essential boundaries
49 DiffusionMultigrid(ParFiniteElementSpaceHierarchy& fespaces,
50 Array<int>& ess_bdr)
51 : GeometricMultigrid(fespaces, ess_bdr), coeff(1.0)
52 {
53 ConstructCoarseOperatorAndSolver(fespaces.GetFESpaceAtLevel(0));
54
55 for (int level = 1; level < fespaces.GetNumLevels(); ++level)
56 {
57 ConstructOperatorAndSmoother(fespaces.GetFESpaceAtLevel(level), level);
58 }
59 }
60
61 virtual ~DiffusionMultigrid()
62 {
63 delete amg;
64 }
65
66private:
67 void ConstructBilinearForm(ParFiniteElementSpace& fespace,
68 bool partial_assembly)
69 {
70 ParBilinearForm* form = new ParBilinearForm(&fespace);
71 if (partial_assembly)
72 {
73 form->SetAssemblyLevel(AssemblyLevel::PARTIAL);
74 }
76 form->Assemble();
77 bfs.Append(form);
78 }
79
80 void ConstructCoarseOperatorAndSolver(ParFiniteElementSpace& coarse_fespace)
81 {
82 ConstructBilinearForm(coarse_fespace, false);
83
84 HypreParMatrix* hypreCoarseMat = new HypreParMatrix();
85 bfs[0]->FormSystemMatrix(*essentialTrueDofs[0], *hypreCoarseMat);
86
87 amg = new HypreBoomerAMG(*hypreCoarseMat);
88 amg->SetPrintLevel(-1);
89
90 CGSolver* pcg = new CGSolver(MPI_COMM_WORLD);
91 pcg->SetPrintLevel(-1);
92 pcg->SetMaxIter(10);
93 pcg->SetRelTol(sqrt(1e-4));
94 pcg->SetAbsTol(0.0);
95 pcg->SetOperator(*hypreCoarseMat);
96 pcg->SetPreconditioner(*amg);
97
98 AddLevel(hypreCoarseMat, pcg, true, true);
99 }
100
101 void ConstructOperatorAndSmoother(ParFiniteElementSpace& fespace, int level)
102 {
103 const Array<int> &ess_tdof_list = *essentialTrueDofs[level];
104 ConstructBilinearForm(fespace, true);
105
106 OperatorPtr opr;
108 bfs.Last()->FormSystemMatrix(ess_tdof_list, opr);
109 opr.SetOperatorOwner(false);
110
111 Vector diag(fespace.GetTrueVSize());
112 bfs.Last()->AssembleDiagonal(diag);
113
114 Solver* smoother = new OperatorChebyshevSmoother(
115 *opr, diag, ess_tdof_list, 2, fespace.GetParMesh()->GetComm());
116
117 AddLevel(opr.Ptr(), smoother, true, true);
118 }
119};
120
121
122int main(int argc, char *argv[])
123{
124 // 1. Initialize MPI and HYPRE.
125 Mpi::Init(argc, argv);
126 int num_procs = Mpi::WorldSize();
127 int myid = Mpi::WorldRank();
128 Hypre::Init();
129
130 // 2. Parse command-line options.
131 const char *mesh_file = "../data/star.mesh";
132 int geometric_refinements = 0;
133 int order_refinements = 2;
134 const char *device_config = "cpu";
135 bool visualization = true;
136
137 OptionsParser args(argc, argv);
138 args.AddOption(&mesh_file, "-m", "--mesh",
139 "Mesh file to use.");
140 args.AddOption(&geometric_refinements, "-gr", "--geometric-refinements",
141 "Number of geometric refinements done prior to order refinements.");
142 args.AddOption(&order_refinements, "-or", "--order-refinements",
143 "Number of order refinements. Finest level in the hierarchy has order 2^{or}.");
144 args.AddOption(&device_config, "-d", "--device",
145 "Device configuration string, see Device::Configure().");
146 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
147 "--no-visualization",
148 "Enable or disable GLVis visualization.");
149 args.Parse();
150 if (!args.Good())
151 {
152 if (myid == 0)
153 {
154 args.PrintUsage(cout);
155 }
156 return 1;
157 }
158 if (myid == 0)
159 {
160 args.PrintOptions(cout);
161 }
162
163 // 3. Enable hardware devices such as GPUs, and programming models such as
164 // CUDA, OCCA, RAJA and OpenMP based on command line options.
165 Device device(device_config);
166 if (myid == 0) { device.Print(); }
167
168 // 4. Read the (serial) mesh from the given mesh file on all processors. We
169 // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
170 // and volume meshes with the same code.
171 Mesh *mesh = new Mesh(mesh_file, 1, 1);
172 int dim = mesh->Dimension();
173
174 // 5. Refine the serial mesh on all processors to increase the resolution. In
175 // this example we do 'ref_levels' of uniform refinement. We choose
176 // 'ref_levels' to be the largest number that gives a final mesh with no
177 // more than 1,000 elements.
178 {
179 int ref_levels =
180 (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
181 for (int l = 0; l < ref_levels; l++)
182 {
183 mesh->UniformRefinement();
184 }
185 }
186
187 // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
188 // this mesh further in parallel to increase the resolution. Once the
189 // parallel mesh is defined, the serial mesh can be deleted.
190 ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
191 delete mesh;
192 {
193 int par_ref_levels = 2;
194 for (int l = 0; l < par_ref_levels; l++)
195 {
196 pmesh->UniformRefinement();
197 }
198 }
199
200 // 7. Define a parallel finite element space hierarchy on the parallel mesh.
201 // Here we use continuous Lagrange finite elements. We start with order 1
202 // on the coarse level and geometrically refine the spaces by the specified
203 // amount. Afterwards, we increase the order of the finite elements by a
204 // factor of 2 for each additional level.
206 ParFiniteElementSpace *coarse_fespace = new ParFiniteElementSpace(pmesh, fec);
207
209 collections.Append(fec);
211 pmesh, coarse_fespace, true, true);
212 for (int level = 0; level < geometric_refinements; ++level)
213 {
214 fespaces->AddUniformlyRefinedLevel();
215 }
216 for (int level = 0; level < order_refinements; ++level)
217 {
218 collections.Append(new H1_FECollection((int)std::pow(2, level+1), dim));
219 fespaces->AddOrderRefinedLevel(collections.Last());
220 }
221
222 HYPRE_BigInt size = fespaces->GetFinestFESpace().GlobalTrueVSize();
223 if (myid == 0)
224 {
225 cout << "Number of finite element unknowns: " << size << endl;
226 }
227
228 // 8. Set up the parallel linear form b(.) which corresponds to the
229 // right-hand side of the FEM linear system, which in this case is
230 // (1,phi_i) where phi_i are the basis functions in fespace.
231 ParLinearForm *b = new ParLinearForm(&fespaces->GetFinestFESpace());
232 ConstantCoefficient one(1.0);
233 b->AddDomainIntegrator(new DomainLFIntegrator(one));
234 b->Assemble();
235
236 // 9. Define the solution vector x as a parallel finite element grid function
237 // corresponding to fespace. Initialize x with initial guess of zero,
238 // which satisfies the boundary conditions.
239 ParGridFunction x(&fespaces->GetFinestFESpace());
240 x = 0.0;
241
242 // 10. Create the multigrid operator using the previously created parallel
243 // FiniteElementSpaceHierarchy and additional boundary information. This
244 // operator is then used to create the MultigridSolver as preconditioner
245 // in the iterative solver.
246 Array<int> ess_bdr(pmesh->bdr_attributes.Max());
247 if (pmesh->bdr_attributes.Size())
248 {
249 ess_bdr = 1;
250 }
251
252 DiffusionMultigrid* M = new DiffusionMultigrid(*fespaces, ess_bdr);
253 M->SetCycleType(Multigrid::CycleType::VCYCLE, 1, 1);
254
255 OperatorPtr A;
256 Vector X, B;
257 M->FormFineLinearSystem(x, *b, A, X, B);
258
259 // 11. Solve the linear system A X = B.
260 CGSolver cg(MPI_COMM_WORLD);
261 cg.SetRelTol(1e-12);
262 cg.SetMaxIter(2000);
263 cg.SetPrintLevel(1);
264 cg.SetOperator(*A);
265 cg.SetPreconditioner(*M);
266 cg.Mult(B, X);
267
268 // 12. Recover the parallel grid function corresponding to X. This is the
269 // local finite element solution on each processor.
270 M->RecoverFineFEMSolution(X, *b, x);
271
272 // 13. Save the refined mesh and the solution in parallel. This output can be
273 // viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
274 {
275 ostringstream mesh_name, sol_name;
276 mesh_name << "mesh." << setfill('0') << setw(6) << myid;
277 sol_name << "sol." << setfill('0') << setw(6) << myid;
278
279 ofstream mesh_ofs(mesh_name.str().c_str());
280 mesh_ofs.precision(8);
281 fespaces->GetFinestFESpace().GetParMesh()->Print(mesh_ofs);
282
283 ofstream sol_ofs(sol_name.str().c_str());
284 sol_ofs.precision(8);
285 x.Save(sol_ofs);
286 }
287
288 // 14. Send the solution by socket to a GLVis server.
289 if (visualization)
290 {
291 char vishost[] = "localhost";
292 int visport = 19916;
293 socketstream sol_sock(vishost, visport);
294 sol_sock << "parallel " << num_procs << " " << myid << "\n";
295 sol_sock.precision(8);
296 sol_sock << "solution\n" << *fespaces->GetFinestFESpace().GetParMesh()
297 << x << flush;
298 }
299
300 // 15. Free the used memory.
301 delete M;
302 delete b;
303 delete fespaces;
304 for (int level = 0; level < collections.Size(); ++level)
305 {
306 delete collections[level];
307 }
308
309 return 0;
310}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
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
T & Last()
Return the last element in the array.
Definition array.hpp:802
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
Conjugate gradient method.
Definition solvers.hpp:513
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition solvers.hpp:526
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the Conjugate Gradient method.
Definition solvers.cpp:718
A coefficient that is constant across space and time.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition device.hpp:123
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition device.cpp:286
Class for domain integration .
Definition lininteg.hpp:109
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
int GetNumLevels() const
Returns the number of levels in the hierarchy.
virtual const FiniteElementSpace & GetFESpaceAtLevel(int level) const
Returns the finite element space at the given level.
Geometric multigrid associated with a hierarchy of finite element spaces.
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)
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
The BoomerAMG solver in hypre.
Definition hypre.hpp:1691
void SetPrintLevel(int print_level)
Definition hypre.hpp:1771
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.hpp:74
void SetRelTol(real_t rtol)
Definition solvers.hpp:209
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition solvers.cpp:173
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition solvers.cpp:71
void SetMaxIter(int max_it)
Definition solvers.hpp:211
void SetAbsTol(real_t atol)
Definition solvers.hpp:210
Mesh data type.
Definition mesh.hpp:56
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:282
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1226
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static int WorldSize()
Return the size of MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
void AddLevel(Operator *op, Solver *smoother, bool ownOperator, bool ownSmoother)
Adds a level to the multigrid operator hierarchy.
Definition multigrid.cpp:87
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
Chebyshev accelerated smoothing with given vector, no matrix necessary.
Definition solvers.hpp:394
Pointer to an Operator of a specified type.
Definition handle.hpp:34
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition handle.hpp:120
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
Definition handle.hpp:132
Operator * Ptr() const
Access the underlying Operator pointer.
Definition handle.hpp:87
@ ANY_TYPE
ID for the base class Operator, i.e. any type.
Definition operator.hpp:285
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set 'var' to receive the value. Enable/disable tags are used to set the bool...
Definition optparser.hpp:82
bool Good() const
Return true if the command line options were parsed successfully.
Class for parallel bilinear form.
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void AddUniformlyRefinedLevel(int dim=1, int ordering=Ordering::byVDIM) override
Adds one level to the hierarchy by uniformly refining the mesh on the previous level.
void AddOrderRefinedLevel(FiniteElementCollection *fec, int dim=1, int ordering=Ordering::byVDIM) override
Adds one level to the hierarchy by using a different finite element order defined through FiniteEleme...
const ParFiniteElementSpace & GetFinestFESpace() const override
Returns the finite element space at the finest level.
Abstract parallel finite element space.
Definition pfespace.hpp:29
HYPRE_BigInt GlobalTrueVSize() const
Definition pfespace.hpp:285
int GetTrueVSize() const override
Return the number of local vector true dofs.
Definition pfespace.hpp:289
ParMesh * GetParMesh() const
Definition pfespace.hpp:277
Class for parallel grid function.
Definition pgridfunc.hpp:33
void Save(std::ostream &out) const override
Class for parallel linear form.
Class for parallel meshes.
Definition pmesh.hpp:34
MPI_Comm GetComm() const
Definition pmesh.hpp:402
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Definition pmesh.cpp:4801
Base class for solvers.
Definition operator.hpp:683
Vector data type.
Definition vector.hpp:80
int dim
Definition ex24.cpp:53
int main()
HYPRE_Int HYPRE_BigInt
real_t b
Definition lissajous.cpp:42
const int visport
const char vishost[]