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