MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ex3p.cpp
Go to the documentation of this file.
1// MFEM Example 3 - Parallel Version
2//
3// Compile with: make ex3p
4//
5// Sample runs: mpirun -np 4 ex3p -m ../data/star.mesh
6// mpirun -np 4 ex3p -m ../data/square-disc.mesh -o 2
7// mpirun -np 4 ex3p -m ../data/beam-tet.mesh
8// mpirun -np 4 ex3p -m ../data/beam-tet.mesh -nc -o 2
9// mpirun -np 4 ex3p -m ../data/beam-hex.mesh
10// mpirun -np 4 ex3p -m ../data/beam-hex.mesh -o 2 -pa
11// mpirun -np 4 ex3p -m ../data/escher.mesh
12// mpirun -np 4 ex3p -m ../data/escher.mesh -o 2
13// mpirun -np 4 ex3p -m ../data/fichera.mesh
14// mpirun -np 4 ex3p -m ../data/fichera-q2.vtk
15// mpirun -np 4 ex3p -m ../data/fichera-q3.mesh
16// mpirun -np 4 ex3p -m ../data/square-disc-nurbs.mesh
17// mpirun -np 4 ex3p -m ../data/beam-hex-nurbs.mesh
18// mpirun -np 4 ex3p -m ../data/amr-quad.mesh -o 2
19// mpirun -np 4 ex3p -m ../data/amr-hex.mesh
20// mpirun -np 4 ex3p -m ../data/ref-prism.mesh -o 1
21// mpirun -np 4 ex3p -m ../data/octahedron.mesh -o 1
22// mpirun -np 4 ex3p -m ../data/star-surf.mesh -o 2
23// mpirun -np 4 ex3p -m ../data/mobius-strip.mesh -o 2 -f 0.1
24// mpirun -np 4 ex3p -m ../data/klein-bottle.mesh -o 2 -f 0.1
25//
26// Device sample runs:
27// mpirun -np 4 ex3p -m ../data/star.mesh -pa -d cuda
28// mpirun -np 4 ex3p -m ../data/star.mesh -no-pa -d cuda
29// mpirun -np 4 ex3p -m ../data/star.mesh -pa -d raja-cuda
30// mpirun -np 4 ex3p -m ../data/star.mesh -pa -d raja-omp
31// mpirun -np 4 ex3p -m ../data/beam-hex.mesh -pa -d cuda
32//
33// Description: This example code solves a simple electromagnetic diffusion
34// problem corresponding to the second order definite Maxwell
35// equation curl curl E + E = f with boundary condition
36// E x n = <given tangential field>. Here, we use a given exact
37// solution E and compute the corresponding r.h.s. f.
38// We discretize with Nedelec finite elements in 2D or 3D.
39//
40// The example demonstrates the use of H(curl) finite element
41// spaces with the curl-curl and the (vector finite element) mass
42// bilinear form, as well as the computation of discretization
43// error when the exact solution is known. Static condensation is
44// also illustrated.
45//
46// We recommend viewing examples 1-2 before viewing this example.
47
48#include "mfem.hpp"
49#include <fstream>
50#include <iostream>
51
52using namespace std;
53using namespace mfem;
54
55// Exact solution, E, and r.h.s., f. See below for implementation.
56void E_exact(const Vector &, Vector &);
57void f_exact(const Vector &, Vector &);
59int dim;
60
61int main(int argc, char *argv[])
62{
63 // 1. Initialize MPI and HYPRE.
64 Mpi::Init(argc, argv);
65 int num_procs = Mpi::WorldSize();
66 int myid = Mpi::WorldRank();
68
69 // 2. Parse command-line options.
70 const char *mesh_file = "../data/beam-tet.mesh";
71 int order = 1;
72 bool static_cond = false;
73 bool pa = false;
74 bool nc = false;
75 const char *device_config = "cpu";
76 bool visualization = true;
77#ifdef MFEM_USE_AMGX
78 bool useAmgX = false;
79#endif
80
81 OptionsParser args(argc, argv);
82 args.AddOption(&mesh_file, "-m", "--mesh",
83 "Mesh file to use.");
84 args.AddOption(&order, "-o", "--order",
85 "Finite element order (polynomial degree).");
86 args.AddOption(&freq, "-f", "--frequency", "Set the frequency for the exact"
87 " solution.");
88 args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
89 "--no-static-condensation", "Enable static condensation.");
90 args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
91 "--no-partial-assembly", "Enable Partial Assembly.");
92 args.AddOption(&nc, "-nc", "--non-conforming", "-c",
93 "--conforming",
94 "Mark the mesh as nonconforming before partitioning.");
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#ifdef MFEM_USE_AMGX
101 args.AddOption(&useAmgX, "-amgx", "--useAmgX", "-no-amgx",
102 "--no-useAmgX",
103 "Enable or disable AmgX in MatrixFreeAMS.");
104#endif
105
106 args.Parse();
107 if (!args.Good())
108 {
109 if (myid == 0)
110 {
111 args.PrintUsage(cout);
112 }
113 return 1;
114 }
115 if (myid == 0)
116 {
117 args.PrintOptions(cout);
118 }
119 kappa = freq * M_PI;
120
121 // 3. Enable hardware devices such as GPUs, and programming models such as
122 // CUDA, OCCA, RAJA and OpenMP based on command line options.
123 Device device(device_config);
124 if (myid == 0) { device.Print(); }
125
126 // 4. Read the (serial) mesh from the given mesh file on all processors. We
127 // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
128 // and volume meshes with the same code.
129 Mesh *mesh = new Mesh(mesh_file, 1, 1);
130 dim = mesh->Dimension();
131 int sdim = mesh->SpaceDimension();
132 if (nc)
133 {
134 // Can set to false to use conformal refinement for simplices.
135 mesh->EnsureNCMesh(true);
136 }
137
138 // 5. Refine the serial mesh on all processors to increase the resolution. In
139 // this example we do 'ref_levels' of uniform refinement. We choose
140 // 'ref_levels' to be the largest number that gives a final mesh with no
141 // more than 1,000 elements.
142 {
143 int ref_levels = (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
144 for (int l = 0; l < ref_levels; l++)
145 {
146 mesh->UniformRefinement();
147 }
148 }
149
150 // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
151 // this mesh further in parallel to increase the resolution. Once the
152 // parallel mesh is defined, the serial mesh can be deleted.
153 ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
154 delete mesh;
155 {
156 int par_ref_levels = 2;
157 for (int l = 0; l < par_ref_levels; l++)
158 {
159 pmesh->UniformRefinement();
160 }
161 }
162
163 // 7. Define a parallel finite element space on the parallel mesh. Here we
164 // use the Nedelec finite elements of the specified order.
166 ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
167 HYPRE_BigInt size = fespace->GlobalTrueVSize();
168 if (myid == 0)
169 {
170 cout << "Number of finite element unknowns: " << size << endl;
171 }
172
173 // 8. Determine the list of true (i.e. parallel conforming) essential
174 // boundary dofs. In this example, the boundary conditions are defined
175 // by marking all the boundary attributes from the mesh as essential
176 // (Dirichlet) and converting them to a list of true dofs.
177 Array<int> ess_tdof_list;
178 Array<int> ess_bdr;
179 if (pmesh->bdr_attributes.Size())
180 {
181 ess_bdr.SetSize(pmesh->bdr_attributes.Max());
182 ess_bdr = 1;
183 fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
184 }
185
186 // 9. Set up the parallel linear form b(.) which corresponds to the
187 // right-hand side of the FEM linear system, which in this case is
188 // (f,phi_i) where f is given by the function f_exact and phi_i are the
189 // basis functions in the finite element fespace.
191 ParLinearForm *b = new ParLinearForm(fespace);
192 b->AddDomainIntegrator(new VectorFEDomainLFIntegrator(f));
193 b->Assemble();
194
195 // 10. Define the solution vector x as a parallel finite element grid function
196 // corresponding to fespace. Initialize x by projecting the exact
197 // solution. Note that only values from the boundary edges will be used
198 // when eliminating the non-homogeneous boundary condition to modify the
199 // r.h.s. vector b.
200 ParGridFunction x(fespace);
203
204 // 11. Set up the parallel bilinear form corresponding to the EM diffusion
205 // operator curl muinv curl + sigma I, by adding the curl-curl and the
206 // mass domain integrators.
207 Coefficient *muinv = new ConstantCoefficient(1.0);
209 ParBilinearForm *a = new ParBilinearForm(fespace);
210 if (pa) { a->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
211 a->AddDomainIntegrator(new CurlCurlIntegrator(*muinv));
212 a->AddDomainIntegrator(new VectorFEMassIntegrator(*sigma));
213
214 // 12. Assemble the parallel bilinear form and the corresponding linear
215 // system, applying any necessary transformations such as: parallel
216 // assembly, eliminating boundary conditions, applying conforming
217 // constraints for non-conforming AMR, static condensation, etc.
218 if (static_cond) { a->EnableStaticCondensation(); }
219 a->Assemble();
220
221 OperatorPtr A;
222 Vector B, X;
223 a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
224
225 // 13. Solve the system AX=B using PCG with an AMS preconditioner.
226 if (pa)
227 {
228#ifdef MFEM_USE_AMGX
229 MatrixFreeAMS ams(*a, *A, *fespace, muinv, sigma, NULL, ess_bdr, useAmgX);
230#else
231 MatrixFreeAMS ams(*a, *A, *fespace, muinv, sigma, NULL, ess_bdr);
232#endif
233 CGSolver cg(MPI_COMM_WORLD);
234 cg.SetRelTol(1e-12);
235 cg.SetMaxIter(1000);
236 cg.SetPrintLevel(1);
237 cg.SetOperator(*A);
238 cg.SetPreconditioner(ams);
239 cg.Mult(B, X);
240 }
241 else
242 {
243 if (myid == 0)
244 {
245 cout << "Size of linear system: "
246 << A.As<HypreParMatrix>()->GetGlobalNumRows() << endl;
247 }
248
249 ParFiniteElementSpace *prec_fespace =
250 (a->StaticCondensationIsEnabled() ? a->SCParFESpace() : fespace);
251 HypreAMS ams(*A.As<HypreParMatrix>(), prec_fespace);
252 HyprePCG pcg(*A.As<HypreParMatrix>());
253 pcg.SetTol(1e-12);
254 pcg.SetMaxIter(500);
255 pcg.SetPrintLevel(2);
256 pcg.SetPreconditioner(ams);
257 pcg.Mult(B, X);
258 }
259
260 // 14. Recover the parallel grid function corresponding to X. This is the
261 // local finite element solution on each processor.
262 a->RecoverFEMSolution(X, *b, x);
263
264 // 15. Compute and print the L^2 norm of the error.
265 {
266 real_t error = x.ComputeL2Error(E);
267 if (myid == 0)
268 {
269 cout << "\n|| E_h - E ||_{L^2} = " << error << '\n' << endl;
270 }
271 }
272
273 // 16. Save the refined mesh and the solution in parallel. This output can
274 // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
275 {
276 ostringstream mesh_name, sol_name;
277 mesh_name << "mesh." << setfill('0') << setw(6) << myid;
278 sol_name << "sol." << setfill('0') << setw(6) << myid;
279
280 ofstream mesh_ofs(mesh_name.str().c_str());
281 mesh_ofs.precision(8);
282 pmesh->Print(mesh_ofs);
283
284 ofstream sol_ofs(sol_name.str().c_str());
285 sol_ofs.precision(8);
286 x.Save(sol_ofs);
287 }
288
289 // 17. Send the solution by socket to a GLVis server.
290 if (visualization)
291 {
292 char vishost[] = "localhost";
293 int visport = 19916;
294 socketstream sol_sock(vishost, visport);
295 sol_sock << "parallel " << num_procs << " " << myid << "\n";
296 sol_sock.precision(8);
297 sol_sock << "solution\n" << *pmesh << x << flush;
298 }
299
300 // 18. Free the used memory.
301 delete a;
302 delete sigma;
303 delete muinv;
304 delete b;
305 delete fespace;
306 delete fec;
307 delete pmesh;
308
309 return 0;
310}
311
312
313void E_exact(const Vector &x, Vector &E)
314{
315 if (dim == 3)
316 {
317 E(0) = sin(kappa * x(1));
318 E(1) = sin(kappa * x(2));
319 E(2) = sin(kappa * x(0));
320 }
321 else
322 {
323 E(0) = sin(kappa * x(1));
324 E(1) = sin(kappa * x(0));
325 if (x.Size() == 3) { E(2) = 0.0; }
326 }
327}
328
329void f_exact(const Vector &x, Vector &f)
330{
331 if (dim == 3)
332 {
333 f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
334 f(1) = (1. + kappa * kappa) * sin(kappa * x(2));
335 f(2) = (1. + kappa * kappa) * sin(kappa * x(0));
336 }
337 else
338 {
339 f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
340 f(1) = (1. + kappa * kappa) * sin(kappa * x(0));
341 if (x.Size() == 3) { f(2) = 0.0; }
342 }
343}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
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
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
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.
Integrator for for Nedelec elements.
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
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 Maxwell Solver in hypre.
Definition hypre.hpp:1845
PCG solver in hypre.
Definition hypre.hpp:1275
void SetPrintLevel(int print_lvl)
Definition hypre.cpp:4156
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition hypre.cpp:4161
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's PCG.
Definition hypre.cpp:4184
void SetMaxIter(int max_iter)
Definition hypre.cpp:4146
void SetTol(real_t tol)
Definition hypre.cpp:4136
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
An auxiliary Maxwell solver for a high-order curl-curl system without high-order assembly.
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
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition mesh.hpp:1163
void EnsureNCMesh(bool simplices_nonconforming=false)
Definition mesh.cpp:10626
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).
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition fe_coll.hpp:465
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
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.
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:285
Class for parallel grid function.
Definition pgridfunc.hpp:33
void Save(std::ostream &out) const override
real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
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:4801
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition lininteg.hpp:347
A general vector function coefficient.
Vector data type.
Definition vector.hpp:80
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
real_t sigma(const Vector &x)
Definition maxwell.cpp:102
real_t kappa
Definition ex3p.cpp:58
int dim
Definition ex3p.cpp:59
void E_exact(const Vector &, Vector &)
Definition ex3p.cpp:313
void f_exact(const Vector &, Vector &)
Definition ex3p.cpp:329
real_t freq
Definition ex3p.cpp:58
int main()
HYPRE_Int HYPRE_BigInt
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
const int visport
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[]