$extrastylesheet
#include <exodusII_io_helper.h>

Classes | |
| class | Conversion |
| class | ElementMaps |
| class | NamesData |
Public Types | |
| enum | ExodusVarType { NODAL = 0, ELEMENTAL = 1, GLOBAL = 2 } |
Public Member Functions | |
| ExodusII_IO_Helper (const ParallelObject &parent, bool v=false, bool run_only_on_proc0=true, bool single_precision=false) | |
| virtual | ~ExodusII_IO_Helper () |
| const char * | get_elem_type () const |
| void | open (const char *filename, bool read_only) |
| void | read_header () |
| void | print_header () |
| void | read_nodes () |
| void | read_node_num_map () |
| void | print_nodes (std::ostream &out=libMesh::out) |
| void | read_block_info () |
| int | get_block_id (int index) |
| std::string | get_block_name (int index) |
| int | get_side_set_id (int index) |
| std::string | get_side_set_name (int index) |
| int | get_node_set_id (int index) |
| std::string | get_node_set_name (int index) |
| void | read_elem_in_block (int block) |
| void | read_elem_num_map () |
| void | read_sideset_info () |
| void | read_nodeset_info () |
| void | read_sideset (int id, int offset) |
| void | read_nodeset (int id) |
| void | close () |
| int | inquire (int req_info, std::string error_msg="") |
| void | read_time_steps () |
| void | read_num_time_steps () |
| void | read_nodal_var_values (std::string nodal_var_name, int time_step) |
| void | read_elemental_var_values (std::string elemental_var_name, int time_step) |
| virtual void | create (std::string filename) |
| virtual void | initialize (std::string title, const MeshBase &mesh, bool use_discontinuous=false) |
| virtual void | write_nodal_coordinates (const MeshBase &mesh, bool use_discontinuous=false) |
| virtual void | write_elements (const MeshBase &mesh, bool use_discontinuous=false) |
| virtual void | write_sidesets (const MeshBase &mesh) |
| virtual void | write_nodesets (const MeshBase &mesh) |
| void | initialize_element_variables (std::vector< std::string > names) |
| void | initialize_nodal_variables (std::vector< std::string > names) |
| void | initialize_global_variables (std::vector< std::string > names) |
| void | write_timestep (int timestep, Real time) |
| void | write_element_values (const MeshBase &mesh, const std::vector< Real > &values, int timestep) |
| void | write_nodal_values (int var_id, const std::vector< Real > &values, int timestep) |
| void | write_information_records (const std::vector< std::string > &records) |
| void | write_global_values (const std::vector< Real > &values, int timestep) |
| void | use_mesh_dimension_instead_of_spatial_dimension (bool val) |
| void | set_coordinate_offset (Point p) |
| std::vector< std::string > | get_complex_names (const std::vector< std::string > &names) const |
| void | message (const std::string &msg) |
| void | message (const std::string &msg, int i) |
| void | read_var_names (ExodusVarType type) |
| const Parallel::Communicator & | comm () const |
| processor_id_type | n_processors () const |
| processor_id_type | processor_id () const |
Public Attributes | |
| int | ex_id |
| int | ex_err |
| int | num_dim |
| int | num_global_vars |
| int | num_nodes |
| int | num_elem |
| int | num_elem_blk |
| int | num_node_sets |
| int | num_side_sets |
| int | num_elem_this_blk |
| int | num_nodes_per_elem |
| int | num_attr |
| int | num_elem_all_sidesets |
| std::vector< int > | block_ids |
| std::vector< int > | connect |
| std::vector< int > | ss_ids |
| std::vector< int > | nodeset_ids |
| std::vector< int > | num_sides_per_set |
| std::vector< int > | num_nodes_per_set |
| std::vector< int > | num_df_per_set |
| std::vector< int > | num_node_df_per_set |
| std::vector< int > | elem_list |
| std::vector< int > | side_list |
| std::vector< int > | node_list |
| std::vector< int > | id_list |
| std::vector< int > | node_num_map |
| std::vector< int > | elem_num_map |
| std::vector< Real > | x |
| std::vector< Real > | y |
| std::vector< Real > | z |
| std::vector< char > | title |
| std::vector< char > | elem_type |
| std::map< int, int > | libmesh_elem_num_to_exodus |
| std::vector< int > | exodus_elem_num_to_libmesh |
| std::map< int, int > | libmesh_node_num_to_exodus |
| std::vector< int > | exodus_node_num_to_libmesh |
| int | num_time_steps |
| std::vector< Real > | time_steps |
| int | num_nodal_vars |
| std::vector< std::string > | nodal_var_names |
| std::vector< Real > | nodal_var_values |
| int | num_elem_vars |
| std::vector< std::string > | elem_var_names |
| std::vector< Real > | elem_var_values |
| std::vector< std::string > | global_var_names |
| std::map< int, std::string > | id_to_block_names |
| std::map< int, std::string > | id_to_ss_names |
| std::map< int, std::string > | id_to_ns_names |
| bool | verbose |
| bool | opened_for_writing |
| bool | opened_for_reading |
| std::string | current_filename |
Protected Attributes | |
| bool | _run_only_on_proc0 |
| bool | _elem_vars_initialized |
| bool | _global_vars_initialized |
| bool | _nodal_vars_initialized |
| bool | _use_mesh_dimension_instead_of_spatial_dimension |
| Point | _coordinate_offset |
| bool | _single_precision |
| const Parallel::Communicator & | _communicator |
Private Member Functions | |
| void | write_var_names (ExodusVarType type, std::vector< std::string > &names) |
| void | check_existing_vars (ExodusVarType type, std::vector< std::string > &names, std::vector< std::string > &names_from_file) |
| void | read_var_names_impl (const char *var_type, int &count, std::vector< std::string > &result) |
| void | write_var_names_impl (const char *var_type, int &count, std::vector< std::string > &names) |
This is the ExodusII_IO_Helper class. This class hides the implementation details of interfacing with the Exodus binary format.
Definition at line 68 of file exodusII_io_helper.h.
Wraps calls to exII::ex_get_var_names() and exII::ex_get_var_param(). The enumeration controls whether nodal, elemental, or global variable names are read and which class members are filled in. NODAL: num_nodal_vars nodal_var_names ELEMENTAL: num_elem_vars elem_var_names GLOBAL: num_global_vars global_var_names
Definition at line 533 of file exodusII_io_helper.h.
| libMesh::ExodusII_IO_Helper::ExodusII_IO_Helper | ( | const ParallelObject & | parent, |
| bool | v = false, |
||
| bool | run_only_on_proc0 = true, |
||
| bool | single_precision = false |
||
| ) |
Constructor. Automatically initializes all the private members of the class. Also allows you to set the verbosity level to v=true (on) or v=false (off). The second argument, if true, tells the class to only perform its actions if running on processor zero. If you initialize this to false, the writing methods will run on all processors instead.
Definition at line 232 of file exodusII_io_helper.C.
References elem_type, and title.
: ParallelObject(parent), ex_id(0), ex_err(0), num_dim(0), num_global_vars(0), num_nodes(0), num_elem(0), num_elem_blk(0), num_node_sets(0), num_side_sets(0), num_elem_this_blk(0), num_nodes_per_elem(0), num_attr(0), num_elem_all_sidesets(0), num_time_steps(0), num_nodal_vars(0), num_elem_vars(0), verbose(v), opened_for_writing(false), opened_for_reading(false), _run_only_on_proc0(run_only_on_proc0), _elem_vars_initialized(false), _global_vars_initialized(false), _nodal_vars_initialized(false), _use_mesh_dimension_instead_of_spatial_dimension(false), _single_precision(single_precision) { title.resize(MAX_LINE_LENGTH+1); elem_type.resize(MAX_STR_LENGTH); }
| libMesh::ExodusII_IO_Helper::~ExodusII_IO_Helper | ( | ) | [virtual] |
| void libMesh::ExodusII_IO_Helper::check_existing_vars | ( | ExodusVarType | type, |
| std::vector< std::string > & | names, | ||
| std::vector< std::string > & | names_from_file | ||
| ) | [private] |
When appending: during initialization, check that variable names in the file match those you attempt to initialize with.
Definition at line 1572 of file exodusII_io_helper.C.
References libMesh::err, libMesh::out, and read_var_names().
Referenced by initialize_element_variables(), initialize_global_variables(), and initialize_nodal_variables().
{
// There may already be global variables in the file (for example,
// if we're appending) and in that case, we
// 1.) Cannot initialize them again.
// 2.) Should check to be sure that the global variable names are the same.
// Fills up names_from_file for us
this->read_var_names(type);
// Both the names of the global variables and their order must match
if (names_from_file != names)
{
libMesh::err << "Error! The Exodus file already contains the variables:" << std::endl;
for (unsigned i=0; i<names_from_file.size(); ++i)
libMesh::out << names_from_file[i] << std::endl;
libMesh::err << "And you asked to write:" << std::endl;
for (unsigned i=0; i<names.size(); ++i)
libMesh::out << names[i] << std::endl;
libmesh_error_msg("Cannot overwrite existing variables in Exodus II file.");
}
}
| void libMesh::ExodusII_IO_Helper::close | ( | ) |
Closes the ExodusII mesh file.
Definition at line 697 of file exodusII_io_helper.C.
References _run_only_on_proc0, ex_err, ex_id, message(), and libMesh::ParallelObject::processor_id().
Referenced by libMesh::ExodusII_IO::~ExodusII_IO(), and libMesh::Nemesis_IO_Helper::~Nemesis_IO_Helper().
{
// Always call close on processor 0.
// If we're running on multiple processors, i.e. as one of several Nemesis files,
// we call close on all processors...
if ((this->processor_id() == 0) || (!_run_only_on_proc0))
{
ex_err = exII::ex_close(ex_id);
EX_CHECK_ERR(ex_err, "Error closing Exodus file.");
message("Exodus file closed successfully.");
}
}
| const Parallel::Communicator& libMesh::ParallelObject::comm | ( | ) | const [inline, inherited] |
Parallel::Communicator object used by this mesh. Definition at line 86 of file parallel_object.h.
References libMesh::ParallelObject::_communicator.
Referenced by libMesh::__libmesh_petsc_diff_solver_monitor(), libMesh::__libmesh_petsc_diff_solver_residual(), libMesh::__libmesh_petsc_snes_residual(), libMesh::MeshRefinement::_coarsen_elements(), libMesh::ExactSolution::_compute_error(), libMesh::MetisPartitioner::_do_partition(), libMesh::ParmetisPartitioner::_do_repartition(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult_add(), libMesh::EquationSystems::_read_impl(), libMesh::MeshRefinement::_refine_elements(), libMesh::ImplicitSystem::add_matrix(), libMesh::System::add_vector(), libMesh::UnstructuredMesh::all_second_order(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assemble_qoi(), libMesh::MeshCommunication::assign_global_indices(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::DofMap::attach_matrix(), libMesh::MeshTools::bounding_box(), libMesh::MeshBase::cache_elem_dims(), libMesh::System::calculate_norm(), libMesh::MeshRefinement::coarsen_elements(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Problem_Interface::computeF(), libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), DMlibMeshFunction(), DMlibMeshSetSystem_libMesh(), libMesh::MeshRefinement::eliminate_unrefined_patches(), libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::ExactErrorEstimator::estimate_error(), libMesh::MeshRefinement::flag_elements_by_elem_fraction(), libMesh::MeshRefinement::flag_elements_by_error_fraction(), libMesh::MeshRefinement::flag_elements_by_nelem_target(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::ImplicitSystem::get_linear_solver(), libMesh::LocationMap< T >::init(), libMesh::TimeSolver::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::EigenSystem::init_data(), libMesh::EigenSystem::init_matrices(), libMesh::ParmetisPartitioner::initialize(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::ParallelMesh::libmesh_assert_valid_parallel_flags(), libMesh::MeshTools::libmesh_assert_valid_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::MeshRefinement::limit_level_mismatch_at_edge(), libMesh::MeshRefinement::limit_level_mismatch_at_node(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshCommunication::make_elems_parallel_consistent(), libMesh::MeshRefinement::make_flags_parallel_consistent(), libMesh::MeshCommunication::make_node_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_nodes_parallel_consistent(), libMesh::MeshRefinement::make_refinement_compatible(), libMesh::FEMSystem::mesh_position_set(), libMesh::MeshSerializer::MeshSerializer(), libMesh::ParallelMesh::n_active_elem(), libMesh::MeshTools::n_active_levels(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::CondensedEigenSystem::n_global_non_condensed_dofs(), libMesh::MeshTools::n_levels(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::MeshTools::n_p_levels(), libMesh::ParallelMesh::parallel_max_elem_id(), libMesh::ParallelMesh::parallel_max_node_id(), libMesh::ParallelMesh::parallel_n_elem(), libMesh::ParallelMesh::parallel_n_nodes(), libMesh::Partitioner::partition(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::petsc_auto_fieldsplit(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshBase::prepare_for_use(), libMesh::System::project_vector(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::MeshBase::recalculate_n_partitions(), libMesh::MeshRefinement::refine_and_coarsen_elements(), libMesh::MeshRefinement::refine_elements(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::MeshBase::subdomain_ids(), libMesh::BoundaryInfo::sync(), libMesh::Parallel::sync_element_data_by_parent_id(), libMesh::Parallel::sync_node_data_by_element_id(), libMesh::MeshRefinement::test_level_one(), libMesh::MeshRefinement::test_unflagged(), libMesh::MeshTools::total_weight(), libMesh::NameBasedIO::write(), libMesh::CheckpointIO::write(), libMesh::XdrIO::write(), libMesh::LegacyXdrIO::write_mesh(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), and libMesh::DivaIO::write_stream().
{ return _communicator; }
| void libMesh::ExodusII_IO_Helper::create | ( | std::string | filename | ) | [virtual] |
Opens an ExodusII mesh file named filename for writing.
Reimplemented in libMesh::Nemesis_IO_Helper.
Definition at line 984 of file exodusII_io_helper.C.
References _run_only_on_proc0, _single_precision, current_filename, ex_id, std::min(), opened_for_writing, libMesh::out, libMesh::ParallelObject::processor_id(), libMesh::Real, and verbose.
Referenced by libMesh::ExodusII_IO::write(), and libMesh::ExodusII_IO::write_nodal_data_common().
{
// If we're processor 0, always create the file.
// If we running on all procs, e.g. as one of several Nemesis files, also
// call create there.
if ((this->processor_id() == 0) || (!_run_only_on_proc0))
{
int
comp_ws = 0,
io_ws = 0;
if(_single_precision)
{
comp_ws = cast_int<int>(sizeof(float));
io_ws = cast_int<int>(sizeof(float));
}
// Fall back on double precision when necessary since ExodusII
// doesn't seem to support long double
else
{
comp_ws = cast_int<int>
(std::min(sizeof(Real), sizeof(double)));
io_ws = cast_int<int>
(std::min(sizeof(Real), sizeof(double)));
}
ex_id = exII::ex_create(filename.c_str(), EX_CLOBBER, &comp_ws, &io_ws);
EX_CHECK_ERR(ex_id, "Error creating ExodusII mesh file.");
if (verbose)
libMesh::out << "File created successfully." << std::endl;
}
opened_for_writing = true;
current_filename = filename;
}
| int libMesh::ExodusII_IO_Helper::get_block_id | ( | int | index | ) |
Get the block number for the given block index.
Definition at line 446 of file exodusII_io_helper.C.
References block_ids.
Referenced by libMesh::ExodusII_IO::read(), and write_element_values().
| std::string libMesh::ExodusII_IO_Helper::get_block_name | ( | int | index | ) |
Get the block name for the given block index if supplied in the mesh file. Otherwise an empty string is returned.
Definition at line 455 of file exodusII_io_helper.C.
References block_ids, and id_to_block_names.
Referenced by libMesh::ExodusII_IO::read().
{
libmesh_assert_less (static_cast<unsigned int>(index), block_ids.size());
return id_to_block_names[block_ids[index]];
}
| std::vector< std::string > libMesh::ExodusII_IO_Helper::get_complex_names | ( | const std::vector< std::string > & | names | ) | const |
Returns a vector with three copies of each element in the provided name vector, starting with r_, i_ and a_ respectively.
Definition at line 1784 of file exodusII_io_helper.C.
Referenced by libMesh::ExodusII_IO::write_element_data(), libMesh::Nemesis_IO::write_global_data(), libMesh::ExodusII_IO::write_global_data(), libMesh::Nemesis_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), and libMesh::ExodusII_IO::write_nodal_data_discontinuous().
{
std::vector<std::string>::const_iterator names_it = names.begin();
std::vector<std::string>::const_iterator names_end = names.end();
std::vector<std::string> complex_names;
// This will loop over all names and create new "complex" names
// (i.e. names that start with r_, i_ or a_
for (; names_it != names_end; ++names_it)
{
// std::cout << "VARIABLE: " << *names_it << std::endl;
std::stringstream name_real, name_imag, name_abs;
name_real << "r_" << *names_it;
name_imag << "i_" << *names_it;
name_abs << "a_" << *names_it;
complex_names.push_back(name_real.str());
complex_names.push_back(name_imag.str());
complex_names.push_back(name_abs.str());
}
return complex_names;
}
| const char * libMesh::ExodusII_IO_Helper::get_elem_type | ( | ) | const |
HEX27. Definition at line 275 of file exodusII_io_helper.C.
References elem_type.
Referenced by libMesh::ExodusII_IO::read().
{
return &elem_type[0];
}
| int libMesh::ExodusII_IO_Helper::get_node_set_id | ( | int | index | ) |
Get the node set id for the given node set index.
Definition at line 482 of file exodusII_io_helper.C.
References nodeset_ids.
{
libmesh_assert_less (static_cast<unsigned int>(index), nodeset_ids.size());
return nodeset_ids[index];
}
| std::string libMesh::ExodusII_IO_Helper::get_node_set_name | ( | int | index | ) |
Get the node set name for the given node set index if supplied in the mesh file. Otherwise an empty string is returned.
Definition at line 491 of file exodusII_io_helper.C.
References id_to_ns_names, and nodeset_ids.
Referenced by libMesh::ExodusII_IO::read().
{
libmesh_assert_less (static_cast<unsigned int>(index), nodeset_ids.size());
return id_to_ns_names[nodeset_ids[index]];
}
| int libMesh::ExodusII_IO_Helper::get_side_set_id | ( | int | index | ) |
Get the side set id for the given side set index.
Definition at line 464 of file exodusII_io_helper.C.
References ss_ids.
Referenced by libMesh::ExodusII_IO::read().
| std::string libMesh::ExodusII_IO_Helper::get_side_set_name | ( | int | index | ) |
Get the side set name for the given side set index if supplied in the mesh file. Otherwise an empty string is returned.
Definition at line 473 of file exodusII_io_helper.C.
References id_to_ss_names, and ss_ids.
Referenced by libMesh::ExodusII_IO::read().
{
libmesh_assert_less (static_cast<unsigned int>(index), ss_ids.size());
return id_to_ss_names[ss_ids[index]];
}
| void libMesh::ExodusII_IO_Helper::initialize | ( | std::string | title, |
| const MeshBase & | mesh, | ||
| bool | use_discontinuous = false |
||
| ) | [virtual] |
Initializes the Exodus file.
Reimplemented in libMesh::Nemesis_IO_Helper.
Definition at line 1024 of file exodusII_io_helper.C.
References _run_only_on_proc0, _use_mesh_dimension_instead_of_spatial_dimension, libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::BoundaryInfo::build_node_boundary_ids(), libMesh::BoundaryInfo::build_side_boundary_ids(), end, libMesh::err, ex_err, ex_id, libMesh::MeshBase::get_boundary_info(), libMesh::DofObject::id(), libMesh::MeshBase::mesh_dimension(), libMesh::MeshBase::n_active_elem(), libMesh::MeshBase::n_elem(), libMesh::MeshBase::n_nodes(), num_dim, num_elem, num_elem_blk, num_node_sets, num_nodes, num_side_sets, libMesh::ParallelObject::processor_id(), libMesh::MeshBase::spatial_dimension(), and libMesh::Elem::subdomain_id().
Referenced by libMesh::ExodusII_IO::write(), and libMesh::ExodusII_IO::write_nodal_data_common().
{
// n_active_elem() is a parallel_only function
unsigned int n_active_elem = mesh.n_active_elem();
if ((_run_only_on_proc0) && (this->processor_id() != 0))
return;
if (_use_mesh_dimension_instead_of_spatial_dimension)
num_dim = mesh.mesh_dimension();
else
num_dim = mesh.spatial_dimension();
num_elem = mesh.n_elem();
if (!use_discontinuous)
num_nodes = mesh.n_nodes();
else
{
MeshBase::const_element_iterator it = mesh.active_elements_begin();
const MeshBase::const_element_iterator end = mesh.active_elements_end();
for (; it!=end; ++it)
num_nodes += (*it)->n_nodes();
}
std::vector<boundary_id_type> unique_side_boundaries;
std::vector<boundary_id_type> unique_node_boundaries;
mesh.get_boundary_info().build_side_boundary_ids(unique_side_boundaries);
mesh.get_boundary_info().build_node_boundary_ids(unique_node_boundaries);
num_side_sets = cast_int<int>(unique_side_boundaries.size());
num_node_sets = cast_int<int>(unique_node_boundaries.size());
//loop through element and map between block and element vector
std::map<subdomain_id_type, std::vector<unsigned int> > subdomain_map;
MeshBase::const_element_iterator it = mesh.active_elements_begin();
const MeshBase::const_element_iterator end = mesh.active_elements_end();
for (; it!=end; ++it)
{
const Elem * elem = *it;
subdomain_id_type cur_subdomain = elem->subdomain_id();
subdomain_map[cur_subdomain].push_back(elem->id());
}
num_elem_blk = cast_int<int>(subdomain_map.size());
if (str_title.size() > MAX_LINE_LENGTH)
{
libMesh::err << "Warning, Exodus files cannot have titles longer than "
<< MAX_LINE_LENGTH
<< " characters. Your title will be truncated."
<< std::endl;
str_title.resize(MAX_LINE_LENGTH);
}
ex_err = exII::ex_put_init(ex_id,
str_title.c_str(),
num_dim,
num_nodes,
n_active_elem,
num_elem_blk,
num_node_sets,
num_side_sets);
EX_CHECK_ERR(ex_err, "Error initializing new Exodus file.");
}
| void libMesh::ExodusII_IO_Helper::initialize_element_variables | ( | std::vector< std::string > | names | ) |
Sets up the nodal variables
Definition at line 1471 of file exodusII_io_helper.C.
References _elem_vars_initialized, _run_only_on_proc0, check_existing_vars(), elem_var_names, ELEMENTAL, ex_err, ex_id, num_elem_blk, num_elem_vars, libMesh::ParallelObject::processor_id(), and write_var_names().
Referenced by libMesh::ExodusII_IO::write_element_data().
{
if ((_run_only_on_proc0) && (this->processor_id() != 0))
return;
// Quick return if there are no element variables to write
if (names.size() == 0)
return;
// Quick return if we have already called this function
if (_elem_vars_initialized)
return;
// Be sure that variables in the file match what we are asking for
if (num_elem_vars > 0)
{
this->check_existing_vars(ELEMENTAL, names, this->elem_var_names);
return;
}
// Set the flag so we can skip this stuff on subsequent calls to
// initialize_element_variables()
_elem_vars_initialized = true;
this->write_var_names(ELEMENTAL, names);
// Form the element variable truth table and send to Exodus.
// This tells which variables are written to which blocks,
// and can dramatically speed up writing element variables
//
// We really should initialize all entries in the truth table to 0
// and then loop over all subdomains, setting their entries to 1
// if a given variable exists on that subdomain. However,
// we don't have that information, and the element variables
// passed to us are padded with zeroes for the blocks where
// they aren't defined. To be consistent with that, fill
// the truth table with ones.
std::vector<int> truth_tab(num_elem_blk*num_elem_vars, 1);
ex_err = exII::ex_put_elem_var_tab(ex_id,
num_elem_blk,
num_elem_vars,
&truth_tab[0]);
EX_CHECK_ERR(ex_err, "Error writing element truth table.");
}
| void libMesh::ExodusII_IO_Helper::initialize_global_variables | ( | std::vector< std::string > | names | ) |
Sets up the global variables
Definition at line 1546 of file exodusII_io_helper.C.
References _global_vars_initialized, _run_only_on_proc0, check_existing_vars(), GLOBAL, global_var_names, num_global_vars, libMesh::ParallelObject::processor_id(), and write_var_names().
Referenced by libMesh::Nemesis_IO::write_global_data(), and libMesh::ExodusII_IO::write_global_data().
{
if ((_run_only_on_proc0) && (this->processor_id() != 0))
return;
// Quick return if there are no global variables to write
if (names.size() == 0)
return;
if (_global_vars_initialized)
return;
// Be sure that variables in the file match what we are asking for
if (num_global_vars > 0)
{
this->check_existing_vars(GLOBAL, names, this->global_var_names);
return;
}
_global_vars_initialized = true;
this->write_var_names(GLOBAL, names);
}
| void libMesh::ExodusII_IO_Helper::initialize_nodal_variables | ( | std::vector< std::string > | names | ) |
Sets up the nodal variables
Definition at line 1518 of file exodusII_io_helper.C.
References _nodal_vars_initialized, _run_only_on_proc0, check_existing_vars(), NODAL, nodal_var_names, num_nodal_vars, libMesh::ParallelObject::processor_id(), and write_var_names().
Referenced by libMesh::Nemesis_IO::write_nodal_data(), and libMesh::ExodusII_IO::write_nodal_data_common().
{
if ((_run_only_on_proc0) && (this->processor_id() != 0))
return;
// Quick return if there are no nodal variables to write
if (names.size() == 0)
return;
// Quick return if we have already called this function
if (_nodal_vars_initialized)
return;
// Be sure that variables in the file match what we are asking for
if (num_nodal_vars > 0)
{
this->check_existing_vars(NODAL, names, this->nodal_var_names);
return;
}
// Set the flag so we can skip the rest of this function on subsequent calls.
_nodal_vars_initialized = true;
this->write_var_names(NODAL, names);
}
| int libMesh::ExodusII_IO_Helper::inquire | ( | int | req_info, |
| std::string | error_msg = "" |
||
| ) |
Generic inquiry, returns the value
Definition at line 712 of file exodusII_io_helper.C.
Referenced by read_num_time_steps(), read_sideset_info(), and write_information_records().
| void libMesh::ExodusII_IO_Helper::message | ( | const std::string & | msg | ) |
Prints the message defined in msg. Can be turned off if verbosity is set to 0.
Definition at line 282 of file exodusII_io_helper.C.
References libMesh::out, and verbose.
Referenced by close(), read_block_info(), read_elem_in_block(), read_elem_num_map(), read_header(), read_node_num_map(), read_nodes(), read_nodeset(), read_nodeset_info(), read_sideset(), and read_sideset_info().
{
if (verbose) libMesh::out << msg << std::endl;
}
| void libMesh::ExodusII_IO_Helper::message | ( | const std::string & | msg, |
| int | i | ||
| ) |
Prints the message defined in msg, and appends the number i to the end of the message. Useful for printing messages in loops. Can be turned off if verbosity is set to 0.
Definition at line 289 of file exodusII_io_helper.C.
References libMesh::out, and verbose.
{
if (verbose) libMesh::out << msg << i << "." << std::endl;
}
| processor_id_type libMesh::ParallelObject::n_processors | ( | ) | const [inline, inherited] |
Definition at line 92 of file parallel_object.h.
References libMesh::ParallelObject::_communicator, and libMesh::Parallel::Communicator::size().
Referenced by libMesh::ParmetisPartitioner::_do_repartition(), libMesh::ParallelMesh::add_elem(), libMesh::ParallelMesh::add_node(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assembly(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::ParallelMesh::assign_unique_ids(), libMesh::AztecLinearSolver< T >::AztecLinearSolver(), libMesh::ParallelMesh::clear(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::UnstructuredMesh::create_pid_mesh(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::EnsightIO::EnsightIO(), libMesh::MeshBase::get_info(), libMesh::EquationSystems::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ParmetisPartitioner::initialize(), libMesh::Nemesis_IO_Helper::initialize(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::DofMap::local_variable_indices(), libMesh::MeshBase::n_active_elem_on_proc(), libMesh::MeshBase::n_elem_on_proc(), libMesh::MeshBase::n_nodes_on_proc(), libMesh::Partitioner::partition(), libMesh::MeshBase::partition(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::PetscLinearSolver< T >::PetscLinearSolver(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshTools::processor_bounding_box(), libMesh::System::project_vector(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::XdrIO::read(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::System::read_serialized_vector(), libMesh::Partitioner::repartition(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::BoundaryInfo::sync(), libMesh::ParallelMesh::update_parallel_id_counts(), libMesh::CheckpointIO::write(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().
{ return cast_int<processor_id_type>(_communicator.size()); }
| void libMesh::ExodusII_IO_Helper::open | ( | const char * | filename, |
| bool | read_only | ||
| ) |
Opens an ExodusII mesh file named filename. If read_only==true, the file will be opened with the EX_READ flag, otherwise it will be opened with the EX_WRITE flag.
Definition at line 296 of file exodusII_io_helper.C.
References current_filename, ex_id, opened_for_reading, opened_for_writing, libMesh::out, libMesh::Real, and verbose.
Referenced by libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), libMesh::Nemesis_IO::write_nodal_data(), and libMesh::ExodusII_IO::write_nodal_data_common().
{
// Version of Exodus you are using
float ex_version = 0.;
// Word size in bytes of the floating point variables used in the
// application program (0, 4, or 8)
int comp_ws = sizeof(Real);
// Word size in bytes of the floating point data as they are stored
// in the ExodusII file. "If this argument is 0, the word size of the
// floating point data already stored in the file is returned"
int io_ws = 0;
ex_id = exII::ex_open(filename,
read_only ? EX_READ : EX_WRITE,
&comp_ws,
&io_ws,
&ex_version);
std::string err_msg = std::string("Error opening ExodusII mesh file: ") + std::string(filename);
EX_CHECK_ERR(ex_id, err_msg);
if (verbose) libMesh::out << "File opened successfully." << std::endl;
if (read_only)
opened_for_reading = true;
else
opened_for_writing = true;
current_filename = std::string(filename);
}
Prints the ExodusII mesh file header, which includes the mesh title, the number of nodes, number of elements, mesh dimension, number of sidesets, and number of nodesets.
Definition at line 360 of file exodusII_io_helper.C.
References num_dim, num_elem, num_elem_blk, num_node_sets, num_nodes, num_side_sets, libMesh::out, title, and verbose.
Referenced by libMesh::Nemesis_IO::read(), and libMesh::ExodusII_IO::read().
{
if (verbose)
libMesh::out << "Title: \t" << &title[0] << std::endl
<< "Mesh Dimension: \t" << num_dim << std::endl
<< "Number of Nodes: \t" << num_nodes << std::endl
<< "Number of elements: \t" << num_elem << std::endl
<< "Number of elt blocks: \t" << num_elem_blk << std::endl
<< "Number of node sets: \t" << num_node_sets << std::endl
<< "Number of side sets: \t" << num_side_sets << std::endl;
}
| void libMesh::ExodusII_IO_Helper::print_nodes | ( | std::ostream & | out = libMesh::out | ) |
| processor_id_type libMesh::ParallelObject::processor_id | ( | ) | const [inline, inherited] |
Definition at line 98 of file parallel_object.h.
References libMesh::ParallelObject::_communicator, and libMesh::Parallel::Communicator::rank().
Referenced by libMesh::MetisPartitioner::_do_partition(), libMesh::EquationSystems::_read_impl(), libMesh::SerialMesh::active_local_elements_begin(), libMesh::ParallelMesh::active_local_elements_begin(), libMesh::SerialMesh::active_local_elements_end(), libMesh::ParallelMesh::active_local_elements_end(), libMesh::SerialMesh::active_local_subdomain_elements_begin(), libMesh::ParallelMesh::active_local_subdomain_elements_begin(), libMesh::SerialMesh::active_local_subdomain_elements_end(), libMesh::ParallelMesh::active_local_subdomain_elements_end(), libMesh::SerialMesh::active_not_local_elements_begin(), libMesh::ParallelMesh::active_not_local_elements_begin(), libMesh::SerialMesh::active_not_local_elements_end(), libMesh::ParallelMesh::active_not_local_elements_end(), libMesh::ParallelMesh::add_elem(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::ParallelMesh::add_node(), libMesh::UnstructuredMesh::all_second_order(), libMesh::FEMSystem::assembly(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::ParallelMesh::assign_unique_ids(), libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::ParmetisPartitioner::build_graph(), libMesh::InfElemBuilder::build_inf_elem(), libMesh::DofMap::build_sparsity(), libMesh::ParallelMesh::clear(), close(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::compute_communication_map_parameters(), libMesh::Nemesis_IO_Helper::compute_internal_and_border_elems_and_internal_nodes(), libMesh::Nemesis_IO_Helper::compute_node_communication_maps(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), create(), libMesh::ParallelMesh::delete_elem(), libMesh::ParallelMesh::delete_node(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::DofMap::end_dof(), libMesh::DofMap::end_old_dof(), libMesh::EnsightIO::EnsightIO(), libMesh::SerialMesh::facelocal_elements_begin(), libMesh::ParallelMesh::facelocal_elements_begin(), libMesh::SerialMesh::facelocal_elements_end(), libMesh::ParallelMesh::facelocal_elements_end(), libMesh::MeshFunction::find_element(), libMesh::UnstructuredMesh::find_neighbors(), libMesh::DofMap::first_dof(), libMesh::DofMap::first_old_dof(), libMesh::Nemesis_IO_Helper::get_cmap_params(), libMesh::Nemesis_IO_Helper::get_eb_info_global(), libMesh::Nemesis_IO_Helper::get_elem_cmap(), libMesh::Nemesis_IO_Helper::get_elem_map(), libMesh::MeshBase::get_info(), libMesh::Nemesis_IO_Helper::get_init_global(), libMesh::Nemesis_IO_Helper::get_init_info(), libMesh::Nemesis_IO_Helper::get_loadbal_param(), libMesh::Nemesis_IO_Helper::get_node_cmap(), libMesh::Nemesis_IO_Helper::get_node_map(), libMesh::Nemesis_IO_Helper::get_ns_param_global(), libMesh::Nemesis_IO_Helper::get_ss_param_global(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ParmetisPartitioner::initialize(), initialize(), initialize_element_variables(), initialize_global_variables(), initialize_nodal_variables(), libMesh::ParallelMesh::insert_elem(), libMesh::SparsityPattern::Build::join(), libMesh::DofMap::last_dof(), libMesh::MeshTools::libmesh_assert_valid_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_procids< Node >(), libMesh::SerialMesh::local_elements_begin(), libMesh::ParallelMesh::local_elements_begin(), libMesh::SerialMesh::local_elements_end(), libMesh::ParallelMesh::local_elements_end(), libMesh::SerialMesh::local_level_elements_begin(), libMesh::ParallelMesh::local_level_elements_begin(), libMesh::SerialMesh::local_level_elements_end(), libMesh::ParallelMesh::local_level_elements_end(), libMesh::SerialMesh::local_nodes_begin(), libMesh::ParallelMesh::local_nodes_begin(), libMesh::SerialMesh::local_nodes_end(), libMesh::ParallelMesh::local_nodes_end(), libMesh::SerialMesh::local_not_level_elements_begin(), libMesh::ParallelMesh::local_not_level_elements_begin(), libMesh::SerialMesh::local_not_level_elements_end(), libMesh::ParallelMesh::local_not_level_elements_end(), libMesh::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_local_elem(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::DofMap::n_local_dofs(), libMesh::System::n_local_dofs(), libMesh::MeshBase::n_local_elem(), libMesh::MeshBase::n_local_nodes(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::SerialMesh::not_local_elements_begin(), libMesh::ParallelMesh::not_local_elements_begin(), libMesh::SerialMesh::not_local_elements_end(), libMesh::ParallelMesh::not_local_elements_end(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::SparsityPattern::Build::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::ParallelMesh::ParallelMesh(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::System::project_vector(), libMesh::Nemesis_IO_Helper::put_cmap_params(), libMesh::Nemesis_IO_Helper::put_elem_cmap(), libMesh::Nemesis_IO_Helper::put_elem_map(), libMesh::Nemesis_IO_Helper::put_loadbal_param(), libMesh::Nemesis_IO_Helper::put_node_cmap(), libMesh::Nemesis_IO_Helper::put_node_map(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::XdrIO::read(), read_elem_num_map(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), read_node_num_map(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::System::read_serialized_data(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::System::read_serialized_vectors(), libMesh::MeshData::read_xdr(), libMesh::SerialMesh::semilocal_elements_begin(), libMesh::ParallelMesh::semilocal_elements_begin(), libMesh::SerialMesh::semilocal_elements_end(), libMesh::ParallelMesh::semilocal_elements_end(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::BoundaryInfo::sync(), libMesh::MeshTools::total_weight(), libMesh::ParallelMesh::update_parallel_id_counts(), libMesh::MeshTools::weight(), libMesh::NameBasedIO::write(), libMesh::ExodusII_IO::write(), libMesh::CheckpointIO::write(), libMesh::XdrIO::write(), libMesh::EquationSystems::write(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::ExodusII_IO::write_element_data(), write_element_values(), write_elements(), libMesh::ExodusII_IO::write_global_data(), write_global_values(), libMesh::System::write_header(), libMesh::ExodusII_IO::write_information_records(), write_information_records(), write_nodal_coordinates(), libMesh::UCDIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), write_nodal_values(), write_nodesets(), libMesh::Nemesis_IO_Helper::write_nodesets(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bc_names(), libMesh::XdrIO::write_serialized_bcs(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::System::write_serialized_data(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), libMesh::System::write_serialized_vector(), libMesh::System::write_serialized_vectors(), write_sidesets(), libMesh::Nemesis_IO_Helper::write_sidesets(), libMesh::ExodusII_IO::write_timestep(), and write_timestep().
{ return cast_int<processor_id_type>(_communicator.rank()); }
Reads information for all of the blocks in the ExodusII mesh file.
Definition at line 419 of file exodusII_io_helper.C.
References block_ids, ex_err, ex_id, id_to_block_names, message(), and num_elem_blk.
Referenced by libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), and libMesh::ExodusII_IO::write_nodal_data_common().
{
block_ids.resize(num_elem_blk);
// Get all element block IDs.
ex_err = exII::ex_get_elem_blk_ids(ex_id,
block_ids.empty() ? NULL : &block_ids[0]);
// Usually, there is only one
// block since there is only
// one type of element.
// However, there could be more.
EX_CHECK_ERR(ex_err, "Error getting block IDs.");
message("All block IDs retrieved successfully.");
char name_buffer[MAX_STR_LENGTH+1];
for (int i=0; i<num_elem_blk; ++i)
{
ex_err = exII::ex_get_name(ex_id, exII::EX_ELEM_BLOCK,
block_ids[i], name_buffer);
EX_CHECK_ERR(ex_err, "Error getting block name.");
id_to_block_names[block_ids[i]] = name_buffer;
}
message("All block names retrieved successfully.");
}
| void libMesh::ExodusII_IO_Helper::read_elem_in_block | ( | int | block | ) |
Reads all of the element connectivity for block block in the ExodusII mesh file.
Definition at line 501 of file exodusII_io_helper.C.
References block_ids, connect, elem_type, ex_err, ex_id, message(), num_attr, num_elem_this_blk, num_nodes_per_elem, libMesh::out, and verbose.
Referenced by libMesh::Nemesis_IO::read(), and libMesh::ExodusII_IO::read().
{
libmesh_assert_less (static_cast<unsigned int>(block), block_ids.size());
ex_err = exII::ex_get_elem_block(ex_id,
block_ids[block],
&elem_type[0],
&num_elem_this_blk,
&num_nodes_per_elem,
&num_attr);
if (verbose)
libMesh::out << "Reading a block of " << num_elem_this_blk
<< " " << &elem_type[0] << "(s)"
<< " having " << num_nodes_per_elem
<< " nodes per element." << std::endl;
EX_CHECK_ERR(ex_err, "Error getting block info.");
message("Info retrieved successfully for block: ", block);
// Read in the connectivity of the elements of this block,
// watching out for the case where we actually have no
// elements in this block (possible with parallel files)
connect.resize(num_nodes_per_elem*num_elem_this_blk);
if (!connect.empty())
{
ex_err = exII::ex_get_elem_conn(ex_id,
block_ids[block],
&connect[0]);
EX_CHECK_ERR(ex_err, "Error reading block connectivity.");
message("Connectivity retrieved successfully for block: ", block);
}
}
Reads the optional node_num_map from the ExodusII mesh file.
Definition at line 541 of file exodusII_io_helper.C.
References elem_num_map, ex_err, ex_id, message(), std::min(), num_elem, libMesh::out, libMesh::ParallelObject::processor_id(), and verbose.
Referenced by libMesh::Nemesis_IO::read(), and libMesh::ExodusII_IO::read().
{
elem_num_map.resize(num_elem);
ex_err = exII::ex_get_elem_num_map (ex_id,
elem_num_map.empty() ? NULL : &elem_num_map[0]);
EX_CHECK_ERR(ex_err, "Error retrieving element number map.");
message("Element numbering map retrieved successfully.");
if (verbose)
{
libMesh::out << "[" << this->processor_id() << "] elem_num_map[i] = ";
for (unsigned int i=0; i<static_cast<unsigned int>(std::min(10, num_elem-1)); ++i)
libMesh::out << elem_num_map[i] << ", ";
libMesh::out << "... " << elem_num_map.back() << std::endl;
}
}
| void libMesh::ExodusII_IO_Helper::read_elemental_var_values | ( | std::string | elemental_var_name, |
| int | time_step | ||
| ) |
Reads elemental values for the variable 'elemental_var_name' at the specified timestep into the 'elem_var_values' array.
Definition at line 910 of file exodusII_io_helper.C.
References block_ids, elem_num_map, elem_var_names, elem_var_values, ELEMENTAL, libMesh::err, ex_err, ex_id, num_elem, num_elem_blk, num_elem_this_blk, and read_var_names().
Referenced by libMesh::ExodusII_IO::copy_elemental_solution().
{
// There is no way to get the whole elemental field from the
// exodus file, so we have to go block by block.
elem_var_values.resize(num_elem);
this->read_var_names(ELEMENTAL);
// See if we can find the variable we are looking for
unsigned var_index = 0;
bool found = false;
// Do a linear search for nodal_var_name in nodal_var_names
for (; var_index<elem_var_names.size(); ++var_index)
{
found = (elem_var_names[var_index] == elemental_var_name);
if (found)
break;
}
if (!found)
{
libMesh::err << "Available variables: " << std::endl;
for (unsigned i=0; i<elem_var_names.size(); ++i)
libMesh::err << elem_var_names[i] << std::endl;
libmesh_error_msg("Unable to locate variable named: " << elemental_var_name);
}
// Sequential index which we can use to look up the element ID in the elem_num_map.
unsigned ex_el_num = 0;
for (unsigned i=0; i<static_cast<unsigned>(num_elem_blk); i++)
{
ex_err = exII::ex_get_elem_block(ex_id,
block_ids[i],
NULL,
&num_elem_this_blk,
NULL,
NULL);
EX_CHECK_ERR(ex_err, "Error getting number of elements in block.");
std::vector<Real> block_elem_var_values(num_elem);
ex_err = exII::ex_get_elem_var(ex_id,
time_step,
var_index+1,
block_ids[i],
num_elem_this_blk,
&block_elem_var_values[0]);
EX_CHECK_ERR(ex_err, "Error getting elemental values.");
for (unsigned j=0; j<static_cast<unsigned>(num_elem_this_blk); j++)
{
// Use the elem_num_map to obtain the ID of this element in the Exodus file,
// and remember to subtract 1 since libmesh is zero-based and Exodus is 1-based.
unsigned mapped_elem_id = this->elem_num_map[ex_el_num] - 1;
// Make sure we can actually write into this location in the elem_var_values vector
if (mapped_elem_id >= elem_var_values.size())
libmesh_error_msg("Error reading elemental variable values in Exodus!");
// Write into the mapped_elem_id entry of the
// elem_var_values vector.
elem_var_values[mapped_elem_id] = block_elem_var_values[j];
// Go to the next sequential element ID.
ex_el_num++;
}
}
}
Reads an ExodusII mesh file header.
Definition at line 330 of file exodusII_io_helper.C.
References ex_err, ex_id, message(), num_dim, num_elem, num_elem_blk, num_elem_vars, num_global_vars, num_nodal_vars, num_node_sets, num_nodes, num_side_sets, read_num_time_steps(), and title.
Referenced by libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), libMesh::Nemesis_IO::write_nodal_data(), and libMesh::ExodusII_IO::write_nodal_data_common().
{
ex_err = exII::ex_get_init(ex_id,
&title[0],
&num_dim,
&num_nodes,
&num_elem,
&num_elem_blk,
&num_node_sets,
&num_side_sets);
EX_CHECK_ERR(ex_err, "Error retrieving header info.");
this->read_num_time_steps();
ex_err = exII::ex_get_var_param(ex_id, "n", &num_nodal_vars);
EX_CHECK_ERR(ex_err, "Error reading number of nodal variables.");
ex_err = exII::ex_get_var_param(ex_id, "e", &num_elem_vars);
EX_CHECK_ERR(ex_err, "Error reading number of elemental variables.");
ex_err = exII::ex_get_var_param(ex_id, "g", &num_global_vars);
EX_CHECK_ERR(ex_err, "Error reading number of global variables.");
message("Exodus header info retrieved successfully.");
}
| void libMesh::ExodusII_IO_Helper::read_nodal_var_values | ( | std::string | nodal_var_name, |
| int | time_step | ||
| ) |
Reads the nodal values for the variable 'nodal_var_name' at the specified time into the 'nodal_var_values' array.
Definition at line 754 of file exodusII_io_helper.C.
References libMesh::err, ex_err, ex_id, NODAL, nodal_var_names, nodal_var_values, num_nodes, and read_var_names().
Referenced by libMesh::ExodusII_IO::copy_nodal_solution().
{
// Read the nodal variable names from file, so we can see if we have the one we're looking for
this->read_var_names(NODAL);
// See if we can find the variable we are looking for
unsigned int var_index = 0;
bool found = false;
// Do a linear search for nodal_var_name in nodal_var_names
for (; var_index<nodal_var_names.size(); ++var_index)
{
found = (nodal_var_names[var_index] == nodal_var_name);
if (found)
break;
}
if (!found)
{
libMesh::err << "Available variables: " << std::endl;
for (unsigned int i=0; i<nodal_var_names.size(); ++i)
libMesh::err << nodal_var_names[i] << std::endl;
libmesh_error_msg("Unable to locate variable named: " << nodal_var_name);
}
// Allocate enough space to store the nodal variable values
nodal_var_values.resize(num_nodes);
// Call the Exodus API to read the nodal variable values
ex_err = exII::ex_get_nodal_var(ex_id,
time_step,
var_index+1,
num_nodes,
&nodal_var_values[0]);
EX_CHECK_ERR(ex_err, "Error reading nodal variable values!");
}
Reads the optional node_num_map from the ExodusII mesh file.
Definition at line 391 of file exodusII_io_helper.C.
References ex_err, ex_id, message(), std::min(), node_num_map, num_nodes, libMesh::out, libMesh::ParallelObject::processor_id(), and verbose.
Referenced by libMesh::Nemesis_IO::read(), and libMesh::ExodusII_IO::read().
{
node_num_map.resize(num_nodes);
ex_err = exII::ex_get_node_num_map (ex_id,
node_num_map.empty() ? NULL : &node_num_map[0]);
EX_CHECK_ERR(ex_err, "Error retrieving nodal number map.");
message("Nodal numbering map retrieved successfully.");
if (verbose)
{
libMesh::out << "[" << this->processor_id() << "] node_num_map[i] = ";
for (unsigned int i=0; i<static_cast<unsigned int>(std::min(10, num_nodes-1)); ++i)
libMesh::out << node_num_map[i] << ", ";
libMesh::out << "... " << node_num_map.back() << std::endl;
}
}
Reads the nodal data (x,y,z coordinates) from the ExodusII mesh file.
Definition at line 374 of file exodusII_io_helper.C.
References ex_err, ex_id, message(), num_nodes, x, y, and z.
Referenced by libMesh::Nemesis_IO::read(), and libMesh::ExodusII_IO::read().
| void libMesh::ExodusII_IO_Helper::read_nodeset | ( | int | id | ) |
Reads information about nodeset id and inserts it into the global nodeset array at the position offset.
Definition at line 667 of file exodusII_io_helper.C.
References ex_err, ex_id, message(), node_list, nodeset_ids, num_node_df_per_set, and num_nodes_per_set.
Referenced by libMesh::Nemesis_IO::read(), and libMesh::ExodusII_IO::read().
{
libmesh_assert_less (static_cast<unsigned int>(id), nodeset_ids.size());
libmesh_assert_less (static_cast<unsigned int>(id), num_nodes_per_set.size());
libmesh_assert_less (static_cast<unsigned int>(id), num_node_df_per_set.size());
ex_err = exII::ex_get_node_set_param(ex_id,
nodeset_ids[id],
&num_nodes_per_set[id],
&num_node_df_per_set[id]);
EX_CHECK_ERR(ex_err, "Error retrieving nodeset parameters.");
message("Parameters retrieved successfully for nodeset: ", id);
node_list.resize(num_nodes_per_set[id]);
// Don't call ex_get_node_set unless there are actually nodes there to get.
// Exodus prints an annoying warning message in DEBUG mode otherwise...
if (num_nodes_per_set[id] > 0)
{
ex_err = exII::ex_get_node_set(ex_id,
nodeset_ids[id],
&node_list[0]);
EX_CHECK_ERR(ex_err, "Error retrieving nodeset data.");
message("Data retrieved successfully for nodeset: ", id);
}
}
Reads information about all of the nodesets in the ExodusII mesh file.
Definition at line 596 of file exodusII_io_helper.C.
References ex_err, ex_id, id_to_ns_names, message(), nodeset_ids, num_node_df_per_set, num_node_sets, and num_nodes_per_set.
Referenced by libMesh::Nemesis_IO::read(), and libMesh::ExodusII_IO::read().
{
nodeset_ids.resize(num_node_sets);
if (num_node_sets > 0)
{
ex_err = exII::ex_get_node_set_ids(ex_id,
&nodeset_ids[0]);
EX_CHECK_ERR(ex_err, "Error retrieving nodeset information.");
message("All nodeset information retrieved successfully.");
// Resize appropriate data structures -- only do this once outnode the loop
num_nodes_per_set.resize(num_node_sets);
num_node_df_per_set.resize(num_node_sets);
}
char name_buffer[MAX_STR_LENGTH+1];
for (int i=0; i<num_node_sets; ++i)
{
ex_err = exII::ex_get_name(ex_id, exII::EX_NODE_SET,
nodeset_ids[i], name_buffer);
EX_CHECK_ERR(ex_err, "Error getting node set name.");
id_to_ns_names[nodeset_ids[i]] = name_buffer;
}
message("All node set names retrieved successfully.");
}
Reads the number of timesteps currently stored in the Exodus file and stores it in the num_time_steps variable.
Definition at line 746 of file exodusII_io_helper.C.
References inquire(), and num_time_steps.
Referenced by libMesh::ExodusII_IO::get_num_time_steps(), read_header(), and read_time_steps().
{
num_time_steps =
this->inquire(exII::EX_INQ_TIME, "Error retrieving number of time steps");
}
| void libMesh::ExodusII_IO_Helper::read_sideset | ( | int | id, |
| int | offset | ||
| ) |
Reads information about sideset id and inserts it into the global sideset array at the position offset.
Definition at line 624 of file exodusII_io_helper.C.
References elem_list, ex_err, ex_id, id_list, message(), num_df_per_set, num_sides_per_set, side_list, and ss_ids.
Referenced by libMesh::Nemesis_IO::read(), and libMesh::ExodusII_IO::read().
{
libmesh_assert_less (static_cast<unsigned int>(id), ss_ids.size());
libmesh_assert_less (static_cast<unsigned int>(id), num_sides_per_set.size());
libmesh_assert_less (static_cast<unsigned int>(id), num_df_per_set.size());
libmesh_assert_less_equal (static_cast<unsigned int>(offset), elem_list.size());
libmesh_assert_less_equal (static_cast<unsigned int>(offset), side_list.size());
ex_err = exII::ex_get_side_set_param(ex_id,
ss_ids[id],
&num_sides_per_set[id],
&num_df_per_set[id]);
EX_CHECK_ERR(ex_err, "Error retrieving sideset parameters.");
message("Parameters retrieved successfully for sideset: ", id);
// It's OK for offset==elem_list.size() as long as num_sides_per_set[id]==0
// because in that case we don't actually read anything...
#ifdef DEBUG
if (static_cast<unsigned int>(offset) == elem_list.size() ||
static_cast<unsigned int>(offset) == side_list.size() )
libmesh_assert_equal_to (num_sides_per_set[id], 0);
#endif
// Don't call ex_get_side_set unless there are actually sides there to get.
// Exodus prints an annoying warning in DEBUG mode otherwise...
if (num_sides_per_set[id] > 0)
{
ex_err = exII::ex_get_side_set(ex_id,
ss_ids[id],
&elem_list[offset],
&side_list[offset]);
EX_CHECK_ERR(ex_err, "Error retrieving sideset data.");
message("Data retrieved successfully for sideset: ", id);
for (int i=0; i<num_sides_per_set[id]; i++)
id_list[i+offset] = ss_ids[id];
}
}
Reads information about all of the sidesets in the ExodusII mesh file.
Definition at line 563 of file exodusII_io_helper.C.
References elem_list, ex_err, ex_id, id_list, id_to_ss_names, inquire(), message(), num_df_per_set, num_elem_all_sidesets, num_side_sets, num_sides_per_set, side_list, and ss_ids.
Referenced by libMesh::Nemesis_IO::read(), and libMesh::ExodusII_IO::read().
{
ss_ids.resize(num_side_sets);
if (num_side_sets > 0)
{
ex_err = exII::ex_get_side_set_ids(ex_id,
&ss_ids[0]);
EX_CHECK_ERR(ex_err, "Error retrieving sideset information.");
message("All sideset information retrieved successfully.");
// Resize appropriate data structures -- only do this once outside the loop
num_sides_per_set.resize(num_side_sets);
num_df_per_set.resize(num_side_sets);
// Inquire about the length of the concatenated side sets element list
num_elem_all_sidesets = inquire(exII::EX_INQ_SS_ELEM_LEN, "Error retrieving length of the concatenated side sets element list!");
elem_list.resize (num_elem_all_sidesets);
side_list.resize (num_elem_all_sidesets);
id_list.resize (num_elem_all_sidesets);
}
char name_buffer[MAX_STR_LENGTH+1];
for (int i=0; i<num_side_sets; ++i)
{
ex_err = exII::ex_get_name(ex_id, exII::EX_SIDE_SET,
ss_ids[i], name_buffer);
EX_CHECK_ERR(ex_err, "Error getting side set name.");
id_to_ss_names[ss_ids[i]] = name_buffer;
}
message("All side set names retrieved successfully.");
}
Reads and stores the timesteps in the 'time_steps' array.
Definition at line 731 of file exodusII_io_helper.C.
References ex_err, ex_id, num_time_steps, read_num_time_steps(), and time_steps.
Referenced by libMesh::ExodusII_IO::get_time_steps().
{
// Make sure we have an up-to-date count of the number of time steps in the file.
this->read_num_time_steps();
if (num_time_steps > 0)
{
time_steps.resize(num_time_steps);
ex_err = exII::ex_get_all_times(ex_id, &time_steps[0]);
EX_CHECK_ERR(ex_err, "Error reading timesteps!");
}
}
| void libMesh::ExodusII_IO_Helper::read_var_names | ( | ExodusVarType | type | ) |
Definition at line 794 of file exodusII_io_helper.C.
References elem_var_names, ELEMENTAL, GLOBAL, global_var_names, NODAL, nodal_var_names, num_elem_vars, num_global_vars, num_nodal_vars, and read_var_names_impl().
Referenced by check_existing_vars(), libMesh::ExodusII_IO::get_elem_var_names(), libMesh::ExodusII_IO::get_nodal_var_names(), read_elemental_var_values(), and read_nodal_var_values().
{
switch (type)
{
case NODAL:
this->read_var_names_impl("n", num_nodal_vars, nodal_var_names);
break;
case ELEMENTAL:
this->read_var_names_impl("e", num_elem_vars, elem_var_names);
break;
case GLOBAL:
this->read_var_names_impl("g", num_global_vars, global_var_names);
break;
default:
libmesh_error_msg("Unrecognized ExodusVarType " << type);
}
}
| void libMesh::ExodusII_IO_Helper::read_var_names_impl | ( | const char * | var_type, |
| int & | count, | ||
| std::vector< std::string > & | result | ||
| ) | [private] |
read_var_names() dispatches to this function.
Definition at line 814 of file exodusII_io_helper.C.
References ex_err, ex_id, libMesh::ExodusII_IO_Helper::NamesData::get_char_star(), libMesh::ExodusII_IO_Helper::NamesData::get_char_star_star(), libMesh::out, and verbose.
Referenced by read_var_names().
{
// First read and store the number of names we have
ex_err = exII::ex_get_var_param(ex_id, var_type, &count);
EX_CHECK_ERR(ex_err, "Error reading number of variables.");
// Do nothing if no variables are detected
if (count == 0)
return;
// Second read the actual names and convert them into a format we can use
NamesData names_table(count, MAX_STR_LENGTH);
ex_err = exII::ex_get_var_names(ex_id,
var_type,
count,
names_table.get_char_star_star()
);
EX_CHECK_ERR(ex_err, "Error reading variable names!");
if (verbose)
{
libMesh::out << "Read the variable(s) from the file:" << std::endl;
for (int i=0; i<count; i++)
libMesh::out << names_table.get_char_star(i) << std::endl;
}
// Allocate enough space for our variable name strings.
result.resize(count);
// Copy the char buffers into strings.
for (int i=0; i<count; i++)
result[i] = names_table.get_char_star(i); // calls string::op=(const char*)
}
Allows you to set a vector that is added to the coordinates of all of the nodes. Effectively, this "moves" the mesh to a particular position
Definition at line 1778 of file exodusII_io_helper.C.
References _coordinate_offset.
Referenced by libMesh::ExodusII_IO::set_coordinate_offset().
{
_coordinate_offset = p;
}
Sets the underlying value of the boolean flag _use_mesh_dimension_instead_of_spatial_dimension. By default, the value of this flag is false.
See the ExodusII_IO class documentation for a detailed description of this flag.
Definition at line 1773 of file exodusII_io_helper.C.
References _use_mesh_dimension_instead_of_spatial_dimension.
Referenced by libMesh::ExodusII_IO::use_mesh_dimension_instead_of_spatial_dimension().
| void libMesh::ExodusII_IO_Helper::write_element_values | ( | const MeshBase & | mesh, |
| const std::vector< Real > & | values, | ||
| int | timestep | ||
| ) |
Writes the vector of values to the element variables.
Definition at line 1615 of file exodusII_io_helper.C.
References _run_only_on_proc0, _single_precision, libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), data, end, ex_err, ex_id, get_block_id(), libMesh::DofObject::id(), libMesh::MeshBase::n_elem(), libMesh::MeshTools::n_elem(), num_elem_vars, libMesh::ParallelObject::processor_id(), and libMesh::Elem::subdomain_id().
Referenced by libMesh::ExodusII_IO::write_element_data().
{
if ((_run_only_on_proc0) && (this->processor_id() != 0))
return;
// Loop over the element blocks and write the data one block at a time
std::map<unsigned int, std::vector<unsigned int> > subdomain_map;
// Ask the file how many element vars it has, store it in the num_elem_vars variable.
ex_err = exII::ex_get_var_param(ex_id, "e", &num_elem_vars);
EX_CHECK_ERR(ex_err, "Error reading number of elemental variables.");
MeshBase::const_element_iterator mesh_it = mesh.active_elements_begin();
const MeshBase::const_element_iterator end = mesh.active_elements_end();
// loop through element and map between block and element vector
for (; mesh_it!=end; ++mesh_it)
{
const Elem * elem = *mesh_it;
subdomain_map[elem->subdomain_id()].push_back(elem->id());
}
// Use mesh.n_elem() to access into the values vector rather than
// the number of elements the Exodus writer thinks the mesh has,
// which may not include inactive elements.
dof_id_type n_elem = mesh.n_elem();
// For each variable, create a 'data' array which holds all the elemental variable
// values *for a given block* on this processor, then write that data vector to file
// before moving onto the next block.
for (unsigned int i=0; i<static_cast<unsigned>(num_elem_vars); ++i)
{
// The size of the subdomain map is the number of blocks.
std::map<unsigned int, std::vector<unsigned int> >::iterator it = subdomain_map.begin();
for (unsigned int j=0; it!=subdomain_map.end(); ++it, ++j)
{
const std::vector<unsigned int> & elem_nums = (*it).second;
const unsigned int num_elems_this_block =
cast_int<unsigned int>(elem_nums.size());
std::vector<Real> data(num_elems_this_block);
for (unsigned int k=0; k<num_elems_this_block; ++k)
data[k] = values[i*n_elem + elem_nums[k]];
if (_single_precision)
{
std::vector<float> cast_data(data.begin(), data.end());
ex_err = exII::ex_put_elem_var(ex_id,
timestep,
i+1,
this->get_block_id(j),
num_elems_this_block,
&cast_data[0]);
}
else
{
ex_err = exII::ex_put_elem_var(ex_id,
timestep,
i+1,
this->get_block_id(j),
num_elems_this_block,
&data[0]);
}
EX_CHECK_ERR(ex_err, "Error writing element values.");
}
}
ex_err = exII::ex_update(ex_id);
EX_CHECK_ERR(ex_err, "Error flushing buffers to file.");
}
| void libMesh::ExodusII_IO_Helper::write_elements | ( | const MeshBase & | mesh, |
| bool | use_discontinuous = false |
||
| ) | [virtual] |
Writes the elements contained in "mesh". FIXME: This only works for Meshes having a single type of element!
Reimplemented in libMesh::Nemesis_IO_Helper.
Definition at line 1202 of file exodusII_io_helper.C.
References _run_only_on_proc0, libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), libMesh::ExodusII_IO_Helper::ElementMaps::assign_conversion(), block_ids, connect, libMesh::MeshBase::elem(), elem_num_map, end, libMesh::Utility::enum_to_string(), ex_err, ex_id, libMesh::ExodusII_IO_Helper::Conversion::exodus_elem_type(), libMesh::ExodusII_IO_Helper::Conversion::get_canonical_type(), libMesh::ExodusII_IO_Helper::NamesData::get_char_star_star(), libMesh::ExodusII_IO_Helper::Conversion::get_inverse_node_map(), libMesh::DofObject::id(), libmesh_elem_num_to_exodus, libMesh::MeshBase::n_active_elem(), libMesh::Elem::n_nodes(), libMesh::Elem::node(), num_elem_blk, num_nodes_per_elem, libMesh::out, libMesh::ParallelObject::processor_id(), libMesh::ExodusII_IO_Helper::NamesData::push_back_entry(), libMesh::Elem::subdomain_id(), libMesh::MeshBase::subdomain_name(), libMesh::Elem::type(), and verbose.
Referenced by libMesh::ExodusII_IO::write(), and libMesh::ExodusII_IO::write_nodal_data_common().
{
// n_active_elem() is a parallel_only function
unsigned int n_active_elem = mesh.n_active_elem();
if ((_run_only_on_proc0) && (this->processor_id() != 0))
return;
// Map from block ID to a vector of element IDs in that block. Element
// IDs are now of type dof_id_type, subdomain IDs are of type subdomain_id_type.
typedef std::map<subdomain_id_type, std::vector<dof_id_type> > subdomain_map_type;
subdomain_map_type subdomain_map;
MeshBase::const_element_iterator mesh_it = mesh.active_elements_begin();
const MeshBase::const_element_iterator end = mesh.active_elements_end();
//loop through element and map between block and element vector
for (; mesh_it!=end; ++mesh_it)
{
const Elem * elem = *mesh_it;
subdomain_map[ elem->subdomain_id() ].push_back(elem->id());
}
// element map vector
num_elem_blk = cast_int<int>(subdomain_map.size());
block_ids.resize(num_elem_blk);
elem_num_map.resize(n_active_elem);
std::vector<int>::iterator curr_elem_map_end = elem_num_map.begin();
// Note: It appears that there is a bug in exodusII::ex_put_name where
// the index returned from the ex_id_lkup is erronously used. For now
// the work around is to use the alternative function ex_put_names, but
// this function requires a char** datastructure.
NamesData names_table(num_elem_blk, MAX_STR_LENGTH);
// This counter is used to fill up the libmesh_elem_num_to_exodus map in the loop below.
unsigned libmesh_elem_num_to_exodus_counter = 0;
// counter indexes into the block_ids vector
unsigned int counter = 0;
// node counter for discontinuous plotting
unsigned int node_counter = 0;
for (subdomain_map_type::iterator it=subdomain_map.begin(); it!=subdomain_map.end(); ++it)
{
block_ids[counter] = (*it).first;
names_table.push_back_entry(mesh.subdomain_name((*it).first));
// Get a reference to a vector of element IDs for this subdomain.
subdomain_map_type::mapped_type& tmp_vec = (*it).second;
ExodusII_IO_Helper::ElementMaps em;
//Use the first element in this block to get representative information.
//Note that Exodus assumes all elements in a block are of the same type!
//We are using that same assumption here!
const ExodusII_IO_Helper::Conversion conv = em.assign_conversion(mesh.elem(tmp_vec[0])->type());
num_nodes_per_elem = mesh.elem(tmp_vec[0])->n_nodes();
ex_err = exII::ex_put_elem_block(ex_id,
(*it).first,
conv.exodus_elem_type().c_str(),
tmp_vec.size(),
num_nodes_per_elem,
/*num_attr=*/0);
EX_CHECK_ERR(ex_err, "Error writing element block.");
connect.resize(tmp_vec.size()*num_nodes_per_elem);
for (unsigned int i=0; i<tmp_vec.size(); i++)
{
unsigned int elem_id = tmp_vec[i];
libmesh_elem_num_to_exodus[elem_id] = ++libmesh_elem_num_to_exodus_counter; // 1-based indexing for Exodus
const Elem* elem = mesh.elem(elem_id);
// We *might* be able to get away with writing mixed element
// types which happen to have the same number of nodes, but
// do we actually *want* to get away with that?
// .) No visualization software would be able to handle it.
// .) There'd be no way for us to read it back in reliably.
// .) Even elements with the same number of nodes may have different connectivities (?)
// This needs to be more than an assert so we don't fail
// with a mysterious segfault while trying to write mixed
// element meshes in optimized mode.
if (elem->type() != conv.get_canonical_type())
libmesh_error_msg("Error: Exodus requires all elements with a given subdomain ID to be the same type.\n" \
<< "Can't write both " \
<< Utility::enum_to_string(elem->type()) \
<< " and " \
<< Utility::enum_to_string(conv.get_canonical_type()) \
<< " in the same block!");
for (unsigned int j=0; j<static_cast<unsigned int>(num_nodes_per_elem); ++j)
{
unsigned connect_index = (i*num_nodes_per_elem)+j;
unsigned elem_node_index = conv.get_inverse_node_map(j); // inverse node map is for writing.
if (verbose)
{
libMesh::out << "Exodus node index: " << j
<< "=LibMesh node index " << elem_node_index << std::endl;
}
// FIXME: We are hard-coding the 1-based node numbering assumption here.
if (!use_discontinuous)
connect[connect_index] = elem->node(elem_node_index)+1;
else
connect[connect_index] = node_counter*num_nodes_per_elem+elem_node_index+1;
}
node_counter++;
}
ex_err = exII::ex_put_elem_conn(ex_id, (*it).first, &connect[0]);
EX_CHECK_ERR(ex_err, "Error writing element connectivities");
// This transform command stores its result in a range that begins at the third argument,
// so this command is adding values to the elem_num_map vector starting from curr_elem_map_end.
curr_elem_map_end = std::transform(tmp_vec.begin(),
tmp_vec.end(),
curr_elem_map_end,
std::bind2nd(std::plus<subdomain_map_type::mapped_type::value_type>(), 1)); // Adds one to each id to make a 1-based exodus file!
// But if we don't want to add one, we just want to put the values
// of tmp_vec into elem_map in the right location, we can use
// std::copy().
// curr_elem_map_end = std::copy(tmp_vec.begin(), tmp_vec.end(), curr_elem_map_end);
counter++;
}
// write out the element number map that we created
ex_err = exII::ex_put_elem_num_map(ex_id, &elem_num_map[0]);
EX_CHECK_ERR(ex_err, "Error writing element map");
// Write out the block names
if (num_elem_blk > 0)
{
ex_err = exII::ex_put_names(ex_id, exII::EX_ELEM_BLOCK, names_table.get_char_star_star());
EX_CHECK_ERR(ex_err, "Error writing element names");
}
}
| void libMesh::ExodusII_IO_Helper::write_global_values | ( | const std::vector< Real > & | values, |
| int | timestep | ||
| ) |
Writes the vector of global variables.
Definition at line 1751 of file exodusII_io_helper.C.
References _run_only_on_proc0, _single_precision, ex_err, ex_id, num_global_vars, and libMesh::ParallelObject::processor_id().
Referenced by libMesh::Nemesis_IO::write_global_data(), and libMesh::ExodusII_IO::write_global_data().
{
if ((_run_only_on_proc0) && (this->processor_id() != 0))
return;
if (_single_precision)
{
std::vector<float> cast_values(values.begin(), values.end());
ex_err = exII::ex_put_glob_vars(ex_id, timestep, num_global_vars, &cast_values[0]);
}
else
{
ex_err = exII::ex_put_glob_vars(ex_id, timestep, num_global_vars, &values[0]);
}
EX_CHECK_ERR(ex_err, "Error writing global values.");
ex_err = exII::ex_update(ex_id);
EX_CHECK_ERR(ex_err, "Error flushing buffers to file.");
}
| void libMesh::ExodusII_IO_Helper::write_information_records | ( | const std::vector< std::string > & | records | ) |
Writes the vector of information records.
Definition at line 1712 of file exodusII_io_helper.C.
References _run_only_on_proc0, libMesh::err, ex_err, ex_id, libMesh::ExodusII_IO_Helper::NamesData::get_char_star_star(), inquire(), libMesh::ParallelObject::processor_id(), and libMesh::ExodusII_IO_Helper::NamesData::push_back_entry().
Referenced by libMesh::Nemesis_IO::write_information_records(), and libMesh::ExodusII_IO::write_information_records().
{
if ((_run_only_on_proc0) && (this->processor_id() != 0))
return;
// There may already be information records in the file (for
// example, if we're appending) and in that case, according to the
// Exodus documentation, writing more information records is not
// supported.
int num_info = inquire(exII::EX_INQ_INFO, "Error retrieving the number of information records from file!");
if (num_info > 0)
{
libMesh::err << "Warning! The Exodus file already contains information records.\n"
<< "Exodus does not support writing additional records in this situation."
<< std::endl;
return;
}
int num_records = cast_int<int>(records.size());
if (num_records > 0)
{
NamesData info(num_records, MAX_LINE_LENGTH);
// If an entry is longer than MAX_LINE_LENGTH characters it's not an error, we just
// write the first MAX_LINE_LENGTH characters to the file.
for (unsigned i=0; i<records.size(); ++i)
info.push_back_entry(records[i]);
ex_err = exII::ex_put_info(ex_id, num_records, info.get_char_star_star());
EX_CHECK_ERR(ex_err, "Error writing global values.");
ex_err = exII::ex_update(ex_id);
EX_CHECK_ERR(ex_err, "Error flushing buffers to file.");
}
}
| void libMesh::ExodusII_IO_Helper::write_nodal_coordinates | ( | const MeshBase & | mesh, |
| bool | use_discontinuous = false |
||
| ) | [virtual] |
Writes the nodal coordinates contained in "mesh"
Reimplemented in libMesh::Nemesis_IO_Helper.
Definition at line 1095 of file exodusII_io_helper.C.
References _coordinate_offset, _run_only_on_proc0, _single_precision, libMesh::MeshBase::active_elements_begin(), libMesh::MeshBase::active_elements_end(), end, ex_err, ex_id, libMesh::DofObject::id(), node_num_map, libMesh::MeshBase::nodes_begin(), libMesh::MeshBase::nodes_end(), num_nodes, libMesh::ParallelObject::processor_id(), x, y, and z.
Referenced by libMesh::ExodusII_IO::write(), and libMesh::ExodusII_IO::write_nodal_data_common().
{
if ((_run_only_on_proc0) && (this->processor_id() != 0))
return;
// Make room in the coordinate vectors
x.resize(num_nodes);
y.resize(num_nodes);
z.resize(num_nodes);
// And in the node_num_map - since the nodes aren't organized in
// blocks, libmesh will always write out the identity map
// here... unless there has been some refinement and coarsening. If
// the nodes aren't renumbered after refinement and coarsening,
// there may be 'holes' in the numbering, so we write out the node
// map just to be on the safe side.
node_num_map.resize(num_nodes);
if (!use_discontinuous)
{
MeshBase::const_node_iterator it = mesh.nodes_begin();
const MeshBase::const_node_iterator end = mesh.nodes_end();
for (unsigned i = 0; it != end; ++it, ++i)
{
const Node* node = *it;
x[i] = (*node)(0) + _coordinate_offset(0);
#if LIBMESH_DIM > 1
y[i]=(*node)(1) + _coordinate_offset(1);
#else
y[i]=0.;
#endif
#if LIBMESH_DIM > 2
z[i]=(*node)(2) + _coordinate_offset(2);
#else
z[i]=0.;
#endif
// Fill in node_num_map entry with the proper (1-based) node id
node_num_map[i] = node->id() + 1;
}
}
else
{
MeshBase::const_element_iterator it = mesh.active_elements_begin();
const MeshBase::const_element_iterator end = mesh.active_elements_end();
unsigned int i = 0;
for (; it!=end; ++it)
for (unsigned int n=0; n<(*it)->n_nodes(); n++)
{
x[i]=(*it)->point(n)(0);
#if LIBMESH_DIM > 1
y[i]=(*it)->point(n)(1);
#else
y[i]=0.;
#endif
#if LIBMESH_DIM > 2
z[i]=(*it)->point(n)(2);
#else
z[i]=0.;
#endif
// Let's skip the node_num_map in the discontinuous case,
// since we're effectively duplicating nodes for the
// sake of discontinuous visualization, so it isn't clear
// how to deal with node_num_map here. It's only optional
// anyway, so no big deal.
i++;
}
}
if (_single_precision)
{
std::vector<float>
x_single(x.begin(), x.end()),
y_single(y.begin(), y.end()),
z_single(z.begin(), z.end());
ex_err = exII::ex_put_coord(ex_id,
x_single.empty() ? NULL : &x_single[0],
y_single.empty() ? NULL : &y_single[0],
z_single.empty() ? NULL : &z_single[0]);
}
else
{
ex_err = exII::ex_put_coord(ex_id,
x.empty() ? NULL : &x[0],
y.empty() ? NULL : &y[0],
z.empty() ? NULL : &z[0]);
}
EX_CHECK_ERR(ex_err, "Error writing coordinates to Exodus file.");
if (!use_discontinuous)
{
// Also write the (1-based) node_num_map to the file.
ex_err = exII::ex_put_node_num_map(ex_id, &node_num_map[0]);
EX_CHECK_ERR(ex_err, "Error writing node_num_map");
}
}
| void libMesh::ExodusII_IO_Helper::write_nodal_values | ( | int | var_id, |
| const std::vector< Real > & | values, | ||
| int | timestep | ||
| ) |
Writes the vector of values to a nodal variable.
Definition at line 1690 of file exodusII_io_helper.C.
References _run_only_on_proc0, _single_precision, ex_err, ex_id, num_nodes, and libMesh::ParallelObject::processor_id().
Referenced by libMesh::ExodusII_IO::write_nodal_data(), and libMesh::ExodusII_IO::write_nodal_data_discontinuous().
{
if ((_run_only_on_proc0) && (this->processor_id() != 0))
return;
if (_single_precision)
{
std::vector<float> cast_values(values.begin(), values.end());
ex_err = exII::ex_put_nodal_var(ex_id, timestep, var_id, num_nodes, &cast_values[0]);
}
else
{
ex_err = exII::ex_put_nodal_var(ex_id, timestep, var_id, num_nodes, &values[0]);
}
EX_CHECK_ERR(ex_err, "Error writing nodal values.");
ex_err = exII::ex_update(ex_id);
EX_CHECK_ERR(ex_err, "Error flushing buffers to file.");
}
| void libMesh::ExodusII_IO_Helper::write_nodesets | ( | const MeshBase & | mesh | ) | [virtual] |
Writes the nodesets contained in "mesh"
Reimplemented in libMesh::Nemesis_IO_Helper.
Definition at line 1423 of file exodusII_io_helper.C.
References _run_only_on_proc0, libMesh::BoundaryInfo::build_node_boundary_ids(), libMesh::BoundaryInfo::build_node_list(), ex_err, ex_id, libMesh::MeshBase::get_boundary_info(), libMesh::BoundaryInfo::get_nodeset_name(), libMesh::ParallelObject::processor_id(), and libMesh::ExodusII_IO_Helper::NamesData::push_back_entry().
Referenced by libMesh::ExodusII_IO::write(), and libMesh::ExodusII_IO::write_nodal_data_common().
{
if ((_run_only_on_proc0) && (this->processor_id() != 0))
return;
std::vector< dof_id_type > nl;
std::vector< boundary_id_type > il;
mesh.get_boundary_info().build_node_list(nl, il);
// Maps from nodeset id to the nodes
std::map<boundary_id_type, std::vector<int> > node;
// Accumulate the vectors to pass into ex_put_node_set
for (unsigned int i=0; i<nl.size(); i++)
node[il[i]].push_back(nl[i]+1);
std::vector<boundary_id_type> node_boundary_ids;
mesh.get_boundary_info().build_node_boundary_ids(node_boundary_ids);
// Write out the nodeset names, but only if there is something to write
if (node_boundary_ids.size() > 0)
{
NamesData names_table(node_boundary_ids.size(), MAX_STR_LENGTH);
for (unsigned int i=0; i<node_boundary_ids.size(); i++)
{
boundary_id_type nodeset_id = node_boundary_ids[i];
int actual_id = nodeset_id;
names_table.push_back_entry(mesh.get_boundary_info().get_nodeset_name(nodeset_id));
ex_err = exII::ex_put_node_set_param(ex_id, actual_id, node[nodeset_id].size(), 0);
EX_CHECK_ERR(ex_err, "Error writing nodeset parameters");
ex_err = exII::ex_put_node_set(ex_id, actual_id, &node[nodeset_id][0]);
EX_CHECK_ERR(ex_err, "Error writing nodesets");
}
// Write out the nodeset names
ex_err = exII::ex_put_names(ex_id, exII::EX_NODE_SET, names_table.get_char_star_star());
EX_CHECK_ERR(ex_err, "Error writing nodeset names");
}
}
| void libMesh::ExodusII_IO_Helper::write_sidesets | ( | const MeshBase & | mesh | ) | [virtual] |
Writes the sidesets contained in "mesh"
We need to build up active elements if AMR is enabled and add them to the exodus sidesets instead of the potentially inactive "parent" elements
Reimplemented in libMesh::Nemesis_IO_Helper.
Definition at line 1351 of file exodusII_io_helper.C.
References _run_only_on_proc0, libMesh::Elem::active_family_tree_by_side(), libMesh::ExodusII_IO_Helper::ElementMaps::assign_conversion(), libMesh::BoundaryInfo::build_side_boundary_ids(), libMesh::BoundaryInfo::build_side_list(), libMesh::MeshBase::elem(), ex_err, ex_id, libMesh::MeshBase::get_boundary_info(), libMesh::ExodusII_IO_Helper::Conversion::get_inverse_side_map(), libMesh::BoundaryInfo::get_sideset_name(), libmesh_elem_num_to_exodus, libMesh::ParallelObject::processor_id(), libMesh::ExodusII_IO_Helper::NamesData::push_back_entry(), and side.
Referenced by libMesh::ExodusII_IO::write(), and libMesh::ExodusII_IO::write_nodal_data_common().
{
if ((_run_only_on_proc0) && (this->processor_id() != 0))
return;
ExodusII_IO_Helper::ElementMaps em;
std::vector< dof_id_type > el;
std::vector< unsigned short int > sl;
std::vector< boundary_id_type > il;
mesh.get_boundary_info().build_side_list(el, sl, il);
// Maps from sideset id to the element and sides
std::map<int, std::vector<int> > elem;
std::map<int, std::vector<int> > side;
// Accumulate the vectors to pass into ex_put_side_set
for (unsigned int i=0; i<el.size(); i++)
{
std::vector<const Elem *> family;
#ifdef LIBMESH_ENABLE_AMR
mesh.elem(el[i])->active_family_tree_by_side(family, sl[i], false);
#else
family.push_back(mesh.elem(el[i]));
#endif
for (unsigned int j=0; j<family.size(); ++j)
{
const ExodusII_IO_Helper::Conversion conv = em.assign_conversion(mesh.elem(family[j]->id())->type());
// Use the libmesh to exodus datastructure map to get the proper sideset IDs
// The datastructure contains the "collapsed" contiguous ids
elem[il[i]].push_back(libmesh_elem_num_to_exodus[family[j]->id()]);
side[il[i]].push_back(conv.get_inverse_side_map(sl[i]));
}
}
std::vector<boundary_id_type> side_boundary_ids;
mesh.get_boundary_info().build_side_boundary_ids(side_boundary_ids);
// Write out the sideset names, but only if there is something to write
if (side_boundary_ids.size() > 0)
{
NamesData names_table(side_boundary_ids.size(), MAX_STR_LENGTH);
for (unsigned int i=0; i<side_boundary_ids.size(); i++)
{
boundary_id_type ss_id = side_boundary_ids[i];
int actual_id = ss_id;
names_table.push_back_entry(mesh.get_boundary_info().get_sideset_name(ss_id));
ex_err = exII::ex_put_side_set_param(ex_id, actual_id, elem[ss_id].size(), 0);
EX_CHECK_ERR(ex_err, "Error writing sideset parameters");
ex_err = exII::ex_put_side_set(ex_id, actual_id, &elem[ss_id][0], &side[ss_id][0]);
EX_CHECK_ERR(ex_err, "Error writing sidesets");
}
ex_err = exII::ex_put_names(ex_id, exII::EX_SIDE_SET, names_table.get_char_star_star());
EX_CHECK_ERR(ex_err, "Error writing sideset names");
}
}
| void libMesh::ExodusII_IO_Helper::write_timestep | ( | int | timestep, |
| Real | time | ||
| ) |
Writes the time for the timestep
Definition at line 1601 of file exodusII_io_helper.C.
References _run_only_on_proc0, ex_err, ex_id, and libMesh::ParallelObject::processor_id().
Referenced by libMesh::Nemesis_IO::write_timestep(), and libMesh::ExodusII_IO::write_timestep().
{
if ((_run_only_on_proc0) && (this->processor_id() != 0))
return;
ex_err = exII::ex_put_time(ex_id, timestep, &time);
EX_CHECK_ERR(ex_err, "Error writing timestep.");
ex_err = exII::ex_update(ex_id);
EX_CHECK_ERR(ex_err, "Error flushing buffers to file.");
}
| void libMesh::ExodusII_IO_Helper::write_var_names | ( | ExodusVarType | type, |
| std::vector< std::string > & | names | ||
| ) | [private] |
Wraps calls to exII::ex_put_var_names() and exII::ex_put_var_param(). The enumeration controls whether nodal, elemental, or global variable names are read and which class members are filled in.
Definition at line 854 of file exodusII_io_helper.C.
References ELEMENTAL, GLOBAL, NODAL, num_elem_vars, num_global_vars, num_nodal_vars, and write_var_names_impl().
Referenced by initialize_element_variables(), initialize_global_variables(), and initialize_nodal_variables().
{
switch (type)
{
case NODAL:
this->write_var_names_impl("n", num_nodal_vars, names);
break;
case ELEMENTAL:
this->write_var_names_impl("e", num_elem_vars, names);
break;
case GLOBAL:
this->write_var_names_impl("g", num_global_vars, names);
break;
default:
libmesh_error_msg("Unrecognized ExodusVarType " << type);
}
}
| void libMesh::ExodusII_IO_Helper::write_var_names_impl | ( | const char * | var_type, |
| int & | count, | ||
| std::vector< std::string > & | names | ||
| ) | [private] |
write_var_names() dispatches to this function.
Definition at line 874 of file exodusII_io_helper.C.
References ex_err, ex_id, libMesh::out, libMesh::ExodusII_IO_Helper::NamesData::push_back_entry(), and verbose.
Referenced by write_var_names().
{
// Update the count variable so that it's available to other parts of the class.
count = cast_int<int>(names.size());
// Write that number of variables to the file.
ex_err = exII::ex_put_var_param(ex_id, var_type, count);
EX_CHECK_ERR(ex_err, "Error setting number of vars.");
if (names.size() > 0)
{
NamesData names_table(names.size(), MAX_STR_LENGTH);
// Store the input names in the format required by Exodus.
for (unsigned i=0; i<names.size(); ++i)
names_table.push_back_entry(names[i]);
if (verbose)
{
libMesh::out << "Writing variable name(s) to file: " << std::endl;
for (unsigned i=0; i<names.size(); ++i)
libMesh::out << names_table.get_char_star(i) << std::endl;
}
ex_err = exII::ex_put_var_names(ex_id,
var_type,
cast_int<int>(names.size()),
names_table.get_char_star_star()
);
EX_CHECK_ERR(ex_err, "Error writing variable names.");
}
}
const Parallel::Communicator& libMesh::ParallelObject::_communicator [protected, inherited] |
Definition at line 104 of file parallel_object.h.
Referenced by libMesh::EquationSystems::build_solution_vector(), libMesh::ParallelObject::comm(), libMesh::EquationSystems::get_solution(), libMesh::ParallelObject::n_processors(), libMesh::ParallelObject::operator=(), and libMesh::ParallelObject::processor_id().
Point libMesh::ExodusII_IO_Helper::_coordinate_offset [protected] |
Definition at line 555 of file exodusII_io_helper.h.
Referenced by set_coordinate_offset(), and write_nodal_coordinates().
bool libMesh::ExodusII_IO_Helper::_elem_vars_initialized [protected] |
Definition at line 541 of file exodusII_io_helper.h.
Referenced by initialize_element_variables().
bool libMesh::ExodusII_IO_Helper::_global_vars_initialized [protected] |
Definition at line 544 of file exodusII_io_helper.h.
Referenced by initialize_global_variables().
bool libMesh::ExodusII_IO_Helper::_nodal_vars_initialized [protected] |
Definition at line 547 of file exodusII_io_helper.h.
Referenced by initialize_nodal_variables().
bool libMesh::ExodusII_IO_Helper::_run_only_on_proc0 [protected] |
Definition at line 538 of file exodusII_io_helper.h.
Referenced by close(), create(), initialize(), initialize_element_variables(), initialize_global_variables(), initialize_nodal_variables(), write_element_values(), write_elements(), write_global_values(), write_information_records(), write_nodal_coordinates(), write_nodal_values(), write_nodesets(), write_sidesets(), and write_timestep().
bool libMesh::ExodusII_IO_Helper::_single_precision [protected] |
Definition at line 558 of file exodusII_io_helper.h.
Referenced by create(), libMesh::Nemesis_IO_Helper::create(), write_element_values(), write_global_values(), write_nodal_coordinates(), libMesh::Nemesis_IO_Helper::write_nodal_coordinates(), and write_nodal_values().
Definition at line 552 of file exodusII_io_helper.h.
Referenced by initialize(), and use_mesh_dimension_instead_of_spatial_dimension().
| std::vector<int> libMesh::ExodusII_IO_Helper::block_ids |
Definition at line 406 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), get_block_id(), get_block_name(), libMesh::Nemesis_IO::read(), read_block_info(), read_elem_in_block(), read_elemental_var_values(), and write_elements().
| std::vector<int> libMesh::ExodusII_IO_Helper::connect |
Definition at line 409 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_elem_in_block(), and write_elements().
| std::string libMesh::ExodusII_IO_Helper::current_filename |
Definition at line 523 of file exodusII_io_helper.h.
Referenced by create(), open(), and libMesh::ExodusII_IO::write_nodal_data_common().
| std::vector<int> libMesh::ExodusII_IO_Helper::elem_list |
Definition at line 430 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_sideset(), and read_sideset_info().
| std::vector<int> libMesh::ExodusII_IO_Helper::elem_num_map |
Definition at line 445 of file exodusII_io_helper.h.
Referenced by libMesh::ExodusII_IO::read(), read_elem_num_map(), read_elemental_var_values(), and write_elements().
| std::vector<char> libMesh::ExodusII_IO_Helper::elem_type |
Definition at line 460 of file exodusII_io_helper.h.
Referenced by ExodusII_IO_Helper(), get_elem_type(), libMesh::Nemesis_IO::read(), and read_elem_in_block().
| std::vector<std::string> libMesh::ExodusII_IO_Helper::elem_var_names |
Definition at line 491 of file exodusII_io_helper.h.
Referenced by libMesh::ExodusII_IO::get_elem_var_names(), initialize_element_variables(), read_elemental_var_values(), and read_var_names().
| std::vector<Real> libMesh::ExodusII_IO_Helper::elem_var_values |
Definition at line 494 of file exodusII_io_helper.h.
Referenced by libMesh::ExodusII_IO::copy_elemental_solution(), and read_elemental_var_values().
Definition at line 370 of file exodusII_io_helper.h.
Referenced by close(), initialize(), initialize_element_variables(), inquire(), read_block_info(), read_elem_in_block(), read_elem_num_map(), read_elemental_var_values(), read_header(), read_nodal_var_values(), read_node_num_map(), read_nodes(), read_nodeset(), read_nodeset_info(), read_sideset(), read_sideset_info(), read_time_steps(), read_var_names_impl(), libMesh::Nemesis_IO::write(), write_element_values(), write_elements(), libMesh::Nemesis_IO_Helper::write_elements(), libMesh::Nemesis_IO_Helper::write_exodus_initialization_info(), write_global_values(), write_information_records(), write_nodal_coordinates(), libMesh::Nemesis_IO_Helper::write_nodal_coordinates(), write_nodal_values(), write_nodesets(), libMesh::Nemesis_IO_Helper::write_nodesets(), write_sidesets(), libMesh::Nemesis_IO_Helper::write_sidesets(), write_timestep(), write_var_names_impl(), and libMesh::Nemesis_IO_Helper::~Nemesis_IO_Helper().
Definition at line 367 of file exodusII_io_helper.h.
Referenced by close(), create(), libMesh::Nemesis_IO_Helper::create(), libMesh::Nemesis_IO_Helper::get_cmap_params(), libMesh::Nemesis_IO_Helper::get_eb_info_global(), libMesh::Nemesis_IO_Helper::get_elem_cmap(), libMesh::Nemesis_IO_Helper::get_elem_map(), libMesh::Nemesis_IO_Helper::get_init_global(), libMesh::Nemesis_IO_Helper::get_init_info(), libMesh::Nemesis_IO_Helper::get_loadbal_param(), libMesh::Nemesis_IO_Helper::get_node_cmap(), libMesh::Nemesis_IO_Helper::get_node_map(), libMesh::Nemesis_IO_Helper::get_ns_param_global(), libMesh::Nemesis_IO_Helper::get_ss_param_global(), initialize(), initialize_element_variables(), inquire(), open(), libMesh::Nemesis_IO_Helper::put_cmap_params(), libMesh::Nemesis_IO_Helper::put_eb_info_global(), libMesh::Nemesis_IO_Helper::put_elem_cmap(), libMesh::Nemesis_IO_Helper::put_elem_map(), libMesh::Nemesis_IO_Helper::put_init_global(), libMesh::Nemesis_IO_Helper::put_init_info(), libMesh::Nemesis_IO_Helper::put_loadbal_param(), libMesh::Nemesis_IO_Helper::put_n_coord(), libMesh::Nemesis_IO_Helper::put_node_cmap(), libMesh::Nemesis_IO_Helper::put_node_map(), libMesh::Nemesis_IO_Helper::put_ns_param_global(), libMesh::Nemesis_IO_Helper::put_ss_param_global(), read_block_info(), read_elem_in_block(), read_elem_num_map(), read_elemental_var_values(), read_header(), read_nodal_var_values(), read_node_num_map(), read_nodes(), read_nodeset(), read_nodeset_info(), read_sideset(), read_sideset_info(), read_time_steps(), read_var_names_impl(), libMesh::Nemesis_IO::write(), write_element_values(), write_elements(), libMesh::Nemesis_IO_Helper::write_elements(), libMesh::Nemesis_IO_Helper::write_exodus_initialization_info(), write_global_values(), write_information_records(), write_nodal_coordinates(), libMesh::Nemesis_IO_Helper::write_nodal_coordinates(), write_nodal_values(), write_nodesets(), libMesh::Nemesis_IO_Helper::write_nodesets(), write_sidesets(), libMesh::Nemesis_IO_Helper::write_sidesets(), write_timestep(), write_var_names_impl(), and libMesh::Nemesis_IO_Helper::~Nemesis_IO_Helper().
| std::vector<int> libMesh::ExodusII_IO_Helper::exodus_elem_num_to_libmesh |
Definition at line 465 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), and libMesh::Nemesis_IO_Helper::write_elements().
| std::vector<int> libMesh::ExodusII_IO_Helper::exodus_node_num_to_libmesh |
Definition at line 470 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), and libMesh::Nemesis_IO_Helper::write_nodal_coordinates().
| std::vector<std::string> libMesh::ExodusII_IO_Helper::global_var_names |
Definition at line 497 of file exodusII_io_helper.h.
Referenced by initialize_global_variables(), and read_var_names().
| std::vector<int> libMesh::ExodusII_IO_Helper::id_list |
Definition at line 439 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_sideset(), and read_sideset_info().
| std::map<int, std::string> libMesh::ExodusII_IO_Helper::id_to_block_names |
Definition at line 500 of file exodusII_io_helper.h.
Referenced by get_block_name(), and read_block_info().
| std::map<int, std::string> libMesh::ExodusII_IO_Helper::id_to_ns_names |
Definition at line 502 of file exodusII_io_helper.h.
Referenced by get_node_set_name(), and read_nodeset_info().
| std::map<int, std::string> libMesh::ExodusII_IO_Helper::id_to_ss_names |
Definition at line 501 of file exodusII_io_helper.h.
Referenced by get_side_set_name(), and read_sideset_info().
| std::map<int, int> libMesh::ExodusII_IO_Helper::libmesh_elem_num_to_exodus |
Definition at line 464 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::Nemesis_IO_Helper::compute_elem_communication_maps(), libMesh::Nemesis_IO_Helper::compute_element_maps(), write_elements(), write_sidesets(), and libMesh::Nemesis_IO_Helper::write_sidesets().
| std::map<int, int> libMesh::ExodusII_IO_Helper::libmesh_node_num_to_exodus |
| std::vector<std::string> libMesh::ExodusII_IO_Helper::nodal_var_names |
Definition at line 482 of file exodusII_io_helper.h.
Referenced by libMesh::ExodusII_IO::get_nodal_var_names(), initialize_nodal_variables(), read_nodal_var_values(), and read_var_names().
| std::vector<Real> libMesh::ExodusII_IO_Helper::nodal_var_values |
Definition at line 485 of file exodusII_io_helper.h.
Referenced by libMesh::ExodusII_IO::copy_nodal_solution(), and read_nodal_var_values().
| std::vector<int> libMesh::ExodusII_IO_Helper::node_list |
Definition at line 436 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), and read_nodeset().
| std::vector<int> libMesh::ExodusII_IO_Helper::node_num_map |
Definition at line 442 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_node_num_map(), and write_nodal_coordinates().
| std::vector<int> libMesh::ExodusII_IO_Helper::nodeset_ids |
Definition at line 415 of file exodusII_io_helper.h.
Referenced by get_node_set_id(), get_node_set_name(), libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_nodeset(), and read_nodeset_info().
Definition at line 400 of file exodusII_io_helper.h.
Referenced by read_elem_in_block().
| std::vector<int> libMesh::ExodusII_IO_Helper::num_df_per_set |
Definition at line 424 of file exodusII_io_helper.h.
Referenced by read_sideset(), and read_sideset_info().
Definition at line 373 of file exodusII_io_helper.h.
Referenced by initialize(), print_header(), read_header(), and libMesh::Nemesis_IO_Helper::write_exodus_initialization_info().
Definition at line 382 of file exodusII_io_helper.h.
Referenced by initialize(), print_header(), libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_elem_num_map(), read_elemental_var_values(), read_header(), and libMesh::Nemesis_IO_Helper::write_exodus_initialization_info().
Definition at line 403 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO::read(), and read_sideset_info().
Definition at line 385 of file exodusII_io_helper.h.
Referenced by initialize(), initialize_element_variables(), print_header(), libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_block_info(), read_elemental_var_values(), read_header(), write_elements(), and libMesh::Nemesis_IO_Helper::write_exodus_initialization_info().
Definition at line 394 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_elem_in_block(), and read_elemental_var_values().
Definition at line 488 of file exodusII_io_helper.h.
Referenced by initialize_element_variables(), read_header(), read_var_names(), write_element_values(), and write_var_names().
Definition at line 376 of file exodusII_io_helper.h.
Referenced by initialize_global_variables(), read_header(), read_var_names(), write_global_values(), and write_var_names().
Definition at line 479 of file exodusII_io_helper.h.
Referenced by initialize_nodal_variables(), read_header(), read_var_names(), and write_var_names().
| std::vector<int> libMesh::ExodusII_IO_Helper::num_node_df_per_set |
Definition at line 427 of file exodusII_io_helper.h.
Referenced by read_nodeset(), and read_nodeset_info().
Definition at line 388 of file exodusII_io_helper.h.
Referenced by initialize(), print_header(), libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_header(), read_nodeset_info(), and libMesh::Nemesis_IO_Helper::write_exodus_initialization_info().
Definition at line 379 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), initialize(), print_header(), print_nodes(), libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_header(), read_nodal_var_values(), read_node_num_map(), read_nodes(), libMesh::Nemesis_IO_Helper::write_exodus_initialization_info(), write_nodal_coordinates(), and write_nodal_values().
Definition at line 397 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_elem_in_block(), write_elements(), and libMesh::Nemesis_IO_Helper::write_elements().
| std::vector<int> libMesh::ExodusII_IO_Helper::num_nodes_per_set |
Definition at line 421 of file exodusII_io_helper.h.
Referenced by read_nodeset(), and read_nodeset_info().
Definition at line 391 of file exodusII_io_helper.h.
Referenced by initialize(), print_header(), libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_header(), read_sideset_info(), and libMesh::Nemesis_IO_Helper::write_exodus_initialization_info().
| std::vector<int> libMesh::ExodusII_IO_Helper::num_sides_per_set |
Definition at line 418 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_sideset(), and read_sideset_info().
Definition at line 473 of file exodusII_io_helper.h.
Referenced by libMesh::ExodusII_IO::get_num_time_steps(), read_num_time_steps(), and read_time_steps().
Definition at line 513 of file exodusII_io_helper.h.
Referenced by libMesh::ExodusII_IO::copy_elemental_solution(), libMesh::ExodusII_IO::copy_nodal_solution(), libMesh::ExodusII_IO::get_num_time_steps(), libMesh::ExodusII_IO::get_time_steps(), and open().
Definition at line 509 of file exodusII_io_helper.h.
Referenced by create(), libMesh::Nemesis_IO_Helper::create(), libMesh::ExodusII_IO::get_num_time_steps(), open(), libMesh::ExodusII_IO::write(), libMesh::ExodusII_IO::write_element_data(), libMesh::Nemesis_IO::write_global_data(), libMesh::ExodusII_IO::write_global_data(), libMesh::Nemesis_IO::write_information_records(), libMesh::ExodusII_IO::write_information_records(), libMesh::Nemesis_IO::write_nodal_data(), and libMesh::ExodusII_IO::write_nodal_data_common().
| std::vector<int> libMesh::ExodusII_IO_Helper::side_list |
Definition at line 433 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_sideset(), and read_sideset_info().
| std::vector<int> libMesh::ExodusII_IO_Helper::ss_ids |
Definition at line 412 of file exodusII_io_helper.h.
Referenced by get_side_set_id(), get_side_set_name(), read_sideset(), and read_sideset_info().
| std::vector<Real> libMesh::ExodusII_IO_Helper::time_steps |
Definition at line 476 of file exodusII_io_helper.h.
Referenced by libMesh::ExodusII_IO::get_time_steps(), and read_time_steps().
| std::vector<char> libMesh::ExodusII_IO_Helper::title |
Definition at line 457 of file exodusII_io_helper.h.
Referenced by ExodusII_IO_Helper(), print_header(), and read_header().
Definition at line 505 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::compute_communication_map_parameters(), libMesh::Nemesis_IO_Helper::compute_internal_and_border_elems_and_internal_nodes(), libMesh::Nemesis_IO_Helper::compute_node_communication_maps(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), create(), libMesh::Nemesis_IO_Helper::create(), libMesh::Nemesis_IO_Helper::get_cmap_params(), libMesh::Nemesis_IO_Helper::get_eb_info_global(), libMesh::Nemesis_IO_Helper::get_elem_cmap(), libMesh::Nemesis_IO_Helper::get_elem_map(), libMesh::Nemesis_IO_Helper::get_init_global(), libMesh::Nemesis_IO_Helper::get_init_info(), libMesh::Nemesis_IO_Helper::get_loadbal_param(), libMesh::Nemesis_IO_Helper::get_node_cmap(), libMesh::Nemesis_IO_Helper::get_node_map(), libMesh::Nemesis_IO_Helper::get_ns_param_global(), libMesh::Nemesis_IO_Helper::get_ss_param_global(), message(), open(), print_header(), libMesh::Nemesis_IO_Helper::put_node_cmap(), read_elem_in_block(), read_elem_num_map(), read_node_num_map(), read_var_names_impl(), libMesh::ExodusII_IO::verbose(), libMesh::Nemesis_IO::verbose(), write_elements(), libMesh::Nemesis_IO_Helper::write_nodesets(), libMesh::Nemesis_IO_Helper::write_sidesets(), and write_var_names_impl().
| std::vector<Real> libMesh::ExodusII_IO_Helper::x |
Definition at line 448 of file exodusII_io_helper.h.
Referenced by print_nodes(), libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_nodes(), write_nodal_coordinates(), and libMesh::Nemesis_IO_Helper::write_nodal_coordinates().
| std::vector<Real> libMesh::ExodusII_IO_Helper::y |
Definition at line 451 of file exodusII_io_helper.h.
Referenced by print_nodes(), libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_nodes(), write_nodal_coordinates(), and libMesh::Nemesis_IO_Helper::write_nodal_coordinates().
| std::vector<Real> libMesh::ExodusII_IO_Helper::z |
Definition at line 454 of file exodusII_io_helper.h.
Referenced by print_nodes(), libMesh::Nemesis_IO::read(), libMesh::ExodusII_IO::read(), read_nodes(), write_nodal_coordinates(), and libMesh::Nemesis_IO_Helper::write_nodal_coordinates().