MFEM  v4.6.0
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  "r) random refinement with a probability\n"
412  "--> " << flush;
413  char sk;
414  cin >> sk;
415  switch (sk)
416  {
417  case 's':
418  mesh->UniformRefinement();
419  // Make sure tet-only meshes are marked for local refinement.
420  mesh->Finalize(true);
421  break;
422  case 'b':
423  mesh->UniformRefinement(1); // ref_algo = 1
424  break;
425  case 'u':
426  case 'g':
427  {
428  cout << "enter refinement factor --> " << flush;
429  int ref_factor;
430  cin >> ref_factor;
431  if (ref_factor <= 1 || ref_factor > 32) { break; }
432  int ref_type = (sk == 'u') ? BasisType::ClosedUniform :
434  *mesh = Mesh::MakeRefined(*mesh, ref_factor, ref_type);
435  break;
436  }
437  case 'l':
438  {
439  Vector pt;
440  Array<int> marked_elements;
441  for (int i = 0; i < mesh->GetNE(); i++)
442  {
443  // check all nodes of the element
445  mesh->GetElementTransformation(i, &T);
446  for (int j = 0; j < T.GetPointMat().Width(); j++)
447  {
448  T.GetPointMat().GetColumnReference(j, pt);
449  if (region(pt) <= region_eps)
450  {
451  marked_elements.Append(i);
452  break;
453  }
454  }
455  }
456  mesh->GeneralRefinement(marked_elements);
457  break;
458  }
459  case 'r':
460  {
461  bool nc_simplices = true;
462  mesh->EnsureNCMesh(nc_simplices);
463  cout << "enter probability --> " << flush;
464  double probability;
465  cin >> probability;
466  if (probability < 0.0 || probability > 1.0) { break; }
467  mesh->RandomRefinement(probability);
468  break;
469  }
470  }
471  print_char = 1;
472  }
473 
474  if (mk == 'c')
475  {
476  int p;
477  cout << "enter new order for mesh curvature --> " << flush;
478  cin >> p;
479  mesh->SetCurvature(p > 0 ? p : -p, p <= 0);
480  print_char = 1;
481  }
482 
483  if (mk == 's')
484  {
485  double factor;
486  cout << "scaling factor ---> " << flush;
487  cin >> factor;
488 
489  GridFunction *nodes = mesh->GetNodes();
490  if (nodes == NULL)
491  {
492  for (int i = 0; i < mesh->GetNV(); i++)
493  {
494  double *v = mesh->GetVertex(i);
495  v[0] *= factor;
496  v[1] *= factor;
497  if (dim == 3)
498  {
499  v[2] *= factor;
500  }
501  }
502  }
503  else
504  {
505  *nodes *= factor;
506  }
507 
508  print_char = 1;
509  }
510 
511  if (mk == 't')
512  {
513  char type;
514  cout << "Choose a transformation:\n"
515  "u) User-defined transform through mesh-explorer::transformation()\n"
516  "k) Kershaw transform\n"<< "---> " << flush;
517  cin >> type;
518  if (type == 'u')
519  {
520  mesh->Transform(transformation);
521  }
522  else if (type == 'k')
523  {
524  cout << "Note: For Kershaw transformation, the input must be "
525  "Cartesian aligned with nx multiple of 6 and "
526  "both ny and nz multiples of 2."
527  "Kershaw transform works for 2D meshes also.\n" << flush;
528 
529  double epsy, epsz = 0.0;
530  cout << "Kershaw transform factor, epsy in (0, 1]) ---> " << flush;
531  cin >> epsy;
532  if (mesh->Dimension() == 3)
533  {
534  cout << "Kershaw transform factor, epsz in (0, 1]) ---> " << flush;
535  cin >> epsz;
536  }
537  common::KershawTransformation kershawT(mesh->Dimension(), epsy, epsz);
538  mesh->Transform(kershawT);
539  }
540  else
541  {
542  MFEM_ABORT("Transformation type not supported.");
543  }
544  print_char = 1;
545  }
546 
547  if (mk == 'j')
548  {
549  double jitter;
550  cout << "jitter factor ---> " << flush;
551  cin >> jitter;
552 
553  GridFunction *nodes = mesh->GetNodes();
554 
555  if (nodes == NULL)
556  {
557  cerr << "The mesh should have nodes, introduce curvature first!\n";
558  }
559  else
560  {
561  FiniteElementSpace *fespace = nodes->FESpace();
562 
563  GridFunction rdm(fespace);
564  rdm.Randomize();
565  rdm -= 0.5; // shift to random values in [-0.5,0.5]
566  rdm *= jitter;
567 
568  // compute minimal local mesh size
569  Vector h0(fespace->GetNDofs());
570  h0 = infinity();
571  {
572  Array<int> dofs;
573  for (int i = 0; i < fespace->GetNE(); i++)
574  {
575  fespace->GetElementDofs(i, dofs);
576  for (int j = 0; j < dofs.Size(); j++)
577  {
578  h0(dofs[j]) = std::min(h0(dofs[j]), mesh->GetElementSize(i));
579  }
580  }
581  }
582 
583  // scale the random values to be of order of the local mesh size
584  for (int i = 0; i < fespace->GetNDofs(); i++)
585  {
586  for (int d = 0; d < dim; d++)
587  {
588  rdm(fespace->DofToVDof(i,d)) *= h0(i);
589  }
590  }
591 
592  char move_bdr = 'n';
593  cout << "move boundary nodes? [y/n] ---> " << flush;
594  cin >> move_bdr;
595 
596  // don't perturb the boundary
597  if (move_bdr == 'n')
598  {
599  Array<int> vdofs;
600  for (int i = 0; i < fespace->GetNBE(); i++)
601  {
602  fespace->GetBdrElementVDofs(i, vdofs);
603  for (int j = 0; j < vdofs.Size(); j++)
604  {
605  rdm(vdofs[j]) = 0.0;
606  }
607  }
608  }
609 
610  *nodes += rdm;
611  }
612 
613  print_char = 1;
614  }
615 
616  if (mk == 'x')
617  {
618  int sd, nz = 0;
619  DenseMatrix J(dim);
620  double min_det_J, max_det_J, min_det_J_z, max_det_J_z;
621  double min_kappa, max_kappa, max_ratio_det_J_z;
622  min_det_J = min_kappa = infinity();
623  max_det_J = max_kappa = max_ratio_det_J_z = -infinity();
624  cout << "subdivision factor ---> " << flush;
625  cin >> sd;
626  Array<int> bad_elems_by_geom(Geometry::NumGeom);
627  bad_elems_by_geom = 0;
628  // Only print so many to keep output compact
629  const int max_to_print = 10;
630  for (int i = 0; i < mesh->GetNE(); i++)
631  {
632  Geometry::Type geom = mesh->GetElementBaseGeometry(i);
634 
635  RefinedGeometry *RefG = GlobGeometryRefiner.Refine(geom, sd, 1);
636  IntegrationRule &ir = RefG->RefPts;
637 
638  min_det_J_z = infinity();
639  max_det_J_z = -infinity();
640  for (int j = 0; j < ir.GetNPoints(); j++)
641  {
642  T->SetIntPoint(&ir.IntPoint(j));
643  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
644 
645  double det_J = J.Det();
646  double kappa =
648 
649  min_det_J_z = fmin(min_det_J_z, det_J);
650  max_det_J_z = fmax(max_det_J_z, det_J);
651 
652  min_kappa = fmin(min_kappa, kappa);
653  max_kappa = fmax(max_kappa, kappa);
654  }
655  max_ratio_det_J_z =
656  fmax(max_ratio_det_J_z, max_det_J_z/min_det_J_z);
657  min_det_J = fmin(min_det_J, min_det_J_z);
658  max_det_J = fmax(max_det_J, max_det_J_z);
659  if (min_det_J_z <= 0.0)
660  {
661  if (nz < max_to_print)
662  {
663  Vector center;
664  mesh->GetElementCenter(i, center);
665  cout << "det(J) < 0 = " << min_det_J_z << " in element "
666  << i << ", centered at: ";
667  center.Print();
668  }
669  nz++;
670  bad_elems_by_geom[geom]++;
671  }
672  }
673  if (nz >= max_to_print)
674  {
675  cout << "det(J) < 0 for " << nz - max_to_print << " more elements "
676  << "not printed.\n";
677  }
678  cout << "\nbad elements = " << nz;
679  if (nz)
680  {
681  cout << " -- ";
682  Mesh::PrintElementsByGeometry(dim, bad_elems_by_geom, cout);
683  }
684  cout << "\nmin det(J) = " << min_det_J
685  << "\nmax det(J) = " << max_det_J
686  << "\nglobal ratio = " << max_det_J/min_det_J
687  << "\nmax el ratio = " << max_ratio_det_J_z
688  << "\nmin kappa = " << min_kappa
689  << "\nmax kappa = " << max_kappa << endl;
690  }
691 
692  if (mk == 'f')
693  {
694  DenseMatrix point_mat(sdim,1);
695  cout << "\npoint in physical space ---> " << flush;
696  for (int i = 0; i < sdim; i++)
697  {
698  cin >> point_mat(i,0);
699  }
700  Array<int> elem_ids;
702 
703  // physical -> reference space
704  mesh->FindPoints(point_mat, elem_ids, ips);
705 
706  cout << "point in reference space:";
707  if (elem_ids[0] == -1)
708  {
709  cout << " NOT FOUND!\n";
710  }
711  else
712  {
713  cout << " element " << elem_ids[0] << ", ip =";
714  cout << " " << ips[0].x;
715  if (sdim > 1)
716  {
717  cout << " " << ips[0].y;
718  if (sdim > 2)
719  {
720  cout << " " << ips[0].z;
721  }
722  }
723  cout << endl;
724  }
725  }
726 
727  if (mk == 'o')
728  {
729  cout << "What type of reordering?\n"
730  "g) Gecko edge-product minimization\n"
731  "h) Hilbert spatial sort\n"
732  "--> " << flush;
733  char rk;
734  cin >> rk;
735 
736  Array<int> ordering, tentative;
737  if (rk == 'h')
738  {
739  mesh->GetHilbertElementOrdering(ordering);
740  mesh->ReorderElements(ordering);
741  }
742  else if (rk == 'g')
743  {
744  int outer, inner, window, period;
745  cout << "Enter number of outer iterations (default 5): " << flush;
746  cin >> outer;
747  cout << "Enter number of inner iterations (default 4): " << flush;
748  cin >> inner;
749  cout << "Enter window size (default 4, beware of exponential cost): "
750  << flush;
751  cin >> window;
752  cout << "Enter period for window size increment (default 2): "
753  << flush;
754  cin >> period;
755 
756  double best_cost = infinity();
757  for (int i = 0; i < outer; i++)
758  {
759  int seed = i+1;
760  double cost = mesh->GetGeckoElementOrdering(
761  tentative, inner, window, period, seed, true);
762 
763  if (cost < best_cost)
764  {
765  ordering = tentative;
766  best_cost = cost;
767  }
768  }
769  cout << "Final cost: " << best_cost << endl;
770 
771  mesh->ReorderElements(ordering);
772  }
773  }
774 
775  // These are most of the cases that open a new GLVis window
776  if (mk == 'm' || mk == 'b' || mk == 'e' || mk == 'v' || mk == 'h' ||
777  mk == 'k' || mk == 'J' || mk == 'p' || mk == 'B' || mk == 'P')
778  {
779  FiniteElementSpace *bdr_attr_fespace = NULL;
780  FiniteElementSpace *attr_fespace =
781  new FiniteElementSpace(mesh, attr_fec);
782  GridFunction bdr_attr;
783  GridFunction attr(attr_fespace);
784 
785  if (mk == 'm')
786  {
787  for (int i = 0; i < mesh->GetNE(); i++)
788  {
789  attr(i) = mesh->GetAttribute(i);
790  }
791  }
792 
793  if (mk == 'P')
794  {
795  for (int i = 0; i < mesh->GetNE(); i++)
796  {
797  attr(i) = partitioning[i] + 1;
798  }
799  }
800 
801  if (mk == 'b' || mk == 'B')
802  {
803  if (dim == 3)
804  {
805  delete bdr_mesh;
806  bdr_mesh = skin_mesh(mesh);
807  bdr_attr_fespace =
808  new FiniteElementSpace(bdr_mesh, bdr_attr_fec);
809  bdr_attr.SetSpace(bdr_attr_fespace);
810  if (mk == 'b')
811  {
812  for (int i = 0; i < bdr_mesh->GetNE(); i++)
813  {
814  bdr_attr(i) = bdr_mesh->GetAttribute(i);
815  }
816  }
817  else if (mk == 'B')
818  {
819  for (int i = 0; i < bdr_mesh->GetNE(); i++)
820  {
821  bdr_attr(i) = bdr_partitioning[i] + 1;
822  }
823  }
824  else
825  {
826  MFEM_WARNING("Unimplemented case.");
827  }
828  }
829  else
830  {
831  MFEM_WARNING("Unsupported mesh dimension.");
832  attr = 1.0;
833  }
834  }
835 
836  if (mk == 'v')
837  {
838  attr = 1.0;
839  }
840 
841  if (mk == 'e')
842  {
843  Array<int> coloring;
844  srand(time(0));
845  double a = double(rand()) / (double(RAND_MAX) + 1.);
846  int el0 = (int)floor(a * mesh->GetNE());
847  cout << "Generating coloring starting with element " << el0+1
848  << " / " << mesh->GetNE() << endl;
849  mesh->GetElementColoring(coloring, el0);
850  for (int i = 0; i < coloring.Size(); i++)
851  {
852  attr(i) = coloring[i];
853  }
854  cout << "Number of colors: " << attr.Max() + 1 << endl;
855  for (int i = 0; i < mesh->GetNE(); i++)
856  {
857  attr(i) = i; // coloring by element number
858  }
859  }
860 
861  if (mk == 'h')
862  {
863  DenseMatrix J(dim);
864  double h_min, h_max;
865  h_min = infinity();
866  h_max = -h_min;
867  for (int i = 0; i < mesh->GetNE(); i++)
868  {
869  Geometry::Type geom = mesh->GetElementBaseGeometry(i);
871  T->SetIntPoint(&Geometries.GetCenter(geom));
872  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
873 
874  attr(i) = J.Det();
875  if (attr(i) < 0.0)
876  {
877  attr(i) = -pow(-attr(i), 1.0/double(dim));
878  }
879  else
880  {
881  attr(i) = pow(attr(i), 1.0/double(dim));
882  }
883  h_min = min(h_min, attr(i));
884  h_max = max(h_max, attr(i));
885  }
886  cout << "h_min = " << h_min << ", h_max = " << h_max << endl;
887  }
888 
889  if (mk == 'k')
890  {
891  DenseMatrix J(dim);
892  for (int i = 0; i < mesh->GetNE(); i++)
893  {
894  Geometry::Type geom = mesh->GetElementBaseGeometry(i);
896  T->SetIntPoint(&Geometries.GetCenter(geom));
897  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
898  attr(i) = J.CalcSingularvalue(0) / J.CalcSingularvalue(dim-1);
899  }
900  }
901 
902  if (mk == 'J')
903  {
904  // The "scaled Jacobian" is the determinant of the Jacobian scaled
905  // by the l2 norms of its columns. It can be used to identify badly
906  // skewed elements, since it takes values between 0 and 1, with 0
907  // corresponding to a flat element, and 1 to orthogonal columns.
908  DenseMatrix J(dim);
909  int sd;
910  cout << "subdivision factor ---> " << flush;
911  cin >> sd;
912  for (int i = 0; i < mesh->GetNE(); i++)
913  {
914  Geometry::Type geom = mesh->GetElementBaseGeometry(i);
916 
917  RefinedGeometry *RefG = GlobGeometryRefiner.Refine(geom, sd, 1);
918  IntegrationRule &ir = RefG->RefPts;
919 
920  // For each element, find the minimal scaled Jacobian in a
921  // lattice of points with the given subdivision factor.
922  attr(i) = infinity();
923  for (int j = 0; j < ir.GetNPoints(); j++)
924  {
925  T->SetIntPoint(&ir.IntPoint(j));
926  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
927 
928  // Jacobian determinant
929  double sJ = J.Det();
930 
931  for (int k = 0; k < J.Width(); k++)
932  {
933  Vector col;
934  J.GetColumnReference(k,col);
935  // Scale by column norms
936  sJ /= col.Norml2();
937  }
938 
939  attr(i) = fmin(sJ, attr(i));
940  }
941  }
942  }
943 
944  if (mk == 'p')
945  {
946  cout << "What type of partitioning?\n"
947  "c) Cartesian\n"
948  "s) Simple 1D split of the element sequence\n"
949  "0) METIS_PartGraphRecursive (sorted neighbor lists)\n"
950  "1) METIS_PartGraphKway (sorted neighbor lists)"
951  " (default)\n"
952  "2) METIS_PartGraphVKway (sorted neighbor lists)\n"
953  "3) METIS_PartGraphRecursive\n"
954  "4) METIS_PartGraphKway\n"
955  "5) METIS_PartGraphVKway\n"
956  "--> " << flush;
957  char pk;
958  cin >> pk;
959  if (pk == 'c')
960  {
961  int nxyz[3];
962  cout << "Enter nx: " << flush;
963  cin >> nxyz[0]; np = nxyz[0];
964  if (mesh->Dimension() > 1)
965  {
966  cout << "Enter ny: " << flush;
967  cin >> nxyz[1]; np *= nxyz[1];
968  if (mesh->Dimension() > 2)
969  {
970  cout << "Enter nz: " << flush;
971  cin >> nxyz[2]; np *= nxyz[2];
972  }
973  }
974  int *part = mesh->CartesianPartitioning(nxyz);
975  partitioning = Array<int>(part, mesh->GetNE());
976  delete [] part;
977  recover_bdr_partitioning(mesh, partitioning, bdr_partitioning);
978  }
979  else if (pk == 's')
980  {
981  cout << "Enter number of processors: " << flush;
982  cin >> np;
983 
984  partitioning.SetSize(mesh->GetNE());
985  for (int i = 0; i < mesh->GetNE(); i++)
986  {
987  partitioning[i] = i * np / mesh->GetNE();
988  }
989  recover_bdr_partitioning(mesh, partitioning, bdr_partitioning);
990  }
991  else
992  {
993  int part_method = pk - '0';
994  if (part_method < 0 || part_method > 5)
995  {
996  continue;
997  }
998  cout << "Enter number of processors: " << flush;
999  cin >> np;
1000  int *part = mesh->GeneratePartitioning(np, part_method);
1001  partitioning = Array<int>(part, mesh->GetNE());
1002  delete [] part;
1003  recover_bdr_partitioning(mesh, partitioning, bdr_partitioning);
1004  }
1005  if (partitioning)
1006  {
1007  const char part_file[] = "partitioning.txt";
1008  ofstream opart(part_file);
1009  opart << "number_of_elements " << mesh->GetNE() << '\n'
1010  << "number_of_processors " << np << '\n';
1011  for (int i = 0; i < mesh->GetNE(); i++)
1012  {
1013  opart << partitioning[i] << '\n';
1014  }
1015  cout << "Partitioning file: " << part_file << endl;
1016 
1017  Array<int> proc_el(np);
1018  proc_el = 0;
1019  for (int i = 0; i < mesh->GetNE(); i++)
1020  {
1021  proc_el[partitioning[i]]++;
1022  }
1023  int min_el = proc_el[0], max_el = proc_el[0];
1024  for (int i = 1; i < np; i++)
1025  {
1026  if (min_el > proc_el[i])
1027  {
1028  min_el = proc_el[i];
1029  }
1030  if (max_el < proc_el[i])
1031  {
1032  max_el = proc_el[i];
1033  }
1034  }
1035  cout << "Partitioning stats:\n"
1036  << " "
1037  << setw(12) << "minimum"
1038  << setw(12) << "average"
1039  << setw(12) << "maximum"
1040  << setw(12) << "total" << '\n';
1041  cout << " elements "
1042  << setw(12) << min_el
1043  << setw(12) << double(mesh->GetNE())/np
1044  << setw(12) << max_el
1045  << setw(12) << mesh->GetNE() << endl;
1046  }
1047  else
1048  {
1049  continue;
1050  }
1051 
1052  for (int i = 0; i < mesh->GetNE(); i++)
1053  {
1054  attr(i) = partitioning[i] + 1;
1055  }
1056  }
1057 
1058  char vishost[] = "localhost";
1059  int visport = 19916;
1060  socketstream sol_sock(vishost, visport);
1061  if (sol_sock.is_open())
1062  {
1063  sol_sock.precision(14);
1064  if (sdim == 2)
1065  {
1066  sol_sock << "fem2d_gf_data_keys\n";
1067  if (mk != 'p')
1068  {
1069  mesh->Print(sol_sock);
1070  }
1071  else
1072  {
1073  // NURBS meshes do not support PrintWithPartitioning
1074  if (mesh->NURBSext)
1075  {
1076  mesh->Print(sol_sock);
1077  for (int i = 0; i < mesh->GetNE(); i++)
1078  {
1079  attr(i) = partitioning[i];
1080  }
1081  }
1082  else
1083  {
1084  mesh->PrintWithPartitioning(partitioning, sol_sock, 1);
1085  }
1086  }
1087  attr.Save(sol_sock);
1088  sol_sock << "RjlmAb***********";
1089  if (mk == 'v')
1090  {
1091  sol_sock << "e";
1092  }
1093  else
1094  {
1095  sol_sock << "\n";
1096  }
1097  }
1098  else
1099  {
1100  sol_sock << "fem3d_gf_data_keys\n";
1101  if (mk == 'v' || mk == 'h' || mk == 'k' || mk == 'J')
1102  {
1103  mesh->Print(sol_sock);
1104  }
1105  else if (mk == 'b' || mk == 'B')
1106  {
1107  bdr_mesh->Print(sol_sock);
1108  bdr_attr.Save(sol_sock);
1109  sol_sock << "mcaaA";
1110  // Switch to a discrete color scale
1111  sol_sock << "pppppp" << "pppppp" << "pppppp";
1112  }
1113  else
1114  {
1115  // NURBS meshes do not support PrintWithPartitioning
1116  if (mesh->NURBSext)
1117  {
1118  mesh->Print(sol_sock);
1119  for (int i = 0; i < mesh->GetNE(); i++)
1120  {
1121  attr(i) = partitioning[i];
1122  }
1123  }
1124  else
1125  {
1126  mesh->PrintWithPartitioning(partitioning, sol_sock);
1127  }
1128  }
1129  if (mk != 'b' && mk != 'B')
1130  {
1131  attr.Save(sol_sock);
1132  sol_sock << "maaA";
1133  if (mk == 'v')
1134  {
1135  sol_sock << "aa";
1136  }
1137  else
1138  {
1139  sol_sock << "\n";
1140  }
1141  }
1142  }
1143  sol_sock << flush;
1144  }
1145  else
1146  {
1147  cout << "Unable to connect to "
1148  << vishost << ':' << visport << endl;
1149  }
1150  delete attr_fespace;
1151  delete bdr_attr_fespace;
1152  }
1153 
1154  if (mk == 'l')
1155  {
1156  // Project and plot the function 'f'
1157  int p;
1158  FiniteElementCollection *fec = NULL;
1159  cout << "Enter projection space order: " << flush;
1160  cin >> p;
1161  if (p >= 1)
1162  {
1163  fec = new H1_FECollection(p, mesh->Dimension(),
1165  }
1166  else
1167  {
1168  fec = new DG_FECollection(-p, mesh->Dimension(),
1170  }
1171  FiniteElementSpace fes(mesh, fec);
1172  GridFunction level(&fes);
1173  FunctionCoefficient coeff(f);
1174  level.ProjectCoefficient(coeff);
1175  char vishost[] = "localhost";
1176  int visport = 19916;
1177  socketstream sol_sock(vishost, visport);
1178  if (sol_sock.is_open())
1179  {
1180  sol_sock.precision(14);
1181  sol_sock << "solution\n" << *mesh << level << flush;
1182  }
1183  else
1184  {
1185  cout << "Unable to connect to "
1186  << vishost << ':' << visport << endl;
1187  }
1188  delete fec;
1189  }
1190 
1191  if (mk == 'S')
1192  {
1193  const char omesh_file[] = "mesh-explorer.mesh";
1194  ofstream omesh(omesh_file);
1195  omesh.precision(14);
1196  mesh->Print(omesh);
1197  cout << "New mesh file: " << omesh_file << endl;
1198  }
1199 
1200  if (mk == 'V')
1201  {
1202  const char omesh_file[] = "mesh-explorer.vtk";
1203  ofstream omesh(omesh_file);
1204  omesh.precision(14);
1205  mesh->PrintVTK(omesh);
1206  cout << "New VTK mesh file: " << omesh_file << endl;
1207  }
1208 
1209  if (mk == 'D')
1210  {
1211  cout << "What type of DataCollection?\n"
1212  "p) ParaView Data Collection\n"
1213  "v) VisIt Data Collection\n"
1214  "--> " << flush;
1215  char dk;
1216  cin >> dk;
1217  if (dk == 'p' || dk == 'P')
1218  {
1219  const char omesh_file[] = "mesh-explorer-paraview";
1220  ParaViewDataCollection dc(omesh_file, mesh);
1221  if (mesh->GetNodes())
1222  {
1223  int order = mesh->GetNodes()->FESpace()->GetMaxElementOrder();
1224  if (order > 1)
1225  {
1226  dc.SetHighOrderOutput(true);
1227  dc.SetLevelsOfDetail(order);
1228  }
1229  }
1230  dc.Save();
1231  cout << "New ParaView mesh file: " << omesh_file << endl;
1232  }
1233  else if (dk == 'v' || dk == 'V')
1234  {
1235  const char omesh_file[] = "mesh-explorer-visit";
1236  VisItDataCollection dc(omesh_file, mesh);
1237  dc.SetPrecision(14);
1238  dc.Save();
1239  cout << "New VisIt mesh file: " << omesh_file << "_000000.mfem_root"
1240  << endl;
1241  }
1242  else
1243  {
1244  cout << "Unrecognized DataCollection type: \"" << dk << "\""
1245  << endl;
1246  }
1247  }
1248 
1249 #ifdef MFEM_USE_ZLIB
1250  if (mk == 'Z')
1251  {
1252  const char omesh_file[] = "mesh-explorer.mesh.gz";
1253  ofgzstream omesh(omesh_file, "zwb9");
1254  omesh.precision(14);
1255  mesh->Print(omesh);
1256  cout << "New mesh file: " << omesh_file << endl;
1257  }
1258 #endif
1259 
1260  }
1261 
1262  delete bdr_attr_fec;
1263  delete attr_fec;
1264  delete bdr_mesh;
1265  delete mesh;
1266  return 0;
1267 }
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:2050
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:605
int * CartesianPartitioning(int nxyz[])
Definition: mesh.cpp:7383
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:253
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:96
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:5630
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:1694
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:331
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:1242
static void PrintElementsByGeometry(int dim, const Array< int > &num_elems_by_geom, std::ostream &out)
Auxiliary method used by PrintCharacteristics().
Definition: mesh.cpp:235
static const int NumGeom
Definition: geom.hpp:42
int Dimension() const
Dimension of the reference space used within the elements.
Definition: mesh.hpp:1020
int * GeneratePartitioning(int nparts, int part_method=1)
Definition: mesh.cpp:7427
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:8329
void Print(std::ostream &out=mfem::out, int width=8) const
Prints vector to stream out.
Definition: vector.cpp:756
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:3793
int AddTriangle(int v1, int v2, int v3, int attr=1)
Definition: mesh.cpp:1680
double Det() const
Definition: densemat.cpp:487
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
void Transform(void(*f)(const Vector &, Vector &))
Definition: mesh.cpp:12045
int GetNDofs() const
Returns number of degrees of freedom. This is the number of Local Degrees of Freedom.
Definition: fespace.hpp:706
double CalcSingularvalue(const int i) const
Return the i-th singular value (decreasing order) of NxN matrix, N=1,2,3.
Definition: densemat.cpp:1267
void PrintError(std::ostream &out) const
Print the error message.
Definition: optparser.cpp:362
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:249
double kappa
Definition: ex24.cpp:54
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:159
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:817
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:1089
Geometry Geometries
Definition: fe.cpp:49
void DeleteAll()
Delete the whole array.
Definition: array.hpp:854
Mesh * skin_mesh(Mesh *mesh)
bool Help() const
Return true if we are flagged to print the help message.
Definition: optparser.hpp:162
std::function< double(const Vector &)> f(double mass_coeff)
Definition: lor_mms.hpp:30
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:1083
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:1190
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:256
int AddVertex(double x, double y=0.0, double z=0.0)
Definition: mesh.cpp:1626
GeometryRefiner GlobGeometryRefiner
Definition: geom.cpp:1773
double GetElementSize(ElementTransformation *T, int type=0)
Definition: mesh.cpp:83
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:10719
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:759
void EnsureNCMesh(bool simplices_nonconforming=false)
Definition: mesh.cpp:9888
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:723
void transformation(const Vector &p, Vector &v)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:10232
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:5635
const double * GetVertex(int i) const
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:1125
virtual const char * Name() const
Definition: fe_coll.hpp:80
void GetHilbertElementOrdering(Array< int > &ordering)
Definition: mesh.cpp:2217
Piecewise-constant discontinuous finite elements in 3D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:1119
void RandomRefinement(double prob, bool aniso=false, int nonconforming=-1, int nc_limit=0)
Refine each element with given probability. Uses GeneralRefinement.
Definition: mesh.cpp:9912
void GetElementCenter(int i, Vector &center)
Definition: mesh.cpp:75
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:941
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:691
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:736
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
Definition: fe_coll.hpp:377
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:1666
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:219
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:2945
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:687
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;. The returned indices are offsets into an ldo...
Definition: fespace.cpp:2731
double Max() const
Returns the maximal element of the vector.
Definition: vector.cpp:925
int main(int argc, char *argv[])
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition: mesh.hpp:1023
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:1086
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:355
DofTransformation * GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for i&#39;th boundary element. The returned indices are offsets int...
Definition: fespace.cpp:297
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition: mesh.cpp:3042
void ReorderElements(const Array< int > &ordering, bool reorder_vertices=true)
Definition: mesh.cpp:2269
int dim
Definition: ex24.cpp:53
void PrintHelp(std::ostream &out) const
Print the help message.
Definition: optparser.cpp:405
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:745
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
Definition: vector.hpp:44
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:2415
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:835
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition: gridfunc.cpp:195
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:12316
int DofToVDof(int dof, int vd, int ndofs=-1) const
Compute a single vdof corresponding to the index dof and the vector index vd.
Definition: fespace.cpp:251
Nodes: x_i = i/(n-1), i=0,...,n-1.
Definition: fe_base.hpp:35
Vector data type.
Definition: vector.hpp:58
virtual void Print(std::ostream &os=mfem::out) const
Definition: mesh.hpp:2011
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:259
void GetElementColoring(Array< int > &colors, int el0=0)
Definition: mesh.cpp:11271
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:8302
virtual void Save() override
RefCoord s[3]
const Element * GetBdrElement(int i) const
Return pointer to the i&#39;th boundary element object.
Definition: mesh.hpp:1158
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3696
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:9820
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:11346
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:6600