MFEM  v4.2.0
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 // Device sample runs:
19 // mpirun -np 4 ex6p -pa -d cuda
20 // mpirun -np 4 ex6p -pa -d occa-cuda
21 // mpirun -np 4 ex6p -pa -d raja-omp
22 // mpirun -np 4 ex6p -pa -d ceed-cpu
23 // * mpirun -np 4 ex6p -pa -d ceed-cuda
24 // mpirun -np 4 ex6p -pa -d ceed-cuda:/gpu/cuda/shared
25 //
26 // Description: This is a version of Example 1 with a simple adaptive mesh
27 // refinement loop. The problem being solved is again the Laplace
28 // equation -Delta u = 1 with homogeneous Dirichlet boundary
29 // conditions. The problem is solved on a sequence of meshes which
30 // are locally refined in a conforming (triangles, tetrahedrons)
31 // or non-conforming (quadrilaterals, hexahedra) manner according
32 // to a simple ZZ error estimator.
33 //
34 // The example demonstrates MFEM's capability to work with both
35 // conforming and nonconforming refinements, in 2D and 3D, on
36 // linear, curved and surface meshes. Interpolation of functions
37 // from coarse to fine meshes, as well as persistent GLVis
38 // visualization are also illustrated.
39 //
40 // We recommend viewing Example 1 before viewing this example.
41 
42 #include "mfem.hpp"
43 #include <fstream>
44 #include <iostream>
45 
46 using namespace std;
47 using namespace mfem;
48 
49 int main(int argc, char *argv[])
50 {
51  // 1. Initialize MPI.
52  int num_procs, myid;
53  MPI_Init(&argc, &argv);
54  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
55  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
56 
57  // 2. Parse command-line options.
58  const char *mesh_file = "../data/star.mesh";
59  int order = 1;
60  bool pa = false;
61  const char *device_config = "cpu";
62  bool visualization = true;
63 
64  OptionsParser args(argc, argv);
65  args.AddOption(&mesh_file, "-m", "--mesh",
66  "Mesh file to use.");
67  args.AddOption(&order, "-o", "--order",
68  "Finite element order (polynomial degree).");
69  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
70  "--no-partial-assembly", "Enable Partial Assembly.");
71  args.AddOption(&device_config, "-d", "--device",
72  "Device configuration string, see Device::Configure().");
73  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
74  "--no-visualization",
75  "Enable or disable GLVis visualization.");
76  args.Parse();
77  if (!args.Good())
78  {
79  if (myid == 0)
80  {
81  args.PrintUsage(cout);
82  }
83  MPI_Finalize();
84  return 1;
85  }
86  if (myid == 0)
87  {
88  args.PrintOptions(cout);
89  }
90 
91  // 3. Enable hardware devices such as GPUs, and programming models such as
92  // CUDA, OCCA, RAJA and OpenMP based on command line options.
93  Device device(device_config);
94  if (myid == 0) { device.Print(); }
95 
96  // 4. Read the (serial) mesh from the given mesh file on all processors. We
97  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
98  // and volume meshes with the same code.
99  Mesh *mesh = new Mesh(mesh_file, 1, 1);
100  int dim = mesh->Dimension();
101  int sdim = mesh->SpaceDimension();
102 
103  // 5. Refine the serial mesh on all processors to increase the resolution.
104  // Also project a NURBS mesh to a piecewise-quadratic curved mesh. Make
105  // sure that the mesh is non-conforming.
106  if (mesh->NURBSext)
107  {
108  mesh->UniformRefinement();
109  mesh->SetCurvature(2);
110  }
111  mesh->EnsureNCMesh();
112 
113  // 6. Define a parallel mesh by partitioning the serial mesh.
114  // Once the parallel mesh is defined, the serial mesh can be deleted.
115  ParMesh pmesh(MPI_COMM_WORLD, *mesh);
116  delete mesh;
117 
118  MFEM_VERIFY(pmesh.bdr_attributes.Size() > 0,
119  "Boundary attributes required in the mesh.");
120  Array<int> ess_bdr(pmesh.bdr_attributes.Max());
121  ess_bdr = 1;
122 
123  // 7. Define a finite element space on the mesh. The polynomial order is
124  // one (linear) by default, but this can be changed on the command line.
125  H1_FECollection fec(order, dim);
126  ParFiniteElementSpace fespace(&pmesh, &fec);
127 
128  // 8. As in Example 1p, we set up bilinear and linear forms corresponding to
129  // the Laplace problem -\Delta u = 1. We don't assemble the discrete
130  // problem yet, this will be done in the main loop.
131  ParBilinearForm a(&fespace);
132  if (pa)
133  {
134  a.SetAssemblyLevel(AssemblyLevel::PARTIAL);
135  a.SetDiagonalPolicy(Operator::DIAG_ONE);
136  }
137  ParLinearForm b(&fespace);
138 
139  ConstantCoefficient one(1.0);
140 
142  a.AddDomainIntegrator(integ);
144 
145  // 9. The solution vector x and the associated finite element grid function
146  // will be maintained over the AMR iterations. We initialize it to zero.
147  ParGridFunction x(&fespace);
148  x = 0;
149 
150  // 10. Connect to GLVis.
151  char vishost[] = "localhost";
152  int visport = 19916;
153 
154  socketstream sout;
155  if (visualization)
156  {
157  sout.open(vishost, visport);
158  if (!sout)
159  {
160  if (myid == 0)
161  {
162  cout << "Unable to connect to GLVis server at "
163  << vishost << ':' << visport << endl;
164  cout << "GLVis visualization disabled.\n";
165  }
166  visualization = false;
167  }
168 
169  sout.precision(8);
170  }
171 
172  // 11. Set up an error estimator. Here we use the Zienkiewicz-Zhu estimator
173  // with L2 projection in the smoothing step to better handle hanging
174  // nodes and parallel partitioning. We need to supply a space for the
175  // discontinuous flux (L2) and a space for the smoothed flux (H(div) is
176  // used here).
177  L2_FECollection flux_fec(order, dim);
178  ParFiniteElementSpace flux_fes(&pmesh, &flux_fec, sdim);
179  RT_FECollection smooth_flux_fec(order-1, dim);
180  ParFiniteElementSpace smooth_flux_fes(&pmesh, &smooth_flux_fec);
181  // Another possible option for the smoothed flux space:
182  // H1_FECollection smooth_flux_fec(order, dim);
183  // ParFiniteElementSpace smooth_flux_fes(&pmesh, &smooth_flux_fec, dim);
184  L2ZienkiewiczZhuEstimator estimator(*integ, x, flux_fes, smooth_flux_fes);
185 
186  // 12. A refiner selects and refines elements based on a refinement strategy.
187  // The strategy here is to refine elements with errors larger than a
188  // fraction of the maximum element error. Other strategies are possible.
189  // The refiner will call the given error estimator.
190  ThresholdRefiner refiner(estimator);
191  refiner.SetTotalErrorFraction(0.7);
192 
193  // 13. The main AMR loop. In each iteration we solve the problem on the
194  // current mesh, visualize the solution, and refine the mesh.
195  const int max_dofs = 100000;
196  for (int it = 0; ; it++)
197  {
198  HYPRE_Int global_dofs = fespace.GlobalTrueVSize();
199  if (myid == 0)
200  {
201  cout << "\nAMR iteration " << it << endl;
202  cout << "Number of unknowns: " << global_dofs << endl;
203  }
204 
205  // 14. Assemble the right-hand side and determine the list of true
206  // (i.e. parallel conforming) essential boundary dofs.
207  Array<int> ess_tdof_list;
208  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
209  b.Assemble();
210 
211  // 15. Assemble the stiffness matrix. Note that MFEM doesn't care at this
212  // point that the mesh is nonconforming and parallel. The FE space is
213  // considered 'cut' along hanging edges/faces, and also across
214  // processor boundaries.
215  a.Assemble();
216 
217  // 16. Create the parallel linear system: eliminate boundary conditions.
218  // The system will be solved for true (unconstrained/unique) DOFs only.
219  OperatorPtr A;
220  Vector B, X;
221 
222  const int copy_interior = 1;
223  a.FormLinearSystem(ess_tdof_list, x, b, A, X, B, copy_interior);
224 
225  // 17. Solve the linear system A X = B.
226  // * With full assembly, use the BoomerAMG preconditioner from hypre.
227  // * With partial assembly, use a diagonal preconditioner.
228  Solver *M = NULL;
229  if (pa)
230  {
231  M = new OperatorJacobiSmoother(a, ess_tdof_list);
232  }
233  else
234  {
235  HypreBoomerAMG *amg = new HypreBoomerAMG;
236  amg->SetPrintLevel(0);
237  M = amg;
238  }
239  CGSolver cg(MPI_COMM_WORLD);
240  cg.SetRelTol(1e-6);
241  cg.SetMaxIter(2000);
242  cg.SetPrintLevel(3); // print the first and the last iterations only
243  cg.SetPreconditioner(*M);
244  cg.SetOperator(*A);
245  cg.Mult(B, X);
246  delete M;
247 
248  // 18. Switch back to the host and extract the parallel grid function
249  // corresponding to the finite element approximation X. This is the
250  // local solution on each processor.
251  a.RecoverFEMSolution(X, b, x);
252 
253  // 19. Send the solution by socket to a GLVis server.
254  if (visualization)
255  {
256  sout << "parallel " << num_procs << " " << myid << "\n";
257  sout << "solution\n" << pmesh << x << flush;
258  }
259 
260  if (global_dofs > max_dofs)
261  {
262  if (myid == 0)
263  {
264  cout << "Reached the maximum number of dofs. Stop." << endl;
265  }
266  break;
267  }
268 
269  // 20. Call the refiner to modify the mesh. The refiner calls the error
270  // estimator to obtain element errors, then it selects elements to be
271  // refined and finally it modifies the mesh. The Stop() method can be
272  // used to determine if a stopping criterion was met.
273  refiner.Apply(pmesh);
274  if (refiner.Stop())
275  {
276  if (myid == 0)
277  {
278  cout << "Stopping criterion satisfied. Stop." << endl;
279  }
280  break;
281  }
282 
283  // 21. Update the finite element space (recalculate the number of DOFs,
284  // etc.) and create a grid function update matrix. Apply the matrix
285  // to any GridFunctions over the space. In this case, the update
286  // matrix is an interpolation matrix so the updated GridFunction will
287  // still represent the same function as before refinement.
288  fespace.Update();
289  x.Update();
290 
291  // 22. Load balance the mesh, and update the space and solution. Currently
292  // available only for nonconforming meshes.
293  if (pmesh.Nonconforming())
294  {
295  pmesh.Rebalance();
296 
297  // Update the space and the GridFunction. This time the update matrix
298  // redistributes the GridFunction among the processors.
299  fespace.Update();
300  x.Update();
301  }
302 
303  // 23. Inform also the bilinear and linear forms that the space has
304  // changed.
305  a.Update();
306  b.Update();
307  }
308 
309  MPI_Finalize();
310  return 0;
311 }
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:775
Conjugate gradient method.
Definition: solvers.hpp:258
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:535
virtual void Update(bool want_transform=true)
Definition: pfespace.cpp:2875
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
bool Stop() const
Check if STOP action is requested, e.g. stopping criterion is satisfied.
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:261
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
Abstract parallel finite element space.
Definition: pfespace.hpp:28
int main(int argc, char *argv[])
Definition: ex1.cpp:66
The L2ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure where t...
Definition: estimators.hpp:210
void Update(ParFiniteElementSpace *pf=NULL)
Update the object according to the given new FE space *pf.
Definition: plinearform.cpp:21
void SetDiagonalPolicy(DiagonalPolicy policy)
Sets diagonal policy used upon construction of the linear system.
virtual void Update(FiniteElementSpace *nfes=NULL)
Update the FiniteElementSpace and delete all data associated with the old one.
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1079
bool Apply(Mesh &mesh)
Perform the mesh operation.
void Rebalance()
Definition: pmesh.cpp:3379
bool Nonconforming() const
Definition: mesh.hpp:1214
Class for parallel linear form.
Definition: plinearform.hpp:26
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:70
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
constexpr char vishost[]
void EnsureNCMesh(bool simplices_nonconforming=false)
Definition: mesh.cpp:8045
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: pgridfunc.cpp:81
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:128
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:8382
void SetPrintLevel(int print_level)
Definition: hypre.hpp:1119
constexpr int visport
Mesh refinement operator using an error threshold.
void SetMaxIter(int max_it)
Definition: solvers.hpp:98
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:4165
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:251
void Assemble(int skip_zeros=1)
Assemble the local matrix.
int Dimension() const
Definition: mesh.hpp:788
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:434
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:272
int SpaceDimension() const
Definition: mesh.hpp:789
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:31
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:201
void SetRelTol(double rtol)
Definition: solvers.hpp:96
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 &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:203
int dim
Definition: ex24.cpp:53
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:304
Class for parallel bilinear form.
int open(const char hostname[], int port)
Open the socket stream on &#39;port&#39; at &#39;hostname&#39;.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:272
Vector data type.
Definition: vector.hpp:51
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:91
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:159
Base class for solvers.
Definition: operator.hpp:634
Class for parallel grid function.
Definition: pgridfunc.hpp:32
void SetTotalErrorFraction(double fraction)
Set the total fraction used in the computation of the threshold. The default value is 1/2...
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:118
Class for parallel meshes.
Definition: pmesh.hpp:32
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:221
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:145