MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
mesh-explorer.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 //
12 // -----------------------------------------------------
13 // Mesh Explorer Miniapp: Explore and manipulate meshes
14 // -----------------------------------------------------
15 //
16 // This miniapp is a handy tool to examine, visualize and manipulate a given
17 // mesh. Some of its features are:
18 //
19 // - visualizing of mesh materials and individual mesh elements
20 // - mesh scaling, randomization, and general transformation
21 // - manipulation of the mesh curvature
22 // - the ability to simulate parallel partitioning
23 // - quantitative and visual reports of mesh quality
24 //
25 // Compile with: make mesh-explorer
26 //
27 // Sample runs: mesh-explorer
28 // mesh-explorer -m ../../data/beam-tri.mesh
29 // mesh-explorer -m ../../data/star-q2.mesh
30 // mesh-explorer -m ../../data/disc-nurbs.mesh
31 // mesh-explorer -m ../../data/escher-p3.mesh
32 // mesh-explorer -m ../../data/mobius-strip.mesh
33 
34 #include "mfem.hpp"
35 #include <fstream>
36 #include <limits>
37 #include <cstdlib>
38 
39 using namespace mfem;
40 using namespace std;
41 
42 // This transformation can be applied to a mesh with the 't' menu option.
43 void transformation(const Vector &p, Vector &v)
44 {
45  // simple shear transformation
46  double s = 0.1;
47 
48  if (p.Size() == 3)
49  {
50  v(0) = p(0) + s*p(1) + s*p(2);
51  v(1) = p(1) + s*p(2) + s*p(0);
52  v(2) = p(2);
53  }
54  else if (p.Size() == 2)
55  {
56  v(0) = p(0) + s*p(1);
57  v(1) = p(1) + s*p(0);
58  }
59  else
60  {
61  v = p;
62  }
63 }
64 
65 // This function is used with the 'r' menu option, sub-option 'l' to refine a
66 // mesh locally in a region, defined by return values <= region_eps.
67 double region_eps = 1e-8;
68 double region(const Vector &p)
69 {
70  const double x = p(0), y = p(1);
71  // here we describe the region: (x <= 1/4) && (y >= 0) && (y <= 1)
72  return std::max(std::max(x - 0.25, -y), y - 1.0);
73 }
74 
75 // The projection of this function can be plotted with the 'l' menu option
76 double f(const Vector &p)
77 {
78  double x = p(0);
79  double y = p.Size() > 1 ? p(1) : 0.0;
80  double z = p.Size() > 2 ? p(2) : 0.0;
81 
82  if (1)
83  {
84  // torus in the xy-plane
85  const double r_big = 2.0;
86  const double r_small = 1.0;
87  return hypot(r_big - hypot(x, y), z) - r_small;
88  }
89  if (0)
90  {
91  // sphere at the origin:
92  const double r = 1.0;
93  return hypot(hypot(x, y), z) - r;
94  }
95 }
96 
97 Mesh *read_par_mesh(int np, const char *mesh_prefix, Array<int>& partitioning,
98  Array<int>& bdr_partitioning)
99 {
100  Mesh *mesh;
101  Array<Mesh *> mesh_array;
102 
103  mesh_array.SetSize(np);
104  for (int p = 0; p < np; p++)
105  {
106  ostringstream fname;
107  fname << mesh_prefix << '.' << setfill('0') << setw(6) << p;
108  ifgzstream meshin(fname.str().c_str());
109  if (!meshin)
110  {
111  cerr << "Can not open mesh file: " << fname.str().c_str()
112  << '!' << endl;
113  for (p--; p >= 0; p--)
114  {
115  delete mesh_array[p];
116  }
117  return NULL;
118  }
119  mesh_array[p] = new Mesh(meshin, 1, 0);
120  // Assign corresponding processor number to element + boundary partitions
121  for (int i = 0; i < mesh_array[p]->GetNE(); i++)
122  {
123  partitioning.Append(p);
124  }
125  for (int i = 0; i < mesh_array[p]->GetNBE(); i++)
126  {
127  bdr_partitioning.Append(p);
128  }
129  }
130  mesh = new Mesh(mesh_array, np);
131 
132  for (int p = 0; p < np; p++)
133  {
134  delete mesh_array[np-1-p];
135  }
136  mesh_array.DeleteAll();
137 
138  return mesh;
139 }
140 
141 // Given a 3D mesh, produce a 2D mesh consisting of its boundary elements.
142 // We guarantee that the skin preserves the boundary index order.
144 {
145  // Determine mapping from vertex to boundary vertex
146  Array<int> v2v(mesh->GetNV());
147  v2v = -1;
148  for (int i = 0; i < mesh->GetNBE(); i++)
149  {
150  Element *el = mesh->GetBdrElement(i);
151  int *v = el->GetVertices();
152  int nv = el->GetNVertices();
153  for (int j = 0; j < nv; j++)
154  {
155  v2v[v[j]] = 0;
156  }
157  }
158  int nbvt = 0;
159  for (int i = 0; i < v2v.Size(); i++)
160  {
161  if (v2v[i] == 0)
162  {
163  v2v[i] = nbvt++;
164  }
165  }
166 
167  // Create a new mesh for the boundary
168  Mesh * bmesh = new Mesh(mesh->Dimension() - 1, nbvt, mesh->GetNBE(),
169  0, mesh->SpaceDimension());
170 
171  // Copy vertices to the boundary mesh
172  nbvt = 0;
173  for (int i = 0; i < v2v.Size(); i++)
174  {
175  if (v2v[i] >= 0)
176  {
177  double *c = mesh->GetVertex(i);
178  bmesh->AddVertex(c);
179  nbvt++;
180  }
181  }
182 
183  // Copy elements to the boundary mesh
184  int bv[4];
185  for (int i = 0; i < mesh->GetNBE(); i++)
186  {
187  Element *el = mesh->GetBdrElement(i);
188  int *v = el->GetVertices();
189  int nv = el->GetNVertices();
190 
191  for (int j = 0; j < nv; j++)
192  {
193  bv[j] = v2v[v[j]];
194  }
195 
196  switch (el->GetGeometryType())
197  {
198  case Geometry::SEGMENT:
199  bmesh->AddSegment(bv, el->GetAttribute());
200  break;
201  case Geometry::TRIANGLE:
202  bmesh->AddTriangle(bv, el->GetAttribute());
203  break;
204  case Geometry::SQUARE:
205  bmesh->AddQuad(bv, el->GetAttribute());
206  break;
207  default:
208  break; // This should not happen
209  }
210 
211  }
212  bmesh->FinalizeTopology();
213 
214  // Copy GridFunction describing nodes if present
215  if (mesh->GetNodes())
216  {
217  FiniteElementSpace *fes = mesh->GetNodes()->FESpace();
218  const FiniteElementCollection *fec = fes->FEColl();
219  if (dynamic_cast<const H1_FECollection*>(fec))
220  {
221  FiniteElementCollection *fec_copy =
223  FiniteElementSpace *fes_copy =
224  new FiniteElementSpace(*fes, bmesh, fec_copy);
225  GridFunction *bdr_nodes = new GridFunction(fes_copy);
226  bdr_nodes->MakeOwner(fec_copy);
227 
228  bmesh->NewNodes(*bdr_nodes, true);
229 
230  Array<int> vdofs;
231  Array<int> bvdofs;
232  Vector v;
233  for (int i=0; i<mesh->GetNBE(); i++)
234  {
235  fes->GetBdrElementVDofs(i, vdofs);
236  mesh->GetNodes()->GetSubVector(vdofs, v);
237 
238  fes_copy->GetElementVDofs(i, bvdofs);
239  bdr_nodes->SetSubVector(bvdofs, v);
240  }
241  }
242  else
243  {
244  cout << "\nDiscontinuous nodes not yet supported" << endl;
245  }
246  }
247 
248  return bmesh;
249 }
250 
251 void recover_bdr_partitioning(const Mesh* mesh, const Array<int>& partitioning,
252  Array<int>& bdr_partitioning)
253 {
254  bdr_partitioning.SetSize(mesh->GetNBE());
255  int info, e;
256  for (int be = 0; be < mesh->GetNBE(); be++)
257  {
258  mesh->GetBdrElementAdjacentElement(be, e, info);
259  bdr_partitioning[be] = partitioning[e];
260  }
261 }
262 
263 int main (int argc, char *argv[])
264 {
265  int np = 0;
266  const char *mesh_file = "../../data/beam-hex.mesh";
267  bool refine = true;
268 
269  OptionsParser args(argc, argv);
270  args.AddOption(&mesh_file, "-m", "--mesh",
271  "Mesh file to visualize.");
272  args.AddOption(&np, "-np", "--num-proc",
273  "Load mesh from multiple processors.");
274  args.AddOption(&refine, "-ref", "--refinement", "-no-ref", "--no-refinement",
275  "Prepare the mesh for refinement or not.");
276  args.Parse();
277  if (!args.Good())
278  {
279  if (!args.Help())
280  {
281  args.PrintError(cout);
282  cout << endl;
283  }
284  cout << "Visualize and manipulate a serial mesh:\n"
285  << " mesh-explorer -m <mesh_file>\n"
286  << "Visualize and manipulate a parallel mesh:\n"
287  << " mesh-explorer -np <#proc> -m <mesh_prefix>\n" << endl
288  << "All Options:\n";
289  args.PrintHelp(cout);
290  return 1;
291  }
292  args.PrintOptions(cout);
293 
294  Mesh *mesh;
295  Mesh *bdr_mesh = NULL;
296 
297  // Helper to distinguish whether we use a parallel or serial mesh.
298  const bool use_par_mesh = np > 0;
299 
300  // Helper for visualizing the partitioning.
301  Array<int> partitioning;
302  Array<int> bdr_partitioning;
303  if (!use_par_mesh)
304  {
305  mesh = new Mesh(mesh_file, 1, refine);
306  partitioning.SetSize(mesh->GetNE());
307  partitioning = 0;
308  bdr_partitioning.SetSize(mesh->GetNBE());
309  bdr_partitioning = 0;
310  }
311  else
312  {
313  mesh = read_par_mesh(np, mesh_file, partitioning, bdr_partitioning);
314  if (mesh == NULL)
315  {
316  return 3;
317  }
318  }
319  int dim = mesh->Dimension();
320  int sdim = mesh->SpaceDimension();
321 
322  FiniteElementCollection *bdr_attr_fec = NULL;
323  FiniteElementCollection *attr_fec;
324  if (dim == 2)
325  {
326  attr_fec = new Const2DFECollection;
327  }
328  else
329  {
330  bdr_attr_fec = new Const2DFECollection;
331  attr_fec = new Const3DFECollection;
332  }
333 
334  int print_char = 1;
335  while (1)
336  {
337  if (print_char)
338  {
339  cout << endl;
340  mesh->PrintCharacteristics();
341  cout << "boundary attribs :";
342  for (int i = 0; i < mesh->bdr_attributes.Size(); i++)
343  {
344  cout << ' ' << mesh->bdr_attributes[i];
345  }
346  cout << '\n' << "material attribs :";
347  for (int i = 0; i < mesh->attributes.Size(); i++)
348  {
349  cout << ' ' << mesh->attributes[i];
350  }
351  cout << endl;
352  cout << "mesh curvature : ";
353  if (mesh->GetNodalFESpace() != NULL)
354  {
355  cout << mesh->GetNodalFESpace()->FEColl()->Name() << endl;
356  }
357  else
358  {
359  cout << "NONE" << endl;
360  }
361  }
362  print_char = 0;
363  cout << endl;
364  cout << "What would you like to do?\n"
365  "r) Refine\n"
366  "c) Change curvature\n"
367  "s) Scale\n"
368  "t) Transform\n"
369  "j) Jitter\n"
370  "v) View mesh\n"
371  "P) View partitioning\n"
372  "m) View materials\n"
373  "b) View boundary\n"
374  "B) View boundary partitioning\n"
375  "e) View elements\n"
376  "h) View element sizes, h\n"
377  "k) View element ratios, kappa\n"
378  "J) View scaled Jacobian\n"
379  "l) Plot a function\n"
380  "x) Print sub-element stats\n"
381  "f) Find physical point in reference space\n"
382  "p) Generate a partitioning\n"
383  "o) Reorder elements\n"
384  "S) Save in MFEM format\n"
385  "V) Save in VTK format (only linear and quadratic meshes)\n"
386  "q) Quit\n"
387 #ifdef MFEM_USE_ZLIB
388  "Z) Save in MFEM format with compression\n"
389 #endif
390  "--> " << flush;
391  char mk;
392  cin >> mk;
393  if (!cin) { break; }
394 
395  if (mk == 'q')
396  {
397  break;
398  }
399 
400  if (mk == 'r')
401  {
402  cout <<
403  "Choose type of refinement:\n"
404  "s) standard refinement with Mesh::UniformRefinement()\n"
405  "b) Mesh::UniformRefinement() (bisection for tet meshes)\n"
406  "u) uniform refinement with a factor\n"
407  "g) non-uniform refinement (Gauss-Lobatto) with a factor\n"
408  "l) refine locally using the region() function\n"
409  "--> " << flush;
410  char sk;
411  cin >> sk;
412  switch (sk)
413  {
414  case 's':
415  mesh->UniformRefinement();
416  // Make sure tet-only meshes are marked for local refinement.
417  mesh->Finalize(true);
418  break;
419  case 'b':
420  mesh->UniformRefinement(1); // ref_algo = 1
421  break;
422  case 'u':
423  case 'g':
424  {
425  cout << "enter refinement factor --> " << flush;
426  int ref_factor;
427  cin >> ref_factor;
428  if (ref_factor <= 1 || ref_factor > 32) { break; }
429  int ref_type = (sk == 'u') ? BasisType::ClosedUniform :
431  *mesh = Mesh::MakeRefined(*mesh, ref_factor, ref_type);
432  break;
433  }
434  case 'l':
435  {
436  Vector pt;
437  Array<int> marked_elements;
438  for (int i = 0; i < mesh->GetNE(); i++)
439  {
440  // check all nodes of the element
442  mesh->GetElementTransformation(i, &T);
443  for (int j = 0; j < T.GetPointMat().Width(); j++)
444  {
445  T.GetPointMat().GetColumnReference(j, pt);
446  if (region(pt) <= region_eps)
447  {
448  marked_elements.Append(i);
449  break;
450  }
451  }
452  }
453  mesh->GeneralRefinement(marked_elements);
454  break;
455  }
456  }
457  print_char = 1;
458  }
459 
460  if (mk == 'c')
461  {
462  int p;
463  cout << "enter new order for mesh curvature --> " << flush;
464  cin >> p;
465  mesh->SetCurvature(p > 0 ? p : -p, p <= 0);
466  print_char = 1;
467  }
468 
469  if (mk == 's')
470  {
471  double factor;
472  cout << "scaling factor ---> " << flush;
473  cin >> factor;
474 
475  GridFunction *nodes = mesh->GetNodes();
476  if (nodes == NULL)
477  {
478  for (int i = 0; i < mesh->GetNV(); i++)
479  {
480  double *v = mesh->GetVertex(i);
481  v[0] *= factor;
482  v[1] *= factor;
483  if (dim == 3)
484  {
485  v[2] *= factor;
486  }
487  }
488  }
489  else
490  {
491  *nodes *= factor;
492  }
493 
494  print_char = 1;
495  }
496 
497  if (mk == 't')
498  {
499  mesh->Transform(transformation);
500  print_char = 1;
501  }
502 
503  if (mk == 'j')
504  {
505  double jitter;
506  cout << "jitter factor ---> " << flush;
507  cin >> jitter;
508 
509  GridFunction *nodes = mesh->GetNodes();
510 
511  if (nodes == NULL)
512  {
513  cerr << "The mesh should have nodes, introduce curvature first!\n";
514  }
515  else
516  {
517  FiniteElementSpace *fespace = nodes->FESpace();
518 
519  GridFunction rdm(fespace);
520  rdm.Randomize();
521  rdm -= 0.5; // shift to random values in [-0.5,0.5]
522  rdm *= jitter;
523 
524  // compute minimal local mesh size
525  Vector h0(fespace->GetNDofs());
526  h0 = infinity();
527  {
528  Array<int> dofs;
529  for (int i = 0; i < fespace->GetNE(); i++)
530  {
531  fespace->GetElementDofs(i, dofs);
532  for (int j = 0; j < dofs.Size(); j++)
533  {
534  h0(dofs[j]) = std::min(h0(dofs[j]), mesh->GetElementSize(i));
535  }
536  }
537  }
538 
539  // scale the random values to be of order of the local mesh size
540  for (int i = 0; i < fespace->GetNDofs(); i++)
541  {
542  for (int d = 0; d < dim; d++)
543  {
544  rdm(fespace->DofToVDof(i,d)) *= h0(i);
545  }
546  }
547 
548  char move_bdr = 'n';
549  cout << "move boundary nodes? [y/n] ---> " << flush;
550  cin >> move_bdr;
551 
552  // don't perturb the boundary
553  if (move_bdr == 'n')
554  {
555  Array<int> vdofs;
556  for (int i = 0; i < fespace->GetNBE(); i++)
557  {
558  fespace->GetBdrElementVDofs(i, vdofs);
559  for (int j = 0; j < vdofs.Size(); j++)
560  {
561  rdm(vdofs[j]) = 0.0;
562  }
563  }
564  }
565 
566  *nodes += rdm;
567  }
568 
569  print_char = 1;
570  }
571 
572  if (mk == 'x')
573  {
574  int sd, nz = 0;
575  DenseMatrix J(dim);
576  double min_det_J, max_det_J, min_det_J_z, max_det_J_z;
577  double min_kappa, max_kappa, max_ratio_det_J_z;
578  min_det_J = min_kappa = infinity();
579  max_det_J = max_kappa = max_ratio_det_J_z = -infinity();
580  cout << "subdivision factor ---> " << flush;
581  cin >> sd;
582  Array<int> bad_elems_by_geom(Geometry::NumGeom);
583  bad_elems_by_geom = 0;
584  // Only print so many to keep output compact
585  const int max_to_print = 10;
586  for (int i = 0; i < mesh->GetNE(); i++)
587  {
588  Geometry::Type geom = mesh->GetElementBaseGeometry(i);
590 
591  RefinedGeometry *RefG = GlobGeometryRefiner.Refine(geom, sd, 1);
592  IntegrationRule &ir = RefG->RefPts;
593 
594  min_det_J_z = infinity();
595  max_det_J_z = -infinity();
596  for (int j = 0; j < ir.GetNPoints(); j++)
597  {
598  T->SetIntPoint(&ir.IntPoint(j));
599  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
600 
601  double det_J = J.Det();
602  double kappa =
603  J.CalcSingularvalue(0) / J.CalcSingularvalue(dim-1);
604 
605  min_det_J_z = fmin(min_det_J_z, det_J);
606  max_det_J_z = fmax(max_det_J_z, det_J);
607 
608  min_kappa = fmin(min_kappa, kappa);
609  max_kappa = fmax(max_kappa, kappa);
610  }
611  max_ratio_det_J_z =
612  fmax(max_ratio_det_J_z, max_det_J_z/min_det_J_z);
613  min_det_J = fmin(min_det_J, min_det_J_z);
614  max_det_J = fmax(max_det_J, max_det_J_z);
615  if (min_det_J_z <= 0.0)
616  {
617  if (nz < max_to_print)
618  {
619  Vector center;
620  mesh->GetElementCenter(i, center);
621  cout << "det(J) < 0 = " << min_det_J_z << " in element "
622  << i << ", centered at: ";
623  center.Print();
624  }
625  nz++;
626  bad_elems_by_geom[geom]++;
627  }
628  }
629  if (nz >= max_to_print)
630  {
631  cout << "det(J) < 0 for " << nz - max_to_print << " more elements "
632  << "not printed.\n";
633  }
634  cout << "\nbad elements = " << nz;
635  if (nz)
636  {
637  cout << " -- ";
638  Mesh::PrintElementsByGeometry(dim, bad_elems_by_geom, cout);
639  }
640  cout << "\nmin det(J) = " << min_det_J
641  << "\nmax det(J) = " << max_det_J
642  << "\nglobal ratio = " << max_det_J/min_det_J
643  << "\nmax el ratio = " << max_ratio_det_J_z
644  << "\nmin kappa = " << min_kappa
645  << "\nmax kappa = " << max_kappa << endl;
646  }
647 
648  if (mk == 'f')
649  {
650  DenseMatrix point_mat(sdim,1);
651  cout << "\npoint in physical space ---> " << flush;
652  for (int i = 0; i < sdim; i++)
653  {
654  cin >> point_mat(i,0);
655  }
656  Array<int> elem_ids;
658 
659  // physical -> reference space
660  mesh->FindPoints(point_mat, elem_ids, ips);
661 
662  cout << "point in reference space:";
663  if (elem_ids[0] == -1)
664  {
665  cout << " NOT FOUND!\n";
666  }
667  else
668  {
669  cout << " element " << elem_ids[0] << ", ip =";
670  cout << " " << ips[0].x;
671  if (sdim > 1)
672  {
673  cout << " " << ips[0].y;
674  if (sdim > 2)
675  {
676  cout << " " << ips[0].z;
677  }
678  }
679  cout << endl;
680  }
681  }
682 
683  if (mk == 'o')
684  {
685  cout << "What type of reordering?\n"
686  "g) Gecko edge-product minimization\n"
687  "h) Hilbert spatial sort\n"
688  "--> " << flush;
689  char rk;
690  cin >> rk;
691 
692  Array<int> ordering, tentative;
693  if (rk == 'h')
694  {
695  mesh->GetHilbertElementOrdering(ordering);
696  mesh->ReorderElements(ordering);
697  }
698  else if (rk == 'g')
699  {
700  int outer, inner, window, period;
701  cout << "Enter number of outer iterations (default 5): " << flush;
702  cin >> outer;
703  cout << "Enter number of inner iterations (default 4): " << flush;
704  cin >> inner;
705  cout << "Enter window size (default 4, beware of exponential cost): "
706  << flush;
707  cin >> window;
708  cout << "Enter period for window size increment (default 2): "
709  << flush;
710  cin >> period;
711 
712  double best_cost = infinity();
713  for (int i = 0; i < outer; i++)
714  {
715  int seed = i+1;
716  double cost = mesh->GetGeckoElementOrdering(
717  tentative, inner, window, period, seed, true);
718 
719  if (cost < best_cost)
720  {
721  ordering = tentative;
722  best_cost = cost;
723  }
724  }
725  cout << "Final cost: " << best_cost << endl;
726 
727  mesh->ReorderElements(ordering);
728  }
729  }
730 
731  // These are most of the cases that open a new GLVis window
732  if (mk == 'm' || mk == 'b' || mk == 'e' || mk == 'v' || mk == 'h' ||
733  mk == 'k' || mk == 'J' || mk == 'p' || mk == 'B' || mk == 'P')
734  {
735  FiniteElementSpace *bdr_attr_fespace = NULL;
736  FiniteElementSpace *attr_fespace =
737  new FiniteElementSpace(mesh, attr_fec);
738  GridFunction bdr_attr;
739  GridFunction attr(attr_fespace);
740 
741  if (mk == 'm')
742  {
743  for (int i = 0; i < mesh->GetNE(); i++)
744  {
745  attr(i) = mesh->GetAttribute(i);
746  }
747  }
748 
749  if (mk == 'P')
750  {
751  for (int i = 0; i < mesh->GetNE(); i++)
752  {
753  attr(i) = partitioning[i] + 1;
754  }
755  }
756 
757  if (mk == 'b' || mk == 'B')
758  {
759  if (dim == 3)
760  {
761  delete bdr_mesh;
762  bdr_mesh = skin_mesh(mesh);
763  bdr_attr_fespace =
764  new FiniteElementSpace(bdr_mesh, bdr_attr_fec);
765  bdr_attr.SetSpace(bdr_attr_fespace);
766  if (mk == 'b')
767  {
768  for (int i = 0; i < bdr_mesh->GetNE(); i++)
769  {
770  bdr_attr(i) = bdr_mesh->GetAttribute(i);
771  }
772  }
773  else if (mk == 'B')
774  {
775  for (int i = 0; i < bdr_mesh->GetNE(); i++)
776  {
777  bdr_attr(i) = bdr_partitioning[i] + 1;
778  }
779  }
780  else
781  {
782  MFEM_WARNING("Unimplemented case.");
783  }
784  }
785  else
786  {
787  MFEM_WARNING("Unsupported mesh dimension.");
788  attr = 1.0;
789  }
790  }
791 
792  if (mk == 'v')
793  {
794  attr = 1.0;
795  }
796 
797  if (mk == 'e')
798  {
799  Array<int> coloring;
800  srand(time(0));
801  double a = double(rand()) / (double(RAND_MAX) + 1.);
802  int el0 = (int)floor(a * mesh->GetNE());
803  cout << "Generating coloring starting with element " << el0+1
804  << " / " << mesh->GetNE() << endl;
805  mesh->GetElementColoring(coloring, el0);
806  for (int i = 0; i < coloring.Size(); i++)
807  {
808  attr(i) = coloring[i];
809  }
810  cout << "Number of colors: " << attr.Max() + 1 << endl;
811  for (int i = 0; i < mesh->GetNE(); i++)
812  {
813  attr(i) = i; // coloring by element number
814  }
815  }
816 
817  if (mk == 'h')
818  {
819  DenseMatrix J(dim);
820  double h_min, h_max;
821  h_min = infinity();
822  h_max = -h_min;
823  for (int i = 0; i < mesh->GetNE(); i++)
824  {
825  Geometry::Type geom = mesh->GetElementBaseGeometry(i);
827  T->SetIntPoint(&Geometries.GetCenter(geom));
828  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
829 
830  attr(i) = J.Det();
831  if (attr(i) < 0.0)
832  {
833  attr(i) = -pow(-attr(i), 1.0/double(dim));
834  }
835  else
836  {
837  attr(i) = pow(attr(i), 1.0/double(dim));
838  }
839  h_min = min(h_min, attr(i));
840  h_max = max(h_max, attr(i));
841  }
842  cout << "h_min = " << h_min << ", h_max = " << h_max << endl;
843  }
844 
845  if (mk == 'k')
846  {
847  DenseMatrix J(dim);
848  for (int i = 0; i < mesh->GetNE(); i++)
849  {
850  Geometry::Type geom = mesh->GetElementBaseGeometry(i);
852  T->SetIntPoint(&Geometries.GetCenter(geom));
853  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
854  attr(i) = J.CalcSingularvalue(0) / J.CalcSingularvalue(dim-1);
855  }
856  }
857 
858  if (mk == 'J')
859  {
860  // The "scaled Jacobian" is the determinant of the Jacobian scaled
861  // by the l2 norms of its columns. It can be used to identify badly
862  // skewed elements, since it takes values between 0 and 1, with 0
863  // corresponding to a flat element, and 1 to orthogonal columns.
864  DenseMatrix J(dim);
865  int sd;
866  cout << "subdivision factor ---> " << flush;
867  cin >> sd;
868  for (int i = 0; i < mesh->GetNE(); i++)
869  {
870  Geometry::Type geom = mesh->GetElementBaseGeometry(i);
872 
873  RefinedGeometry *RefG = GlobGeometryRefiner.Refine(geom, sd, 1);
874  IntegrationRule &ir = RefG->RefPts;
875 
876  // For each element, find the minimal scaled Jacobian in a
877  // lattice of points with the given subdivision factor.
878  attr(i) = infinity();
879  for (int j = 0; j < ir.GetNPoints(); j++)
880  {
881  T->SetIntPoint(&ir.IntPoint(j));
882  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
883 
884  // Jacobian determinant
885  double sJ = J.Det();
886 
887  for (int k = 0; k < J.Width(); k++)
888  {
889  Vector col;
890  J.GetColumnReference(k,col);
891  // Scale by column norms
892  sJ /= col.Norml2();
893  }
894 
895  attr(i) = fmin(sJ, attr(i));
896  }
897  }
898  }
899 
900  if (mk == 'p')
901  {
902  cout << "What type of partitioning?\n"
903  "c) Cartesian\n"
904  "s) Simple 1D split of the element sequence\n"
905  "0) METIS_PartGraphRecursive (sorted neighbor lists)\n"
906  "1) METIS_PartGraphKway (sorted neighbor lists)"
907  " (default)\n"
908  "2) METIS_PartGraphVKway (sorted neighbor lists)\n"
909  "3) METIS_PartGraphRecursive\n"
910  "4) METIS_PartGraphKway\n"
911  "5) METIS_PartGraphVKway\n"
912  "--> " << flush;
913  char pk;
914  cin >> pk;
915  if (pk == 'c')
916  {
917  int nxyz[3];
918  cout << "Enter nx: " << flush;
919  cin >> nxyz[0]; np = nxyz[0];
920  if (mesh->Dimension() > 1)
921  {
922  cout << "Enter ny: " << flush;
923  cin >> nxyz[1]; np *= nxyz[1];
924  if (mesh->Dimension() > 2)
925  {
926  cout << "Enter nz: " << flush;
927  cin >> nxyz[2]; np *= nxyz[2];
928  }
929  }
930  partitioning = Array<int>(mesh->CartesianPartitioning(nxyz), mesh->GetNE());
931  recover_bdr_partitioning(mesh, partitioning, bdr_partitioning);
932  }
933  else if (pk == 's')
934  {
935  cout << "Enter number of processors: " << flush;
936  cin >> np;
937 
938  partitioning.SetSize(mesh->GetNE());
939  for (int i = 0; i < mesh->GetNE(); i++)
940  {
941  partitioning[i] = i * np / mesh->GetNE();
942  }
943  recover_bdr_partitioning(mesh, partitioning, bdr_partitioning);
944  }
945  else
946  {
947  int part_method = pk - '0';
948  if (part_method < 0 || part_method > 5)
949  {
950  continue;
951  }
952  cout << "Enter number of processors: " << flush;
953  cin >> np;
954  partitioning = Array<int>(mesh->GeneratePartitioning(np, part_method),
955  mesh->GetNE());
956  recover_bdr_partitioning(mesh, partitioning, bdr_partitioning);
957  }
958  if (partitioning)
959  {
960  const char part_file[] = "partitioning.txt";
961  ofstream opart(part_file);
962  opart << "number_of_elements " << mesh->GetNE() << '\n'
963  << "number_of_processors " << np << '\n';
964  for (int i = 0; i < mesh->GetNE(); i++)
965  {
966  opart << partitioning[i] << '\n';
967  }
968  cout << "Partitioning file: " << part_file << endl;
969 
970  Array<int> proc_el(np);
971  proc_el = 0;
972  for (int i = 0; i < mesh->GetNE(); i++)
973  {
974  proc_el[partitioning[i]]++;
975  }
976  int min_el = proc_el[0], max_el = proc_el[0];
977  for (int i = 1; i < np; i++)
978  {
979  if (min_el > proc_el[i])
980  {
981  min_el = proc_el[i];
982  }
983  if (max_el < proc_el[i])
984  {
985  max_el = proc_el[i];
986  }
987  }
988  cout << "Partitioning stats:\n"
989  << " "
990  << setw(12) << "minimum"
991  << setw(12) << "average"
992  << setw(12) << "maximum"
993  << setw(12) << "total" << '\n';
994  cout << " elements "
995  << setw(12) << min_el
996  << setw(12) << double(mesh->GetNE())/np
997  << setw(12) << max_el
998  << setw(12) << mesh->GetNE() << endl;
999  }
1000  else
1001  {
1002  continue;
1003  }
1004 
1005  for (int i = 0; i < mesh->GetNE(); i++)
1006  {
1007  attr(i) = partitioning[i] + 1;
1008  }
1009  }
1010 
1011  char vishost[] = "localhost";
1012  int visport = 19916;
1013  socketstream sol_sock(vishost, visport);
1014  if (sol_sock.is_open())
1015  {
1016  sol_sock.precision(14);
1017  if (sdim == 2)
1018  {
1019  sol_sock << "fem2d_gf_data_keys\n";
1020  if (mk != 'p')
1021  {
1022  mesh->Print(sol_sock);
1023  }
1024  else
1025  {
1026  // NURBS meshes do not support PrintWithPartitioning
1027  if (mesh->NURBSext)
1028  {
1029  mesh->Print(sol_sock);
1030  for (int i = 0; i < mesh->GetNE(); i++)
1031  {
1032  attr(i) = partitioning[i];
1033  }
1034  }
1035  else
1036  {
1037  mesh->PrintWithPartitioning(partitioning, sol_sock, 1);
1038  }
1039  }
1040  attr.Save(sol_sock);
1041  sol_sock << "RjlmAb***********";
1042  if (mk == 'v')
1043  {
1044  sol_sock << "e";
1045  }
1046  else
1047  {
1048  sol_sock << "\n";
1049  }
1050  }
1051  else
1052  {
1053  sol_sock << "fem3d_gf_data_keys\n";
1054  if (mk == 'v' || mk == 'h' || mk == 'k' || mk == 'J')
1055  {
1056  mesh->Print(sol_sock);
1057  }
1058  else if (mk == 'b' || mk == 'B')
1059  {
1060  bdr_mesh->Print(sol_sock);
1061  bdr_attr.Save(sol_sock);
1062  sol_sock << "mcaaA";
1063  // Switch to a discrete color scale
1064  sol_sock << "pppppp" << "pppppp" << "pppppp";
1065  }
1066  else
1067  {
1068  // NURBS meshes do not support PrintWithPartitioning
1069  if (mesh->NURBSext)
1070  {
1071  mesh->Print(sol_sock);
1072  for (int i = 0; i < mesh->GetNE(); i++)
1073  {
1074  attr(i) = partitioning[i];
1075  }
1076  }
1077  else
1078  {
1079  mesh->PrintWithPartitioning(partitioning, sol_sock);
1080  }
1081  }
1082  if (mk != 'b' && mk != 'B')
1083  {
1084  attr.Save(sol_sock);
1085  sol_sock << "maaA";
1086  if (mk == 'v')
1087  {
1088  sol_sock << "aa";
1089  }
1090  else
1091  {
1092  sol_sock << "\n";
1093  }
1094  }
1095  }
1096  sol_sock << flush;
1097  }
1098  else
1099  {
1100  cout << "Unable to connect to "
1101  << vishost << ':' << visport << endl;
1102  }
1103  delete attr_fespace;
1104  delete bdr_attr_fespace;
1105  }
1106 
1107  if (mk == 'l')
1108  {
1109  // Project and plot the function 'f'
1110  int p;
1111  FiniteElementCollection *fec = NULL;
1112  cout << "Enter projection space order: " << flush;
1113  cin >> p;
1114  if (p >= 1)
1115  {
1116  fec = new H1_FECollection(p, mesh->Dimension(),
1118  }
1119  else
1120  {
1121  fec = new DG_FECollection(-p, mesh->Dimension(),
1123  }
1124  FiniteElementSpace fes(mesh, fec);
1125  GridFunction level(&fes);
1126  FunctionCoefficient coeff(f);
1127  level.ProjectCoefficient(coeff);
1128  char vishost[] = "localhost";
1129  int visport = 19916;
1130  socketstream sol_sock(vishost, visport);
1131  if (sol_sock.is_open())
1132  {
1133  sol_sock.precision(14);
1134  sol_sock << "solution\n" << *mesh << level << flush;
1135  }
1136  else
1137  {
1138  cout << "Unable to connect to "
1139  << vishost << ':' << visport << endl;
1140  }
1141  delete fec;
1142  }
1143 
1144  if (mk == 'S')
1145  {
1146  const char omesh_file[] = "mesh-explorer.mesh";
1147  ofstream omesh(omesh_file);
1148  omesh.precision(14);
1149  mesh->Print(omesh);
1150  cout << "New mesh file: " << omesh_file << endl;
1151  }
1152 
1153  if (mk == 'V')
1154  {
1155  const char omesh_file[] = "mesh-explorer.vtk";
1156  ofstream omesh(omesh_file);
1157  omesh.precision(14);
1158  mesh->PrintVTK(omesh);
1159  cout << "New VTK mesh file: " << omesh_file << endl;
1160  }
1161 
1162 #ifdef MFEM_USE_ZLIB
1163  if (mk == 'Z')
1164  {
1165  const char omesh_file[] = "mesh-explorer.mesh.gz";
1166  ofgzstream omesh(omesh_file, "zwb9");
1167  omesh.precision(14);
1168  mesh->Print(omesh);
1169  cout << "New mesh file: " << omesh_file << endl;
1170  }
1171 #endif
1172 
1173  }
1174 
1175  delete bdr_attr_fec;
1176  delete attr_fec;
1177  delete bdr_mesh;
1178  delete mesh;
1179  return 0;
1180 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
const char vishost[]
bool Help() const
Return true if we are flagged to print the help message.
Definition: optparser.hpp:153
Geometry::Type GetGeometryType() const
Definition: element.hpp:52
double GetGeckoElementOrdering(Array< int > &ordering, int iterations=4, int window=4, int period=2, int seed=0, bool verbose=false, double time_limit=0)
Definition: mesh.cpp:2018
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:560
int * CartesianPartitioning(int nxyz[])
Definition: mesh.cpp:6936
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:560
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
const double * GetVertex(int i) const
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:1005
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:251
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1662
virtual void GetVertices(Array< int > &v) const =0
Returns element&#39;s vertices.
static void PrintElementsByGeometry(int dim, const Array< int > &num_elems_by_geom, std::ostream &out)
Auxiliary method used by PrintCharacteristics().
Definition: mesh.cpp:228
static const int NumGeom
Definition: geom.hpp:42
void JacToPerfJac(int GeomType, const DenseMatrix &J, DenseMatrix &PJ) const
Definition: geom.cpp:867
int * GeneratePartitioning(int nparts, int part_method=1)
Definition: mesh.cpp:6975
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:931
double Det() const
Definition: densemat.cpp:436
void NewNodes(GridFunction &nodes, bool make_owner=false)
Replace the internal node GridFunction with the given GridFunction.
Definition: mesh.cpp:7865
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:792
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:93
static Mesh MakeRefined(Mesh &orig_mesh, int ref_factor, int ref_type)
Create a refined (by any factor) version of orig_mesh.
Definition: mesh.cpp:3721
int AddTriangle(int v1, int v2, int v3, int attr=1)
Definition: mesh.cpp:1648
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
void Transform(void(*f)(const Vector &, Vector &))
Definition: mesh.cpp:11551
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:928
void PrintHelp(std::ostream &out) const
Print the help message.
Definition: optparser.cpp:398
void PrintCharacteristics(Vector *Vh=NULL, Vector *Vk=NULL, std::ostream &os=mfem::out)
Compute and print mesh characteristics such as number of vertices, number of elements, number of boundary elements, minimal and maximal element sizes, minimal and maximal element aspect ratios, etc.
Definition: mesh.cpp:242
double kappa
Definition: ex24.cpp:54
void MakeOwner(FiniteElementCollection *fec_)
Make the GridFunction the owner of fec and fes.
Definition: gridfunc.hpp:124
void Randomize(int seed=0)
Set random values in the vector.
Definition: vector.cpp:772
const IntegrationPoint & GetCenter(int GeomType)
Return the center of the given Geometry::Type, GeomType.
Definition: geom.hpp:71
Mesh * read_par_mesh(int np, const char *mesh_prefix, Array< int > &partitioning, Array< int > &bdr_partitioning)
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:1062
void DeleteAll()
Delete the whole array.
Definition: array.hpp:846
Mesh * skin_mesh(Mesh *mesh)
DofTransformation * GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:297
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3619
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:590
Geometry Geometries
Definition: fe.cpp:49
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
int AddVertex(double x, double y=0.0, double z=0.0)
Definition: mesh.cpp:1601
void PrintWithPartitioning(int *partitioning, std::ostream &os, int elem_attr=0) const
Prints the mesh with boundary elements given by the boundary of the subdomains, so that the boundary ...
Definition: mesh.cpp:10852
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
double GetElementSize(ElementTransformation *T, int type=0)
Definition: mesh.cpp:76
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:599
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
void PrintVTK(std::ostream &os)
Definition: mesh.cpp:10224
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:751
void transformation(const Vector &p, Vector &v)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9737
double region_eps
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:119
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:5312
GeometryRefiner GlobGeometryRefiner
Definition: geom.cpp:1773
int GetAttribute() const
Return element&#39;s attribute.
Definition: element.hpp:55
void GetBdrElementAdjacentElement(int bdr_el, int &el, int &info) const
For the given boundary element, bdr_el, return its adjacent element and its info, i...
Definition: mesh.cpp:6175
void GetHilbertElementOrdering(Array< int > &ordering)
Definition: mesh.cpp:2185
Piecewise-constant discontinuous finite elements in 3D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:1071
void GetElementCenter(int i, Vector &center)
Definition: mesh.cpp:68
IntegrationRule RefPts
Definition: geom.hpp:282
RefinedGeometry * Refine(Geometry::Type Geom, int Times, int ETimes=1)
Definition: geom.cpp:1099
const int visport
Piecewise-constant discontinuous finite elements in 2D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:893
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:652
int Dimension() const
Definition: mesh.hpp:999
virtual DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element &#39;elem&#39;.
Definition: fespace.cpp:2680
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
Definition: fe_coll.hpp:334
FDualNumber< tbase > pow(const FDualNumber< tbase > &a, const FDualNumber< tbase > &b)
pow([dual number],[dual number])
Definition: fdual.hpp:543
int SpaceDimension() const
Definition: mesh.hpp:1000
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:270
bool is_open()
True if the socketstream is open, false otherwise.
int AddSegment(int v1, int v2, int attr=1)
Definition: mesh.cpp:1634
void Print(std::ostream &out=mfem::out, int width=8) const
Prints vector to stream out.
Definition: vector.cpp:711
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
void FinalizeTopology(bool generate_bdr=true)
Finalize the construction of the secondary topology (connectivity) data of a Mesh.
Definition: mesh.cpp:2878
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:679
void GetColumnReference(int c, Vector &col)
Definition: densemat.hpp:285
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
Definition: fe_coll.cpp:51
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:272
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:925
void PrintError(std::ostream &out) const
Print the error message.
Definition: optparser.cpp:355
virtual const char * Name() const
Definition: fe_coll.hpp:61
A standard isoparametric element transformation.
Definition: eltrans.hpp:361
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:348
double Max() const
Returns the maximal element of the vector.
Definition: vector.cpp:882
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:5307
virtual void Print(std::ostream &os=mfem::out) const
Definition: mesh.hpp:1646
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition: mesh.cpp:2971
void ReorderElements(const Array< int > &ordering, bool reorder_vertices=true)
Definition: mesh.cpp:2237
int dim
Definition: ex24.cpp:53
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:324
double CalcSingularvalue(const int i) const
Return the i-th singular value (decreasing order) of NxN matrix, N=1,2,3.
Definition: densemat.cpp:1207
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:46
virtual int GetNVertices() const =0
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: gridfunc.cpp:2395
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition: gridfunc.cpp:193
A general function coefficient.
virtual int FindPoints(DenseMatrix &point_mat, Array< int > &elem_ids, Array< IntegrationPoint > &ips, bool warn=true, InverseElementTransformation *inv_trans=NULL)
Find the ids of the elements that contain the given points, and their corresponding reference coordin...
Definition: mesh.cpp:11820
Vector data type.
Definition: vector.hpp:60
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:577
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7841
Nodes: x_i = i/(n-1), i=0,...,n-1.
Definition: fe_base.hpp:34
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
void GetElementColoring(Array< int > &colors, int el0=0)
Definition: mesh.cpp:10777
RefCoord s[3]
const DenseMatrix & GetPointMat() const
Return the stored point matrix.
Definition: eltrans.hpp:404
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:9326
double region(const Vector &p)
Abstract data type element.
Definition: element.hpp:28
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:1455
int main()
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:268
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:1037
void recover_bdr_partitioning(const Mesh *mesh, const Array< int > &partitioning, Array< int > &bdr_partitioning)
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150