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