MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex29p.cpp
Go to the documentation of this file.
1 // MFEM Example 29 - Parallel Version
2 //
3 // Compile with: make ex29p
4 //
5 // Sample runs: mpirun -np 4 ex29p
6 // mpirun -np 4 ex29p -sc
7 // mpirun -np 4 ex29p -mt 3 -o 3 -sc
8 // mpirun -np 4 ex29p -mt 3 -rs 1 -o 4 -sc
9 //
10 // Description: This example code demonstrates the use of MFEM to define a
11 // finite element discretization of a PDE on a 2 dimensional
12 // surface embedded in a 3 dimensional domain. In this case we
13 // solve the Laplace problem -Div(sigma Grad u) = 1, with
14 // homogeneous Dirichlet boundary conditions, where sigma is an
15 // anisotropic diffusion constant defined as a 3x3 matrix
16 // coefficient.
17 //
18 // This example demonstrates the use of finite element integrators
19 // on 2D domains with 3D coefficients.
20 //
21 // We recommend viewing examples 1 and 7 before viewing this
22 // example.
23 
24 #include "mfem.hpp"
25 #include <fstream>
26 #include <iostream>
27 
28 using namespace std;
29 using namespace mfem;
30 
31 Mesh * GetMesh(int type);
32 
33 void trans(const Vector &x, Vector &r);
34 
35 void sigmaFunc(const Vector &x, DenseMatrix &s);
36 
37 double uExact(const Vector &x)
38 {
39  return (0.25 * (2.0 + x[0]) - x[2]) * (x[2] + 0.25 * (2.0 + x[0]));
40 }
41 
42 void duExact(const Vector &x, Vector &du)
43 {
44  du.SetSize(3);
45  du[0] = 0.125 * (2.0 + x[0]) * x[1] * x[1];
46  du[1] = -0.125 * (2.0 + x[0]) * x[0] * x[1];
47  du[2] = -2.0 * x[2];
48 }
49 
50 void fluxExact(const Vector &x, Vector &f)
51 {
52  f.SetSize(3);
53 
54  DenseMatrix s(3);
55  sigmaFunc(x, s);
56 
57  Vector du(3);
58  duExact(x, du);
59 
60  s.Mult(du, f);
61  f *= -1.0;
62 }
63 
64 int main(int argc, char *argv[])
65 {
66  // 1. Initialize MPI.
67  MPI_Session mpi(argc, argv);
68  int num_procs = mpi.WorldSize();
69  int myid = mpi.WorldRank();
70 
71  // 2. Parse command-line options.
72  int order = 3;
73  int mesh_type = 4; // Default to Quadrilateral mesh
74  int mesh_order = 3;
75  int ser_ref_levels = 2;
76  int par_ref_levels = 1;
77  bool static_cond = false;
78  bool visualization = true;
79 
80  OptionsParser args(argc, argv);
81  args.AddOption(&mesh_type, "-mt", "--mesh-type",
82  "Mesh type: 3 - Triangular, 4 - Quadrilateral.");
83  args.AddOption(&mesh_order, "-mo", "--mesh-order",
84  "Geometric order of the curved mesh.");
85  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
86  "Number of times to refine the mesh uniformly in serial.");
87  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
88  "Number of times to refine the mesh uniformly in parallel.");
89  args.AddOption(&order, "-o", "--order",
90  "Finite element order (polynomial degree) or -1 for"
91  " isoparametric space.");
92  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
93  "--no-static-condensation", "Enable static condensation.");
94  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
95  "--no-visualization",
96  "Enable or disable GLVis visualization.");
97  args.ParseCheck();
98 
99  // 3. Construct a quadrilateral or triangular mesh with the topology of a
100  // cylindrical surface.
101  Mesh *mesh = GetMesh(mesh_type);
102  int dim = mesh->Dimension();
103 
104  // 4. Refine the mesh to increase the resolution. In this example we do
105  // 'ser_ref_levels' of uniform refinement.
106  for (int l = 0; l < ser_ref_levels; l++)
107  {
108  mesh->UniformRefinement();
109  }
110 
111  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
112  // this mesh further in parallel to increase the resolution. Once the
113  // parallel mesh is defined, the serial mesh can be deleted.
114  ParMesh pmesh(MPI_COMM_WORLD, *mesh);
115  delete mesh;
116  for (int l = 0; l < par_ref_levels; l++)
117  {
118  pmesh.UniformRefinement();
119  }
120 
121  // 6. Transform the mesh so that it has a more interesting geometry.
122  pmesh.SetCurvature(mesh_order);
123  pmesh.Transform(trans);
124 
125  // 7. Define a finite element space on the mesh. Here we use continuous
126  // Lagrange finite elements of the specified order.
127  H1_FECollection fec(order, dim);
128  ParFiniteElementSpace fespace(&pmesh, &fec);
129  HYPRE_Int total_num_dofs = fespace.GlobalTrueVSize();
130  if (mpi.Root()) { cout << "Number of unknowns: " << total_num_dofs << endl; }
131 
132  // 8. Determine the list of true (i.e. conforming) essential boundary dofs.
133  // In this example, the boundary conditions are defined by marking all
134  // the boundary attributes from the mesh as essential (Dirichlet) and
135  // converting them to a list of true dofs.
136  Array<int> ess_tdof_list;
137  if (pmesh.bdr_attributes.Size())
138  {
139  Array<int> ess_bdr(pmesh.bdr_attributes.Max());
140  ess_bdr = 1;
141  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
142  }
143 
144  // 9. Set up the linear form b(.) which corresponds to the right-hand side of
145  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
146  // the basis functions in the finite element fespace.
147  ParLinearForm b(&fespace);
148  ConstantCoefficient one(1.0);
150  b.Assemble();
151 
152  // 10. Define the solution vector x as a finite element grid function
153  // corresponding to fespace. Initialize x with initial guess of zero,
154  // which satisfies the boundary conditions.
155  ParGridFunction x(&fespace);
156  x = 0.0;
157 
158  // 11. Set up the bilinear form a(.,.) on the finite element space
159  // corresponding to the Laplacian operator -Delta, by adding the
160  // Diffusion domain integrator.
161  ParBilinearForm a(&fespace);
163  BilinearFormIntegrator *integ = new DiffusionIntegrator(sigma);
164  a.AddDomainIntegrator(integ);
165 
166  // 12. Assemble the bilinear form and the corresponding linear system,
167  // applying any necessary transformations such as: eliminating boundary
168  // conditions, applying conforming constraints for non-conforming AMR,
169  // static condensation, etc.
170  if (static_cond) { a.EnableStaticCondensation(); }
171  a.Assemble();
172 
173  OperatorPtr A;
174  Vector B, X;
175  a.FormLinearSystem(ess_tdof_list, x, b, A, X, B);
176 
177  if (myid == 0)
178  {
179  cout << "Size of linear system: "
180  << A.As<HypreParMatrix>()->GetGlobalNumRows() << endl;
181  }
182 
183  // 13. Define and apply a parallel PCG solver for A X = B with the BoomerAMG
184  // preconditioner from hypre.
185  HypreBoomerAMG *amg = new HypreBoomerAMG;
186  CGSolver cg(MPI_COMM_WORLD);
187  cg.SetRelTol(1e-12);
188  cg.SetMaxIter(2000);
189  cg.SetPrintLevel(1);
190  cg.SetPreconditioner(*amg);
191  cg.SetOperator(*A);
192  cg.Mult(B, X);
193  delete amg;
194 
195  // 14. Recover the solution as a finite element grid function.
196  a.RecoverFEMSolution(X, b, x);
197 
198  // 15. Compute error in the solution and its flux
200  double err = x.ComputeL2Error(uCoef);
201 
202  if (myid == 0) { cout << "|u - u_h|_2 = " << err << endl; }
203 
204  ParFiniteElementSpace flux_fespace(&pmesh, &fec, 3);
205  ParGridFunction flux(&flux_fespace);
206  x.ComputeFlux(*integ, flux); flux *= -1.0;
207 
209  double flux_err = flux.ComputeL2Error(fluxCoef);
210 
211  if (myid == 0) { cout << "|f - f_h|_2 = " << flux_err << endl; }
212 
213  // 16. Save the refined mesh and the solution. This output can be viewed
214  // later using GLVis: "glvis -np <np> -m mesh -g sol".
215  {
216  ostringstream mesh_name, sol_name, flux_name;
217  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
218  sol_name << "sol." << setfill('0') << setw(6) << myid;
219  flux_name << "flux." << setfill('0') << setw(6) << myid;
220 
221  ofstream mesh_ofs(mesh_name.str().c_str());
222  mesh_ofs.precision(8);
223  pmesh.Print(mesh_ofs);
224 
225  ofstream sol_ofs(sol_name.str().c_str());
226  sol_ofs.precision(8);
227  x.Save(sol_ofs);
228 
229  ofstream flux_ofs(flux_name.str().c_str());
230  flux_ofs.precision(8);
231  flux.Save(flux_ofs);
232  }
233 
234  // 17. Send the solution by socket to a GLVis server.
235  if (visualization)
236  {
237  char vishost[] = "localhost";
238  int visport = 19916;
239  socketstream sol_sock(vishost, visport);
240  sol_sock << "parallel " << num_procs << " " << myid << "\n";
241  sol_sock.precision(8);
242  sol_sock << "solution\n" << pmesh << x
243  << "window_title 'Solution'\n" << flush;
244 
245  socketstream flux_sock(vishost, visport);
246  flux_sock << "parallel " << num_procs << " " << myid << "\n";
247  flux_sock.precision(8);
248  flux_sock << "solution\n" << pmesh << flux
249  << "keys vvv\n"
250  << "window_geometry 402 0 400 350\n"
251  << "window_title 'Flux'\n" << flush;
252  }
253 
254  return 0;
255 }
256 
257 // Defines a mesh consisting of four flat rectangular surfaces connected to form
258 // a loop.
259 Mesh * GetMesh(int type)
260 {
261  Mesh * mesh = NULL;
262 
263  if (type == 3)
264  {
265  mesh = new Mesh(2, 12, 16, 8, 3);
266 
267  mesh->AddVertex(-1.0, -1.0, 0.0);
268  mesh->AddVertex( 1.0, -1.0, 0.0);
269  mesh->AddVertex( 1.0, 1.0, 0.0);
270  mesh->AddVertex(-1.0, 1.0, 0.0);
271  mesh->AddVertex(-1.0, -1.0, 1.0);
272  mesh->AddVertex( 1.0, -1.0, 1.0);
273  mesh->AddVertex( 1.0, 1.0, 1.0);
274  mesh->AddVertex(-1.0, 1.0, 1.0);
275  mesh->AddVertex( 0.0, -1.0, 0.5);
276  mesh->AddVertex( 1.0, 0.0, 0.5);
277  mesh->AddVertex( 0.0, 1.0, 0.5);
278  mesh->AddVertex(-1.0, 0.0, 0.5);
279 
280  mesh->AddTriangle(0, 1, 8);
281  mesh->AddTriangle(1, 5, 8);
282  mesh->AddTriangle(5, 4, 8);
283  mesh->AddTriangle(4, 0, 8);
284  mesh->AddTriangle(1, 2, 9);
285  mesh->AddTriangle(2, 6, 9);
286  mesh->AddTriangle(6, 5, 9);
287  mesh->AddTriangle(5, 1, 9);
288  mesh->AddTriangle(2, 3, 10);
289  mesh->AddTriangle(3, 7, 10);
290  mesh->AddTriangle(7, 6, 10);
291  mesh->AddTriangle(6, 2, 10);
292  mesh->AddTriangle(3, 0, 11);
293  mesh->AddTriangle(0, 4, 11);
294  mesh->AddTriangle(4, 7, 11);
295  mesh->AddTriangle(7, 3, 11);
296 
297  mesh->AddBdrSegment(0, 1, 1);
298  mesh->AddBdrSegment(1, 2, 1);
299  mesh->AddBdrSegment(2, 3, 1);
300  mesh->AddBdrSegment(3, 0, 1);
301  mesh->AddBdrSegment(5, 4, 2);
302  mesh->AddBdrSegment(6, 5, 2);
303  mesh->AddBdrSegment(7, 6, 2);
304  mesh->AddBdrSegment(4, 7, 2);
305  }
306  else if (type == 4)
307  {
308  mesh = new Mesh(2, 8, 4, 8, 3);
309 
310  mesh->AddVertex(-1.0, -1.0, 0.0);
311  mesh->AddVertex( 1.0, -1.0, 0.0);
312  mesh->AddVertex( 1.0, 1.0, 0.0);
313  mesh->AddVertex(-1.0, 1.0, 0.0);
314  mesh->AddVertex(-1.0, -1.0, 1.0);
315  mesh->AddVertex( 1.0, -1.0, 1.0);
316  mesh->AddVertex( 1.0, 1.0, 1.0);
317  mesh->AddVertex(-1.0, 1.0, 1.0);
318 
319  mesh->AddQuad(0, 1, 5, 4);
320  mesh->AddQuad(1, 2, 6, 5);
321  mesh->AddQuad(2, 3, 7, 6);
322  mesh->AddQuad(3, 0, 4, 7);
323 
324  mesh->AddBdrSegment(0, 1, 1);
325  mesh->AddBdrSegment(1, 2, 1);
326  mesh->AddBdrSegment(2, 3, 1);
327  mesh->AddBdrSegment(3, 0, 1);
328  mesh->AddBdrSegment(5, 4, 2);
329  mesh->AddBdrSegment(6, 5, 2);
330  mesh->AddBdrSegment(7, 6, 2);
331  mesh->AddBdrSegment(4, 7, 2);
332  }
333  else
334  {
335  MFEM_ABORT("Unrecognized mesh type " << type << "!");
336  }
337  mesh->FinalizeTopology();
338 
339  return mesh;
340 }
341 
342 // Transforms the four-sided loop into a curved cylinder with skewed top and
343 // base.
344 void trans(const Vector &x, Vector &r)
345 {
346  r.SetSize(3);
347 
348  double tol = 1e-6;
349  double theta = 0.0;
350  if (fabs(x[1] + 1.0) < tol)
351  {
352  theta = 0.25 * M_PI * (x[0] - 2.0);
353  }
354  else if (fabs(x[0] - 1.0) < tol)
355  {
356  theta = 0.25 * M_PI * x[1];
357  }
358  else if (fabs(x[1] - 1.0) < tol)
359  {
360  theta = 0.25 * M_PI * (2.0 - x[0]);
361  }
362  else if (fabs(x[0] + 1.0) < tol)
363  {
364  theta = 0.25 * M_PI * (4.0 - x[1]);
365  }
366  else
367  {
368  cerr << "side not recognized "
369  << x[0] << " " << x[1] << " " << x[2] << endl;
370  }
371 
372  r[0] = cos(theta);
373  r[1] = sin(theta);
374  r[2] = 0.25 * (2.0 * x[2] - 1.0) * (r[0] + 2.0);
375 }
376 
377 // Anisotropic diffusion coefficient
378 void sigmaFunc(const Vector &x, DenseMatrix &s)
379 {
380  s.SetSize(3);
381  double a = 17.0 - 2.0 * x[0] * (1.0 + x[0]);
382  s(0,0) = 0.5 + x[0] * x[0] * (8.0 / a - 0.5);
383  s(0,1) = x[0] * x[1] * (8.0 / a - 0.5);
384  s(0,2) = 0.0;
385  s(1,0) = s(0,1);
386  s(1,1) = 0.5 * x[0] * x[0] + 8.0 * x[1] * x[1] / a;
387  s(1,2) = 0.0;
388  s(2,0) = 0.0;
389  s(2,1) = 0.0;
390  s(2,2) = a / 32.0;
391 }
int WorldSize() const
Return MPI_COMM_WORLD&#39;s size.
int Size() const
Return the logical size of the array.
Definition: array.hpp:134
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:97
void duExact(const Vector &x, Vector &du)
Definition: ex29.cpp:42
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:790
Conjugate gradient method.
Definition: solvers.hpp:316
void trans(const Vector &u, Vector &x)
Definition: ex27.cpp:421
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:99
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(...
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: pmesh.cpp:1958
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:616
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:275
void sigmaFunc(const Vector &x, DenseMatrix &s)
Definition: ex29.cpp:337
int AddTriangle(int v1, int v2, int v3, int attr=1)
Definition: mesh.cpp:1310
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
void Transform(void(*f)(const Vector &, Vector &))
Definition: mesh.cpp:11008
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:841
Abstract parallel finite element space.
Definition: pfespace.hpp:28
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1387
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
A simple convenience class that calls MPI_Init() at construction and MPI_Finalize() at destruction...
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:71
A matrix coefficient with an optional scalar coefficient multiplier q. The matrix function can either...
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
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.
bool Root() const
Return true if WorldRank() == 0.
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
int Dimension() const
Definition: mesh.hpp:911
A general vector function coefficient.
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:34
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
int WorldRank() const
Return MPI_COMM_WORLD&#39;s rank.
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
Definition: globals.hpp:71
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, bool wcoef=true, int subdomain=-1)
Definition: pgridfunc.cpp:1068
void FinalizeTopology(bool generate_bdr=true)
Finalize the construction of the secondary topology (connectivity) data of a Mesh.
Definition: mesh.cpp:2507
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
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:273
double uExact(const Vector &x)
Definition: ex29.cpp:37
Mesh * GetMesh(int type)
Definition: ex29.cpp:218
int dim
Definition: ex24.cpp:53
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
Class for parallel bilinear form.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:330
A general function coefficient.
Vector data type.
Definition: vector.hpp:60
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:175
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:92
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
RefCoord s[3]
Class for parallel grid function.
Definition: pgridfunc.hpp:32
void fluxExact(const Vector &x, Vector &f)
Definition: ex29.cpp:50
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:277
Class for parallel meshes.
Definition: pmesh.hpp:32
int main()
void ParseCheck(std::ostream &out=mfem::out)
Definition: optparser.cpp:251
void EnableStaticCondensation()
Enable the use of static condensation. For details see the description for class StaticCondensation i...
double f(const Vector &p)
double sigma(const Vector &x)
Definition: maxwell.cpp:102