MFEM  v4.6.0
Finite element discretization library
ex28p.cpp
Go to the documentation of this file.
1 // MFEM Example 28 - Parallel Version
2 //
3 // Compile with: make ex28p
4 //
5 // Sample runs: ex28p
6 // ex28p --visit-datafiles
7 // ex28p --order 4
8 // ex28p --penalty 1e+5
9 //
10 // mpirun -np 4 ex28p
11 // mpirun -np 4 ex28p --penalty 1e+5
12 //
13 // Description: Demonstrates a sliding boundary condition in an elasticity
14 // problem. A trapezoid, roughly as pictured below, is pushed
15 // from the right into a rigid notch. Normal displacement is
16 // restricted, but tangential movement is allowed, so the
17 // trapezoid compresses into the notch.
18 //
19 // /-------+
20 // normal constrained --->/ | <--- boundary force (2)
21 // boundary (4) /---------+
22 // ^
23 // |
24 // normal constrained boundary (1)
25 //
26 // This example demonstrates the use of the ConstrainedSolver
27 // framework.
28 //
29 // We recommend viewing Example 2 before viewing this example.
30 
31 #include "mfem.hpp"
32 #include <fstream>
33 #include <iostream>
34 
35 using namespace std;
36 using namespace mfem;
37 
38 // Return a mesh with a single element with vertices (0, 0), (1, 0), (1, 1),
39 // (offset, 1) to demonstrate boundary conditions on a surface that is not
40 // axis-aligned.
41 Mesh * build_trapezoid_mesh(double offset)
42 {
43  MFEM_VERIFY(offset < 0.9, "offset is too large!");
44 
45  const int dimension = 2;
46  const int nvt = 4; // vertices
47  const int nbe = 4; // num boundary elements
48  Mesh * mesh = new Mesh(dimension, nvt, 1, nbe);
49 
50  // vertices
51  double vc[dimension];
52  vc[0] = 0.0; vc[1] = 0.0;
53  mesh->AddVertex(vc);
54  vc[0] = 1.0; vc[1] = 0.0;
55  mesh->AddVertex(vc);
56  vc[0] = offset; vc[1] = 1.0;
57  mesh->AddVertex(vc);
58  vc[0] = 1.0; vc[1] = 1.0;
59  mesh->AddVertex(vc);
60 
61  // element
62  Array<int> vert(4);
63  vert[0] = 0; vert[1] = 1; vert[2] = 3; vert[3] = 2;
64  mesh->AddQuad(vert, 1);
65 
66  // boundary
67  Array<int> sv(2);
68  sv[0] = 0; sv[1] = 1;
69  mesh->AddBdrSegment(sv, 1);
70  sv[0] = 1; sv[1] = 3;
71  mesh->AddBdrSegment(sv, 2);
72  sv[0] = 2; sv[1] = 3;
73  mesh->AddBdrSegment(sv, 3);
74  sv[0] = 0; sv[1] = 2;
75  mesh->AddBdrSegment(sv, 4);
76 
77  mesh->FinalizeQuadMesh(1, 0, true);
78 
79  return mesh;
80 }
81 
82 int main(int argc, char *argv[])
83 {
84 #ifdef HYPRE_USING_GPU
85  cout << "\nAs of mfem-4.3 and hypre-2.22.0 (July 2021) this example\n"
86  << "is NOT supported with the GPU version of hypre.\n\n";
87  return 242;
88 #endif
89 
90  // 1. Initialize MPI and HYPRE.
91  Mpi::Init(argc, argv);
92  int num_procs = Mpi::WorldSize();
93  int myid = Mpi::WorldRank();
94  Hypre::Init();
95 
96  // 2. Parse command-line options.
97  int order = 1;
98  bool visualization = 1;
99  bool reorder_space = false;
100  double offset = 0.3;
101  bool visit = false;
102  double penalty = 0.0;
103 
104  OptionsParser args(argc, argv);
105  args.AddOption(&order, "-o", "--order",
106  "Finite element order (polynomial degree).");
107  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
108  "--no-visualization",
109  "Enable or disable GLVis visualization.");
110  args.AddOption(&reorder_space, "-nodes", "--by-nodes", "-vdim", "--by-vdim",
111  "Use byNODES ordering of vector space instead of byVDIM");
112  args.AddOption(&offset, "--offset", "--offset",
113  "How much to offset the trapezoid.");
114  args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
115  "--no-visit-datafiles",
116  "Save data files for VisIt (visit.llnl.gov) visualization.");
117  args.AddOption(&penalty, "-p", "--penalty",
118  "Penalty parameter; 0 means use elimination solver.");
119  args.Parse();
120  if (!args.Good())
121  {
122  if (myid == 0)
123  {
124  args.PrintUsage(cout);
125  }
126  return 1;
127  }
128  if (myid == 0)
129  {
130  args.PrintOptions(cout);
131  }
132 
133  // 3. Build a trapezoidal mesh with a single quadrilateral element, where
134  // 'offset' determines how far off it is from a rectangle.
135  Mesh *mesh = build_trapezoid_mesh(offset);
136  int dim = mesh->Dimension();
137 
138  // 4. 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 =
144  (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
145  for (int l = 0; l < ref_levels; l++)
146  {
147  mesh->UniformRefinement();
148  }
149  }
150 
151  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
152  // this mesh further in parallel to increase the resolution. Once the
153  // parallel mesh is defined, the serial mesh can be deleted.
154  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
155  delete mesh;
156  {
157  int par_ref_levels = 1;
158  for (int l = 0; l < par_ref_levels; l++)
159  {
160  pmesh->UniformRefinement();
161  }
162  }
163 
164  // 6. Define a parallel finite element space on the parallel mesh. Here we
165  // use vector finite elements, i.e. dim copies of a scalar finite element
166  // space. We use the ordering by vector dimension (the last argument of
167  // the FiniteElementSpace constructor) which is expected in the systems
168  // version of BoomerAMG preconditioner. For NURBS meshes, we use the
169  // (degree elevated) NURBS space associated with the mesh nodes.
171  ParFiniteElementSpace *fespace;
172  const bool use_nodal_fespace = pmesh->NURBSext;
173  if (use_nodal_fespace)
174  {
175  fec = NULL;
176  fespace = (ParFiniteElementSpace *)pmesh->GetNodes()->FESpace();
177  }
178  else
179  {
180  fec = new H1_FECollection(order, dim);
181  if (reorder_space)
182  {
183  fespace = new ParFiniteElementSpace(pmesh, fec, dim, Ordering::byNODES);
184  }
185  else
186  {
187  fespace = new ParFiniteElementSpace(pmesh, fec, dim, Ordering::byVDIM);
188  }
189  }
190  HYPRE_BigInt size = fespace->GlobalTrueVSize();
191  if (myid == 0)
192  {
193  cout << "Number of finite element unknowns: " << size << endl
194  << "Assembling matrix and r.h.s... " << flush;
195  }
196 
197  // 7. Determine the list of true (i.e. parallel conforming) essential
198  // boundary dofs. In this example, there are no essential boundary
199  // conditions in the usual sense, but we leave the machinery here for
200  // users to modify if they wish.
201  Array<int> ess_tdof_list, ess_bdr(pmesh->bdr_attributes.Max());
202  ess_bdr = 0;
203  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
204 
205  // 8. Set up the parallel linear form b(.) which corresponds to the
206  // right-hand side of the FEM linear system. In this case, b_i equals the
207  // boundary integral of f*phi_i where f represents a "pull down" force on
208  // the Neumann part of the boundary and phi_i are the basis functions in
209  // the finite element fespace. The force is defined by the object f, which
210  // is a vector of Coefficient objects. The fact that f is non-zero on
211  // boundary attribute 2 is indicated by the use of piece-wise constants
212  // coefficient for its last component.
214  for (int i = 0; i < dim-1; i++)
215  {
216  f.Set(i, new ConstantCoefficient(0.0));
217  }
218 
219  // 9. Put a leftward force on the right side of the trapezoid
220  {
221  Vector push_force(pmesh->bdr_attributes.Max());
222  push_force = 0.0;
223  push_force(1) = -5.0e-2; // index 1 attribute 2
224  f.Set(0, new PWConstCoefficient(push_force));
225  }
226 
227  ParLinearForm *b = new ParLinearForm(fespace);
228  b->AddBoundaryIntegrator(new VectorBoundaryLFIntegrator(f));
229  b->Assemble();
230 
231  // 10. Define the solution vector x as a parallel finite element grid
232  // function corresponding to fespace. Initialize x with initial guess of
233  // zero, which satisfies the boundary conditions.
234  ParGridFunction x(fespace);
235  x = 0.0;
236 
237  // 11. Set up the parallel bilinear form a(.,.) on the finite element space
238  // corresponding to the linear elasticity integrator with piece-wise
239  // constants coefficient lambda and mu. We use constant coefficients,
240  // but see ex2 for how to set up piecewise constant coefficients based
241  // on attribute.
242  Vector lambda(pmesh->attributes.Max());
243  lambda = 1.0;
244  PWConstCoefficient lambda_func(lambda);
245  Vector mu(pmesh->attributes.Max());
246  mu = 1.0;
247  PWConstCoefficient mu_func(mu);
248  ParBilinearForm *a = new ParBilinearForm(fespace);
249  a->AddDomainIntegrator(new ElasticityIntegrator(lambda_func, mu_func));
250 
251  // 12. Assemble the parallel bilinear form and the corresponding linear
252  // system, applying any necessary transformations such as: parallel
253  // assembly, eliminating boundary conditions, applying conforming
254  // constraints for non-conforming AMR, etc.
255  a->Assemble();
256 
257  HypreParMatrix A;
258  Vector B, X;
259  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
260  if (myid == 0)
261  {
262  cout << "done." << endl;
263  cout << "Size of linear system: " << A.GetGlobalNumRows() << endl;
264  }
265 
266  // 13. Set up constraint matrix to constrain normal displacement (but
267  // allow tangential displacement) on specified boundaries.
268  Array<int> constraint_atts(2);
269  constraint_atts[0] = 1; // attribute 1 bottom
270  constraint_atts[1] = 4; // attribute 4 left side
271  Array<int> constraint_rowstarts;
272  SparseMatrix* local_constraints =
273  ParBuildNormalConstraints(*fespace, constraint_atts,
274  constraint_rowstarts);
275 
276  // 14. Define and apply a parallel PCG solver for the constrained system
277  // where the normal boundary constraints have been separately eliminated
278  // from the system.
279  ConstrainedSolver * solver;
280  if (penalty == 0.0)
281  {
282  solver = new EliminationCGSolver(A, *local_constraints,
283  constraint_rowstarts, dim,
284  reorder_space);
285  }
286  else
287  {
288  solver = new PenaltyPCGSolver(A, *local_constraints, penalty,
289  dim, reorder_space);
290  }
291 
292  solver->SetRelTol(1e-8);
293  solver->SetMaxIter(500);
294  solver->SetPrintLevel(1);
295  solver->Mult(B, X);
296 
297  // 15. Recover the parallel grid function corresponding to X. This is the
298  // local finite element solution on each processor.
299  a->RecoverFEMSolution(X, *b, x);
300 
301  // 16. For non-NURBS meshes, make the mesh curved based on the finite element
302  // space. This means that we define the mesh elements through a fespace
303  // based transformation of the reference element. This allows us to save
304  // the displaced mesh as a curved mesh when using high-order finite
305  // element displacement field. We assume that the initial mesh (read from
306  // the file) is not higher order curved mesh compared to the chosen FE
307  // space.
308  if (!use_nodal_fespace)
309  {
310  pmesh->SetNodalFESpace(fespace);
311  }
312 
313  GridFunction *nodes = pmesh->GetNodes();
314  *nodes += x;
315 
316  // 17. Save the refined mesh and the solution in VisIt format.
317  if (visit)
318  {
319  VisItDataCollection visit_dc(MPI_COMM_WORLD, "ex28p", pmesh);
320  visit_dc.SetLevelsOfDetail(4);
321  visit_dc.RegisterField("displacement", &x);
322  visit_dc.Save();
323  }
324 
325  // 18. Save in parallel the displaced mesh and the inverted solution (which
326  // gives the backward displacements to the original grid). This output
327  // can be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
328  {
329  x *= -1; // sign convention for GLVis displacements
330 
331  ostringstream mesh_name, sol_name;
332  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
333  sol_name << "sol." << setfill('0') << setw(6) << myid;
334 
335  ofstream mesh_ofs(mesh_name.str().c_str());
336  mesh_ofs.precision(8);
337  pmesh->Print(mesh_ofs);
338 
339  ofstream sol_ofs(sol_name.str().c_str());
340  sol_ofs.precision(8);
341  x.Save(sol_ofs);
342  }
343 
344  // 19. Send the above data by socket to a GLVis server. Use the "n" and "b"
345  // keys in GLVis to visualize the displacements.
346  if (visualization)
347  {
348  char vishost[] = "localhost";
349  int visport = 19916;
350  socketstream sol_sock(vishost, visport);
351  sol_sock << "parallel " << num_procs << " " << myid << "\n";
352  sol_sock.precision(8);
353  sol_sock << "solution\n" << *pmesh << x << flush;
354  }
355 
356  // 20. Free the used memory.
357  delete local_constraints;
358  delete solver;
359  delete a;
360  delete b;
361  if (fec)
362  {
363  delete fespace;
364  delete fec;
365  }
366  delete pmesh;
367 
368  // HYPRE_Finalize();
369 
370  return 0;
371 }
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:1031
int visport
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
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1694
A coefficient that is constant across space and time.
Definition: coefficient.hpp:84
SparseMatrix * ParBuildNormalConstraints(ParFiniteElementSpace &fespace, Array< int > &constrained_att, Array< int > &constraint_rowstarts)
Parallel wrapper for BuildNormalConstraints.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:331
int Dimension() const
Dimension of the reference space used within the elements.
Definition: mesh.hpp:1020
int main(int argc, char *argv[])
Definition: ex28p.cpp:82
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:462
An abstract class to solve the constrained system subject to the constraint .
Definition: constraints.hpp:53
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition: array.cpp:68
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:159
Abstract parallel finite element space.
Definition: pfespace.hpp:28
STL namespace.
virtual void Mult(const Vector &f, Vector &x) const override
Solve for given .
int AddVertex(double x, double y=0.0, double z=0.0)
Definition: mesh.cpp:1626
Class for parallel linear form.
Definition: plinearform.hpp:26
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
char vishost[]
Data type sparse matrix.
Definition: sparsemat.hpp:50
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:10232
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
Definition: solvers.hpp:201
int AddBdrSegment(int v1, int v2, int attr=1)
Definition: mesh.cpp:1843
void SetLevelsOfDetail(int levels_of_detail)
Set VisIt parameter: default levels of detail for the MultiresControl.
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:281
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:275
void SetRelTol(double rtol)
Definition: solvers.hpp:199
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
virtual void Save() override
Save the collection and a VisIt root file.
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
HYPRE_Int HYPRE_BigInt
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:909
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:1086
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:277
A piecewise constant coefficient with the constants keyed off the element attribute numbers...
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a quadrilateral Mesh.
Definition: mesh.cpp:1988
double mu
Definition: ex25.cpp:140
int dim
Definition: ex24.cpp:53
Mesh * build_trapezoid_mesh(double offset)
Definition: ex28p.cpp:41
constexpr int dimension
This example only works in 3D. Kernels for 2D are not implemented.
Definition: hooke.cpp:45
Class for parallel bilinear form.
Vector data type.
Definition: vector.hpp:58
void SetNodalFESpace(FiniteElementSpace *nfes) override
Definition: pmesh.cpp:2076
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:259
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:8302
void Print(std::ostream &out=mfem::out) const override
Definition: pmesh.cpp:4825
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:343
Class for parallel meshes.
Definition: pmesh.hpp:32
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:635
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:273
double f(const Vector &p)
virtual void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.