MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex27.cpp
Go to the documentation of this file.
1 // MFEM Example 27 - Serial Version
2 //
3 // Compile with: make ex27
4 //
5 // Sample runs: ex27
6 // ex27 -dg
7 // ex27 -dg -dbc 8 -nbc -2
8 // ex27 -rbc-a 1 -rbc-b 8
9 //
10 // Description: This example code demonstrates the use of MFEM to define a
11 // simple finite element discretization of the Laplace problem
12 // -Delta u = 0 with a variety of boundary conditions.
13 //
14 // Specifically, we discretize using a FE space of the specified
15 // order using a continuous or discontinuous space. We then apply
16 // Dirichlet, Neumann (both homogeneous and inhomogeneous), Robin,
17 // and Periodic boundary conditions on different portions of a
18 // predefined mesh.
19 //
20 // The predefined mesh consists of a rectangle with two holes
21 // removed (see below). The narrow ends of the mesh are connected
22 // to form a Periodic boundary condition. The lower edge (tagged
23 // with attribute 1) receives an inhomogeneous Neumann boundary
24 // condition. A Robin boundary condition is applied to upper edge
25 // (attribute 2). The circular hole on the left (attribute 3)
26 // enforces a Dirichlet boundary condition. Finally, a natural
27 // boundary condition, or homogeneous Neumann BC, is applied to
28 // the circular hole on the right (attribute 4).
29 //
30 // Attribute 3 ^ y Attribute 2
31 // \ | /
32 // +-----------+-----------+
33 // | \_ | _ |
34 // | / \ | / \ |
35 // <--+---+---+---+---+---+---+--> x
36 // | \_/ | \_/ |
37 // | | \ |
38 // +-----------+-----------+ (hole radii are
39 // / | \ adjustable)
40 // Attribute 1 v Attribute 4
41 //
42 // The boundary conditions are defined as (where u is the solution
43 // field):
44 //
45 // Dirichlet: u = d
46 // Neumann: n.Grad(u) = g
47 // Robin: n.Grad(u) + a u = b
48 //
49 // The user can adjust the values of 'd', 'g', 'a', and 'b' with
50 // command line options.
51 //
52 // This example highlights the differing implementations of
53 // boundary conditions with continuous and discontinuous Galerkin
54 // formulations of the Laplace problem.
55 //
56 // We recommend viewing Examples 1 and 14 before viewing this
57 // example.
58 
59 #include "mfem.hpp"
60 #include <fstream>
61 #include <iostream>
62 
63 using namespace std;
64 using namespace mfem;
65 
66 static double a_ = 0.2;
67 
68 // Normal to hole with boundary attribute 4
69 void n4Vec(const Vector &x, Vector &n) { n = x; n[0] -= 0.5; n /= -n.Norml2(); }
70 
71 Mesh * GenerateSerialMesh(int ref);
72 
73 // Compute the average value of alpha*n.Grad(sol) + beta*sol over the boundary
74 // attributes marked in bdr_marker. Also computes the L2 norm of
75 // alpha*n.Grad(sol) + beta*sol - gamma over the same boundary.
76 double IntegrateBC(const GridFunction &sol, const Array<int> &bdr_marker,
77  double alpha, double beta, double gamma,
78  double &err);
79 
80 int main(int argc, char *argv[])
81 {
82  // 1. Parse command-line options.
83  int ser_ref_levels = 2;
84  int order = 1;
85  double sigma = -1.0;
86  double kappa = -1.0;
87  bool h1 = true;
88  bool visualization = true;
89 
90  double mat_val = 1.0;
91  double dbc_val = 0.0;
92  double nbc_val = 1.0;
93  double rbc_a_val = 1.0; // du/dn + a * u = b
94  double rbc_b_val = 1.0;
95 
96  OptionsParser args(argc, argv);
97  args.AddOption(&h1, "-h1", "--continuous", "-dg", "--discontinuous",
98  "Select continuous \"H1\" or discontinuous \"DG\" basis.");
99  args.AddOption(&order, "-o", "--order",
100  "Finite element order (polynomial degree) or -1 for"
101  " isoparametric space.");
102  args.AddOption(&sigma, "-s", "--sigma",
103  "One of the two DG penalty parameters, typically +1/-1."
104  " See the documentation of class DGDiffusionIntegrator.");
105  args.AddOption(&kappa, "-k", "--kappa",
106  "One of the two DG penalty parameters, should be positive."
107  " Negative values are replaced with (order+1)^2.");
108  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
109  "Number of times to refine the mesh uniformly in serial.");
110  args.AddOption(&mat_val, "-mat", "--material-value",
111  "Constant value for material coefficient "
112  "in the Laplace operator.");
113  args.AddOption(&dbc_val, "-dbc", "--dirichlet-value",
114  "Constant value for Dirichlet Boundary Condition.");
115  args.AddOption(&nbc_val, "-nbc", "--neumann-value",
116  "Constant value for Neumann Boundary Condition.");
117  args.AddOption(&rbc_a_val, "-rbc-a", "--robin-a-value",
118  "Constant 'a' value for Robin Boundary Condition: "
119  "du/dn + a * u = b.");
120  args.AddOption(&rbc_b_val, "-rbc-b", "--robin-b-value",
121  "Constant 'b' value for Robin Boundary Condition: "
122  "du/dn + a * u = b.");
123  args.AddOption(&a_, "-a", "--radius",
124  "Radius of holes in the mesh.");
125  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
126  "--no-visualization",
127  "Enable or disable GLVis visualization.");
128  args.Parse();
129  if (!args.Good())
130  {
131  args.PrintUsage(mfem::out);
132  return 1;
133  }
134  if (kappa < 0 && !h1)
135  {
136  kappa = (order+1)*(order+1);
137  }
138  args.PrintOptions(mfem::out);
139 
140  if (a_ < 0.01)
141  {
142  mfem::out << "Hole radius too small, resetting to 0.01.\n";
143  a_ = 0.01;
144  }
145  if (a_ > 0.49)
146  {
147  mfem::out << "Hole radius too large, resetting to 0.49.\n";
148  a_ = 0.49;
149  }
150 
151  // 2. Construct the (serial) mesh and refine it if requested.
152  Mesh *mesh = GenerateSerialMesh(ser_ref_levels);
153  int dim = mesh->Dimension();
154 
155  // 3. Define a finite element space on the serial mesh. Here we use either
156  // continuous Lagrange finite elements or discontinuous Galerkin finite
157  // elements of the specified order.
159  h1 ? (FiniteElementCollection*)new H1_FECollection(order, dim) :
160  (FiniteElementCollection*)new DG_FECollection(order, dim);
161  FiniteElementSpace fespace(mesh, fec);
162  int size = fespace.GetTrueVSize();
163  mfem::out << "Number of finite element unknowns: " << size << endl;
164 
165  // 4. Create "marker arrays" to define the portions of boundary associated
166  // with each type of boundary condition. These arrays have an entry
167  // corresponding to each boundary attribute. Placing a '1' in entry i
168  // marks attribute i+1 as being active, '0' is inactive.
169  Array<int> nbc_bdr(mesh->bdr_attributes.Max());
170  Array<int> rbc_bdr(mesh->bdr_attributes.Max());
171  Array<int> dbc_bdr(mesh->bdr_attributes.Max());
172 
173  nbc_bdr = 0; nbc_bdr[0] = 1;
174  rbc_bdr = 0; rbc_bdr[1] = 1;
175  dbc_bdr = 0; dbc_bdr[2] = 1;
176 
177  Array<int> ess_tdof_list(0);
178  if (h1 && mesh->bdr_attributes.Size())
179  {
180  // For a continuous basis the linear system must be modified to enforce an
181  // essential (Dirichlet) boundary condition. In the DG case this is not
182  // necessary as the boundary condition will only be enforced weakly.
183  fespace.GetEssentialTrueDofs(dbc_bdr, ess_tdof_list);
184  }
185 
186  // 5. Setup the various coefficients needed for the Laplace operator and the
187  // various boundary conditions. In general these coefficients could be
188  // functions of position but here we use only constants.
189  ConstantCoefficient matCoef(mat_val);
190  ConstantCoefficient dbcCoef(dbc_val);
191  ConstantCoefficient nbcCoef(nbc_val);
192  ConstantCoefficient rbcACoef(rbc_a_val);
193  ConstantCoefficient rbcBCoef(rbc_b_val);
194 
195  // Since the n.Grad(u) terms arise by integrating -Div(m Grad(u)) by parts we
196  // must introduce the coefficient 'm' into the boundary conditions.
197  // Therefore, in the case of the Neumann BC, we actually enforce m n.Grad(u)
198  // = m g rather than simply n.Grad(u) = g.
199  ProductCoefficient m_nbcCoef(matCoef, nbcCoef);
200  ProductCoefficient m_rbcACoef(matCoef, rbcACoef);
201  ProductCoefficient m_rbcBCoef(matCoef, rbcBCoef);
202 
203  // 6. Define the solution vector u as a finite element grid function
204  // corresponding to fespace. Initialize u with initial guess of zero.
205  GridFunction u(&fespace);
206  u = 0.0;
207 
208  // 7. Set up the bilinear form a(.,.) on the finite element space
209  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
210  // domain integrator.
211  BilinearForm a(&fespace);
213  if (h1)
214  {
215  // Add a Mass integrator on the Robin boundary
216  a.AddBoundaryIntegrator(new MassIntegrator(m_rbcACoef), rbc_bdr);
217  }
218  else
219  {
220  // Add the interfacial portion of the Laplace operator
222  sigma, kappa));
223 
224  // Counteract the n.Grad(u) term on the Dirichlet portion of the boundary
225  a.AddBdrFaceIntegrator(new DGDiffusionIntegrator(matCoef, sigma, kappa),
226  dbc_bdr);
227 
228  // Augment the n.Grad(u) term with a*u on the Robin portion of boundary
230  rbc_bdr);
231  }
232  a.Assemble();
233 
234  // 8. Assemble the linear form for the right hand side vector.
235  LinearForm b(&fespace);
236 
237  if (h1)
238  {
239  // Set the Dirichlet values in the solution vector
240  u.ProjectBdrCoefficient(dbcCoef, dbc_bdr);
241 
242  // Add the desired value for n.Grad(u) on the Neumann boundary
243  b.AddBoundaryIntegrator(new BoundaryLFIntegrator(m_nbcCoef), nbc_bdr);
244 
245  // Add the desired value for n.Grad(u) + a*u on the Robin boundary
246  b.AddBoundaryIntegrator(new BoundaryLFIntegrator(m_rbcBCoef), rbc_bdr);
247  }
248  else
249  {
250  // Add the desired value for the Dirichlet boundary
251  b.AddBdrFaceIntegrator(new DGDirichletLFIntegrator(dbcCoef, matCoef,
252  sigma, kappa),
253  dbc_bdr);
254 
255  // Add the desired value for n.Grad(u) on the Neumann boundary
257  nbc_bdr);
258 
259  // Add the desired value for n.Grad(u) + a*u on the Robin boundary
260  b.AddBdrFaceIntegrator(new BoundaryLFIntegrator(m_rbcBCoef),
261  rbc_bdr);
262  }
263  b.Assemble();
264 
265  // 9. Construct the linear system.
266  OperatorPtr A;
267  Vector B, X;
268  a.FormLinearSystem(ess_tdof_list, u, b, A, X, B);
269 
270 #ifndef MFEM_USE_SUITESPARSE
271  // 10. Define a simple symmetric Gauss-Seidel preconditioner and use it to
272  // solve the system AX=B with PCG in the symmetric case, and GMRES in the
273  // non-symmetric one.
274  {
275  GSSmoother M((SparseMatrix&)(*A));
276  if (sigma == -1.0)
277  {
278  PCG(*A, M, B, X, 1, 500, 1e-12, 0.0);
279  }
280  else
281  {
282  GMRES(*A, M, B, X, 1, 500, 10, 1e-12, 0.0);
283  }
284  }
285 #else
286  // 11. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the
287  // system.
288  UMFPackSolver umf_solver;
289  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
290  umf_solver.SetOperator(*A);
291  umf_solver.Mult(B, X);
292 #endif
293 
294  // 12. Recover the grid function corresponding to U. This is the local finite
295  // element solution.
296  a.RecoverFEMSolution(X, b, u);
297 
298  // 13. Build a mass matrix to help solve for n.Grad(u) where 'n' is a surface
299  // normal.
300  BilinearForm m(&fespace);
302  m.Assemble();
303 
304  ess_tdof_list.SetSize(0);
305  OperatorPtr M;
306  m.FormSystemMatrix(ess_tdof_list, M);
307 
308  // 14. Compute the various boundary integrals.
309  mfem::out << endl
310  << "Verifying boundary conditions" << endl
311  << "=============================" << endl;
312  {
313  // Integrate the solution on the Dirichlet boundary and compare to the
314  // expected value.
315  double err, avg = IntegrateBC(u, dbc_bdr, 0.0, 1.0, dbc_val, err);
316 
317  bool hom_dbc = (dbc_val == 0.0);
318  err /= hom_dbc ? 1.0 : fabs(dbc_val);
319  mfem::out << "Average of solution on Gamma_dbc:\t"
320  << avg << ", \t"
321  << (hom_dbc ? "absolute" : "relative")
322  << " error " << err << endl;
323  }
324  {
325  // Integrate n.Grad(u) on the inhomogeneous Neumann boundary and compare
326  // to the expected value.
327  double err, avg = IntegrateBC(u, nbc_bdr, 1.0, 0.0, nbc_val, err);
328 
329  bool hom_nbc = (nbc_val == 0.0);
330  err /= hom_nbc ? 1.0 : fabs(nbc_val);
331  mfem::out << "Average of n.Grad(u) on Gamma_nbc:\t"
332  << avg << ", \t"
333  << (hom_nbc ? "absolute" : "relative")
334  << " error " << err << endl;
335  }
336  {
337  // Integrate n.Grad(u) on the homogeneous Neumann boundary and compare to
338  // the expected value of zero.
339  Array<int> nbc0_bdr(mesh->bdr_attributes.Max());
340  nbc0_bdr = 0;
341  nbc0_bdr[3] = 1;
342 
343  double err, avg = IntegrateBC(u, nbc0_bdr, 1.0, 0.0, 0.0, err);
344 
345  bool hom_nbc = true;
346  mfem::out << "Average of n.Grad(u) on Gamma_nbc0:\t"
347  << avg << ", \t"
348  << (hom_nbc ? "absolute" : "relative")
349  << " error " << err << endl;
350  }
351  {
352  // Integrate n.Grad(u) + a * u on the Robin boundary and compare to the
353  // expected value.
354  double err, avg = IntegrateBC(u, rbc_bdr, 1.0, rbc_a_val, rbc_b_val, err);
355 
356  bool hom_rbc = (rbc_b_val == 0.0);
357  err /= hom_rbc ? 1.0 : fabs(rbc_b_val);
358  mfem::out << "Average of n.Grad(u)+a*u on Gamma_rbc:\t"
359  << avg << ", \t"
360  << (hom_rbc ? "absolute" : "relative")
361  << " error " << err << endl;
362  }
363 
364  // 15. Save the refined mesh and the solution. This output can be viewed
365  // later using GLVis: "glvis -m refined.mesh -g sol.gf".
366  {
367  ofstream mesh_ofs("refined.mesh");
368  mesh_ofs.precision(8);
369  mesh->Print(mesh_ofs);
370  ofstream sol_ofs("sol.gf");
371  sol_ofs.precision(8);
372  u.Save(sol_ofs);
373  }
374 
375  // 16. Send the solution by socket to a GLVis server.
376  if (visualization)
377  {
378  string title_str = h1 ? "H1" : "DG";
379  char vishost[] = "localhost";
380  int visport = 19916;
381  socketstream sol_sock(vishost, visport);
382  sol_sock.precision(8);
383  sol_sock << "solution\n" << *mesh << u
384  << "window_title '" << title_str << " Solution'"
385  << " keys 'mmc'" << flush;
386  }
387 
388  // 17. Free the used memory.
389  delete fec;
390  delete mesh;
391 
392  return 0;
393 }
394 
395 void quad_trans(double u, double v, double &x, double &y, bool log = false)
396 {
397  double a = a_; // Radius of disc
398 
399  double d = 4.0 * a * (M_SQRT2 - 2.0 * a) * (1.0 - 2.0 * v);
400 
401  double v0 = (1.0 + M_SQRT2) * (M_SQRT2 * a - 2.0 * v) *
402  ((4.0 - 3 * M_SQRT2) * a +
403  (8.0 * (M_SQRT2 - 1.0) * a - 2.0) * v) / d;
404 
405  double r = 2.0 * ((M_SQRT2 - 1.0) * a * a * (1.0 - 4.0 *v) +
406  2.0 * (1.0 + M_SQRT2 *
407  (1.0 + 2.0 * (2.0 * a - M_SQRT2 - 1.0) * a)) * v * v
408  ) / d;
409 
410  double t = asin(v / r) * u / v;
411  if (log)
412  {
413  mfem::out << "u, v, r, v0, t "
414  << u << " " << v << " " << r << " " << v0 << " " << t
415  << endl;
416  }
417  x = r * sin(t);
418  y = r * cos(t) - v0;
419 }
420 
421 void trans(const Vector &u, Vector &x)
422 {
423  double tol = 1e-4;
424 
425  if (u[1] > 0.5 - tol || u[1] < -0.5 + tol)
426  {
427  x = u;
428  return;
429  }
430  if (u[0] > 1.0 - tol || u[0] < -1.0 + tol || fabs(u[0]) < tol)
431  {
432  x = u;
433  return;
434  }
435 
436  if (u[0] > 0.0)
437  {
438  if (u[1] > fabs(u[0] - 0.5))
439  {
440  quad_trans(u[0] - 0.5, u[1], x[0], x[1]);
441  x[0] += 0.5;
442  return;
443  }
444  if (u[1] < -fabs(u[0] - 0.5))
445  {
446  quad_trans(u[0] - 0.5, -u[1], x[0], x[1]);
447  x[0] += 0.5;
448  x[1] *= -1.0;
449  return;
450  }
451  if (u[0] - 0.5 > fabs(u[1]))
452  {
453  quad_trans(u[1], u[0] - 0.5, x[1], x[0]);
454  x[0] += 0.5;
455  return;
456  }
457  if (u[0] - 0.5 < -fabs(u[1]))
458  {
459  quad_trans(u[1], 0.5 - u[0], x[1], x[0]);
460  x[0] *= -1.0;
461  x[0] += 0.5;
462  return;
463  }
464  }
465  else
466  {
467  if (u[1] > fabs(u[0] + 0.5))
468  {
469  quad_trans(u[0] + 0.5, u[1], x[0], x[1]);
470  x[0] -= 0.5;
471  return;
472  }
473  if (u[1] < -fabs(u[0] + 0.5))
474  {
475  quad_trans(u[0] + 0.5, -u[1], x[0], x[1]);
476  x[0] -= 0.5;
477  x[1] *= -1.0;
478  return;
479  }
480  if (u[0] + 0.5 > fabs(u[1]))
481  {
482  quad_trans(u[1], u[0] + 0.5, x[1], x[0]);
483  x[0] -= 0.5;
484  return;
485  }
486  if (u[0] + 0.5 < -fabs(u[1]))
487  {
488  quad_trans(u[1], -0.5 - u[0], x[1], x[0]);
489  x[0] *= -1.0;
490  x[0] -= 0.5;
491  return;
492  }
493  }
494  x = u;
495 }
496 
498 {
499  Mesh * mesh = new Mesh(2, 29, 16, 24, 2);
500 
501  int vi[4];
502 
503  for (int i=0; i<2; i++)
504  {
505  int o = 13 * i;
506  vi[0] = o + 0; vi[1] = o + 3; vi[2] = o + 4; vi[3] = o + 1;
507  mesh->AddQuad(vi);
508 
509  vi[0] = o + 1; vi[1] = o + 4; vi[2] = o + 5; vi[3] = o + 2;
510  mesh->AddQuad(vi);
511 
512  vi[0] = o + 5; vi[1] = o + 8; vi[2] = o + 9; vi[3] = o + 2;
513  mesh->AddQuad(vi);
514 
515  vi[0] = o + 8; vi[1] = o + 12; vi[2] = o + 15; vi[3] = o + 9;
516  mesh->AddQuad(vi);
517 
518  vi[0] = o + 11; vi[1] = o + 14; vi[2] = o + 15; vi[3] = o + 12;
519  mesh->AddQuad(vi);
520 
521  vi[0] = o + 10; vi[1] = o + 13; vi[2] = o + 14; vi[3] = o + 11;
522  mesh->AddQuad(vi);
523 
524  vi[0] = o + 6; vi[1] = o + 13; vi[2] = o + 10; vi[3] = o + 7;
525  mesh->AddQuad(vi);
526 
527  vi[0] = o + 0; vi[1] = o + 6; vi[2] = o + 7; vi[3] = o + 3;
528  mesh->AddQuad(vi);
529  }
530 
531  vi[0] = 0; vi[1] = 6; mesh->AddBdrSegment(vi, 1);
532  vi[0] = 6; vi[1] = 13; mesh->AddBdrSegment(vi, 1);
533  vi[0] = 13; vi[1] = 19; mesh->AddBdrSegment(vi, 1);
534  vi[0] = 19; vi[1] = 26; mesh->AddBdrSegment(vi, 1);
535 
536  vi[0] = 28; vi[1] = 22; mesh->AddBdrSegment(vi, 2);
537  vi[0] = 22; vi[1] = 15; mesh->AddBdrSegment(vi, 2);
538  vi[0] = 15; vi[1] = 9; mesh->AddBdrSegment(vi, 2);
539  vi[0] = 9; vi[1] = 2; mesh->AddBdrSegment(vi, 2);
540 
541  for (int i=0; i<2; i++)
542  {
543  int o = 13 * i;
544  vi[0] = o + 7; vi[1] = o + 3; mesh->AddBdrSegment(vi, 3 + i);
545  vi[0] = o + 10; vi[1] = o + 7; mesh->AddBdrSegment(vi, 3 + i);
546  vi[0] = o + 11; vi[1] = o + 10; mesh->AddBdrSegment(vi, 3 + i);
547  vi[0] = o + 12; vi[1] = o + 11; mesh->AddBdrSegment(vi, 3 + i);
548  vi[0] = o + 8; vi[1] = o + 12; mesh->AddBdrSegment(vi, 3 + i);
549  vi[0] = o + 5; vi[1] = o + 8; mesh->AddBdrSegment(vi, 3 + i);
550  vi[0] = o + 4; vi[1] = o + 5; mesh->AddBdrSegment(vi, 3 + i);
551  vi[0] = o + 3; vi[1] = o + 4; mesh->AddBdrSegment(vi, 3 + i);
552  }
553 
554  double d[2];
555  double a = a_ / M_SQRT2;
556 
557  d[0] = -1.0; d[1] = -0.5; mesh->AddVertex(d);
558  d[0] = -1.0; d[1] = 0.0; mesh->AddVertex(d);
559  d[0] = -1.0; d[1] = 0.5; mesh->AddVertex(d);
560 
561  d[0] = -0.5 - a; d[1] = -a; mesh->AddVertex(d);
562  d[0] = -0.5 - a; d[1] = 0.0; mesh->AddVertex(d);
563  d[0] = -0.5 - a; d[1] = a; mesh->AddVertex(d);
564 
565  d[0] = -0.5; d[1] = -0.5; mesh->AddVertex(d);
566  d[0] = -0.5; d[1] = -a; mesh->AddVertex(d);
567  d[0] = -0.5; d[1] = a; mesh->AddVertex(d);
568  d[0] = -0.5; d[1] = 0.5; mesh->AddVertex(d);
569 
570  d[0] = -0.5 + a; d[1] = -a; mesh->AddVertex(d);
571  d[0] = -0.5 + a; d[1] = 0.0; mesh->AddVertex(d);
572  d[0] = -0.5 + a; d[1] = a; mesh->AddVertex(d);
573 
574  d[0] = 0.0; d[1] = -0.5; mesh->AddVertex(d);
575  d[0] = 0.0; d[1] = 0.0; mesh->AddVertex(d);
576  d[0] = 0.0; d[1] = 0.5; mesh->AddVertex(d);
577 
578  d[0] = 0.5 - a; d[1] = -a; mesh->AddVertex(d);
579  d[0] = 0.5 - a; d[1] = 0.0; mesh->AddVertex(d);
580  d[0] = 0.5 - a; d[1] = a; mesh->AddVertex(d);
581 
582  d[0] = 0.5; d[1] = -0.5; mesh->AddVertex(d);
583  d[0] = 0.5; d[1] = -a; mesh->AddVertex(d);
584  d[0] = 0.5; d[1] = a; mesh->AddVertex(d);
585  d[0] = 0.5; d[1] = 0.5; mesh->AddVertex(d);
586 
587  d[0] = 0.5 + a; d[1] = -a; mesh->AddVertex(d);
588  d[0] = 0.5 + a; d[1] = 0.0; mesh->AddVertex(d);
589  d[0] = 0.5 + a; d[1] = a; mesh->AddVertex(d);
590 
591  d[0] = 1.0; d[1] = -0.5; mesh->AddVertex(d);
592  d[0] = 1.0; d[1] = 0.0; mesh->AddVertex(d);
593  d[0] = 1.0; d[1] = 0.5; mesh->AddVertex(d);
594 
595  mesh->FinalizeTopology();
596 
597  mesh->SetCurvature(1, true);
598 
599  // Stitch the ends of the stack together
600  {
601  Array<int> v2v(mesh->GetNV());
602  for (int i = 0; i < v2v.Size() - 3; i++)
603  {
604  v2v[i] = i;
605  }
606  // identify vertices on the narrow ends of the rectangle
607  v2v[v2v.Size() - 3] = 0;
608  v2v[v2v.Size() - 2] = 1;
609  v2v[v2v.Size() - 1] = 2;
610 
611  // renumber elements
612  for (int i = 0; i < mesh->GetNE(); i++)
613  {
614  Element *el = mesh->GetElement(i);
615  int *v = el->GetVertices();
616  int nv = el->GetNVertices();
617  for (int j = 0; j < nv; j++)
618  {
619  v[j] = v2v[v[j]];
620  }
621  }
622  // renumber boundary elements
623  for (int i = 0; i < mesh->GetNBE(); i++)
624  {
625  Element *el = mesh->GetBdrElement(i);
626  int *v = el->GetVertices();
627  int nv = el->GetNVertices();
628  for (int j = 0; j < nv; j++)
629  {
630  v[j] = v2v[v[j]];
631  }
632  }
633  mesh->RemoveUnusedVertices();
634  mesh->RemoveInternalBoundaries();
635  }
636  mesh->SetCurvature(3, true);
637 
638  for (int l = 0; l < ref; l++)
639  {
640  mesh->UniformRefinement();
641  }
642 
643  mesh->Transform(trans);
644 
645  return mesh;
646 }
647 
648 double IntegrateBC(const GridFunction &x, const Array<int> &bdr,
649  double alpha, double beta, double gamma,
650  double &err)
651 {
652  double nrm = 0.0;
653  double avg = 0.0;
654  err = 0.0;
655 
656  const bool a_is_zero = alpha == 0.0;
657  const bool b_is_zero = beta == 0.0;
658 
659  const FiniteElementSpace &fes = *x.FESpace();
660  MFEM_ASSERT(fes.GetVDim() == 1, "");
661  Mesh &mesh = *fes.GetMesh();
662  Vector shape, loc_dofs, w_nor;
663  DenseMatrix dshape;
664  Array<int> dof_ids;
665  for (int i = 0; i < mesh.GetNBE(); i++)
666  {
667  if (bdr[mesh.GetBdrAttribute(i)-1] == 0) { continue; }
668 
669  FaceElementTransformations *FTr = mesh.GetBdrFaceTransformations(i);
670  if (FTr == nullptr) { continue; }
671 
672  const FiniteElement &fe = *fes.GetFE(FTr->Elem1No);
673  MFEM_ASSERT(fe.GetMapType() == FiniteElement::VALUE, "");
674  const int int_order = 2*fe.GetOrder() + 3;
675  const IntegrationRule &ir = IntRules.Get(FTr->FaceGeom, int_order);
676 
677  fes.GetElementDofs(FTr->Elem1No, dof_ids);
678  x.GetSubVector(dof_ids, loc_dofs);
679  if (!a_is_zero)
680  {
681  const int sdim = FTr->Face->GetSpaceDim();
682  w_nor.SetSize(sdim);
683  dshape.SetSize(fe.GetDof(), sdim);
684  }
685  if (!b_is_zero)
686  {
687  shape.SetSize(fe.GetDof());
688  }
689  for (int j = 0; j < ir.GetNPoints(); j++)
690  {
691  const IntegrationPoint &ip = ir.IntPoint(j);
692  IntegrationPoint eip;
693  FTr->Loc1.Transform(ip, eip);
694  FTr->Face->SetIntPoint(&ip);
695  double face_weight = FTr->Face->Weight();
696  double val = 0.0;
697  if (!a_is_zero)
698  {
699  FTr->Elem1->SetIntPoint(&eip);
700  fe.CalcPhysDShape(*FTr->Elem1, dshape);
701  CalcOrtho(FTr->Face->Jacobian(), w_nor);
702  val += alpha * dshape.InnerProduct(w_nor, loc_dofs) / face_weight;
703  }
704  if (!b_is_zero)
705  {
706  fe.CalcShape(eip, shape);
707  val += beta * (shape * loc_dofs);
708  }
709 
710  // Measure the length of the boundary
711  nrm += ip.weight * face_weight;
712 
713  // Integrate alpha * n.Grad(x) + beta * x
714  avg += val * ip.weight * face_weight;
715 
716  // Integrate |alpha * n.Grad(x) + beta * x - gamma|^2
717  val -= gamma;
718  err += (val*val) * ip.weight * face_weight;
719  }
720  }
721 
722  // Normalize by the length of the boundary
723  if (std::abs(nrm) > 0.0)
724  {
725  err /= nrm;
726  avg /= nrm;
727  }
728 
729  // Compute l2 norm of the error in the boundary condition (negative
730  // quadrature weights may produce negative 'err')
731  err = (err >= 0.0) ? sqrt(err) : -sqrt(-err);
732 
733  // Return the average value of alpha * n.Grad(x) + beta * x
734  return avg;
735 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:245
Abstract class for all finite elements.
Definition: fe.hpp:235
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1234
void trans(const Vector &u, Vector &x)
Definition: ex27.cpp:421
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
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:1293
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:915
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
virtual void GetVertices(Array< int > &v) const =0
Returns element&#39;s vertices.
void n4Vec(const Vector &x, Vector &n)
Definition: ex27.cpp:69
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:459
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:740
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:744
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
Scalar coefficient defined as the product of two scalar coefficients or a scalar and a scalar coeffic...
double InnerProduct(const double *x, const double *y) const
Compute y^t A x.
Definition: densemat.cpp:299
void quad_trans(double u, double v, double &x, double &y, bool log=false)
Definition: ex27.cpp:395
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:85
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:467
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition: vector.cpp:495
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 RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Recover the solution of a linear system formed with FormLinearSystem().
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe.hpp:319
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Class for boundary integration L(v) := (g, v)
Definition: lininteg.hpp:150
void Transform(void(*f)(const Vector &, Vector &))
Definition: mesh.cpp:10292
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Definition: fespace.cpp:415
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
double kappa
Definition: ex24.cpp:54
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
Definition: fe.hpp:341
int main(int argc, char *argv[])
Definition: ex1.cpp:66
void CalcOrtho(const DenseMatrix &J, Vector &n)
Definition: densemat.cpp:2305
void RemoveInternalBoundaries()
Definition: mesh.cpp:10448
Data type for Gauss-Seidel smoother of sparse matrix.
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:731
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3417
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:248
int AddVertex(double x, double y=0.0, double z=0.0)
Definition: mesh.cpp:1232
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:46
constexpr char vishost[]
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:314
IntegrationPointTransformation Loc1
Definition: eltrans.hpp:511
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:8382
constexpr int visport
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:111
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:4165
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
virtual void GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:1696
int AddBdrSegment(int v1, int v2, int attr=1)
Definition: mesh.cpp:1409
double Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:123
const Element * GetElement(int i) const
Definition: mesh.hpp:819
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
Definition: solvers.cpp:733
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:66
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:403
void CalcPhysDShape(ElementTransformation &Trans, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in physical space at the poi...
Definition: fe.cpp:201
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:389
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:595
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:53
int Dimension() const
Definition: mesh.hpp:788
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:434
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
Definition: fe_coll.hpp:269
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:201
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:742
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:315
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
void RemoveUnusedVertices()
Remove unused vertices and rebuild mesh connectivity.
Definition: mesh.cpp:10341
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
void FinalizeTopology(bool generate_bdr=true)
Finalize the construction of the secondary topology (connectivity) data of a Mesh.
Definition: mesh.cpp:2469
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
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:654
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition: eltrans.cpp:533
virtual int GetSpaceDim() const =0
Get the dimension of the target (physical) space.
double a
Definition: lissajous.cpp:41
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:734
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
Class for integration point with weight.
Definition: intrules.hpp:25
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi)
Adds new Boundary Integrator. Assumes ownership of bfi.
int dim
Definition: ex24.cpp:53
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.
Mesh * GenerateSerialMesh(int ref)
Definition: ex27.cpp:497
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:1798
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:304
ElementTransformation * Elem1
Definition: eltrans.hpp:509
ElementTransformation * Face
Definition: eltrans.hpp:510
virtual int GetNVertices() const =0
const double alpha
Definition: ex15.cpp:336
Vector data type.
Definition: vector.hpp:51
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:159
Vector with associated FE space and LinearFormIntegrators.
Definition: linearform.hpp:23
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:90
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:378
double IntegrateBC(const GridFunction &sol, const Array< int > &bdr_marker, double alpha, double beta, double gamma, double &err)
Definition: ex27.cpp:648
Abstract data type element.
Definition: element.hpp:28
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:2818
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:823
int GMRES(const Operator &A, Vector &x, const Vector &b, Solver &M, int &max_iter, int m, double &tol, double atol, int printit)
GMRES method. (tolerances are squared)
Definition: solvers.cpp:1156
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator. Assumes ownership of bfi.
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the bo...
Definition: gridfunc.hpp:399
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:2721
double sigma(const Vector &x)
Definition: maxwell.cpp:102
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:145