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