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