MFEM  v3.3
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex6p.cpp
Go to the documentation of this file.
1 // MFEM Example 6 - Parallel Version
2 //
3 // Compile with: make ex6p
4 //
5 // Sample runs: mpirun -np 4 ex6p -m ../data/square-disc.mesh -o 1
6 // mpirun -np 4 ex6p -m ../data/square-disc.mesh -o 2
7 // mpirun -np 4 ex6p -m ../data/square-disc-nurbs.mesh -o 2
8 // mpirun -np 4 ex6p -m ../data/star.mesh -o 3
9 // mpirun -np 4 ex6p -m ../data/escher.mesh -o 2
10 // mpirun -np 4 ex6p -m ../data/fichera.mesh -o 2
11 // mpirun -np 4 ex6p -m ../data/disc-nurbs.mesh -o 2
12 // mpirun -np 4 ex6p -m ../data/ball-nurbs.mesh
13 // mpirun -np 4 ex6p -m ../data/pipe-nurbs.mesh
14 // mpirun -np 4 ex6p -m ../data/star-surf.mesh -o 2
15 // mpirun -np 4 ex6p -m ../data/square-disc-surf.mesh -o 2
16 // mpirun -np 4 ex6p -m ../data/amr-quad.mesh
17 //
18 // Description: This is a version of Example 1 with a simple adaptive mesh
19 // refinement loop. The problem being solved is again the Laplace
20 // equation -Delta u = 1 with homogeneous Dirichlet boundary
21 // conditions. The problem is solved on a sequence of meshes which
22 // are locally refined in a conforming (triangles, tetrahedrons)
23 // or non-conforming (quadrilateral, hexahedrons) manner according
24 // to a simple ZZ error estimator.
25 //
26 // The example demonstrates MFEM's capability to work with both
27 // conforming and nonconforming refinements, in 2D and 3D, on
28 // linear, curved and surface meshes. Interpolation of functions
29 // from coarse to fine meshes, as well as persistent GLVis
30 // visualization are also illustrated.
31 //
32 // We recommend viewing Example 1 before viewing this example.
33 
34 #include "mfem.hpp"
35 #include <fstream>
36 #include <iostream>
37 
38 using namespace std;
39 using namespace mfem;
40 
41 int main(int argc, char *argv[])
42 {
43  // 1. Initialize MPI.
44  int num_procs, myid;
45  MPI_Init(&argc, &argv);
46  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
47  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
48 
49  // 2. Parse command-line options.
50  const char *mesh_file = "../data/star.mesh";
51  int order = 1;
52  bool visualization = true;
53 
54  OptionsParser args(argc, argv);
55  args.AddOption(&mesh_file, "-m", "--mesh",
56  "Mesh file to use.");
57  args.AddOption(&order, "-o", "--order",
58  "Finite element order (polynomial degree).");
59  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
60  "--no-visualization",
61  "Enable or disable GLVis visualization.");
62  args.Parse();
63  if (!args.Good())
64  {
65  if (myid == 0)
66  {
67  args.PrintUsage(cout);
68  }
69  MPI_Finalize();
70  return 1;
71  }
72  if (myid == 0)
73  {
74  args.PrintOptions(cout);
75  }
76 
77  // 3. Read the (serial) mesh from the given mesh file on all processors. We
78  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
79  // and volume meshes with the same code.
80  Mesh *mesh = new Mesh(mesh_file, 1, 1);
81  int dim = mesh->Dimension();
82  int sdim = mesh->SpaceDimension();
83 
84  // 4. Refine the serial mesh on all processors to increase the resolution.
85  // Also project a NURBS mesh to a piecewise-quadratic curved mesh. Make
86  // sure that the mesh is non-conforming.
87  if (mesh->NURBSext)
88  {
89  mesh->UniformRefinement();
90  mesh->SetCurvature(2);
91  }
92  mesh->EnsureNCMesh();
93 
94  // 5. Define a parallel mesh by partitioning the serial mesh.
95  // Once the parallel mesh is defined, the serial mesh can be deleted.
96  ParMesh pmesh(MPI_COMM_WORLD, *mesh);
97  delete mesh;
98 
99  MFEM_VERIFY(pmesh.bdr_attributes.Size() > 0,
100  "Boundary attributes required in the mesh.");
101  Array<int> ess_bdr(pmesh.bdr_attributes.Max());
102  ess_bdr = 1;
103 
104  // 6. Define a finite element space on the mesh. The polynomial order is
105  // one (linear) by default, but this can be changed on the command line.
106  H1_FECollection fec(order, dim);
107  ParFiniteElementSpace fespace(&pmesh, &fec);
108 
109  // 7. As in Example 1p, we set up bilinear and linear forms corresponding to
110  // the Laplace problem -\Delta u = 1. We don't assemble the discrete
111  // problem yet, this will be done in the main loop.
112  ParBilinearForm a(&fespace);
113  ParLinearForm b(&fespace);
114 
115  ConstantCoefficient one(1.0);
116 
118  a.AddDomainIntegrator(integ);
120 
121  // 8. The solution vector x and the associated finite element grid function
122  // will be maintained over the AMR iterations. We initialize it to zero.
123  ParGridFunction x(&fespace);
124  x = 0;
125 
126  // 9. Connect to GLVis.
127  char vishost[] = "localhost";
128  int visport = 19916;
129 
130  socketstream sout;
131  if (visualization)
132  {
133  sout.open(vishost, visport);
134  if (!sout)
135  {
136  if (myid == 0)
137  {
138  cout << "Unable to connect to GLVis server at "
139  << vishost << ':' << visport << endl;
140  cout << "GLVis visualization disabled.\n";
141  }
142  visualization = false;
143  }
144 
145  sout.precision(8);
146  }
147 
148  // 10. Set up an error estimator. Here we use the Zienkiewicz-Zhu estimator
149  // with L2 projection in the smoothing step to better handle hanging
150  // nodes and parallel partitioning. We need to supply a space for the
151  // discontinuous flux (L2) and a space for the smoothed flux (H(div) is
152  // used here).
153  L2_FECollection flux_fec(order, dim);
154  ParFiniteElementSpace flux_fes(&pmesh, &flux_fec, sdim);
155  RT_FECollection smooth_flux_fec(order-1, dim);
156  ParFiniteElementSpace smooth_flux_fes(&pmesh, &smooth_flux_fec);
157  // Another possible option for the smoothed flux space:
158  // H1_FECollection smooth_flux_fec(order, dim);
159  // ParFiniteElementSpace smooth_flux_fes(&pmesh, &smooth_flux_fec, dim);
160  L2ZienkiewiczZhuEstimator estimator(*integ, x, flux_fes, smooth_flux_fes);
161 
162  // 11. A refiner selects and refines elements based on a refinement strategy.
163  // The strategy here is to refine elements with errors larger than a
164  // fraction of the maximum element error. Other strategies are possible.
165  // The refiner will call the given error estimator.
166  ThresholdRefiner refiner(estimator);
167  refiner.SetTotalErrorFraction(0.7);
168 
169  // 12. The main AMR loop. In each iteration we solve the problem on the
170  // current mesh, visualize the solution, and refine the mesh.
171  const int max_dofs = 100000;
172  for (int it = 0; ; it++)
173  {
174  HYPRE_Int global_dofs = fespace.GlobalTrueVSize();
175  if (myid == 0)
176  {
177  cout << "\nAMR iteration " << it << endl;
178  cout << "Number of unknowns: " << global_dofs << endl;
179  }
180 
181  // 13. Assemble the stiffness matrix and the right-hand side. Note that
182  // MFEM doesn't care at this point that the mesh is nonconforming
183  // and parallel. The FE space is considered 'cut' along hanging
184  // edges/faces, and also across processor boundaries.
185  a.Assemble();
186  b.Assemble();
187 
188  // 14. Create the parallel linear system: eliminate boundary conditions,
189  // constrain hanging nodes and nodes across processor boundaries.
190  // The system will be solved for true (unconstrained/unique) DOFs only.
191  Array<int> ess_tdof_list;
192  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
193 
194  HypreParMatrix A;
195  Vector B, X;
196  const int copy_interior = 1;
197  a.FormLinearSystem(ess_tdof_list, x, b, A, X, B, copy_interior);
198 
199  // 15. Define and apply a parallel PCG solver for AX=B with the BoomerAMG
200  // preconditioner from hypre.
201  HypreBoomerAMG amg;
202  amg.SetPrintLevel(0);
203  CGSolver pcg(A.GetComm());
204  pcg.SetPreconditioner(amg);
205  pcg.SetOperator(A);
206  pcg.SetRelTol(1e-6);
207  pcg.SetMaxIter(200);
208  pcg.SetPrintLevel(3); // print the first and the last iterations only
209  pcg.Mult(B, X);
210 
211  // 16. Extract the parallel grid function corresponding to the finite element
212  // approximation X. This is the local solution on each processor.
213  a.RecoverFEMSolution(X, b, x);
214 
215  // 17. Send the solution by socket to a GLVis server.
216  if (visualization)
217  {
218  sout << "parallel " << num_procs << " " << myid << "\n";
219  sout << "solution\n" << pmesh << x << flush;
220  }
221 
222  if (global_dofs > max_dofs)
223  {
224  if (myid == 0)
225  {
226  cout << "Reached the maximum number of dofs. Stop." << endl;
227  }
228  break;
229  }
230 
231  // 18. Call the refiner to modify the mesh. The refiner calls the error
232  // estimator to obtain element errors, then it selects elements to be
233  // refined and finally it modifies the mesh. The Stop() method can be
234  // used to determine if a stopping criterion was met.
235  refiner.Apply(pmesh);
236  if (refiner.Stop())
237  {
238  if (myid == 0)
239  {
240  cout << "Stopping criterion satisfied. Stop." << endl;
241  }
242  break;
243  }
244 
245  // 19. Update the finite element space (recalculate the number of DOFs,
246  // etc.) and create a grid function update matrix. Apply the matrix
247  // to any GridFunctions over the space. In this case, the update
248  // matrix is an interpolation matrix so the updated GridFunction will
249  // still represent the same function as before refinement.
250  fespace.Update();
251  x.Update();
252 
253  // 20. Load balance the mesh, and update the space and solution. Currently
254  // available only for nonconforming meshes.
255  if (pmesh.Nonconforming())
256  {
257  pmesh.Rebalance();
258 
259  // Update the space and the GridFunction. This time the update matrix
260  // redistributes the GridFunction among the processors.
261  fespace.Update();
262  x.Update();
263  }
264 
265  // 21. Inform also the bilinear and linear forms that the space has
266  // changed.
267  a.Update();
268  b.Update();
269  }
270 
271  MPI_Finalize();
272  return 0;
273 }
int Size() const
Logical size of the array.
Definition: array.hpp:109
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:47
Conjugate gradient method.
Definition: solvers.hpp:111
MPI_Comm GetComm() const
MPI communicator.
Definition: hypre.hpp:295
Subclass constant coefficient.
Definition: coefficient.hpp:57
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
virtual void Update(bool want_transform=true)
Definition: pfespace.cpp:2198
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:42
bool Stop() const
Check if STOP action is requested, e.g. stopping criterion is satisfied.
Abstract parallel finite element space.
Definition: pfespace.hpp:28
int main(int argc, char *argv[])
The L2ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure where t...
Definition: estimators.hpp:183
void Update(ParFiniteElementSpace *pf=NULL)
Definition: plinearform.cpp:21
virtual void Update(FiniteElementSpace *nfes=NULL)
The BoomerAMG solver in hypre.
Definition: hypre.hpp:770
bool Apply(Mesh &mesh)
Perform the mesh operation.
void Rebalance()
Load balance the mesh. NC meshes only.
Definition: pmesh.cpp:2713
bool Nonconforming() const
Definition: mesh.hpp:968
Class for parallel linear form.
Definition: plinearform.hpp:26
int dim
Definition: ex3.cpp:47
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6684
void SetPrintLevel(int print_level)
Definition: hypre.hpp:811
Mesh refinement operator using an error threshold.
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3278
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:108
void Assemble(int skip_zeros=1)
Assemble the local matrix.
int Dimension() const
Definition: mesh.hpp:611
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:239
int SpaceDimension() const
Definition: mesh.hpp:612
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:22
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:142
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list)
Definition: pfespace.cpp:521
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)
Definition: optparser.hpp:74
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:144
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
Class for parallel bilinear form.
void EnsureNCMesh(bool triangles_nonconforming=false)
Definition: mesh.cpp:6350
int open(const char hostname[], int port)
Vector data type.
Definition: vector.hpp:36
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:92
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:146
Class for parallel grid function.
Definition: pgridfunc.hpp:31
void SetTotalErrorFraction(double fraction)
Set the total fraction used in the computation of the threshold. The default value is 1/2...
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:174
Class for parallel meshes.
Definition: pmesh.hpp:28
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:195
bool Good() const
Definition: optparser.hpp:120