MFEM  v4.5.2
Finite element discretization library
mesh-explorer.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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 "../common/mfem-common.hpp"
36 #include <fstream>
37 #include <limits>
38 #include <cstdlib>
39 
40 using namespace mfem;
41 using namespace std;
42 
43 // This transformation can be applied to a mesh with the 't' menu option.
44 void transformation(const Vector &p, Vector &v)
45 {
46  // simple shear transformation
47  double s = 0.1;
48 
49  if (p.Size() == 3)
50  {
51  v(0) = p(0) + s*p(1) + s*p(2);
52  v(1) = p(1) + s*p(2) + s*p(0);
53  v(2) = p(2);
54  }
55  else if (p.Size() == 2)
56  {
57  v(0) = p(0) + s*p(1);
58  v(1) = p(1) + s*p(0);
59  }
60  else
61  {
62  v = p;
63  }
64 }
65 
66 // This function is used with the 'r' menu option, sub-option 'l' to refine a
67 // mesh locally in a region, defined by return values <= region_eps.
68 double region_eps = 1e-8;
69 double region(const Vector &p)
70 {
71  const double x = p(0), y = p(1);
72  // here we describe the region: (x <= 1/4) && (y >= 0) && (y <= 1)
73  return std::max(std::max(x - 0.25, -y), y - 1.0);
74 }
75 
76 // The projection of this function can be plotted with the 'l' menu option
77 double f(const Vector &p)
78 {
79  double x = p(0);
80  double y = p.Size() > 1 ? p(1) : 0.0;
81  double z = p.Size() > 2 ? p(2) : 0.0;
82 
83  if (1)
84  {
85  // torus in the xy-plane
86  const double r_big = 2.0;
87  const double r_small = 1.0;
88  return hypot(r_big - hypot(x, y), z) - r_small;
89  }
90  if (0)
91  {
92  // sphere at the origin:
93  const double r = 1.0;
94  return hypot(hypot(x, y), z) - r;
95  }
96 }
97 
98 Mesh *read_par_mesh(int np, const char *mesh_prefix, Array<int>& partitioning,
99  Array<int>& bdr_partitioning)
100 {
101  Mesh *mesh;
102  Array<Mesh *> mesh_array;
103 
104  mesh_array.SetSize(np);
105  for (int p = 0; p < np; p++)
106  {
107  ostringstream fname;
108  fname << mesh_prefix << '.' << setfill('0') << setw(6) << p;
109  ifgzstream meshin(fname.str().c_str());
110  if (!meshin)
111  {
112  cerr << "Can not open mesh file: " << fname.str().c_str()
113  << '!' << endl;
114  for (p--; p >= 0; p--)
115  {
116  delete mesh_array[p];
117  }
118  return NULL;
119  }
120  mesh_array[p] = new Mesh(meshin, 1, 0);
121  // Assign corresponding processor number to element + boundary partitions
122  for (int i = 0; i < mesh_array[p]->GetNE(); i++)
123  {
124  partitioning.Append(p);
125  }
126  for (int i = 0; i < mesh_array[p]->GetNBE(); i++)
127  {
128  bdr_partitioning.Append(p);
129  }
130  }
131  mesh = new Mesh(mesh_array, np);
132 
133  for (int p = 0; p < np; p++)
134  {
135  delete mesh_array[np-1-p];
136  }
137  mesh_array.DeleteAll();
138 
139  return mesh;
140 }
141 
142 // Given a 3D mesh, produce a 2D mesh consisting of its boundary elements.
143 // We guarantee that the skin preserves the boundary index order.
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 void recover_bdr_partitioning(const Mesh* mesh, const Array<int>& partitioning,
253  Array<int>& bdr_partitioning)
254 {
255  bdr_partitioning.SetSize(mesh->GetNBE());
256  int info, e;
257  for (int be = 0; be < mesh->GetNBE(); be++)
258  {
259  mesh->GetBdrElementAdjacentElement(be, e, info);
260  bdr_partitioning[be] = partitioning[e];
261  }
262 }
263 
264 int main (int argc, char *argv[])
265 {
266  int np = 0;
267  const char *mesh_file = "../../data/beam-hex.mesh";
268  bool refine = true;
269 
270  OptionsParser args(argc, argv);
271  args.AddOption(&mesh_file, "-m", "--mesh",
272  "Mesh file to visualize.");
273  args.AddOption(&np, "-np", "--num-proc",
274  "Load mesh from multiple processors.");
275  args.AddOption(&refine, "-ref", "--refinement", "-no-ref", "--no-refinement",
276  "Prepare the mesh for refinement or not.");
277  args.Parse();
278  if (!args.Good())
279  {
280  if (!args.Help())
281  {
282  args.PrintError(cout);
283  cout << endl;
284  }
285  cout << "Visualize and manipulate a serial mesh:\n"
286  << " mesh-explorer -m <mesh_file>\n"
287  << "Visualize and manipulate a parallel mesh:\n"
288  << " mesh-explorer -np <#proc> -m <mesh_prefix>\n" << endl
289  << "All Options:\n";
290  args.PrintHelp(cout);
291  return 1;
292  }
293  args.PrintOptions(cout);
294 
295  Mesh *mesh;
296  Mesh *bdr_mesh = NULL;
297 
298  // Helper to distinguish whether we use a parallel or serial mesh.
299  const bool use_par_mesh = np > 0;
300 
301  // Helper for visualizing the partitioning.
302  Array<int> partitioning;
303  Array<int> bdr_partitioning;
304  if (!use_par_mesh)
305  {
306  mesh = new Mesh(mesh_file, 1, refine);
307  partitioning.SetSize(mesh->GetNE());
308  partitioning = 0;
309  bdr_partitioning.SetSize(mesh->GetNBE());
310  bdr_partitioning = 0;
311  }
312  else
313  {
314  mesh = read_par_mesh(np, mesh_file, partitioning, bdr_partitioning);
315  if (mesh == NULL)
316  {
317  return 3;
318  }
319  }
320  int dim = mesh->Dimension();
321  int sdim = mesh->SpaceDimension();
322 
323  FiniteElementCollection *bdr_attr_fec = NULL;
324  FiniteElementCollection *attr_fec;
325  if (dim == 2)
326  {
327  attr_fec = new Const2DFECollection;
328  }
329  else
330  {
331  bdr_attr_fec = new Const2DFECollection;
332  attr_fec = new Const3DFECollection;
333  }
334 
335  int print_char = 1;
336  while (1)
337  {
338  if (print_char)
339  {
340  cout << endl;
341  mesh->PrintCharacteristics();
342  cout << "boundary attribs :";
343  for (int i = 0; i < mesh->bdr_attributes.Size(); i++)
344  {
345  cout << ' ' << mesh->bdr_attributes[i];
346  }
347  cout << '\n' << "material attribs :";
348  for (int i = 0; i < mesh->attributes.Size(); i++)
349  {
350  cout << ' ' << mesh->attributes[i];
351  }
352  cout << endl;
353  cout << "mesh curvature : ";
354  if (mesh->GetNodalFESpace() != NULL)
355  {
356  cout << mesh->GetNodalFESpace()->FEColl()->Name() << endl;
357  }
358  else
359  {
360  cout << "NONE" << endl;
361  }
362  }
363  print_char = 0;
364  cout << endl;
365  cout << "What would you like to do?\n"
366  "r) Refine\n"
367  "c) Change curvature\n"
368  "s) Scale\n"
369  "t) Transform\n"
370  "j) Jitter\n"
371  "v) View mesh\n"
372  "P) View partitioning\n"
373  "m) View materials\n"
374  "b) View boundary\n"
375  "B) View boundary partitioning\n"
376  "e) View elements\n"
377  "h) View element sizes, h\n"
378  "k) View element ratios, kappa\n"
379  "J) View scaled Jacobian\n"
380  "l) Plot a function\n"
381  "x) Print sub-element stats\n"
382  "f) Find physical point in reference space\n"
383  "p) Generate a partitioning\n"
384  "o) Reorder elements\n"
385  "S) Save in MFEM format\n"
386  "V) Save in VTK format (only linear and quadratic meshes)\n"
387  "D) Save as a DataCollection\n"
388  "q) Quit\n"
389 #ifdef MFEM_USE_ZLIB
390  "Z) Save in MFEM format with compression\n"
391 #endif
392  "--> " << flush;
393  char mk;
394  cin >> mk;
395  if (!cin) { break; }
396 
397  if (mk == 'q')
398  {
399  break;
400  }
401 
402  if (mk == 'r')
403  {
404  cout <<
405  "Choose type of refinement:\n"
406  "s) standard refinement with Mesh::UniformRefinement()\n"
407  "b) Mesh::UniformRefinement() (bisection for tet meshes)\n"
408  "u) uniform refinement with a factor\n"
409  "g) non-uniform refinement (Gauss-Lobatto) with a factor\n"
410  "l) refine locally using the region() function\n"
411  "--> " << flush;
412  char sk;
413  cin >> sk;
414  switch (sk)
415  {
416  case 's':
417  mesh->UniformRefinement();
418  // Make sure tet-only meshes are marked for local refinement.
419  mesh->Finalize(true);
420  break;
421  case 'b':
422  mesh->UniformRefinement(1); // ref_algo = 1
423  break;
424  case 'u':
425  case 'g':
426  {
427  cout << "enter refinement factor --> " << flush;
428  int ref_factor;
429  cin >> ref_factor;
430  if (ref_factor <= 1 || ref_factor > 32) { break; }
431  int ref_type = (sk == 'u') ? BasisType::ClosedUniform :
433  *mesh = Mesh::MakeRefined(*mesh, ref_factor, ref_type);
434  break;
435  }
436  case 'l':
437  {
438  Vector pt;
439  Array<int> marked_elements;
440  for (int i = 0; i < mesh->GetNE(); i++)
441  {
442  // check all nodes of the element
444  mesh->GetElementTransformation(i, &T);
445  for (int j = 0; j < T.GetPointMat().Width(); j++)
446  {
447  T.GetPointMat().GetColumnReference(j, pt);
448  if (region(pt) <= region_eps)
449  {
450  marked_elements.Append(i);
451  break;
452  }
453  }
454  }
455  mesh->GeneralRefinement(marked_elements);
456  break;
457  }
458  }
459  print_char = 1;
460  }
461 
462  if (mk == 'c')
463  {
464  int p;
465  cout << "enter new order for mesh curvature --> " << flush;
466  cin >> p;
467  mesh->SetCurvature(p > 0 ? p : -p, p <= 0);
468  print_char = 1;
469  }
470 
471  if (mk == 's')
472  {
473  double factor;
474  cout << "scaling factor ---> " << flush;
475  cin >> factor;
476 
477  GridFunction *nodes = mesh->GetNodes();
478  if (nodes == NULL)
479  {
480  for (int i = 0; i < mesh->GetNV(); i++)
481  {
482  double *v = mesh->GetVertex(i);
483  v[0] *= factor;
484  v[1] *= factor;
485  if (dim == 3)
486  {
487  v[2] *= factor;
488  }
489  }
490  }
491  else
492  {
493  *nodes *= factor;
494  }
495 
496  print_char = 1;
497  }
498 
499  if (mk == 't')
500  {
501  char type;
502  cout << "Choose a transformation:\n"
503  "u) User-defined transform through mesh-explorer::transformation()\n"
504  "k) Kershaw transform\n"<< "---> " << flush;
505  cin >> type;
506  if (type == 'u')
507  {
508  mesh->Transform(transformation);
509  }
510  else if (type == 'k')
511  {
512  cout << "Note: For Kershaw transformation, the input must be "
513  "Cartesian aligned with nx multiple of 6 and "
514  "both ny and nz multiples of 2."
515  "Kershaw transform works for 2D meshes also.\n" << flush;
516 
517  double epsy, epsz = 0.0;
518  cout << "Kershaw transform factor, epsy in (0, 1]) ---> " << flush;
519  cin >> epsy;
520  if (mesh->Dimension() == 3)
521  {
522  cout << "Kershaw transform factor, epsz in (0, 1]) ---> " << flush;
523  cin >> epsz;
524  }
525  common::KershawTransformation kershawT(mesh->Dimension(), epsy, epsz);
526  mesh->Transform(kershawT);
527  }
528  else
529  {
530  MFEM_ABORT("Transformation type not supported.");
531  }
532  print_char = 1;
533  }
534 
535  if (mk == 'j')
536  {
537  double jitter;
538  cout << "jitter factor ---> " << flush;
539  cin >> jitter;
540 
541  GridFunction *nodes = mesh->GetNodes();
542 
543  if (nodes == NULL)
544  {
545  cerr << "The mesh should have nodes, introduce curvature first!\n";
546  }
547  else
548  {
549  FiniteElementSpace *fespace = nodes->FESpace();
550 
551  GridFunction rdm(fespace);
552  rdm.Randomize();
553  rdm -= 0.5; // shift to random values in [-0.5,0.5]
554  rdm *= jitter;
555 
556  // compute minimal local mesh size
557  Vector h0(fespace->GetNDofs());
558  h0 = infinity();
559  {
560  Array<int> dofs;
561  for (int i = 0; i < fespace->GetNE(); i++)
562  {
563  fespace->GetElementDofs(i, dofs);
564  for (int j = 0; j < dofs.Size(); j++)
565  {
566  h0(dofs[j]) = std::min(h0(dofs[j]), mesh->GetElementSize(i));
567  }
568  }
569  }
570 
571  // scale the random values to be of order of the local mesh size
572  for (int i = 0; i < fespace->GetNDofs(); i++)
573  {
574  for (int d = 0; d < dim; d++)
575  {
576  rdm(fespace->DofToVDof(i,d)) *= h0(i);
577  }
578  }
579 
580  char move_bdr = 'n';
581  cout << "move boundary nodes? [y/n] ---> " << flush;
582  cin >> move_bdr;
583 
584  // don't perturb the boundary
585  if (move_bdr == 'n')
586  {
587  Array<int> vdofs;
588  for (int i = 0; i < fespace->GetNBE(); i++)
589  {
590  fespace->GetBdrElementVDofs(i, vdofs);
591  for (int j = 0; j < vdofs.Size(); j++)
592  {
593  rdm(vdofs[j]) = 0.0;
594  }
595  }
596  }
597 
598  *nodes += rdm;
599  }
600 
601  print_char = 1;
602  }
603 
604  if (mk == 'x')
605  {
606  int sd, nz = 0;
607  DenseMatrix J(dim);
608  double min_det_J, max_det_J, min_det_J_z, max_det_J_z;
609  double min_kappa, max_kappa, max_ratio_det_J_z;
610  min_det_J = min_kappa = infinity();
611  max_det_J = max_kappa = max_ratio_det_J_z = -infinity();
612  cout << "subdivision factor ---> " << flush;
613  cin >> sd;
614  Array<int> bad_elems_by_geom(Geometry::NumGeom);
615  bad_elems_by_geom = 0;
616  // Only print so many to keep output compact
617  const int max_to_print = 10;
618  for (int i = 0; i < mesh->GetNE(); i++)
619  {
620  Geometry::Type geom = mesh->GetElementBaseGeometry(i);
622 
623  RefinedGeometry *RefG = GlobGeometryRefiner.Refine(geom, sd, 1);
624  IntegrationRule &ir = RefG->RefPts;
625 
626  min_det_J_z = infinity();
627  max_det_J_z = -infinity();
628  for (int j = 0; j < ir.GetNPoints(); j++)
629  {
630  T->SetIntPoint(&ir.IntPoint(j));
631  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
632 
633  double det_J = J.Det();
634  double kappa =
636 
637  min_det_J_z = fmin(min_det_J_z, det_J);
638  max_det_J_z = fmax(max_det_J_z, det_J);
639 
640  min_kappa = fmin(min_kappa, kappa);
641  max_kappa = fmax(max_kappa, kappa);
642  }
643  max_ratio_det_J_z =
644  fmax(max_ratio_det_J_z, max_det_J_z/min_det_J_z);
645  min_det_J = fmin(min_det_J, min_det_J_z);
646  max_det_J = fmax(max_det_J, max_det_J_z);
647  if (min_det_J_z <= 0.0)
648  {
649  if (nz < max_to_print)
650  {
651  Vector center;
652  mesh->GetElementCenter(i, center);
653  cout << "det(J) < 0 = " << min_det_J_z << " in element "
654  << i << ", centered at: ";
655  center.Print();
656  }
657  nz++;
658  bad_elems_by_geom[geom]++;
659  }
660  }
661  if (nz >= max_to_print)
662  {
663  cout << "det(J) < 0 for " << nz - max_to_print << " more elements "
664  << "not printed.\n";
665  }
666  cout << "\nbad elements = " << nz;
667  if (nz)
668  {
669  cout << " -- ";
670  Mesh::PrintElementsByGeometry(dim, bad_elems_by_geom, cout);
671  }
672  cout << "\nmin det(J) = " << min_det_J
673  << "\nmax det(J) = " << max_det_J
674  << "\nglobal ratio = " << max_det_J/min_det_J
675  << "\nmax el ratio = " << max_ratio_det_J_z
676  << "\nmin kappa = " << min_kappa
677  << "\nmax kappa = " << max_kappa << endl;
678  }
679 
680  if (mk == 'f')
681  {
682  DenseMatrix point_mat(sdim,1);
683  cout << "\npoint in physical space ---> " << flush;
684  for (int i = 0; i < sdim; i++)
685  {
686  cin >> point_mat(i,0);
687  }
688  Array<int> elem_ids;
690 
691  // physical -> reference space
692  mesh->FindPoints(point_mat, elem_ids, ips);
693 
694  cout << "point in reference space:";
695  if (elem_ids[0] == -1)
696  {
697  cout << " NOT FOUND!\n";
698  }
699  else
700  {
701  cout << " element " << elem_ids[0] << ", ip =";
702  cout << " " << ips[0].x;
703  if (sdim > 1)
704  {
705  cout << " " << ips[0].y;
706  if (sdim > 2)
707  {
708  cout << " " << ips[0].z;
709  }
710  }
711  cout << endl;
712  }
713  }
714 
715  if (mk == 'o')
716  {
717  cout << "What type of reordering?\n"
718  "g) Gecko edge-product minimization\n"
719  "h) Hilbert spatial sort\n"
720  "--> " << flush;
721  char rk;
722  cin >> rk;
723 
724  Array<int> ordering, tentative;
725  if (rk == 'h')
726  {
727  mesh->GetHilbertElementOrdering(ordering);
728  mesh->ReorderElements(ordering);
729  }
730  else if (rk == 'g')
731  {
732  int outer, inner, window, period;
733  cout << "Enter number of outer iterations (default 5): " << flush;
734  cin >> outer;
735  cout << "Enter number of inner iterations (default 4): " << flush;
736  cin >> inner;
737  cout << "Enter window size (default 4, beware of exponential cost): "
738  << flush;
739  cin >> window;
740  cout << "Enter period for window size increment (default 2): "
741  << flush;
742  cin >> period;
743 
744  double best_cost = infinity();
745  for (int i = 0; i < outer; i++)
746  {
747  int seed = i+1;
748  double cost = mesh->GetGeckoElementOrdering(
749  tentative, inner, window, period, seed, true);
750 
751  if (cost < best_cost)
752  {
753  ordering = tentative;
754  best_cost = cost;
755  }
756  }
757  cout << "Final cost: " << best_cost << endl;
758 
759  mesh->ReorderElements(ordering);
760  }
761  }
762 
763  // These are most of the cases that open a new GLVis window
764  if (mk == 'm' || mk == 'b' || mk == 'e' || mk == 'v' || mk == 'h' ||
765  mk == 'k' || mk == 'J' || mk == 'p' || mk == 'B' || mk == 'P')
766  {
767  FiniteElementSpace *bdr_attr_fespace = NULL;
768  FiniteElementSpace *attr_fespace =
769  new FiniteElementSpace(mesh, attr_fec);
770  GridFunction bdr_attr;
771  GridFunction attr(attr_fespace);
772 
773  if (mk == 'm')
774  {
775  for (int i = 0; i < mesh->GetNE(); i++)
776  {
777  attr(i) = mesh->GetAttribute(i);
778  }
779  }
780 
781  if (mk == 'P')
782  {
783  for (int i = 0; i < mesh->GetNE(); i++)
784  {
785  attr(i) = partitioning[i] + 1;
786  }
787  }
788 
789  if (mk == 'b' || mk == 'B')
790  {
791  if (dim == 3)
792  {
793  delete bdr_mesh;
794  bdr_mesh = skin_mesh(mesh);
795  bdr_attr_fespace =
796  new FiniteElementSpace(bdr_mesh, bdr_attr_fec);
797  bdr_attr.SetSpace(bdr_attr_fespace);
798  if (mk == 'b')
799  {
800  for (int i = 0; i < bdr_mesh->GetNE(); i++)
801  {
802  bdr_attr(i) = bdr_mesh->GetAttribute(i);
803  }
804  }
805  else if (mk == 'B')
806  {
807  for (int i = 0; i < bdr_mesh->GetNE(); i++)
808  {
809  bdr_attr(i) = bdr_partitioning[i] + 1;
810  }
811  }
812  else
813  {
814  MFEM_WARNING("Unimplemented case.");
815  }
816  }
817  else
818  {
819  MFEM_WARNING("Unsupported mesh dimension.");
820  attr = 1.0;
821  }
822  }
823 
824  if (mk == 'v')
825  {
826  attr = 1.0;
827  }
828 
829  if (mk == 'e')
830  {
831  Array<int> coloring;
832  srand(time(0));
833  double a = double(rand()) / (double(RAND_MAX) + 1.);
834  int el0 = (int)floor(a * mesh->GetNE());
835  cout << "Generating coloring starting with element " << el0+1
836  << " / " << mesh->GetNE() << endl;
837  mesh->GetElementColoring(coloring, el0);
838  for (int i = 0; i < coloring.Size(); i++)
839  {
840  attr(i) = coloring[i];
841  }
842  cout << "Number of colors: " << attr.Max() + 1 << endl;
843  for (int i = 0; i < mesh->GetNE(); i++)
844  {
845  attr(i) = i; // coloring by element number
846  }
847  }
848 
849  if (mk == 'h')
850  {
851  DenseMatrix J(dim);
852  double h_min, h_max;
853  h_min = infinity();
854  h_max = -h_min;
855  for (int i = 0; i < mesh->GetNE(); i++)
856  {
857  Geometry::Type geom = mesh->GetElementBaseGeometry(i);
859  T->SetIntPoint(&Geometries.GetCenter(geom));
860  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
861 
862  attr(i) = J.Det();
863  if (attr(i) < 0.0)
864  {
865  attr(i) = -pow(-attr(i), 1.0/double(dim));
866  }
867  else
868  {
869  attr(i) = pow(attr(i), 1.0/double(dim));
870  }
871  h_min = min(h_min, attr(i));
872  h_max = max(h_max, attr(i));
873  }
874  cout << "h_min = " << h_min << ", h_max = " << h_max << endl;
875  }
876 
877  if (mk == 'k')
878  {
879  DenseMatrix J(dim);
880  for (int i = 0; i < mesh->GetNE(); i++)
881  {
882  Geometry::Type geom = mesh->GetElementBaseGeometry(i);
884  T->SetIntPoint(&Geometries.GetCenter(geom));
885  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
886  attr(i) = J.CalcSingularvalue(0) / J.CalcSingularvalue(dim-1);
887  }
888  }
889 
890  if (mk == 'J')
891  {
892  // The "scaled Jacobian" is the determinant of the Jacobian scaled
893  // by the l2 norms of its columns. It can be used to identify badly
894  // skewed elements, since it takes values between 0 and 1, with 0
895  // corresponding to a flat element, and 1 to orthogonal columns.
896  DenseMatrix J(dim);
897  int sd;
898  cout << "subdivision factor ---> " << flush;
899  cin >> sd;
900  for (int i = 0; i < mesh->GetNE(); i++)
901  {
902  Geometry::Type geom = mesh->GetElementBaseGeometry(i);
904 
905  RefinedGeometry *RefG = GlobGeometryRefiner.Refine(geom, sd, 1);
906  IntegrationRule &ir = RefG->RefPts;
907 
908  // For each element, find the minimal scaled Jacobian in a
909  // lattice of points with the given subdivision factor.
910  attr(i) = infinity();
911  for (int j = 0; j < ir.GetNPoints(); j++)
912  {
913  T->SetIntPoint(&ir.IntPoint(j));
914  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
915 
916  // Jacobian determinant
917  double sJ = J.Det();
918 
919  for (int k = 0; k < J.Width(); k++)
920  {
921  Vector col;
922  J.GetColumnReference(k,col);
923  // Scale by column norms
924  sJ /= col.Norml2();
925  }
926 
927  attr(i) = fmin(sJ, attr(i));
928  }
929  }
930  }
931 
932  if (mk == 'p')
933  {
934  cout << "What type of partitioning?\n"
935  "c) Cartesian\n"
936  "s) Simple 1D split of the element sequence\n"
937  "0) METIS_PartGraphRecursive (sorted neighbor lists)\n"
938  "1) METIS_PartGraphKway (sorted neighbor lists)"
939  " (default)\n"
940  "2) METIS_PartGraphVKway (sorted neighbor lists)\n"
941  "3) METIS_PartGraphRecursive\n"
942  "4) METIS_PartGraphKway\n"
943  "5) METIS_PartGraphVKway\n"
944  "--> " << flush;
945  char pk;
946  cin >> pk;
947  if (pk == 'c')
948  {
949  int nxyz[3];
950  cout << "Enter nx: " << flush;
951  cin >> nxyz[0]; np = nxyz[0];
952  if (mesh->Dimension() > 1)
953  {
954  cout << "Enter ny: " << flush;
955  cin >> nxyz[1]; np *= nxyz[1];
956  if (mesh->Dimension() > 2)
957  {
958  cout << "Enter nz: " << flush;
959  cin >> nxyz[2]; np *= nxyz[2];
960  }
961  }
962  int *part = mesh->CartesianPartitioning(nxyz);
963  partitioning = Array<int>(part, mesh->GetNE());
964  delete [] part;
965  recover_bdr_partitioning(mesh, partitioning, bdr_partitioning);
966  }
967  else if (pk == 's')
968  {
969  cout << "Enter number of processors: " << flush;
970  cin >> np;
971 
972  partitioning.SetSize(mesh->GetNE());
973  for (int i = 0; i < mesh->GetNE(); i++)
974  {
975  partitioning[i] = i * np / mesh->GetNE();
976  }
977  recover_bdr_partitioning(mesh, partitioning, bdr_partitioning);
978  }
979  else
980  {
981  int part_method = pk - '0';
982  if (part_method < 0 || part_method > 5)
983  {
984  continue;
985  }
986  cout << "Enter number of processors: " << flush;
987  cin >> np;
988  int *part = mesh->GeneratePartitioning(np, part_method);
989  partitioning = Array<int>(part, mesh->GetNE());
990  delete [] part;
991  recover_bdr_partitioning(mesh, partitioning, bdr_partitioning);
992  }
993  if (partitioning)
994  {
995  const char part_file[] = "partitioning.txt";
996  ofstream opart(part_file);
997  opart << "number_of_elements " << mesh->GetNE() << '\n'
998  << "number_of_processors " << np << '\n';
999  for (int i = 0; i < mesh->GetNE(); i++)
1000  {
1001  opart << partitioning[i] << '\n';
1002  }
1003  cout << "Partitioning file: " << part_file << endl;
1004 
1005  Array<int> proc_el(np);
1006  proc_el = 0;
1007  for (int i = 0; i < mesh->GetNE(); i++)
1008  {
1009  proc_el[partitioning[i]]++;
1010  }
1011  int min_el = proc_el[0], max_el = proc_el[0];
1012  for (int i = 1; i < np; i++)
1013  {
1014  if (min_el > proc_el[i])
1015  {
1016  min_el = proc_el[i];
1017  }
1018  if (max_el < proc_el[i])
1019  {
1020  max_el = proc_el[i];
1021  }
1022  }
1023  cout << "Partitioning stats:\n"
1024  << " "
1025  << setw(12) << "minimum"
1026  << setw(12) << "average"
1027  << setw(12) << "maximum"
1028  << setw(12) << "total" << '\n';
1029  cout << " elements "
1030  << setw(12) << min_el
1031  << setw(12) << double(mesh->GetNE())/np
1032  << setw(12) << max_el
1033  << setw(12) << mesh->GetNE() << endl;
1034  }
1035  else
1036  {
1037  continue;
1038  }
1039 
1040  for (int i = 0; i < mesh->GetNE(); i++)
1041  {
1042  attr(i) = partitioning[i] + 1;
1043  }
1044  }
1045 
1046  char vishost[] = "localhost";
1047  int visport = 19916;
1048  socketstream sol_sock(vishost, visport);
1049  if (sol_sock.is_open())
1050  {
1051  sol_sock.precision(14);
1052  if (sdim == 2)
1053  {
1054  sol_sock << "fem2d_gf_data_keys\n";
1055  if (mk != 'p')
1056  {
1057  mesh->Print(sol_sock);
1058  }
1059  else
1060  {
1061  // NURBS meshes do not support PrintWithPartitioning
1062  if (mesh->NURBSext)
1063  {
1064  mesh->Print(sol_sock);
1065  for (int i = 0; i < mesh->GetNE(); i++)
1066  {
1067  attr(i) = partitioning[i];
1068  }
1069  }
1070  else
1071  {
1072  mesh->PrintWithPartitioning(partitioning, sol_sock, 1);
1073  }
1074  }
1075  attr.Save(sol_sock);
1076  sol_sock << "RjlmAb***********";
1077  if (mk == 'v')
1078  {
1079  sol_sock << "e";
1080  }
1081  else
1082  {
1083  sol_sock << "\n";
1084  }
1085  }
1086  else
1087  {
1088  sol_sock << "fem3d_gf_data_keys\n";
1089  if (mk == 'v' || mk == 'h' || mk == 'k' || mk == 'J')
1090  {
1091  mesh->Print(sol_sock);
1092  }
1093  else if (mk == 'b' || mk == 'B')
1094  {
1095  bdr_mesh->Print(sol_sock);
1096  bdr_attr.Save(sol_sock);
1097  sol_sock << "mcaaA";
1098  // Switch to a discrete color scale
1099  sol_sock << "pppppp" << "pppppp" << "pppppp";
1100  }
1101  else
1102  {
1103  // NURBS meshes do not support PrintWithPartitioning
1104  if (mesh->NURBSext)
1105  {
1106  mesh->Print(sol_sock);
1107  for (int i = 0; i < mesh->GetNE(); i++)
1108  {
1109  attr(i) = partitioning[i];
1110  }
1111  }
1112  else
1113  {
1114  mesh->PrintWithPartitioning(partitioning, sol_sock);
1115  }
1116  }
1117  if (mk != 'b' && mk != 'B')
1118  {
1119  attr.Save(sol_sock);
1120  sol_sock << "maaA";
1121  if (mk == 'v')
1122  {
1123  sol_sock << "aa";
1124  }
1125  else
1126  {
1127  sol_sock << "\n";
1128  }
1129  }
1130  }
1131  sol_sock << flush;
1132  }
1133  else
1134  {
1135  cout << "Unable to connect to "
1136  << vishost << ':' << visport << endl;
1137  }
1138  delete attr_fespace;
1139  delete bdr_attr_fespace;
1140  }
1141 
1142  if (mk == 'l')
1143  {
1144  // Project and plot the function 'f'
1145  int p;
1146  FiniteElementCollection *fec = NULL;
1147  cout << "Enter projection space order: " << flush;
1148  cin >> p;
1149  if (p >= 1)
1150  {
1151  fec = new H1_FECollection(p, mesh->Dimension(),
1153  }
1154  else
1155  {
1156  fec = new DG_FECollection(-p, mesh->Dimension(),
1158  }
1159  FiniteElementSpace fes(mesh, fec);
1160  GridFunction level(&fes);
1161  FunctionCoefficient coeff(f);
1162  level.ProjectCoefficient(coeff);
1163  char vishost[] = "localhost";
1164  int visport = 19916;
1165  socketstream sol_sock(vishost, visport);
1166  if (sol_sock.is_open())
1167  {
1168  sol_sock.precision(14);
1169  sol_sock << "solution\n" << *mesh << level << flush;
1170  }
1171  else
1172  {
1173  cout << "Unable to connect to "
1174  << vishost << ':' << visport << endl;
1175  }
1176  delete fec;
1177  }
1178 
1179  if (mk == 'S')
1180  {
1181  const char omesh_file[] = "mesh-explorer.mesh";
1182  ofstream omesh(omesh_file);
1183  omesh.precision(14);
1184  mesh->Print(omesh);
1185  cout << "New mesh file: " << omesh_file << endl;
1186  }
1187 
1188  if (mk == 'V')
1189  {
1190  const char omesh_file[] = "mesh-explorer.vtk";
1191  ofstream omesh(omesh_file);
1192  omesh.precision(14);
1193  mesh->PrintVTK(omesh);
1194  cout << "New VTK mesh file: " << omesh_file << endl;
1195  }
1196 
1197  if (mk == 'D')
1198  {
1199  cout << "What type of DataCollection?\n"
1200  "p) ParaView Data Collection\n"
1201  "v) VisIt Data Collection\n"
1202  "--> " << flush;
1203  char dk;
1204  cin >> dk;
1205  if (dk == 'p' || dk == 'P')
1206  {
1207  const char omesh_file[] = "mesh-explorer-paraview";
1208  ParaViewDataCollection dc(omesh_file, mesh);
1209  if (mesh->GetNodes())
1210  {
1211  int order = mesh->GetNodes()->FESpace()->GetMaxElementOrder();
1212  if (order > 1)
1213  {
1214  dc.SetHighOrderOutput(true);
1215  dc.SetLevelsOfDetail(order);
1216  }
1217  }
1218  dc.Save();
1219  cout << "New ParaView mesh file: " << omesh_file << endl;
1220  }
1221  else if (dk == 'v' || dk == 'V')
1222  {
1223  const char omesh_file[] = "mesh-explorer-visit";
1224  VisItDataCollection dc(omesh_file, mesh);
1225  dc.SetPrecision(14);
1226  dc.Save();
1227  cout << "New VisIt mesh file: " << omesh_file << "_000000.mfem_root"
1228  << endl;
1229  }
1230  else
1231  {
1232  cout << "Unrecognized DataCollection type: \"" << dk << "\""
1233  << endl;
1234  }
1235  }
1236 
1237 #ifdef MFEM_USE_ZLIB
1238  if (mk == 'Z')
1239  {
1240  const char omesh_file[] = "mesh-explorer.mesh.gz";
1241  ofgzstream omesh(omesh_file, "zwb9");
1242  omesh.precision(14);
1243  mesh->Print(omesh);
1244  cout << "New mesh file: " << omesh_file << endl;
1245  }
1246 #endif
1247 
1248  }
1249 
1250  delete bdr_attr_fec;
1251  delete attr_fec;
1252  delete bdr_mesh;
1253  delete mesh;
1254  return 0;
1255 }
const char vishost[]
void SetPrecision(int prec)
Set the precision (number of digits) used for the text output of doubles.
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:2042
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:574
int * CartesianPartitioning(int nxyz[])
Definition: mesh.cpp:7030
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:5351
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:1686
virtual void GetVertices(Array< int > &v) const =0
Returns element&#39;s vertices.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:324
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:1108
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
int Dimension() const
Definition: mesh.hpp:1047
int * GeneratePartitioning(int nparts, int part_method=1)
Definition: mesh.cpp:7074
Helper class for ParaView visualization data.
void NewNodes(GridFunction &nodes, bool make_owner=false)
Replace the internal node GridFunction with the given GridFunction.
Definition: mesh.cpp:7976
void Print(std::ostream &out=mfem::out, int width=8) const
Prints vector to stream out.
Definition: vector.cpp:725
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
int GetAttribute() const
Return element&#39;s attribute.
Definition: element.hpp:55
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:3746
int AddTriangle(int v1, int v2, int v3, int attr=1)
Definition: mesh.cpp:1672
double Det() const
Definition: densemat.cpp:488
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
void Transform(void(*f)(const Vector &, Vector &))
Definition: mesh.cpp:11691
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:584
double CalcSingularvalue(const int i) const
Return the i-th singular value (decreasing order) of NxN matrix, N=1,2,3.
Definition: densemat.cpp:1268
void PrintError(std::ostream &out) const
Print the error message.
Definition: optparser.cpp:355
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
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150
void MakeOwner(FiniteElementCollection *fec_)
Make the GridFunction the owner of fec and fes.
Definition: gridfunc.hpp:122
STL namespace.
void Randomize(int seed=0)
Set random values in the vector.
Definition: vector.cpp:786
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)
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:939
Geometry Geometries
Definition: fe.cpp:49
void DeleteAll()
Delete the whole array.
Definition: array.hpp:851
Mesh * skin_mesh(Mesh *mesh)
bool Help() const
Return true if we are flagged to print the help message.
Definition: optparser.hpp:153
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:933
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:1544
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:1618
GeometryRefiner GlobGeometryRefiner
Definition: geom.cpp:1773
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
double GetElementSize(ElementTransformation *T, int type=0)
Definition: mesh.cpp:76
Geometry::Type GetGeometryType() const
Definition: element.hpp:52
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:10365
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:756
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:601
void transformation(const Vector &p, Vector &v)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9878
double region_eps
Data collection with VisIt I/O routines.
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
Definition: mesh.cpp:5356
const double * GetVertex(int i) const
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:1053
virtual const char * Name() const
Definition: fe_coll.hpp:73
void GetHilbertElementOrdering(Array< int > &ordering)
Definition: mesh.cpp:2209
Piecewise-constant discontinuous finite elements in 3D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:1107
void GetElementCenter(int i, Vector &center)
Definition: mesh.cpp:68
IntegrationRule RefPts
Definition: geom.hpp:314
void SetHighOrderOutput(bool high_order_output_)
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:929
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:683
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:614
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
Definition: fe_coll.hpp:370
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:275
bool is_open()
True if the socketstream is open, false otherwise.
int AddSegment(int v1, int v2, int attr=1)
Definition: mesh.cpp:1658
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
virtual void Save() override
Save the collection and a VisIt root file.
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:2902
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:684
void GetColumnReference(int c, Vector &col)
Definition: densemat.hpp:312
void JacToPerfJac(int GeomType, const DenseMatrix &J, DenseMatrix &PJ) const
Definition: geom.cpp:867
virtual DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element &#39;elem&#39;.
Definition: fespace.cpp:2680
double Max() const
Returns the maximal element of the vector.
Definition: vector.cpp:894
int main(int argc, char *argv[])
int SpaceDimension() const
Definition: mesh.hpp:1048
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
Definition: fe_coll.cpp:116
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:936
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:277
A standard isoparametric element transformation.
Definition: eltrans.hpp:361
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:348
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 Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition: mesh.cpp:2995
void ReorderElements(const Array< int > &ordering, bool reorder_vertices=true)
Definition: mesh.cpp:2261
int dim
Definition: ex24.cpp:53
void PrintHelp(std::ostream &out) const
Print the help message.
Definition: optparser.cpp:398
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:623
double infinity()
Define a shortcut for std::numeric_limits<double>::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:2396
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:804
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition: gridfunc.cpp:196
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
void SetLevelsOfDetail(int levels_of_detail_)
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:11962
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:251
Nodes: x_i = i/(n-1), i=0,...,n-1.
Definition: fe_base.hpp:34
Vector data type.
Definition: vector.hpp:60
virtual void Print(std::ostream &os=mfem::out) const
Definition: mesh.hpp:1749
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:252
void GetElementColoring(Array< int > &colors, int el0=0)
Definition: mesh.cpp:10917
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7949
virtual void Save() override
RefCoord s[3]
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:1085
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3673
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:9467
double region(const Vector &p)
Abstract data type element.
Definition: element.hpp:28
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:273
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:10992
void recover_bdr_partitioning(const Mesh *mesh, const Array< int > &partitioning, Array< int > &bdr_partitioning)
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:6247