MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
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_CUDA
85  cout << "\nAs of mfem-4.3 and hypre-2.22.0 (July 2021) this example\n"
86  << "is NOT supported with the CUDA version of hypre.\n\n";
87  return 255;
88 #endif
89 
90  // 1. Initialize MPI.
91  int num_procs, myid;
92  MPI_Init(&argc, &argv);
93  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
94  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
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  MPI_Finalize();
127  return 1;
128  }
129  if (myid == 0)
130  {
131  args.PrintOptions(cout);
132  }
133 
134  // 3. Build a trapezoidal mesh with a single quadrilateral element, where
135  // 'offset' determines how far off it is from a rectangle.
136  Mesh *mesh = build_trapezoid_mesh(offset);
137  int dim = mesh->Dimension();
138 
139  // 4. Refine the serial mesh on all processors to increase the resolution. In
140  // this example we do 'ref_levels' of uniform refinement. We choose
141  // 'ref_levels' to be the largest number that gives a final mesh with no
142  // more than 1,000 elements.
143  {
144  int ref_levels =
145  (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
146  for (int l = 0; l < ref_levels; l++)
147  {
148  mesh->UniformRefinement();
149  }
150  }
151 
152  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
153  // this mesh further in parallel to increase the resolution. Once the
154  // parallel mesh is defined, the serial mesh can be deleted.
155  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
156  delete mesh;
157  {
158  int par_ref_levels = 1;
159  for (int l = 0; l < par_ref_levels; l++)
160  {
161  pmesh->UniformRefinement();
162  }
163  }
164 
165  // 6. Define a parallel finite element space on the parallel mesh. Here we
166  // use vector finite elements, i.e. dim copies of a scalar finite element
167  // space. We use the ordering by vector dimension (the last argument of
168  // the FiniteElementSpace constructor) which is expected in the systems
169  // version of BoomerAMG preconditioner. For NURBS meshes, we use the
170  // (degree elevated) NURBS space associated with the mesh nodes.
172  ParFiniteElementSpace *fespace;
173  const bool use_nodal_fespace = pmesh->NURBSext;
174  if (use_nodal_fespace)
175  {
176  fec = NULL;
177  fespace = (ParFiniteElementSpace *)pmesh->GetNodes()->FESpace();
178  }
179  else
180  {
181  fec = new H1_FECollection(order, dim);
182  if (reorder_space)
183  {
184  fespace = new ParFiniteElementSpace(pmesh, fec, dim, Ordering::byNODES);
185  }
186  else
187  {
188  fespace = new ParFiniteElementSpace(pmesh, fec, dim, Ordering::byVDIM);
189  }
190  }
191  HYPRE_BigInt size = fespace->GlobalTrueVSize();
192  if (myid == 0)
193  {
194  cout << "Number of finite element unknowns: " << size << endl
195  << "Assembling matrix and r.h.s... " << flush;
196  }
197 
198  // 7. Determine the list of true (i.e. parallel conforming) essential
199  // boundary dofs. In this example, there are no essential boundary
200  // conditions in the usual sense, but we leave the machinery here for
201  // users to modify if they wish.
202  Array<int> ess_tdof_list, ess_bdr(pmesh->bdr_attributes.Max());
203  ess_bdr = 0;
204  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
205 
206  // 8. Set up the parallel linear form b(.) which corresponds to the
207  // right-hand side of the FEM linear system. In this case, b_i equals the
208  // boundary integral of f*phi_i where f represents a "pull down" force on
209  // the Neumann part of the boundary and phi_i are the basis functions in
210  // the finite element fespace. The force is defined by the object f, which
211  // is a vector of Coefficient objects. The fact that f is non-zero on
212  // boundary attribute 2 is indicated by the use of piece-wise constants
213  // coefficient for its last component.
215  for (int i = 0; i < dim-1; i++)
216  {
217  f.Set(i, new ConstantCoefficient(0.0));
218  }
219 
220  // 9. Put a leftward force on the right side of the trapezoid
221  {
222  Vector push_force(pmesh->bdr_attributes.Max());
223  push_force = 0.0;
224  push_force(1) = -5.0e-2; // index 1 attribute 2
225  f.Set(0, new PWConstCoefficient(push_force));
226  }
227 
228  ParLinearForm *b = new ParLinearForm(fespace);
230  b->Assemble();
231 
232  // 10. Define the solution vector x as a parallel finite element grid
233  // function corresponding to fespace. Initialize x with initial guess of
234  // zero, which satisfies the boundary conditions.
235  ParGridFunction x(fespace);
236  x = 0.0;
237 
238  // 11. Set up the parallel bilinear form a(.,.) on the finite element space
239  // corresponding to the linear elasticity integrator with piece-wise
240  // constants coefficient lambda and mu. We use constant coefficients,
241  // but see ex2 for how to set up piecewise constant coefficients based
242  // on attribute.
243  Vector lambda(pmesh->attributes.Max());
244  lambda = 1.0;
245  PWConstCoefficient lambda_func(lambda);
246  Vector mu(pmesh->attributes.Max());
247  mu = 1.0;
248  PWConstCoefficient mu_func(mu);
249  ParBilinearForm *a = new ParBilinearForm(fespace);
250  a->AddDomainIntegrator(new ElasticityIntegrator(lambda_func, mu_func));
251 
252  // 12. Assemble the parallel bilinear form and the corresponding linear
253  // system, applying any necessary transformations such as: parallel
254  // assembly, eliminating boundary conditions, applying conforming
255  // constraints for non-conforming AMR, etc.
256  a->Assemble();
257 
258  HypreParMatrix A;
259  Vector B, X;
260  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
261  if (myid == 0)
262  {
263  cout << "done." << endl;
264  cout << "Size of linear system: " << A.GetGlobalNumRows() << endl;
265  }
266 
267  // 13. Set up constraint matrix to constrain normal displacement (but
268  // allow tangential displacement) on specified boundaries.
269  Array<int> constraint_atts(2);
270  constraint_atts[0] = 1; // attribute 1 bottom
271  constraint_atts[1] = 4; // attribute 4 left side
272  Array<int> constraint_rowstarts;
273  SparseMatrix* local_constraints =
274  ParBuildNormalConstraints(*fespace, constraint_atts,
275  constraint_rowstarts);
276 
277  // 14. Define and apply a parallel PCG solver for the constrained system
278  // where the normal boundary constraints have been separately eliminated
279  // from the system.
280  ConstrainedSolver * solver;
281  if (penalty == 0.0)
282  {
283  solver = new EliminationCGSolver(A, *local_constraints,
284  constraint_rowstarts, dim,
285  reorder_space);
286  }
287  else
288  {
289  solver = new PenaltyPCGSolver(A, *local_constraints, penalty,
290  dim, reorder_space);
291  }
292 
293  solver->SetRelTol(1e-8);
294  solver->SetMaxIter(500);
295  solver->SetPrintLevel(1);
296  solver->Mult(B, X);
297 
298  // 15. Recover the parallel grid function corresponding to X. This is the
299  // local finite element solution on each processor.
300  a->RecoverFEMSolution(X, *b, x);
301 
302  // 16. For non-NURBS meshes, make the mesh curved based on the finite element
303  // space. This means that we define the mesh elements through a fespace
304  // based transformation of the reference element. This allows us to save
305  // the displaced mesh as a curved mesh when using high-order finite
306  // element displacement field. We assume that the initial mesh (read from
307  // the file) is not higher order curved mesh compared to the chosen FE
308  // space.
309  if (!use_nodal_fespace)
310  {
311  pmesh->SetNodalFESpace(fespace);
312  }
313 
314  GridFunction *nodes = pmesh->GetNodes();
315  *nodes += x;
316 
317  // 17. Save the refined mesh and the solution in VisIt format.
318  if (visit)
319  {
320  VisItDataCollection visit_dc(MPI_COMM_WORLD, "ex28p", pmesh);
321  visit_dc.SetLevelsOfDetail(4);
322  visit_dc.RegisterField("displacement", &x);
323  visit_dc.Save();
324  }
325 
326  // 18. Save in parallel the displaced mesh and the inverted solution (which
327  // gives the backward displacements to the original grid). This output
328  // can be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
329  {
330  x *= -1; // sign convention for GLVis displacements
331 
332  ostringstream mesh_name, sol_name;
333  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
334  sol_name << "sol." << setfill('0') << setw(6) << myid;
335 
336  ofstream mesh_ofs(mesh_name.str().c_str());
337  mesh_ofs.precision(8);
338  pmesh->Print(mesh_ofs);
339 
340  ofstream sol_ofs(sol_name.str().c_str());
341  sol_ofs.precision(8);
342  x.Save(sol_ofs);
343  }
344 
345  // 19. Send the above data by socket to a GLVis server. Use the "n" and "b"
346  // keys in GLVis to visualize the displacements.
347  if (visualization)
348  {
349  char vishost[] = "localhost";
350  int visport = 19916;
351  socketstream sol_sock(vishost, visport);
352  sol_sock << "parallel " << num_procs << " " << myid << "\n";
353  sol_sock.precision(8);
354  sol_sock << "solution\n" << *pmesh << x << flush;
355  }
356 
357  // 20. Free the used memory.
358  delete local_constraints;
359  delete solver;
360  delete a;
361  delete b;
362  if (fec)
363  {
364  delete fespace;
365  delete fec;
366  }
367  delete pmesh;
368 
369  // HYPRE_Finalize();
370  MPI_Finalize();
371 
372  return 0;
373 }
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:790
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:1324
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(...
SparseMatrix * ParBuildNormalConstraints(ParFiniteElementSpace &fespace, Array< int > &constrained_att, Array< int > &constraint_rowstarts)
Parallel wrapper for BuildNormalConstraints.
An abstract class to solve the constrained system subject to the constraint .
Definition: constraints.hpp:53
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:275
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:846
virtual void Save()
Save the collection and a VisIt root file.
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:841
Abstract parallel finite element space.
Definition: pfespace.hpp:28
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:1263
Class for parallel linear form.
Definition: plinearform.hpp:26
void SetPrintLevel(int print_lvl)
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:150
Data type sparse matrix.
Definition: sparsemat.hpp:41
constexpr char vishost[]
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9143
constexpr int visport
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
Definition: solvers.hpp:100
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
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:4382
int AddBdrSegment(int v1, int v2, int attr=1)
Definition: mesh.cpp:1440
void SetLevelsOfDetail(int levels_of_detail)
Set VisIt parameter: default levels of detail for the MultiresControl.
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:4843
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:70
int Dimension() const
Definition: mesh.hpp:911
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:204
void SetRelTol(double rtol)
Definition: solvers.hpp:98
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
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 RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:206
A piecewise constant coefficient with the constants keyed off the element attribute numbers...
Definition: coefficient.hpp:94
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:573
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a quadrilateral Mesh.
Definition: mesh.cpp:1585
double mu
Definition: ex25.cpp:139
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:327
Mesh * build_trapezoid_mesh(double offset)
Definition: ex28.cpp:38
Class for parallel bilinear form.
Vector data type.
Definition: vector.hpp:60
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7343
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:277
Class for parallel meshes.
Definition: pmesh.hpp:32
int main()
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:202
double f(const Vector &p)
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150