MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ex6p.cpp
Go to the documentation of this file.
1// MFEM Example 6 - Parallel Version
2// PETSc Modification
3//
4// Compile with: make ex6p
5//
6// Sample runs:
7// mpirun -np 4 ex6p -m ../../data/amr-quad.mesh
8// mpirun -np 4 ex6p -m ../../data/amr-quad.mesh -nonoverlapping
9//
10// Description: This is a version of Example 1 with a simple adaptive mesh
11// refinement loop. The problem being solved is again the Laplace
12// equation -Delta u = 1 with homogeneous Dirichlet boundary
13// conditions. The problem is solved on a sequence of meshes which
14// are locally refined in a conforming (triangles, tetrahedrons)
15// or non-conforming (quadrilaterals, hexahedra) manner according
16// to a simple ZZ error estimator.
17//
18// The example demonstrates MFEM's capability to work with both
19// conforming and nonconforming refinements, in 2D and 3D, on
20// linear, curved and surface meshes. Interpolation of functions
21// from coarse to fine meshes, as well as persistent GLVis
22// visualization are also illustrated.
23//
24// PETSc assembly timings can be benchmarked if requested by
25// command line.
26//
27// We recommend viewing Example 1 before viewing this example.
28
29#include "mfem.hpp"
30#include <fstream>
31#include <iostream>
32
33#ifndef MFEM_USE_PETSC
34#error This example requires that MFEM is built with MFEM_USE_PETSC=YES
35#endif
36
37using namespace std;
38using namespace mfem;
39
40int main(int argc, char *argv[])
41{
42 // 1. Initialize MPI and HYPRE.
43 Mpi::Init(argc, argv);
44 int num_procs = Mpi::WorldSize();
45 int myid = Mpi::WorldRank();
47
48 // 2. Parse command-line options.
49 const char *mesh_file = "../../data/star.mesh";
50 int order = 1;
51 bool visualization = true;
52 int max_dofs = 100000;
53 bool use_petsc = true;
54 const char *petscrc_file = "";
55 bool use_nonoverlapping = false;
56
57 OptionsParser args(argc, argv);
58 args.AddOption(&mesh_file, "-m", "--mesh",
59 "Mesh file to use.");
60 args.AddOption(&order, "-o", "--order",
61 "Finite element order (polynomial degree).");
62 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
63 "--no-visualization",
64 "Enable or disable GLVis visualization.");
65 args.AddOption(&max_dofs, "-md", "--max_dofs",
66 "Maximum number of dofs.");
67 args.AddOption(&use_petsc, "-usepetsc", "--usepetsc", "-no-petsc",
68 "--no-petsc",
69 "Use or not PETSc to solve the linear system.");
70 args.AddOption(&petscrc_file, "-petscopts", "--petscopts",
71 "PetscOptions file to use.");
72 args.AddOption(&use_nonoverlapping, "-nonoverlapping", "--nonoverlapping",
73 "-no-nonoverlapping", "--no-nonoverlapping",
74 "Use or not the block diagonal PETSc's matrix format "
75 "for non-overlapping domain decomposition.");
76 args.Parse();
77 if (!args.Good())
78 {
79 if (myid == 0)
80 {
81 args.PrintUsage(cout);
82 }
83 return 1;
84 }
85 if (myid == 0)
86 {
87 args.PrintOptions(cout);
88 }
89 // 2b. We initialize PETSc
90 if (use_petsc) { MFEMInitializePetsc(NULL,NULL,petscrc_file,NULL); }
91
92 // 3. Read the (serial) mesh from the given mesh file on all processors. We
93 // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
94 // and volume meshes with the same code.
95 Mesh *mesh = new Mesh(mesh_file, 1, 1);
96 int dim = mesh->Dimension();
97 int sdim = mesh->SpaceDimension();
98
99 // 4. Refine the serial mesh on all processors to increase the resolution.
100 // Also project a NURBS mesh to a piecewise-quadratic curved mesh. Make
101 // sure that the mesh is non-conforming.
102 if (mesh->NURBSext)
103 {
104 mesh->UniformRefinement();
105 mesh->SetCurvature(2);
106 }
107 mesh->EnsureNCMesh();
108
109 // 5. Define a parallel mesh by partitioning the serial mesh.
110 // Once the parallel mesh is defined, the serial mesh can be deleted.
111 ParMesh pmesh(MPI_COMM_WORLD, *mesh);
112 delete mesh;
113
114 MFEM_VERIFY(pmesh.bdr_attributes.Size() > 0,
115 "Boundary attributes required in the mesh.");
116 Array<int> ess_bdr(pmesh.bdr_attributes.Max());
117 ess_bdr = 1;
118
119 // 6. Define a finite element space on the mesh. The polynomial order is
120 // one (linear) by default, but this can be changed on the command line.
121 H1_FECollection fec(order, dim);
122 ParFiniteElementSpace fespace(&pmesh, &fec);
123
124 // 7. As in Example 1p, we set up bilinear and linear forms corresponding to
125 // the Laplace problem -\Delta u = 1. We don't assemble the discrete
126 // problem yet, this will be done in the main loop.
127 ParBilinearForm a(&fespace);
128 ParLinearForm b(&fespace);
129
130 ConstantCoefficient one(1.0);
131
133 a.AddDomainIntegrator(integ);
134 b.AddDomainIntegrator(new DomainLFIntegrator(one));
135
136 // 8. The solution vector x and the associated finite element grid function
137 // will be maintained over the AMR iterations. We initialize it to zero.
138 ParGridFunction x(&fespace);
139 x = 0;
140
141 // 9. Connect to GLVis.
142 char vishost[] = "localhost";
143 int visport = 19916;
144
145 socketstream sout;
146 if (visualization)
147 {
148 sout.open(vishost, visport);
149 if (!sout)
150 {
151 if (myid == 0)
152 {
153 cout << "Unable to connect to GLVis server at "
154 << vishost << ':' << visport << endl;
155 cout << "GLVis visualization disabled.\n";
156 }
157 visualization = false;
158 }
159
160 sout.precision(8);
161 }
162
163 // 10. Set up an error estimator. Here we use the Zienkiewicz-Zhu estimator
164 // with L2 projection in the smoothing step to better handle hanging
165 // nodes and parallel partitioning. We need to supply a space for the
166 // discontinuous flux (L2) and a space for the smoothed flux (H(div) is
167 // used here).
168 L2_FECollection flux_fec(order, dim);
169 ParFiniteElementSpace flux_fes(&pmesh, &flux_fec, sdim);
170 RT_FECollection smooth_flux_fec(order-1, dim);
171 ParFiniteElementSpace smooth_flux_fes(&pmesh, &smooth_flux_fec);
172 // Another possible option for the smoothed flux space:
173 // H1_FECollection smooth_flux_fec(order, dim);
174 // ParFiniteElementSpace smooth_flux_fes(&pmesh, &smooth_flux_fec, dim);
175 L2ZienkiewiczZhuEstimator estimator(*integ, x, flux_fes, smooth_flux_fes);
176
177 // 11. A refiner selects and refines elements based on a refinement strategy.
178 // The strategy here is to refine elements with errors larger than a
179 // fraction of the maximum element error. Other strategies are possible.
180 // The refiner will call the given error estimator.
181 ThresholdRefiner refiner(estimator);
182 refiner.SetTotalErrorFraction(0.7);
183
184 // 12. The main AMR loop. In each iteration we solve the problem on the
185 // current mesh, visualize the solution, and refine the mesh.
186 for (int it = 0; ; it++)
187 {
188 HYPRE_BigInt global_dofs = fespace.GlobalTrueVSize();
189 if (myid == 0)
190 {
191 cout << "\nAMR iteration " << it << endl;
192 cout << "Number of unknowns: " << global_dofs << endl;
193 }
194
195 // 13. Assemble the stiffness matrix and the right-hand side. Note that
196 // MFEM doesn't care at this point that the mesh is nonconforming
197 // and parallel. The FE space is considered 'cut' along hanging
198 // edges/faces, and also across processor boundaries.
199 a.Assemble();
200 b.Assemble();
201
202 // 14. Create the parallel linear system: eliminate boundary conditions,
203 // constrain hanging nodes and nodes across processor boundaries.
204 // The system will be solved for true (unconstrained/unique) DOFs only.
205 Array<int> ess_tdof_list;
206 fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
207 real_t time;
208 const int copy_interior = 1;
209
210 if (use_petsc)
211 {
212 a.SetOperatorType(use_nonoverlapping ?
215 Vector pX,pB;
216 MPI_Barrier(MPI_COMM_WORLD);
217 time = -MPI_Wtime();
218 a.FormLinearSystem(ess_tdof_list, x, b, pA, pX, pB, copy_interior);
219 MPI_Barrier(MPI_COMM_WORLD);
220 time += MPI_Wtime();
221 if (myid == 0) { cout << "PETSc assembly timing : " << time << endl; }
222 }
223
224 a.Assemble();
225 b.Assemble();
226 a.SetOperatorType(Operator::Hypre_ParCSR);
228 Vector B, X;
229 MPI_Barrier(MPI_COMM_WORLD);
230 time = -MPI_Wtime();
231 a.FormLinearSystem(ess_tdof_list, x, b, A, X, B, copy_interior);
232 MPI_Barrier(MPI_COMM_WORLD);
233 time += MPI_Wtime();
234 if (myid == 0) { cout << "HYPRE assembly timing : " << time << endl; }
235
236 // 15. Define and apply a parallel PCG solver for AX=B with the BoomerAMG
237 // preconditioner from hypre.
238 HypreBoomerAMG amg;
239 amg.SetPrintLevel(0);
240 CGSolver pcg(A.GetComm());
241 pcg.SetPreconditioner(amg);
242 pcg.SetOperator(A);
243 pcg.SetRelTol(1e-6);
244 pcg.SetMaxIter(200);
245 pcg.SetPrintLevel(3); // print the first and the last iterations only
246 pcg.Mult(B, X);
247
248 // 16. Extract the parallel grid function corresponding to the finite element
249 // approximation X. This is the local solution on each processor.
250 a.RecoverFEMSolution(X, b, x);
251
252 // 17. Send the solution by socket to a GLVis server.
253 if (visualization)
254 {
255 sout << "parallel " << num_procs << " " << myid << "\n";
256 sout << "solution\n" << pmesh << x << flush;
257 }
258
259 if (global_dofs > max_dofs)
260 {
261 if (myid == 0)
262 {
263 cout << "Reached the maximum number of dofs. Stop." << endl;
264 }
265 // we need to call Update here to delete any internal PETSc object that
266 // have been created by the ParBilinearForm; otherwise, these objects
267 // will be destroyed at the end of the main scope, when PETSc has been
268 // already finalized.
269 a.Update();
270 b.Update();
271 break;
272 }
273
274 // 18. Call the refiner to modify the mesh. The refiner calls the error
275 // estimator to obtain element errors, then it selects elements to be
276 // refined and finally it modifies the mesh. The Stop() method can be
277 // used to determine if a stopping criterion was met.
278 refiner.Apply(pmesh);
279 if (refiner.Stop())
280 {
281 if (myid == 0)
282 {
283 cout << "Stopping criterion satisfied. Stop." << endl;
284 }
285 a.Update();
286 b.Update();
287 break;
288 }
289
290 // 19. Update the finite element space (recalculate the number of DOFs,
291 // etc.) and create a grid function update matrix. Apply the matrix
292 // to any GridFunctions over the space. In this case, the update
293 // matrix is an interpolation matrix so the updated GridFunction will
294 // still represent the same function as before refinement.
295 fespace.Update();
296 x.Update();
297
298 // 20. Load balance the mesh, and update the space and solution. Currently
299 // available only for nonconforming meshes.
300 if (pmesh.Nonconforming())
301 {
302 pmesh.Rebalance();
303
304 // Update the space and the GridFunction. This time the update matrix
305 // redistributes the GridFunction among the processors.
306 fespace.Update();
307 x.Update();
308 }
309
310 // 21. Inform also the bilinear and linear forms that the space has
311 // changed.
312 a.Update();
313 b.Update();
314 }
315
316 // We finalize PETSc
317 if (use_petsc) { MFEMFinalizePetsc(); }
318
319 return 0;
320}
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
Abstract base class BilinearFormIntegrator.
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
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
MPI_Comm GetComm() const
MPI communicator.
Definition hypre.hpp:578
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
The L2ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure where t...
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
bool Apply(Mesh &mesh)
Perform the mesh operation.
bool Stop() const
Check if STOP action is requested, e.g. stopping criterion is satisfied.
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
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition mesh.hpp:290
bool Nonconforming() const
Definition mesh.hpp:2229
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
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
Definition mesh.cpp:6211
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).
@ PETSC_MATIS
ID for class PetscParMatrix, MATIS format.
Definition operator.hpp:289
@ Hypre_ParCSR
ID for class HypreParMatrix.
Definition operator.hpp:287
@ PETSC_MATAIJ
ID for class PetscParMatrix, MATAIJ format.
Definition operator.hpp:288
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
void Update(bool want_transform=true) override
Class for parallel grid function.
Definition pgridfunc.hpp:33
void Update() override
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition pgridfunc.cpp:90
Class for parallel linear form.
Class for parallel meshes.
Definition pmesh.hpp:34
void Rebalance()
Definition pmesh.cpp:4009
Wrapper for PETSc's matrix class.
Definition petsc.hpp:319
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition fe_coll.hpp:386
Mesh refinement operator using an error threshold.
void SetTotalErrorFraction(real_t fraction)
Set the total fraction used in the computation of the threshold. The default value is 1/2.
Vector data type.
Definition vector.hpp:80
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
int dim
Definition ex24.cpp:53
int main()
HYPRE_Int HYPRE_BigInt
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
const int visport
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
Definition petsc.cpp:201
void MFEMFinalizePetsc()
Definition petsc.cpp:241
float real_t
Definition config.hpp:43
const char vishost[]