#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 | read_qa_records () |
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, std::map< dof_id_type, Real > &elem_var_value_map) |
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) |
virtual void | initialize_element_variables (std::vector< std::string > names, const std::vector< std::set< subdomain_id_type >> &vars_active_subdomains) |
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, const std::vector< std::set< subdomain_id_type >> &vars_active_subdomains) |
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 | read_global_values (std::vector< Real > &values, int timestep) |
void | use_mesh_dimension_instead_of_spatial_dimension (bool val) |
void | write_as_dimension (unsigned dim) |
void | set_coordinate_offset (Point p) |
std::vector< std::string > | get_complex_names (const std::vector< std::string > &names) const |
std::vector< std::set< subdomain_id_type > > | get_complex_vars_active_subdomains (const std::vector< std::set< subdomain_id_type >> &vars_active_subdomains) 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 Member Functions | |
void | check_existing_vars (ExodusVarType type, std::vector< std::string > &names, std::vector< std::string > &names_from_file) |
void | write_var_names (ExodusVarType type, std::vector< std::string > &names) |
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 |
unsigned | _write_as_dimension |
Point | _coordinate_offset |
bool | _single_precision |
const Parallel::Communicator & | _communicator |
Private Member Functions | |
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 96 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
Enumerator | |
---|---|
NODAL | |
ELEMENTAL | |
GLOBAL |
Definition at line 601 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 270 of file exodusII_io_helper.C.
References elem_type, and title.
|
virtual |
|
protected |
When appending: during initialization, check that variable names in the file match those you attempt to initialize with.
Definition at line 1814 of file exodusII_io_helper.C.
References libMesh::err, libMesh::Quality::name(), libMesh::out, and read_var_names().
Referenced by initialize_element_variables(), libMesh::Nemesis_IO_Helper::initialize_element_variables(), initialize_global_variables(), and initialize_nodal_variables().
void libMesh::ExodusII_IO_Helper::close | ( | ) |
Closes the ExodusII
mesh file.
Definition at line 795 of file exodusII_io_helper.C.
References _run_only_on_proc0, ex_close(), ex_err, ex_id, message(), opened_for_reading, opened_for_writing, and libMesh::ParallelObject::processor_id().
Referenced by libMesh::Nemesis_IO_Helper::~Nemesis_IO_Helper().
|
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(), libMesh::ExodusII_IO::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().
|
virtual |
Opens an ExodusII
mesh file named filename
for writing.
Reimplemented in libMesh::Nemesis_IO_Helper.
Definition at line 1079 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.
int libMesh::ExodusII_IO_Helper::get_block_id | ( | int | index | ) |
Get the block number for the given block index.
Definition at line 544 of file exodusII_io_helper.C.
References block_ids.
Referenced by 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 553 of file exodusII_io_helper.C.
References block_ids, and id_to_block_names.
std::vector< std::string > libMesh::ExodusII_IO_Helper::get_complex_names | ( | const std::vector< std::string > & | names | ) | const |
Definition at line 2064 of file exodusII_io_helper.C.
std::vector< std::set< subdomain_id_type > > libMesh::ExodusII_IO_Helper::get_complex_vars_active_subdomains | ( | const std::vector< std::set< subdomain_id_type >> & | vars_active_subdomains | ) | const |
returns a "tripled" copy of vars_active_subdomains
, which is necessary in the complex-valued case.
Definition at line 2090 of file exodusII_io_helper.C.
const char * libMesh::ExodusII_IO_Helper::get_elem_type | ( | ) | const |
HEX27
. Definition at line 314 of file exodusII_io_helper.C.
References elem_type.
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 580 of file exodusII_io_helper.C.
References nodeset_ids.
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 589 of file exodusII_io_helper.C.
References id_to_ns_names, and nodeset_ids.
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 562 of file exodusII_io_helper.C.
References ss_ids.
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 571 of file exodusII_io_helper.C.
References id_to_ss_names, and ss_ids.
|
virtual |
Initializes the Exodus file.
Reimplemented in libMesh::Nemesis_IO_Helper.
Definition at line 1119 of file exodusII_io_helper.C.
References _run_only_on_proc0, _use_mesh_dimension_instead_of_spatial_dimension, _write_as_dimension, libMesh::MeshBase::active_element_ptr_range(), libMesh::BoundaryInfo::build_node_boundary_ids(), libMesh::BoundaryInfo::build_shellface_boundary_ids(), libMesh::BoundaryInfo::build_side_boundary_ids(), libMesh::err, ex_err, ex_id, ex_put_init(), libMesh::MeshBase::get_boundary_info(), mesh, libMesh::MeshBase::mesh_dimension(), libMesh::MeshBase::n_active_elem(), libMesh::MeshBase::n_elem(), libMesh::MeshBase::nodes_begin(), libMesh::MeshBase::nodes_end(), num_dim, num_elem, num_elem_blk, num_node_sets, num_nodes, num_side_sets, libMesh::ParallelObject::processor_id(), and libMesh::MeshBase::spatial_dimension().
|
virtual |
Sets up the nodal variables
Reimplemented in libMesh::Nemesis_IO_Helper.
Definition at line 1695 of file exodusII_io_helper.C.
References _elem_vars_initialized, _run_only_on_proc0, block_ids, check_existing_vars(), elem_var_names, ELEMENTAL, ex_err, ex_id, ex_put_elem_var_tab(), libMesh::index_range(), num_elem_blk, num_elem_vars, libMesh::ParallelObject::processor_id(), and write_var_names().
void libMesh::ExodusII_IO_Helper::initialize_global_variables | ( | std::vector< std::string > | names | ) |
Sets up the global variables
Definition at line 1788 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().
void libMesh::ExodusII_IO_Helper::initialize_nodal_variables | ( | std::vector< std::string > | names | ) |
Sets up the nodal variables
Definition at line 1760 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().
int libMesh::ExodusII_IO_Helper::inquire | ( | int | req_info, |
std::string | error_msg = "" |
||
) |
Definition at line 814 of file exodusII_io_helper.C.
References ex_err, ex_id, and ex_inquire().
Referenced by read_num_time_steps(), read_qa_records(), 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 321 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().
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 328 of file exodusII_io_helper.C.
References libMesh::out, and verbose.
|
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().
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 335 of file exodusII_io_helper.C.
References _single_precision, current_filename, ex_id, std::min(), opened_for_reading, opened_for_writing, libMesh::out, libMesh::Real, and verbose.
void libMesh::ExodusII_IO_Helper::print_header | ( | ) |
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 455 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.
void libMesh::ExodusII_IO_Helper::print_nodes | ( | std::ostream & | out = libMesh::out | ) |
Prints the nodal information, by default to libMesh::out
.
Definition at line 509 of file exodusII_io_helper.C.
References num_nodes, x, y, and z.
|
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(), 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(), 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(), initialize(), initialize_element_variables(), initialize_global_variables(), 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(), read_elem_num_map(), read_global_values(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), 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_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(), 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(), libMesh::Nemesis_IO_Helper::write_nodesets(), 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(), write_sidesets(), libMesh::ExodusII_IO::write_timestep(), write_timestep(), and libMesh::ExodusII_IO::write_timestep_discontinuous().
void libMesh::ExodusII_IO_Helper::read_block_info | ( | ) |
Reads information for all of the blocks in the ExodusII
mesh file.
Definition at line 517 of file exodusII_io_helper.C.
References block_ids, EX_ELEM_BLOCK, ex_err, ex_get_elem_blk_ids(), ex_get_name(), ex_id, id_to_block_names, message(), and num_elem_blk.
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 599 of file exodusII_io_helper.C.
References block_ids, connect, elem_type, ex_err, ex_get_elem_block(), ex_get_elem_conn(), ex_id, message(), num_attr, num_elem_this_blk, num_nodes_per_elem, libMesh::out, and verbose.
void libMesh::ExodusII_IO_Helper::read_elem_num_map | ( | ) |
Reads the optional node_num_map
from the ExodusII
mesh file.
Definition at line 639 of file exodusII_io_helper.C.
References elem_num_map, ex_err, ex_get_elem_num_map(), ex_id, message(), std::min(), num_elem, libMesh::out, libMesh::ParallelObject::processor_id(), and verbose.
void libMesh::ExodusII_IO_Helper::read_elemental_var_values | ( | std::string | elemental_var_name, |
int | time_step, | ||
std::map< dof_id_type, Real > & | elem_var_value_map | ||
) |
Reads elemental values for the variable 'elemental_var_name' at the specified timestep into the 'elem_var_value_map' which is passed in.
Definition at line 1012 of file exodusII_io_helper.C.
References block_ids, elem_num_map, elem_var_names, ELEMENTAL, libMesh::err, ex_err, ex_get_elem_block(), ex_get_elem_var(), ex_id, num_elem, num_elem_blk, num_elem_this_blk, and read_var_names().
void libMesh::ExodusII_IO_Helper::read_global_values | ( | std::vector< Real > & | values, |
int | timestep | ||
) |
Reads the vector of global variables.
Definition at line 2031 of file exodusII_io_helper.C.
References _run_only_on_proc0, ex_err, ex_get_glob_vars(), ex_id, num_global_vars, and libMesh::ParallelObject::processor_id().
void libMesh::ExodusII_IO_Helper::read_header | ( | ) |
Reads an ExodusII
mesh file header.
Definition at line 375 of file exodusII_io_helper.C.
References ex_err, ex_get_init(), ex_get_var_param(), 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.
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 856 of file exodusII_io_helper.C.
References libMesh::err, ex_err, ex_get_nodal_var(), ex_id, NODAL, nodal_var_names, nodal_var_values, num_nodes, and read_var_names().
void libMesh::ExodusII_IO_Helper::read_node_num_map | ( | ) |
Reads the optional node_num_map
from the ExodusII
mesh file.
Definition at line 489 of file exodusII_io_helper.C.
References ex_err, ex_get_node_num_map(), ex_id, message(), std::min(), node_num_map, num_nodes, libMesh::out, libMesh::ParallelObject::processor_id(), and verbose.
void libMesh::ExodusII_IO_Helper::read_nodes | ( | ) |
Reads the nodal data (x,y,z coordinates) from the ExodusII
mesh file.
Definition at line 469 of file exodusII_io_helper.C.
References ex_err, ex_get_coord(), ex_id, message(), num_nodes, x, y, and z.
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 765 of file exodusII_io_helper.C.
References ex_err, ex_get_node_set(), ex_get_node_set_param(), ex_id, message(), node_list, nodeset_ids, num_node_df_per_set, and num_nodes_per_set.
void libMesh::ExodusII_IO_Helper::read_nodeset_info | ( | ) |
Reads information about all of the nodesets in the ExodusII
mesh file.
Definition at line 694 of file exodusII_io_helper.C.
References ex_err, ex_get_name(), ex_get_node_set_ids(), ex_id, EX_NODE_SET, id_to_ns_names, message(), nodeset_ids, num_node_df_per_set, num_node_sets, and num_nodes_per_set.
void libMesh::ExodusII_IO_Helper::read_num_time_steps | ( | ) |
Reads the number of timesteps currently stored in the Exodus file and stores it in the num_time_steps variable.
Definition at line 848 of file exodusII_io_helper.C.
References EX_INQ_TIME, inquire(), and num_time_steps.
Referenced by read_header(), and read_time_steps().
void libMesh::ExodusII_IO_Helper::read_qa_records | ( | ) |
Reads the QA records from an ExodusII file. We can use this to detect when e.g. CUBIT 14+ was used to generate a Mesh file, and work around certain known bugs in that version.
Definition at line 405 of file exodusII_io_helper.C.
References ex_err, ex_get_qa(), ex_id, EX_INQ_QA, inquire(), libMesh::out, and verbose.
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 722 of file exodusII_io_helper.C.
References elem_list, ex_err, ex_get_side_set(), ex_get_side_set_param(), ex_id, id_list, message(), num_df_per_set, num_sides_per_set, side_list, and ss_ids.
void libMesh::ExodusII_IO_Helper::read_sideset_info | ( | ) |
Reads information about all of the sidesets in the ExodusII
mesh file.
Definition at line 661 of file exodusII_io_helper.C.
References elem_list, ex_err, ex_get_name(), ex_get_side_set_ids(), ex_id, EX_INQ_SS_ELEM_LEN, EX_SIDE_SET, 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.
void libMesh::ExodusII_IO_Helper::read_time_steps | ( | ) |
Reads and stores the timesteps in the 'time_steps' array.
Definition at line 833 of file exodusII_io_helper.C.
References ex_err, ex_get_all_times(), ex_id, num_time_steps, read_num_time_steps(), and time_steps.
void libMesh::ExodusII_IO_Helper::read_var_names | ( | ExodusVarType | type | ) |
Definition at line 896 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(), read_elemental_var_values(), and read_nodal_var_values().
|
private |
read_var_names() dispatches to this function.
Definition at line 916 of file exodusII_io_helper.C.
References ex_err, ex_get_var_names(), ex_get_var_param(), 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().
void libMesh::ExodusII_IO_Helper::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 2058 of file exodusII_io_helper.C.
References _coordinate_offset.
void libMesh::ExodusII_IO_Helper::use_mesh_dimension_instead_of_spatial_dimension | ( | bool | val | ) |
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 2044 of file exodusII_io_helper.C.
References _use_mesh_dimension_instead_of_spatial_dimension.
void libMesh::ExodusII_IO_Helper::write_as_dimension | ( | unsigned | dim | ) |
Sets the value of _write_as_dimension.
This directly controls 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 2051 of file exodusII_io_helper.C.
References _write_as_dimension.
void libMesh::ExodusII_IO_Helper::write_element_values | ( | const MeshBase & | mesh, |
const std::vector< Real > & | values, | ||
int | timestep, | ||
const std::vector< std::set< subdomain_id_type >> & | vars_active_subdomains | ||
) |
Writes the vector of values to the element variables.
Definition at line 1865 of file exodusII_io_helper.C.
References _run_only_on_proc0, _single_precision, libMesh::MeshBase::active_element_ptr_range(), data, ex_err, ex_get_var_param(), ex_id, ex_put_elem_var(), ex_update(), get_block_id(), mesh, libMesh::MeshBase::n_elem(), libMesh::MeshTools::n_elem(), num_elem_vars, and libMesh::ParallelObject::processor_id().
|
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 1324 of file exodusII_io_helper.C.
References _run_only_on_proc0, libMesh::MeshBase::active_element_ptr_range(), libMesh::ExodusII_IO_Helper::ElementMaps::assign_conversion(), block_ids, connect, elem_num_map, libMesh::MeshBase::elem_ref(), libMesh::Utility::enum_to_string(), EX_ELEM_BLOCK, ex_err, ex_id, ex_put_concat_elem_block(), ex_put_elem_conn(), ex_put_elem_num_map(), ex_put_names(), 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::index_range(), libmesh_elem_num_to_exodus, libmesh_node_num_to_exodus, mesh, libMesh::MeshBase::n_active_elem(), libMesh::Elem::n_nodes(), libMesh::Elem::node_id(), num_elem_blk, num_nodes_per_elem, libMesh::out, libMesh::ParallelObject::processor_id(), libMesh::ExodusII_IO_Helper::NamesData::push_back_entry(), libMesh::MeshBase::subdomain_name(), libMesh::Elem::type(), and verbose.
void libMesh::ExodusII_IO_Helper::write_global_values | ( | const std::vector< Real > & | values, |
int | timestep | ||
) |
Writes the vector of global variables.
Definition at line 2009 of file exodusII_io_helper.C.
References _run_only_on_proc0, _single_precision, ex_err, ex_id, ex_put_glob_vars(), ex_update(), num_global_vars, and libMesh::ParallelObject::processor_id().
void libMesh::ExodusII_IO_Helper::write_information_records | ( | const std::vector< std::string > & | records | ) |
Writes the vector of information records.
Definition at line 1970 of file exodusII_io_helper.C.
References _run_only_on_proc0, libMesh::err, ex_err, ex_id, EX_INQ_INFO, ex_put_info(), ex_update(), libMesh::ExodusII_IO_Helper::NamesData::get_char_star_star(), inquire(), libMesh::ParallelObject::processor_id(), and libMesh::ExodusII_IO_Helper::NamesData::push_back_entry().
|
virtual |
Writes the nodal coordinates contained in "mesh"
Reimplemented in libMesh::Nemesis_IO_Helper.
Definition at line 1208 of file exodusII_io_helper.C.
References _coordinate_offset, _run_only_on_proc0, _single_precision, libMesh::MeshBase::active_element_ptr_range(), ex_err, ex_id, ex_put_coord(), ex_put_node_num_map(), libMesh::DofObject::id(), libmesh_node_num_to_exodus, mesh, node_num_map, libMesh::MeshBase::node_ptr_range(), num_nodes, libMesh::ParallelObject::processor_id(), x, y, and z.
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 1948 of file exodusII_io_helper.C.
References _run_only_on_proc0, _single_precision, ex_err, ex_id, ex_put_nodal_var(), ex_update(), num_nodes, and libMesh::ParallelObject::processor_id().
Referenced by libMesh::Nemesis_IO_Helper::write_nodal_solution().
|
virtual |
Writes the nodesets contained in "mesh"
Reimplemented in libMesh::Nemesis_IO_Helper.
Definition at line 1653 of file exodusII_io_helper.C.
References _run_only_on_proc0, libMesh::BoundaryInfo::build_node_boundary_ids(), libMesh::BoundaryInfo::build_node_list(), data, ex_err, ex_id, EX_NODE_SET, ex_put_names(), ex_put_node_set(), ex_put_node_set_param(), libMesh::MeshBase::get_boundary_info(), libMesh::BoundaryInfo::get_nodeset_name(), mesh, libMesh::ParallelObject::processor_id(), and libMesh::ExodusII_IO_Helper::NamesData::push_back_entry().
|
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
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 1544 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_shellface_boundary_ids(), libMesh::BoundaryInfo::build_shellface_list(), libMesh::BoundaryInfo::build_side_boundary_ids(), libMesh::BoundaryInfo::build_side_list(), libMesh::MeshBase::elem_ptr(), libMesh::MeshBase::elem_ref(), ex_err, ex_id, ex_put_names(), ex_put_sets(), EX_SIDE_SET, libMesh::MeshBase::get_boundary_info(), libMesh::ExodusII_IO_Helper::Conversion::get_inverse_shellface_map(), libMesh::ExodusII_IO_Helper::Conversion::get_inverse_side_map(), libMesh::BoundaryInfo::get_sideset_name(), libMesh::index_range(), libmesh_elem_num_to_exodus, mesh, libMesh::ParallelObject::processor_id(), libMesh::ExodusII_IO_Helper::NamesData::push_back_entry(), and side.
void libMesh::ExodusII_IO_Helper::write_timestep | ( | int | timestep, |
Real | time | ||
) |
Writes the time for the timestep
Definition at line 1843 of file exodusII_io_helper.C.
References _run_only_on_proc0, _single_precision, ex_err, ex_id, ex_put_time(), ex_update(), and libMesh::ParallelObject::processor_id().
|
protected |
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 956 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(), libMesh::Nemesis_IO_Helper::initialize_element_variables(), initialize_global_variables(), and initialize_nodal_variables().
|
private |
write_var_names() dispatches to this function.
Definition at line 976 of file exodusII_io_helper.C.
References ex_err, ex_id, ex_put_var_names(), ex_put_var_param(), libMesh::ExodusII_IO_Helper::NamesData::get_char_star(), libMesh::ExodusII_IO_Helper::NamesData::get_char_star_star(), libMesh::out, libMesh::ExodusII_IO_Helper::NamesData::push_back_entry(), and verbose.
Referenced by write_var_names().
|
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().
|
protected |
Definition at line 640 of file exodusII_io_helper.h.
Referenced by set_coordinate_offset(), and write_nodal_coordinates().
|
protected |
Definition at line 622 of file exodusII_io_helper.h.
Referenced by initialize_element_variables(), and libMesh::Nemesis_IO_Helper::initialize_element_variables().
|
protected |
Definition at line 625 of file exodusII_io_helper.h.
Referenced by initialize_global_variables().
|
protected |
Definition at line 628 of file exodusII_io_helper.h.
Referenced by initialize_nodal_variables().
|
protected |
Definition at line 619 of file exodusII_io_helper.h.
Referenced by close(), create(), initialize(), initialize_element_variables(), initialize_global_variables(), initialize_nodal_variables(), read_global_values(), write_element_values(), write_elements(), write_global_values(), write_information_records(), write_nodal_coordinates(), write_nodal_values(), write_nodesets(), write_sidesets(), and write_timestep().
|
protected |
Definition at line 643 of file exodusII_io_helper.h.
Referenced by create(), libMesh::Nemesis_IO_Helper::create(), open(), write_element_values(), write_global_values(), libMesh::Nemesis_IO_Helper::write_nodal_coordinates(), write_nodal_coordinates(), write_nodal_values(), and write_timestep().
|
protected |
Definition at line 633 of file exodusII_io_helper.h.
Referenced by initialize(), and use_mesh_dimension_instead_of_spatial_dimension().
|
protected |
Definition at line 637 of file exodusII_io_helper.h.
Referenced by initialize(), and write_as_dimension().
std::vector<int> libMesh::ExodusII_IO_Helper::block_ids |
Definition at line 474 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), get_block_id(), get_block_name(), initialize_element_variables(), 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 477 of file exodusII_io_helper.h.
Referenced by read_elem_in_block(), and write_elements().
std::string libMesh::ExodusII_IO_Helper::current_filename |
Definition at line 591 of file exodusII_io_helper.h.
std::vector<int> libMesh::ExodusII_IO_Helper::elem_list |
Definition at line 498 of file exodusII_io_helper.h.
Referenced by read_sideset(), and read_sideset_info().
std::vector<int> libMesh::ExodusII_IO_Helper::elem_num_map |
Definition at line 513 of file exodusII_io_helper.h.
Referenced by read_elem_num_map(), read_elemental_var_values(), and write_elements().
std::vector<char> libMesh::ExodusII_IO_Helper::elem_type |
Definition at line 528 of file exodusII_io_helper.h.
Referenced by ExodusII_IO_Helper(), get_elem_type(), and read_elem_in_block().
std::vector<std::string> libMesh::ExodusII_IO_Helper::elem_var_names |
Definition at line 559 of file exodusII_io_helper.h.
Referenced by initialize_element_variables(), libMesh::Nemesis_IO_Helper::initialize_element_variables(), read_elemental_var_values(), and read_var_names().
std::vector<Real> libMesh::ExodusII_IO_Helper::elem_var_values |
Definition at line 562 of file exodusII_io_helper.h.
int libMesh::ExodusII_IO_Helper::ex_err |
Definition at line 438 of file exodusII_io_helper.h.
Referenced by close(), initialize(), initialize_element_variables(), libMesh::Nemesis_IO_Helper::initialize_element_variables(), inquire(), read_block_info(), read_elem_in_block(), read_elem_num_map(), read_elemental_var_values(), read_global_values(), read_header(), read_nodal_var_values(), read_node_num_map(), read_nodes(), read_nodeset(), read_nodeset_info(), read_qa_records(), read_sideset(), read_sideset_info(), read_time_steps(), read_var_names_impl(), write_element_values(), libMesh::Nemesis_IO_Helper::write_element_values(), libMesh::Nemesis_IO_Helper::write_elements(), write_elements(), libMesh::Nemesis_IO_Helper::write_exodus_initialization_info(), write_global_values(), write_information_records(), libMesh::Nemesis_IO_Helper::write_nodal_coordinates(), write_nodal_coordinates(), write_nodal_values(), libMesh::Nemesis_IO_Helper::write_nodesets(), write_nodesets(), libMesh::Nemesis_IO_Helper::write_sidesets(), write_sidesets(), write_timestep(), write_var_names_impl(), and libMesh::Nemesis_IO_Helper::~Nemesis_IO_Helper().
int libMesh::ExodusII_IO_Helper::ex_id |
Definition at line 435 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(), libMesh::Nemesis_IO_Helper::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_global_values(), read_header(), read_nodal_var_values(), read_node_num_map(), read_nodes(), read_nodeset(), read_nodeset_info(), read_qa_records(), read_sideset(), read_sideset_info(), read_time_steps(), read_var_names_impl(), write_element_values(), libMesh::Nemesis_IO_Helper::write_element_values(), libMesh::Nemesis_IO_Helper::write_elements(), write_elements(), libMesh::Nemesis_IO_Helper::write_exodus_initialization_info(), write_global_values(), write_information_records(), libMesh::Nemesis_IO_Helper::write_nodal_coordinates(), write_nodal_coordinates(), write_nodal_values(), libMesh::Nemesis_IO_Helper::write_nodesets(), write_nodesets(), libMesh::Nemesis_IO_Helper::write_sidesets(), 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 533 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 538 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::Nemesis_IO_Helper::write_nodal_coordinates(), and libMesh::Nemesis_IO_Helper::write_nodal_solution().
std::vector<std::string> libMesh::ExodusII_IO_Helper::global_var_names |
Definition at line 565 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 507 of file exodusII_io_helper.h.
Referenced by read_sideset(), and read_sideset_info().
std::map<int, std::string> libMesh::ExodusII_IO_Helper::id_to_block_names |
Definition at line 568 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 570 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 569 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 532 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(), libMesh::Nemesis_IO_Helper::write_sidesets(), and write_sidesets().
std::map<int, int> libMesh::ExodusII_IO_Helper::libmesh_node_num_to_exodus |
Definition at line 537 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::Nemesis_IO_Helper::compute_node_communication_maps(), libMesh::Nemesis_IO_Helper::compute_node_maps(), write_elements(), write_nodal_coordinates(), and libMesh::Nemesis_IO_Helper::write_nodesets().
std::vector<std::string> libMesh::ExodusII_IO_Helper::nodal_var_names |
Definition at line 550 of file exodusII_io_helper.h.
Referenced by initialize_nodal_variables(), read_nodal_var_values(), and read_var_names().
std::vector<Real> libMesh::ExodusII_IO_Helper::nodal_var_values |
Definition at line 553 of file exodusII_io_helper.h.
Referenced by read_nodal_var_values().
std::vector<int> libMesh::ExodusII_IO_Helper::node_list |
Definition at line 504 of file exodusII_io_helper.h.
Referenced by read_nodeset().
std::vector<int> libMesh::ExodusII_IO_Helper::node_num_map |
Definition at line 510 of file exodusII_io_helper.h.
Referenced by read_node_num_map(), and write_nodal_coordinates().
std::vector<int> libMesh::ExodusII_IO_Helper::nodeset_ids |
Definition at line 483 of file exodusII_io_helper.h.
Referenced by get_node_set_id(), get_node_set_name(), read_nodeset(), and read_nodeset_info().
int libMesh::ExodusII_IO_Helper::num_attr |
Definition at line 468 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 492 of file exodusII_io_helper.h.
Referenced by read_sideset(), and read_sideset_info().
int libMesh::ExodusII_IO_Helper::num_dim |
Definition at line 441 of file exodusII_io_helper.h.
Referenced by initialize(), print_header(), read_header(), and libMesh::Nemesis_IO_Helper::write_exodus_initialization_info().
int libMesh::ExodusII_IO_Helper::num_elem |
Definition at line 450 of file exodusII_io_helper.h.
Referenced by initialize(), print_header(), read_elem_num_map(), read_elemental_var_values(), read_header(), and libMesh::Nemesis_IO_Helper::write_exodus_initialization_info().
int libMesh::ExodusII_IO_Helper::num_elem_all_sidesets |
Definition at line 471 of file exodusII_io_helper.h.
Referenced by read_sideset_info().
int libMesh::ExodusII_IO_Helper::num_elem_blk |
Definition at line 453 of file exodusII_io_helper.h.
Referenced by initialize(), initialize_element_variables(), print_header(), read_block_info(), read_elemental_var_values(), read_header(), write_elements(), and libMesh::Nemesis_IO_Helper::write_exodus_initialization_info().
int libMesh::ExodusII_IO_Helper::num_elem_this_blk |
Definition at line 462 of file exodusII_io_helper.h.
Referenced by read_elem_in_block(), and read_elemental_var_values().
int libMesh::ExodusII_IO_Helper::num_elem_vars |
Definition at line 556 of file exodusII_io_helper.h.
Referenced by initialize_element_variables(), libMesh::Nemesis_IO_Helper::initialize_element_variables(), read_header(), read_var_names(), write_element_values(), and write_var_names().
int libMesh::ExodusII_IO_Helper::num_global_vars |
Definition at line 444 of file exodusII_io_helper.h.
Referenced by initialize_global_variables(), read_global_values(), read_header(), read_var_names(), write_global_values(), and write_var_names().
int libMesh::ExodusII_IO_Helper::num_nodal_vars |
Definition at line 547 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 495 of file exodusII_io_helper.h.
Referenced by read_nodeset(), and read_nodeset_info().
int libMesh::ExodusII_IO_Helper::num_node_sets |
Definition at line 456 of file exodusII_io_helper.h.
Referenced by initialize(), print_header(), read_header(), read_nodeset_info(), and libMesh::Nemesis_IO_Helper::write_exodus_initialization_info().
int libMesh::ExodusII_IO_Helper::num_nodes |
Definition at line 447 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), initialize(), print_header(), print_nodes(), read_header(), read_nodal_var_values(), read_node_num_map(), read_nodes(), libMesh::Nemesis_IO_Helper::write_exodus_initialization_info(), write_nodal_coordinates(), libMesh::Nemesis_IO_Helper::write_nodal_solution(), and write_nodal_values().
int libMesh::ExodusII_IO_Helper::num_nodes_per_elem |
Definition at line 465 of file exodusII_io_helper.h.
Referenced by libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), read_elem_in_block(), libMesh::Nemesis_IO_Helper::write_elements(), and write_elements().
std::vector<int> libMesh::ExodusII_IO_Helper::num_nodes_per_set |
Definition at line 489 of file exodusII_io_helper.h.
Referenced by read_nodeset(), and read_nodeset_info().
int libMesh::ExodusII_IO_Helper::num_side_sets |
Definition at line 459 of file exodusII_io_helper.h.
Referenced by initialize(), print_header(), 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 486 of file exodusII_io_helper.h.
Referenced by read_sideset(), and read_sideset_info().
int libMesh::ExodusII_IO_Helper::num_time_steps |
Definition at line 541 of file exodusII_io_helper.h.
Referenced by read_num_time_steps(), and read_time_steps().
bool libMesh::ExodusII_IO_Helper::opened_for_reading |
Definition at line 581 of file exodusII_io_helper.h.
bool libMesh::ExodusII_IO_Helper::opened_for_writing |
Definition at line 577 of file exodusII_io_helper.h.
Referenced by close(), create(), libMesh::Nemesis_IO_Helper::create(), open(), and libMesh::Nemesis_IO_Helper::~Nemesis_IO_Helper().
std::vector<int> libMesh::ExodusII_IO_Helper::side_list |
Definition at line 501 of file exodusII_io_helper.h.
Referenced by read_sideset(), and read_sideset_info().
std::vector<int> libMesh::ExodusII_IO_Helper::ss_ids |
Definition at line 480 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 544 of file exodusII_io_helper.h.
Referenced by read_time_steps().
std::vector<char> libMesh::ExodusII_IO_Helper::title |
Definition at line 525 of file exodusII_io_helper.h.
Referenced by ExodusII_IO_Helper(), print_header(), and read_header().
bool libMesh::ExodusII_IO_Helper::verbose |
Definition at line 573 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_qa_records(), read_var_names_impl(), 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 516 of file exodusII_io_helper.h.
Referenced by print_nodes(), read_nodes(), libMesh::Nemesis_IO_Helper::write_nodal_coordinates(), and write_nodal_coordinates().
std::vector<Real> libMesh::ExodusII_IO_Helper::y |
Definition at line 519 of file exodusII_io_helper.h.
Referenced by print_nodes(), read_nodes(), libMesh::Nemesis_IO_Helper::write_nodal_coordinates(), and write_nodal_coordinates().
std::vector<Real> libMesh::ExodusII_IO_Helper::z |
Definition at line 522 of file exodusII_io_helper.h.
Referenced by print_nodes(), read_nodes(), libMesh::Nemesis_IO_Helper::write_nodal_coordinates(), and write_nodal_coordinates().