MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
ex4p.cpp
Go to the documentation of this file.
1// MFEM Example 4 - Parallel Version
2//
3// Compile with: make ex4p
4//
5// Sample runs: mpirun -np 4 ex4p -m ../data/square-disc.mesh
6// mpirun -np 4 ex4p -m ../data/star.mesh
7// mpirun -np 4 ex4p -m ../data/beam-tet.mesh
8// mpirun -np 4 ex4p -m ../data/beam-hex.mesh
9// mpirun -np 4 ex4p -m ../data/beam-hex.mesh -o 2 -pa
10// mpirun -np 4 ex4p -m ../data/escher.mesh -o 2 -sc
11// mpirun -np 4 ex4p -m ../data/fichera.mesh -o 2 -hb
12// mpirun -np 4 ex4p -m ../data/fichera-q2.vtk
13// mpirun -np 4 ex4p -m ../data/fichera-q3.mesh -o 2 -sc
14// mpirun -np 4 ex4p -m ../data/square-disc-nurbs.mesh -o 3
15// mpirun -np 4 ex4p -m ../data/beam-hex-nurbs.mesh -o 3
16// mpirun -np 4 ex4p -m ../data/periodic-square.mesh -no-bc
17// mpirun -np 4 ex4p -m ../data/periodic-cube.mesh -no-bc
18// mpirun -np 4 ex4p -m ../data/amr-quad.mesh
19// mpirun -np 3 ex4p -m ../data/amr-quad.mesh -o 2 -hb
20// mpirun -np 4 ex4p -m ../data/amr-hex.mesh -o 2 -sc
21// mpirun -np 4 ex4p -m ../data/amr-hex.mesh -o 2 -hb
22// mpirun -np 4 ex4p -m ../data/ref-prism.mesh -o 1
23// mpirun -np 4 ex4p -m ../data/octahedron.mesh -o 1
24// mpirun -np 4 ex4p -m ../data/star-surf.mesh -o 3 -hb
25//
26// Device sample runs:
27// mpirun -np 4 ex4p -m ../data/star.mesh -pa -d cuda
28// mpirun -np 4 ex4p -m ../data/star.mesh -pa -d raja-cuda
29// mpirun -np 4 ex4p -m ../data/star.mesh -pa -d raja-omp
30// mpirun -np 4 ex4p -m ../data/beam-hex.mesh -pa -d cuda
31//
32// Description: This example code solves a simple 2D/3D H(div) diffusion
33// problem corresponding to the second order definite equation
34// -grad(alpha div F) + beta F = f with boundary condition F dot n
35// = <given normal field>. Here, we use a given exact solution F
36// and compute the corresponding r.h.s. f. We discretize with
37// Raviart-Thomas finite elements.
38//
39// The example demonstrates the use of H(div) finite element
40// spaces with the grad-div and H(div) vector finite element mass
41// bilinear form, as well as the computation of discretization
42// error when the exact solution is known. Bilinear form
43// hybridization and static condensation are also illustrated.
44//
45// We recommend viewing examples 1-3 before viewing this example.
46
47#include "mfem.hpp"
48#include <fstream>
49#include <iostream>
50
51using namespace std;
52using namespace mfem;
53
54// Exact solution, F, and r.h.s., f. See below for implementation.
55void F_exact(const Vector &, Vector &);
56void f_exact(const Vector &, Vector &);
58
59int main(int argc, char *argv[])
60{
61 // 1. Initialize MPI and HYPRE.
62 Mpi::Init(argc, argv);
63 int num_procs = Mpi::WorldSize();
64 int myid = Mpi::WorldRank();
66
67 // 2. Parse command-line options.
68 const char *mesh_file = "../data/star.mesh";
69 int order = 1;
70 bool set_bc = true;
71 bool static_cond = false;
72 bool hybridization = false;
73 bool pa = false;
74 bool ea = false;
75 const char *device_config = "cpu";
76 bool visualization = 1;
77
78 OptionsParser args(argc, argv);
79 args.AddOption(&mesh_file, "-m", "--mesh",
80 "Mesh file to use.");
81 args.AddOption(&order, "-o", "--order",
82 "Finite element order (polynomial degree).");
83 args.AddOption(&set_bc, "-bc", "--impose-bc", "-no-bc", "--dont-impose-bc",
84 "Impose or not essential boundary conditions.");
85 args.AddOption(&freq, "-f", "--frequency", "Set the frequency for the exact"
86 " solution.");
87 args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
88 "--no-static-condensation", "Enable static condensation.");
89 args.AddOption(&hybridization, "-hb", "--hybridization", "-no-hb",
90 "--no-hybridization", "Enable hybridization.");
91 args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
92 "--no-partial-assembly", "Enable Partial Assembly.");
93 args.AddOption(&ea, "-ea", "--element-assembly", "-no-ea",
94 "--no-element-assembly", "Enable Element Assembly.");
95 args.AddOption(&device_config, "-d", "--device",
96 "Device configuration string, see Device::Configure().");
97 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
98 "--no-visualization",
99 "Enable or disable GLVis visualization.");
100 args.ParseCheck();
101 kappa = freq * M_PI;
102
103 // 3. Enable hardware devices such as GPUs, and programming models such as
104 // CUDA, OCCA, RAJA and OpenMP based on command line options.
105 Device device(device_config);
106 if (myid == 0) { device.Print(); }
107
108 // 4. Read the (serial) mesh from the given mesh file on all processors. We
109 // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
110 // and volume, as well as periodic meshes with the same code.
111 Mesh *mesh = new Mesh(mesh_file, 1, 1);
112 int dim = mesh->Dimension();
113 int sdim = mesh->SpaceDimension();
114
115 // 5. Refine the serial mesh on all processors to increase the resolution. In
116 // this example we do 'ref_levels' of uniform refinement. We choose
117 // 'ref_levels' to be the largest number that gives a final mesh with no
118 // more than 1,000 elements.
119 {
120 int ref_levels =
121 (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
122 for (int l = 0; l < ref_levels; l++)
123 {
124 mesh->UniformRefinement();
125 }
126 }
127
128 // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
129 // this mesh further in parallel to increase the resolution. Once the
130 // parallel mesh is defined, the serial mesh can be deleted.
131 ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
132 delete mesh;
133 {
134 int par_ref_levels = 2;
135 for (int l = 0; l < par_ref_levels; l++)
136 {
137 pmesh->UniformRefinement();
138 }
139 }
140
141 // 7. Define a parallel finite element space on the parallel mesh. Here we
142 // use the Raviart-Thomas finite elements of the specified order.
143 FiniteElementCollection *fec = new RT_FECollection(order-1, dim);
144 ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
145 HYPRE_BigInt size = fespace->GlobalTrueVSize();
146 if (myid == 0)
147 {
148 cout << "Number of finite element unknowns: " << size << endl;
149 }
150
151 // 8. Determine the list of true (i.e. parallel conforming) essential
152 // boundary dofs. In this example, the boundary conditions are defined
153 // by marking all the boundary attributes from the mesh as essential
154 // (Dirichlet) and converting them to a list of true dofs.
155 Array<int> ess_tdof_list;
156 if (pmesh->bdr_attributes.Size())
157 {
158 Array<int> ess_bdr(pmesh->bdr_attributes.Max());
159 ess_bdr = set_bc ? 1 : 0;
160 fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
161 }
162
163 // 9. Set up the parallel linear form b(.) which corresponds to the
164 // right-hand side of the FEM linear system, which in this case is
165 // (f,phi_i) where f is given by the function f_exact and phi_i are the
166 // basis functions in the finite element fespace.
168 ParLinearForm *b = new ParLinearForm(fespace);
169 b->AddDomainIntegrator(new VectorFEDomainLFIntegrator(f));
170 b->Assemble();
171
172 // 10. Define the solution vector x as a parallel finite element grid function
173 // corresponding to fespace. Initialize x by projecting the exact
174 // solution. Note that only values from the boundary faces will be used
175 // when eliminating the non-homogeneous boundary condition to modify the
176 // r.h.s. vector b.
177 ParGridFunction x(fespace);
180
181 // 11. Set up the parallel bilinear form corresponding to the H(div)
182 // diffusion operator grad alpha div + beta I, by adding the div-div and
183 // the mass domain integrators.
186 ParBilinearForm *a = new ParBilinearForm(fespace);
187 if (pa) { a->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
188 if (ea) { a->SetAssemblyLevel(AssemblyLevel::ELEMENT); }
189 a->AddDomainIntegrator(new DivDivIntegrator(*alpha));
190 a->AddDomainIntegrator(new VectorFEMassIntegrator(*beta));
191
192 // 12. Assemble the parallel bilinear form and the corresponding linear
193 // system, applying any necessary transformations such as: parallel
194 // assembly, eliminating boundary conditions, applying conforming
195 // constraints for non-conforming AMR, static condensation,
196 // hybridization, etc.
197 FiniteElementCollection *hfec = NULL;
198 ParFiniteElementSpace *hfes = NULL;
199 if (static_cond)
200 {
201 a->EnableStaticCondensation();
202 }
203 else if (hybridization)
204 {
205 hfec = new DG_Interface_FECollection(order-1, dim);
206 hfes = new ParFiniteElementSpace(pmesh, hfec);
207 a->EnableHybridization(hfes, new NormalTraceJumpIntegrator(),
208 ess_tdof_list);
209 }
210 a->Assemble();
211
212 OperatorPtr A;
213 Vector B, X;
214 a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
215
216 if (myid == 0 && !pa)
217 {
218 cout << "Size of linear system: "
219 << A.As<HypreParMatrix>()->GetGlobalNumRows() << endl;
220 }
221
222 // 13. Define and apply a parallel PCG solver for A X = B with the 2D AMS or
223 // the 3D ADS preconditioners from hypre. If using hybridization, the
224 // system is preconditioned with hypre's BoomerAMG. In the partial
225 // assembly case, use Jacobi preconditioning.
226 Solver *prec = NULL;
227 CGSolver *pcg = new CGSolver(MPI_COMM_WORLD);
228 pcg->SetOperator(*A);
229 pcg->SetRelTol(1e-12);
230 pcg->SetMaxIter(2000);
231 pcg->SetPrintLevel(1);
232 if (hybridization) { prec = new HypreBoomerAMG(*A.As<HypreParMatrix>()); }
233 else if (pa) { prec = new OperatorJacobiSmoother(*a, ess_tdof_list); }
234 else
235 {
236 ParFiniteElementSpace *prec_fespace =
237 (a->StaticCondensationIsEnabled() ? a->SCParFESpace() : fespace);
238 if (dim == 2) { prec = new HypreAMS(*A.As<HypreParMatrix>(), prec_fespace); }
239 else { prec = new HypreADS(*A.As<HypreParMatrix>(), prec_fespace); }
240 }
241 pcg->SetPreconditioner(*prec);
242 pcg->Mult(B, X);
243
244 // 14. Recover the parallel grid function corresponding to X. This is the
245 // local finite element solution on each processor.
246 a->RecoverFEMSolution(X, *b, x);
247
248 // 15. Compute and print the L^2 norm of the error.
249 {
250 real_t error = x.ComputeL2Error(F);
251 if (myid == 0)
252 {
253 cout << "\n|| F_h - F ||_{L^2} = " << error << '\n' << endl;
254 }
255 }
256
257 // 16. Save the refined mesh and the solution in parallel. This output can
258 // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
259 {
260 ostringstream mesh_name, sol_name;
261 mesh_name << "mesh." << setfill('0') << setw(6) << myid;
262 sol_name << "sol." << setfill('0') << setw(6) << myid;
263
264 ofstream mesh_ofs(mesh_name.str().c_str());
265 mesh_ofs.precision(8);
266 pmesh->Print(mesh_ofs);
267
268 ofstream sol_ofs(sol_name.str().c_str());
269 sol_ofs.precision(8);
270 x.Save(sol_ofs);
271 }
272
273 // 17. Send the solution by socket to a GLVis server.
274 if (visualization)
275 {
276 char vishost[] = "localhost";
277 int visport = 19916;
278 socketstream sol_sock(vishost, visport);
279 sol_sock << "parallel " << num_procs << " " << myid << "\n";
280 sol_sock.precision(8);
281 sol_sock << "solution\n" << *pmesh << x << flush;
282 }
283
284 // 18. Free the used memory.
285 delete pcg;
286 delete prec;
287 delete hfes;
288 delete hfec;
289 delete a;
290 delete alpha;
291 delete beta;
292 delete b;
293 delete fespace;
294 delete fec;
295 delete pmesh;
296
297 return 0;
298}
299
300
301// The exact solution (for non-surface meshes)
302void F_exact(const Vector &p, Vector &F)
303{
304 int dim = p.Size();
305
306 real_t x = p(0);
307 real_t y = p(1);
308 // real_t z = (dim == 3) ? p(2) : 0.0; // Uncomment if F is changed to depend on z
309
310 F(0) = cos(kappa*x)*sin(kappa*y);
311 F(1) = cos(kappa*y)*sin(kappa*x);
312 if (dim == 3)
313 {
314 F(2) = 0.0;
315 }
316}
317
318// The right hand side
319void f_exact(const Vector &p, Vector &f)
320{
321 int dim = p.Size();
322
323 real_t x = p(0);
324 real_t y = p(1);
325 // real_t z = (dim == 3) ? p(2) : 0.0; // Uncomment if f is changed to depend on z
326
327 real_t temp = 1 + 2*kappa*kappa;
328
329 f(0) = temp*cos(kappa*x)*sin(kappa*y);
330 f(1) = temp*cos(kappa*y)*sin(kappa*x);
331 if (dim == 3)
332 {
333 f(2) = 0;
334 }
335}
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:147
Conjugate gradient method.
Definition solvers.hpp:538
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the Conjugate Gradient method.
Definition solvers.cpp:751
void SetOperator(const Operator &op) override
Set/update the solver for the given operator.
Definition solvers.hpp:551
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
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:297
for Raviart-Thomas elements
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
The Auxiliary-space Divergence Solver in hypre.
Definition hypre.hpp:1948
The Auxiliary-space Maxwell Solver in hypre.
Definition hypre.hpp:1871
The BoomerAMG solver in hypre.
Definition hypre.hpp:1717
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:408
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.cpp:33
void SetRelTol(real_t rtol)
Definition solvers.hpp:229
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition solvers.cpp:174
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition solvers.cpp:72
void SetMaxIter(int max_it)
Definition solvers.hpp:231
Mesh data type.
Definition mesh.hpp:64
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:290
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1282
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1216
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition mesh.hpp:1219
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:11295
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).
Pointer to an Operator of a specified type.
Definition handle.hpp:34
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get().
Definition handle.hpp:104
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition solvers.hpp:333
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
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
HYPRE_BigInt GlobalTrueVSize() const
Definition pfespace.hpp:346
Class for parallel grid function.
Definition pgridfunc.hpp:50
void Save(std::ostream &out) const override
real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
Returns ||u_ex - u_h||_L2 in parallel for H1 or L2 elements.
void ProjectCoefficient(Coefficient &coeff) override
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Class for parallel linear form.
Class for parallel meshes.
Definition pmesh.hpp:34
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Definition pmesh.cpp:4800
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition fe_coll.hpp:403
Base class for solvers.
Definition operator.hpp:780
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition lininteg.hpp:344
A general vector function coefficient.
Vector data type.
Definition vector.hpp:82
int Size() const
Returns the size of the vector.
Definition vector.hpp:226
Vector beta
const real_t alpha
Definition ex15.cpp:369
int dim
Definition ex24.cpp:53
real_t kappa
Definition ex4p.cpp:57
void f_exact(const Vector &, Vector &)
Definition ex4p.cpp:319
real_t freq
Definition ex4p.cpp:57
void F_exact(const Vector &, Vector &)
Definition ex4p.cpp:302
int main()
HYPRE_Int HYPRE_BigInt
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
float real_t
Definition config.hpp:43
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
const char vishost[]
real_t p(const Vector &x, real_t t)