MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
mesh-explorer.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 //
12 // -----------------------------------------------------
13 // Mesh Explorer Miniapp: Explore and manipulate meshes
14 // -----------------------------------------------------
15 //
16 // This miniapp is a handy tool to examine, visualize and manipulate a given
17 // mesh. Some of its features are:
18 //
19 // - visualizing of mesh materials and individual mesh elements
20 // - mesh scaling, randomization, and general transformation
21 // - manipulation of the mesh curvature
22 // - the ability to simulate parallel partitioning
23 // - quantitative and visual reports of mesh quality
24 //
25 // Compile with: make mesh-explorer
26 //
27 // Sample runs: mesh-explorer
28 // mesh-explorer -m ../../data/beam-tri.mesh
29 // mesh-explorer -m ../../data/star-q2.mesh
30 // mesh-explorer -m ../../data/disc-nurbs.mesh
31 // mesh-explorer -m ../../data/escher-p3.mesh
32 // mesh-explorer -m ../../data/mobius-strip.mesh
33 
34 #include "mfem.hpp"
35 #include "../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  "q) Quit\n"
388 #ifdef MFEM_USE_ZLIB
389  "Z) Save in MFEM format with compression\n"
390 #endif
391  "--> " << flush;
392  char mk;
393  cin >> mk;
394  if (!cin) { break; }
395 
396  if (mk == 'q')
397  {
398  break;
399  }
400 
401  if (mk == 'r')
402  {
403  cout <<
404  "Choose type of refinement:\n"
405  "s) standard refinement with Mesh::UniformRefinement()\n"
406  "b) Mesh::UniformRefinement() (bisection for tet meshes)\n"
407  "u) uniform refinement with a factor\n"
408  "g) non-uniform refinement (Gauss-Lobatto) with a factor\n"
409  "l) refine locally using the region() function\n"
410  "--> " << flush;
411  char sk;
412  cin >> sk;
413  switch (sk)
414  {
415  case 's':
416  mesh->UniformRefinement();
417  // Make sure tet-only meshes are marked for local refinement.
418  mesh->Finalize(true);
419  break;
420  case 'b':
421  mesh->UniformRefinement(1); // ref_algo = 1
422  break;
423  case 'u':
424  case 'g':
425  {
426  cout << "enter refinement factor --> " << flush;
427  int ref_factor;
428  cin >> ref_factor;
429  if (ref_factor <= 1 || ref_factor > 32) { break; }
430  int ref_type = (sk == 'u') ? BasisType::ClosedUniform :
432  *mesh = Mesh::MakeRefined(*mesh, ref_factor, ref_type);
433  break;
434  }
435  case 'l':
436  {
437  Vector pt;
438  Array<int> marked_elements;
439  for (int i = 0; i < mesh->GetNE(); i++)
440  {
441  // check all nodes of the element
443  mesh->GetElementTransformation(i, &T);
444  for (int j = 0; j < T.GetPointMat().Width(); j++)
445  {
446  T.GetPointMat().GetColumnReference(j, pt);
447  if (region(pt) <= region_eps)
448  {
449  marked_elements.Append(i);
450  break;
451  }
452  }
453  }
454  mesh->GeneralRefinement(marked_elements);
455  break;
456  }
457  }
458  print_char = 1;
459  }
460 
461  if (mk == 'c')
462  {
463  int p;
464  cout << "enter new order for mesh curvature --> " << flush;
465  cin >> p;
466  mesh->SetCurvature(p > 0 ? p : -p, p <= 0);
467  print_char = 1;
468  }
469 
470  if (mk == 's')
471  {
472  double factor;
473  cout << "scaling factor ---> " << flush;
474  cin >> factor;
475 
476  GridFunction *nodes = mesh->GetNodes();
477  if (nodes == NULL)
478  {
479  for (int i = 0; i < mesh->GetNV(); i++)
480  {
481  double *v = mesh->GetVertex(i);
482  v[0] *= factor;
483  v[1] *= factor;
484  if (dim == 3)
485  {
486  v[2] *= factor;
487  }
488  }
489  }
490  else
491  {
492  *nodes *= factor;
493  }
494 
495  print_char = 1;
496  }
497 
498  if (mk == 't')
499  {
500  char type;
501  cout << "Choose a transformation:\n"
502  "u) User-defined transform through mesh-explorer::transformation()\n"
503  "k) Kershaw transform\n"<< "---> " << flush;
504  cin >> type;
505  if (type == 'u')
506  {
507  mesh->Transform(transformation);
508  }
509  else if (type == 'k')
510  {
511  cout << "Note: For Kershaw transformation, the input must be "
512  "Cartesian aligned with nx multiple of 6 and "
513  "both ny and nz multiples of 2."
514  "Kershaw transform works for 2D meshes also.\n" << flush;
515 
516  double epsy, epsz = 0.0;
517  cout << "Kershaw transform factor, epsy in (0, 1]) ---> " << flush;
518  cin >> epsy;
519  if (mesh->Dimension() == 3)
520  {
521  cout << "Kershaw transform factor, epsz in (0, 1]) ---> " << flush;
522  cin >> epsz;
523  }
524  common::KershawTransformation kershawT(mesh->Dimension(), epsy, epsz);
525  mesh->Transform(kershawT);
526  }
527  else
528  {
529  MFEM_ABORT("Transformation type not supported.");
530  }
531  print_char = 1;
532  }
533 
534  if (mk == 'j')
535  {
536  double jitter;
537  cout << "jitter factor ---> " << flush;
538  cin >> jitter;
539 
540  GridFunction *nodes = mesh->GetNodes();
541 
542  if (nodes == NULL)
543  {
544  cerr << "The mesh should have nodes, introduce curvature first!\n";
545  }
546  else
547  {
548  FiniteElementSpace *fespace = nodes->FESpace();
549 
550  GridFunction rdm(fespace);
551  rdm.Randomize();
552  rdm -= 0.5; // shift to random values in [-0.5,0.5]
553  rdm *= jitter;
554 
555  // compute minimal local mesh size
556  Vector h0(fespace->GetNDofs());
557  h0 = infinity();
558  {
559  Array<int> dofs;
560  for (int i = 0; i < fespace->GetNE(); i++)
561  {
562  fespace->GetElementDofs(i, dofs);
563  for (int j = 0; j < dofs.Size(); j++)
564  {
565  h0(dofs[j]) = std::min(h0(dofs[j]), mesh->GetElementSize(i));
566  }
567  }
568  }
569 
570  // scale the random values to be of order of the local mesh size
571  for (int i = 0; i < fespace->GetNDofs(); i++)
572  {
573  for (int d = 0; d < dim; d++)
574  {
575  rdm(fespace->DofToVDof(i,d)) *= h0(i);
576  }
577  }
578 
579  char move_bdr = 'n';
580  cout << "move boundary nodes? [y/n] ---> " << flush;
581  cin >> move_bdr;
582 
583  // don't perturb the boundary
584  if (move_bdr == 'n')
585  {
586  Array<int> vdofs;
587  for (int i = 0; i < fespace->GetNBE(); i++)
588  {
589  fespace->GetBdrElementVDofs(i, vdofs);
590  for (int j = 0; j < vdofs.Size(); j++)
591  {
592  rdm(vdofs[j]) = 0.0;
593  }
594  }
595  }
596 
597  *nodes += rdm;
598  }
599 
600  print_char = 1;
601  }
602 
603  if (mk == 'x')
604  {
605  int sd, nz = 0;
606  DenseMatrix J(dim);
607  double min_det_J, max_det_J, min_det_J_z, max_det_J_z;
608  double min_kappa, max_kappa, max_ratio_det_J_z;
609  min_det_J = min_kappa = infinity();
610  max_det_J = max_kappa = max_ratio_det_J_z = -infinity();
611  cout << "subdivision factor ---> " << flush;
612  cin >> sd;
613  Array<int> bad_elems_by_geom(Geometry::NumGeom);
614  bad_elems_by_geom = 0;
615  // Only print so many to keep output compact
616  const int max_to_print = 10;
617  for (int i = 0; i < mesh->GetNE(); i++)
618  {
619  Geometry::Type geom = mesh->GetElementBaseGeometry(i);
621 
622  RefinedGeometry *RefG = GlobGeometryRefiner.Refine(geom, sd, 1);
623  IntegrationRule &ir = RefG->RefPts;
624 
625  min_det_J_z = infinity();
626  max_det_J_z = -infinity();
627  for (int j = 0; j < ir.GetNPoints(); j++)
628  {
629  T->SetIntPoint(&ir.IntPoint(j));
630  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
631 
632  double det_J = J.Det();
633  double kappa =
634  J.CalcSingularvalue(0) / J.CalcSingularvalue(dim-1);
635 
636  min_det_J_z = fmin(min_det_J_z, det_J);
637  max_det_J_z = fmax(max_det_J_z, det_J);
638 
639  min_kappa = fmin(min_kappa, kappa);
640  max_kappa = fmax(max_kappa, kappa);
641  }
642  max_ratio_det_J_z =
643  fmax(max_ratio_det_J_z, max_det_J_z/min_det_J_z);
644  min_det_J = fmin(min_det_J, min_det_J_z);
645  max_det_J = fmax(max_det_J, max_det_J_z);
646  if (min_det_J_z <= 0.0)
647  {
648  if (nz < max_to_print)
649  {
650  Vector center;
651  mesh->GetElementCenter(i, center);
652  cout << "det(J) < 0 = " << min_det_J_z << " in element "
653  << i << ", centered at: ";
654  center.Print();
655  }
656  nz++;
657  bad_elems_by_geom[geom]++;
658  }
659  }
660  if (nz >= max_to_print)
661  {
662  cout << "det(J) < 0 for " << nz - max_to_print << " more elements "
663  << "not printed.\n";
664  }
665  cout << "\nbad elements = " << nz;
666  if (nz)
667  {
668  cout << " -- ";
669  Mesh::PrintElementsByGeometry(dim, bad_elems_by_geom, cout);
670  }
671  cout << "\nmin det(J) = " << min_det_J
672  << "\nmax det(J) = " << max_det_J
673  << "\nglobal ratio = " << max_det_J/min_det_J
674  << "\nmax el ratio = " << max_ratio_det_J_z
675  << "\nmin kappa = " << min_kappa
676  << "\nmax kappa = " << max_kappa << endl;
677  }
678 
679  if (mk == 'f')
680  {
681  DenseMatrix point_mat(sdim,1);
682  cout << "\npoint in physical space ---> " << flush;
683  for (int i = 0; i < sdim; i++)
684  {
685  cin >> point_mat(i,0);
686  }
687  Array<int> elem_ids;
689 
690  // physical -> reference space
691  mesh->FindPoints(point_mat, elem_ids, ips);
692 
693  cout << "point in reference space:";
694  if (elem_ids[0] == -1)
695  {
696  cout << " NOT FOUND!\n";
697  }
698  else
699  {
700  cout << " element " << elem_ids[0] << ", ip =";
701  cout << " " << ips[0].x;
702  if (sdim > 1)
703  {
704  cout << " " << ips[0].y;
705  if (sdim > 2)
706  {
707  cout << " " << ips[0].z;
708  }
709  }
710  cout << endl;
711  }
712  }
713 
714  if (mk == 'o')
715  {
716  cout << "What type of reordering?\n"
717  "g) Gecko edge-product minimization\n"
718  "h) Hilbert spatial sort\n"
719  "--> " << flush;
720  char rk;
721  cin >> rk;
722 
723  Array<int> ordering, tentative;
724  if (rk == 'h')
725  {
726  mesh->GetHilbertElementOrdering(ordering);
727  mesh->ReorderElements(ordering);
728  }
729  else if (rk == 'g')
730  {
731  int outer, inner, window, period;
732  cout << "Enter number of outer iterations (default 5): " << flush;
733  cin >> outer;
734  cout << "Enter number of inner iterations (default 4): " << flush;
735  cin >> inner;
736  cout << "Enter window size (default 4, beware of exponential cost): "
737  << flush;
738  cin >> window;
739  cout << "Enter period for window size increment (default 2): "
740  << flush;
741  cin >> period;
742 
743  double best_cost = infinity();
744  for (int i = 0; i < outer; i++)
745  {
746  int seed = i+1;
747  double cost = mesh->GetGeckoElementOrdering(
748  tentative, inner, window, period, seed, true);
749 
750  if (cost < best_cost)
751  {
752  ordering = tentative;
753  best_cost = cost;
754  }
755  }
756  cout << "Final cost: " << best_cost << endl;
757 
758  mesh->ReorderElements(ordering);
759  }
760  }
761 
762  // These are most of the cases that open a new GLVis window
763  if (mk == 'm' || mk == 'b' || mk == 'e' || mk == 'v' || mk == 'h' ||
764  mk == 'k' || mk == 'J' || mk == 'p' || mk == 'B' || mk == 'P')
765  {
766  FiniteElementSpace *bdr_attr_fespace = NULL;
767  FiniteElementSpace *attr_fespace =
768  new FiniteElementSpace(mesh, attr_fec);
769  GridFunction bdr_attr;
770  GridFunction attr(attr_fespace);
771 
772  if (mk == 'm')
773  {
774  for (int i = 0; i < mesh->GetNE(); i++)
775  {
776  attr(i) = mesh->GetAttribute(i);
777  }
778  }
779 
780  if (mk == 'P')
781  {
782  for (int i = 0; i < mesh->GetNE(); i++)
783  {
784  attr(i) = partitioning[i] + 1;
785  }
786  }
787 
788  if (mk == 'b' || mk == 'B')
789  {
790  if (dim == 3)
791  {
792  delete bdr_mesh;
793  bdr_mesh = skin_mesh(mesh);
794  bdr_attr_fespace =
795  new FiniteElementSpace(bdr_mesh, bdr_attr_fec);
796  bdr_attr.SetSpace(bdr_attr_fespace);
797  if (mk == 'b')
798  {
799  for (int i = 0; i < bdr_mesh->GetNE(); i++)
800  {
801  bdr_attr(i) = bdr_mesh->GetAttribute(i);
802  }
803  }
804  else if (mk == 'B')
805  {
806  for (int i = 0; i < bdr_mesh->GetNE(); i++)
807  {
808  bdr_attr(i) = bdr_partitioning[i] + 1;
809  }
810  }
811  else
812  {
813  MFEM_WARNING("Unimplemented case.");
814  }
815  }
816  else
817  {
818  MFEM_WARNING("Unsupported mesh dimension.");
819  attr = 1.0;
820  }
821  }
822 
823  if (mk == 'v')
824  {
825  attr = 1.0;
826  }
827 
828  if (mk == 'e')
829  {
830  Array<int> coloring;
831  srand(time(0));
832  double a = double(rand()) / (double(RAND_MAX) + 1.);
833  int el0 = (int)floor(a * mesh->GetNE());
834  cout << "Generating coloring starting with element " << el0+1
835  << " / " << mesh->GetNE() << endl;
836  mesh->GetElementColoring(coloring, el0);
837  for (int i = 0; i < coloring.Size(); i++)
838  {
839  attr(i) = coloring[i];
840  }
841  cout << "Number of colors: " << attr.Max() + 1 << endl;
842  for (int i = 0; i < mesh->GetNE(); i++)
843  {
844  attr(i) = i; // coloring by element number
845  }
846  }
847 
848  if (mk == 'h')
849  {
850  DenseMatrix J(dim);
851  double h_min, h_max;
852  h_min = infinity();
853  h_max = -h_min;
854  for (int i = 0; i < mesh->GetNE(); i++)
855  {
856  Geometry::Type geom = mesh->GetElementBaseGeometry(i);
858  T->SetIntPoint(&Geometries.GetCenter(geom));
859  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
860 
861  attr(i) = J.Det();
862  if (attr(i) < 0.0)
863  {
864  attr(i) = -pow(-attr(i), 1.0/double(dim));
865  }
866  else
867  {
868  attr(i) = pow(attr(i), 1.0/double(dim));
869  }
870  h_min = min(h_min, attr(i));
871  h_max = max(h_max, attr(i));
872  }
873  cout << "h_min = " << h_min << ", h_max = " << h_max << endl;
874  }
875 
876  if (mk == 'k')
877  {
878  DenseMatrix J(dim);
879  for (int i = 0; i < mesh->GetNE(); i++)
880  {
881  Geometry::Type geom = mesh->GetElementBaseGeometry(i);
883  T->SetIntPoint(&Geometries.GetCenter(geom));
884  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
885  attr(i) = J.CalcSingularvalue(0) / J.CalcSingularvalue(dim-1);
886  }
887  }
888 
889  if (mk == 'J')
890  {
891  // The "scaled Jacobian" is the determinant of the Jacobian scaled
892  // by the l2 norms of its columns. It can be used to identify badly
893  // skewed elements, since it takes values between 0 and 1, with 0
894  // corresponding to a flat element, and 1 to orthogonal columns.
895  DenseMatrix J(dim);
896  int sd;
897  cout << "subdivision factor ---> " << flush;
898  cin >> sd;
899  for (int i = 0; i < mesh->GetNE(); i++)
900  {
901  Geometry::Type geom = mesh->GetElementBaseGeometry(i);
903 
904  RefinedGeometry *RefG = GlobGeometryRefiner.Refine(geom, sd, 1);
905  IntegrationRule &ir = RefG->RefPts;
906 
907  // For each element, find the minimal scaled Jacobian in a
908  // lattice of points with the given subdivision factor.
909  attr(i) = infinity();
910  for (int j = 0; j < ir.GetNPoints(); j++)
911  {
912  T->SetIntPoint(&ir.IntPoint(j));
913  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
914 
915  // Jacobian determinant
916  double sJ = J.Det();
917 
918  for (int k = 0; k < J.Width(); k++)
919  {
920  Vector col;
921  J.GetColumnReference(k,col);
922  // Scale by column norms
923  sJ /= col.Norml2();
924  }
925 
926  attr(i) = fmin(sJ, attr(i));
927  }
928  }
929  }
930 
931  if (mk == 'p')
932  {
933  cout << "What type of partitioning?\n"
934  "c) Cartesian\n"
935  "s) Simple 1D split of the element sequence\n"
936  "0) METIS_PartGraphRecursive (sorted neighbor lists)\n"
937  "1) METIS_PartGraphKway (sorted neighbor lists)"
938  " (default)\n"
939  "2) METIS_PartGraphVKway (sorted neighbor lists)\n"
940  "3) METIS_PartGraphRecursive\n"
941  "4) METIS_PartGraphKway\n"
942  "5) METIS_PartGraphVKway\n"
943  "--> " << flush;
944  char pk;
945  cin >> pk;
946  if (pk == 'c')
947  {
948  int nxyz[3];
949  cout << "Enter nx: " << flush;
950  cin >> nxyz[0]; np = nxyz[0];
951  if (mesh->Dimension() > 1)
952  {
953  cout << "Enter ny: " << flush;
954  cin >> nxyz[1]; np *= nxyz[1];
955  if (mesh->Dimension() > 2)
956  {
957  cout << "Enter nz: " << flush;
958  cin >> nxyz[2]; np *= nxyz[2];
959  }
960  }
961  partitioning = Array<int>(mesh->CartesianPartitioning(nxyz), mesh->GetNE());
962  recover_bdr_partitioning(mesh, partitioning, bdr_partitioning);
963  }
964  else if (pk == 's')
965  {
966  cout << "Enter number of processors: " << flush;
967  cin >> np;
968 
969  partitioning.SetSize(mesh->GetNE());
970  for (int i = 0; i < mesh->GetNE(); i++)
971  {
972  partitioning[i] = i * np / mesh->GetNE();
973  }
974  recover_bdr_partitioning(mesh, partitioning, bdr_partitioning);
975  }
976  else
977  {
978  int part_method = pk - '0';
979  if (part_method < 0 || part_method > 5)
980  {
981  continue;
982  }
983  cout << "Enter number of processors: " << flush;
984  cin >> np;
985  partitioning = Array<int>(mesh->GeneratePartitioning(np, part_method),
986  mesh->GetNE());
987  recover_bdr_partitioning(mesh, partitioning, bdr_partitioning);
988  }
989  if (partitioning)
990  {
991  const char part_file[] = "partitioning.txt";
992  ofstream opart(part_file);
993  opart << "number_of_elements " << mesh->GetNE() << '\n'
994  << "number_of_processors " << np << '\n';
995  for (int i = 0; i < mesh->GetNE(); i++)
996  {
997  opart << partitioning[i] << '\n';
998  }
999  cout << "Partitioning file: " << part_file << endl;
1000 
1001  Array<int> proc_el(np);
1002  proc_el = 0;
1003  for (int i = 0; i < mesh->GetNE(); i++)
1004  {
1005  proc_el[partitioning[i]]++;
1006  }
1007  int min_el = proc_el[0], max_el = proc_el[0];
1008  for (int i = 1; i < np; i++)
1009  {
1010  if (min_el > proc_el[i])
1011  {
1012  min_el = proc_el[i];
1013  }
1014  if (max_el < proc_el[i])
1015  {
1016  max_el = proc_el[i];
1017  }
1018  }
1019  cout << "Partitioning stats:\n"
1020  << " "
1021  << setw(12) << "minimum"
1022  << setw(12) << "average"
1023  << setw(12) << "maximum"
1024  << setw(12) << "total" << '\n';
1025  cout << " elements "
1026  << setw(12) << min_el
1027  << setw(12) << double(mesh->GetNE())/np
1028  << setw(12) << max_el
1029  << setw(12) << mesh->GetNE() << endl;
1030  }
1031  else
1032  {
1033  continue;
1034  }
1035 
1036  for (int i = 0; i < mesh->GetNE(); i++)
1037  {
1038  attr(i) = partitioning[i] + 1;
1039  }
1040  }
1041 
1042  char vishost[] = "localhost";
1043  int visport = 19916;
1044  socketstream sol_sock(vishost, visport);
1045  if (sol_sock.is_open())
1046  {
1047  sol_sock.precision(14);
1048  if (sdim == 2)
1049  {
1050  sol_sock << "fem2d_gf_data_keys\n";
1051  if (mk != 'p')
1052  {
1053  mesh->Print(sol_sock);
1054  }
1055  else
1056  {
1057  // NURBS meshes do not support PrintWithPartitioning
1058  if (mesh->NURBSext)
1059  {
1060  mesh->Print(sol_sock);
1061  for (int i = 0; i < mesh->GetNE(); i++)
1062  {
1063  attr(i) = partitioning[i];
1064  }
1065  }
1066  else
1067  {
1068  mesh->PrintWithPartitioning(partitioning, sol_sock, 1);
1069  }
1070  }
1071  attr.Save(sol_sock);
1072  sol_sock << "RjlmAb***********";
1073  if (mk == 'v')
1074  {
1075  sol_sock << "e";
1076  }
1077  else
1078  {
1079  sol_sock << "\n";
1080  }
1081  }
1082  else
1083  {
1084  sol_sock << "fem3d_gf_data_keys\n";
1085  if (mk == 'v' || mk == 'h' || mk == 'k' || mk == 'J')
1086  {
1087  mesh->Print(sol_sock);
1088  }
1089  else if (mk == 'b' || mk == 'B')
1090  {
1091  bdr_mesh->Print(sol_sock);
1092  bdr_attr.Save(sol_sock);
1093  sol_sock << "mcaaA";
1094  // Switch to a discrete color scale
1095  sol_sock << "pppppp" << "pppppp" << "pppppp";
1096  }
1097  else
1098  {
1099  // NURBS meshes do not support PrintWithPartitioning
1100  if (mesh->NURBSext)
1101  {
1102  mesh->Print(sol_sock);
1103  for (int i = 0; i < mesh->GetNE(); i++)
1104  {
1105  attr(i) = partitioning[i];
1106  }
1107  }
1108  else
1109  {
1110  mesh->PrintWithPartitioning(partitioning, sol_sock);
1111  }
1112  }
1113  if (mk != 'b' && mk != 'B')
1114  {
1115  attr.Save(sol_sock);
1116  sol_sock << "maaA";
1117  if (mk == 'v')
1118  {
1119  sol_sock << "aa";
1120  }
1121  else
1122  {
1123  sol_sock << "\n";
1124  }
1125  }
1126  }
1127  sol_sock << flush;
1128  }
1129  else
1130  {
1131  cout << "Unable to connect to "
1132  << vishost << ':' << visport << endl;
1133  }
1134  delete attr_fespace;
1135  delete bdr_attr_fespace;
1136  }
1137 
1138  if (mk == 'l')
1139  {
1140  // Project and plot the function 'f'
1141  int p;
1142  FiniteElementCollection *fec = NULL;
1143  cout << "Enter projection space order: " << flush;
1144  cin >> p;
1145  if (p >= 1)
1146  {
1147  fec = new H1_FECollection(p, mesh->Dimension(),
1149  }
1150  else
1151  {
1152  fec = new DG_FECollection(-p, mesh->Dimension(),
1154  }
1155  FiniteElementSpace fes(mesh, fec);
1156  GridFunction level(&fes);
1157  FunctionCoefficient coeff(f);
1158  level.ProjectCoefficient(coeff);
1159  char vishost[] = "localhost";
1160  int visport = 19916;
1161  socketstream sol_sock(vishost, visport);
1162  if (sol_sock.is_open())
1163  {
1164  sol_sock.precision(14);
1165  sol_sock << "solution\n" << *mesh << level << flush;
1166  }
1167  else
1168  {
1169  cout << "Unable to connect to "
1170  << vishost << ':' << visport << endl;
1171  }
1172  delete fec;
1173  }
1174 
1175  if (mk == 'S')
1176  {
1177  const char omesh_file[] = "mesh-explorer.mesh";
1178  ofstream omesh(omesh_file);
1179  omesh.precision(14);
1180  mesh->Print(omesh);
1181  cout << "New mesh file: " << omesh_file << endl;
1182  }
1183 
1184  if (mk == 'V')
1185  {
1186  const char omesh_file[] = "mesh-explorer.vtk";
1187  ofstream omesh(omesh_file);
1188  omesh.precision(14);
1189  mesh->PrintVTK(omesh);
1190  cout << "New VTK mesh file: " << omesh_file << endl;
1191  }
1192 
1193 #ifdef MFEM_USE_ZLIB
1194  if (mk == 'Z')
1195  {
1196  const char omesh_file[] = "mesh-explorer.mesh.gz";
1197  ofgzstream omesh(omesh_file, "zwb9");
1198  omesh.precision(14);
1199  mesh->Print(omesh);
1200  cout << "New mesh file: " << omesh_file << endl;
1201  }
1202 #endif
1203 
1204  }
1205 
1206  delete bdr_attr_fec;
1207  delete attr_fec;
1208  delete bdr_mesh;
1209  delete mesh;
1210  return 0;
1211 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
const char vishost[]
bool Help() const
Return true if we are flagged to print the help message.
Definition: optparser.hpp:153
Geometry::Type GetGeometryType() const
Definition: element.hpp:52
double GetGeckoElementOrdering(Array< int > &ordering, int iterations=4, int window=4, int period=2, int seed=0, bool verbose=false, double time_limit=0)
Definition: mesh.cpp:2030
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:573
int * CartesianPartitioning(int nxyz[])
Definition: mesh.cpp:6989
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:584
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:1012
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:251
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1674
virtual void GetVertices(Array< int > &v) const =0
Returns element&#39;s vertices.
static void PrintElementsByGeometry(int dim, const Array< int > &num_elems_by_geom, std::ostream &out)
Auxiliary method used by PrintCharacteristics().
Definition: mesh.cpp:228
static const int NumGeom
Definition: geom.hpp:42
void JacToPerfJac(int GeomType, const DenseMatrix &J, DenseMatrix &PJ) const
Definition: geom.cpp:867
int * GeneratePartitioning(int nparts, int part_method=1)
Definition: mesh.cpp:7033
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:926
double Det() const
Definition: densemat.cpp:449
void NewNodes(GridFunction &nodes, bool make_owner=false)
Replace the internal node GridFunction with the given GridFunction.
Definition: mesh.cpp:7935
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:73
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:812
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:93
static Mesh MakeRefined(Mesh &orig_mesh, int ref_factor, int ref_type)
Create a refined (by any factor) version of orig_mesh.
Definition: mesh.cpp:3733
int AddTriangle(int v1, int v2, int v3, int attr=1)
Definition: mesh.cpp:1660
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:200
void Transform(void(*f)(const Vector &, Vector &))
Definition: mesh.cpp:11629
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:923
void PrintHelp(std::ostream &out) const
Print the help message.
Definition: optparser.cpp:398
void PrintCharacteristics(Vector *Vh=NULL, Vector *Vk=NULL, std::ostream &os=mfem::out)
Compute and print mesh characteristics such as number of vertices, number of elements, number of boundary elements, minimal and maximal element sizes, minimal and maximal element aspect ratios, etc.
Definition: mesh.cpp:242
double kappa
Definition: ex24.cpp:54
void MakeOwner(FiniteElementCollection *fec_)
Make the GridFunction the owner of fec and fes.
Definition: gridfunc.hpp:124
void Randomize(int seed=0)
Set random values in the vector.
Definition: vector.cpp:785
const IntegrationPoint & GetCenter(int GeomType)
Return the center of the given Geometry::Type, GeomType.
Definition: geom.hpp:71
Mesh * read_par_mesh(int np, const char *mesh_prefix, Array< int > &partitioning, Array< int > &bdr_partitioning)
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:1069
void DeleteAll()
Delete the whole array.
Definition: array.hpp:846
Mesh * skin_mesh(Mesh *mesh)
DofTransformation * GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:297
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3676
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:614
Geometry Geometries
Definition: fe.cpp:49
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
int AddVertex(double x, double y=0.0, double z=0.0)
Definition: mesh.cpp:1613
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:10930
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
double GetElementSize(ElementTransformation *T, int type=0)
Definition: mesh.cpp:76
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:623
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:10303
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:751
void transformation(const Vector &p, Vector &v)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9816
double region_eps
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:119
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:5343
GeometryRefiner GlobGeometryRefiner
Definition: geom.cpp:1773
int GetAttribute() const
Return element&#39;s attribute.
Definition: element.hpp:55
void GetBdrElementAdjacentElement(int bdr_el, int &el, int &info) const
For the given boundary element, bdr_el, return its adjacent element and its info, i...
Definition: mesh.cpp:6206
void GetHilbertElementOrdering(Array< int > &ordering)
Definition: mesh.cpp:2197
Piecewise-constant discontinuous finite elements in 3D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:1075
void GetElementCenter(int i, Vector &center)
Definition: mesh.cpp:68
IntegrationRule RefPts
Definition: geom.hpp:282
RefinedGeometry * Refine(Geometry::Type Geom, int Times, int ETimes=1)
Definition: geom.cpp:1099
const int visport
Piecewise-constant discontinuous finite elements in 2D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:897
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:668
int Dimension() const
Definition: mesh.hpp:1006
virtual DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element &#39;elem&#39;.
Definition: fespace.cpp:2678
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
Definition: fe_coll.hpp:338
int SpaceDimension() const
Definition: mesh.hpp:1007
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:270
bool is_open()
True if the socketstream is open, false otherwise.
int AddSegment(int v1, int v2, int attr=1)
Definition: mesh.cpp:1646
void Print(std::ostream &out=mfem::out, int width=8) const
Prints vector to stream out.
Definition: vector.cpp:724
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
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:2890
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:679
void GetColumnReference(int c, Vector &col)
Definition: densemat.hpp:289
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
Definition: fe_coll.cpp:51
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:272
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:920
void PrintError(std::ostream &out) const
Print the error message.
Definition: optparser.cpp:355
virtual const char * Name() const
Definition: fe_coll.hpp:65
A standard isoparametric element transformation.
Definition: eltrans.hpp:361
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:348
double Max() const
Returns the maximal element of the vector.
Definition: vector.cpp:902
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:5338
virtual void Print(std::ostream &os=mfem::out) const
Definition: mesh.hpp:1671
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition: mesh.cpp:2983
void ReorderElements(const Array< int > &ordering, bool reorder_vertices=true)
Definition: mesh.cpp:2249
int dim
Definition: ex24.cpp:53
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:324
double CalcSingularvalue(const int i) const
Return the i-th singular value (decreasing order) of NxN matrix, N=1,2,3.
Definition: densemat.cpp:1229
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:46
virtual int GetNVertices() const =0
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: gridfunc.cpp:2399
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition: gridfunc.cpp:196
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:11898
Nodes: x_i = i/(n-1), i=0,...,n-1.
Definition: fe_base.hpp:34
Vector data type.
Definition: vector.hpp:60
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:601
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7908
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:220
void GetElementColoring(Array< int > &colors, int el0=0)
Definition: mesh.cpp:10855
RefCoord s[3]
const DenseMatrix & GetPointMat() const
Return the stored point matrix.
Definition: eltrans.hpp:404
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:9405
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:1480
int main()
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:268
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:1044
void recover_bdr_partitioning(const Mesh *mesh, const Array< int > &partitioning, Array< int > &bdr_partitioning)
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150