MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex21p.cpp
Go to the documentation of this file.
1 // MFEM Example 21
2 //
3 // Compile with: make ex21p
4 //
5 // Sample runs: mpirun -np 4 ex21p
6 // mpirun -np 4 ex21p -o 3
7 // mpirun -np 4 ex21p -m ../data/beam-quad.mesh
8 // mpirun -np 4 ex21p -m ../data/beam-quad.mesh -o 3
9 // mpirun -np 4 ex21p -m ../data/beam-tet.mesh
10 // mpirun -np 4 ex21p -m ../data/beam-tet.mesh -o 2
11 // mpirun -np 4 ex21p -m ../data/beam-hex.mesh
12 // mpirun -np 4 ex21p -m ../data/beam-hex.mesh -o 2
13 //
14 // Description: This is a version of Example 2p with a simple adaptive mesh
15 // refinement loop. The problem being solved is again the linear
16 // elasticity describing a multi-material cantilever beam.
17 // The problem is solved on a sequence of meshes which
18 // are locally refined in a conforming (triangles, tetrahedrons)
19 // or non-conforming (quadrilaterals, hexahedra) manner according
20 // to a simple ZZ error estimator.
21 //
22 // The example demonstrates MFEM's capability to work with both
23 // conforming and nonconforming refinements, in 2D and 3D, on
24 // linear and curved meshes. Interpolation of functions from
25 // coarse to fine meshes, as well as persistent GLVis
26 // visualization are also illustrated.
27 //
28 // We recommend viewing Examples 2p and 6p before viewing this
29 // example.
30 
31 #include "mfem.hpp"
32 #include <fstream>
33 #include <iostream>
34 
35 using namespace std;
36 using namespace mfem;
37 
38 int main(int argc, char *argv[])
39 {
40  // 0. Initialize MPI.
41  int num_procs, myid;
42  MPI_Init(&argc, &argv);
43  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
44  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
45 
46  // 1. Parse command-line options.
47  const char *mesh_file = "../data/beam-tri.mesh";
48  int serial_ref_levels = 0;
49  int order = 1;
50  bool static_cond = false;
51  bool visualization = 1;
52 
53  OptionsParser args(argc, argv);
54  args.AddOption(&mesh_file, "-m", "--mesh",
55  "Mesh file to use.");
56  args.AddOption(&serial_ref_levels, "-rs", "--refine-serial",
57  "Number of uniform serial refinements (before parallel"
58  " partitioning)");
59  args.AddOption(&order, "-o", "--order",
60  "Finite element order (polynomial degree).");
61  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
62  "--no-static-condensation", "Enable static condensation.");
63  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
64  "--no-visualization",
65  "Enable or disable GLVis visualization.");
66  args.Parse();
67  if (!args.Good())
68  {
69  if (myid == 0)
70  {
71  args.PrintUsage(cout);
72  }
73  MPI_Finalize();
74  return 1;
75  }
76  if (myid == 0)
77  {
78  args.PrintOptions(cout);
79  }
80 
81  // 2. Read the mesh from the given mesh file. We can handle triangular,
82  // quadrilateral, tetrahedral, and hexahedral meshes with the same code.
83  Mesh mesh(mesh_file, 1, 1);
84  int dim = mesh.Dimension();
85  MFEM_VERIFY(mesh.SpaceDimension() == dim, "invalid mesh");
86 
87  if (mesh.attributes.Max() < 2 || mesh.bdr_attributes.Max() < 2)
88  {
89  cerr << "\nInput mesh should have at least two materials and "
90  << "two boundary attributes! (See schematic in ex2.cpp)\n"
91  << endl;
92  MPI_Finalize();
93  return 3;
94  }
95 
96  // 3. Refine the mesh before parallel partitioning. Since a NURBS mesh can
97  // currently only be refined uniformly, we need to convert it to a
98  // piecewise-polynomial curved mesh. First we refine the NURBS mesh a bit
99  // more and then project the curvature to quadratic Nodes.
100  if (mesh.NURBSext && serial_ref_levels == 0)
101  {
102  serial_ref_levels = 2;
103  }
104  for (int i = 0; i < serial_ref_levels; i++)
105  {
106  mesh.UniformRefinement();
107  }
108  if (mesh.NURBSext)
109  {
110  mesh.SetCurvature(2);
111  }
112  mesh.EnsureNCMesh();
113 
114  ParMesh pmesh(MPI_COMM_WORLD, mesh);
115  mesh.Clear();
116 
117  // 4. Define a finite element space on the mesh. The polynomial order is
118  // one (linear) by default, but this can be changed on the command line.
119  H1_FECollection fec(order, dim);
120  ParFiniteElementSpace fespace(&pmesh, &fec, dim);
121 
122  // 5. As in Example 2, we set up the linear form b(.) which corresponds to
123  // the right-hand side of the FEM linear system. In this case, b_i equals
124  // the boundary integral of f*phi_i where f represents a "pull down"
125  // force on the Neumann part of the boundary and phi_i are the basis
126  // functions in the finite element fespace. The force is defined by the
127  // VectorArrayCoefficient object f, which is a vector of Coefficient
128  // objects. The fact that f is non-zero on boundary attribute 2 is
129  // indicated by the use of piece-wise constants coefficient for its last
130  // component. We don't assemble the discrete problem yet, this will be
131  // done in the main loop.
133  for (int i = 0; i < dim-1; i++)
134  {
135  f.Set(i, new ConstantCoefficient(0.0));
136  }
137  {
138  Vector pull_force(pmesh.bdr_attributes.Max());
139  pull_force = 0.0;
140  pull_force(1) = -1.0e-2;
141  f.Set(dim-1, new PWConstCoefficient(pull_force));
142  }
143 
144  ParLinearForm b(&fespace);
146 
147  // 6. Set up the bilinear form a(.,.) on the finite element space
148  // corresponding to the linear elasticity integrator with piece-wise
149  // constants coefficient lambda and mu.
150  Vector lambda(pmesh.attributes.Max());
151  lambda = 1.0;
152  lambda(0) = lambda(1)*50;
153  PWConstCoefficient lambda_func(lambda);
154  Vector mu(pmesh.attributes.Max());
155  mu = 1.0;
156  mu(0) = mu(1)*50;
157  PWConstCoefficient mu_func(mu);
158 
159  ParBilinearForm a(&fespace);
160  BilinearFormIntegrator *integ =
161  new ElasticityIntegrator(lambda_func,mu_func);
162  a.AddDomainIntegrator(integ);
163  if (static_cond) { a.EnableStaticCondensation(); }
164 
165  // 7. The solution vector x and the associated finite element grid function
166  // will be maintained over the AMR iterations. We initialize it to zero.
167  Vector zero_vec(dim);
168  zero_vec = 0.0;
169  VectorConstantCoefficient zero_vec_coeff(zero_vec);
170  ParGridFunction x(&fespace);
171  x = 0.0;
172 
173  // 8. Determine the list of true (i.e. conforming) essential boundary dofs.
174  // In this example, the boundary conditions are defined by marking only
175  // boundary attribute 1 from the mesh as essential and converting it to a
176  // list of true dofs. The conversion to true dofs will be done in the
177  // main loop.
178  Array<int> ess_bdr(pmesh.bdr_attributes.Max());
179  ess_bdr = 0;
180  ess_bdr[0] = 1;
181 
182  // 9. GLVis visualization.
183  char vishost[] = "localhost";
184  int visport = 19916;
185  socketstream sol_sock;
186 
187  // 10. Set up an error estimator. Here we use the Zienkiewicz-Zhu estimator
188  // that uses the ComputeElementFlux method of the ElasticityIntegrator to
189  // recover a smoothed flux (stress) that is subtracted from the element
190  // flux to get an error indicator. We need to supply the space for the
191  // smoothed flux: an (H1)^tdim (i.e., vector-valued) space is used here.
192  // Here, tdim represents the number of components for a symmetric (dim x
193  // dim) tensor.
194  const int tdim = dim*(dim+1)/2;
195  L2_FECollection flux_fec(order, dim);
196  ParFiniteElementSpace flux_fespace(&pmesh, &flux_fec, tdim);
197  ParFiniteElementSpace smooth_flux_fespace(&pmesh, &fec, tdim);
198  L2ZienkiewiczZhuEstimator estimator(*integ, x, flux_fespace,
199  smooth_flux_fespace);
200 
201  // 11. A refiner selects and refines elements based on a refinement strategy.
202  // The strategy here is to refine elements with errors larger than a
203  // fraction of the maximum element error. Other strategies are possible.
204  // The refiner will call the given error estimator.
205  ThresholdRefiner refiner(estimator);
206  refiner.SetTotalErrorFraction(0.7);
207 
208  // 12. The main AMR loop. In each iteration we solve the problem on the
209  // current mesh, visualize the solution, and refine the mesh.
210  const int max_dofs = 50000;
211  const int max_amr_itr = 20;
212  for (int it = 0; it <= max_amr_itr; it++)
213  {
214  HYPRE_Int global_dofs = fespace.GlobalTrueVSize();
215  if (myid == 0)
216  {
217  cout << "\nAMR iteration " << it << endl;
218  cout << "Number of unknowns: " << global_dofs << endl;
219  }
220 
221  // 13. Assemble the stiffness matrix and the right-hand side.
222  a.Assemble();
223  b.Assemble();
224 
225  // 14. Set Dirichlet boundary values in the GridFunction x.
226  // Determine the list of Dirichlet true DOFs in the linear system.
227  Array<int> ess_tdof_list;
228  x.ProjectBdrCoefficient(zero_vec_coeff, ess_bdr);
229  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
230 
231  // 15. Create the linear system: eliminate boundary conditions, constrain
232  // hanging nodes and possibly apply other transformations. The system
233  // will be solved for true (unconstrained) DOFs only.
234 
235  HypreParMatrix A;
236  Vector B, X;
237  const int copy_interior = 1;
238  a.FormLinearSystem(ess_tdof_list, x, b, A, X, B, copy_interior);
239 
240  // 16. Define and apply a parallel PCG solver for AX=B with the BoomerAMG
241  // preconditioner from hypre.
242  HypreBoomerAMG amg;
243  amg.SetPrintLevel(0);
244  // amg.SetSystemsOptions(dim); // optional
245  CGSolver pcg(A.GetComm());
246  pcg.SetPreconditioner(amg);
247  pcg.SetOperator(A);
248  pcg.SetRelTol(1e-6);
249  pcg.SetMaxIter(500);
250  pcg.SetPrintLevel(3); // print the first and the last iterations only
251  pcg.Mult(B, X);
252 
253  // 17. After solving the linear system, reconstruct the solution as a
254  // finite element GridFunction. Constrained nodes are interpolated
255  // from true DOFs (it may therefore happen that x.Size() >= X.Size()).
256  a.RecoverFEMSolution(X, b, x);
257 
258  // 18. Send solution by socket to the GLVis server.
259  if (visualization && it == 0)
260  {
261  sol_sock.open(vishost, visport);
262  sol_sock.precision(8);
263  }
264  if (visualization && sol_sock.good())
265  {
266  GridFunction nodes(&fespace), *nodes_p = &nodes;
267  pmesh.GetNodes(nodes);
268  nodes += x;
269  int own_nodes = 0;
270  pmesh.SwapNodes(nodes_p, own_nodes);
271  x.Neg(); // visualize the backward displacement
272  sol_sock << "parallel " << num_procs << ' ' << myid << '\n';
273  sol_sock << "solution\n" << pmesh << x << flush;
274  x.Neg();
275  pmesh.SwapNodes(nodes_p, own_nodes);
276  if (it == 0)
277  {
278  sol_sock << "keys '" << ((dim == 2) ? "Rjl" : "") << "m'" << endl;
279  }
280  sol_sock << "window_title 'AMR iteration: " << it << "'\n"
281  << "pause" << endl;
282  if (myid == 0)
283  {
284  cout << "Visualization paused. "
285  "Press <space> in the GLVis window to continue." << endl;
286  }
287  }
288 
289  if (global_dofs > max_dofs)
290  {
291  if (myid == 0)
292  {
293  cout << "Reached the maximum number of dofs. Stop." << endl;
294  }
295  break;
296  }
297 
298  // 19. Call the refiner to modify the mesh. The refiner calls the error
299  // estimator to obtain element errors, then it selects elements to be
300  // refined and finally it modifies the mesh. The Stop() method can be
301  // used to determine if a stopping criterion was met.
302  refiner.Apply(pmesh);
303  if (refiner.Stop())
304  {
305  if (myid == 0)
306  {
307  cout << "Stopping criterion satisfied. Stop." << endl;
308  }
309  break;
310  }
311 
312  // 20. Update the space to reflect the new state of the mesh. Also,
313  // interpolate the solution x so that it lies in the new space but
314  // represents the same function. This saves solver iterations later
315  // since we'll have a good initial guess of x in the next step.
316  // Internally, FiniteElementSpace::Update() calculates an
317  // interpolation matrix which is then used by GridFunction::Update().
318  fespace.Update();
319  x.Update();
320 
321  // 21. Load balance the mesh, and update the space and solution. Currently
322  // available only for nonconforming meshes.
323  if (pmesh.Nonconforming())
324  {
325  pmesh.Rebalance();
326 
327  // Update the space and the GridFunction. This time the update matrix
328  // redistributes the GridFunction among the processors.
329  fespace.Update();
330  x.Update();
331  }
332 
333  // 21. Inform also the bilinear and linear forms that the space has
334  // changed.
335  a.Update();
336  b.Update();
337  }
338 
339  {
340  ostringstream mref_name, mesh_name, sol_name;
341  mref_name << "ex21p_reference_mesh." << setfill('0') << setw(6) << myid;
342  mesh_name << "ex21p_deformed_mesh." << setfill('0') << setw(6) << myid;
343  sol_name << "ex21p_displacement." << setfill('0') << setw(6) << myid;
344 
345  ofstream mesh_ref_out(mref_name.str().c_str());
346  mesh_ref_out.precision(16);
347  pmesh.Print(mesh_ref_out);
348 
349  ofstream mesh_out(mesh_name.str().c_str());
350  mesh_out.precision(16);
351  GridFunction nodes(&fespace), *nodes_p = &nodes;
352  pmesh.GetNodes(nodes);
353  nodes += x;
354  int own_nodes = 0;
355  pmesh.SwapNodes(nodes_p, own_nodes);
356  pmesh.Print(mesh_out);
357  pmesh.SwapNodes(nodes_p, own_nodes);
358 
359  ofstream x_out(sol_name.str().c_str());
360  x_out.precision(16);
361  x.Save(x_out);
362  }
363 
364  MPI_Finalize();
365  return 0;
366 }
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
MPI_Comm GetComm() const
MPI communicator.
Definition: hypre.hpp:330
Vector coefficient defined by an array of scalar coefficients. Coefficients that are not set will eva...
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
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(...
Vector coefficient that is constant in space and time.
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Definition: mesh.cpp:6641
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
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:819
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[])
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
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 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
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.
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.
void Set(int i, Coefficient *c, bool own=true)
Sets coefficient in the vector.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4169
int Dimension() const
Definition: mesh.hpp:788
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:434
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 ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr)
Definition: pgridfunc.cpp:582
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
A piecewise constant coefficient with the constants keyed off the element attribute numbers...
Definition: coefficient.hpp:94
double mu
Definition: ex25.cpp:135
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;.
void Clear()
Clear the contents of the Mesh.
Definition: mesh.hpp:719
Vector data type.
Definition: vector.hpp:51
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:6603
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
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...
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:181
Class for parallel meshes.
Definition: pmesh.hpp:32
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:199
void EnableStaticCondensation()
Enable the use of static condensation. For details see the description for class StaticCondensation i...
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:221
double f(const Vector &p)
void Neg()
(*this) = -(*this)
Definition: vector.cpp:253
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:145