$extrastylesheet
#include <tecplot_io.h>

Public Member Functions | |
| TecplotIO (const MeshBase &, const bool binary=false, const double time=0., const int strand_offset=0) | |
| virtual void | write (const std::string &) |
| virtual void | write_nodal_data (const std::string &, const std::vector< Number > &, const std::vector< std::string > &) |
| bool & | binary () |
| double & | time () |
| int & | strand_offset () |
| std::string & | zone_title () |
| bool & | ascii_append () |
| virtual void | write_equation_systems (const std::string &, const EquationSystems &, const std::set< std::string > *system_names=NULL) |
| unsigned int & | ascii_precision () |
Protected Member Functions | |
| const MeshBase & | mesh () const |
Protected Attributes | |
| const bool | _is_parallel_format |
Private Member Functions | |
| void | write_ascii (const std::string &, const std::vector< Number > *=NULL, const std::vector< std::string > *=NULL) |
| void | write_binary (const std::string &, const std::vector< Number > *=NULL, const std::vector< std::string > *=NULL) |
| unsigned | elem_dimension () |
Private Attributes | |
| bool | _binary |
| double | _time |
| int | _strand_offset |
| std::string | _zone_title |
| bool | _ascii_append |
| std::set< subdomain_id_type > | _subdomain_ids |
This class implements writing meshes in the Tecplot format.
Definition at line 47 of file tecplot_io.h.
| libMesh::TecplotIO::TecplotIO | ( | const MeshBase & | mesh_in, |
| const bool | binary = false, |
||
| const double | time = 0., |
||
| const int | strand_offset = 0 |
||
| ) | [explicit] |
Constructor. Takes a reference to a constant mesh object. This constructor will only allow us to write the mesh. The optional parameter binary can be used to switch between ASCII (false, the default) or binary (true) output files.
Definition at line 120 of file tecplot_io.C.
References _subdomain_ids, and libMesh::MeshBase::subdomain_ids().
: MeshOutput<MeshBase> (mesh_in), _binary (binary_in), _time (time_in), _strand_offset (strand_offset_in), _zone_title ("zone"), _ascii_append(false) { // Gather a list of subdomain ids in the mesh. // We must do this now, while we have every // processor's attention // (some of the write methods only execute on processor 0). mesh_in.subdomain_ids (_subdomain_ids); }
| bool & libMesh::TecplotIO::ascii_append | ( | ) |
Set to true to write multiple solutions to a single file (ASCII only). Tecplot will read multiple zones in a single file, but currently you have to repeat the mesh information each time.
Definition at line 167 of file tecplot_io.C.
References _ascii_append.
{
return _ascii_append;
}
| 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 write_ascii(), libMesh::GMVIO::write_ascii_new_impl(), and libMesh::GMVIO::write_ascii_old_impl().
| bool & libMesh::TecplotIO::binary | ( | ) |
Flag indicating whether or not to write a binary file (if the tecio.a library was found by configure).
Definition at line 140 of file tecplot_io.C.
References _binary.
Referenced by write(), and write_nodal_data().
{
return _binary;
}
| unsigned libMesh::TecplotIO::elem_dimension | ( | ) | [private] |
Determines the logical spatial dimension of the elements in the Mesh. Ex: A 1D edge element living in 3D is a logically one-dimensional element as far as Tecplot is concerned. Throws an error if mixed-dimension element types are found, since I'm not sure how to handle that case currently.
Definition at line 205 of file tecplot_io.C.
References libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), end, and libMesh::MeshOutput< MeshBase >::mesh().
Referenced by write_ascii(), and write_binary().
{
// Get a constant reference to the mesh.
const MeshBase& the_mesh = MeshOutput<MeshBase>::mesh();
std::vector<unsigned> elem_dims(3);
// Loop over all the elements and mark the proper dimension entry in
// the elem_dims vector.
MeshBase::const_element_iterator it = the_mesh.active_elements_begin();
const MeshBase::const_element_iterator end = the_mesh.active_elements_end();
for ( ; it != end; ++it)
elem_dims[(*it)->dim() - 1] = 1;
// Detect and disallow (for now) the writing of mixed dimension meshes.
if (std::count(elem_dims.begin(), elem_dims.end(), 1) > 1)
libmesh_error_msg("Error, cannot write Mesh with mixed element dimensions to Tecplot file!");
if (elem_dims[0])
return 1;
else if (elem_dims[1])
return 2;
else if (elem_dims[2])
return 3;
else
libmesh_error_msg("No 1, 2, or 3D elements detected!");
}
| const MeshBase & libMesh::MeshOutput< MeshBase >::mesh | ( | ) | const [protected, inherited] |
Returns the object as a read-only reference.
Referenced by elem_dimension(), libMesh::FroIO::write(), libMesh::DivaIO::write(), write(), libMesh::PostscriptIO::write(), libMesh::MEDITIO::write(), libMesh::EnsightIO::write(), libMesh::MEDITIO::write_ascii(), write_ascii(), write_binary(), libMesh::EnsightIO::write_geometry_ascii(), write_nodal_data(), libMesh::MEDITIO::write_nodal_data(), libMesh::EnsightIO::write_scalar_ascii(), libMesh::GnuPlotIO::write_solution(), libMesh::DivaIO::write_stream(), and libMesh::EnsightIO::write_vector_ascii().
| int & libMesh::TecplotIO::strand_offset | ( | ) |
Strand offset for this file. Each mesh block will be written to (strand_id=block_id+1+strand_offset). Written to newer binary formats that are time-aware, defaults to 0.
Definition at line 154 of file tecplot_io.C.
References _strand_offset.
Referenced by write_binary().
{
return _strand_offset;
}
| double & libMesh::TecplotIO::time | ( | ) |
Solution time for transient data. Written to newer binary formats that are time-aware.
Definition at line 147 of file tecplot_io.C.
References _time.
{
return _time;
}
| void libMesh::TecplotIO::write | ( | const std::string & | fname | ) | [virtual] |
This method implements writing a mesh to a specified file.
Implements libMesh::MeshOutput< MeshBase >.
Definition at line 173 of file tecplot_io.C.
References binary(), libMesh::MeshOutput< MeshBase >::mesh(), libMesh::processor_id(), write_ascii(), and write_binary().
Referenced by libMesh::NameBasedIO::write().
{
if (this->mesh().processor_id() == 0)
{
if (this->binary())
this->write_binary (fname);
else
this->write_ascii (fname);
}
}
| void libMesh::TecplotIO::write_ascii | ( | const std::string & | fname, |
| const std::vector< Number > * | v = NULL, |
||
| const std::vector< std::string > * | solution_names = NULL |
||
| ) | [private] |
This method implements writing a mesh with nodal data to a specified file where the nodal data and variable names are optionally provided. This will write an ASCII file.
Definition at line 235 of file tecplot_io.C.
References _ascii_append, _time, std::abs(), libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::MeshOutput< MeshBase >::ascii_precision(), elem_dimension(), end, libMesh::MeshOutput< MeshBase >::mesh(), libMesh::MeshBase::n_active_sub_elem(), libMesh::MeshBase::n_nodes(), n_vars, libMesh::out, libMesh::MeshBase::point(), libMesh::processor_id(), libMesh::TECPLOT, and libMesh::TypeVector< T >::write_unformatted().
Referenced by write(), write_binary(), and write_nodal_data().
{
// Should only do this on processor 0!
libmesh_assert_equal_to (this->mesh().processor_id(), 0);
// Create an output stream, possibly in append mode.
std::ofstream out_stream(fname.c_str(), _ascii_append ? std::ofstream::app : std::ofstream::out);
// Make sure it opened correctly
if (!out_stream.good())
libmesh_file_error(fname.c_str());
// Get a constant reference to the mesh.
const MeshBase& the_mesh = MeshOutput<MeshBase>::mesh();
// Write header to stream
{
{
// TODO: We used to print out the SVN revision here when we did keyword expansions...
out_stream << "# For a description of the Tecplot format see the Tecplot User's guide.\n"
<< "#\n";
}
out_stream << "Variables=x,y,z";
if (solution_names != NULL)
for (unsigned int n=0; n<solution_names->size(); n++)
{
#ifdef LIBMESH_USE_REAL_NUMBERS
// Write variable names for real variables
out_stream << "," << (*solution_names)[n];
#else
// Write variable names for complex variables
out_stream << "," << "r_" << (*solution_names)[n]
<< "," << "i_" << (*solution_names)[n]
<< "," << "a_" << (*solution_names)[n];
#endif
}
out_stream << '\n';
out_stream << "Zone f=fepoint, n=" << the_mesh.n_nodes() << ", e=" << the_mesh.n_active_sub_elem();
// We cannot choose the element type simply based on the mesh
// dimension... there might be 1D elements living in a 3D mesh.
// So look at the elements which are actually in the Mesh, and
// choose either "lineseg", "quadrilateral", or "brick" depending
// on if the elements are 1, 2, or 3D.
// Write the element type we've determined to the header.
out_stream << ", et=";
switch (this->elem_dimension())
{
case 1:
out_stream << "lineseg";
break;
case 2:
out_stream << "quadrilateral";
break;
case 3:
out_stream << "brick";
break;
default:
libmesh_error_msg("Unsupported element dimension: " << this->elem_dimension());
}
// Output the time in the header
out_stream << ", t=\"T " << _time << "\"";
// Use default mesh color = black
out_stream << ", c=black\n";
} // finished writing header
for (unsigned int i=0; i<the_mesh.n_nodes(); i++)
{
// Print the point without a newline
the_mesh.point(i).write_unformatted(out_stream, false);
if ((v != NULL) && (solution_names != NULL))
{
const std::size_t n_vars = solution_names->size();
for (std::size_t c=0; c<n_vars; c++)
{
#ifdef LIBMESH_USE_REAL_NUMBERS
// Write real data
out_stream << std::setprecision(this->ascii_precision())
<< (*v)[i*n_vars + c] << " ";
#else
// Write complex data
out_stream << std::setprecision(this->ascii_precision())
<< (*v)[i*n_vars + c].real() << " "
<< (*v)[i*n_vars + c].imag() << " "
<< std::abs((*v)[i*n_vars + c]) << " ";
#endif
}
}
// Write a new line after the data for this node
out_stream << '\n';
}
MeshBase::const_element_iterator it = the_mesh.active_elements_begin();
const MeshBase::const_element_iterator end = the_mesh.active_elements_end();
for ( ; it != end; ++it)
(*it)->write_connectivity(out_stream, TECPLOT);
}
| void libMesh::TecplotIO::write_binary | ( | const std::string & | fname, |
| const std::vector< Number > * | vec = NULL, |
||
| const std::vector< std::string > * | solution_names = NULL |
||
| ) | [private] |
This method implements writing a mesh with nodal data to a specified file where the nodal data and variable names are optionally provided. This will write a binary file if the tecio.a library was found at compile time, otherwise a warning message will be printed and an ASCII file will be created.
Definition at line 357 of file tecplot_io.C.
References _subdomain_ids, _time, std::abs(), libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::MeshBase::active_subdomain_elements_begin(), libMesh::MeshBase::active_subdomain_elements_end(), elem_dimension(), end, libMesh::err, libMesh::ierr, std::max(), libMesh::MeshOutput< MeshBase >::mesh(), libMesh::MeshBase::mesh_dimension(), libMesh::MeshBase::n_active_sub_elem(), libMesh::MeshBase::n_nodes(), n_vars, libMesh::Quality::name(), libMesh::MeshBase::point(), libMesh::processor_id(), strand_offset(), libMesh::MeshBase::subdomain_name(), libMesh::TECPLOT, write_ascii(), and zone_title().
Referenced by write(), and write_nodal_data().
{
//-----------------------------------------------------------
// Call the ASCII output function if configure did not detect
// the Tecplot binary API
#ifndef LIBMESH_HAVE_TECPLOT_API
libMesh::err << "WARNING: Tecplot Binary files require the Tecplot API." << std::endl
<< "Continuing with ASCII output."
<< std::endl;
if (this->mesh().processor_id() == 0)
this->write_ascii (fname, vec, solution_names);
return;
//------------------------------------------------------------
// New binary formats, time aware and whatnot
#elif defined(LIBMESH_HAVE_TECPLOT_API_112)
// Get a constant reference to the mesh.
const MeshBase& the_mesh = MeshOutput<MeshBase>::mesh();
// Required variables
std::string tecplot_variable_names;
int
ierr = 0,
file_type = 0, // full
is_double = 0,
#ifdef DEBUG
tec_debug = 1,
#else
tec_debug = 0,
#endif
cell_type = -1,
nn_per_elem = -1;
switch (this->elem_dimension())
{
case 1:
cell_type = 1; // FELINESEG
nn_per_elem = 2;
break;
case 2:
cell_type = 3; // FEQUADRILATERAL
nn_per_elem = 4;
break;
case 3:
cell_type = 5; // FEBRICK
nn_per_elem = 8;
break;
default:
libmesh_error_msg("Unsupported element dimension: " << this->elem_dimension());
}
// Build a string containing all the variable names to pass to Tecplot
{
tecplot_variable_names += "x, y, z";
if (solution_names != NULL)
{
for (unsigned int name=0; name<solution_names->size(); name++)
{
#ifdef LIBMESH_USE_REAL_NUMBERS
tecplot_variable_names += ", ";
tecplot_variable_names += (*solution_names)[name];
#else
tecplot_variable_names += ", ";
tecplot_variable_names += "r_";
tecplot_variable_names += (*solution_names)[name];
tecplot_variable_names += ", ";
tecplot_variable_names += "i_";
tecplot_variable_names += (*solution_names)[name];
tecplot_variable_names += ", ";
tecplot_variable_names += "a_";
tecplot_variable_names += (*solution_names)[name];
#endif
}
}
}
// Instantiate a TecplotMacros interface. In 2D the most nodes per
// face should be 4, in 3D it's 8.
TecplotMacros tm(the_mesh.n_nodes(),
#ifdef LIBMESH_USE_REAL_NUMBERS
(3 + ((solution_names == NULL) ? 0 :
cast_int<unsigned int>(solution_names->size()))),
#else
(3 + 3*((solution_names == NULL) ? 0 :
cast_int<unsigned int>(solution_names->size()))),
#endif
the_mesh.n_active_sub_elem(),
nn_per_elem
);
// Copy the nodes and data to the TecplotMacros class. Note that we store
// everything as a float here since the eye doesn't require a double to
// understand what is going on
for (unsigned int v=0; v<the_mesh.n_nodes(); v++)
{
tm.nd(0,v) = static_cast<float>(the_mesh.point(v)(0));
tm.nd(1,v) = static_cast<float>(the_mesh.point(v)(1));
tm.nd(2,v) = static_cast<float>(the_mesh.point(v)(2));
if ((vec != NULL) &&
(solution_names != NULL))
{
const std::size_t n_vars = solution_names->size();
for (std::size_t c=0; c<n_vars; c++)
{
#ifdef LIBMESH_USE_REAL_NUMBERS
tm.nd((3+c),v) = static_cast<float>((*vec)[v*n_vars + c]);
#else
tm.nd((3+3*c),v) = static_cast<float>((*vec)[v*n_vars + c].real());
tm.nd((3+3*c+1),v) = static_cast<float>((*vec)[v*n_vars + c].imag());
tm.nd((3+3*c+2),v) = static_cast<float>(std::abs((*vec)[v*n_vars + c]));
#endif
}
}
}
// Initialize the file
ierr = TECINI112 (NULL,
const_cast<char*>(tecplot_variable_names.c_str()),
const_cast<char*>(fname.c_str()),
const_cast<char*>("."),
&file_type,
&tec_debug,
&is_double);
if (ierr)
libmesh_file_error(fname);
// A zone for each subdomain
bool firstzone=true;
for (std::set<subdomain_id_type>::const_iterator sbd_it=_subdomain_ids.begin();
sbd_it!=_subdomain_ids.end(); ++sbd_it)
{
// Copy the connectivity for this subdomain
{
MeshBase::const_element_iterator it = the_mesh.active_subdomain_elements_begin (*sbd_it);
const MeshBase::const_element_iterator end = the_mesh.active_subdomain_elements_end (*sbd_it);
unsigned int n_subcells_in_subdomain=0;
for (; it != end; ++it)
n_subcells_in_subdomain += (*it)->n_sub_elem();
// update the connectivty array to include only the elements in this subdomain
tm.set_n_cells (n_subcells_in_subdomain);
unsigned int te = 0;
for (it = the_mesh.active_subdomain_elements_begin (*sbd_it);
it != end; ++it)
{
std::vector<dof_id_type> conn;
for (unsigned int se=0; se<(*it)->n_sub_elem(); se++)
{
(*it)->connectivity(se, TECPLOT, conn);
for (unsigned int node=0; node<conn.size(); node++)
tm.cd(node,te) = conn[node];
te++;
}
}
}
// Ready to call the Tecplot API for this subdomain
{
int
num_nodes = static_cast<int>(the_mesh.n_nodes()),
num_cells = static_cast<int>(tm.n_cells),
num_faces = 0,
i_cell_max = 0,
j_cell_max = 0,
k_cell_max = 0,
strand_id = std::max(*sbd_it,static_cast<subdomain_id_type>(1)) + this->strand_offset(),
parent_zone = 0,
is_block = 1,
num_face_connect = 0,
face_neighbor_mode = 0,
tot_num_face_nodes = 0,
num_connect_boundary_faces = 0,
tot_num_boundary_connect = 0,
share_connect_from_zone=0;
std::vector<int>
passive_var_list (tm.n_vars, 0),
share_var_from_zone (tm.n_vars, 1); // We only write data for the first zone, all other
// zones will share from this one.
// get the subdomain name from libMesh, if there is one.
std::string subdomain_name = the_mesh.subdomain_name(*sbd_it);
std::ostringstream zone_name;
zone_name << this->zone_title();
// We will title this
// "{zone_title()}_{subdomain_name}", or
// "{zone_title()}_{subdomain_id}", or
// "{zone_title()}"
if (subdomain_name.size())
{
zone_name << "_";
zone_name << subdomain_name;
}
else if (_subdomain_ids.size() > 1)
{
zone_name << "_";
zone_name << *sbd_it;
}
ierr = TECZNE112 (const_cast<char*>(zone_name.str().c_str()),
&cell_type,
&num_nodes,
&num_cells,
&num_faces,
&i_cell_max,
&j_cell_max,
&k_cell_max,
&_time,
&strand_id,
&parent_zone,
&is_block,
&num_face_connect,
&face_neighbor_mode,
&tot_num_face_nodes,
&num_connect_boundary_faces,
&tot_num_boundary_connect,
&passive_var_list[0],
NULL, // = all are node centered
(firstzone) ? NULL : &share_var_from_zone[0],
&share_connect_from_zone);
if (ierr)
libmesh_file_error(fname);
// Write *all* the data for the first zone, then share it with the others
if (firstzone)
{
int total = cast_int<int>
#ifdef LIBMESH_USE_REAL_NUMBERS
((3 + ((solution_names == NULL) ? 0 : solution_names->size()))*num_nodes);
#else
((3 + 3*((solution_names == NULL) ? 0 : solution_names->size()))*num_nodes);
#endif
ierr = TECDAT112 (&total,
&tm.nodalData[0],
&is_double);
if (ierr)
libmesh_file_error(fname);
}
// Write the connectivity
ierr = TECNOD112 (&tm.connData[0]);
if (ierr)
libmesh_file_error(fname);
}
firstzone = false;
}
// Done, close the file.
ierr = TECEND112 ();
if (ierr)
libmesh_file_error(fname);
//------------------------------------------------------------
// Legacy binary format
#else
// Get a constant reference to the mesh.
const MeshBase& the_mesh = MeshOutput<MeshBase>::mesh();
// Tecplot binary output only good for dim=2,3
if (the_mesh.mesh_dimension() == 1)
{
this->write_ascii (fname, vec, solution_names);
return;
}
// Required variables
std::string tecplot_variable_names;
int is_double = 0,
tec_debug = 0,
cell_type = ((the_mesh.mesh_dimension()==2) ? (1) : (3));
// Build a string containing all the variable names to pass to Tecplot
{
tecplot_variable_names += "x, y, z";
if (solution_names != NULL)
{
for (unsigned int name=0; name<solution_names->size(); name++)
{
#ifdef LIBMESH_USE_REAL_NUMBERS
tecplot_variable_names += ", ";
tecplot_variable_names += (*solution_names)[name];
#else
tecplot_variable_names += ", ";
tecplot_variable_names += "r_";
tecplot_variable_names += (*solution_names)[name];
tecplot_variable_names += ", ";
tecplot_variable_names += "i_";
tecplot_variable_names += (*solution_names)[name];
tecplot_variable_names += ", ";
tecplot_variable_names += "a_";
tecplot_variable_names += (*solution_names)[name];
#endif
}
}
}
// Instantiate a TecplotMacros interface. In 2D the most nodes per
// face should be 4, in 3D it's 8.
TecplotMacros tm(cast_int<unsigned int>(the_mesh.n_nodes()),
cast_int<unsigned int>
#ifdef LIBMESH_USE_REAL_NUMBERS
(3 + ((solution_names == NULL) ? 0 : solution_names->size())),
#else
(3 + 3*((solution_names == NULL) ? 0 : solution_names->size())),
#endif
cast_int<unsigned int>
(the_mesh.n_active_sub_elem()),
((the_mesh.mesh_dimension() == 2) ? 4 : 8)
);
// Copy the nodes and data to the TecplotMacros class. Note that we store
// everything as a float here since the eye doesn't require a double to
// understand what is going on
for (unsigned int v=0; v<the_mesh.n_nodes(); v++)
{
tm.nd(0,v) = static_cast<float>(the_mesh.point(v)(0));
tm.nd(1,v) = static_cast<float>(the_mesh.point(v)(1));
tm.nd(2,v) = static_cast<float>(the_mesh.point(v)(2));
if ((vec != NULL) &&
(solution_names != NULL))
{
const unsigned int n_vars = solution_names->size();
for (unsigned int c=0; c<n_vars; c++)
{
#ifdef LIBMESH_USE_REAL_NUMBERS
tm.nd((3+c),v) = static_cast<float>((*vec)[v*n_vars + c]);
#else
tm.nd((3+3*c),v) = static_cast<float>((*vec)[v*n_vars + c].real());
tm.nd((3+3*c+1),v) = static_cast<float>((*vec)[v*n_vars + c].imag());
tm.nd((3+3*c+2),v) = static_cast<float>(std::abs((*vec)[v*n_vars + c]));
#endif
}
}
}
// Copy the connectivity
{
unsigned int te = 0;
MeshBase::const_element_iterator it = the_mesh.active_elements_begin();
const MeshBase::const_element_iterator end = the_mesh.active_elements_end();
for ( ; it != end; ++it)
{
std::vector<dof_id_type> conn;
for (unsigned int se=0; se<(*it)->n_sub_elem(); se++)
{
(*it)->connectivity(se, TECPLOT, conn);
for (unsigned int node=0; node<conn.size(); node++)
tm.cd(node,te) = conn[node];
te++;
}
}
}
// Ready to call the Tecplot API
{
int ierr = 0,
num_nodes = static_cast<int>(the_mesh.n_nodes()),
num_cells = static_cast<int>(the_mesh.n_active_sub_elem());
ierr = TECINI (NULL,
(char*) tecplot_variable_names.c_str(),
(char*) fname.c_str(),
(char*) ".",
&tec_debug,
&is_double);
if (ierr)
libmesh_file_error(fname);
ierr = TECZNE (NULL,
&num_nodes,
&num_cells,
&cell_type,
(char*) "FEBLOCK",
NULL);
if (ierr)
libmesh_file_error(fname);
int total =
#ifdef LIBMESH_USE_REAL_NUMBERS
((3 + ((solution_names == NULL) ? 0 : solution_names->size()))*num_nodes);
#else
((3 + 3*((solution_names == NULL) ? 0 : solution_names->size()))*num_nodes);
#endif
ierr = TECDAT (&total,
&tm.nodalData[0],
&is_double);
if (ierr)
libmesh_file_error(fname);
ierr = TECNOD (&tm.connData[0]);
if (ierr)
libmesh_file_error(fname);
ierr = TECEND ();
if (ierr)
libmesh_file_error(fname);
}
#endif
}
| virtual void libMesh::MeshOutput< MeshBase >::write_equation_systems | ( | const std::string & | , |
| const EquationSystems & | , | ||
| const std::set< std::string > * | system_names = NULL |
||
| ) | [virtual, inherited] |
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().
| void libMesh::TecplotIO::write_nodal_data | ( | const std::string & | fname, |
| const std::vector< Number > & | soln, | ||
| const std::vector< std::string > & | names | ||
| ) | [virtual] |
This method implements writing a mesh with nodal data to a specified file where the nodal data and variable names are provided.
Reimplemented from libMesh::MeshOutput< MeshBase >.
Definition at line 186 of file tecplot_io.C.
References binary(), libMesh::MeshOutput< MeshBase >::mesh(), libMesh::processor_id(), libMesh::START_LOG(), write_ascii(), and write_binary().
Referenced by libMesh::NameBasedIO::write_nodal_data().
{
START_LOG("write_nodal_data()", "TecplotIO");
if (this->mesh().processor_id() == 0)
{
if (this->binary())
this->write_binary (fname, &soln, &names);
else
this->write_ascii (fname, &soln, &names);
}
STOP_LOG("write_nodal_data()", "TecplotIO");
}
| std::string & libMesh::TecplotIO::zone_title | ( | ) |
The zone title to write.
Definition at line 161 of file tecplot_io.C.
References _zone_title.
Referenced by write_binary().
{
return _zone_title;
}
bool libMesh::TecplotIO::_ascii_append [private] |
If true, when writing in ASCII format, open the file in std::ofstream::app mode.
Definition at line 165 of file tecplot_io.h.
Referenced by ascii_append(), and write_ascii().
bool libMesh::TecplotIO::_binary [private] |
const bool libMesh::MeshOutput< MeshBase >::_is_parallel_format [protected, inherited] |
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 126 of file mesh_output.h.
Referenced by libMesh::FroIO::write(), libMesh::DivaIO::write(), libMesh::PostscriptIO::write(), and libMesh::EnsightIO::write().
int libMesh::TecplotIO::_strand_offset [private] |
Offset for Tecplot's STRANDID.
Definition at line 154 of file tecplot_io.h.
Referenced by strand_offset().
std::set<subdomain_id_type> libMesh::TecplotIO::_subdomain_ids [private] |
The subdomains in the mesh.
Definition at line 170 of file tecplot_io.h.
Referenced by TecplotIO(), and write_binary().
double libMesh::TecplotIO::_time [private] |
Solution time.
Definition at line 149 of file tecplot_io.h.
Referenced by time(), write_ascii(), and write_binary().
std::string libMesh::TecplotIO::_zone_title [private] |