MFEM  v4.2.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-2020, 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)
98 {
99  Mesh *mesh;
100  Array<Mesh *> mesh_array;
101 
102  mesh_array.SetSize(np);
103  for (int p = 0; p < np; p++)
104  {
105  ostringstream fname;
106  fname << mesh_prefix << '.' << setfill('0') << setw(6) << p;
107  ifgzstream meshin(fname.str().c_str());
108  if (!meshin)
109  {
110  cerr << "Can not open mesh file: " << fname.str().c_str()
111  << '!' << endl;
112  for (p--; p >= 0; p--)
113  {
114  delete mesh_array[p];
115  }
116  return NULL;
117  }
118  mesh_array[p] = new Mesh(meshin, 1, 0);
119  // set element and boundary attributes to be the processor number + 1
120  if (1)
121  {
122  for (int i = 0; i < mesh_array[p]->GetNE(); i++)
123  {
124  mesh_array[p]->GetElement(i)->SetAttribute(p+1);
125  }
126  for (int i = 0; i < mesh_array[p]->GetNBE(); i++)
127  {
128  mesh_array[p]->GetBdrElement(i)->SetAttribute(p+1);
129  }
130  }
131  }
132  mesh = new Mesh(mesh_array, np);
133 
134  for (int p = 0; p < np; p++)
135  {
136  delete mesh_array[np-1-p];
137  }
138  mesh_array.DeleteAll();
139 
140  return mesh;
141 }
142 
143 // Given a 3D mesh, produce a 2D mesh consisting of its boundary elements.
145 {
146  // Determine mapping from vertex to boundary vertex
147  Array<int> v2v(mesh->GetNV());
148  v2v = -1;
149  for (int i = 0; i < mesh->GetNBE(); i++)
150  {
151  Element *el = mesh->GetBdrElement(i);
152  int *v = el->GetVertices();
153  int nv = el->GetNVertices();
154  for (int j = 0; j < nv; j++)
155  {
156  v2v[v[j]] = 0;
157  }
158  }
159  int nbvt = 0;
160  for (int i = 0; i < v2v.Size(); i++)
161  {
162  if (v2v[i] == 0)
163  {
164  v2v[i] = nbvt++;
165  }
166  }
167 
168  // Create a new mesh for the boundary
169  Mesh * bmesh = new Mesh(mesh->Dimension() - 1, nbvt, mesh->GetNBE(),
170  0, mesh->SpaceDimension());
171 
172  // Copy vertices to the boundary mesh
173  nbvt = 0;
174  for (int i = 0; i < v2v.Size(); i++)
175  {
176  if (v2v[i] >= 0)
177  {
178  double *c = mesh->GetVertex(i);
179  bmesh->AddVertex(c);
180  nbvt++;
181  }
182  }
183 
184  // Copy elements to the boundary mesh
185  int bv[4];
186  for (int i = 0; i < mesh->GetNBE(); i++)
187  {
188  Element *el = mesh->GetBdrElement(i);
189  int *v = el->GetVertices();
190  int nv = el->GetNVertices();
191 
192  for (int j = 0; j < nv; j++)
193  {
194  bv[j] = v2v[v[j]];
195  }
196 
197  switch (el->GetGeometryType())
198  {
199  case Geometry::SEGMENT:
200  bmesh->AddSegment(bv, el->GetAttribute());
201  break;
202  case Geometry::TRIANGLE:
203  bmesh->AddTriangle(bv, el->GetAttribute());
204  break;
205  case Geometry::SQUARE:
206  bmesh->AddQuad(bv, el->GetAttribute());
207  break;
208  default:
209  break; /// This should not happen
210  }
211 
212  }
213  bmesh->FinalizeTopology();
214 
215  // Copy GridFunction describing nodes if present
216  if (mesh->GetNodes())
217  {
218  FiniteElementSpace *fes = mesh->GetNodes()->FESpace();
219  const FiniteElementCollection *fec = fes->FEColl();
220  if (dynamic_cast<const H1_FECollection*>(fec))
221  {
222  FiniteElementCollection *fec_copy =
224  FiniteElementSpace *fes_copy =
225  new FiniteElementSpace(*fes, bmesh, fec_copy);
226  GridFunction *bdr_nodes = new GridFunction(fes_copy);
227  bdr_nodes->MakeOwner(fec_copy);
228 
229  bmesh->NewNodes(*bdr_nodes, true);
230 
231  Array<int> vdofs;
232  Array<int> bvdofs;
233  Vector v;
234  for (int i=0; i<mesh->GetNBE(); i++)
235  {
236  fes->GetBdrElementVDofs(i, vdofs);
237  mesh->GetNodes()->GetSubVector(vdofs, v);
238 
239  fes_copy->GetElementVDofs(i, bvdofs);
240  bdr_nodes->SetSubVector(bvdofs, v);
241  }
242  }
243  else
244  {
245  cout << "\nDiscontinuous nodes not yet supported" << endl;
246  }
247  }
248 
249  return bmesh;
250 }
251 
252 int main (int argc, char *argv[])
253 {
254  int np = 0;
255  const char *mesh_file = "../../data/beam-hex.mesh";
256  bool refine = true;
257 
258  OptionsParser args(argc, argv);
259  args.AddOption(&mesh_file, "-m", "--mesh",
260  "Mesh file to visualize.");
261  args.AddOption(&np, "-np", "--num-proc",
262  "Load mesh from multiple processors.");
263  args.AddOption(&refine, "-ref", "--refinement", "-no-ref", "--no-refinement",
264  "Prepare the mesh for refinement or not.");
265  args.Parse();
266  if (!args.Good())
267  {
268  if (!args.Help())
269  {
270  args.PrintError(cout);
271  cout << endl;
272  }
273  cout << "Visualize and manipulate a serial mesh:\n"
274  << " mesh-explorer -m <mesh_file>\n"
275  << "Visualize and manipulate a parallel mesh:\n"
276  << " mesh-explorer -np <#proc> -m <mesh_prefix>\n" << endl
277  << "All Options:\n";
278  args.PrintHelp(cout);
279  return 1;
280  }
281  args.PrintOptions(cout);
282 
283  Mesh *mesh;
284  Mesh *bdr_mesh = NULL;
285  if (np <= 0)
286  {
287  mesh = new Mesh(mesh_file, 1, refine);
288  }
289  else
290  {
291  mesh = read_par_mesh(np, mesh_file);
292  if (mesh == NULL)
293  {
294  return 3;
295  }
296  }
297  int dim = mesh->Dimension();
298  int sdim = mesh->SpaceDimension();
299 
300  FiniteElementCollection *bdr_attr_fec = NULL;
301  FiniteElementCollection *attr_fec;
302  if (dim == 2)
303  {
304  attr_fec = new Const2DFECollection;
305  }
306  else
307  {
308  bdr_attr_fec = new Const2DFECollection;
309  attr_fec = new Const3DFECollection;
310  }
311 
312  int print_char = 1;
313  while (1)
314  {
315  if (print_char)
316  {
317  cout << endl;
318  mesh->PrintCharacteristics();
319  cout << "boundary attribs :";
320  for (int i = 0; i < mesh->bdr_attributes.Size(); i++)
321  {
322  cout << ' ' << mesh->bdr_attributes[i];
323  }
324  cout << '\n' << "material attribs :";
325  for (int i = 0; i < mesh->attributes.Size(); i++)
326  {
327  cout << ' ' << mesh->attributes[i];
328  }
329  cout << endl;
330  cout << "mesh curvature : ";
331  if (mesh->GetNodalFESpace() != NULL)
332  {
333  cout << mesh->GetNodalFESpace()->FEColl()->Name() << endl;
334  }
335  else
336  {
337  cout << "NONE" << endl;
338  }
339  }
340  print_char = 0;
341  cout << endl;
342  cout << "What would you like to do?\n"
343  "r) Refine\n"
344  "c) Change curvature\n"
345  "s) Scale\n"
346  "t) Transform\n"
347  "j) Jitter\n"
348  "v) View\n"
349  "m) View materials\n"
350  "b) View boundary\n"
351  "e) View elements\n"
352  "h) View element sizes, h\n"
353  "k) View element ratios, kappa\n"
354  "l) Plot a function\n"
355  "x) Print sub-element stats\n"
356  "f) Find physical point in reference space\n"
357  "p) Generate a partitioning\n"
358  "o) Reorder elements\n"
359  "S) Save in MFEM format\n"
360  "V) Save in VTK format (only linear and quadratic meshes)\n"
361  "q) Quit\n"
362 #ifdef MFEM_USE_ZLIB
363  "Z) Save in MFEM format with compression\n"
364 #endif
365  "--> " << flush;
366  char mk;
367  cin >> mk;
368  if (!cin) { break; }
369 
370  if (mk == 'q')
371  {
372  break;
373  }
374 
375  if (mk == 'r')
376  {
377  cout <<
378  "Choose type of refinement:\n"
379  "s) standard refinement with Mesh::UniformRefinement()\n"
380  "b) Mesh::UniformRefinement() (bisection for tet meshes)\n"
381  "u) uniform refinement with a factor\n"
382  "g) non-uniform refinement (Gauss-Lobatto) with a factor\n"
383  "l) refine locally using the region() function\n"
384  "--> " << flush;
385  char sk;
386  cin >> sk;
387  switch (sk)
388  {
389  case 's':
390  mesh->UniformRefinement();
391  // Make sure tet-only meshes are marked for local refinement.
392  mesh->Finalize(true);
393  break;
394  case 'b':
395  mesh->UniformRefinement(1); // ref_algo = 1
396  break;
397  case 'u':
398  case 'g':
399  {
400  cout << "enter refinement factor --> " << flush;
401  int ref_factor;
402  cin >> ref_factor;
403  if (ref_factor <= 1 || ref_factor > 32) { break; }
404  int ref_type = (sk == 'u') ? BasisType::ClosedUniform :
406  Mesh *rmesh = new Mesh(mesh, ref_factor, ref_type);
407  delete mesh;
408  mesh = rmesh;
409  break;
410  }
411  case 'l':
412  {
413  Vector pt;
414  Array<int> marked_elements;
415  for (int i = 0; i < mesh->GetNE(); i++)
416  {
417  // check all nodes of the element
419  mesh->GetElementTransformation(i, &T);
420  for (int j = 0; j < T.GetPointMat().Width(); j++)
421  {
422  T.GetPointMat().GetColumnReference(j, pt);
423  if (region(pt) <= region_eps)
424  {
425  marked_elements.Append(i);
426  break;
427  }
428  }
429  }
430  mesh->GeneralRefinement(marked_elements);
431  break;
432  }
433  }
434  print_char = 1;
435  }
436 
437  if (mk == 'c')
438  {
439  int p;
440  cout << "enter new order for mesh curvature --> " << flush;
441  cin >> p;
442  mesh->SetCurvature(p > 0 ? p : -p, p <= 0);
443  print_char = 1;
444  }
445 
446  if (mk == 's')
447  {
448  double factor;
449  cout << "scaling factor ---> " << flush;
450  cin >> factor;
451 
452  GridFunction *nodes = mesh->GetNodes();
453  if (nodes == NULL)
454  {
455  for (int i = 0; i < mesh->GetNV(); i++)
456  {
457  double *v = mesh->GetVertex(i);
458  v[0] *= factor;
459  v[1] *= factor;
460  if (dim == 3)
461  {
462  v[2] *= factor;
463  }
464  }
465  }
466  else
467  {
468  *nodes *= factor;
469  }
470 
471  print_char = 1;
472  }
473 
474  if (mk == 't')
475  {
476  mesh->Transform(transformation);
477  print_char = 1;
478  }
479 
480  if (mk == 'j')
481  {
482  double jitter;
483  cout << "jitter factor ---> " << flush;
484  cin >> jitter;
485 
486  GridFunction *nodes = mesh->GetNodes();
487 
488  if (nodes == NULL)
489  {
490  cerr << "The mesh should have nodes, introduce curvature first!\n";
491  }
492  else
493  {
494  FiniteElementSpace *fespace = nodes->FESpace();
495 
496  GridFunction rdm(fespace);
497  rdm.Randomize();
498  rdm -= 0.5; // shift to random values in [-0.5,0.5]
499  rdm *= jitter;
500 
501  // compute minimal local mesh size
502  Vector h0(fespace->GetNDofs());
503  h0 = infinity();
504  {
505  Array<int> dofs;
506  for (int i = 0; i < fespace->GetNE(); i++)
507  {
508  fespace->GetElementDofs(i, dofs);
509  for (int j = 0; j < dofs.Size(); j++)
510  {
511  h0(dofs[j]) = std::min(h0(dofs[j]), mesh->GetElementSize(i));
512  }
513  }
514  }
515 
516  // scale the random values to be of order of the local mesh size
517  for (int i = 0; i < fespace->GetNDofs(); i++)
518  {
519  for (int d = 0; d < dim; d++)
520  {
521  rdm(fespace->DofToVDof(i,d)) *= h0(i);
522  }
523  }
524 
525  char move_bdr = 'n';
526  cout << "move boundary nodes? [y/n] ---> " << flush;
527  cin >> move_bdr;
528 
529  // don't perturb the boundary
530  if (move_bdr == 'n')
531  {
532  Array<int> vdofs;
533  for (int i = 0; i < fespace->GetNBE(); i++)
534  {
535  fespace->GetBdrElementVDofs(i, vdofs);
536  for (int j = 0; j < vdofs.Size(); j++)
537  {
538  rdm(vdofs[j]) = 0.0;
539  }
540  }
541  }
542 
543  *nodes += rdm;
544  }
545 
546  print_char = 1;
547  }
548 
549  if (mk == 'x')
550  {
551  int sd, nz = 0;
552  DenseMatrix J(dim);
553  double min_det_J, max_det_J, min_det_J_z, max_det_J_z;
554  double min_kappa, max_kappa, max_ratio_det_J_z;
555  min_det_J = min_kappa = infinity();
556  max_det_J = max_kappa = max_ratio_det_J_z = -infinity();
557  cout << "subdivision factor ---> " << flush;
558  cin >> sd;
559  Array<int> bad_elems_by_geom(Geometry::NumGeom);
560  bad_elems_by_geom = 0;
561  for (int i = 0; i < mesh->GetNE(); i++)
562  {
565 
566  RefinedGeometry *RefG = GlobGeometryRefiner.Refine(geom, sd, 1);
567  IntegrationRule &ir = RefG->RefPts;
568 
569  min_det_J_z = infinity();
570  max_det_J_z = -infinity();
571  for (int j = 0; j < ir.GetNPoints(); j++)
572  {
573  T->SetIntPoint(&ir.IntPoint(j));
574  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
575 
576  double det_J = J.Det();
577  double kappa =
578  J.CalcSingularvalue(0) / J.CalcSingularvalue(dim-1);
579 
580  min_det_J_z = fmin(min_det_J_z, det_J);
581  max_det_J_z = fmax(max_det_J_z, det_J);
582 
583  min_kappa = fmin(min_kappa, kappa);
584  max_kappa = fmax(max_kappa, kappa);
585  }
586  max_ratio_det_J_z =
587  fmax(max_ratio_det_J_z, max_det_J_z/min_det_J_z);
588  min_det_J = fmin(min_det_J, min_det_J_z);
589  max_det_J = fmax(max_det_J, max_det_J_z);
590  if (min_det_J_z <= 0.0)
591  {
592  nz++;
593  bad_elems_by_geom[geom]++;
594  }
595  }
596  cout << "\nbad elements = " << nz;
597  if (nz)
598  {
599  cout << " -- ";
600  Mesh::PrintElementsByGeometry(dim, bad_elems_by_geom, cout);
601  }
602  cout << "\nmin det(J) = " << min_det_J
603  << "\nmax det(J) = " << max_det_J
604  << "\nglobal ratio = " << max_det_J/min_det_J
605  << "\nmax el ratio = " << max_ratio_det_J_z
606  << "\nmin kappa = " << min_kappa
607  << "\nmax kappa = " << max_kappa << endl;
608  }
609 
610  if (mk == 'f')
611  {
612  DenseMatrix point_mat(sdim,1);
613  cout << "\npoint in physical space ---> " << flush;
614  for (int i = 0; i < sdim; i++)
615  {
616  cin >> point_mat(i,0);
617  }
618  Array<int> elem_ids;
620 
621  // physical -> reference space
622  mesh->FindPoints(point_mat, elem_ids, ips);
623 
624  cout << "point in reference space:";
625  if (elem_ids[0] == -1)
626  {
627  cout << " NOT FOUND!\n";
628  }
629  else
630  {
631  cout << " element " << elem_ids[0] << ", ip =";
632  cout << " " << ips[0].x;
633  if (sdim > 1)
634  {
635  cout << " " << ips[0].y;
636  if (sdim > 2)
637  {
638  cout << " " << ips[0].z;
639  }
640  }
641  cout << endl;
642  }
643  }
644 
645  if (mk == 'o')
646  {
647  cout << "What type of reordering?\n"
648  "g) Gecko edge-product minimization\n"
649  "h) Hilbert spatial sort\n"
650  "--> " << flush;
651  char rk;
652  cin >> rk;
653 
654  Array<int> ordering, tentative;
655  if (rk == 'h')
656  {
657  mesh->GetHilbertElementOrdering(ordering);
658  mesh->ReorderElements(ordering);
659  }
660  else if (rk == 'g')
661  {
662  int outer, inner, window, period;
663  cout << "Enter number of outer iterations (default 5): " << flush;
664  cin >> outer;
665  cout << "Enter number of inner iterations (default 4): " << flush;
666  cin >> inner;
667  cout << "Enter window size (default 4, beware of exponential cost): "
668  << flush;
669  cin >> window;
670  cout << "Enter period for window size increment (default 2): "
671  << flush;
672  cin >> period;
673 
674  double best_cost = infinity();
675  for (int i = 0; i < outer; i++)
676  {
677  int seed = i+1;
678  double cost = mesh->GetGeckoElementOrdering(
679  tentative, inner, window, period, seed, true);
680 
681  if (cost < best_cost)
682  {
683  ordering = tentative;
684  best_cost = cost;
685  }
686  }
687  cout << "Final cost: " << best_cost << endl;
688 
689  mesh->ReorderElements(ordering);
690  }
691  }
692 
693  // These are most of the cases that open a new GLVis window
694  if (mk == 'm' || mk == 'b' || mk == 'e' || mk == 'v' || mk == 'h' ||
695  mk == 'k' || mk == 'p')
696  {
697  Array<int> bdr_part;
698  Array<int> part(mesh->GetNE());
699  FiniteElementSpace *bdr_attr_fespace = NULL;
700  FiniteElementSpace *attr_fespace =
701  new FiniteElementSpace(mesh, attr_fec);
702  GridFunction bdr_attr;
703  GridFunction attr(attr_fespace);
704 
705  if (mk == 'm')
706  {
707  for (int i = 0; i < mesh->GetNE(); i++)
708  {
709  part[i] = (attr(i) = mesh->GetAttribute(i)) - 1;
710  }
711  }
712 
713  if (mk == 'b')
714  {
715  if (dim == 3)
716  {
717  delete bdr_mesh;
718  bdr_mesh = skin_mesh(mesh);
719  bdr_attr_fespace =
720  new FiniteElementSpace(bdr_mesh, bdr_attr_fec);
721  bdr_part.SetSize(bdr_mesh->GetNE());
722  bdr_attr.SetSpace(bdr_attr_fespace);
723  for (int i = 0; i < bdr_mesh->GetNE(); i++)
724  {
725  bdr_part[i] = (bdr_attr(i) = bdr_mesh->GetAttribute(i)) - 1;
726  }
727  }
728  else
729  {
730  attr = 1.0;
731  }
732  }
733 
734  if (mk == 'v')
735  {
736  attr = 1.0;
737  }
738 
739  if (mk == 'e')
740  {
741  Array<int> coloring;
742  srand(time(0));
743  double a = double(rand()) / (double(RAND_MAX) + 1.);
744  int el0 = (int)floor(a * mesh->GetNE());
745  cout << "Generating coloring starting with element " << el0+1
746  << " / " << mesh->GetNE() << endl;
747  mesh->GetElementColoring(coloring, el0);
748  for (int i = 0; i < coloring.Size(); i++)
749  {
750  attr(i) = coloring[i];
751  }
752  cout << "Number of colors: " << attr.Max() + 1 << endl;
753  for (int i = 0; i < mesh->GetNE(); i++)
754  {
755  // part[i] = i; // checkerboard element coloring
756  attr(i) = part[i] = i; // coloring by element number
757  }
758  }
759 
760  if (mk == 'h')
761  {
762  DenseMatrix J(dim);
763  double h_min, h_max;
764  h_min = infinity();
765  h_max = -h_min;
766  for (int i = 0; i < mesh->GetNE(); i++)
767  {
768  int geom = mesh->GetElementBaseGeometry(i);
770  T->SetIntPoint(&Geometries.GetCenter(geom));
771  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
772 
773  attr(i) = J.Det();
774  if (attr(i) < 0.0)
775  {
776  attr(i) = -pow(-attr(i), 1.0/double(dim));
777  }
778  else
779  {
780  attr(i) = pow(attr(i), 1.0/double(dim));
781  }
782  h_min = min(h_min, attr(i));
783  h_max = max(h_max, attr(i));
784  }
785  cout << "h_min = " << h_min << ", h_max = " << h_max << endl;
786  }
787 
788  if (mk == 'k')
789  {
790  DenseMatrix J(dim);
791  for (int i = 0; i < mesh->GetNE(); i++)
792  {
793  int geom = mesh->GetElementBaseGeometry(i);
795  T->SetIntPoint(&Geometries.GetCenter(geom));
796  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
797  attr(i) = J.CalcSingularvalue(0) / J.CalcSingularvalue(dim-1);
798  }
799  }
800 
801  if (mk == 'p')
802  {
803  int *partitioning = NULL, np;
804  cout << "What type of partitioning?\n"
805  "c) Cartesian\n"
806  "s) Simple 1D split of the element sequence\n"
807  "0) METIS_PartGraphRecursive (sorted neighbor lists)\n"
808  "1) METIS_PartGraphKway (sorted neighbor lists)"
809  " (default)\n"
810  "2) METIS_PartGraphVKway (sorted neighbor lists)\n"
811  "3) METIS_PartGraphRecursive\n"
812  "4) METIS_PartGraphKway\n"
813  "5) METIS_PartGraphVKway\n"
814  "--> " << flush;
815  char pk;
816  cin >> pk;
817  if (pk == 'c')
818  {
819  int nxyz[3];
820  cout << "Enter nx: " << flush;
821  cin >> nxyz[0]; np = nxyz[0];
822  if (mesh->Dimension() > 1)
823  {
824  cout << "Enter ny: " << flush;
825  cin >> nxyz[1]; np *= nxyz[1];
826  if (mesh->Dimension() > 2)
827  {
828  cout << "Enter nz: " << flush;
829  cin >> nxyz[2]; np *= nxyz[2];
830  }
831  }
832  partitioning = mesh->CartesianPartitioning(nxyz);
833  }
834  else if (pk == 's')
835  {
836  cout << "Enter number of processors: " << flush;
837  cin >> np;
838 
839  partitioning = new int[mesh->GetNE()];
840  for (int i = 0; i < mesh->GetNE(); i++)
841  {
842  partitioning[i] = i * np / mesh->GetNE();
843  }
844  }
845  else
846  {
847  int part_method = pk - '0';
848  if (part_method < 0 || part_method > 5)
849  {
850  continue;
851  }
852  cout << "Enter number of processors: " << flush;
853  cin >> np;
854  partitioning = mesh->GeneratePartitioning(np, part_method);
855  }
856  if (partitioning)
857  {
858  const char part_file[] = "partitioning.txt";
859  ofstream opart(part_file);
860  opart << "number_of_elements " << mesh->GetNE() << '\n'
861  << "number_of_processors " << np << '\n';
862  for (int i = 0; i < mesh->GetNE(); i++)
863  {
864  opart << partitioning[i] << '\n';
865  }
866  cout << "Partitioning file: " << part_file << endl;
867 
868  Array<int> proc_el(np);
869  proc_el = 0;
870  for (int i = 0; i < mesh->GetNE(); i++)
871  {
872  proc_el[partitioning[i]]++;
873  }
874  int min_el = proc_el[0], max_el = proc_el[0];
875  for (int i = 1; i < np; i++)
876  {
877  if (min_el > proc_el[i])
878  {
879  min_el = proc_el[i];
880  }
881  if (max_el < proc_el[i])
882  {
883  max_el = proc_el[i];
884  }
885  }
886  cout << "Partitioning stats:\n"
887  << " "
888  << setw(12) << "minimum"
889  << setw(12) << "average"
890  << setw(12) << "maximum"
891  << setw(12) << "total" << '\n';
892  cout << " elements "
893  << setw(12) << min_el
894  << setw(12) << double(mesh->GetNE())/np
895  << setw(12) << max_el
896  << setw(12) << mesh->GetNE() << endl;
897  }
898  else
899  {
900  continue;
901  }
902 
903  for (int i = 0; i < mesh->GetNE(); i++)
904  {
905  attr(i) = part[i] = partitioning[i];
906  }
907  delete [] partitioning;
908  }
909 
910  char vishost[] = "localhost";
911  int visport = 19916;
912  socketstream sol_sock(vishost, visport);
913  if (sol_sock.is_open())
914  {
915  sol_sock.precision(14);
916  if (sdim == 2)
917  {
918  sol_sock << "fem2d_gf_data_keys\n";
919  if (mk != 'p')
920  {
921  mesh->Print(sol_sock);
922  }
923  else
924  {
925  // NURBS meshes do not support PrintWithPartitioning
926  if (mesh->NURBSext)
927  {
928  mesh->Print(sol_sock);
929  for (int i = 0; i < mesh->GetNE(); i++)
930  {
931  attr(i) = part[i];
932  }
933  }
934  else
935  {
936  mesh->PrintWithPartitioning(part, sol_sock, 1);
937  }
938  }
939  attr.Save(sol_sock);
940  sol_sock << "RjlmAb***********";
941  if (mk == 'v')
942  {
943  sol_sock << "e";
944  }
945  else
946  {
947  sol_sock << "\n";
948  }
949  }
950  else
951  {
952  sol_sock << "fem3d_gf_data_keys\n";
953  if (mk == 'v' || mk == 'h' || mk == 'k')
954  {
955  mesh->Print(sol_sock);
956  }
957  else if (mk == 'b')
958  {
959  bdr_mesh->Print(sol_sock);
960  bdr_attr.Save(sol_sock);
961  sol_sock << "mcaaA";
962  // Switch to a discrete color scale
963  sol_sock << "pppppp" << "pppppp" << "pppppp";
964  }
965  else
966  {
967  // NURBS meshes do not support PrintWithPartitioning
968  if (mesh->NURBSext)
969  {
970  mesh->Print(sol_sock);
971  for (int i = 0; i < mesh->GetNE(); i++)
972  {
973  attr(i) = part[i];
974  }
975  }
976  else
977  {
978  mesh->PrintWithPartitioning(part, sol_sock);
979  }
980  }
981  if (mk != 'b')
982  {
983  attr.Save(sol_sock);
984  sol_sock << "maaA";
985  if (mk == 'v')
986  {
987  sol_sock << "aa";
988  }
989  else
990  {
991  sol_sock << "\n";
992  }
993  }
994  }
995  sol_sock << flush;
996  }
997  else
998  {
999  cout << "Unable to connect to "
1000  << vishost << ':' << visport << endl;
1001  }
1002  delete attr_fespace;
1003  delete bdr_attr_fespace;
1004  }
1005 
1006  if (mk == 'l')
1007  {
1008  // Project and plot the function 'f'
1009  int p;
1010  FiniteElementCollection *fec = NULL;
1011  cout << "Enter projection space order: " << flush;
1012  cin >> p;
1013  if (p >= 1)
1014  {
1015  fec = new H1_FECollection(p, mesh->Dimension(),
1017  }
1018  else
1019  {
1020  fec = new DG_FECollection(-p, mesh->Dimension(),
1022  }
1023  FiniteElementSpace fes(mesh, fec);
1024  GridFunction level(&fes);
1025  FunctionCoefficient coeff(f);
1026  level.ProjectCoefficient(coeff);
1027  char vishost[] = "localhost";
1028  int visport = 19916;
1029  socketstream sol_sock(vishost, visport);
1030  if (sol_sock.is_open())
1031  {
1032  sol_sock.precision(14);
1033  sol_sock << "solution\n" << *mesh << level << flush;
1034  }
1035  else
1036  {
1037  cout << "Unable to connect to "
1038  << vishost << ':' << visport << endl;
1039  }
1040  delete fec;
1041  }
1042 
1043  if (mk == 'S')
1044  {
1045  const char mesh_file[] = "mesh-explorer.mesh";
1046  ofstream omesh(mesh_file);
1047  omesh.precision(14);
1048  mesh->Print(omesh);
1049  cout << "New mesh file: " << mesh_file << endl;
1050  }
1051 
1052  if (mk == 'V')
1053  {
1054  const char mesh_file[] = "mesh-explorer.vtk";
1055  ofstream omesh(mesh_file);
1056  omesh.precision(14);
1057  mesh->PrintVTK(omesh);
1058  cout << "New VTK mesh file: " << mesh_file << endl;
1059  }
1060 
1061 #ifdef MFEM_USE_ZLIB
1062  if (mk == 'Z')
1063  {
1064  const char mesh_file[] = "mesh-explorer.mesh.gz";
1065  ofgzstream omesh(mesh_file, "zwb9");
1066  omesh.precision(14);
1067  mesh->Print(omesh);
1068  cout << "New mesh file: " << mesh_file << endl;
1069  }
1070 #endif
1071 
1072  }
1073 
1074  delete bdr_attr_fec;
1075  delete attr_fec;
1076  delete bdr_mesh;
1077  delete mesh;
1078  return 0;
1079 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:245
bool Help() const
Return true if we are flagged to print the help message.
Definition: optparser.hpp:148
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:1609
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:521
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1234
int * CartesianPartitioning(int nxyz[])
Definition: mesh.cpp:5693
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:397
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:794
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:144
Mesh * read_par_mesh(int np, const char *mesh_prefix)
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1293
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:217
static const int NumGeom
Definition: geom.hpp:41
void JacToPerfJac(int GeomType, const DenseMatrix &J, DenseMatrix &PJ) const
Definition: geom.cpp:734
int * GeneratePartitioning(int nparts, int part_method=1)
Definition: mesh.cpp:5732
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:740
void PrintWithPartitioning(int *partitioning, std::ostream &out, 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:9594
void MakeOwner(FiniteElementCollection *_fec)
Make the GridFunction the owner of fec and fes.
Definition: gridfunc.hpp:115
const Geometry::Type geom
Definition: ex1.cpp:40
double Det() const
Definition: densemat.cpp:451
void NewNodes(GridFunction &nodes, bool make_owner=false)
Replace the internal node GridFunction with the given GridFunction.
Definition: mesh.cpp:6627
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:71
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:85
int AddTriangle(int v1, int v2, int v3, int attr=1)
Definition: mesh.cpp:1279
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:160
void Transform(void(*f)(const Vector &, Vector &))
Definition: mesh.cpp:10292
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
void PrintHelp(std::ostream &out) const
Print the help message.
Definition: optparser.cpp:378
double kappa
Definition: ex24.cpp:54
void Randomize(int seed=0)
Set random values in the vector.
Definition: vector.cpp:725
const IntegrationPoint & GetCenter(int GeomType)
Return the center of the given Geometry::Type, GeomType.
Definition: geom.hpp:70
int main(int argc, char *argv[])
Definition: ex1.cpp:66
Nodes: x_i = i/(n-1), i=0,...,n-1.
Definition: fe.hpp:37
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:834
void DeleteAll()
Delete the whole array.
Definition: array.hpp:821
Mesh * skin_mesh(Mesh *mesh)
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3417
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:427
Geometry Geometries
Definition: fe.cpp:12850
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
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:436
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:726
constexpr char vishost[]
void transformation(const Vector &p, Vector &v)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:8382
double region_eps
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 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
GeometryRefiner GlobGeometryRefiner
Definition: geom.cpp:1520
int GetAttribute() const
Return element&#39;s attribute.
Definition: element.hpp:55
void GetHilbertElementOrdering(Array< int > &ordering)
Definition: mesh.cpp:1776
Piecewise-constant discontinuous finite elements in 3D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:836
IntegrationRule RefPts
Definition: geom.hpp:243
RefinedGeometry * Refine(Geometry::Type Geom, int Times, int ETimes=1)
Definition: geom.cpp:942
Piecewise-constant discontinuous finite elements in 2D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:658
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:595
int Dimension() const
Definition: mesh.hpp:788
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
Definition: fe_coll.hpp:269
double GetElementSize(int i, int type=0)
Get the size of the i-th element relative to the perfect reference element.
Definition: mesh.cpp:74
int SpaceDimension() const
Definition: mesh.hpp:789
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:201
bool is_open()
True if the socketstream is open, false otherwise.
void PrintVTK(std::ostream &out)
Definition: mesh.cpp:8853
int AddSegment(int v1, int v2, int attr=1)
Definition: mesh.cpp:1265
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
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:2469
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:654
void GetColumnReference(int c, Vector &col)
Definition: densemat.hpp:270
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
Definition: fe_coll.cpp:47
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:203
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:734
void PrintError(std::ostream &out) const
Print the error message.
Definition: optparser.cpp:335
virtual const char * Name() const
Definition: fe_coll.hpp:61
A standard isoparametric element transformation.
Definition: eltrans.hpp:348
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:336
double Max() const
Returns the maximal element of the vector.
Definition: vector.cpp:832
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:4160
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition: mesh.cpp:2552
void ReorderElements(const Array< int > &ordering, bool reorder_vertices=true)
Definition: mesh.cpp:1828
int dim
Definition: ex24.cpp:53
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:304
double CalcSingularvalue(const int i) const
Return the i-th singular value (decreasing order) of NxN matrix, N=1,2,3.
Definition: densemat.cpp:1226
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:45
virtual int GetNVertices() const =0
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:2252
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition: gridfunc.cpp:180
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:10561
Vector data type.
Definition: vector.hpp:51
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:414
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:6603
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:159
void GetElementColoring(Array< int > &colors, int el0=0)
Definition: mesh.cpp:9519
void GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:179
const DenseMatrix & GetPointMat() const
Return the stored point matrix.
Definition: eltrans.hpp:391
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:7977
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:1047
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:199
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:823
void PrintCharacteristics(Vector *Vh=NULL, Vector *Vk=NULL, std::ostream &out=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:231
double f(const Vector &p)
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:145