libMesh
Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
libMesh::TetGenIO Class Reference

This class implements reading and writing meshes in the TetGen format. More...

#include <tetgen_io.h>

Inheritance diagram for libMesh::TetGenIO:
[legend]

Public Member Functions

 TetGenIO (MeshBase &mesh)
 Constructor. More...
 
 TetGenIO (const MeshBase &mesh)
 Constructor. More...
 
virtual void read (const std::string &) libmesh_override
 This method implements reading a mesh from a specified file in TetGen format. More...
 
virtual void write (const std::string &) libmesh_override
 This method implements writing a mesh to a specified ".poly" file. More...
 
virtual void write_equation_systems (const std::string &, const EquationSystems &, const std::set< std::string > *system_names=libmesh_nullptr)
 This method implements writing a mesh with data to a specified file where the data is taken from the EquationSystems object. More...
 
virtual void write_nodal_data (const std::string &, const std::vector< Number > &, const std::vector< std::string > &)
 This method implements writing a mesh with nodal data to a specified file where the nodal data and variable names are provided. More...
 
virtual void write_nodal_data (const std::string &, const NumericVector< Number > &, const std::vector< std::string > &)
 This method should be overridden by "parallel" output formats for writing nodal data. More...
 
unsigned intascii_precision ()
 Return/set the precision to use when writing ASCII files. More...
 

Public Attributes

std::vector< std::vector< Real > > node_attributes
 Data structure to hold node attributes read in from file. More...
 
std::vector< std::vector< Real > > element_attributes
 Data structure to hold element attributes read in from file. More...
 

Protected Member Functions

MeshBasemesh ()
 
void set_n_partitions (unsigned int n_parts)
 Sets the number of partitions in the mesh. More...
 
void skip_comment_lines (std::istream &in, const char comment_start)
 Reads input from in, skipping all the lines that start with the character comment_start. More...
 
const MeshBasemesh () const
 

Protected Attributes

std::vector< bool > elems_of_dimension
 A vector of bools describing what dimension elements have been encountered when reading a mesh. More...
 
const bool _is_parallel_format
 Flag specifying whether this format is parallel-capable. More...
 
const bool _serial_only_needed_on_proc_0
 Flag specifying whether this format can be written by only serializing the mesh to processor zero. More...
 

Private Member Functions

void read_nodes_and_elem (std::istream &node_stream, std::istream &ele_stream)
 Reads a mesh (nodes & elements) from the file provided through node_stream and ele_stream. More...
 
void node_in (std::istream &node_stream)
 Method reads nodes from node_stream and stores them in vector<Node *> nodes in the order they come in. More...
 
void element_in (std::istream &ele_stream)
 Method reads elements and stores them in vector<Elem *> elements in the same order as they come in. More...
 

Private Attributes

std::map< dof_id_type, dof_id_type_assign_nodes
 stores new positions of nodes. More...
 
dof_id_type _num_nodes
 total number of nodes. More...
 
dof_id_type _num_elements
 total number of elements. More...
 

Detailed Description

This class implements reading and writing meshes in the TetGen format.

Format description: cf. TetGen home page.

Author
Benjamin S. Kirk
Date
2004

Definition at line 47 of file tetgen_io.h.

Constructor & Destructor Documentation

libMesh::TetGenIO::TetGenIO ( MeshBase mesh)
explicit

Constructor.

Takes a writable reference to a mesh object. This is the constructor required to read a mesh.

Definition at line 151 of file tetgen_io.h.

151  :
152  MeshInput<MeshBase> (mesh),
153  MeshOutput<MeshBase>(mesh)
154 {
155 }
libMesh::TetGenIO::TetGenIO ( const MeshBase mesh)
explicit

Constructor.

Takes a read-only reference to a mesh object. This is the constructor required to write a mesh.

Definition at line 160 of file tetgen_io.h.

160  :
161  MeshOutput<MeshBase>(mesh)
162 {
163 }

Member Function Documentation

unsigned int& libMesh::MeshOutput< MeshBase >::ascii_precision ( )
inherited

Return/set the precision to use when writing ASCII files.

By default we use numeric_limits<Real>::digits10 + 2, which should be enough to write out to ASCII and get the exact same Real back when reading in.

Referenced by libMesh::TecplotIO::write_ascii(), libMesh::GMVIO::write_ascii_new_impl(), and libMesh::GMVIO::write_ascii_old_impl().

void libMesh::TetGenIO::element_in ( std::istream &  ele_stream)
private

Method reads elements and stores them in vector<Elem *> elements in the same order as they come in.

Within TetGenMeshInterface, element labels are ignored.

Definition at line 170 of file tetgen_io.C.

References _assign_nodes, _num_elements, libMesh::libmesh_assert(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshInput< MT >::mesh(), and n_nodes.

Referenced by read_nodes_and_elem().

171 {
172  // Check input buffer
173  libmesh_assert (ele_stream.good());
174 
175  // Get a reference to the mesh
176  MeshBase & mesh = MeshInput<MeshBase>::mesh();
177 
178  // Read the elements from the ele_stream (*.ele file).
179  unsigned int element_lab=0, n_nodes=0, region_attribute=0;
180 
181  ele_stream >> _num_elements // Read the number of tetrahedrons from the stream.
182  >> n_nodes // Read the number of nodes per tetrahedron from the stream (defaults to 4).
183  >> region_attribute; // Read the number of attributes from stream.
184 
185  // According to the Tetgen docs for .ele files:
186  // http://wias-berlin.de/software/tetgen/1.5/doc/manual/manual006.html#ff_ele
187  // region_attribute can either 0 or 1, and specifies whether, for
188  // each tetrahedron, there is an extra integer specifying which
189  // region it belongs to. Normally, this id matches a value in a
190  // corresponding .poly or .smesh file, but here we simply use it to
191  // set the subdomain_id of the element in question.
192  if (region_attribute > 1)
193  libmesh_error_msg("Invalid region_attribute " << region_attribute << " specified in .ele file.");
194 
195  // Vector that assigns element nodes to their correct position.
196  // TetGen is normally 0-based
197  // (right now this is strictly not necessary since it is the identity map,
198  // but in the future TetGen could change their numbering scheme.)
199  static const unsigned int assign_elm_nodes[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9};
200 
201  for (dof_id_type i=0; i<_num_elements; i++)
202  {
203  libmesh_assert (ele_stream.good());
204 
205  // TetGen only supports Tet4 and Tet10 elements.
206  Elem * elem;
207 
208  if (n_nodes==4)
209  elem = new Tet4;
210 
211  else if (n_nodes==10)
212  elem = new Tet10;
213 
214  else
215  libmesh_error_msg("Elements with " << n_nodes << " nodes are not supported in the LibMesh tetgen module.");
216 
217  elem->set_id(i);
218 
219  mesh.add_elem (elem);
220 
221  libmesh_assert(elem);
222  libmesh_assert_equal_to (elem->n_nodes(), n_nodes);
223 
224  // The first number on the line is the tetrahedron number. We
225  // have previously ignored this, preferring to set our own ids,
226  // but this could be changed to respect the Tetgen numbering if
227  // desired.
228  ele_stream >> element_lab;
229 
230  // Read node labels
231  for (dof_id_type j=0; j<n_nodes; j++)
232  {
233  dof_id_type node_label;
234  ele_stream >> node_label;
235 
236  // Assign node to element
237  elem->set_node(assign_elm_nodes[j]) =
238  mesh.node_ptr(_assign_nodes[node_label]);
239  }
240 
241  // Read the region attribute (if present) and use it to set the subdomain id.
242  if (region_attribute)
243  {
244  unsigned int region;
245  ele_stream >> region;
246 
247  // Make sure that the id we read can be successfully cast to
248  // an integral value of type subdomain_id_type.
249  elem->subdomain_id() = cast_int<subdomain_id_type>(region);
250  }
251  }
252 }
dof_id_type _num_elements
total number of elements.
Definition: tetgen_io.h:143
std::map< dof_id_type, dof_id_type > _assign_nodes
stores new positions of nodes.
Definition: tetgen_io.h:133
libmesh_assert(j)
const dof_id_type n_nodes
Definition: tecplot_io.C:67
uint8_t dof_id_type
Definition: id_types.h:64
MeshBase & libMesh::MeshInput< MeshBase >::mesh ( )
protectedinherited
Returns
The object as a writable reference.

Referenced by libMesh::GMVIO::_read_one_cell(), element_in(), libMesh::UNVIO::elements_in(), libMesh::UNVIO::elements_out(), libMesh::UNVIO::groups_in(), node_in(), libMesh::UNVIO::nodes_in(), libMesh::UNVIO::nodes_out(), libMesh::GMVIO::read(), libMesh::ExodusII_IO::read(), libMesh::CheckpointIO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read_bcs(), libMesh::CheckpointIO::read_connectivity(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::UCDIO::read_implementation(), libMesh::UNVIO::read_implementation(), libMesh::GmshIO::read_mesh(), libMesh::CheckpointIO::read_nodes(), libMesh::CheckpointIO::read_nodesets(), libMesh::CheckpointIO::read_remote_elem(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::OFFIO::read_stream(), libMesh::MatlabIO::read_stream(), libMesh::CheckpointIO::read_subdomain_names(), libMesh::UCDIO::UCDIO(), libMesh::VTKIO::VTKIO(), write(), libMesh::ExodusII_IO::write(), libMesh::CheckpointIO::write(), libMesh::XdrIO::write(), libMesh::GMVIO::write_ascii_new_impl(), libMesh::GMVIO::write_ascii_old_impl(), libMesh::CheckpointIO::write_bcs(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::ExodusII_IO::write_element_data(), libMesh::UCDIO::write_implementation(), libMesh::GmshIO::write_mesh(), libMesh::UCDIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_common(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), libMesh::CheckpointIO::write_nodesets(), libMesh::XdrIO::write_parallel(), libMesh::GmshIO::write_post(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), and libMesh::CheckpointIO::write_subdomain_names().

const MeshBase & libMesh::MeshOutput< MeshBase >::mesh ( ) const
protectedinherited
void libMesh::TetGenIO::node_in ( std::istream &  node_stream)
private

Method reads nodes from node_stream and stores them in vector<Node *> nodes in the order they come in.

The original node labels are being stored in the map _assign_nodes in order to assign the elements to the right nodes later.

Definition at line 112 of file tetgen_io.C.

References _assign_nodes, _num_nodes, libMesh::libmesh_assert(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshInput< MT >::mesh(), node_attributes, and libMesh::Real.

Referenced by read_nodes_and_elem().

113 {
114  // Check input buffer
115  libmesh_assert (node_stream.good());
116 
117  // Get a reference to the mesh
118  MeshBase & mesh = MeshInput<MeshBase>::mesh();
119 
120  unsigned int dimension=0, nAttributes=0, BoundaryMarkers=0;
121 
122  node_stream >> _num_nodes // Read the number of nodes from the stream
123  >> dimension // Read the dimension from the stream
124  >> nAttributes // Read the number of attributes from stream
125  >> BoundaryMarkers; // Read if or not boundary markers are included in *.node (0 or 1)
126 
127  // Read the nodal coordinates from the node_stream (*.node file).
128  unsigned int node_lab=0;
129  Point xyz;
130  Real dummy;
131 
132  // If present, make room for node attributes to be stored.
133  this->node_attributes.resize(nAttributes);
134  for (unsigned i=0; i<nAttributes; ++i)
135  this->node_attributes[i].resize(_num_nodes);
136 
137 
138  for (unsigned int i=0; i<_num_nodes; i++)
139  {
140  // Check input buffer
141  libmesh_assert (node_stream.good());
142 
143  node_stream >> node_lab // node number
144  >> xyz(0) // x-coordinate value
145  >> xyz(1) // y-coordinate value
146  >> xyz(2); // z-coordinate value
147 
148  // Read and store attributes from the stream.
149  for (unsigned int j=0; j<nAttributes; j++)
150  node_stream >> node_attributes[j][i];
151 
152  // Read (and discard) boundary marker if BoundaryMarker=1.
153  // TODO: should we store this somehow?
154  if (BoundaryMarkers == 1)
155  node_stream >> dummy;
156 
157  // Store the new position of the node under its label.
158  //_assign_nodes.insert (std::make_pair(node_lab,i));
159  _assign_nodes[node_lab] = i;
160 
161  // Add this point to the Mesh.
162  mesh.add_point(xyz, i);
163  }
164 }
std::map< dof_id_type, dof_id_type > _assign_nodes
stores new positions of nodes.
Definition: tetgen_io.h:133
libmesh_assert(j)
dof_id_type _num_nodes
total number of nodes.
Definition: tetgen_io.h:138
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::vector< std::vector< Real > > node_attributes
Data structure to hold node attributes read in from file.
Definition: tetgen_io.h:82
void libMesh::TetGenIO::read ( const std::string &  name)
virtual

This method implements reading a mesh from a specified file in TetGen format.

Implements libMesh::MeshInput< MeshBase >.

Definition at line 33 of file tetgen_io.C.

References libMesh::MeshInput< MT >::mesh(), libMesh::Quality::name(), libMesh::out, libMesh::processor_id(), read_nodes_and_elem(), and libMesh::MeshInput< MeshBase >::skip_comment_lines().

Referenced by libMesh::NameBasedIO::read().

34 {
35  // This is a serial-only process for now;
36  // the Mesh should be read on processor 0 and
37  // broadcast later
38  libmesh_assert_equal_to (MeshOutput<MeshBase>::mesh().processor_id(), 0);
39 
40  std::string name_node, name_ele, dummy;
41 
42  // tetgen only works in 3D
43  MeshInput<MeshBase>::mesh().set_mesh_dimension(3);
44 
45 #if LIBMESH_DIM < 3
46  libmesh_error_msg("Cannot open dimension 3 mesh file when configured without 3D support.");
47 #endif
48 
49  // Check name for *.node or *.ele extension.
50  // Set std::istream for node_stream and ele_stream.
51  //
52  if (name.rfind(".node") < name.size())
53  {
54  name_node = name;
55  dummy = name;
56  std::size_t position = dummy.rfind(".node");
57  name_ele = dummy.replace(position, 5, ".ele");
58  }
59  else if (name.rfind(".ele") < name.size())
60  {
61  name_ele = name;
62  dummy = name;
63  std::size_t position = dummy.rfind(".ele");
64  name_node = dummy.replace(position, 4, ".node");
65  }
66  else
67  libmesh_error_msg("ERROR: Unrecognized file name: " << name);
68 
69 
70 
71  // Set the streams from which to read in
72  std::ifstream node_stream (name_node.c_str());
73  std::ifstream ele_stream (name_ele.c_str());
74 
75  if (!node_stream.good() || !ele_stream.good())
76  libmesh_error_msg("Error while opening either " \
77  << name_node \
78  << " or " \
79  << name_ele);
80 
81  libMesh::out<< "TetGenIO found the tetgen files to read " <<std::endl;
82 
83  // Skip the comment lines at the beginning
84  this->skip_comment_lines (node_stream, '#');
85  this->skip_comment_lines (ele_stream, '#');
86 
87  // Read the nodes and elements from the streams
88  this->read_nodes_and_elem (node_stream, ele_stream);
89  libMesh::out<< "TetGenIO read in nodes and elements " <<std::endl;
90 }
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:39
const MeshBase & mesh() const
void skip_comment_lines(std::istream &in, const char comment_start)
Reads input from in, skipping all the lines that start with the character comment_start.
OStreamProxy out
void read_nodes_and_elem(std::istream &node_stream, std::istream &ele_stream)
Reads a mesh (nodes & elements) from the file provided through node_stream and ele_stream.
Definition: tetgen_io.C:94
processor_id_type processor_id()
Definition: libmesh_base.h:96
void libMesh::TetGenIO::read_nodes_and_elem ( std::istream &  node_stream,
std::istream &  ele_stream 
)
private

Reads a mesh (nodes & elements) from the file provided through node_stream and ele_stream.

Definition at line 94 of file tetgen_io.C.

References _assign_nodes, _num_elements, _num_nodes, element_in(), and node_in().

Referenced by read().

96 {
97  _num_nodes = 0;
98  _num_elements = 0;
99 
100  // Read all the datasets.
101  this->node_in (node_stream);
102  this->element_in (ele_stream);
103 
104  // some more clean-up
105  _assign_nodes.clear();
106 }
dof_id_type _num_elements
total number of elements.
Definition: tetgen_io.h:143
std::map< dof_id_type, dof_id_type > _assign_nodes
stores new positions of nodes.
Definition: tetgen_io.h:133
dof_id_type _num_nodes
total number of nodes.
Definition: tetgen_io.h:138
void node_in(std::istream &node_stream)
Method reads nodes from node_stream and stores them in vector<Node *> nodes in the order they come in...
Definition: tetgen_io.C:112
void element_in(std::istream &ele_stream)
Method reads elements and stores them in vector<Elem *> elements in the same order as they come in...
Definition: tetgen_io.C:170
void libMesh::MeshInput< MeshBase >::set_n_partitions ( unsigned int  n_parts)
protectedinherited

Sets the number of partitions in the mesh.

Typically this gets done by the partitioner, but some parallel file formats begin "pre-partitioned".

Definition at line 91 of file mesh_input.h.

References libMesh::MeshInput< MT >::mesh().

Referenced by libMesh::Nemesis_IO::read(), and libMesh::XdrIO::read_header().

91 { this->mesh().set_n_partitions() = n_parts; }
unsigned int & set_n_partitions()
Definition: mesh_base.h:1329
void libMesh::MeshInput< MeshBase >::skip_comment_lines ( std::istream &  in,
const char  comment_start 
)
protectedinherited

Reads input from in, skipping all the lines that start with the character comment_start.

Referenced by read(), and libMesh::UCDIO::read_implementation().

void libMesh::TetGenIO::write ( const std::string &  fname)
virtual

This method implements writing a mesh to a specified ".poly" file.

".poly" files defines so called Piecewise Linear Complex (PLC).

Implements libMesh::MeshOutput< MeshBase >.

Definition at line 259 of file tetgen_io.C.

References libMesh::MeshBase::active_element_ptr_range(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), libMesh::MeshBase::n_elem(), libMesh::MeshBase::n_nodes(), and libMesh::MeshBase::point().

Referenced by libMesh::NameBasedIO::write().

260 {
261  // libmesh_assert three dimensions (should be extended later)
262  libmesh_assert_equal_to (MeshOutput<MeshBase>::mesh().mesh_dimension(), 3);
263 
264  if (!(fname.rfind(".poly") < fname.size()))
265  libmesh_error_msg("ERROR: Unrecognized file name: " << fname);
266 
267  // Open the output file stream
268  std::ofstream out_stream (fname.c_str());
269 
270  // Make sure it opened correctly
271  if (!out_stream.good())
272  libmesh_file_error(fname.c_str());
273 
274  // Get a reference to the mesh
275  const MeshBase & mesh = MeshOutput<MeshBase>::mesh();
276 
277  // Begin interfacing with the .poly file
278  {
279  // header:
280  out_stream << "# poly file output generated by libmesh\n"
281  << mesh.n_nodes() << " 3 0 0\n";
282 
283  // write the nodes:
284  for (dof_id_type v=0; v<mesh.n_nodes(); v++)
285  out_stream << v << " "
286  << mesh.point(v)(0) << " "
287  << mesh.point(v)(1) << " "
288  << mesh.point(v)(2) << "\n";
289  }
290 
291  {
292  // write the connectivity:
293  out_stream << "# Facets:\n"
294  << mesh.n_elem() << " 0\n";
295 
296  for (const auto & elem : mesh.active_element_ptr_range())
297  out_stream << "1\n3 " // no. of facet polygons
298  // << elem->n_nodes() << " "
299  << elem->node_id(0) << " "
300  << elem->node_id(1) << " "
301  << elem->node_id(2) << "\n";
302  }
303 
304  // end of the file
305  out_stream << "0\n"; // no holes output!
306  out_stream << "\n\n# end of file\n";
307 }
const MeshBase & mesh() const
uint8_t dof_id_type
Definition: id_types.h:64
virtual void libMesh::MeshOutput< MeshBase >::write_equation_systems ( const std::string &  ,
const EquationSystems ,
const std::set< std::string > *  system_names = libmesh_nullptr 
)
virtualinherited

This method implements writing a mesh with data to a specified file where the data is taken from the EquationSystems object.

Reimplemented in libMesh::NameBasedIO.

Referenced by libMesh::Nemesis_IO::write_timestep(), and libMesh::ExodusII_IO::write_timestep().

virtual void libMesh::MeshOutput< MeshBase >::write_nodal_data ( const std::string &  ,
const std::vector< Number > &  ,
const std::vector< std::string > &   
)
virtualinherited

This method implements writing a mesh with nodal data to a specified file where the nodal data and variable names are provided.

Reimplemented in libMesh::ExodusII_IO, libMesh::NameBasedIO, libMesh::GmshIO, libMesh::Nemesis_IO, libMesh::VTKIO, libMesh::UCDIO, libMesh::GMVIO, libMesh::MEDITIO, libMesh::GnuPlotIO, and libMesh::TecplotIO.

Definition at line 96 of file mesh_output.h.

References libMesh::MeshOutput< MT >::ascii_precision(), libMesh::MeshOutput< MT >::mesh(), and libMesh::MeshOutput< MT >::write_nodal_data().

99  { libmesh_not_implemented(); }
virtual void libMesh::MeshOutput< MeshBase >::write_nodal_data ( const std::string &  ,
const NumericVector< Number > &  ,
const std::vector< std::string > &   
)
virtualinherited

This method should be overridden by "parallel" output formats for writing nodal data.

Instead of getting a localized copy of the nodal solution vector, it is passed a NumericVector of type=PARALLEL which is in node-major order i.e. (u0,v0,w0, u1,v1,w1, u2,v2,w2, u3,v3,w3, ...) and contains n_nodes*n_vars total entries. Then, it is up to the individual I/O class to extract the required solution values from this vector and write them in parallel.

If not implemented, localizes the parallel vector into a std::vector and calls the other version of this function.

Reimplemented in libMesh::Nemesis_IO.

Member Data Documentation

std::map<dof_id_type,dof_id_type> libMesh::TetGenIO::_assign_nodes
private

stores new positions of nodes.

Used when reading.

Definition at line 133 of file tetgen_io.h.

Referenced by element_in(), node_in(), and read_nodes_and_elem().

const bool libMesh::MeshOutput< MeshBase >::_is_parallel_format
protectedinherited

Flag specifying whether this format is parallel-capable.

If this is false (default) I/O is only permitted when the mesh has been serialized.

Definition at line 141 of file mesh_output.h.

Referenced by libMesh::FroIO::write(), libMesh::PostscriptIO::write(), and libMesh::EnsightIO::write().

dof_id_type libMesh::TetGenIO::_num_elements
private

total number of elements.

Primarily used when reading.

Definition at line 143 of file tetgen_io.h.

Referenced by element_in(), and read_nodes_and_elem().

dof_id_type libMesh::TetGenIO::_num_nodes
private

total number of nodes.

Primarily used when reading.

Definition at line 138 of file tetgen_io.h.

Referenced by node_in(), and read_nodes_and_elem().

const bool libMesh::MeshOutput< MeshBase >::_serial_only_needed_on_proc_0
protectedinherited

Flag specifying whether this format can be written by only serializing the mesh to processor zero.

If this is false (default) the mesh will be serialized to all processors

Definition at line 150 of file mesh_output.h.

std::vector<std::vector<Real> > libMesh::TetGenIO::element_attributes

Data structure to hold element attributes read in from file.

Note
This vector is no longer filled or used for anything. If region attributes are present in the .ele file, they are used to set the subdomain ids of the elements as they are created.
Deprecated:
This member, since it was originally a part of the public interface, remains for now, but will be removed some time in the near future.

Definition at line 95 of file tetgen_io.h.

std::vector<bool> libMesh::MeshInput< MeshBase >::elems_of_dimension
protectedinherited
std::vector<std::vector<Real> > libMesh::TetGenIO::node_attributes

Data structure to hold node attributes read in from file.

What you do with these is up to you!

Definition at line 82 of file tetgen_io.h.

Referenced by node_in().


The documentation for this class was generated from the following files: