libMesh::FroIO Class Reference
#include <fro_io.h>

Public Member Functions | |
| FroIO (const MeshBase &) | |
| virtual void | write (const std::string &) |
| virtual void | write_equation_systems (const std::string &, const EquationSystems &, const std::set< std::string > *system_names=NULL) |
| virtual void | write_nodal_data (const std::string &, const std::vector< Number > &, const std::vector< std::string > &) |
| unsigned int & | ascii_precision () |
Protected Member Functions | |
| const MeshBase & | mesh () const |
Protected Attributes | |
| const bool | _is_parallel_format |
Detailed Description
This class implements writing meshes in the .fro format used by the MIT ACDL. Valid only for triangular meshes.
Definition at line 46 of file fro_io.h.
Constructor & Destructor Documentation
| libMesh::FroIO::FroIO | ( | const MeshBase & | mesh_in | ) | [inline, explicit] |
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().
| const MeshBase & libMesh::MeshOutput< MeshBase >::mesh | ( | ) | const [protected, inherited] |
Returns the object as a read-only reference.
Referenced by libMesh::PostscriptIO::write(), write(), libMesh::EnsightIO::write(), libMesh::DivaIO::write(), libMesh::TecplotIO::write_ascii(), libMesh::MEDITIO::write_ascii(), libMesh::TecplotIO::write_binary(), libMesh::EnsightIO::write_geometry_ascii(), libMesh::EnsightIO::write_scalar_ascii(), libMesh::GnuPlotIO::write_solution(), libMesh::DivaIO::write_stream(), and libMesh::EnsightIO::write_vector_ascii().
| void libMesh::FroIO::write | ( | const std::string & | fname | ) | [virtual] |
This method implements writing a mesh to a specified file.
Implements libMesh::MeshOutput< MeshBase >.
Definition at line 41 of file fro_io.C.
References libMesh::MeshOutput< MeshBase >::_is_parallel_format, libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), bc_id, libMesh::MeshBase::boundary_info, libMesh::Elem::build_side(), libMesh::MeshBase::elem(), end, libMesh::err, libMesh::MeshOutput< MeshBase >::mesh(), libMesh::MeshBase::n_elem(), libMesh::MeshBase::n_nodes(), libMesh::MeshBase::point(), libMesh::processor_id(), side, and libMeshEnums::TRI3.
00042 { 00043 // We may need to gather a ParallelMesh to output it, making that 00044 // const qualifier in our constructor a dirty lie 00045 MeshSerializer serialize(const_cast<MeshBase&>(this->mesh()), !_is_parallel_format); 00046 00047 if (libMesh::processor_id() == 0) 00048 { 00049 // Open the output file stream 00050 std::ofstream out_stream (fname.c_str()); 00051 libmesh_assert (out_stream.good()); 00052 00053 // Make sure it opened correctly 00054 if (!out_stream.good()) 00055 libmesh_file_error(fname.c_str()); 00056 00057 // Get a reference to the mesh 00058 const MeshBase& the_mesh = MeshOutput<MeshBase>::mesh(); 00059 00060 // Write the header 00061 out_stream << the_mesh.n_elem() << " " 00062 << the_mesh.n_nodes() << " " 00063 << "0 0 " 00064 << the_mesh.boundary_info->n_boundary_ids() << " 1\n"; 00065 00066 // Write the nodes -- 1-based! 00067 for (unsigned int n=0; n<the_mesh.n_nodes(); n++) 00068 out_stream << n+1 << " \t" 00069 << std::scientific 00070 << std::setprecision(12) 00071 << the_mesh.point(n)(0) << " \t" 00072 << the_mesh.point(n)(1) << " \t" 00073 << 0. << '\n'; 00074 00075 // Write the elements -- 1-based! 00076 MeshBase::const_element_iterator it = the_mesh.active_elements_begin(); 00077 const MeshBase::const_element_iterator end = the_mesh.active_elements_end(); 00078 00079 for (unsigned int e=0 ; it != end; ++it) 00080 { 00081 // .fro likes TRI3's 00082 if ((*it)->type() != TRI3) 00083 { 00084 libMesh::err << "ERROR: .fro format only valid for triangles!\n" 00085 << " writing of " << fname << " aborted.\n" 00086 << std::endl; 00087 libmesh_error(); 00088 } 00089 00090 out_stream << ++e << " \t"; 00091 00092 for (unsigned int n=0; n<(*it)->n_nodes(); n++) 00093 out_stream << (*it)->node(n)+1 << " \t"; 00094 00095 // // LHS -> RHS Mapping, for inverted triangles 00096 // out_stream << (*it)->node(0)+1 << " \t"; 00097 // out_stream << (*it)->node(2)+1 << " \t"; 00098 // out_stream << (*it)->node(1)+1 << " \t"; 00099 00100 out_stream << "1\n"; 00101 } 00102 00103 // Write BCs. 00104 { 00105 const std::set<boundary_id_type>& bc_ids = 00106 the_mesh.boundary_info->get_boundary_ids(); 00107 00108 std::vector<dof_id_type> el; 00109 std::vector<unsigned short int> sl; 00110 std::vector<boundary_id_type> il; 00111 00112 the_mesh.boundary_info->build_side_list (el, sl, il); 00113 00114 00115 // Map the boundary ids into [1,n_bc_ids], 00116 // treat them one at a time. 00117 boundary_id_type bc_id=0; 00118 for (std::set<boundary_id_type>::const_iterator id = bc_ids.begin(); 00119 id != bc_ids.end(); ++id) 00120 { 00121 std::deque<dof_id_type> node_list; 00122 00123 std::map<dof_id_type, dof_id_type> 00124 forward_edges, backward_edges; 00125 00126 // Get all sides on this element with the relevant BC id. 00127 for (std::size_t e=0; e<el.size(); e++) 00128 if (il[e] == *id) 00129 { 00130 // need to build up node_list as a sorted array of edge nodes... 00131 // for the following: 00132 // a---b---c---d---e 00133 // node_list [ a b c d e]; 00134 // 00135 // the issue is just how to get this out of the elem/side based data structure. 00136 // the approach is to build up 'chain links' like this: 00137 // a---b b---c c---d d---e 00138 // and piece them together. 00139 // 00140 // so, for an arbitray edge n0---n1, we build the 00141 // "forward_edges" map n0-->n1 00142 // "backward_edges" map n1-->n0 00143 // and then start with one chain link, and add on... 00144 // 00145 AutoPtr<Elem> side = the_mesh.elem(el[e])->build_side(sl[e]); 00146 00147 const dof_id_type 00148 n0 = side->node(0), 00149 n1 = side->node(1); 00150 00151 // insert into forward-edge set 00152 forward_edges.insert (std::make_pair(n0, n1)); 00153 00154 // insert into backward-edge set 00155 backward_edges.insert (std::make_pair(n1, n0)); 00156 00157 // go ahead and add one edge to the list -- this will give us the beginning of a 00158 // chain to work from! 00159 if (node_list.empty()) 00160 { 00161 node_list.push_front(n0); 00162 node_list.push_back (n1); 00163 } 00164 } 00165 00166 // we now have the node_list with one edge, the forward_edges, and the backward_edges 00167 // the node_list will be filled when (node_list.size() == (n_edges+1)) 00168 // until that is the case simply add on to the beginning and end of the node_list, 00169 // building up a chain of ordered nodes... 00170 const std::size_t n_edges = forward_edges.size(); 00171 00172 while (node_list.size() != (n_edges+1)) 00173 { 00174 const dof_id_type 00175 front_node = node_list.front(), 00176 back_node = node_list.back(); 00177 00178 // look for front_pair in the backward_edges list 00179 { 00180 std::map<dof_id_type, dof_id_type>::iterator 00181 pos = backward_edges.find(front_node); 00182 00183 if (pos != backward_edges.end()) 00184 { 00185 node_list.push_front(pos->second); 00186 00187 backward_edges.erase(pos); 00188 } 00189 } 00190 00191 // look for back_pair in the forward_edges list 00192 { 00193 std::map<dof_id_type, dof_id_type>::iterator 00194 pos = forward_edges.find(back_node); 00195 00196 if (pos != forward_edges.end()) 00197 { 00198 node_list.push_back(pos->second); 00199 00200 forward_edges.erase(pos); 00201 } 00202 } 00203 00204 // libMesh::out << "node_list.size()=" << node_list.size() 00205 // << ", n_edges+1=" << n_edges+1 << std::endl; 00206 } 00207 00208 00209 out_stream << ++bc_id << " " << node_list.size() << '\n'; 00210 00211 std::deque<dof_id_type>::iterator pos = node_list.begin(); 00212 for ( ; pos != node_list.end(); ++pos) 00213 out_stream << *pos+1 << " \t0\n"; 00214 } 00215 } 00216 } 00217 }
| 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.
Referenced by libMesh::Nemesis_IO::write_timestep().
| virtual void libMesh::MeshOutput< MeshBase >::write_nodal_data | ( | const std::string & | , | |
| const std::vector< Number > & | , | |||
| const std::vector< std::string > & | ||||
| ) | [inline, virtual, inherited] |
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::GmshIO, libMesh::GMVIO, libMesh::GnuPlotIO, libMesh::MEDITIO, libMesh::Nemesis_IO, libMesh::TecplotIO, libMesh::UCDIO, and libMesh::VTKIO.
Definition at line 98 of file mesh_output.h.
Member Data Documentation
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::PostscriptIO::write(), write(), libMesh::EnsightIO::write(), and libMesh::DivaIO::write().
The documentation for this class was generated from the following files:
Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:26 UTC
Hosted By: