34 #include "libmesh/getpot.h" 
   35 #include "libmesh/gmv_io.h" 
   36 #include "libmesh/exodusII_io.h" 
   37 #include "libmesh/mesh.h" 
   38 #include "libmesh/nemesis_io.h" 
   39 #include "libmesh/tecplot_io.h" 
   40 #include "libmesh/vtk_io.h" 
   43 #include <sys/errno.h> 
   45 #include <sys/types.h> 
   51                                 const libMesh::Parallel::Communicator &comm )
 
   52     : _vis_output_file_prefix( input(
"vis-options/vis_output_file_prefix", 
"unknown" ) )
 
   54     unsigned int num_formats = input.vector_variable_size(
"vis-options/output_format");
 
   57     if( num_formats == 0 )
 
   60         std::string mesh_class = input(
"mesh-options/mesh_class", 
"default");
 
   62         if (mesh_class == 
"parallel")
 
   64         else if (mesh_class == 
"serial")
 
   66         else if (mesh_class == 
"default")
 
   69             libMesh::Mesh testdefault(comm);
 
   70             testdefault.delete_remote_elements();
 
   71             if (testdefault.is_serial())
 
   78             std::cerr << 
" Visualization::Visualization:" 
   79                       << 
" mesh-options/mesh_class had invalid value " << mesh_class
 
   85     for( 
unsigned int i = 0; i < num_formats; i++ )
 
   87         _output_format.push_back( input(
"vis-options/output_format", 
"DIE", i ) );
 
  106     ( SharedPtr<libMesh::EquationSystems> equation_system,
 
  107       const unsigned int time_step,
 
  108       const libMesh::Real time )
 
  110     std::stringstream suffix;
 
  114     std::string filename = this->_vis_output_file_prefix;
 
  115     filename+=
"."+suffix.str();
 
  117     this->dump_visualization( equation_system, filename, time );
 
  130     (SharedPtr<libMesh::EquationSystems> equation_system,
 
  132      const libMesh::ParameterVector & params)
 
  134     this->output_residual_sensitivities
 
  135       ( equation_system, system, params, 0, 0.0 );
 
  145     (SharedPtr<libMesh::EquationSystems> equation_system,
 
  147      const libMesh::ParameterVector & params)
 
  149     this->output_solution_sensitivities
 
  150       ( equation_system, system, params, 0, 0.0 );
 
  154     ( SharedPtr<libMesh::EquationSystems> equation_system,
 
  155       const std::string& filename_prefix,
 
  156       const libMesh::Real time )
 
  158     libMesh::MeshBase& mesh = equation_system->get_mesh();
 
  160     if( this->_vis_output_file_prefix == 
"unknown" )
 
  163         std::cout << 
" WARNING in Visualization::dump_visualization :" 
  164                   << 
" using 'unknown' as file prefix since it was not set " 
  170     if (!mesh.comm().rank())
 
  171       for (std::size_t pos = this->_vis_output_file_prefix.find(
'/');
 
  172            pos != std::string::npos;
 
  173            pos = this->_vis_output_file_prefix.find(
'/',++pos))
 
  174         if (mkdir(this->_vis_output_file_prefix.substr(0,pos).c_str(),
 
  175                   0777) != 0 && errno != EEXIST)
 
  176           libmesh_file_error(this->_vis_output_file_prefix.substr(0,pos));
 
  178     for( std::vector<std::string>::const_iterator format = _output_format.begin();
 
  179          format != _output_format.end();
 
  183         if ((*format) == 
"tecplot" ||
 
  186             std::string filename = filename_prefix+
".dat";
 
  187             libMesh::TecplotIO(mesh,
false).write_equation_systems( filename,
 
  190         else if ((*format) == 
"tecplot_binary" ||
 
  193             std::string filename = filename_prefix+
".plt";
 
  194             libMesh::TecplotIO(mesh,
true).write_equation_systems( filename,
 
  197         else if ((*format) == 
"gmv")
 
  199             std::string filename = filename_prefix+
".gmv";
 
  200             libMesh::GMVIO(mesh).write_equation_systems( filename,
 
  203         else if ((*format) == 
"pvtu")
 
  205             std::string filename = filename_prefix+
".pvtu";
 
  206             libMesh::VTKIO(mesh).write_equation_systems( filename,
 
  209         else if ((*format) == 
"ExodusII")
 
  211             std::string filename = filename_prefix+
".exo";
 
  216             libMesh::ExodusII_IO(mesh).write_timestep
 
  217               ( filename, *equation_system, 1, time );
 
  219         else if ((*format) == 
"Nemesis")
 
  221             std::string filename = filename_prefix+
".nem";
 
  226             libMesh::Nemesis_IO(mesh).write_timestep
 
  227               ( filename, *equation_system, 1, time );
 
  229         else if ((*format).find(
"xda") != std::string::npos ||
 
  230                  (*format).find(
"xdr") != std::string::npos)
 
  232             std::string filename = filename_prefix+
"."+(*format);
 
  233             const bool binary = ((*format).find(
"xdr") != std::string::npos);
 
  234             equation_system->write( filename,
 
  235                                     binary ? GRINSEnums::ENCODE : GRINSEnums::WRITE,
 
  236                                     libMesh::EquationSystems::WRITE_DATA |
 
  237                                     libMesh::EquationSystems::WRITE_ADDITIONAL_DATA );
 
  239         else if ((*format) == 
"mesh_only" )
 
  241             std::string filename = filename_prefix+
"_mesh.xda";
 
  242             equation_system->get_mesh().write( filename );
 
void output_residual(SharedPtr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system)
 
void dump_visualization(SharedPtr< libMesh::EquationSystems > equation_system, const std::string &filename_prefix, const libMesh::Real time)
 
void output_adjoint(SharedPtr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system)
 
void output_residual_sensitivities(SharedPtr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system, const libMesh::ParameterVector ¶ms)
 
std::vector< std::string > _output_format
 
std::string _vis_output_file_prefix
 
Visualization(const GetPot &input, const libMesh::Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
 
Interface with libMesh for solving Multiphysics problems. 
 
void output(SharedPtr< libMesh::EquationSystems > equation_system)
 
void output_solution_sensitivities(SharedPtr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system, const libMesh::ParameterVector ¶ms)