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