$extrastylesheet
#include <ensight_io.h>

Classes | |
| struct | Scalars |
| struct | SystemVars |
| struct | Vectors |
Public Member Functions | |
| EnsightIO (const std::string &filename, const EquationSystems &eq) | |
| ~EnsightIO () | |
| void | add_vector (const std::string &system, const std::string &vec_description, const std::string &u, const std::string &v) |
| void | add_vector (const std::string &system, const std::string &vec_description, const std::string &u, const std::string &v, const std::string &w) |
| void | add_scalar (const std::string &system, const std::string &scalar_description, const std::string &s) |
| virtual void | write (const std::string &name) |
| void | write (const double time=0) |
| bool & | has_mesh_refinement () |
| 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 |
Private Types | |
| typedef std::map< std::string, SystemVars > | SystemsVarsMap |
| typedef std::map< std::string, SystemVars >::iterator | SystemsVarsMapIterator |
| typedef std::pair< std::string, SystemVars > | SystemsVarsValue |
Private Member Functions | |
| void | write_ascii (const double time=0) |
| void | write_scalar_ascii (const std::string &sys, const std::string &var) |
| void | write_vector_ascii (const std::string &sys, const std::vector< std::string > &vec, const std::string &var_name) |
| void | write_solution_ascii () |
| void | write_geometry_ascii () |
| void | write_case () |
| void | elem_type_to_string (ElemType, char *) |
Private Attributes | |
| std::string | _ensight_file_name |
| std::vector< double > | _time_steps |
| SystemsVarsMap | _systems_vars_map |
| const EquationSystems & | _equation_systems |
This class implements writing meshes and solutions in Ensight's Gold format. author Camata
Definition at line 43 of file ensight_io.h.
typedef std::map<std::string, SystemVars> libMesh::EnsightIO::SystemsVarsMap [private] |
Definition at line 111 of file ensight_io.h.
typedef std::map<std::string, SystemVars>::iterator libMesh::EnsightIO::SystemsVarsMapIterator [private] |
Definition at line 112 of file ensight_io.h.
typedef std::pair<std::string, SystemVars> libMesh::EnsightIO::SystemsVarsValue [private] |
Definition at line 113 of file ensight_io.h.
| libMesh::EnsightIO::EnsightIO | ( | const std::string & | filename, |
| const EquationSystems & | eq | ||
| ) |
Constructor.
Definition at line 39 of file ensight_io.C.
References _ensight_file_name, _equation_systems, libMesh::ParallelObject::n_processors(), and libMesh::ParallelObject::processor_id().
: MeshOutput<MeshBase> (eq.get_mesh()), _equation_systems(eq) { if (_equation_systems.n_processors() == 1) _ensight_file_name = filename; else { std::stringstream tmp_file; tmp_file << filename << "_rank" << _equation_systems.processor_id(); _ensight_file_name = tmp_file.str(); } }
Definition at line 56 of file ensight_io.C.
{}
| void libMesh::EnsightIO::add_scalar | ( | const std::string & | system, |
| const std::string & | scalar_description, | ||
| const std::string & | s | ||
| ) |
add scalar: tell the EnsightIO interface that the variable s is a scalar
Definition at line 96 of file ensight_io.C.
References _equation_systems, _systems_vars_map, libMesh::EnsightIO::Scalars::description, libMesh::EquationSystems::get_system(), libMesh::EquationSystems::has_system(), and libMesh::libmesh_assert().
{
libmesh_assert(_equation_systems.has_system(system_name));
libmesh_assert(_equation_systems.get_system(system_name).has_variable(s));
Scalars scl;
scl.description = scl_description;
scl.scalar_name = s;
_systems_vars_map[system_name].EnsightScalars.push_back(scl);
}
| void libMesh::EnsightIO::add_vector | ( | const std::string & | system, |
| const std::string & | vec_description, | ||
| const std::string & | u, | ||
| const std::string & | v | ||
| ) |
add 2D vector: Tell the EnsightIO interface that the variables u and v are a vector. Note that u and v should be the same variables defined in the system.
Definition at line 61 of file ensight_io.C.
References _equation_systems, _systems_vars_map, libMesh::EnsightIO::Vectors::description, libMesh::EquationSystems::get_system(), libMesh::EquationSystems::has_system(), and libMesh::libmesh_assert().
{
libmesh_assert (_equation_systems.has_system(system_name));
libmesh_assert (_equation_systems.get_system(system_name).has_variable(u));
libmesh_assert (_equation_systems.get_system(system_name).has_variable(v));
Vectors vec;
vec.description = vec_description;
vec.components.push_back(u);
vec.components.push_back(v);
_systems_vars_map[system_name].EnsightVectors.push_back(vec);
}
| void libMesh::EnsightIO::add_vector | ( | const std::string & | system, |
| const std::string & | vec_description, | ||
| const std::string & | u, | ||
| const std::string & | v, | ||
| const std::string & | w | ||
| ) |
add 3D vector: tell the EnsightIO interface that the variables u, v and w are vector components
Definition at line 78 of file ensight_io.C.
References _equation_systems, _systems_vars_map, libMesh::EnsightIO::Vectors::description, libMesh::EquationSystems::get_system(), libMesh::EquationSystems::has_system(), and libMesh::libmesh_assert().
{
libmesh_assert(_equation_systems.has_system(system_name));
libmesh_assert(_equation_systems.get_system(system_name).has_variable(u));
libmesh_assert(_equation_systems.get_system(system_name).has_variable(v));
libmesh_assert(_equation_systems.get_system(system_name).has_variable(w));
Vectors vec;
vec.description = vec_name;
vec.components.push_back(u);
vec.components.push_back(v);
vec.components.push_back(w);
_systems_vars_map[system_name].EnsightVectors.push_back(vec);
}
| 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::EnsightIO::elem_type_to_string | ( | ElemType | type, |
| char * | buffer | ||
| ) | [private] |
Definition at line 556 of file ensight_io.C.
References libMesh::EDGE2, libMesh::EDGE3, libMesh::HEX20, libMesh::HEX27, libMesh::HEX8, libMesh::PYRAMID5, libMesh::QUAD4, libMesh::QUAD8, libMesh::QUAD9, libMesh::TET10, libMesh::TET4, libMesh::TRI3, and libMesh::TRI6.
Referenced by write_geometry_ascii().
{
switch(type){
case EDGE2:
std::strcpy(buffer,"bar2");
break;
case EDGE3:
std::strcpy(buffer,"bar3");
break;
case QUAD4:
std::strcpy(buffer,"quad4");
break;
case QUAD8:
std::strcpy(buffer,"quad8");
break;
case QUAD9:
libmesh_error_msg("QUAD9: element not supported!");
break;
case TRI3:
std::strcpy(buffer,"tria3");
break;
case TRI6:
std::strcpy(buffer,"tria6");
break;
case TET4:
std::strcpy(buffer,"tetra4");
break;
case TET10:
std::strcpy(buffer,"tetra10");
break;
case HEX8:
std::strcpy(buffer,"hexa8");
break;
case HEX20:
std::strcpy(buffer,"hexa20");
break;
case HEX27:
libmesh_error_msg("HEX27: element not supported!");
break;
case PYRAMID5:
std::strcpy(buffer,"pyramid5");
break;
default:
break;
}
}
| bool& libMesh::EnsightIO::has_mesh_refinement | ( | ) |
| const MeshBase & libMesh::MeshOutput< MeshBase >::mesh | ( | ) | const [protected, inherited] |
Returns the object as a read-only reference.
Referenced by libMesh::TecplotIO::elem_dimension(), libMesh::FroIO::write(), libMesh::DivaIO::write(), libMesh::TecplotIO::write(), libMesh::PostscriptIO::write(), libMesh::MEDITIO::write(), write(), libMesh::MEDITIO::write_ascii(), libMesh::TecplotIO::write_ascii(), libMesh::TecplotIO::write_binary(), write_geometry_ascii(), libMesh::TecplotIO::write_nodal_data(), libMesh::MEDITIO::write_nodal_data(), write_scalar_ascii(), libMesh::GnuPlotIO::write_solution(), libMesh::DivaIO::write_stream(), and write_vector_ascii().
| void libMesh::EnsightIO::write | ( | const std::string & | name | ) | [virtual] |
write solution
Implements libMesh::MeshOutput< MeshBase >.
Definition at line 113 of file ensight_io.C.
References _ensight_file_name, libMesh::MeshOutput< MeshBase >::_is_parallel_format, libMesh::MeshOutput< MeshBase >::mesh(), and libMesh::Quality::name().
{
// We may need to gather a ParallelMesh to output it, making that
// const qualifier in our constructor a dirty lie
MeshSerializer serialize(const_cast<MeshBase&>(this->mesh()), !_is_parallel_format);
_ensight_file_name = name;
this->write();
}
| void libMesh::EnsightIO::write | ( | const double | time = 0 | ) |
write solution
Definition at line 125 of file ensight_io.C.
References write_ascii(), and write_case().
{
this->write_ascii(time);
this->write_case();
}
| void libMesh::EnsightIO::write_ascii | ( | const double | time = 0 | ) | [private] |
Definition at line 133 of file ensight_io.C.
References _time_steps, write_geometry_ascii(), and write_solution_ascii().
Referenced by write().
{
_time_steps.push_back(time);
this->write_geometry_ascii();
this->write_solution_ascii();
}
| void libMesh::EnsightIO::write_case | ( | ) | [private] |
Definition at line 269 of file ensight_io.C.
References _ensight_file_name, _systems_vars_map, _time_steps, libMesh::EnsightIO::Vectors::description, libMesh::EnsightIO::Scalars::description, libMesh::EnsightIO::Scalars::scalar_name, and libMesh::sys.
Referenced by write().
{
std::stringstream case_file, geo_file;
case_file << _ensight_file_name << ".case";
FILE* fout = fopen(case_file.str().c_str(),"w");
fprintf(fout,"FORMAT\n");
fprintf(fout,"type: ensight gold\n\n");
fprintf(fout,"GEOMETRY\n");
geo_file << _ensight_file_name << ".geo";
fprintf(fout,"model: 1 %s***\n",geo_file.str().c_str());
SystemsVarsMapIterator sys = _systems_vars_map.begin();
const SystemsVarsMapIterator sys_end = _systems_vars_map.end();
// Write Variable per node section
if ( sys != sys_end )
fprintf(fout,"\n\nVARIABLE\n");
for (; sys != sys_end; ++sys)
{
SystemsVarsValue value = *sys;
for (unsigned int i=0; i < value.second.EnsightScalars.size(); i++)
{
std::stringstream scl_file;
Scalars scalar = value.second.EnsightScalars[i];
scl_file << _ensight_file_name
<< "_" << scalar.scalar_name
<< ".scl";
fprintf(fout,"scalar per node: 1 %s %s***\n",scalar.description.c_str(), scl_file.str().c_str());
}
for (unsigned int i=0; i < value.second.EnsightVectors.size(); i++)
{
std::stringstream vec_file;
Vectors vec = value.second.EnsightVectors[i];
vec_file<<_ensight_file_name<<"_"<<vec.description<<".vec";
fprintf(fout,"vector per node: 1 %s %s***\n",vec.description.c_str(), vec_file.str().c_str());
}
// Write time step section
if( _time_steps.size() != 0)
{
fprintf(fout,"\n\nTIME\n");
fprintf(fout,"time set: 1\n");
fprintf(fout,"number of steps: %10d\n", static_cast<int>(_time_steps.size()));
fprintf(fout,"filename start number: %10d\n", 0);
fprintf(fout,"filename increment: %10d\n", 1);
fprintf(fout,"time values:\n");
for (unsigned int i = 0; i < _time_steps.size(); i++)
fprintf(fout,"%12.5e\n", _time_steps[i]);
}
}
fclose(fout);
}
| 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::EnsightIO::write_geometry_ascii | ( | ) | [private] |
Definition at line 143 of file ensight_io.C.
References _ensight_file_name, _time_steps, libMesh::MeshBase::active_local_elements_begin(), libMesh::MeshBase::active_local_elements_end(), elem_type_to_string(), libMesh::HEX27, libMesh::MeshOutput< MeshBase >::mesh(), libMesh::Elem::n_nodes(), libMesh::Elem::node(), libMesh::Elem::point(), libMesh::QUAD9, and libMesh::Elem::type().
Referenced by write_ascii().
{
std::ostringstream file;
file << _ensight_file_name << ".geo";
file << std::setw(3)
<< std::setprecision(0)
<< std::setfill('0')
<< std::right
<< _time_steps.size()-1;
FILE* fout = fopen(file.str().c_str(),"w");
char buffer[80];
fprintf(fout,"EnSight Gold Geometry File Format\n");
fprintf(fout,"Generated by \n");
fprintf(fout,"node id off\n");
fprintf(fout,"element id given\n");
fprintf(fout,"part\n");
fprintf(fout,"%10d\n",1);
fprintf(fout,"uns-elements\n");
fprintf(fout,"coordinates\n");
// mapping between nodal index and your coordinates
std::map<int, Point> mesh_nodes_map;
typedef std::map <int, Point>::iterator mesh_nodes_iterator;
typedef std::pair<int, Point> mesh_node_value;
// Mapping between global and local indices
std::map <int, int> ensight_node_index;
// Grouping elements of the same type
std::map<ElemType, std::vector<const Elem*> > ensight_parts_map;
typedef std::map<ElemType, std::vector<const Elem*> >::iterator ensight_parts_iterator;
typedef std::pair<ElemType, std::vector<const Elem*> > ensight_parts_value;
const MeshBase& the_mesh = MeshOutput<MeshBase>::mesh();
MeshBase::const_element_iterator el = the_mesh.active_local_elements_begin();
const MeshBase::const_element_iterator end_el = the_mesh.active_local_elements_end();
for ( ; el != end_el ; ++el)
{
const Elem* elem = *el;
ensight_parts_map[elem->type()].push_back(elem);
for (unsigned int i = 0; i < elem->n_nodes(); i++)
mesh_nodes_map[elem->node(i)] = elem->point(i);
}
// Write number of local points
fprintf(fout,"%10d\n",static_cast<int>(mesh_nodes_map.size()));
mesh_nodes_iterator no_it = mesh_nodes_map.begin();
const mesh_nodes_iterator no_end_it = mesh_nodes_map.end();
// write x
for(int i = 1; no_it != no_end_it; ++no_it, i++)
{
const mesh_node_value pn = *no_it;
fprintf(fout,"%12.5e\n",static_cast<double>(pn.second(0)));
ensight_node_index[pn.first] = i;
}
// write y
no_it = mesh_nodes_map.begin();
for(; no_it != no_end_it; ++no_it)
{
const mesh_node_value pn = *no_it;
fprintf(fout,"%12.5e\n",static_cast<double>(pn.second(1)));
}
// write z
no_it = mesh_nodes_map.begin();
for(; no_it != no_end_it; ++no_it)
{
const mesh_node_value pn = *no_it;
fprintf(fout,"%12.5e\n",static_cast<double>(pn.second(2)));
}
ensight_parts_iterator parts_it = ensight_parts_map.begin();
const ensight_parts_iterator end_parts_it = ensight_parts_map.end();
// Write parts
for (; parts_it != end_parts_it; ++parts_it)
{
ensight_parts_value kvp = *parts_it;
// Write element type
elem_type_to_string(kvp.first,buffer);
fprintf(fout,"\n%s\n", buffer);
std::vector<const Elem*> elem_ref = kvp.second;
// Write number of element
fprintf(fout,"%10d\n",static_cast<int>(elem_ref.size()));
// Write element id
for (unsigned int i = 0; i < elem_ref.size(); i++)
fprintf(fout,"%10lu\n",static_cast<unsigned long>(elem_ref[i]->id()));
// Write connectivity
for (unsigned int i = 0; i < elem_ref.size(); i++)
{
for (unsigned int j = 0; j < elem_ref[i]->n_nodes(); j++) {
// tests!
if(kvp.first == QUAD9 && i==4)
continue;
// tests!
if(kvp.first == HEX27 && (i==4 || i ==10 || i == 12 ||
i == 13 || i ==14 || i == 16 || i == 22))
continue;
fprintf(fout,"%10d",ensight_node_index[elem_ref[i]->node(j)]);
}
fprintf(fout,"\n");
}
}
fclose(fout);
}
| 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::GMVIO, libMesh::NameBasedIO, libMesh::GmshIO, libMesh::Nemesis_IO, libMesh::VTKIO, libMesh::UCDIO, libMesh::MEDITIO, libMesh::GnuPlotIO, and libMesh::TecplotIO.
Definition at line 98 of file mesh_output.h.
{ libmesh_not_implemented(); }
| void libMesh::EnsightIO::write_scalar_ascii | ( | const std::string & | sys, |
| const std::string & | var | ||
| ) | [private] |
Definition at line 355 of file ensight_io.C.
References _ensight_file_name, _equation_systems, _time_steps, libMesh::MeshBase::active_local_elements_begin(), libMesh::MeshBase::active_local_elements_end(), libMesh::System::current_solution(), libMesh::DofMap::dof_indices(), libMesh::err, libMesh::System::get_dof_map(), libMesh::EquationSystems::get_system(), libMesh::libmesh_real(), libMesh::MeshOutput< MeshBase >::mesh(), libMesh::MeshBase::mesh_dimension(), libMesh::Elem::n_nodes(), libMesh::FEInterface::nodal_soln(), libMesh::Elem::node(), libMesh::System::variable_number(), and libMesh::System::variable_type().
Referenced by write_solution_ascii().
{
std::ostringstream scl_file;
scl_file << _ensight_file_name << "_" << var_name << ".scl";
scl_file << std::setw(3)
<< std::setprecision(0)
<< std::setfill('0')
<< std::right
<< _time_steps.size()-1;
FILE * fout = fopen(scl_file.str().c_str(),"w");
fprintf(fout,"Per node scalar value\n");
fprintf(fout,"part\n");
fprintf(fout,"%10d\n",1);
fprintf(fout,"coordinates\n");
const MeshBase& the_mesh = MeshOutput<MeshBase>::mesh();
const unsigned int dim = the_mesh.mesh_dimension();
const System &system = _equation_systems.get_system(sys);
const DofMap& dof_map = system.get_dof_map();
int var = system.variable_number(var_name);
std::vector<dof_id_type> dof_indices;
std::vector<dof_id_type> dof_indices_scl;
// Now we will loop over all the elements in the mesh.
MeshBase::const_element_iterator el = the_mesh.active_local_elements_begin();
const MeshBase::const_element_iterator end_el = the_mesh.active_local_elements_end();
typedef std::map<int,Real> map_local_soln;
typedef map_local_soln::iterator local_soln_iterator;
map_local_soln local_soln;
std::vector<Number> elem_soln;
std::vector<Number> nodal_soln;
for ( ; el != end_el ; ++el){
const Elem* elem = *el;
const FEType& fe_type = system.variable_type(var);
dof_map.dof_indices (elem, dof_indices);
dof_map.dof_indices (elem, dof_indices_scl, var);
elem_soln.resize(dof_indices_scl.size());
for (unsigned int i = 0; i < dof_indices_scl.size(); i++)
elem_soln[i] = system.current_solution(dof_indices_scl[i]);
FEInterface::nodal_soln (dim,fe_type, elem, elem_soln, nodal_soln);
libmesh_assert_equal_to (nodal_soln.size(), elem->n_nodes());
#ifdef LIBMESH_USE_COMPLEX_NUMBERS
libMesh::err << "Complex-valued Ensight output not yet supported" << std::endl;
libmesh_not_implemented();
#endif
for (unsigned int n=0; n<elem->n_nodes(); n++)
local_soln[elem->node(n)] = libmesh_real(nodal_soln[n]);
}
local_soln_iterator sol = local_soln.begin();
const local_soln_iterator sol_end = local_soln.end();
for(; sol != sol_end; ++sol)
fprintf(fout,"%12.5e\n",static_cast<double>((*sol).second));
fclose(fout);
}
| void libMesh::EnsightIO::write_solution_ascii | ( | ) | [private] |
Definition at line 333 of file ensight_io.C.
References _systems_vars_map, libMesh::sys, write_scalar_ascii(), and write_vector_ascii().
Referenced by write_ascii().
{
SystemsVarsMapIterator sys = _systems_vars_map.begin();
const SystemsVarsMapIterator sys_end = _systems_vars_map.end();
for (; sys != sys_end; ++sys)
{
SystemsVarsValue value = *sys;
for (unsigned int i = 0; i < value.second.EnsightScalars.size(); i++)
this->write_scalar_ascii(value.first,
value.second.EnsightScalars[i].scalar_name);
for (unsigned int i = 0; i < value.second.EnsightVectors.size(); i++)
this->write_vector_ascii(value.first,
value.second.EnsightVectors[i].components,
value.second.EnsightVectors[i].description);
}
}
| void libMesh::EnsightIO::write_vector_ascii | ( | const std::string & | sys, |
| const std::vector< std::string > & | vec, | ||
| const std::string & | var_name | ||
| ) | [private] |
Definition at line 439 of file ensight_io.C.
References _ensight_file_name, _equation_systems, _time_steps, libMesh::MeshBase::active_local_elements_begin(), libMesh::MeshBase::active_local_elements_end(), libMesh::System::current_solution(), libMesh::DofMap::dof_indices(), libMesh::err, libMesh::System::get_dof_map(), libMesh::EquationSystems::get_system(), libMesh::libmesh_real(), libMesh::MeshOutput< MeshBase >::mesh(), libMesh::MeshBase::mesh_dimension(), libMesh::Elem::n_nodes(), libMesh::FEInterface::nodal_soln(), libMesh::Elem::node(), libMesh::System::variable_number(), and libMesh::System::variable_type().
Referenced by write_solution_ascii().
{
std::ostringstream vec_file;
vec_file<<_ensight_file_name<<"_"<<var_name<<".vec";
vec_file << std::setw(3)
<< std::setprecision(0)
<< std::setfill('0')
<< std::right
<< _time_steps.size()-1;
FILE * fout = fopen(vec_file.str().c_str(),"w");
fprintf(fout,"Per vector per value\n");
fprintf(fout,"part\n");
fprintf(fout,"%10d\n",1);
fprintf(fout,"coordinates\n");
// Get a constant reference to the mesh object.
const MeshBase& the_mesh = MeshOutput<MeshBase>::mesh();
// The dimension that we are running
const unsigned int dim = the_mesh.mesh_dimension();
const System &system = _equation_systems.get_system(sys);
const DofMap& dof_map = system.get_dof_map();
const unsigned int u_var = system.variable_number(vec[0]);
const unsigned int v_var = system.variable_number(vec[1]);
const unsigned int w_var = (dim==3) ? system.variable_number(vec[2]) : 0;
std::vector<dof_id_type> dof_indices;
std::vector<dof_id_type> dof_indices_u;
std::vector<dof_id_type> dof_indices_v;
std::vector<dof_id_type> dof_indices_w;
// Now we will loop over all the elements in the mesh.
MeshBase::const_element_iterator el = the_mesh.active_local_elements_begin();
const MeshBase::const_element_iterator end_el = the_mesh.active_local_elements_end();
typedef std::map<int,std::vector<Real> > map_local_soln;
typedef map_local_soln::iterator local_soln_iterator;
map_local_soln local_soln;
for ( ; el != end_el ; ++el){
const Elem* elem = *el;
const FEType& fe_type = system.variable_type(u_var);
dof_map.dof_indices (elem, dof_indices);
dof_map.dof_indices (elem, dof_indices_u,u_var);
dof_map.dof_indices (elem, dof_indices_v,v_var);
if(dim==3) dof_map.dof_indices (elem, dof_indices,w_var);
std::vector<Number> elem_soln_u;
std::vector<Number> elem_soln_v;
std::vector<Number> elem_soln_w;
std::vector<Number> nodal_soln_u;
std::vector<Number> nodal_soln_v;
std::vector<Number> nodal_soln_w;
elem_soln_u.resize(dof_indices_u.size());
elem_soln_v.resize(dof_indices_v.size());
if(dim == 3) elem_soln_w.resize(dof_indices_w.size());
for (unsigned int i = 0; i < dof_indices_u.size(); i++)
{
elem_soln_u[i] = system.current_solution(dof_indices_u[i]);
elem_soln_v[i] = system.current_solution(dof_indices_v[i]);
if(dim==3) elem_soln_w[i] = system.current_solution(dof_indices_w[i]);
}
FEInterface::nodal_soln (dim,fe_type,elem,elem_soln_u,nodal_soln_u);
FEInterface::nodal_soln (dim,fe_type,elem,elem_soln_v,nodal_soln_v);
if(dim == 3) FEInterface::nodal_soln (dim,fe_type,elem,elem_soln_w,nodal_soln_w);
libmesh_assert_equal_to (nodal_soln_u.size(), elem->n_nodes());
libmesh_assert_equal_to (nodal_soln_v.size(), elem->n_nodes());
#ifdef LIBMESH_ENABLE_COMPLEX
libMesh::err << "Complex-valued Ensight output not yet supported" << std::endl;
libmesh_not_implemented()
#endif
for (unsigned int n=0; n<elem->n_nodes(); n++)
{
std::vector<Real> node_vec(3);
node_vec[0]= libmesh_real(nodal_soln_u[n]);
node_vec[1]= libmesh_real(nodal_soln_v[n]);
node_vec[2]=0.0;
if(dim==3) node_vec[2]= libmesh_real(nodal_soln_w[n]);
local_soln[elem->node(n)] = node_vec;
}
}
local_soln_iterator sol = local_soln.begin();
const local_soln_iterator sol_end = local_soln.end();
for(; sol != sol_end; ++sol)
fprintf(fout,"%12.5e\n",static_cast<double>((*sol).second[0]));
sol = local_soln.begin();
for(; sol != sol_end; ++sol)
fprintf(fout,"%12.5e\n",static_cast<double>((*sol).second[1]));
sol = local_soln.begin();
for(; sol != sol_end; ++sol)
fprintf(fout,"%12.5e\n",static_cast<double>((*sol).second[2]));
fclose(fout);
}
std::string libMesh::EnsightIO::_ensight_file_name [private] |
Definition at line 128 of file ensight_io.h.
Referenced by EnsightIO(), write(), write_case(), write_geometry_ascii(), write_scalar_ascii(), and write_vector_ascii().
const EquationSystems& libMesh::EnsightIO::_equation_systems [private] |
Definition at line 131 of file ensight_io.h.
Referenced by add_scalar(), add_vector(), EnsightIO(), write_scalar_ascii(), and write_vector_ascii().
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 write().
Definition at line 130 of file ensight_io.h.
Referenced by add_scalar(), add_vector(), write_case(), and write_solution_ascii().
std::vector<double> libMesh::EnsightIO::_time_steps [private] |
Definition at line 129 of file ensight_io.h.
Referenced by write_ascii(), write_case(), write_geometry_ascii(), write_scalar_ascii(), and write_vector_ascii().