MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ex0p.cpp
Go to the documentation of this file.
1// MFEM Example 0 - Parallel Version
2//
3// Compile with: make ex0p
4//
5// Sample runs: mpirun -np 4 ex0p
6// mpirun -np 4 ex0p -m ../data/fichera.mesh
7// mpirun -np 4 ex0p -m ../data/square-disc.mesh -o 2
8//
9// Description: This example code demonstrates the most basic parallel usage of
10// MFEM to define a simple finite element discretization of the
11// Laplace problem -Delta u = 1 with zero Dirichlet boundary
12// conditions. General 2D/3D serial mesh files and finite element
13// polynomial degrees can be specified by command line options.
14
15#include "mfem.hpp"
16#include <fstream>
17#include <iostream>
18
19using namespace std;
20using namespace mfem;
21
22int main(int argc, char *argv[])
23{
24 // 1. Initialize MPI and HYPRE.
25 Mpi::Init(argc, argv);
27
28 // 2. Parse command line options.
29 string mesh_file = "../data/star.mesh";
30 int order = 1;
31
32 OptionsParser args(argc, argv);
33 args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use.");
34 args.AddOption(&order, "-o", "--order", "Finite element polynomial degree");
35 args.ParseCheck();
36
37 // 3. Read the serial mesh from the given mesh file.
38 Mesh serial_mesh(mesh_file);
39
40 // 4. Define a parallel mesh by a partitioning of the serial mesh. Refine
41 // this mesh once in parallel to increase the resolution.
42 ParMesh mesh(MPI_COMM_WORLD, serial_mesh);
43 serial_mesh.Clear(); // the serial mesh is no longer needed
44 mesh.UniformRefinement();
45
46 // 5. Define a finite element space on the mesh. Here we use H1 continuous
47 // high-order Lagrange finite elements of the given order.
48 H1_FECollection fec(order, mesh.Dimension());
49 ParFiniteElementSpace fespace(&mesh, &fec);
50 HYPRE_BigInt total_num_dofs = fespace.GlobalTrueVSize();
51 if (Mpi::Root())
52 {
53 cout << "Number of unknowns: " << total_num_dofs << endl;
54 }
55
56 // 6. Extract the list of all the boundary DOFs. These will be marked as
57 // Dirichlet in order to enforce zero boundary conditions.
58 Array<int> boundary_dofs;
59 fespace.GetBoundaryTrueDofs(boundary_dofs);
60
61 // 7. Define the solution x as a finite element grid function in fespace. Set
62 // the initial guess to zero, which also sets the boundary conditions.
63 ParGridFunction x(&fespace);
64 x = 0.0;
65
66 // 8. Set up the linear form b(.) corresponding to the right-hand side.
67 ConstantCoefficient one(1.0);
68 ParLinearForm b(&fespace);
69 b.AddDomainIntegrator(new DomainLFIntegrator(one));
70 b.Assemble();
71
72 // 9. Set up the bilinear form a(.,.) corresponding to the -Delta operator.
73 ParBilinearForm a(&fespace);
74 a.AddDomainIntegrator(new DiffusionIntegrator);
75 a.Assemble();
76
77 // 10. Form the linear system A X = B. This includes eliminating boundary
78 // conditions, applying AMR constraints, parallel assembly, etc.
80 Vector B, X;
81 a.FormLinearSystem(boundary_dofs, x, b, A, X, B);
82
83 // 11. Solve the system using PCG with hypre's BoomerAMG preconditioner.
84 HypreBoomerAMG M(A);
85 CGSolver cg(MPI_COMM_WORLD);
86 cg.SetRelTol(1e-12);
87 cg.SetMaxIter(2000);
88 cg.SetPrintLevel(1);
90 cg.SetOperator(A);
91 cg.Mult(B, X);
92
93 // 12. Recover the solution x as a grid function and save to file. The output
94 // can be viewed using GLVis as follows: "glvis -np <np> -m mesh -g sol"
95 a.RecoverFEMSolution(X, b, x);
96 x.Save("sol");
97 mesh.Save("mesh");
98
99 return 0;
100}
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.
Class for domain integration .
Definition lininteg.hpp:109
void GetBoundaryTrueDofs(Array< int > &boundary_dofs, int component=-1)
Get a list of all boundary true dofs, boundary_dofs. For spaces with 'vdim' > 1, the 'component' para...
Definition fespace.cpp:632
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
The BoomerAMG solver in hypre.
Definition hypre.hpp:1691
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
Mesh data type.
Definition mesh.hpp:56
void Clear()
Clear the contents of the Mesh.
Definition mesh.hpp:730
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 bool Root()
Return true if the rank in MPI_COMM_WORLD is zero.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
void ParseCheck(std::ostream &out=mfem::out)
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
Class for parallel bilinear form.
Abstract parallel finite element space.
Definition pfespace.hpp:29
HYPRE_BigInt GlobalTrueVSize() const
Definition pfespace.hpp:285
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
void Save(const std::string &fname, int precision=16) const override
Definition pmesh.cpp:4940
Vector data type.
Definition vector.hpp:80
int main()
HYPRE_Int HYPRE_BigInt
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41