MFEM  v3.3
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 void transformation(const Vector &p, Vector &v)
43 {
44  // simple shear transformation
45  double s = 0.1;
46  v(0) = p(0) + s*p(1) + s*p(2);
47  v(1) = p(1) + s*p(2) + s*p(0);
48  v(2) = p(2);
49 }
50 
51 Mesh *read_par_mesh(int np, const char *mesh_prefix)
52 {
53  Mesh *mesh;
54  Array<Mesh *> mesh_array;
55 
56  mesh_array.SetSize(np);
57  for (int p = 0; p < np; p++)
58  {
59  ostringstream fname;
60  fname << mesh_prefix << '.' << setfill('0') << setw(6) << p;
61  ifgzstream meshin(fname.str().c_str());
62  if (!meshin)
63  {
64  cerr << "Can not open mesh file: " << fname.str().c_str()
65  << '!' << endl;
66  for (p--; p >= 0; p--)
67  {
68  delete mesh_array[p];
69  }
70  return NULL;
71  }
72  mesh_array[p] = new Mesh(meshin, 1, 0);
73  // set element and boundary attributes to be the processor number + 1
74  if (1)
75  {
76  for (int i = 0; i < mesh_array[p]->GetNE(); i++)
77  {
78  mesh_array[p]->GetElement(i)->SetAttribute(p+1);
79  }
80  for (int i = 0; i < mesh_array[p]->GetNBE(); i++)
81  {
82  mesh_array[p]->GetBdrElement(i)->SetAttribute(p+1);
83  }
84  }
85  }
86  mesh = new Mesh(mesh_array, np);
87 
88  for (int p = 0; p < np; p++)
89  {
90  delete mesh_array[np-1-p];
91  }
92  mesh_array.DeleteAll();
93 
94  return mesh;
95 }
96 
97 int main (int argc, char *argv[])
98 {
99  int np = 0;
100  const char *mesh_file = "../../data/beam-hex.mesh";
101 
102  OptionsParser args(argc, argv);
103  args.AddOption(&mesh_file, "-m", "--mesh",
104  "Mesh file to visualize.");
105  args.AddOption(&np, "-np", "--num-proc",
106  "Load mesh from multiple processors.");
107  args.Parse();
108  if (!args.Good())
109  {
110  if (!args.Help())
111  {
112  args.PrintError(cout);
113  cout << endl;
114  }
115  cout << "Visualize and manipulate a serial mesh:\n"
116  << " mesh-explorer -m <mesh_file>\n"
117  << "Visualize and manipulate a parallel mesh:\n"
118  << " mesh-explorer -np <#proc> -m <mesh_prefix>\n" << endl
119  << "All Options:\n";
120  args.PrintHelp(cout);
121  return 1;
122  }
123  args.PrintOptions(cout);
124 
125  Mesh *mesh;
126  if (np <= 0)
127  {
128  mesh = new Mesh(mesh_file, 1, 1);
129  }
130  else
131  {
132  mesh = read_par_mesh(np, mesh_file);
133  if (mesh == NULL)
134  {
135  return 3;
136  }
137  }
138  int dim = mesh->Dimension();
139  int sdim = mesh->SpaceDimension();
140 
141  FiniteElementCollection *attr_fec;
142  if (dim == 2)
143  {
144  attr_fec = new Const2DFECollection;
145  }
146  else
147  {
148  attr_fec = new Const3DFECollection;
149  }
150 
151  int print_char = 1;
152  while (1)
153  {
154  if (print_char)
155  {
156  cout << endl;
157  mesh->PrintCharacteristics();
158  cout << "boundary attribs :";
159  for (int i = 0; i < mesh->bdr_attributes.Size(); i++)
160  {
161  cout << ' ' << mesh->bdr_attributes[i];
162  }
163  cout << '\n' << "material attribs :";
164  for (int i = 0; i < mesh->attributes.Size(); i++)
165  {
166  cout << ' ' << mesh->attributes[i];
167  }
168  cout << endl;
169  cout << "mesh curvature : ";
170  if (mesh->GetNodalFESpace() != NULL)
171  {
172  cout << mesh->GetNodalFESpace()->FEColl()->Name() << endl;
173  }
174  else
175  {
176  cout << "NONE" << endl;
177  }
178  }
179  print_char = 0;
180  cout << endl;
181  cout << "What would you like to do?\n"
182  "q) Quit\n"
183  "r) Refine\n"
184  "c) Change curvature\n"
185  "s) Scale\n"
186  "t) Transform\n"
187  "j) Jitter\n"
188  "v) View\n"
189  "m) View materials\n"
190  "b) View boundary\n"
191  "e) View elements\n"
192  "h) View element sizes, h\n"
193  "k) View element ratios, kappa\n"
194  "x) Print sub-element stats\n"
195  "p) Generate a partitioning\n"
196  "S) Save\n"
197  "--> " << flush;
198  char mk;
199  cin >> mk;
200 
201  if (mk == 'q')
202  {
203  break;
204  }
205 
206  if (mk == 'r')
207  {
208  cout <<
209  "Choose type of refinement:\n"
210  "s) standard refinement with Mesh::UniformRefinement()\n"
211  "u) uniform refinement with a factor\n"
212  "g) non-uniform refinement (Gauss-Lobatto) with a factor\n"
213  "--> " << flush;
214  char sk;
215  cin >> sk;
216  switch (sk)
217  {
218  case 's':
219  mesh->UniformRefinement();
220  break;
221  case 'u':
222  case 'g':
223  {
224  cout << "enter refinement factor --> " << flush;
225  int ref_factor;
226  cin >> ref_factor;
227  if (ref_factor <= 1 || ref_factor > 32) { break; }
228  int ref_type = (sk == 'u') ? BasisType::ClosedUniform :
230  Mesh *rmesh = new Mesh(mesh, ref_factor, ref_type);
231  delete mesh;
232  mesh = rmesh;
233  break;
234  }
235  }
236  print_char = 1;
237  }
238 
239  if (mk == 'c')
240  {
241  int p;
242  cout << "enter new order for mesh curvature --> " << flush;
243  cin >> p;
244  mesh->SetCurvature(p);
245  print_char = 1;
246  }
247 
248  if (mk == 's')
249  {
250  double factor;
251  cout << "scaling factor ---> " << flush;
252  cin >> factor;
253 
254  GridFunction *nodes = mesh->GetNodes();
255  if (nodes == NULL)
256  {
257  for (int i = 0; i < mesh->GetNV(); i++)
258  {
259  double *v = mesh->GetVertex(i);
260  v[0] *= factor;
261  v[1] *= factor;
262  if (dim == 3)
263  {
264  v[2] *= factor;
265  }
266  }
267  }
268  else
269  {
270  *nodes *= factor;
271  }
272 
273  print_char = 1;
274  }
275 
276  if (mk == 't')
277  {
278  mesh->Transform(transformation);
279  print_char = 1;
280  }
281 
282  if (mk == 'j')
283  {
284  double jitter;
285  cout << "jitter factor ---> " << flush;
286  cin >> jitter;
287 
288  GridFunction *nodes = mesh->GetNodes();
289 
290  if (nodes == NULL)
291  {
292  cerr << "The mesh should have nodes, introduce curvature first!\n";
293  }
294  else
295  {
296  FiniteElementSpace *fespace = nodes->FESpace();
297 
298  GridFunction rdm(fespace);
299  rdm.Randomize();
300  rdm -= 0.5; // shift to random values in [-0.5,0.5]
301  rdm *= jitter;
302 
303  // compute minimal local mesh size
304  Vector h0(fespace->GetNDofs());
305  h0 = std::numeric_limits<double>::infinity();
306  {
307  Array<int> dofs;
308  for (int i = 0; i < fespace->GetNE(); i++)
309  {
310  fespace->GetElementDofs(i, dofs);
311  for (int j = 0; j < dofs.Size(); j++)
312  {
313  h0(dofs[j]) = std::min(h0(dofs[j]), mesh->GetElementSize(i));
314  }
315  }
316  }
317 
318  // scale the random values to be of order of the local mesh size
319  for (int i = 0; i < fespace->GetNDofs(); i++)
320  {
321  for (int d = 0; d < dim; d++)
322  {
323  rdm(fespace->DofToVDof(i,d)) *= h0(i);
324  }
325  }
326 
327  int bdr = 0;
328  // don't perturb the boundary
329  if (!bdr)
330  {
331  Array<int> vdofs;
332  for (int i = 0; i < fespace->GetNBE(); i++)
333  {
334  fespace->GetBdrElementVDofs(i, vdofs);
335  for (int j = 0; j < vdofs.Size(); j++)
336  {
337  rdm(vdofs[j]) = 0.0;
338  }
339  }
340  }
341 
342  *nodes += rdm;
343  }
344 
345  print_char = 1;
346  }
347 
348  if (mk == 'x')
349  {
350  int sd, nz = 0;
351  DenseMatrix J(dim);
352  double min_det_J, max_det_J, min_det_J_z, max_det_J_z;
353  double min_kappa, max_kappa, max_ratio_det_J_z;
354  min_det_J = min_kappa = numeric_limits<double>::infinity();
355  max_det_J = max_kappa = max_ratio_det_J_z = -min_det_J;
356  cout << "subdivision factor ---> " << flush;
357  cin >> sd;
358  for (int i = 0; i < mesh->GetNE(); i++)
359  {
360  int geom = mesh->GetElementBaseGeometry(i);
362 
363  RefinedGeometry *RefG = GlobGeometryRefiner.Refine(geom, sd, 1);
364  IntegrationRule &ir = RefG->RefPts;
365 
366  min_det_J_z = numeric_limits<double>::infinity();
367  max_det_J_z = -min_det_J_z;
368  for (int j = 0; j < ir.GetNPoints(); j++)
369  {
370  T->SetIntPoint(&ir.IntPoint(j));
371  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
372 
373  double det_J = J.Det();
374  double kappa =
375  J.CalcSingularvalue(0) / J.CalcSingularvalue(dim-1);
376 
377  min_det_J_z = fmin(min_det_J_z, det_J);
378  max_det_J_z = fmax(max_det_J_z, det_J);
379 
380  min_kappa = fmin(min_kappa, kappa);
381  max_kappa = fmax(max_kappa, kappa);
382  }
383  max_ratio_det_J_z =
384  fmax(max_ratio_det_J_z, max_det_J_z/min_det_J_z);
385  min_det_J = fmin(min_det_J, min_det_J_z);
386  max_det_J = fmax(max_det_J, max_det_J_z);
387  if (min_det_J_z <= 0.0)
388  {
389  nz++;
390  }
391  }
392  cout
393  << "\nbad elements = " << nz
394  << "\nmin det(J) = " << min_det_J
395  << "\nmax det(J) = " << max_det_J
396  << "\nglobal ratio = " << max_det_J/min_det_J
397  << "\nmax el ratio = " << max_ratio_det_J_z
398  << "\nmin kappa = " << min_kappa
399  << "\nmax kappa = " << max_kappa << endl;
400  }
401 
402  if (mk == 'm' || mk == 'b' || mk == 'e' || mk == 'v' || mk == 'h' ||
403  mk == 'k' || mk == 'p')
404  {
405  Array<int> part(mesh->GetNE());
406  FiniteElementSpace *attr_fespace =
407  new FiniteElementSpace(mesh, attr_fec);
408  GridFunction attr(attr_fespace);
409 
410  if (mk == 'm')
411  for (int i = 0; i < mesh->GetNE(); i++)
412  {
413  part[i] = (attr(i) = mesh->GetAttribute(i)) - 1;
414  }
415 
416  if (mk == 'b' || mk == 'v')
417  {
418  attr = 1.0;
419  }
420 
421  if (mk == 'e')
422  {
423  Array<int> coloring;
424  srand(time(0));
425  double a = double(rand()) / (double(RAND_MAX) + 1.);
426  int el0 = (int)floor(a * mesh->GetNE());
427  cout << "Generating coloring starting with element " << el0+1
428  << " / " << mesh->GetNE() << endl;
429  mesh->GetElementColoring(coloring, el0);
430  for (int i = 0; i < coloring.Size(); i++)
431  {
432  attr(i) = coloring[i];
433  }
434  cout << "Number of colors: " << attr.Max() + 1 << endl;
435  for (int i = 0; i < mesh->GetNE(); i++)
436  {
437  // part[i] = i; // checkerboard element coloring
438  attr(i) = part[i] = i; // coloring by element number
439  }
440  }
441 
442  if (mk == 'h')
443  {
444  DenseMatrix J(dim);
445  double h_min, h_max;
446  h_min = numeric_limits<double>::infinity();
447  h_max = -h_min;
448  for (int i = 0; i < mesh->GetNE(); i++)
449  {
450  int geom = mesh->GetElementBaseGeometry(i);
452  T->SetIntPoint(&Geometries.GetCenter(geom));
453  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
454 
455  attr(i) = J.Det();
456  if (attr(i) < 0.0)
457  {
458  attr(i) = -pow(-attr(i), 1.0/double(dim));
459  }
460  else
461  {
462  attr(i) = pow(attr(i), 1.0/double(dim));
463  }
464  h_min = min(h_min, attr(i));
465  h_max = max(h_max, attr(i));
466  }
467  cout << "h_min = " << h_min << ", h_max = " << h_max << endl;
468  }
469 
470  if (mk == 'k')
471  {
472  DenseMatrix J(dim);
473  for (int i = 0; i < mesh->GetNE(); i++)
474  {
475  int geom = mesh->GetElementBaseGeometry(i);
477  T->SetIntPoint(&Geometries.GetCenter(geom));
478  Geometries.JacToPerfJac(geom, T->Jacobian(), J);
479  attr(i) = J.CalcSingularvalue(0) / J.CalcSingularvalue(dim-1);
480  }
481  }
482 
483  if (mk == 'p')
484  {
485  int *partitioning = NULL, n;
486  cout << "What type of partitioning?\n"
487  "c) Cartesian\n"
488  "0) METIS_PartGraphRecursive (sorted neighbor lists)\n"
489  "1) METIS_PartGraphKway (sorted neighbor lists)\n"
490  "2) METIS_PartGraphVKway (sorted neighbor lists)\n"
491  "3) METIS_PartGraphRecursive\n"
492  "4) METIS_PartGraphKway\n"
493  "5) METIS_PartGraphVKway\n"
494  "--> " << flush;
495  char pk;
496  cin >> pk;
497  if (pk == 'c')
498  {
499  int nxyz[3];
500  cout << "Enter nx: " << flush;
501  cin >> nxyz[0]; n = nxyz[0];
502  if (mesh->Dimension() > 1)
503  {
504  cout << "Enter ny: " << flush;
505  cin >> nxyz[1]; n *= nxyz[1];
506  if (mesh->Dimension() > 2)
507  {
508  cout << "Enter nz: " << flush;
509  cin >> nxyz[2]; n *= nxyz[2];
510  }
511  }
512  partitioning = mesh->CartesianPartitioning(nxyz);
513  }
514  else
515  {
516  int part_method = pk - '0';
517  if (part_method < 0 || part_method > 5)
518  {
519  continue;
520  }
521  cout << "Enter number of processors: " << flush;
522  cin >> n;
523  partitioning = mesh->GeneratePartitioning(n, part_method);
524  }
525  if (partitioning)
526  {
527  const char part_file[] = "partitioning.txt";
528  ofstream opart(part_file);
529  opart << "number_of_elements " << mesh->GetNE() << '\n'
530  << "number_of_processors " << n << '\n';
531  for (int i = 0; i < mesh->GetNE(); i++)
532  {
533  opart << partitioning[i] << '\n';
534  }
535  cout << "Partitioning file: " << part_file << endl;
536 
537  Array<int> proc_el(n);
538  proc_el = 0;
539  for (int i = 0; i < mesh->GetNE(); i++)
540  {
541  proc_el[partitioning[i]]++;
542  }
543  int min_el = proc_el[0], max_el = proc_el[0];
544  for (int i = 1; i < n; i++)
545  {
546  if (min_el > proc_el[i])
547  {
548  min_el = proc_el[i];
549  }
550  if (max_el < proc_el[i])
551  {
552  max_el = proc_el[i];
553  }
554  }
555  cout << "Partitioning stats:\n"
556  << " "
557  << setw(12) << "minimum"
558  << setw(12) << "average"
559  << setw(12) << "maximum"
560  << setw(12) << "total" << '\n';
561  cout << " elements "
562  << setw(12) << min_el
563  << setw(12) << double(mesh->GetNE())/n
564  << setw(12) << max_el
565  << setw(12) << mesh->GetNE() << endl;
566  }
567  else
568  {
569  continue;
570  }
571 
572  for (int i = 0; i < mesh->GetNE(); i++)
573  {
574  attr(i) = part[i] = partitioning[i];
575  }
576  delete [] partitioning;
577  }
578 
579  char vishost[] = "localhost";
580  int visport = 19916;
581  socketstream sol_sock(vishost, visport);
582  if (sol_sock.is_open())
583  {
584  sol_sock.precision(14);
585  if (sdim == 2)
586  {
587  sol_sock << "fem2d_gf_data_keys\n";
588  if (mk != 'p')
589  {
590  mesh->Print(sol_sock);
591  }
592  else
593  {
594  // NURBS meshes do not support PrintWithPartitioning
595  if (mesh->NURBSext)
596  {
597  mesh->Print(sol_sock);
598  for (int i = 0; i < mesh->GetNE(); i++)
599  {
600  attr(i) = part[i];
601  }
602  }
603  else
604  {
605  mesh->PrintWithPartitioning(part, sol_sock, 1);
606  }
607  }
608  attr.Save(sol_sock);
609  sol_sock << "RjlmAb***********";
610  if (mk == 'v')
611  {
612  sol_sock << "e";
613  }
614  else
615  {
616  sol_sock << "\n";
617  }
618  }
619  else
620  {
621  sol_sock << "fem3d_gf_data_keys\n";
622  if (mk == 'b' || mk == 'v' || mk == 'h' || mk == 'k')
623  {
624  mesh->Print(sol_sock);
625  }
626  else
627  {
628  // NURBS meshes do not support PrintWithPartitioning
629  if (mesh->NURBSext)
630  {
631  mesh->Print(sol_sock);
632  for (int i = 0; i < mesh->GetNE(); i++)
633  {
634  attr(i) = part[i];
635  }
636  }
637  else
638  {
639  mesh->PrintWithPartitioning(part, sol_sock);
640  }
641  }
642  attr.Save(sol_sock);
643  sol_sock << "maaA";
644  if (mk == 'v')
645  {
646  sol_sock << "aa";
647  }
648  else
649  {
650  sol_sock << "\n";
651  }
652  }
653  sol_sock << flush;
654  }
655  else
656  {
657  cout << "Unable to connect to "
658  << vishost << ':' << visport << endl;
659  }
660  delete attr_fespace;
661  }
662 
663  if (mk == 'S')
664  {
665  const char mesh_file[] = "mesh-explorer.mesh";
666  ofstream omesh(mesh_file);
667  // Save gzip-ed mesh, requires MFEM_USE_GZSTREAM = YES
668  // ofgzstream omesh(mesh_file, "zwb9");
669  omesh.precision(14);
670  mesh->Print(omesh);
671  cout << "New mesh file: " << mesh_file << endl;
672  }
673  }
674 
675  delete attr_fec;
676  delete mesh;
677 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:222
bool Help() const
Definition: optparser.hpp:121
int Size() const
Logical size of the array.
Definition: array.hpp:109
int * CartesianPartitioning(int nxyz[])
Definition: mesh.cpp:4541
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:161
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:83
const double * GetVertex(int i) const
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:617
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:104
Mesh * read_par_mesh(int np, const char *mesh_prefix)
void JacToPerfJac(int GeomType, const DenseMatrix &J, DenseMatrix &PJ) const
Definition: geom.cpp:428
int * GeneratePartitioning(int nparts, int part_method=1)
Definition: mesh.cpp:4586
void PrintWithPartitioning(int *partitioning, std::ostream &out, int elem_attr=0) const
Definition: mesh.cpp:7538
const Geometry::Type geom
double Det() const
Calculates the determinant of the matrix (for 2x2 or 3x3 matrices)
Definition: densemat.cpp:416
void SetIntPoint(const IntegrationPoint *ip)
Definition: eltrans.hpp:51
RefinedGeometry * Refine(int Geom, int Times, int ETimes=1)
Definition: geom.cpp:566
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
void Transform(void(*f)(const Vector &, Vector &))
Definition: mesh.cpp:8228
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:584
void PrintHelp(std::ostream &out) const
Definition: optparser.cpp:378
void Randomize(int seed=0)
Set random values in the vector.
Definition: vector.cpp:648
const IntegrationPoint & GetCenter(int GeomType)
Definition: geom.hpp:57
int main(int argc, char *argv[])
void DeleteAll()
Delete whole array.
Definition: array.hpp:481
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2185
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:184
Geometry Geometries
Definition: geom.cpp:544
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:225
int dim
Definition: ex3.cpp:47
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:190
void transformation(const Vector &p, Vector &v)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6684
const DenseMatrix & Jacobian()
Definition: eltrans.hpp:64
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3278
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:1035
GeometryRefiner GlobGeometryRefiner
Definition: geom.cpp:1022
IntegrationRule RefPts
Definition: geom.hpp:190
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:293
int Dimension() const
Definition: mesh.hpp:611
double GetElementSize(int i, int type=0)
Definition: mesh.cpp:44
int SpaceDimension() const
Definition: mesh.hpp:612
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:142
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:74
void PrintCharacteristics(Vector *Vh=NULL, Vector *Vk=NULL, std::ostream &out=std::cout)
Definition: mesh.cpp:151
Nodes: x_i = i/(n-1), i=0,...,n-1.
Definition: fe_coll.hpp:32
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:349
int GetElementBaseGeometry(int i=0) const
Definition: mesh.hpp:654
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:144
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:581
void PrintError(std::ostream &out) const
Definition: optparser.cpp:335
virtual const char * Name() const
Definition: fe_coll.hpp:117
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:238
double Max() const
Returns the maximal element of the vector.
Definition: vector.cpp:718
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:3273
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:1650
double kappa
Definition: ex3.cpp:46
Vector data type.
Definition: vector.hpp:36
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:177
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5368
void GetElementColoring(Array< int > &colors, int el0=0)
Definition: mesh.cpp:7463
void GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:139
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:821
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:140
virtual void Print(std::ostream &out=std::cout) const
Definition: mesh.hpp:988
bool Good() const
Definition: optparser.hpp:120