MFEM  v3.3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex17p.cpp
Go to the documentation of this file.
1 // MFEM Example 17 - Parallel Version
2 //
3 // Compile with: make ex17p
4 //
5 // Sample runs:
6 //
7 // mpirun -np 4 ex17p -m ../data/beam-tri.mesh
8 // mpirun -np 4 ex17p -m ../data/beam-quad.mesh
9 // mpirun -np 4 ex17p -m ../data/beam-tet.mesh
10 // mpirun -np 4 ex17p -m ../data/beam-hex.mesh
11 // mpirun -np 4 ex17p -m ../data/beam-quad.mesh -rs 2 -rp 2 -o 3 -elast
12 // mpirun -np 4 ex17p -m ../data/beam-quad.mesh -rs 2 -rp 3 -o 2 -a 1 -k 1
13 // mpirun -np 4 ex17p -m ../data/beam-hex.mesh -rs 2 -rp 1 -o 2
14 //
15 // Description: This example code solves a simple linear elasticity problem
16 // describing a multi-material cantilever beam using symmetric or
17 // non-symmetric discontinuous Galerkin (DG) formulation.
18 //
19 // Specifically, we approximate the weak form of -div(sigma(u))=0
20 // where sigma(u)=lambda*div(u)*I+mu*(grad*u+u*grad) is the stress
21 // tensor corresponding to displacement field u, and lambda and mu
22 // are the material Lame constants. The boundary conditions are
23 // Dirichlet, u=u_D on the fixed part of the boundary, namely
24 // boundary attributes 1 and 2; on the rest of the boundary we use
25 // sigma(u).n=0 b.c. The geometry of the domain is assumed to be
26 // as follows:
27 //
28 // +----------+----------+
29 // boundary --->| material | material |<--- boundary
30 // attribute 1 | 1 | 2 | attribute 2
31 // (fixed) +----------+----------+ (fixed, nonzero)
32 //
33 // The example demonstrates the use of high-order DG vector finite
34 // element spaces with the linear DG elasticity bilinear form,
35 // meshes with curved elements, and the definition of piece-wise
36 // constant and function vector-coefficient objects. The use of
37 // non-homogeneous Dirichlet b.c. imposed weakly, is also
38 // illustrated.
39 //
40 // We recommend viewing examples 2p and 14p before viewing this
41 // example.
42 
43 #include "mfem.hpp"
44 #include <fstream>
45 #include <iostream>
46 
47 using namespace std;
48 using namespace mfem;
49 
50 // Initial displacement, used for Dirichlet boundary conditions on boundary
51 // attributes 1 and 2.
52 void InitDisplacement(const Vector &x, Vector &u);
53 
54 // A Coefficient for computing the components of the stress.
55 class StressCoefficient : public Coefficient
56 {
57 protected:
58  Coefficient &lambda, &mu;
59  GridFunction *u; // displacement
60  int si, sj; // component of the stress to evaluate, 0 <= si,sj < dim
61 
62  DenseMatrix grad; // auxiliary matrix, used in Eval
63 
64 public:
65  StressCoefficient(Coefficient &lambda_, Coefficient &mu_)
66  : lambda(lambda_), mu(mu_), u(NULL), si(0), sj(0) { }
67 
68  void SetDisplacement(GridFunction &u_) { u = &u_; }
69  void SetComponent(int i, int j) { si = i; sj = j; }
70 
71  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
72 };
73 
74 // Simple GLVis visualization manager.
75 class VisMan : public iostream
76 {
77 protected:
78  const char *host;
79  int port;
81  int sid; // active socket, index inside 'sock'.
82 
83  int win_x, win_y, win_w, win_h;
84  int win_stride_x, win_stride_y, win_nx;
85 
86 public:
87  VisMan(const char *vishost, const int visport);
88  void NewWindow();
89  void CloseConnection();
90  void PositionWindow();
91  virtual ~VisMan();
92 };
93 
94 // Manipulators for the GLVis visualization manager.
95 void new_window (VisMan &v) { v.NewWindow(); }
96 void position_window (VisMan &v) { v.PositionWindow(); }
97 void close_connection(VisMan &v) { v.CloseConnection(); }
98 ostream &operator<<(ostream &v, void (*f)(VisMan&));
99 
100 int main(int argc, char *argv[])
101 {
102  MPI_Session mpi(argc, argv);
103 
104  // 1. Define and parse command-line options.
105  const char *mesh_file = "../data/beam-tri.mesh";
106  int ser_ref_levels = -1;
107  int par_ref_levels = 1;
108  int order = 1;
109  double alpha = -1.0;
110  double kappa = -1.0;
111  bool amg_elast = false;
112  bool visualization = 1;
113 
114  OptionsParser args(argc, argv);
115  args.AddOption(&mesh_file, "-m", "--mesh",
116  "Mesh file to use.");
117  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
118  "Number of times to refine the mesh uniformly before parallel"
119  " partitioning, -1 for auto.");
120  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
121  "Number of times to refine the mesh uniformly after parallel"
122  " partitioning.");
123  args.AddOption(&order, "-o", "--order",
124  "Finite element order (polynomial degree).");
125  args.AddOption(&alpha, "-a", "--alpha",
126  "One of the two DG penalty parameters, typically +1/-1."
127  " See the documentation of class DGElasticityIntegrator.");
128  args.AddOption(&kappa, "-k", "--kappa",
129  "One of the two DG penalty parameters, should be positive."
130  " Negative values are replaced with (order+1)^2.");
131  args.AddOption(&amg_elast, "-elast", "--amg-for-elasticity", "-sys",
132  "--amg-for-systems",
133  "Use the special AMG elasticity solver (GM/LN approaches), "
134  "or standard AMG for systems (unknown approach).");
135  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
136  "--no-visualization",
137  "Enable or disable GLVis visualization.");
138  args.Parse();
139  if (!args.Good())
140  {
141  if (mpi.Root()) { args.PrintUsage(cout); }
142  return 1;
143  }
144  if (kappa < 0)
145  {
146  kappa = (order+1)*(order+1);
147  }
148  if (mpi.Root()) { args.PrintOptions(cout); }
149 
150  // 2. Read the mesh from the given mesh file.
151  Mesh mesh(mesh_file, 1, 1);
152  int dim = mesh.Dimension();
153 
154  if (mesh.attributes.Max() < 2 || mesh.bdr_attributes.Max() < 2)
155  {
156  if (mpi.Root())
157  {
158  cerr << "\nInput mesh should have at least two materials and "
159  << "two boundary attributes! (See schematic in ex17p.cpp)\n"
160  << endl;
161  }
162  return 3;
163  }
164 
165  // 3. Refine the mesh to increase the resolution.
166  if (ser_ref_levels < 0)
167  {
168  ser_ref_levels = (int)floor(log(5000./mesh.GetNE())/log(2.)/dim);
169  }
170  for (int l = 0; l < ser_ref_levels; l++)
171  {
172  mesh.UniformRefinement();
173  }
174  // Since NURBS meshes do not support DG integrators, we convert them to
175  // regular polynomial mesh of the specified (solution) order.
176  if (mesh.NURBSext) { mesh.SetCurvature(order); }
177 
178  ParMesh pmesh(MPI_COMM_WORLD, mesh);
179  mesh.Clear();
180 
181  for (int l = 0; l < par_ref_levels; l++)
182  {
183  pmesh.UniformRefinement();
184  }
185 
186  // 4. Define a DG vector finite element space on the mesh. Here, we use
187  // Gauss-Lobatto nodal basis because it gives rise to a sparser matrix
188  // compared to the default Gauss-Legendre nodal basis.
189  DG_FECollection fec(order, dim, BasisType::GaussLobatto);
190  ParFiniteElementSpace fespace(&pmesh, &fec, dim, Ordering::byVDIM);
191 
192  HYPRE_Int glob_size = fespace.GlobalTrueVSize();
193  if (mpi.Root())
194  {
195  cout << "Number of finite element unknowns: " << glob_size
196  << "\nAssembling: " << flush;
197  }
198 
199  // 5. In this example, the Dirichlet boundary conditions are defined by
200  // marking boundary attributes 1 and 2 in the marker Array 'dir_bdr'.
201  // These b.c. are imposed weakly, by adding the appropriate boundary
202  // integrators over the marked 'dir_bdr' to the bilinear and linear forms.
203  // With this DG formulation, there are no essential boundary conditions.
204  Array<int> ess_tdof_list; // no essential b.c. (empty list)
205  Array<int> dir_bdr(pmesh.bdr_attributes.Max());
206  dir_bdr = 0;
207  dir_bdr[0] = 1; // boundary attribute 1 is Dirichlet
208  dir_bdr[1] = 1; // boundary attribute 2 is Dirichlet
209 
210  // 6. Define the DG solution vector 'x' as a finite element grid function
211  // corresponding to fespace. Initialize 'x' using the 'InitDisplacement'
212  // function.
213  ParGridFunction x(&fespace);
215  x.ProjectCoefficient(init_x);
216 
217  // 7. Set up the Lame constants for the two materials. They are defined as
218  // piece-wise (with respect to the element attributes) constant
219  // coefficients, i.e. type PWConstCoefficient.
220  Vector lambda(pmesh.attributes.Max());
221  lambda = 1.0; // Set lambda = 1 for all element attributes.
222  lambda(0) = 50.0; // Set lambda = 50 for element attribute 1.
223  PWConstCoefficient lambda_c(lambda);
224  Vector mu(pmesh.attributes.Max());
225  mu = 1.0; // Set mu = 1 for all element attributes.
226  mu(0) = 50.0; // Set mu = 50 for element attribute 1.
227  PWConstCoefficient mu_c(mu);
228 
229  // 8. Set up the linear form b(.) which corresponds to the right-hand side of
230  // the FEM linear system. In this example, the linear form b(.) consists
231  // only of the terms responsible for imposing weakly the Dirichlet
232  // boundary conditions, over the attributes marked in 'dir_bdr'. The
233  // values for the Dirichlet boundary condition are taken from the
234  // VectorFunctionCoefficient 'x_init' which in turn is based on the
235  // function 'InitDisplacement'.
236  ParLinearForm b(&fespace);
237  if (mpi.Root()) { cout << "r.h.s. ... " << flush; }
240  init_x, lambda_c, mu_c, alpha, kappa), dir_bdr);
241  b.Assemble();
242 
243  // 9. Set up the bilinear form a(.,.) on the DG finite element space
244  // corresponding to the linear elasticity integrator with coefficients
245  // lambda and mu as defined above. The additional interior face integrator
246  // ensures the weak continuity of the displacement field. The additional
247  // boundary face integrator works together with the boundary integrator
248  // added to the linear form b(.) to impose weakly the Dirichlet boundary
249  // conditions.
250  ParBilinearForm a(&fespace);
251  a.AddDomainIntegrator(new ElasticityIntegrator(lambda_c, mu_c));
253  new DGElasticityIntegrator(lambda_c, mu_c, alpha, kappa));
255  new DGElasticityIntegrator(lambda_c, mu_c, alpha, kappa), dir_bdr);
256 
257  // 10. Assemble the bilinear form and the corresponding linear system.
258  if (mpi.Root()) { cout << "matrix ... " << flush; }
259  a.Assemble();
260 
261  HypreParMatrix A;
262  Vector B, X;
263  a.FormLinearSystem(ess_tdof_list, x, b, A, X, B);
264  if (mpi.Root()) { cout << "done." << endl; }
265 
266  // 11. Define a simple symmetric Gauss-Seidel preconditioner and use it to
267  // solve the system Ax=b with PCG for the symmetric formulation, or GMRES
268  // for the non-symmetric.
269  const double rtol = 1e-6;
270  HypreBoomerAMG amg(A);
271  if (amg_elast)
272  {
273  amg.SetElasticityOptions(&fespace);
274  }
275  else
276  {
277  amg.SetSystemsOptions(dim);
278  }
279  CGSolver pcg(A.GetComm());
280  GMRESSolver gmres(A.GetComm());
281  gmres.SetKDim(50);
282  IterativeSolver &ipcg = pcg, &igmres = gmres;
283  IterativeSolver &solver = (alpha == -1.0) ? ipcg : igmres;
284  solver.SetRelTol(rtol);
285  solver.SetMaxIter(500);
286  solver.SetPrintLevel(1);
287  solver.SetOperator(A);
288  solver.SetPreconditioner(amg);
289  solver.Mult(B, X);
290 
291  // 12. Recover the solution as a finite element grid function 'x'.
292  a.RecoverFEMSolution(X, b, x);
293 
294  // 13. Use the DG solution space as the mesh nodal space. This allows us to
295  // save the displaced mesh as a curved DG mesh.
296  pmesh.SetNodalFESpace(&fespace);
297 
298  Vector reference_nodes;
299  if (visualization) { reference_nodes = *pmesh.GetNodes(); }
300 
301  // 14. Save the displaced mesh and minus the solution (which gives the
302  // backward displacements to the reference mesh). This output can be
303  // viewed later using GLVis: "glvis -m displaced.mesh -g sol.gf".
304  {
305  *pmesh.GetNodes() += x;
306  x.Neg(); // x = -x
307 
308  ostringstream mesh_name, sol_name;
309  mesh_name << "mesh." << setfill('0') << setw(6) << mpi.WorldRank();
310  sol_name << "sol." << setfill('0') << setw(6) << mpi.WorldRank();
311 
312  ofstream mesh_ofs(mesh_name.str().c_str());
313  mesh_ofs.precision(8);
314  mesh_ofs << pmesh;
315 
316  ofstream sol_ofs(sol_name.str().c_str());
317  sol_ofs.precision(8);
318  sol_ofs << x;
319  }
320 
321  // 15. Visualization: send data by socket to a GLVis server.
322  if (visualization)
323  {
324  char vishost[] = "localhost";
325  int visport = 19916;
326  VisMan vis(vishost, visport);
327  const char *glvis_keys = (dim < 3) ? "Rjlc" : "c";
328 
329  // Visualize the deformed configuration.
330  vis << new_window << setprecision(8)
331  << "parallel " << pmesh.GetNRanks() << ' ' << pmesh.GetMyRank()
332  << '\n'
333  << "solution\n" << pmesh << x << flush
334  << "keys " << glvis_keys << endl
335  << "window_title 'Deformed configuration'" << endl
336  << "plot_caption 'Backward displacement'" << endl
338 
339  // Visualize the stress components.
340  const char *c = "xyz";
341  ParFiniteElementSpace scalar_dg_space(&pmesh, &fec);
342  ParGridFunction stress(&scalar_dg_space);
343  StressCoefficient stress_c(lambda_c, mu_c);
344  *pmesh.GetNodes() = reference_nodes;
345  x.Neg(); // x = -x
346  stress_c.SetDisplacement(x);
347  for (int si = 0; si < dim; si++)
348  {
349  for (int sj = si; sj < dim; sj++)
350  {
351  stress_c.SetComponent(si, sj);
352  stress.ProjectCoefficient(stress_c);
353 
354  MPI_Barrier(MPI_COMM_WORLD);
355  vis << new_window << setprecision(8)
356  << "parallel " << pmesh.GetNRanks() << ' ' << pmesh.GetMyRank()
357  << '\n'
358  << "solution\n" << pmesh << stress << flush
359  << "keys " << glvis_keys << endl
360  << "window_title |Stress " << c[si] << c[sj] << '|' << endl
362  }
363  }
364  }
365 
366  return 0;
367 }
368 
369 
370 void InitDisplacement(const Vector &x, Vector &u)
371 {
372  u = 0.0;
373  u(u.Size()-1) = -0.2*x(0);
374 }
375 
376 
378  const IntegrationPoint &ip)
379 {
380  MFEM_ASSERT(u != NULL, "displacement field is not set");
381 
382  double L = lambda.Eval(T, ip);
383  double M = mu.Eval(T, ip);
384  u->GetVectorGradient(T, grad);
385  if (si == sj)
386  {
387  double div_u = grad.Trace();
388  return L*div_u + 2*M*grad(si,si);
389  }
390  else
391  {
392  return M*(grad(si,sj) + grad(sj,si));
393  }
394 }
395 
396 
397 VisMan::VisMan(const char *vishost, const int visport)
398  : iostream(0),
399  host(vishost), port(visport), sid(0)
400 {
401  win_x = 0;
402  win_y = 0;
403  win_w = 400; // window width
404  win_h = 350; // window height
405  win_stride_x = win_w;
406  win_stride_y = win_h + 20;
407  win_nx = 4; // number of windows in a row
408 }
409 
410 void VisMan::NewWindow()
411 {
412  sock.Append(new socketstream(host, port));
413  sid = sock.Size()-1;
414  iostream::rdbuf(sock[sid]->rdbuf());
415 }
416 
417 void VisMan::CloseConnection()
418 {
419  if (sid < sock.Size())
420  {
421  delete sock[sid];
422  sock[sid] = NULL;
423  iostream::rdbuf(0);
424  }
425 }
426 
427 void VisMan::PositionWindow()
428 {
429  *this << "window_geometry "
430  << win_x + win_stride_x*(sid%win_nx) << ' '
431  << win_y + win_stride_y*(sid/win_nx) << ' '
432  << win_w << ' ' << win_h << endl;
433 }
434 
435 VisMan::~VisMan()
436 {
437  for (int i = sock.Size()-1; i >= 0; i--)
438  {
439  delete sock[i];
440  }
441 }
442 
443 
444 ostream &operator<<(ostream &v, void (*f)(VisMan&))
445 {
446  VisMan *vp = dynamic_cast<VisMan*>(&v);
447  if (vp) { (*f)(*vp); }
448  return v;
449 }
std::ostream & operator<<(std::ostream &out, const Mesh &mesh)
Definition: mesh.cpp:8533
double Eval(ElementTransformation &T, const IntegrationPoint &ip, double t)
Definition: coefficient.hpp:45
Conjugate gradient method.
Definition: solvers.hpp:111
MPI_Comm GetComm() const
MPI communicator.
Definition: hypre.hpp:298
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:51
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:113
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:614
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
Abstract parallel finite element space.
Definition: pfespace.hpp:31
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:275
void position_window(VisMan &v)
Definition: ex17.cpp:96
int GetNRanks() const
Definition: pmesh.hpp:118
void SetSystemsOptions(int dim)
Definition: hypre.cpp:2502
The BoomerAMG solver in hypre.
Definition: hypre.hpp:783
Class for parallel linear form.
Definition: plinearform.hpp:26
A simple convenience class that calls MPI_Init() at construction and MPI_Finalize() at destruction...
int dim
Definition: ex3.cpp:47
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:72
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:99
void new_window(VisMan &v)
Definition: ex17.cpp:95
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6718
void SetMaxIter(int max_it)
Definition: solvers.hpp:63
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3315
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:108
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:178
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void SetKDim(int dim)
Definition: solvers.hpp:155
bool Root() const
Return true if WorldRank() == 0.
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator.
Definition: linearform.cpp:38
void SetElasticityOptions(ParFiniteElementSpace *fespace)
Definition: hypre.cpp:2584
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:3298
int Dimension() const
Definition: mesh.hpp:641
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
int main()
void close_connection(VisMan &v)
Definition: ex17.cpp:97
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:172
int GetMyRank() const
Definition: pmesh.hpp:119
int WorldRank() const
Return MPI_COMM_WORLD&#39;s rank.
void SetRelTol(double rtol)
Definition: solvers.hpp:61
Abstract base class for iterative solver.
Definition: solvers.hpp:32
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
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
GMRES method.
Definition: solvers.hpp:143
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:174
class for piecewise constant coefficient
Definition: coefficient.hpp:72
Class for integration point with weight.
Definition: intrules.hpp:25
void InitDisplacement(const Vector &x, Vector &u)
Definition: ex17.cpp:332
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new interior Face Integrator.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
Class for parallel bilinear form.
void Clear()
Clear the contents of the Mesh.
Definition: mesh.hpp:599
const double alpha
Definition: ex15.cpp:337
double kappa
Definition: ex3.cpp:46
Vector data type.
Definition: vector.hpp:41
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5401
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:93
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:175
Class for parallel meshes.
Definition: pmesh.hpp:29
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:170
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator.
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:195
void Neg()
(*this) = -(*this)
Definition: vector.cpp:256
bool Good() const
Definition: optparser.hpp:120