Handles reading and writing of Exodus binary files. More...
#include <exodusII_io.h>
Public Member Functions | |
ExodusII_IO (MeshBase &mesh, bool single_precision=false) | |
virtual | ~ExodusII_IO () |
virtual void | read (const std::string &name) override |
virtual void | write (const std::string &fname) override |
void | verbose (bool set_verbosity) |
const std::vector< Real > & | get_time_steps () |
int | get_num_time_steps () |
void | copy_nodal_solution (System &system, std::string var_name, unsigned int timestep=1) |
void | copy_nodal_solution (System &system, std::string system_var_name, std::string exodus_var_name, unsigned int timestep=1) |
void | copy_elemental_solution (System &system, std::string system_var_name, std::string exodus_var_name, unsigned int timestep=1) |
void | read_elemental_variable (std::string elemental_var_name, unsigned int timestep, std::map< unsigned int, Real > &unique_id_to_value_map) |
void | read_global_variable (std::vector< std::string > global_var_names, unsigned int timestep, std::vector< Real > &global_values) |
void | write_discontinuous_exodusII (const std::string &name, const EquationSystems &es, const std::set< std::string > *system_names=nullptr) |
void | write_timestep_discontinuous (const std::string &fname, const EquationSystems &es, const int timestep, const Real time, const std::set< std::string > *system_names=nullptr) |
void | write_element_data (const EquationSystems &es) |
virtual void | write_nodal_data (const std::string &, const std::vector< Number > &, const std::vector< std::string > &) override |
void | write_nodal_data_discontinuous (const std::string &, const std::vector< Number > &, const std::vector< std::string > &) override |
void | write_global_data (const std::vector< Number > &, const std::vector< std::string > &) |
void | write_information_records (const std::vector< std::string > &) |
void | write_timestep (const std::string &fname, const EquationSystems &es, const int timestep, const Real time, const std::set< std::string > *system_names=nullptr) |
void | set_output_variables (const std::vector< std::string > &output_variables, bool allow_empty=true) |
void | use_mesh_dimension_instead_of_spatial_dimension (bool val) |
void | write_as_dimension (unsigned dim) |
void | set_coordinate_offset (Point p) |
void | append (bool val) |
const std::vector< std::string > & | get_elem_var_names () |
const std::vector< std::string > & | get_nodal_var_names () |
ExodusII_IO_Helper & | get_exio_helper () |
void | write_nodal_data_common (std::string fname, const std::vector< std::string > &names, bool continuous=true) |
virtual void | write_equation_systems (const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr) |
virtual void | write_discontinuous_equation_systems (const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr) |
virtual void | write_nodal_data (const std::string &, const NumericVector< Number > &, const std::vector< std::string > &) |
unsigned int & | ascii_precision () |
const Parallel::Communicator & | comm () const |
processor_id_type | n_processors () const |
processor_id_type | processor_id () const |
Protected Member Functions | |
MeshBase & | mesh () |
void | set_n_partitions (unsigned int n_parts) |
void | skip_comment_lines (std::istream &in, const char comment_start) |
const MeshBase & | mesh () const |
Protected Attributes | |
std::vector< bool > | elems_of_dimension |
const bool | _is_parallel_format |
const bool | _serial_only_needed_on_proc_0 |
const Parallel::Communicator & | _communicator |
Private Attributes | |
std::unique_ptr< ExodusII_IO_Helper > | exio_helper |
int | _timestep |
bool | _verbose |
bool | _append |
std::vector< std::string > | _output_variables |
bool | _allow_empty_variables |
Handles reading and writing of Exodus binary files.
The ExodusII_IO
class implements reading meshes in the ExodusII
file format from Sandia National Labs. By default, LibMesh expects ExodusII files to have a ".exd" or ".e" file extension.
Definition at line 52 of file exodusII_io.h.
|
explicit |
Constructor. Takes a writable reference to a mesh object. This is the constructor required to read a mesh.
Definition at line 45 of file exodusII_io.C.
|
virtual |
Destructor.
Definition at line 133 of file exodusII_io.C.
References exio_helper.
void libMesh::ExodusII_IO::append | ( | bool | val | ) |
If true, this flag will cause the ExodusII_IO object to attempt to open an existing file for writing, rather than creating a new file. Obviously this will only work if the file already exists.
Definition at line 450 of file exodusII_io.C.
References _append.
|
inlineinherited |
Return/set the precision to use when writing ASCII files.
By default we use numeric_limits<Real>::digits10 + 2, which should be enough to write out to ASCII and get the exact same Real back when reading in.
Definition at line 244 of file mesh_output.h.
Referenced by libMesh::TecplotIO::write_ascii(), libMesh::GMVIO::write_ascii_new_impl(), and libMesh::GMVIO::write_ascii_old_impl().
|
inlineinherited |
Parallel::Communicator
object used by this mesh. Definition at line 89 of file parallel_object.h.
References libMesh::ParallelObject::_communicator.
Referenced by libMesh::__libmesh_petsc_diff_solver_jacobian(), libMesh::__libmesh_petsc_diff_solver_monitor(), libMesh::__libmesh_petsc_diff_solver_residual(), libMesh::__libmesh_tao_equality_constraints(), libMesh::__libmesh_tao_equality_constraints_jacobian(), libMesh::__libmesh_tao_gradient(), libMesh::__libmesh_tao_hessian(), libMesh::__libmesh_tao_inequality_constraints(), libMesh::__libmesh_tao_inequality_constraints_jacobian(), libMesh::__libmesh_tao_objective(), libMesh::MeshRefinement::_coarsen_elements(), libMesh::ExactSolution::_compute_error(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::BoundaryInfo::_find_id_maps(), 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::MeshRefinement::_smooth_flags(), libMesh::PetscDMWrapper::add_dofs_helper(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::ImplicitSystem::add_matrix(), libMesh::System::add_vector(), libMesh::UnstructuredMesh::all_second_order(), libMesh::MeshTools::Modification::all_tri(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assemble_qoi(), libMesh::MeshCommunication::assign_global_indices(), libMesh::DofMap::attach_matrix(), libMesh::MeshTools::Generation::build_extrusion(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::PetscDMWrapper::build_section(), libMesh::PetscDMWrapper::build_sf(), libMesh::MeshBase::cache_elem_dims(), libMesh::System::calculate_norm(), libMesh::DofMap::check_dirichlet_bcid_consistency(), libMesh::PetscDMWrapper::check_section_n_dofs(), 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(), copy_elemental_solution(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::MeshTools::create_bounding_box(), libMesh::MeshTools::create_nodal_bounding_box(), libMesh::MeshRefinement::create_parent_error_vector(), libMesh::MeshTools::create_processor_bounding_box(), libMesh::MeshTools::create_subdomain_bounding_box(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), DMlibMeshFunction(), DMlibMeshJacobian(), DMlibMeshSetSystem_libMesh(), DMVariableBounds_libMesh(), libMesh::MeshRefinement::eliminate_unrefined_patches(), libMesh::EpetraVector< T >::EpetraVector(), 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::DofMap::get_info(), libMesh::ImplicitSystem::get_linear_solver(), libMesh::LocationMap< T >::init(), libMesh::TimeSolver::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_and_attach_petscdm(), libMesh::EigenSystem::init_data(), libMesh::EigenSystem::init_matrices(), libMesh::OptimizationSystem::initialize_equality_constraints_storage(), libMesh::OptimizationSystem::initialize_inequality_constraints_storage(), libMesh::MeshTools::libmesh_assert_consistent_distributed(), libMesh::MeshTools::libmesh_assert_consistent_distributed_nodes(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_new_node_procids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_flags(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_p_levels(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::MeshTools::libmesh_assert_valid_unique_ids(), libMesh::libmesh_petsc_snes_fd_residual(), libMesh::libmesh_petsc_snes_jacobian(), libMesh::libmesh_petsc_snes_mffd_residual(), libMesh::libmesh_petsc_snes_postcheck(), libMesh::libmesh_petsc_snes_residual(), libMesh::libmesh_petsc_snes_residual_helper(), libMesh::MeshRefinement::limit_level_mismatch_at_edge(), libMesh::MeshRefinement::limit_level_mismatch_at_node(), libMesh::MeshRefinement::limit_overrefined_boundary(), libMesh::MeshRefinement::limit_underrefined_boundary(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshCommunication::make_elems_parallel_consistent(), libMesh::MeshRefinement::make_flags_parallel_consistent(), libMesh::MeshCommunication::make_new_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_new_nodes_parallel_consistent(), libMesh::MeshCommunication::make_node_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_unique_ids_parallel_consistent(), libMesh::MeshCommunication::make_nodes_parallel_consistent(), libMesh::MeshCommunication::make_p_levels_parallel_consistent(), libMesh::MeshRefinement::make_refinement_compatible(), libMesh::FEMSystem::mesh_position_set(), libMesh::DistributedMesh::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::BoundaryInfo::n_shellface_conds(), libMesh::DistributedMesh::parallel_max_elem_id(), libMesh::DistributedMesh::parallel_max_node_id(), libMesh::ReplicatedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_n_elem(), libMesh::DistributedMesh::parallel_n_nodes(), libMesh::SparsityPattern::Build::parallel_sync(), libMesh::MeshTools::paranoid_n_levels(), libMesh::petsc_auto_fieldsplit(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshBase::prepare_for_use(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), 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_helper(), 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::DistributedMesh::renumber_dof_objects(), libMesh::CheckpointIO::select_split_config(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::PetscDiffSolver::setup_petsc_data(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::split_mesh(), libMesh::MeshBase::subdomain_ids(), libMesh::BoundaryInfo::sync(), libMesh::MeshRefinement::test_level_one(), libMesh::MeshRefinement::test_unflagged(), libMesh::MeshTools::total_weight(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().
void libMesh::ExodusII_IO::copy_elemental_solution | ( | System & | system, |
std::string | system_var_name, | ||
std::string | exodus_var_name, | ||
unsigned int | timestep = 1 |
||
) |
If we read in a elemental solution while reading in a mesh, we can attempt to copy that elemental solution into an EquationSystems object.
Definition at line 513 of file exodusII_io.C.
References libMesh::ParallelObject::comm(), libMesh::CONSTANT, libMesh::DofObject::dof_number(), end, exio_helper, libMesh::MeshInput< MT >::mesh(), libMesh::MONOMIAL, libMesh::DofObject::n_comp(), libMesh::System::number(), libMesh::Parallel::Communicator::rank(), libMesh::System::solution, libMesh::System::update(), libMesh::System::variable_number(), and libMesh::System::variable_type().
void libMesh::ExodusII_IO::copy_nodal_solution | ( | System & | system, |
std::string | var_name, | ||
unsigned int | timestep = 1 |
||
) |
Backward compatibility version of function that takes a single variable name.
Definition at line 78 of file exodusII_io.C.
void libMesh::ExodusII_IO::copy_nodal_solution | ( | System & | system, |
std::string | system_var_name, | ||
std::string | exodus_var_name, | ||
unsigned int | timestep = 1 |
||
) |
If we read in a nodal solution while reading in a mesh, we can attempt to copy that nodal solution into an EquationSystems object.
Definition at line 479 of file exodusII_io.C.
References libMesh::DofObject::dof_number(), exio_helper, libMesh::MeshInput< MT >::mesh(), libMesh::DofObject::n_comp(), libMesh::System::number(), libMesh::System::solution, libMesh::System::update(), and libMesh::System::variable_number().
const std::vector< std::string > & libMesh::ExodusII_IO::get_elem_var_names | ( | ) |
Return list of the elemental variable names
Definition at line 1032 of file exodusII_io.C.
References libMesh::ExodusII_IO_Helper::ELEMENTAL, and exio_helper.
ExodusII_IO_Helper & libMesh::ExodusII_IO::get_exio_helper | ( | ) |
Return a reference to the ExodusII_IO_Helper object.
Definition at line 1038 of file exodusII_io.C.
References exio_helper.
const std::vector< std::string > & libMesh::ExodusII_IO::get_nodal_var_names | ( | ) |
Return list of the nodal variable names
Definition at line 1026 of file exodusII_io.C.
References exio_helper, and libMesh::ExodusII_IO_Helper::NODAL.
int libMesh::ExodusII_IO::get_num_time_steps | ( | ) |
Knowing the number of time steps currently stored in the file is sometimes necessary when appending, so we can know where to start writing new data. Throws an error if the file is not currently open for reading or writing.
Definition at line 468 of file exodusII_io.C.
References exio_helper.
const std::vector< Real > & libMesh::ExodusII_IO::get_time_steps | ( | ) |
Definition at line 457 of file exodusII_io.C.
References exio_helper.
|
inlineprotectedinherited |
Definition at line 169 of file mesh_input.h.
Referenced by libMesh::GMVIO::_read_one_cell(), libMesh::VTKIO::cells_to_vtk(), libMesh::TetGenIO::element_in(), libMesh::UNVIO::elements_in(), libMesh::UNVIO::elements_out(), libMesh::UNVIO::groups_in(), libMesh::TetGenIO::node_in(), libMesh::UNVIO::nodes_in(), libMesh::UNVIO::nodes_out(), libMesh::VTKIO::nodes_to_vtk(), libMesh::Nemesis_IO::prepare_to_write_nodal_data(), libMesh::Nemesis_IO::read(), read(), libMesh::GMVIO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read(), libMesh::VTKIO::read(), libMesh::CheckpointIO::read_bcs(), libMesh::CheckpointIO::read_connectivity(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::UCDIO::read_implementation(), libMesh::UNVIO::read_implementation(), libMesh::GmshIO::read_mesh(), libMesh::CheckpointIO::read_nodes(), libMesh::CheckpointIO::read_nodesets(), libMesh::CheckpointIO::read_remote_elem(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::OFFIO::read_stream(), libMesh::MatlabIO::read_stream(), libMesh::CheckpointIO::read_subdomain_names(), libMesh::TetGenIO::write(), libMesh::Nemesis_IO::write(), write(), libMesh::XdrIO::write(), libMesh::CheckpointIO::write(), libMesh::GMVIO::write_ascii_new_impl(), libMesh::GMVIO::write_ascii_old_impl(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::Nemesis_IO::write_element_data(), write_element_data(), libMesh::UCDIO::write_header(), libMesh::UCDIO::write_implementation(), libMesh::UCDIO::write_interior_elems(), libMesh::GmshIO::write_mesh(), libMesh::VTKIO::write_nodal_data(), libMesh::UCDIO::write_nodal_data(), write_nodal_data(), write_nodal_data_common(), write_nodal_data_discontinuous(), libMesh::UCDIO::write_nodes(), libMesh::CheckpointIO::write_nodesets(), libMesh::XdrIO::write_parallel(), libMesh::GmshIO::write_post(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), libMesh::UCDIO::write_soln(), and libMesh::CheckpointIO::write_subdomain_names().
|
inlineprotectedinherited |
Definition at line 234 of file mesh_output.h.
Referenced by libMesh::FroIO::write(), libMesh::TecplotIO::write(), libMesh::MEDITIO::write(), libMesh::PostscriptIO::write(), libMesh::EnsightIO::write(), libMesh::TecplotIO::write_ascii(), libMesh::TecplotIO::write_binary(), libMesh::TecplotIO::write_nodal_data(), libMesh::MEDITIO::write_nodal_data(), and libMesh::GnuPlotIO::write_solution().
|
inlineinherited |
Definition at line 95 of file parallel_object.h.
References libMesh::ParallelObject::_communicator, and libMesh::Parallel::Communicator::size().
Referenced by libMesh::BoundaryInfo::_find_id_maps(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::DistributedMesh::add_elem(), libMesh::DistributedMesh::add_node(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assembly(), libMesh::AztecLinearSolver< T >::AztecLinearSolver(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::DistributedMesh::clear(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::UnstructuredMesh::create_pid_mesh(), libMesh::MeshTools::create_processor_bounding_box(), 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::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_and_attach_petscdm(), libMesh::Nemesis_IO_Helper::initialize(), libMesh::DistributedMesh::insert_elem(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_new_node_procids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_elem_on_proc(), libMesh::MeshBase::n_elem_on_proc(), libMesh::MeshBase::n_nodes_on_proc(), libMesh::MeshBase::partition(), libMesh::PetscLinearSolver< T >::PetscLinearSolver(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::CheckpointIO::read_connectivity(), libMesh::XdrIO::read_header(), libMesh::CheckpointIO::read_nodes(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::System::read_serialized_vector(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().
|
inlineinherited |
Definition at line 101 of file parallel_object.h.
References libMesh::ParallelObject::_communicator, and libMesh::Parallel::Communicator::rank().
Referenced by libMesh::BoundaryInfo::_find_id_maps(), libMesh::EquationSystems::_read_impl(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::DistributedMesh::add_elem(), libMesh::BoundaryInfo::add_elements(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::DistributedMesh::add_node(), libMesh::UnstructuredMesh::all_second_order(), libMesh::MeshTools::Modification::all_tri(), libMesh::FEMSystem::assembly(), libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::InfElemBuilder::build_inf_elem(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::DistributedMesh::clear(), libMesh::ExodusII_IO_Helper::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(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::ExodusII_IO_Helper::create(), libMesh::DistributedMesh::delete_elem(), libMesh::DistributedMesh::delete_node(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::DistributedMesh::DistributedMesh(), libMesh::DofMap::end_dof(), libMesh::DofMap::end_old_dof(), libMesh::EnsightIO::EnsightIO(), libMesh::MeshFunction::find_element(), libMesh::MeshFunction::find_elements(), 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::DofMap::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::SparsityPattern::Build::handle_vi_vj(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::ExodusII_IO_Helper::initialize_element_variables(), libMesh::ExodusII_IO_Helper::initialize_global_variables(), libMesh::ExodusII_IO_Helper::initialize_nodal_variables(), libMesh::DistributedMesh::insert_elem(), libMesh::DofMap::is_evaluable(), libMesh::SparsityPattern::Build::join(), libMesh::DofMap::last_dof(), libMesh::MeshTools::libmesh_assert_consistent_distributed(), libMesh::MeshTools::libmesh_assert_consistent_distributed_nodes(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), 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::BoundaryInfo::n_shellface_conds(), libMesh::SparsityPattern::Build::operator()(), libMesh::DistributedMesh::own_node(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), 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::XdrIO::read(), libMesh::CheckpointIO::read(), libMesh::ExodusII_IO_Helper::read_elem_num_map(), libMesh::ExodusII_IO_Helper::read_global_values(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::ExodusII_IO_Helper::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_helper(), 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::DistributedMesh::renumber_dof_objects(), libMesh::CheckpointIO::select_split_config(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::MeshTools::total_weight(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::Parallel::Packing< T >::unpack(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::CheckpointIO::write(), libMesh::EquationSystems::write(), libMesh::GMVIO::write_discontinuous_gmv(), write_element_data(), libMesh::ExodusII_IO_Helper::write_element_values(), libMesh::ExodusII_IO_Helper::write_elements(), write_global_data(), libMesh::ExodusII_IO_Helper::write_global_values(), libMesh::System::write_header(), write_information_records(), libMesh::ExodusII_IO_Helper::write_information_records(), libMesh::ExodusII_IO_Helper::write_nodal_coordinates(), libMesh::UCDIO::write_nodal_data(), write_nodal_data(), write_nodal_data_discontinuous(), libMesh::ExodusII_IO_Helper::write_nodal_values(), libMesh::Nemesis_IO_Helper::write_nodesets(), libMesh::ExodusII_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_helper(), 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(), libMesh::Nemesis_IO_Helper::write_sidesets(), libMesh::ExodusII_IO_Helper::write_sidesets(), write_timestep(), libMesh::ExodusII_IO_Helper::write_timestep(), and write_timestep_discontinuous().
|
overridevirtual |
This method implements reading a mesh from a specified file. Open the file named name
and read the mesh in Sandia National Lab's ExodusII format. This is the method to use for reading in meshes generated by cubit. Works in 2D for TRIs
, TRI6s
, QUAD
s, and QUAD9s
. Works in 3D for TET4s
, TET10s
, HEX8s
, and HEX27s
.
Implements libMesh::MeshInput< MeshBase >.
Definition at line 140 of file exodusII_io.C.
References libMesh::MeshBase::add_elem(), libMesh::BoundaryInfo::add_node(), libMesh::MeshBase::add_point(), libMesh::BoundaryInfo::add_shellface(), libMesh::BoundaryInfo::add_side(), libMesh::ExodusII_IO_Helper::ElementMaps::assign_conversion(), libMesh::Elem::build(), libMesh::MeshBase::clear(), libMesh::Elem::dim(), libMesh::MeshBase::elem_ref(), libMesh::MeshInput< MeshBase >::elems_of_dimension, libMesh::Utility::enum_to_string(), exio_helper, libMesh::MeshBase::get_boundary_info(), libMesh::ExodusII_IO_Helper::Conversion::get_shellface_index_offset(), libMesh::ExodusII_IO_Helper::Conversion::get_side_map(), libMesh::DofObject::id(), libMesh::index_range(), libMesh::ExodusII_IO_Helper::Conversion::invalid_id, libMesh::MeshInput< MT >::mesh(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::mesh_dimension(), libMesh::MeshBase::node_ptr(), libMesh::BoundaryInfo::nodeset_name(), libMesh::MeshBase::reserve_elem(), libMesh::MeshBase::reserve_nodes(), libMesh::DofObject::set_id(), libMesh::MeshBase::set_mesh_dimension(), libMesh::Elem::set_node(), libMesh::BoundaryInfo::sideset_name(), libMesh::Elem::subdomain_id(), and libMesh::MeshBase::subdomain_name().
Referenced by libMesh::NameBasedIO::read(), and libMesh::Nemesis_IO::read().
void libMesh::ExodusII_IO::read_elemental_variable | ( | std::string | elemental_var_name, |
unsigned int | timestep, | ||
std::map< unsigned int, Real > & | unique_id_to_value_map | ||
) |
Given an elemental variable and a time step, returns a mapping from the elements (top parent) unique IDs to the value of the elemental variable at the corresponding time step index. Note that this function MUST only be called before renumbering! This function is essentially a wrapper for read_elemental_var_values from the exodus helper (which is not accessible outside this class).
elemental_var_name | Name of an elemental variable |
timestep | The corresponding time step index |
unique_id_to_value_map | The map to be filled |
Definition at line 555 of file exodusII_io.C.
References exio_helper, libMesh::MeshInput< MT >::mesh(), libMesh::Elem::top_parent(), and libMesh::DofObject::unique_id().
void libMesh::ExodusII_IO::read_global_variable | ( | std::vector< std::string > | global_var_names, |
unsigned int | timestep, | ||
std::vector< Real > & | global_values | ||
) |
Given a vector of global variables and a time step, returns the values of the global variable at the corresponding time step index.
global_var_names | Vector of names of global variables |
timestep | The corresponding time step index |
global_values | The vector to be filled |
Definition at line 570 of file exodusII_io.C.
References exio_helper, and libMesh::ExodusII_IO_Helper::GLOBAL.
void libMesh::ExodusII_IO::set_coordinate_offset | ( | Point | p | ) |
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 442 of file exodusII_io.C.
References exio_helper.
|
inlineprotectedinherited |
Sets the number of partitions in the mesh. Typically this gets done by the partitioner, but some parallel file formats begin "pre-partitioned".
Definition at line 91 of file mesh_input.h.
References libMesh::MeshInput< MT >::mesh().
Referenced by libMesh::Nemesis_IO::read(), and libMesh::XdrIO::read_header().
void libMesh::ExodusII_IO::set_output_variables | ( | const std::vector< std::string > & | output_variables, |
bool | allow_empty = true |
||
) |
Sets the list of variable names to be included in the output. This is optional. If this is never called then all variables will be present. If this is called and an empty vector is supplied no variables will be output. Setting the allow_empty = false will result in empty vectors supplied here to also be populated with all variables.
Definition at line 68 of file exodusII_io.C.
References _allow_empty_variables, and _output_variables.
|
protectedinherited |
Reads input from in
, skipping all the lines that start with the character comment_start
.
Definition at line 179 of file mesh_input.h.
Referenced by libMesh::TetGenIO::read(), and libMesh::UCDIO::read_implementation().
void libMesh::ExodusII_IO::use_mesh_dimension_instead_of_spatial_dimension | ( | bool | val | ) |
In the general case, meshes containing 2D elements can be manifolds living in 3D space, thus by default we write all meshes with the Exodus dimension set to LIBMESH_DIM = mesh.spatial_dimension().
In certain cases, however, the user may know his 2D mesh actually lives in the z=0 plane, and therefore wants to write a truly 2D Exodus mesh. In such a case, he should call this function with val=true.
Definition at line 428 of file exodusII_io.C.
References exio_helper.
void libMesh::ExodusII_IO::verbose | ( | bool | set_verbosity | ) |
Set the flag indicating if we should be verbose.
Definition at line 418 of file exodusII_io.C.
References _verbose, and exio_helper.
|
overridevirtual |
This method implements writing a mesh to a specified file.
Implements libMesh::MeshOutput< MeshBase >.
Definition at line 879 of file exodusII_io.C.
References _append, _verbose, exio_helper, libMesh::MeshBase::get_boundary_info(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), and libMesh::BoundaryInfo::n_edge_conds().
Referenced by libMesh::ErrorVector::plot_error(), and libMesh::NameBasedIO::write().
void libMesh::ExodusII_IO::write_as_dimension | ( | unsigned | dim | ) |
Directly control the num_dim which is written to the Exodus file. If non-zero, this value supersedes all other dimensions, including: 1.) MeshBase::spatial_dimension() 2.) MeshBase::mesh_dimension() 3.) Any value passed to use_mesh_dimension_instead_of_spatial_dimension() This is useful/necessary for working around a bug in Paraview which prevents the "Plot Over Line" filter from working on 1D meshes.
Definition at line 435 of file exodusII_io.C.
References exio_helper.
|
virtualinherited |
This method implements writing a mesh with discontinuous data to a specified file where the data is taken from the EquationSystems
object.
Definition at line 92 of file mesh_output.C.
References libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::EquationSystems::build_variable_names(), libMesh::EquationSystems::get_mesh(), and libMesh::out.
Referenced by write_timestep_discontinuous().
void libMesh::ExodusII_IO::write_discontinuous_exodusII | ( | const std::string & | name, |
const EquationSystems & | es, | ||
const std::set< std::string > * | system_names = nullptr |
||
) |
Writes a exodusII file with discontinuous data
Definition at line 89 of file exodusII_io.C.
References libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::EquationSystems::build_variable_names(), and write_nodal_data_discontinuous().
void libMesh::ExodusII_IO::write_element_data | ( | const EquationSystems & | es | ) |
Write out element solution.
Definition at line 599 of file exodusII_io.C.
References _output_variables, _timestep, std::abs(), libMesh::EquationSystems::build_elemental_solution_vector(), libMesh::EquationSystems::build_variable_names(), libMesh::CONSTANT, exio_helper, libMesh::EquationSystems::get_vars_active_subdomains(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), libMesh::MONOMIAL, libMesh::ParallelObject::processor_id(), and value.
Referenced by libMesh::ErrorVector::plot_error().
|
virtualinherited |
This method implements writing a mesh with data to a specified file where the data is taken from the EquationSystems
object.
Reimplemented in libMesh::NameBasedIO.
Definition at line 31 of file mesh_output.C.
References libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::EquationSystems::build_solution_vector(), libMesh::EquationSystems::build_variable_names(), libMesh::EquationSystems::get_mesh(), and libMesh::out.
Referenced by libMesh::Nemesis_IO::write_timestep(), and write_timestep().
void libMesh::ExodusII_IO::write_global_data | ( | const std::vector< Number > & | soln, |
const std::vector< std::string > & | names | ||
) |
Write out global variables.
Definition at line 809 of file exodusII_io.C.
References _timestep, std::abs(), exio_helper, libMesh::ParallelObject::processor_id(), and value.
void libMesh::ExodusII_IO::write_information_records | ( | const std::vector< std::string > & | records | ) |
Write out information records.
Definition at line 796 of file exodusII_io.C.
References exio_helper, and libMesh::ParallelObject::processor_id().
|
virtualinherited |
This method should be overridden by "parallel" output formats for writing nodal data. Instead of getting a localized copy of the nodal solution vector, it is passed a NumericVector of type=PARALLEL which is in node-major order i.e. (u0,v0,w0, u1,v1,w1, u2,v2,w2, u3,v3,w3, ...) and contains n_nodes*n_vars total entries. Then, it is up to the individual I/O class to extract the required solution values from this vector and write them in parallel.
If not implemented, localizes the parallel vector into a std::vector and calls the other version of this function.
Reimplemented in libMesh::Nemesis_IO.
Definition at line 150 of file mesh_output.C.
References libMesh::NumericVector< T >::localize().
|
overridevirtual |
Write out a nodal solution.
Reimplemented from libMesh::MeshOutput< MeshBase >.
Definition at line 700 of file exodusII_io.C.
References _allow_empty_variables, _output_variables, _timestep, std::abs(), exio_helper, libMesh::MeshTools::Generation::Private::idx(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), libMesh::MeshBase::n_nodes(), libMesh::MeshBase::node_ptr_range(), libMesh::ParallelObject::processor_id(), and write_nodal_data_common().
Referenced by libMesh::NameBasedIO::write_nodal_data().
void libMesh::ExodusII_IO::write_nodal_data_common | ( | std::string | fname, |
const std::vector< std::string > & | names, | ||
bool | continuous = true |
||
) |
This function factors out a bunch of code which is common to the write_nodal_data() and write_nodal_data_discontinuous() functions
Definition at line 970 of file exodusII_io.C.
References _append, _verbose, exio_helper, libMesh::MeshBase::get_boundary_info(), libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), and libMesh::BoundaryInfo::n_edge_conds().
Referenced by write_nodal_data(), and write_nodal_data_discontinuous().
|
overridevirtual |
Write out a discontinuous nodal solution.
Reimplemented from libMesh::MeshOutput< MeshBase >.
Definition at line 912 of file exodusII_io.C.
References _timestep, std::abs(), libMesh::MeshBase::active_element_ptr_range(), exio_helper, libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshOutput< MT >::mesh(), libMesh::ParallelObject::processor_id(), and write_nodal_data_common().
Referenced by write_discontinuous_exodusII().
void libMesh::ExodusII_IO::write_timestep | ( | const std::string & | fname, |
const EquationSystems & | es, | ||
const int | timestep, | ||
const Real | time, | ||
const std::set< std::string > * | system_names = nullptr |
||
) |
Writes out the solution at a specific timestep.
fname | Name of the file to write to |
es | EquationSystems object which contains the solution vector. |
timestep | The timestep to write out, should be 1 indexed. |
time | The current simulation time. |
system_names | Optional list of systems to write solutions for. |
Definition at line 862 of file exodusII_io.C.
References _timestep, exio_helper, libMesh::ParallelObject::processor_id(), and libMesh::MeshOutput< MeshBase >::write_equation_systems().
void libMesh::ExodusII_IO::write_timestep_discontinuous | ( | const std::string & | fname, |
const EquationSystems & | es, | ||
const int | timestep, | ||
const Real | time, | ||
const std::set< std::string > * | system_names = nullptr |
||
) |
Writes a discontinuous solution at a specific timestep
fname | Name of the file to be written |
es | EquationSystems object which contains the solution vector |
timestep | The timestep to write out. (should be 1 indexed) |
time | The current simulation time |
system_names | Optional list of systems to write solutions for. |
Definition at line 104 of file exodusII_io.C.
References _timestep, exio_helper, libMesh::ParallelObject::processor_id(), and libMesh::MeshOutput< MeshBase >::write_discontinuous_equation_systems().
|
private |
If true, _output_variables is allowed to remain empty. If false, if _output_variables is empty it will be populated with a complete list of all variables By default, calling set_output_variables() sets this flag to true, but it provides an override.
Definition at line 354 of file exodusII_io.h.
Referenced by set_output_variables(), and write_nodal_data().
|
private |
Default false. If true, files will be opened with EX_WRITE rather than created from scratch when writing.
Definition at line 340 of file exodusII_io.h.
Referenced by append(), write(), and write_nodal_data_common().
|
protectedinherited |
Definition at line 107 of file parallel_object.h.
Referenced by libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::ParallelObject::comm(), libMesh::ParallelObject::n_processors(), libMesh::ParallelObject::operator=(), and libMesh::ParallelObject::processor_id().
|
protectedinherited |
Flag specifying whether this format is parallel-capable. If this is false (default) I/O is only permitted when the mesh has been serialized.
Definition at line 159 of file mesh_output.h.
Referenced by libMesh::FroIO::write(), libMesh::PostscriptIO::write(), and libMesh::EnsightIO::write().
|
private |
The names of the variables to be output. If this is empty then all variables are output.
Definition at line 347 of file exodusII_io.h.
Referenced by set_output_variables(), write_element_data(), and write_nodal_data().
|
protectedinherited |
Flag specifying whether this format can be written by only serializing the mesh to processor zero
If this is false (default) the mesh will be serialized to all processors
Definition at line 168 of file mesh_output.h.
|
private |
Stores the current value of the timestep when calling ExodusII_IO::write_timestep().
Definition at line 329 of file exodusII_io.h.
Referenced by write_element_data(), write_global_data(), write_nodal_data(), write_nodal_data_discontinuous(), write_timestep(), and write_timestep_discontinuous().
|
private |
should we be verbose?
Definition at line 334 of file exodusII_io.h.
Referenced by verbose(), write(), and write_nodal_data_common().
|
protectedinherited |
A vector of bools describing what dimension elements have been encountered when reading a mesh.
Definition at line 97 of file mesh_input.h.
Referenced by libMesh::GMVIO::_read_one_cell(), libMesh::UNVIO::elements_in(), libMesh::UNVIO::max_elem_dimension_seen(), libMesh::AbaqusIO::max_elem_dimension_seen(), libMesh::AbaqusIO::read(), libMesh::Nemesis_IO::read(), read(), libMesh::GMVIO::read(), libMesh::VTKIO::read(), libMesh::AbaqusIO::read_elements(), libMesh::UCDIO::read_implementation(), libMesh::UNVIO::read_implementation(), and libMesh::XdrIO::read_serialized_connectivity().
|
private |
Only attempt to instantiate an ExodusII helper class if the Exodus API is defined. This class will have no functionality when LIBMESH_HAVE_EXODUS_API is not defined.
Definition at line 323 of file exodusII_io.h.
Referenced by copy_elemental_solution(), copy_nodal_solution(), get_elem_var_names(), get_exio_helper(), get_nodal_var_names(), get_num_time_steps(), get_time_steps(), read(), read_elemental_variable(), read_global_variable(), set_coordinate_offset(), use_mesh_dimension_instead_of_spatial_dimension(), verbose(), write(), write_as_dimension(), write_element_data(), write_global_data(), write_information_records(), write_nodal_data(), write_nodal_data_common(), write_nodal_data_discontinuous(), write_timestep(), write_timestep_discontinuous(), and ~ExodusII_IO().