libMesh::EnsightIO Class Reference

#include <ensight_io.h>

Inheritance diagram for libMesh::EnsightIO:

List of all members.

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 MeshBasemesh () 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

Detailed Description

This class implements writing meshes and solutions in Ensight's Gold format. author Camata

Definition at line 43 of file ensight_io.h.


Member Typedef Documentation

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.


Constructor & Destructor Documentation

libMesh::EnsightIO::EnsightIO ( const std::string &  filename,
const EquationSystems eq 
)

Constructor.

Definition at line 39 of file ensight_io.C.

References _ensight_file_name, libMesh::n_processors(), and libMesh::processor_id().

00039                                                                           :
00040   MeshOutput<MeshBase> (eq.get_mesh()),
00041   _equation_systems(eq)
00042 {
00043 
00044   if (libMesh::n_processors() == 1)
00045     _ensight_file_name = filename;
00046   else
00047     {
00048       std::stringstream tmp_file;
00049       tmp_file << filename << "_rank" << libMesh::processor_id();
00050       _ensight_file_name = tmp_file.str();
00051     }
00052 }

libMesh::EnsightIO::~EnsightIO (  ) 

Definition at line 56 of file ensight_io.C.

00057 {}


Member Function Documentation

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::EquationSystems::get_system(), and libMesh::EquationSystems::has_system().

00098 {
00099   libmesh_assert(_equation_systems.has_system(system_name));
00100   libmesh_assert(_equation_systems.get_system(system_name).has_variable(s));
00101 
00102   Scalars scl;
00103   scl.description = scl_description;
00104   scl.scalar_name = s;
00105 
00106   _systems_vars_map[system_name].EnsightScalars.push_back(scl);
00107 }

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::EquationSystems::get_system(), and libMesh::EquationSystems::has_system().

00080 {
00081   libmesh_assert(_equation_systems.has_system(system_name));
00082   libmesh_assert(_equation_systems.get_system(system_name).has_variable(u));
00083   libmesh_assert(_equation_systems.get_system(system_name).has_variable(v));
00084   libmesh_assert(_equation_systems.get_system(system_name).has_variable(w));
00085 
00086   Vectors vec;
00087   vec.description = vec_name;
00088   vec.components.push_back(u);
00089   vec.components.push_back(v);
00090   vec.components.push_back(w);
00091   _systems_vars_map[system_name].EnsightVectors.push_back(vec);
00092 }

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::EquationSystems::get_system(), and libMesh::EquationSystems::has_system().

00063 {
00064   libmesh_assert (_equation_systems.has_system(system_name));
00065   libmesh_assert (_equation_systems.get_system(system_name).has_variable(u));
00066   libmesh_assert (_equation_systems.get_system(system_name).has_variable(v));
00067 
00068   Vectors vec;
00069   vec.description = vec_description;
00070   vec.components.push_back(u);
00071   vec.components.push_back(v);
00072 
00073   _systems_vars_map[system_name].EnsightVectors.push_back(vec);
00074 }

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 libMeshEnums::EDGE2, libMeshEnums::EDGE3, libMeshEnums::HEX20, libMeshEnums::HEX27, libMeshEnums::HEX8, libMesh::out, libMeshEnums::PYRAMID5, libMeshEnums::QUAD4, libMeshEnums::QUAD8, libMeshEnums::QUAD9, libMeshEnums::TET10, libMeshEnums::TET4, libMeshEnums::TRI3, and libMeshEnums::TRI6.

Referenced by write_geometry_ascii().

00557 {
00558   switch(type){
00559   case EDGE2:
00560     std::strcpy(buffer,"bar2");
00561     break;
00562   case EDGE3:
00563     std::strcpy(buffer,"bar3");
00564     break;
00565   case QUAD4:
00566     std::strcpy(buffer,"quad4");
00567     break;
00568   case QUAD8:
00569     std::strcpy(buffer,"quad8");
00570     break;
00571   case QUAD9:
00572     libMesh::out<<"QUAD9: element not supported!"<<std::endl;
00573     libmesh_error();
00574     break;
00575 
00576   case TRI3:
00577     std::strcpy(buffer,"tria3");
00578     break;
00579   case TRI6:
00580     std::strcpy(buffer,"tria6");
00581     break;
00582   case TET4:
00583     std::strcpy(buffer,"tetra4");
00584     break;
00585   case TET10:
00586     std::strcpy(buffer,"tetra10");
00587     break;
00588   case HEX8:
00589     std::strcpy(buffer,"hexa8");
00590     break;
00591   case HEX20:
00592     std::strcpy(buffer,"hexa20");
00593     break;
00594   case HEX27:
00595     libMesh::out<<"HEX27: element not supported!"<<std::endl;
00596     libmesh_error();
00597     break;
00598   case PYRAMID5:
00599     std::strcpy(buffer,"pyramid5");
00600     break;
00601   default:
00602     break;
00603   }
00604 }

bool& libMesh::EnsightIO::has_mesh_refinement (  ) 
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().

00126 {
00127   this->write_ascii(time);
00128   this->write_case();
00129 }

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, and libMesh::MeshOutput< MeshBase >::mesh().

00114 {
00115   // We may need to gather a ParallelMesh to output it, making that
00116   // const qualifier in our constructor a dirty lie
00117   MeshSerializer serialize(const_cast<MeshBase&>(this->mesh()), !_is_parallel_format);
00118 
00119   _ensight_file_name = name;
00120   this->write();
00121 }

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().

00134 {
00135   _time_steps.push_back(time);
00136 
00137   this->write_geometry_ascii();
00138   this->write_solution_ascii();
00139 }

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, and libMesh::EnsightIO::Scalars::scalar_name.

Referenced by write().

00270 {
00271   std::stringstream case_file, geo_file;
00272   case_file << _ensight_file_name << ".case";
00273 
00274   FILE* fout = fopen(case_file.str().c_str(),"w");
00275   fprintf(fout,"FORMAT\n");
00276   fprintf(fout,"type:  ensight gold\n\n");
00277   fprintf(fout,"GEOMETRY\n");
00278 
00279   geo_file << _ensight_file_name << ".geo";
00280 
00281 
00282   fprintf(fout,"model:            1     %s***\n",geo_file.str().c_str());
00283 
00284   SystemsVarsMapIterator       sys      = _systems_vars_map.begin();
00285   const SystemsVarsMapIterator sys_end  = _systems_vars_map.end();
00286 
00287   // Write Variable per node section
00288   if ( sys != sys_end )
00289     fprintf(fout,"\n\nVARIABLE\n");
00290 
00291   for (; sys != sys_end; ++sys)
00292     {
00293       SystemsVarsValue value = *sys;
00294 
00295       for (unsigned int i=0; i < value.second.EnsightScalars.size(); i++)
00296         {
00297           std::stringstream scl_file;
00298           Scalars scalar = value.second.EnsightScalars[i];
00299           scl_file << _ensight_file_name
00300                    << "_" << scalar.scalar_name
00301                    << ".scl";
00302 
00303           fprintf(fout,"scalar per node:   1  %s %s***\n",scalar.description.c_str(), scl_file.str().c_str());
00304         }
00305 
00306       for (unsigned int i=0; i < value.second.EnsightVectors.size(); i++)
00307         {
00308           std::stringstream vec_file;
00309           Vectors vec = value.second.EnsightVectors[i];
00310           vec_file<<_ensight_file_name<<"_"<<vec.description<<".vec";
00311 
00312           fprintf(fout,"vector per node:      1    %s %s***\n",vec.description.c_str(), vec_file.str().c_str());
00313         }
00314 
00315       // Write time step section
00316       if( _time_steps.size() != 0)
00317         {
00318           fprintf(fout,"\n\nTIME\n");
00319           fprintf(fout,"time set:             1\n");
00320           fprintf(fout,"number of steps:   %10d\n", static_cast<int>(_time_steps.size()));
00321           fprintf(fout,"filename start number:   %10d\n", 0);
00322           fprintf(fout,"filename increment:  %10d\n", 1);
00323           fprintf(fout,"time values:\n");
00324           for (unsigned int i = 0; i < _time_steps.size(); i++)
00325             fprintf(fout,"%12.5e\n", _time_steps[i]);
00326         }
00327     }
00328   fclose(fout);
00329 }

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().

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(), libMeshEnums::HEX27, libMesh::MeshOutput< MeshBase >::mesh(), libMesh::Elem::n_nodes(), libMesh::Elem::node(), libMesh::Elem::point(), libMeshEnums::QUAD9, and libMesh::Elem::type().

Referenced by write_ascii().

00144 {
00145   std::ostringstream file;
00146   file << _ensight_file_name << ".geo";
00147 
00148   file << std::setw(3)
00149        << std::setprecision(0)
00150        << std::setfill('0')
00151        << std::right
00152        << _time_steps.size()-1;
00153 
00154   FILE* fout = fopen(file.str().c_str(),"w");
00155 
00156   char buffer[80];
00157 
00158   fprintf(fout,"EnSight Gold Geometry File Format\n");
00159   fprintf(fout,"Generated by \n");
00160   fprintf(fout,"node id off\n");
00161   fprintf(fout,"element id given\n");
00162   fprintf(fout,"part\n");
00163   fprintf(fout,"%10d\n",1);
00164   fprintf(fout,"uns-elements\n");
00165   fprintf(fout,"coordinates\n");
00166 
00167   // mapping between nodal index and your coordinates
00168   std::map<int, Point>                     mesh_nodes_map;
00169   typedef std::map <int, Point>::iterator  mesh_nodes_iterator;
00170   typedef std::pair<int, Point>            mesh_node_value;
00171 
00172   // Mapping between global and local indices
00173   std::map <int, int>       ensight_node_index;
00174 
00175   // Grouping elements of the same type
00176   std::map<ElemType, std::vector<const Elem*> >                    ensight_parts_map;
00177   typedef std::map<ElemType, std::vector<const Elem*> >::iterator  ensight_parts_iterator;
00178   typedef std::pair<ElemType, std::vector<const Elem*> >           ensight_parts_value;
00179 
00180   const MeshBase& the_mesh = MeshOutput<MeshBase>::mesh();
00181 
00182   MeshBase::const_element_iterator       el     = the_mesh.active_local_elements_begin();
00183   const MeshBase::const_element_iterator end_el = the_mesh.active_local_elements_end();
00184 
00185   for ( ; el != end_el ; ++el)
00186     {
00187       const Elem* elem = *el;
00188       ensight_parts_map[elem->type()].push_back(elem);
00189 
00190       for (unsigned int i = 0; i < elem->n_nodes(); i++)
00191         mesh_nodes_map[elem->node(i)] = elem->point(i);
00192     }
00193 
00194   // Write number of local points
00195   fprintf(fout,"%10d\n",static_cast<int>(mesh_nodes_map.size()));
00196 
00197   mesh_nodes_iterator           no_it = mesh_nodes_map.begin();
00198   const mesh_nodes_iterator no_end_it = mesh_nodes_map.end();
00199 
00200   // write x
00201   for(int i = 1; no_it != no_end_it; ++no_it, i++)
00202     {
00203       const mesh_node_value pn = *no_it;
00204       fprintf(fout,"%12.5e\n",static_cast<double>(pn.second(0)));
00205       ensight_node_index[pn.first] = i;
00206     }
00207 
00208   // write y
00209   no_it = mesh_nodes_map.begin();
00210   for(; no_it != no_end_it; ++no_it)
00211     {
00212       const mesh_node_value pn = *no_it;
00213       fprintf(fout,"%12.5e\n",static_cast<double>(pn.second(1)));
00214     }
00215 
00216   // write z
00217   no_it = mesh_nodes_map.begin();
00218   for(; no_it != no_end_it; ++no_it)
00219     {
00220       const mesh_node_value pn = *no_it;
00221       fprintf(fout,"%12.5e\n",static_cast<double>(pn.second(2)));
00222     }
00223 
00224   ensight_parts_iterator            parts_it  =  ensight_parts_map.begin();
00225   const ensight_parts_iterator  end_parts_it  =  ensight_parts_map.end();
00226 
00227   // Write parts
00228   for (; parts_it != end_parts_it; ++parts_it)
00229     {
00230       ensight_parts_value kvp = *parts_it;
00231 
00232       // Write element type
00233       elem_type_to_string(kvp.first,buffer);
00234       fprintf(fout,"\n%s\n", buffer);
00235 
00236       std::vector<const Elem*> elem_ref  = kvp.second;
00237 
00238       // Write number of element
00239       fprintf(fout,"%10d\n",static_cast<int>(elem_ref.size()));
00240 
00241       // Write element id
00242       for (unsigned int i = 0; i < elem_ref.size(); i++)
00243         fprintf(fout,"%10lu\n",static_cast<unsigned long>(elem_ref[i]->id()));
00244 
00245       // Write connectivity
00246       for (unsigned int i = 0; i < elem_ref.size(); i++)
00247         {
00248           for (unsigned int j = 0; j < elem_ref[i]->n_nodes(); j++) {
00249             // tests!
00250             if(kvp.first == QUAD9 && i==4)
00251               continue;
00252             // tests!
00253             if(kvp.first == HEX27 && (i==4    || i ==10 || i == 12 ||
00254                                       i == 13 || i ==14 || i == 16 || i == 22))
00255               continue;
00256 
00257             fprintf(fout,"%10d",ensight_node_index[elem_ref[i]->node(j)]);
00258           }
00259           fprintf(fout,"\n");
00260         }
00261     }
00262   fclose(fout);
00263 }

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.

00101   { libmesh_error(); }

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().

00356 {
00357   std::ostringstream scl_file;
00358   scl_file << _ensight_file_name << "_" << var_name << ".scl";
00359 
00360   scl_file << std::setw(3)
00361            << std::setprecision(0)
00362            << std::setfill('0')
00363            << std::right
00364            << _time_steps.size()-1;
00365 
00366   FILE * fout = fopen(scl_file.str().c_str(),"w");
00367 
00368   fprintf(fout,"Per node scalar value\n");
00369   fprintf(fout,"part\n");
00370   fprintf(fout,"%10d\n",1);
00371   fprintf(fout,"coordinates\n");
00372 
00373   const MeshBase& the_mesh = MeshOutput<MeshBase>::mesh();
00374 
00375   const unsigned int dim = the_mesh.mesh_dimension();
00376 
00377   const System &system = _equation_systems.get_system(sys);
00378 
00379   const DofMap& dof_map = system.get_dof_map();
00380 
00381 
00382   int var = system.variable_number(var_name);
00383 
00384 
00385   std::vector<dof_id_type> dof_indices;
00386   std::vector<dof_id_type> dof_indices_scl;
00387 
00388   // Now we will loop over all the elements in the mesh.
00389 
00390   MeshBase::const_element_iterator       el     = the_mesh.active_local_elements_begin();
00391   const MeshBase::const_element_iterator end_el = the_mesh.active_local_elements_end();
00392 
00393   typedef std::map<int,Real> map_local_soln;
00394   typedef map_local_soln::iterator local_soln_iterator;
00395 
00396   map_local_soln local_soln;
00397 
00398   std::vector<Number>       elem_soln;
00399   std::vector<Number>       nodal_soln;
00400 
00401   for ( ; el != end_el ; ++el){
00402 
00403     const Elem* elem = *el;
00404 
00405     const FEType& fe_type    = system.variable_type(var);
00406 
00407     dof_map.dof_indices (elem, dof_indices);
00408     dof_map.dof_indices (elem, dof_indices_scl, var);
00409 
00410     elem_soln.resize(dof_indices_scl.size());
00411 
00412     for (unsigned int i = 0; i < dof_indices_scl.size(); i++)
00413       elem_soln[i] = system.current_solution(dof_indices_scl[i]);
00414 
00415     FEInterface::nodal_soln (dim,fe_type, elem, elem_soln, nodal_soln);
00416 
00417     libmesh_assert_equal_to (nodal_soln.size(), elem->n_nodes());
00418 
00419 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
00420     libMesh::err << "Complex-valued Ensight output not yet supported" << std::endl;
00421     libmesh_not_implemented();
00422 #endif
00423 
00424     for (unsigned int n=0; n<elem->n_nodes(); n++)
00425       local_soln[elem->node(n)] = libmesh_real(nodal_soln[n]);
00426 
00427   }
00428 
00429   local_soln_iterator sol = local_soln.begin();
00430   const local_soln_iterator sol_end = local_soln.end();
00431   for(; sol != sol_end; ++sol)
00432     fprintf(fout,"%12.5e\n",static_cast<double>((*sol).second));
00433 
00434   fclose(fout);
00435 
00436 }

void libMesh::EnsightIO::write_solution_ascii (  )  [private]

Definition at line 333 of file ensight_io.C.

References _systems_vars_map, write_scalar_ascii(), and write_vector_ascii().

Referenced by write_ascii().

00334 {
00335 
00336   SystemsVarsMapIterator       sys     = _systems_vars_map.begin();
00337   const SystemsVarsMapIterator sys_end = _systems_vars_map.end();
00338 
00339   for (; sys != sys_end; ++sys)
00340     {
00341       SystemsVarsValue value = *sys;
00342 
00343       for (unsigned int i = 0; i < value.second.EnsightScalars.size(); i++)
00344         this->write_scalar_ascii(value.first,
00345                                  value.second.EnsightScalars[i].scalar_name);
00346 
00347       for (unsigned int i = 0; i < value.second.EnsightVectors.size(); i++)
00348         this->write_vector_ascii(value.first,
00349                                  value.second.EnsightVectors[i].components,
00350                                  value.second.EnsightVectors[i].description);
00351     }
00352 }

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().

00440 {
00441   std::ostringstream vec_file;
00442   vec_file<<_ensight_file_name<<"_"<<var_name<<".vec";
00443 
00444   vec_file << std::setw(3)
00445            << std::setprecision(0)
00446            << std::setfill('0')
00447            << std::right
00448            << _time_steps.size()-1;
00449 
00450   FILE * fout = fopen(vec_file.str().c_str(),"w");
00451   fprintf(fout,"Per vector per value\n");
00452   fprintf(fout,"part\n");
00453   fprintf(fout,"%10d\n",1);
00454   fprintf(fout,"coordinates\n");
00455 
00456   // Get a constant reference to the mesh object.
00457   const MeshBase& the_mesh = MeshOutput<MeshBase>::mesh();
00458 
00459   // The dimension that we are running
00460   const unsigned int dim = the_mesh.mesh_dimension();
00461 
00462   const System &system = _equation_systems.get_system(sys);
00463 
00464   const DofMap& dof_map = system.get_dof_map();
00465 
00466   const unsigned int u_var = system.variable_number(vec[0]);
00467   const unsigned int v_var = system.variable_number(vec[1]);
00468   const unsigned int w_var = (dim==3) ? system.variable_number(vec[2]) : 0;
00469 
00470   std::vector<dof_id_type> dof_indices;
00471   std::vector<dof_id_type> dof_indices_u;
00472   std::vector<dof_id_type> dof_indices_v;
00473   std::vector<dof_id_type> dof_indices_w;
00474 
00475   // Now we will loop over all the elements in the mesh.
00476   MeshBase::const_element_iterator       el     = the_mesh.active_local_elements_begin();
00477   const MeshBase::const_element_iterator end_el = the_mesh.active_local_elements_end();
00478 
00479   typedef std::map<int,std::vector<Real> > map_local_soln;
00480   typedef map_local_soln::iterator  local_soln_iterator;
00481 
00482   map_local_soln local_soln;
00483 
00484   for ( ; el != end_el ; ++el){
00485 
00486     const Elem* elem = *el;
00487 
00488     const FEType& fe_type    = system.variable_type(u_var);
00489 
00490     dof_map.dof_indices (elem, dof_indices);
00491     dof_map.dof_indices (elem, dof_indices_u,u_var);
00492     dof_map.dof_indices (elem, dof_indices_v,v_var);
00493     if(dim==3)  dof_map.dof_indices (elem, dof_indices,w_var);
00494 
00495 
00496     std::vector<Number>       elem_soln_u;
00497     std::vector<Number>       elem_soln_v;
00498     std::vector<Number>       elem_soln_w;
00499 
00500     std::vector<Number>       nodal_soln_u;
00501     std::vector<Number>       nodal_soln_v;
00502     std::vector<Number>       nodal_soln_w;
00503 
00504     elem_soln_u.resize(dof_indices_u.size());
00505     elem_soln_v.resize(dof_indices_v.size());
00506     if(dim == 3) elem_soln_w.resize(dof_indices_w.size());
00507 
00508     for (unsigned int i = 0; i < dof_indices_u.size(); i++)
00509       {
00510         elem_soln_u[i] = system.current_solution(dof_indices_u[i]);
00511         elem_soln_v[i] = system.current_solution(dof_indices_v[i]);
00512         if(dim==3) elem_soln_w[i] = system.current_solution(dof_indices_w[i]);
00513       }
00514 
00515     FEInterface::nodal_soln (dim,fe_type,elem,elem_soln_u,nodal_soln_u);
00516     FEInterface::nodal_soln (dim,fe_type,elem,elem_soln_v,nodal_soln_v);
00517     if(dim == 3) FEInterface::nodal_soln (dim,fe_type,elem,elem_soln_w,nodal_soln_w);
00518 
00519 
00520     libmesh_assert_equal_to (nodal_soln_u.size(), elem->n_nodes());
00521     libmesh_assert_equal_to (nodal_soln_v.size(), elem->n_nodes());
00522 
00523 #ifdef LIBMESH_ENABLE_COMPLEX
00524     libMesh::err << "Complex-valued Ensight output not yet supported" << std::endl;
00525     libmesh_not_implemented()
00526 #endif
00527 
00528     for (unsigned int n=0; n<elem->n_nodes(); n++)
00529       {
00530         std::vector<Real> node_vec(3);
00531         node_vec[0]= libmesh_real(nodal_soln_u[n]);
00532         node_vec[1]= libmesh_real(nodal_soln_v[n]);
00533         node_vec[2]=0.0;
00534         if(dim==3) node_vec[2]= libmesh_real(nodal_soln_w[n]);
00535         local_soln[elem->node(n)] = node_vec;
00536       }
00537 
00538   }
00539 
00540   local_soln_iterator sol = local_soln.begin();
00541   const local_soln_iterator sol_end = local_soln.end();
00542 
00543   for(; sol != sol_end; ++sol)
00544     fprintf(fout,"%12.5e\n",static_cast<double>((*sol).second[0]));
00545   sol = local_soln.begin();
00546   for(; sol != sol_end; ++sol)
00547     fprintf(fout,"%12.5e\n",static_cast<double>((*sol).second[1]));
00548   sol = local_soln.begin();
00549   for(; sol != sol_end; ++sol)
00550     fprintf(fout,"%12.5e\n",static_cast<double>((*sol).second[2]));
00551 
00552   fclose(fout);
00553 
00554 }


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(), libMesh::FroIO::write(), write(), and libMesh::DivaIO::write().

std::vector<double> libMesh::EnsightIO::_time_steps [private]

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

Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:17 UTC

Hosted By:
SourceForge.net Logo