vtk_io.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 // C++ includes
20 #include <fstream>
21 
22 // Local includes
23 #include "libmesh/libmesh_config.h"
24 #include "libmesh/vtk_io.h"
25 #include "libmesh/mesh_base.h"
27 #include "libmesh/numeric_vector.h"
28 #include "libmesh/system.h"
29 #include "libmesh/node.h"
30 #include "libmesh/elem.h"
32 
33 #ifdef LIBMESH_HAVE_VTK
34 
35 // I get a lot of "warning: extra ';' inside a class [-Wextra-semi]" from clang
36 // on VTK header files.
38 
39 #include "vtkXMLUnstructuredGridReader.h"
40 #include "vtkXMLUnstructuredGridWriter.h"
41 #include "vtkXMLPUnstructuredGridWriter.h"
42 #include "vtkUnstructuredGrid.h"
43 #include "vtkIntArray.h"
44 #include "vtkCellArray.h"
45 #include "vtkCellData.h"
46 #include "vtkConfigure.h"
47 #include "vtkDoubleArray.h"
48 #include "vtkGenericCell.h"
49 #include "vtkPointData.h"
50 #include "vtkPoints.h"
51 #include "vtkSmartPointer.h"
52 
54 
55 // A convenient macro for comparing VTK versions. Returns 1 if the
56 // current VTK version is < major.minor.subminor and zero otherwise.
57 //
58 // It relies on the VTK version numbers detected during configure. Note that if
59 // LIBMESH_HAVE_VTK is not defined, none of the LIBMESH_DETECTED_VTK_VERSION_* variables will
60 // be defined either.
61 #define VTK_VERSION_LESS_THAN(major,minor,subminor) \
62  ((LIBMESH_DETECTED_VTK_VERSION_MAJOR < (major) || \
63  (LIBMESH_DETECTED_VTK_VERSION_MAJOR == (major) && (LIBMESH_DETECTED_VTK_VERSION_MINOR < (minor) || \
64  (LIBMESH_DETECTED_VTK_VERSION_MINOR == (minor) && \
65  LIBMESH_DETECTED_VTK_VERSION_SUBMINOR < (subminor))))) ? 1 : 0)
66 
67 #endif // LIBMESH_HAVE_VTK
68 
69 
70 
71 namespace libMesh
72 {
73 
74 // Constructor for reading
76  MeshInput<MeshBase> (mesh, /*is_parallel_format=*/true),
77  MeshOutput<MeshBase>(mesh, /*is_parallel_format=*/true)
78 #ifdef LIBMESH_HAVE_VTK
79  ,_compress(false)
80 #endif
81 {
82 }
83 
84 
85 
86 // Constructor for writing
87 VTKIO::VTKIO (const MeshBase & mesh) :
88  MeshOutput<MeshBase>(mesh, /*is_parallel_format=*/true)
89 #ifdef LIBMESH_HAVE_VTK
90  ,_compress(false)
91 #endif
92 {
93 }
94 
95 
96 
97 // Output the mesh without solutions to a .pvtu file
98 void VTKIO::write (const std::string & name)
99 {
100  std::vector<Number> soln;
101  std::vector<std::string> names;
102  this->write_nodal_data(name, soln, names);
103 }
104 
105 
106 
107 // The rest of the file is wrapped in ifdef LIBMESH_HAVE_VTK except for
108 // a couple of "stub" functions at the bottom.
109 #ifdef LIBMESH_HAVE_VTK
110 
111 // Initialize the static _element_maps struct.
112 VTKIO::ElementMaps VTKIO::_element_maps = VTKIO::build_element_maps();
113 
114 // Static function which constructs the ElementMaps object.
116 {
117  // Object to be filled up
118  ElementMaps em;
119 
120  em.associate(EDGE2, VTK_LINE);
121  em.associate(EDGE3, VTK_QUADRATIC_EDGE);
122  em.associate(TRI3, VTK_TRIANGLE);
123  em.associate(TRI6, VTK_QUADRATIC_TRIANGLE);
124  em.associate(QUAD4, VTK_QUAD);
125  em.associate(QUAD8, VTK_QUADRATIC_QUAD);
126  em.associate(TET4, VTK_TETRA);
127  em.associate(TET10, VTK_QUADRATIC_TETRA);
128  em.associate(HEX8, VTK_HEXAHEDRON);
129  em.associate(HEX20, VTK_QUADRATIC_HEXAHEDRON);
130  em.associate(HEX27, VTK_TRIQUADRATIC_HEXAHEDRON);
131  em.associate(PRISM6, VTK_WEDGE);
132  em.associate(PRISM15, VTK_QUADRATIC_WEDGE);
133  em.associate(PRISM18, VTK_BIQUADRATIC_QUADRATIC_WEDGE);
134  em.associate(PYRAMID5, VTK_PYRAMID);
135 
136  // VTK_BIQUADRATIC_QUAD has been around since VTK 5.0
137 #if VTK_MAJOR_VERSION > 5 || (VTK_MAJOR_VERSION == 5 && VTK_MINOR_VERSION > 0)
138  em.associate(QUAD9, VTK_BIQUADRATIC_QUAD);
139 #endif
140 
141  // TRI3SUBDIVISION is for writing only
142  em.writing_map[TRI3SUBDIVISION] = VTK_TRIANGLE;
143 
144  return em;
145 }
146 
147 
148 
149 void VTKIO::read (const std::string & name)
150 {
151  // This is a serial-only process for now;
152  // the Mesh should be read on processor 0 and
153  // broadcast later
154  libmesh_assert_equal_to (MeshOutput<MeshBase>::mesh().processor_id(), 0);
155 
156  // Keep track of what kinds of elements this file contains
157  elems_of_dimension.clear();
158  elems_of_dimension.resize(4, false);
159 
160  // Use a typedef, because these names are just crazy
161  typedef vtkSmartPointer<vtkXMLUnstructuredGridReader> MyReader;
162  MyReader reader = MyReader::New();
163 
164  // Pass the filename along to the reader
165  reader->SetFileName(name.c_str());
166 
167  // Force reading
168  reader->Update();
169 
170  // read in the grid
171  _vtk_grid = reader->GetOutput();
172 
173  // Get a reference to the mesh
175 
176  // Clear out any pre-existing data from the Mesh
177  mesh.clear();
178 
179  // Get the number of points from the _vtk_grid object
180  const unsigned int vtk_num_points = static_cast<unsigned int>(_vtk_grid->GetNumberOfPoints());
181 
182  // always numbered nicely so we can loop like this
183  for (unsigned int i=0; i<vtk_num_points; ++i)
184  {
185  // add to the id map
186  // and add the actual point
187  double pnt[3];
188  _vtk_grid->GetPoint(static_cast<vtkIdType>(i), pnt);
189  Point xyz(pnt[0], pnt[1], pnt[2]);
190  mesh.add_point(xyz, i);
191  }
192 
193  // Get the number of cells from the _vtk_grid object
194  const unsigned int vtk_num_cells = static_cast<unsigned int>(_vtk_grid->GetNumberOfCells());
195 
196  vtkSmartPointer<vtkGenericCell> cell = vtkSmartPointer<vtkGenericCell>::New();
197  for (unsigned int i=0; i<vtk_num_cells; ++i)
198  {
199  _vtk_grid->GetCell(i, cell);
200 
201  // Get the libMesh element type corresponding to this VTK element type.
202  ElemType libmesh_elem_type = _element_maps.find(cell->GetCellType());
203  Elem * elem = Elem::build(libmesh_elem_type).release();
204 
205  // get the straightforward numbering from the VTK cells
206  for (unsigned int j=0; j<elem->n_nodes(); ++j)
207  elem->set_node(j) =
208  mesh.node_ptr(cast_int<dof_id_type>(cell->GetPointId(j)));
209 
210  // then get the connectivity
211  std::vector<dof_id_type> conn;
212  elem->connectivity(0, VTK, conn);
213 
214  // then reshuffle the nodes according to the connectivity, this
215  // two-time-assign would evade the definition of the vtk_mapping
216  for (unsigned int j=0,
217  n_conn = cast_int<unsigned int>(conn.size());
218  j != n_conn; ++j)
219  elem->set_node(j) = mesh.node_ptr(conn[j]);
220 
221  elem->set_id(i);
222 
223  elems_of_dimension[elem->dim()] = true;
224 
225  mesh.add_elem(elem);
226  } // end loop over VTK cells
227 
228  // Set the mesh dimension to the largest encountered for an element
229  for (unsigned char i=0; i!=4; ++i)
230  if (elems_of_dimension[i])
232 
233 #if LIBMESH_DIM < 3
234  if (mesh.mesh_dimension() > LIBMESH_DIM)
235  libmesh_error_msg("Cannot open dimension " \
236  << mesh.mesh_dimension() \
237  << " mesh file when configured without " \
238  << mesh.mesh_dimension() \
239  << "D support.");
240 #endif // LIBMESH_DIM < 3
241 }
242 
243 
244 
245 void VTKIO::write_nodal_data (const std::string & fname,
246  const std::vector<Number> & soln,
247  const std::vector<std::string> & names)
248 {
250 
251  // Warn that the .pvtu file extension should be used. Paraview
252  // recognizes this, and it works in both serial and parallel. Only
253  // warn about this once.
254  if (fname.substr(fname.rfind("."), fname.size()) != ".pvtu")
255  libmesh_do_once(libMesh::err << "The .pvtu extension should be used when writing VTK files in libMesh.");
256 
257  // If there are variable names being written, the solution vector
258  // should not be empty, it should have been broadcast to all
259  // processors by the MeshOutput base class, since VTK is a parallel
260  // format. Verify this before going further.
261  if (!names.empty() && soln.empty())
262  libmesh_error_msg("Empty soln vector in VTKIO::write_nodal_data().");
263 
264  // we only use Unstructured grids
265  _vtk_grid = vtkSmartPointer<vtkUnstructuredGrid>::New();
266  vtkSmartPointer<vtkXMLPUnstructuredGridWriter> writer = vtkSmartPointer<vtkXMLPUnstructuredGridWriter>::New();
267 
268  // add nodes to the grid and update _local_node_map
269  _local_node_map.clear();
270  this->nodes_to_vtk();
271 
272  // add cells to the grid
273  this->cells_to_vtk();
274 
275  // add nodal solutions to the grid, if solutions are given
276  if (names.size() > 0)
277  {
278  std::size_t num_vars = names.size();
279  dof_id_type num_nodes = mesh.n_nodes();
280 
281  for (std::size_t variable=0; variable<num_vars; ++variable)
282  {
283  vtkSmartPointer<vtkDoubleArray> data = vtkSmartPointer<vtkDoubleArray>::New();
284  data->SetName(names[variable].c_str());
285 
286  // number of local and ghost nodes
287  data->SetNumberOfValues(_local_node_map.size());
288 
289  // loop over all nodes and get the solution for the current
290  // variable, if the node is in the current partition
291  for (dof_id_type k=0; k<num_nodes; ++k)
292  {
293  std::map<dof_id_type, dof_id_type>::iterator local_node_it = _local_node_map.find(k);
294  if (local_node_it == _local_node_map.end())
295  continue; // not a local node
296 
297 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
298  libmesh_do_once (libMesh::err << "Only writing the real part for complex numbers!\n"
299  << "if you need this support contact " << LIBMESH_PACKAGE_BUGREPORT
300  << std::endl);
301  data->SetValue(local_node_it->second, soln[k*num_vars + variable].real());
302 #else
303  data->SetValue(local_node_it->second, soln[k*num_vars + variable]);
304 #endif
305  }
306  _vtk_grid->GetPointData()->AddArray(data);
307  }
308  }
309 
310  // Tell the writer how many partitions exist and on which processor
311  // we are currently
312  writer->SetNumberOfPieces(MeshOutput<MeshBase>::mesh().n_processors());
313  writer->SetStartPiece(MeshOutput<MeshBase>::mesh().processor_id());
314  writer->SetEndPiece(MeshOutput<MeshBase>::mesh().processor_id());
315 
316  // partitions overlap by one node
317  // FIXME: According to this document
318  // http://paraview.org/Wiki/images/5/51/SC07_tut107_ParaView_Handouts.pdf
319  // the ghosts are cells rather than nodes.
320  writer->SetGhostLevel(1);
321 
322  // VTK 6 replaces SetInput() with SetInputData(). See
323  // http://www.vtk.org/Wiki/VTK/VTK_6_Migration/Replacement_of_SetInput
324  // for the full explanation.
325 #if VTK_VERSION_LESS_THAN(6,0,0)
326  writer->SetInput(_vtk_grid);
327 #else
328  writer->SetInputData(_vtk_grid);
329 #endif
330 
331  writer->SetFileName(fname.c_str());
332  writer->SetDataModeToAscii();
333 
334  // compress the output, if desired (switches also to binary)
335  if (this->_compress)
336  {
337 #if !VTK_VERSION_LESS_THAN(5,6,0)
338  writer->SetCompressorTypeToZLib();
339 #else
340  libmesh_do_once(libMesh::err << "Compression not implemented with old VTK libs!" << std::endl;);
341 #endif
342  }
343 
344  writer->Write();
345 
346 }
347 
348 
349 
350 vtkUnstructuredGrid * VTKIO::get_vtk_grid()
351 {
352  return _vtk_grid;
353 }
354 
355 
356 
358 {
359  this->_compress = b;
360 }
361 
362 
363 
365 {
367 
368  // containers for points and coordinates of points
369  vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
370  vtkSmartPointer<vtkDoubleArray> pcoords = vtkSmartPointer<vtkDoubleArray>::New();
371  // if this grid is to be used in VTK then the dimension of the points should be 3
372  pcoords->SetNumberOfComponents(LIBMESH_DIM);
373  pcoords->Allocate(3*mesh.n_local_nodes());
374  points->SetNumberOfPoints(mesh.n_local_nodes()); // it seems that it needs this to prevent a segfault
375 
376  unsigned int local_node_counter = 0;
377 
378  for (const auto & node_ptr : mesh.local_node_ptr_range())
379  {
380  const Node & node = *node_ptr;
381 
382  double pnt[3] = {0, 0, 0};
383  for (unsigned int i=0; i<LIBMESH_DIM; ++i)
384  pnt[i] = node(i);
385 
386  // Fill mapping between global and local node numbers
387  _local_node_map[node.id()] = local_node_counter;
388 
389  // add point
390 #if VTK_VERSION_LESS_THAN(7,1,0)
391  pcoords->InsertNextTupleValue(pnt);
392 #else
393  pcoords->InsertNextTuple(pnt);
394 #endif
395  ++local_node_counter;
396  }
397 
398  // add coordinates to points
399  points->SetData(pcoords);
400 
401  // add points to grid
402  _vtk_grid->SetPoints(points);
403 }
404 
405 
406 
408 {
410 
411  vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
412  vtkSmartPointer<vtkIdList> pts = vtkSmartPointer<vtkIdList>::New();
413 
414  std::vector<int> types(mesh.n_active_local_elem());
415 
416  vtkSmartPointer<vtkIntArray> elem_id = vtkSmartPointer<vtkIntArray>::New();
417  elem_id->SetName("libmesh_elem_id");
418  elem_id->SetNumberOfComponents(1);
419 
420  vtkSmartPointer<vtkIntArray> subdomain_id = vtkSmartPointer<vtkIntArray>::New();
421  subdomain_id->SetName("subdomain_id");
422  subdomain_id->SetNumberOfComponents(1);
423 
424  vtkSmartPointer<vtkIntArray> elem_proc_id = vtkSmartPointer<vtkIntArray>::New();
425  elem_proc_id->SetName("processor_id");
426  elem_proc_id->SetNumberOfComponents(1);
427 
428  unsigned active_element_counter = 0;
429  for (const auto & elem : mesh.active_local_element_ptr_range())
430  {
431  pts->SetNumberOfIds(elem->n_nodes());
432 
433  // get the connectivity for this element
434  std::vector<dof_id_type> conn;
435  elem->connectivity(0, VTK, conn);
436 
437  for (unsigned int i=0,
438  n_conn = cast_int<unsigned int>(conn.size());
439  i != n_conn; ++i)
440  {
441  // If the node ID is not found in the _local_node_map, we'll
442  // add it to the _vtk_grid. NOTE[JWP]: none of the examples
443  // I have actually enters this section of code...
444  if (_local_node_map.find(conn[i]) == _local_node_map.end())
445  {
446  dof_id_type global_node_id = elem->node_id(i);
447 
448  const Point & the_node = mesh.point(global_node_id);
449 
450  // InsertNextPoint accepts either a double or float array of length 3.
451  double pt[3] = {0., 0., 0.};
452  for (unsigned int d=0; d<LIBMESH_DIM; ++d)
453  pt[d] = the_node(d);
454 
455  // Insert the point into the _vtk_grid
456  vtkIdType local = _vtk_grid->GetPoints()->InsertNextPoint(pt);
457 
458  // Update the _local_node_map with the ID returned by VTK
459  _local_node_map[global_node_id] =
460  cast_int<dof_id_type>(local);
461  }
462 
463  // Otherwise, the node ID was found in the _local_node_map, so
464  // insert it into the vtkIdList.
465  pts->InsertId(i, _local_node_map[conn[i]]);
466  }
467 
468  vtkIdType vtkcellid = cells->InsertNextCell(pts);
469  types[active_element_counter] = cast_int<int>(_element_maps.find(elem->type()));
470 
471  elem_id->InsertTuple1(vtkcellid, elem->id());
472  subdomain_id->InsertTuple1(vtkcellid, elem->subdomain_id());
473  elem_proc_id->InsertTuple1(vtkcellid, elem->processor_id());
474  ++active_element_counter;
475  } // end loop over active elements
476 
477  _vtk_grid->SetCells(types.data(), cells);
478  _vtk_grid->GetCellData()->AddArray(elem_id);
479  _vtk_grid->GetCellData()->AddArray(subdomain_id);
480  _vtk_grid->GetCellData()->AddArray(elem_proc_id);
481 }
482 
483 
484 
492 // void VTKIO::system_vectors_to_vtk(const EquationSystems & es,
493 // vtkUnstructuredGrid *& grid)
494 // {
495 // if (MeshOutput<MeshBase>::mesh().processor_id() == 0)
496 // {
497 // std::map<std::string, std::vector<Number>> vecs;
498 // for (unsigned int i=0; i<es.n_systems(); ++i)
499 // {
500 // const System & sys = es.get_system(i);
501 // System::const_vectors_iterator v_end = sys.vectors_end();
502 // System::const_vectors_iterator it = sys.vectors_begin();
503 // for (; it!= v_end; ++it)
504 // {
505 // // for all vectors on this system
506 // std::vector<Number> values;
507 // // libMesh::out<<"it "<<it->first<<std::endl;
508 //
509 // it->second->localize_to_one(values, 0);
510 // // libMesh::out<<"finish localize"<<std::endl;
511 // vecs[it->first] = values;
512 // }
513 // }
514 //
515 // std::map<std::string, std::vector<Number>>::iterator it = vecs.begin();
516 //
517 // for (; it!=vecs.end(); ++it)
518 // {
519 // vtkSmartPointer<vtkDoubleArray> data = vtkSmartPointer<vtkDoubleArray>::New();
520 // data->SetName(it->first.c_str());
521 // libmesh_assert_equal_to (it->second.size(), es.get_mesh().n_nodes());
522 // data->SetNumberOfValues(it->second.size());
523 //
524 // for (std::size_t i=0; i<it->second.size(); ++i)
525 // {
526 // #ifdef LIBMESH_USE_COMPLEX_NUMBERS
527 // libmesh_do_once (libMesh::err << "Only writing the real part for complex numbers!\n"
528 // << "if you need this support contact " << LIBMESH_PACKAGE_BUGREPORT
529 // << std::endl);
530 // data->SetValue(i, it->second[i].real());
531 // #else
532 // data->SetValue(i, it->second[i]);
533 // #endif
534 //
535 // }
536 // grid->GetPointData()->AddArray(data);
537 // }
538 // }
539 // }
540 
541 
542 
543 #else // !LIBMESH_HAVE_VTK
544 
545 void VTKIO::read (const std::string & name)
546 {
547  libmesh_error_msg("Cannot read VTK file: " << name \
548  << "\nYou must have VTK installed and correctly configured to read VTK meshes.");
549 }
550 
551 
552 
553 void VTKIO::write_nodal_data (const std::string & fname,
554  const std::vector<Number> &,
555  const std::vector<std::string> &)
556 {
557  libmesh_error_msg("Cannot write VTK file: " << fname \
558  << "\nYou must have VTK installed and correctly configured to read VTK meshes.");
559 }
560 
561 
562 #endif // LIBMESH_HAVE_VTK
563 
564 
565 
566 } // namespace libMesh
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
void set_compression(bool b)
Definition: vtk_io.C:357
const MT & mesh() const
Definition: mesh_output.h:234
std::map< dof_id_type, dof_id_type > _local_node_map
Definition: vtk_io.h:161
virtual Node *& set_node(const unsigned int i)
Definition: elem.h:2024
A geometric point in (x,y,z) space associated with a DOF.
Definition: node.h:52
virtual SimpleRange< node_iterator > local_node_ptr_range()=0
std::vector< bool > elems_of_dimension
Definition: mesh_input.h:97
vtkIdType find(ElemType libmesh_type)
Definition: vtk_io.h:178
vtkSmartPointer< vtkUnstructuredGrid > _vtk_grid
Definition: vtk_io.h:151
The base class for all geometric element types.
Definition: elem.h:100
MeshBase & mesh
bool _compress
Definition: vtk_io.h:156
virtual void write(const std::string &) override
dof_id_type n_local_nodes() const
Definition: mesh_base.h:286
static ElementMaps _element_maps
Definition: vtk_io.h:207
virtual Node * add_point(const Point &p, const dof_id_type id=DofObject::invalid_id, const processor_id_type proc_id=DofObject::invalid_processor_id)=0
Base class for Mesh.
Definition: mesh_base.h:77
dof_id_type & set_id()
Definition: dof_object.h:664
virtual void read(const std::string &) override
Definition: vtk_io.C:149
void nodes_to_vtk()
Definition: vtk_io.C:364
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
dof_id_type n_active_local_elem() const
Definition: mesh_base.h:403
dof_id_type id() const
Definition: dof_object.h:655
virtual unsigned int n_nodes() const =0
virtual Elem * add_elem(Elem *e)=0
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
Definition: elem.C:245
virtual void write_nodal_data(const std::string &, const std::vector< Number > &, const std::vector< std::string > &) override
Definition: vtk_io.C:245
static ElementMaps build_element_maps()
Definition: vtk_io.C:115
vtkUnstructuredGrid * get_vtk_grid()
Definition: vtk_io.C:350
OStreamProxy err(std::cerr)
void set_mesh_dimension(unsigned char d)
Definition: mesh_base.h:213
std::map< ElemType, vtkIdType > writing_map
Definition: vtk_io.h:199
virtual void clear()
Definition: mesh_base.C:260
VTKIO(MeshBase &mesh)
Definition: vtk_io.C:75
void associate(ElemType libmesh_type, vtkIdType vtk_type)
Definition: vtk_io.h:171
virtual unsigned short dim() const =0
void cells_to_vtk()
Definition: vtk_io.C:407
unsigned int mesh_dimension() const
Definition: mesh_base.C:126
virtual void connectivity(const unsigned int sc, const IOPackage iop, std::vector< dof_id_type > &conn) const =0
virtual const Point & point(const dof_id_type i) const =0
IterBase * data
virtual const Node * node_ptr(const dof_id_type i) const =0
A geometric point in (x,y,z) space.
Definition: point.h:38
virtual dof_id_type n_nodes() const =0
uint8_t dof_id_type
Definition: id_types.h:64