MFEM  v4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
mesh-explorer.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 //
12 // -----------------------------------------------------
13 // Mesh Explorer Miniapp: Explore and manipulate meshes
14 // -----------------------------------------------------
15 //
16 // This miniapp is a handy tool to examine, visualize and manipulate a given
17 // mesh. Some of its features are:
18 //
19 // - visualizing of mesh materials and individual mesh elements
20 // - mesh scaling, randomization, and general transformation
21 // - manipulation of the mesh curvature
22 // - the ability to simulate parallel partitioning
23 // - quantitative and visual reports of mesh quality
24 //
25 // Compile with: make mesh-explorer
26 //
27 // Sample runs: mesh-explorer
28 // mesh-explorer -m ../../data/beam-tri.mesh
29 // mesh-explorer -m ../../data/star-q2.mesh
30 // mesh-explorer -m ../../data/disc-nurbs.mesh
31 // mesh-explorer -m ../../data/escher-p3.mesh
32 // mesh-explorer -m ../../data/mobius-strip.mesh
33 
34 #include "mfem.hpp"
35 #include <fstream>
36 #include <limits>
37 #include <cstdlib>
38 
39 using namespace mfem;
40 using namespace std;
41 
42 // This tranformation can be applied to a mesh with the 't' menu option.
43 void transformation(const Vector &p, Vector &v)
44 {
45  // simple shear transformation
46  double s = 0.1;
47 
48  if (p.Size() == 3)
49  {
50  v(0) = p(0) + s*p(1) + s*p(2);
51  v(1) = p(1) + s*p(2) + s*p(0);
52  v(2) = p(2);
53  }
54  else if (p.Size() == 2)
55  {
56  v(0) = p(0) + s*p(1);
57  v(1) = p(1) + s*p(0);
58  }
59  else
60  {
61  v = p;
62  }
63 }
64 
65 // This function is used with the 'r' menu option, sub-option 'l' to refine a
66 // mesh locally in a region, defined by return values <= region_eps.
67 double region_eps = 1e-8;
68 double region(const Vector &p)
69 {
70  const double x = p(0), y = p(1);
71  // here we describe the region: (x <= 1/4) && (y >= 0) && (y <= 1)
72  return std::max(std::max(x - 0.25, -y), y - 1.0);
73 }
74 
75 Mesh *read_par_mesh(int np, const char *mesh_prefix)
76 {
77  Mesh *mesh;
78  Array<Mesh *> mesh_array;
79 
80  mesh_array.SetSize(np);
81  for (int p = 0; p < np; p++)
82  {
83  ostringstream fname;
84  fname << mesh_prefix << '.' << setfill('0') << setw(6) << p;
85  ifgzstream meshin(fname.str().c_str());
86  if (!meshin)
87  {
88  cerr << "Can not open mesh file: " << fname.str().c_str()
89  << '!' << endl;
90  for (p--; p >= 0; p--)
91  {
92  delete mesh_array[p];
93  }
94  return NULL;
95  }
96  mesh_array[p] = new Mesh(meshin, 1, 0);
97  // set element and boundary attributes to be the processor number + 1
98  if (1)
99  {
100  for (int i = 0; i < mesh_array[p]->GetNE(); i++)
101  {
102  mesh_array[p]->GetElement(i)->SetAttribute(p+1);
103  }
104  for (int i = 0; i < mesh_array[p]->GetNBE(); i++)
105  {
106  mesh_array[p]->GetBdrElement(i)->SetAttribute(p+1);
107  }
108  }
109  }
110  mesh = new Mesh(mesh_array, np);
111 
112  for (int p = 0; p < np; p++)
113  {
114  delete mesh_array[np-1-p];
115  }
116  mesh_array.DeleteAll();
117 
118  return mesh;
119 }
120 
121 int main (int argc, char *argv[])
122 {
123  int np = 0;
124  const char *mesh_file = "../../data/beam-hex.mesh";
125  bool refine = true;
126 
127  OptionsParser args(argc, argv);
128  args.AddOption(&mesh_file, "-m", "--mesh",
129  "Mesh file to visualize.");
130  args.AddOption(&np, "-np", "--num-proc",
131  "Load mesh from multiple processors.");
132  args.AddOption(&refine, "-ref", "--refinement", "-no-ref", "--no-refinement",
133  "Prepare the mesh for refinement or not.");
134  args.Parse();
135  if (!args.Good())
136  {
137  if (!args.Help())
138  {
139  args.PrintError(cout);
140  cout << endl;
141  }
142  cout << "Visualize and manipulate a serial mesh:\n"
143  << " mesh-explorer -m <mesh_file>\n"
144  << "Visualize and manipulate a parallel mesh:\n"
145  << " mesh-explorer -np <#proc> -m <mesh_prefix>\n" << endl
146  << "All Options:\n";
147  args.PrintHelp(cout);
148  return 1;
149  }
150  args.PrintOptions(cout);
151 
152  Mesh *mesh;
153  if (np <= 0)
154  {
155  mesh = new Mesh(mesh_file, 1, refine);
156  }
157  else
158  {
159  mesh = read_par_mesh(np, mesh_file);
160  if (mesh == NULL)
161  {
162  return 3;
163  }
164  }
165  int dim = mesh->Dimension();
166  int sdim = mesh->SpaceDimension();
167 
168  FiniteElementCollection *attr_fec;
169  if (dim == 2)
170  {
171  attr_fec = new Const2DFECollection;
172  }
173  else
174  {
175  attr_fec = new Const3DFECollection;
176  }
177 
178  int print_char = 1;
179  while (1)
180  {
181  if (print_char)
182  {
183  cout << endl;
184  mesh->PrintCharacteristics();
185  cout << "boundary attribs :";
186  for (int i = 0; i < mesh->bdr_attributes.Size(); i++)
187  {
188  cout << ' ' << mesh->bdr_attributes[i];
189  }
190  cout << '\n' << "material attribs :";
191  for (int i = 0; i < mesh->attributes.Size(); i++)
192  {
193  cout << ' ' << mesh->attributes[i];
194  }
195  cout << endl;
196  cout << "mesh curvature : ";
197  if (mesh->GetNodalFESpace() != NULL)
198  {
199  cout << mesh->GetNodalFESpace()->FEColl()->Name() << endl;
200  }
201  else
202  {
203  cout << "NONE" << endl;
204  }
205  }
206  print_char = 0;
207  cout << endl;
208  cout << "What would you like to do?\n"
209  "q) Quit\n"
210  "r) Refine\n"
211  "c) Change curvature\n"
212  "s) Scale\n"
213  "t) Transform\n"
214  "j) Jitter\n"
215  "v) View\n"
216  "m) View materials\n"
217  "b) View boundary\n"
218  "e) View elements\n"
219  "h) View element sizes, h\n"
220  "k) View element ratios, kappa\n"
221  "x) Print sub-element stats\n"
222  "f) Find physical point in reference space\n"
223  "p) Generate a partitioning\n"
224  "S) Save in MFEM format\n"
225  "V) Save in VTK format (only linear and quadratic meshes)\n"
226 #ifdef MFEM_USE_GZSTREAM
227  "Z) Save in MFEM format with compression\n"
228 #endif
229  "--> " << flush;
230  char mk;
231  cin >> mk;
232  if (!cin) { break; }
233 
234  if (mk == 'q')
235  {
236  break;
237  }
238 
239  if (mk == 'r')
240  {
241  cout <<
242  "Choose type of refinement:\n"
243  "s) standard refinement with Mesh::UniformRefinement()\n"
244  "b) Mesh::UniformRefinement() (bisection for tet meshes)\n"
245  "u) uniform refinement with a factor\n"
246  "g) non-uniform refinement (Gauss-Lobatto) with a factor\n"
247  "l) refine locally using the region() function\n"
248  "--> " << flush;
249  char sk;
250  cin >> sk;
251  switch (sk)
252  {
253  case 's':
254  mesh->UniformRefinement();
255  // Make sure tet-only meshes are marked for local refinement.
256  mesh->Finalize(true);
257  break;
258  case 'b':
259  mesh->UniformRefinement(1); // ref_algo = 1
260  break;
261  case 'u':
262  case 'g':
263  {
264  cout << "enter refinement factor --> " << flush;
265  int ref_factor;
266  cin >> ref_factor;
267  if (ref_factor <= 1 || ref_factor > 32) { break; }
268  int ref_type = (sk == 'u') ? BasisType::ClosedUniform :
270  Mesh *rmesh = new Mesh(mesh, ref_factor, ref_type);
271  delete mesh;
272  mesh = rmesh;
273  break;
274  }
275  case 'l':
276  {
277  Vector pt;
278  Array<int> marked_elements;
279  for (int i = 0; i < mesh->GetNE(); i++)
280  {
281  // check all nodes of the element
283  mesh->GetElementTransformation(i, &T);
284  for (int j = 0; j < T.GetPointMat().Width(); j++)
285  {
286  T.GetPointMat().GetColumnReference(j, pt);
287  if (region(pt) <= region_eps)
288  {
289  marked_elements.Append(i);
290  break;
291  }
292  }
293  }
294  mesh->GeneralRefinement(marked_elements);
295  break;
296  }
297  }
298  print_char = 1;
299  }
300 
301  if (mk == 'c')
302  {
303  int p;
304  cout << "enter new order for mesh curvature --> " << flush;
305  cin >> p;
306  mesh->SetCurvature(p > 0 ? p : -p, p <= 0);
307  print_char = 1;
308  }
309 
310  if (mk == 's')
311  {
312  double factor;
313  cout << "scaling factor ---> " << flush;
314  cin >> factor;
315 
316  GridFunction *nodes = mesh->GetNodes();
317  if (nodes == NULL)
318  {
319  for (int i = 0; i < mesh->GetNV(); i++)
320  {
321  double *v = mesh->GetVertex(i);
322  v[0] *= factor;
323  v[1] *= factor;
324  if (dim == 3)
325  {
326  v[2] *= factor;
327  }
328  }
329  }
330  else
331  {
332  *nodes *= factor;
333  }
334 
335  print_char = 1;
336  }
337 
338  if (mk == 't')
339  {
340  mesh->Transform(transformation);
341  print_char = 1;
342  }
343 
344  if (mk == 'j')
345  {
346  double jitter;
347  cout << "jitter factor ---> " << flush;
348  cin >> jitter;
349 
350  GridFunction *nodes = mesh->GetNodes();
351 
352  if (nodes == NULL)
353  {
354  cerr << "The mesh should have nodes, introduce curvature first!\n";
355  }
356  else
357  {
358  FiniteElementSpace *fespace = nodes->FESpace();
359 
360  GridFunction rdm(fespace);
361  rdm.Randomize();
362  rdm -= 0.5; // shift to random values in [-0.5,0.5]
363  rdm *= jitter;
364 
365  // compute minimal local mesh size
366  Vector h0(fespace->GetNDofs());
367  h0 = infinity();
368  {
369  Array<int> dofs;
370  for (int i = 0; i < fespace->GetNE(); i++)
371  {
372  fespace->GetElementDofs(i, dofs);
373  for (int j = 0; j < dofs.Size(); j++)
374  {
375  h0(dofs[j]) = std::min(h0(dofs[j]), mesh->GetElementSize(i));
376  }
377  }
378  }
379 
380  // scale the random values to be of order of the local mesh size
381  for (int i = 0; i < fespace->GetNDofs(); i++)
382  {
383  for (int d = 0; d < dim; d++)
384  {
385  rdm(fespace->DofToVDof(i,d)) *= h0(i);
386  }
387  }
388 
389  char move_bdr = 'n';
390  cout << "move boundary nodes? [y/n] ---> " << flush;
391  cin >> move_bdr;
392 
393  // don't perturb the boundary
394  if (move_bdr == 'n')
395  {
396  Array<int> vdofs;
397  for (int i = 0; i < fespace->GetNBE(); i++)
398  {
399  fespace->GetBdrElementVDofs(i, vdofs);
400  for (int j = 0; j < vdofs.Size(); j++)
401  {
402  rdm(vdofs[j]) = 0.0;
403  }
404  }
405  }
406 
407  *nodes += rdm;
408  }
409 
410  print_char = 1;
411  }
412 
413  if (mk == 'x')
414  {
415  int sd, nz = 0;
416  DenseMatrix J(dim);
417  double min_det_J, max_det_J, min_det_J_z, max_det_J_z;
418  double min_kappa, max_kappa, max_ratio_det_J_z;
419  min_det_J = min_kappa = infinity();
420  max_det_J = max_kappa = max_ratio_det_J_z = -infinity();
421  cout << "subdivision factor ---> " << flush;
422  cin >> sd;
423  Array<int> bad_elems_by_geom(Geometry::NumGeom);
424  bad_elems_by_geom = 0;
425  for (int i = 0; i < mesh->GetNE(); i++)
426  {
429 
430  RefinedGeometry *RefG = GlobGeometryRefiner.Refine(geom, sd, 1);
431  IntegrationRule &ir = RefG->RefPts;
432 
433  min_det_J_z = infinity();
434  max_det_J_z = -infinity();
435  for (int j = 0; j < ir.GetNPoints(); j++)
436  {
437  T->SetIntPoint(&ir.IntPoint(j));
438  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
439 
440  double det_J = J.Det();
441  double kappa =
442  J.CalcSingularvalue(0) / J.CalcSingularvalue(dim-1);
443 
444  min_det_J_z = fmin(min_det_J_z, det_J);
445  max_det_J_z = fmax(max_det_J_z, det_J);
446 
447  min_kappa = fmin(min_kappa, kappa);
448  max_kappa = fmax(max_kappa, kappa);
449  }
450  max_ratio_det_J_z =
451  fmax(max_ratio_det_J_z, max_det_J_z/min_det_J_z);
452  min_det_J = fmin(min_det_J, min_det_J_z);
453  max_det_J = fmax(max_det_J, max_det_J_z);
454  if (min_det_J_z <= 0.0)
455  {
456  nz++;
457  bad_elems_by_geom[geom]++;
458  }
459  }
460  cout << "\nbad elements = " << nz;
461  if (nz)
462  {
463  cout << " -- ";
464  Mesh::PrintElementsByGeometry(dim, bad_elems_by_geom, cout);
465  }
466  cout << "\nmin det(J) = " << min_det_J
467  << "\nmax det(J) = " << max_det_J
468  << "\nglobal ratio = " << max_det_J/min_det_J
469  << "\nmax el ratio = " << max_ratio_det_J_z
470  << "\nmin kappa = " << min_kappa
471  << "\nmax kappa = " << max_kappa << endl;
472  }
473 
474  if (mk == 'f')
475  {
476  DenseMatrix point_mat(sdim,1);
477  cout << "\npoint in physical space ---> " << flush;
478  for (int i = 0; i < sdim; i++)
479  {
480  cin >> point_mat(i,0);
481  }
482  Array<int> elem_ids;
484 
485  // physical -> reference space
486  mesh->FindPoints(point_mat, elem_ids, ips);
487 
488  cout << "point in reference space:";
489  if (elem_ids[0] == -1)
490  {
491  cout << " NOT FOUND!\n";
492  }
493  else
494  {
495  cout << " element " << elem_ids[0] << ", ip =";
496  cout << " " << ips[0].x;
497  if (sdim > 1)
498  {
499  cout << " " << ips[0].y;
500  if (sdim > 2)
501  {
502  cout << " " << ips[0].z;
503  }
504  }
505  cout << endl;
506  }
507  }
508 
509  // These are the cases that open a new GLVis window
510  if (mk == 'm' || mk == 'b' || mk == 'e' || mk == 'v' || mk == 'h' ||
511  mk == 'k' || mk == 'p')
512  {
513  Array<int> part(mesh->GetNE());
514  FiniteElementSpace *attr_fespace =
515  new FiniteElementSpace(mesh, attr_fec);
516  GridFunction attr(attr_fespace);
517 
518  if (mk == 'm')
519  for (int i = 0; i < mesh->GetNE(); i++)
520  {
521  part[i] = (attr(i) = mesh->GetAttribute(i)) - 1;
522  }
523 
524  if (mk == 'b' || mk == 'v')
525  {
526  attr = 1.0;
527  }
528 
529  if (mk == 'e')
530  {
531  Array<int> coloring;
532  srand(time(0));
533  double a = double(rand()) / (double(RAND_MAX) + 1.);
534  int el0 = (int)floor(a * mesh->GetNE());
535  cout << "Generating coloring starting with element " << el0+1
536  << " / " << mesh->GetNE() << endl;
537  mesh->GetElementColoring(coloring, el0);
538  for (int i = 0; i < coloring.Size(); i++)
539  {
540  attr(i) = coloring[i];
541  }
542  cout << "Number of colors: " << attr.Max() + 1 << endl;
543  for (int i = 0; i < mesh->GetNE(); i++)
544  {
545  // part[i] = i; // checkerboard element coloring
546  attr(i) = part[i] = i; // coloring by element number
547  }
548  }
549 
550  if (mk == 'h')
551  {
552  DenseMatrix J(dim);
553  double h_min, h_max;
554  h_min = infinity();
555  h_max = -h_min;
556  for (int i = 0; i < mesh->GetNE(); i++)
557  {
558  int geom = mesh->GetElementBaseGeometry(i);
560  T->SetIntPoint(&Geometries.GetCenter(geom));
561  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
562 
563  attr(i) = J.Det();
564  if (attr(i) < 0.0)
565  {
566  attr(i) = -pow(-attr(i), 1.0/double(dim));
567  }
568  else
569  {
570  attr(i) = pow(attr(i), 1.0/double(dim));
571  }
572  h_min = min(h_min, attr(i));
573  h_max = max(h_max, attr(i));
574  }
575  cout << "h_min = " << h_min << ", h_max = " << h_max << endl;
576  }
577 
578  if (mk == 'k')
579  {
580  DenseMatrix J(dim);
581  for (int i = 0; i < mesh->GetNE(); i++)
582  {
583  int geom = mesh->GetElementBaseGeometry(i);
585  T->SetIntPoint(&Geometries.GetCenter(geom));
586  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
587  attr(i) = J.CalcSingularvalue(0) / J.CalcSingularvalue(dim-1);
588  }
589  }
590 
591  if (mk == 'p')
592  {
593  int *partitioning = NULL, n;
594  cout << "What type of partitioning?\n"
595  "c) Cartesian\n"
596  "0) METIS_PartGraphRecursive (sorted neighbor lists)\n"
597  "1) METIS_PartGraphKway (sorted neighbor lists)"
598  " (default)\n"
599  "2) METIS_PartGraphVKway (sorted neighbor lists)\n"
600  "3) METIS_PartGraphRecursive\n"
601  "4) METIS_PartGraphKway\n"
602  "5) METIS_PartGraphVKway\n"
603  "--> " << flush;
604  char pk;
605  cin >> pk;
606  if (pk == 'c')
607  {
608  int nxyz[3];
609  cout << "Enter nx: " << flush;
610  cin >> nxyz[0]; n = nxyz[0];
611  if (mesh->Dimension() > 1)
612  {
613  cout << "Enter ny: " << flush;
614  cin >> nxyz[1]; n *= nxyz[1];
615  if (mesh->Dimension() > 2)
616  {
617  cout << "Enter nz: " << flush;
618  cin >> nxyz[2]; n *= nxyz[2];
619  }
620  }
621  partitioning = mesh->CartesianPartitioning(nxyz);
622  }
623  else
624  {
625  int part_method = pk - '0';
626  if (part_method < 0 || part_method > 5)
627  {
628  continue;
629  }
630  cout << "Enter number of processors: " << flush;
631  cin >> n;
632  partitioning = mesh->GeneratePartitioning(n, part_method);
633  }
634  if (partitioning)
635  {
636  const char part_file[] = "partitioning.txt";
637  ofstream opart(part_file);
638  opart << "number_of_elements " << mesh->GetNE() << '\n'
639  << "number_of_processors " << n << '\n';
640  for (int i = 0; i < mesh->GetNE(); i++)
641  {
642  opart << partitioning[i] << '\n';
643  }
644  cout << "Partitioning file: " << part_file << endl;
645 
646  Array<int> proc_el(n);
647  proc_el = 0;
648  for (int i = 0; i < mesh->GetNE(); i++)
649  {
650  proc_el[partitioning[i]]++;
651  }
652  int min_el = proc_el[0], max_el = proc_el[0];
653  for (int i = 1; i < n; i++)
654  {
655  if (min_el > proc_el[i])
656  {
657  min_el = proc_el[i];
658  }
659  if (max_el < proc_el[i])
660  {
661  max_el = proc_el[i];
662  }
663  }
664  cout << "Partitioning stats:\n"
665  << " "
666  << setw(12) << "minimum"
667  << setw(12) << "average"
668  << setw(12) << "maximum"
669  << setw(12) << "total" << '\n';
670  cout << " elements "
671  << setw(12) << min_el
672  << setw(12) << double(mesh->GetNE())/n
673  << setw(12) << max_el
674  << setw(12) << mesh->GetNE() << endl;
675  }
676  else
677  {
678  continue;
679  }
680 
681  for (int i = 0; i < mesh->GetNE(); i++)
682  {
683  attr(i) = part[i] = partitioning[i];
684  }
685  delete [] partitioning;
686  }
687 
688  char vishost[] = "localhost";
689  int visport = 19916;
690  socketstream sol_sock(vishost, visport);
691  if (sol_sock.is_open())
692  {
693  sol_sock.precision(14);
694  if (sdim == 2)
695  {
696  sol_sock << "fem2d_gf_data_keys\n";
697  if (mk != 'p')
698  {
699  mesh->Print(sol_sock);
700  }
701  else
702  {
703  // NURBS meshes do not support PrintWithPartitioning
704  if (mesh->NURBSext)
705  {
706  mesh->Print(sol_sock);
707  for (int i = 0; i < mesh->GetNE(); i++)
708  {
709  attr(i) = part[i];
710  }
711  }
712  else
713  {
714  mesh->PrintWithPartitioning(part, sol_sock, 1);
715  }
716  }
717  attr.Save(sol_sock);
718  sol_sock << "RjlmAb***********";
719  if (mk == 'v')
720  {
721  sol_sock << "e";
722  }
723  else
724  {
725  sol_sock << "\n";
726  }
727  }
728  else
729  {
730  sol_sock << "fem3d_gf_data_keys\n";
731  if (mk == 'b' || mk == 'v' || mk == 'h' || mk == 'k')
732  {
733  mesh->Print(sol_sock);
734  }
735  else
736  {
737  // NURBS meshes do not support PrintWithPartitioning
738  if (mesh->NURBSext)
739  {
740  mesh->Print(sol_sock);
741  for (int i = 0; i < mesh->GetNE(); i++)
742  {
743  attr(i) = part[i];
744  }
745  }
746  else
747  {
748  mesh->PrintWithPartitioning(part, sol_sock);
749  }
750  }
751  attr.Save(sol_sock);
752  sol_sock << "maaA";
753  if (mk == 'v')
754  {
755  sol_sock << "aa";
756  }
757  else
758  {
759  sol_sock << "\n";
760  }
761  }
762  sol_sock << flush;
763  }
764  else
765  {
766  cout << "Unable to connect to "
767  << vishost << ':' << visport << endl;
768  }
769  delete attr_fespace;
770  }
771 
772  if (mk == 'S')
773  {
774  const char mesh_file[] = "mesh-explorer.mesh";
775  ofstream omesh(mesh_file);
776  omesh.precision(14);
777  mesh->Print(omesh);
778  cout << "New mesh file: " << mesh_file << endl;
779  }
780 
781  if (mk == 'V')
782  {
783  const char mesh_file[] = "mesh-explorer.vtk";
784  ofstream omesh(mesh_file);
785  omesh.precision(14);
786  mesh->PrintVTK(omesh);
787  cout << "New VTK mesh file: " << mesh_file << endl;
788  }
789 
790 #ifdef MFEM_USE_GZSTREAM
791  if (mk == 'Z')
792  {
793  const char mesh_file[] = "mesh-explorer.mesh.gz";
794  ofgzstream omesh(mesh_file, "zwb9");
795  omesh.precision(14);
796  mesh->Print(omesh);
797  cout << "New mesh file: " << mesh_file << endl;
798  }
799 #endif
800 
801  }
802 
803  delete attr_fec;
804  delete mesh;
805  return 0;
806 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:237
bool Help() const
Definition: optparser.hpp:123
int Size() const
Logical size of the array.
Definition: array.hpp:118
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1156
int * CartesianPartitioning(int nxyz[])
Definition: mesh.cpp:4977
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:344
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:85
DenseMatrix & GetPointMat()
Read and write access to the underlying point matrix describing the transformation.
Definition: eltrans.hpp:313
const double * GetVertex(int i) const
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:719
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:143
Mesh * read_par_mesh(int np, const char *mesh_prefix)
static void PrintElementsByGeometry(int dim, const Array< int > &num_elems_by_geom, std::ostream &out)
Auxiliary method used by PrintCharacteristics().
Definition: mesh.cpp:217
static const int NumGeom
Definition: geom.hpp:41
void JacToPerfJac(int GeomType, const DenseMatrix &J, DenseMatrix &PJ) const
Definition: geom.cpp:737
int * GeneratePartitioning(int nparts, int part_method=1)
Definition: mesh.cpp:5016
void PrintWithPartitioning(int *partitioning, std::ostream &out, int elem_attr=0) const
Prints the mesh with boundary elements given by the boundary of the subdomains, so that the boundary ...
Definition: mesh.cpp:8482
const Geometry::Type geom
Definition: ex1.cpp:40
double Det() const
Definition: densemat.cpp:472
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:42
void SetIntPoint(const IntegrationPoint *ip)
Definition: eltrans.hpp:53
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:150
void Transform(void(*f)(const Vector &, Vector &))
Definition: mesh.cpp:9180
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:676
void PrintHelp(std::ostream &out) const
Definition: optparser.cpp:378
void Randomize(int seed=0)
Set random values in the vector.
Definition: vector.cpp:690
const IntegrationPoint & GetCenter(int GeomType)
Return the center of the given Geometry::Type, GeomType.
Definition: geom.hpp:70
int main(int argc, char *argv[])
Definition: ex1.cpp:58
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:759
void DeleteAll()
Delete whole array.
Definition: array.hpp:785
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2627
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:374
Geometry Geometries
Definition: fe.cpp:11972
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:240
int dim
Definition: ex3.cpp:48
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:383
int Append(const T &el)
Append element to array, resize if necessary.
Definition: array.hpp:690
void transformation(const Vector &p, Vector &v)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7591
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:68
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3644
virtual void GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:1439
GeometryRefiner GlobGeometryRefiner
Definition: geom.cpp:1395
IntegrationRule RefPts
Definition: geom.hpp:239
RefinedGeometry * Refine(Geometry::Type Geom, int Times, int ETimes=1)
Definition: geom.cpp:925
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:430
int Dimension() const
Definition: mesh.hpp:713
double GetElementSize(int i, int type=0)
Get the size of the i-th element relative to the perfect reference element.
Definition: mesh.cpp:74
int SpaceDimension() const
Definition: mesh.hpp:714
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:179
void PrintVTK(std::ostream &out)
Definition: mesh.cpp:8060
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:85
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)
Definition: optparser.hpp:76
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:618
void GetColumnReference(int c, Vector &col)
Definition: densemat.hpp:262
Nodes: x_i = i/(n-1), i=0,...,n-1.
Definition: fe.hpp:37
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:181
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:673
void PrintError(std::ostream &out) const
Definition: optparser.cpp:335
virtual const char * Name() const
Definition: fe_coll.hpp:53
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:336
double Max() const
Returns the maximal element of the vector.
Definition: vector.cpp:816
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:3639
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition: mesh.cpp:2098
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
double CalcSingularvalue(const int i) const
Return the i-th singular value (decreasing order) of NxN matrix, N=1,2,3.
Definition: densemat.cpp:1842
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:42
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:9449
double kappa
Definition: ex3.cpp:47
Vector data type.
Definition: vector.hpp:48
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:361
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5837
void GetElementColoring(Array< int > &colors, int el0=0)
Definition: mesh.cpp:8407
void GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:178
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:7185
double region(const Vector &p)
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:970
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:177
void PrintCharacteristics(Vector *Vh=NULL, Vector *Vk=NULL, std::ostream &out=mfem::out)
Compute and print mesh characteristics such as number of vertices, number of elements, number of boundary elements, minimal and maximal element sizes, minimal and maximal element aspect ratios, etc.
Definition: mesh.cpp:231
bool Good() const
Definition: optparser.hpp:122