#include <Range.hpp>
Classes | |
class | const_iterator |
a const iterator which iterates over an Range More... | |
class | const_pair_iterator |
class | const_reverse_iterator |
a const reverse iterator which iterates over an Range More... | |
class | pair_iterator |
used to iterate over sub-ranges of a range More... | |
struct | PairNode |
Public Types | |
typedef const_iterator | iterator |
typedef const_reverse_iterator | reverse_iterator |
typedef EntityHandle | value_type |
for short hand notation, lets typedef the container class that holds the ranges More... | |
Public Member Functions | |
Range & | operator-= (const Range &rhs) |
just like subtract, but as an operator More... | |
Range () | |
default constructor More... | |
Range (const Range ©) | |
copy constructor More... | |
Range (EntityHandle val1, EntityHandle val2) | |
another constructor that takes an initial range More... | |
Range & | operator= (const Range ©) |
operator= More... | |
~Range () | |
destructor More... | |
const_iterator | begin () const |
return the beginning const iterator of this range More... | |
const_reverse_iterator | rbegin () const |
return the beginning const reverse iterator of this range More... | |
const_iterator | end () const |
return the ending const iterator for this range More... | |
const_reverse_iterator | rend () const |
return the ending const reverse iterator for this range More... | |
size_t | size () const |
return the number of values this Ranges represents More... | |
size_t | psize () const |
return the number of range pairs in the list More... | |
bool | empty () const |
return whether empty or not always use "if(!Ranges::empty())" instead of "if(Ranges::size())" More... | |
iterator | insert (iterator hint, EntityHandle val) |
iterator | insert (EntityHandle val) |
insert an item into the list and return the iterator for the inserted item More... | |
iterator | insert (iterator hint, EntityHandle first, EntityHandle last) |
insert a range of items into this list and return the iterator for the first inserted item More... | |
iterator | insert (EntityHandle val1, EntityHandle val2) |
insert a range of items into this list and return the iterator for the first inserted item More... | |
template<typename T > | |
iterator | insert_list (T begin_iter, T end_iter) |
template<class T > | |
iterator | insert (typename T::const_iterator begin_iter, typename T::const_iterator end_iter) |
template<typename T > | |
iterator | insert (const T *begin_iter, const T *end_iter) |
iterator | erase (iterator iter) |
remove an item from this list and return an iterator to the next item More... | |
iterator | erase (iterator iter1, iterator iter2) |
remove a range of items from the list More... | |
iterator | erase (EntityHandle val) |
erases a value from this container More... | |
const EntityHandle & | front () const |
get first entity in range More... | |
const EntityHandle & | back () const |
get last entity in range More... | |
EntityHandle | pop_front () |
remove first entity from range More... | |
EntityHandle | pop_back () |
remove last entity from range More... | |
const_iterator | find (EntityHandle val) const |
find an item int the list and return an iterator at that value More... | |
const_iterator | lower_bound (EntityHandle val) const |
const_iterator | upper_bound (EntityHandle val) const |
const_iterator | lower_bound (EntityType type) const |
const_iterator | upper_bound (EntityType type) const |
std::pair< const_iterator, const_iterator > | equal_range (EntityType type) const |
const_iterator | lower_bound (EntityType type, const_iterator first) const |
const_iterator | upper_bound (EntityType type, const_iterator first) const |
bool | all_of_type (EntityType type) const |
True if all entities in range are of passed type (also true if range is empty) More... | |
bool | all_of_dimension (int dimension) const |
True if all entities in range are of passed dimension (also true if range is empty) More... | |
unsigned | num_of_type (EntityType type) const |
unsigned | num_of_dimension (int dim) const |
void | clear () |
clears the contents of the list More... | |
const std::string | str_rep (const char *indent_prefix=NULL) const |
for debugging More... | |
void | print (const char *indent_prefix=NULL) const |
void | print (std::ostream &s, const char *indent_prefix=NULL) const |
unsigned long | get_memory_use () const |
double | compactness () const |
void | insert (Range::const_iterator begin, Range::const_iterator end) |
void | merge (const Range &range) |
merges this Range with another range More... | |
void | merge (Range::const_iterator beginr, Range::const_iterator endr) |
merge a subset of some other range More... | |
void | swap (Range &range) |
swap the contents of this range with another one More... | |
void | sanity_check () const |
check for internal consistency More... | |
bool | contains (const Range &other) const |
Check if this range is a non-strict superset of some other range. More... | |
Range | subset_by_type (EntityType t) const |
return a subset of this range, by type More... | |
Range | subset_by_dimension (int dim) const |
return a subset of this range, by dimension More... | |
EntityHandle | operator[] (EntityID index) const |
int | index (EntityHandle handle) const |
pair_iterator | pair_begin () |
pair_iterator | pair_end () |
const_pair_iterator | const_pair_begin () const |
const_pair_iterator | const_pair_end () const |
const_pair_iterator | pair_begin () const |
const_pair_iterator | pair_end () const |
template<typename Iterator > | |
Range::iterator | insert_list (Iterator begin_iter, Iterator end_iter) |
Static Public Member Functions | |
static const_iterator | lower_bound (const_iterator first, const_iterator last, EntityHandle val) |
return an iterator to the first value >= val More... | |
static const_iterator | upper_bound (const_iterator first, const_iterator last, EntityHandle val) |
Protected Member Functions | |
void | delete_pair_node (PairNode *dead_node) |
if dead_node is not mHead, remove it from the list and free it's memory. More... | |
Protected Attributes | |
PairNode | mHead |
the head of the list that contains pairs that represent the ranges this list is sorted and unique at all times More... | |
Friends | |
Range | intersect (const Range &, const Range &) |
intersect two ranges, placing the results in the return range More... | |
Range | subtract (const Range &, const Range &) |
subtract range2 from this, placing the results in the return range More... | |
the class Range
typedef const_iterator moab::Range::iterator |
typedef EntityHandle moab::Range::value_type |
|
inline |
moab::Range::Range | ( | const Range & | copy | ) |
copy constructor
Definition at line 179 of file Range.cpp.
References alloc_pair(), mHead, moab::Range::PairNode::mNext, and moab::Range::PairNode::mPrev.
moab::Range::Range | ( | EntityHandle | val1, |
EntityHandle | val2 | ||
) |
another constructor that takes an initial range
Definition at line 172 of file Range.cpp.
References alloc_pair(), mHead, moab::Range::PairNode::mNext, and moab::Range::PairNode::mPrev.
|
inline |
bool moab::Range::all_of_dimension | ( | int | dimension | ) | const |
True if all entities in range are of passed dimension (also true if range is empty)
Definition at line 884 of file Range.cpp.
References back(), moab::CN::Dimension(), empty(), front(), and moab::TYPE_FROM_HANDLE().
Referenced by moab::Skinner::find_skin(), moab::Skinner::find_skin_vertices_1D(), moab::Skinner::find_skin_vertices_2D(), moab::Skinner::find_skin_vertices_3D(), and moab::Core::get_vertices().
bool moab::Range::all_of_type | ( | EntityType | type | ) | const |
True if all entities in range are of passed type (also true if range is empty)
Definition at line 879 of file Range.cpp.
References back(), empty(), front(), and moab::TYPE_FROM_HANDLE().
Referenced by moab::WriteNCDF::gather_mesh_information(), moab::ReorderTool::handle_order_from_sets_and_adj(), iMOAB_GetBlockInfo(), moab::ReorderTool::int_order_from_sets_and_adj(), moab::OrientedBoxTreeTool::join_trees(), main(), moab::GeomQueryTool::measure_area(), moab::GeomQueryTool::measure_volume(), and moab::WriteNCDF::write_elementblocks().
|
inline |
get last entity in range
Definition at line 912 of file Range.hpp.
References mHead, and moab::Range::PairNode::mPrev.
Referenced by all_of_dimension(), all_of_type(), moab::WriteNCDF::gather_mesh_information(), moab::WriteUtil::gather_nodes_from_elements(), moab::ReadHDF5::get_partition(), moab::AdaptiveKDTree::merge_leaf(), moab::BSPTree::merge_leaf(), pop_back(), moab::ReadHDF5::read_elems(), moab::ReadHDF5::read_sparse_tag_indices(), and moab::GeomTopoTool::resize_rootSets().
|
inline |
return the beginning const iterator of this range
Definition at line 860 of file Range.hpp.
References mHead, and moab::Range::PairNode::mNext.
Referenced by moab::GeomTopoTool::A_is_in_B(), moab::SpatialLocator::add_elems(), moab::ParallelComm::add_verts(), moab::SharedSetData::append_local_handles(), moab::SmoothFace::area(), moab::IntxAreaUtils::area_on_sphere(), moab::WriteSTL::ascii_write_triangles(), MetisPartitioner::assemble_graph(), ZoltanPartitioner::assemble_graph(), MetisPartitioner::assemble_taggedents_graph(), MetisPartitioner::assemble_taggedsets_graph(), moab::ScdInterface::assign_global_ids(), moab::WriteUtil::assign_ids(), moab::DualTool::atomic_pillow(), moab::WriteSTL::binary_write_triangles(), moab::FBEngine::boundary_mesh_edges_on_face(), moab::FBEngine::boundary_nodes_on_face(), moab::box_from_axes(), moab::MeshGeneration::BrickInstance(), SphereDecomp::build_hexes(), moab::FBEngine::chain_able_edge(), moab::FBEngine::chain_two_edges(), moab::Core::check_adjacencies(), moab::ParallelComm::check_clean_iface(), moab::DualTool::check_dual_adjs(), moab::DualTool::check_dual_equiv_edges(), moab::AEntityFactory::check_equiv_entities(), moab::NCHelperMPAS::check_existing_mesh(), moab::GeomTopoTool::check_model(), moab::ParallelComm::check_sent_ents(), moab::Skinner::classify_2d_boundary(), moab::Intx2Mesh::clean(), moab::ParallelComm::clean_shared_tags(), moab::Core::clear_meshset(), moab::OrientedBoxTreeTool::closest_to_location(), moab::closest_to_triangles(), moab::HalfFacetRep::collect_and_compare(), moab::MeshTopoUtil::common_entity(), moab::WriteHDF5Parallel::communicate_shared_set_data(), moab::WriteHDF5Parallel::communicate_shared_set_ids(), moab::HiReconstruction::compute_average_vertex_normals_surf(), moab::HiReconstruction::compute_average_vertex_tangents_curve(), moab::SmoothFace::compute_control_points_on_edges(), moab::NestedRefine::compute_coordinates(), compute_dual_mesh(), moab::SmoothFace::compute_internal_control_points_on_facets(), SphereDecomp::compute_nodes(), moab::ParCommGraph::compute_partition(), moab::SmoothFace::compute_tangents_for_each_edge(), moab::TempestOnlineMap::ComputeAdjacencyRelations(), moab::TempestRemapper::ComputeOverlapMesh(), moab::DualTool::construct_dual(), moab::DualTool::construct_dual_cells(), moab::DualTool::construct_dual_edges(), moab::DualTool::construct_dual_faces(), moab::DualTool::construct_dual_vertices(), moab::BVHTree::construct_element_vec(), moab::NestedRefine::construct_hm_1D(), moab::NestedRefine::construct_hm_2D(), moab::NestedRefine::construct_hm_3D(), moab::DualTool::construct_hp_parent_child(), moab::GeomTopoTool::construct_obb_tree(), moab::GeomTopoTool::construct_obb_trees(), moab::GeomTopoTool::construct_vertex_ranges(), moab::TempestRemapper::convert_mesh_to_tempest_private(), moab::Tqdcfr::convert_nodesets_sidesets(), moab::TempestRemapper::convert_overlap_mesh_sorted_by_source(), moab::ReadHDF5::convert_range_to_handle(), moab::HigherOrderFactory::convert_sequence(), moab::SpectralMeshTool::convert_to_coarse(), moab::copy_sorted_file_ids(), moab::WriteHDF5::count_adjacencies(), moab::WriteHDF5::count_set_size(), moab::ReadCCMIO::create_cell_from_faces(), moab::AEntityFactory::create_explicit_adjs(), moab::NCHelperGCRM::create_gather_set_vertices(), moab::NCHelperMPAS::create_gather_set_vertices(), moab::NestedRefine::create_hm_storage_single_level(), moab::ParallelComm::create_iface_pc_links(), moab::ReadABAQUS::create_instance_of_part(), moab::NCHelperESMF::create_local_cells(), moab::NCHelperMPAS::create_local_cells(), moab::NCHelperGCRM::create_local_edges(), moab::NCHelperMPAS::create_local_edges(), moab::NCHelperESMF::create_local_vertices(), moab::NCHelperGCRM::create_local_vertices(), moab::NCHelperMPAS::create_local_vertices(), moab::ScdNCHelper::create_mesh(), moab::NCHelperDomain::create_mesh(), moab::NCHelperHOMME::create_mesh(), moab::NCHelperScrip::create_mesh(), moab::NCHelperESMF::create_padded_local_cells(), moab::NCHelperGCRM::create_padded_local_cells(), moab::NCHelperMPAS::create_padded_local_cells(), moab::ScdNCHelper::create_quad_coordinate_tag(), moab::ReadGmsh::create_sets(), moab::SpectralMeshTool::create_spectral_elems(), moab::AEntityFactory::create_vert_elem_adjacencies(), moab::Intx2Mesh::createTags(), moab::IntxUtils::deep_copy_set_with_quads(), moab::DualTool::delete_dual_entities(), moab::ParallelComm::delete_entities(), moab::ReadHDF5::delete_non_side_elements(), moab::ReadParallel::delete_nonlocal_entities(), moab::GeomTopoTool::delete_obb_tree(), moab::HalfFacetRep::determine_border_vertices(), moab::HalfFacetRep::determine_incident_halfedges(), moab::HalfFacetRep::determine_incident_halffaces(), moab::HalfFacetRep::determine_incident_halfverts(), moab::HalfFacetRep::determine_sibling_halfedges(), moab::HalfFacetRep::determine_sibling_halffaces(), moab::HalfFacetRep::determine_sibling_halfverts(), moab::Intx2Mesh::DetermineOrderedNeighbors(), do_test_mode(), dot_children(), dot_contained(), dot_nodes(), dot_write_id_nodes(), moab::SmoothFace::DumpModelControlPoints(), moab::IntxUtils::enforce_convexity(), equal_range(), moab::NestedRefine::estimate_hm_storage(), moab::Core::estimated_memory_use(), moab::WriteHDF5Parallel::exchange_file_ids(), moab::DualTool::face_shrink(), moab::ParallelComm::filter_pstatus(), moab::Intx2Mesh::filterByMask(), moab::Tree::find_all_trees(), moab::ScdInterface::find_boxes(), moab::VarLenDenseTag::find_entities_with_value(), moab::DenseTag::find_entities_with_value(), moab::SparseTag::find_entities_with_value(), moab::VarLenSparseTag::find_entities_with_value(), moab::ParallelComm::find_existing_entity(), moab::Skinner::find_geometric_skin(), moab::Skinner::find_inferred_edges(), moab::find_map_values(), moab::HalfFacetRep::find_matching_halfedge(), moab::HalfFacetRep::find_matching_halfface(), moab::HalfFacetRep::find_matching_implicit_edge_in_cell(), moab::MergeMesh::find_merged_to(), moab::ReadHDF5::find_sets_containing(), moab::Skinner::find_skin_noadj(), moab::Skinner::find_skin_vertices_1D(), moab::Skinner::find_skin_vertices_2D(), moab::Skinner::find_skin_vertices_3D(), moab::find_tag_values(), moab::HalfFacetRep::find_total_edges_2d(), moab::HalfFacetRep::find_total_edges_faces_3d(), moab::FBEngine::find_vertex_set_for_node(), moab::GeomQueryTool::find_volume_slow(), moab::Intx2Mesh::FindMaxEdgesInSet(), moab::IntxUtils::fix_degenerate_quads(), moab::DualTool::foc_delete_dual(), moab::DualTool::foc_get_ents(), moab::ParCommGraph::form_tuples_to_migrate_mesh(), moab::DualTool::fs_get_quad_loops(), moab::DualTool::fs_get_quads(), moab::DualTool::fsr_get_fourth_quad(), moab::ParallelComm::gather_data(), moab::WriteSLAC::gather_interior_exterior(), moab::WriteVtk::gather_mesh(), moab::WriteHDF5::gather_mesh_info(), moab::WriteNCDF::gather_mesh_information(), moab::WriteSLAC::gather_mesh_information(), moab::WriteTemplate::gather_mesh_information(), moab::WriteUtil::gather_nodes_from_elements(), moab::ReadUtil::gather_related_ents(), gather_set_stats(), moab::GeomTopoTool::generate_implicit_complement(), moab::TempestRemapper::GenerateMeshMetadata(), moab::GeomTopoTool::geometrize_surface_set(), moab::ScdBox::get_adj_edge_or_face(), moab::Core::get_adjacencies(), moab::get_adjacencies_intersection(), moab::get_adjacencies_union(), moab::MeshTopoUtil::get_average_position(), moab::DualTool::get_cell_points(), moab::Core::get_connectivity(), moab::NestedRefine::get_connectivity(), moab::Core::get_connectivity_by_type(), moab::GeomTopoTool::get_ct_children_by_dimension(), moab::AEntityFactory::get_down_adjacency_elements_poly(), moab::DualTool::get_dual_entities(), moab::WriteUtil::get_element_connect(), moab::TypeSequenceManager::get_entities(), moab::MeshSet::get_entities_by_dimension(), moab::Core::get_entities_by_handle(), moab::MeshSet::get_entities_by_type(), moab::ParallelComm::get_ghosted_entities(), moab::Coupler::get_gl_points_on_elements(), moab::DualTool::get_graphics_points(), moab::HalfFacetRep::get_half_facet_in_comp(), moab::ParallelComm::get_iface_entities(), moab::ParallelComm::get_interface_procs(), moab::ParallelComm::get_interface_sets(), moab::ParallelData::get_interface_sets(), moab::NestedRefine::get_lid_inci_child(), moab::ParallelComm::get_local_handles(), moab::MeshTopoUtil::get_manifold(), get_max_volume(), moab::AEntityFactory::get_memory_use(), moab::WriteCCMIO::get_neuset_elems(), moab::WriteSLAC::get_neuset_elems(), moab::WriteTemplate::get_neuset_elems(), moab::ReorderTool::get_new_handles(), moab::MeshSet::get_non_set_entities(), moab::HiReconstruction::get_normals_surf(), moab::DualTool::get_opposite_verts(), moab::ParallelComm::get_part_entities(), moab::ReadHDF5::get_partition(), moab::ParallelComm::get_pstatus_entities(), moab::ParallelComm::get_remote_handles(), moab::ReorderTool::get_reordered_handles(), moab::ReadABAQUS::get_set_by_name(), moab::WriteCCMIO::get_sets(), moab::HalfFacetRep::get_sibling_map(), moab::WriteNCDF::get_sideset_elems(), get_signed_volume(), moab::WriteHDF5::get_tag_data_length(), moab::ReadHDF5::get_tagged_entities(), moab::HiReconstruction::get_tangents_curve(), moab::MeshSetSequence::get_type(), moab::HalfFacetRep::get_up_adjacencies_edg_3d(), moab::HalfFacetRep::get_up_adjacencies_face_3d(), moab::HalfFacetRep::get_up_adjacencies_vert_2d(), moab::HalfFacetRep::get_up_adjacencies_vert_3d(), moab::AEntityFactory::get_up_adjacency_elements(), moab::WriteSLAC::get_valid_sides(), moab::WriteTemplate::get_valid_sides(), moab::WriteNCDF::get_valid_sides(), moab::Core::get_vertex_coordinates(), moab::NestedRefine::get_vertex_duplicates(), moab::IntxUtils::global_gnomonic_projection(), hcFilter(), iMOAB_GetElementOwnership(), iMOAB_GetMeshInfo(), iMOAB_GetPointerToSurfaceBC(), iMOAB_GetPointerToVertexBC(), iMOAB_GetVertexOwnership(), iMOAB_GetVisibleElementsInfo(), iMOAB_SetDoubleTagStorageWithGid(), ZoltanPartitioner::include_closure(), index(), moab::SmoothFace::init_gradient(), moab::HalfFacetRep::init_surface(), moab::HalfFacetRep::init_volume(), moab::WriteNCDF::initialize_exodus_file(), moab::ReadHDF5::insert_in_id_map(), insert_list(), moab::ReorderTool::int_order_from_sets_and_adj(), moab::intersect(), moab::AdaptiveKDTree::intersect_children_with_elems(), moab::Intx2Mesh::intersect_meshes(), moab::Intx2Mesh::intersect_meshes_kdtree(), moab::DualTool::is_blind(), moab::NestedRefine::is_cell_on_boundary(), moab::NestedRefine::is_vertex_on_boundary(), moab::OrientedBoxTreeTool::join_trees(), laplacianFilter(), moab::RayIntersectSets::leaf(), moab::ParallelComm::list_entities(), moab::Core::list_entities(), moab::ReadParallel::load_file(), moab::ReadABAQUS::load_file(), moab::ReadSms::load_file_impl(), moab::WriteGMV::local_write_mesh(), moab::Coupler::locate_points(), lower_bound(), main(), make_tris_from_quads(), moab::WriteDamsel::map_sparse_tags(), moab::GeomQueryTool::measure_area(), moab::GeomQueryTool::measure_volume(), merge(), moab::AEntityFactory::merge_adjust_adjacencies(), merge_duplicate_vertices(), moab::MergeMesh::merge_higher_dimensions(), merge_input_surfs(), moab::MergeMesh::merge_using_integer_tag(), moab::Coupler::nat_param(), moab::DualTool::next_loop_vertex(), obbstat_write(), obbvis_create(), operator-=(), operator[](), moab::DualTool::order_chord(), moab::GeomTopoTool::other_entity(), moab::ParallelComm::pack_entities(), moab::ParallelComm::pack_entity_seq(), moab::ParallelComm::pack_sets(), ZoltanPartitioner::partition_owned_cells(), perform_laplacian_smoothing(), perform_lloyd_relaxation(), moab::LloydSmoother::perform_smooth(), moab::GeomQueryTool::point_in_volume_slow(), moab::IntxAreaUtils::positive_orientation(), moab::AdaptiveKDTree::print(), moab::Core::print(), moab::TreeNodePrinter::print_contents(), moab::WriteHDF5Parallel::print_set_sharing_data(), moab::HalfFacetRep::print_tags(), moab::print_type_sets(), moab::ReadDamsel::process_entity_tags(), quads_to_tris(), moab::AdaptiveKDTree::ray_intersect_triangles(), moab::OrientedBoxTreeTool::ray_intersect_triangles(), moab::ReadHDF5VarLen::read_data(), moab::ReadHDF5::read_dense_tag(), McnpData::read_mcnpfile(), moab::ReadHDF5::read_node_adj_elems(), moab::Tqdcfr::read_nodes(), moab::ReadNCDF::read_nodesets(), moab::ReadHDF5VarLen::read_offsets(), moab::ScdNCHelper::read_scd_variables_to_nonset_allocate(), moab::ReadHDF5::read_set_data(), moab::ReadHDF5::read_sets(), moab::ReadNCDF::read_sidesets(), moab::ReadHDF5::read_sparse_tag(), moab::ReadHDF5::read_sparse_tag_indices(), moab::ReadHDF5::read_tag_values_partial(), moab::NCHelperMPAS::read_ucd_variables_to_nonset(), moab::NCHelperGCRM::read_ucd_variables_to_nonset_allocate(), moab::NCHelperHOMME::read_ucd_variables_to_nonset_allocate(), moab::NCHelperMPAS::read_ucd_variables_to_nonset_allocate(), moab::ParCommGraph::receive_mesh(), moab::ParCommGraph::receive_tag_values(), moab::HiReconstruction::reconstruct3D_curve_geom(), moab::HiReconstruction::reconstruct3D_surf_geom(), moab::MeshSetSequence::recursive_get_sets(), moab::ParallelComm::reduce_tags(), remove_entities_from_sets(), remove_from_vector(), moab::IntxUtils::remove_padded_vertices(), moab::WriteHDF5Parallel::remove_remote_entities(), moab::NestedRefine::reorder_indices(), moab::ReorderTool::reorder_tag_data(), moab::HalfFacetRep::resize_hf_maps(), moab::ParallelComm::resolve_shared_ents(), moab::ParallelComm::resolve_shared_sets(), moab::GeomTopoTool::restore_topology_from_adjacency(), moab::GeomTopoTool::restore_topology_from_geometric_inclusion(), moab::DualTool::rev_atomic_pillow(), moab::DualTool::rev_face_shrink(), moab::IntxUtils::ScaleToRadius(), moab::ScdBox::ScdBox(), moab::BitPage::search(), moab::ParCommGraph::send_tag_values(), moab::FBEngine::separate(), moab::GeomTopoTool::separate_by_dimension(), moab::Core::set_coords(), moab::ReadHDF5Dataset::set_file_ids(), moab::HalfFacetRep::set_sibling_map(), moab::ParallelComm::settle_intersection_points(), moab::Core::side_element(), moab::FBEngine::smooth_new_intx_points(), moab::OrientedBoxTreeTool::sphere_intersect_triangles(), moab::AdaptiveKDTree::sphere_intersect_triangles(), moab::FBEngine::split_bedge_at_new_mesh_node(), moab::FBEngine::split_boundary(), moab::FBEngine::split_edge_at_mesh_node(), moab::MeshTopoUtil::split_entities_manifold(), moab::MeshTopoUtil::split_entity_nonmanifold(), moab::FBEngine::split_internal_edge(), moab::ParCommGraph::split_owned_range(), moab::FBEngine::split_quads(), moab::FBEngine::split_surface_with_direction(), moab::MeshTopoUtil::star_entities_nonmanifold(), moab::MeshTopoUtil::star_next_entity(), moab::NestedRefine::subdivide_cells(), moab::NestedRefine::subdivide_tets(), subset_by_dimension(), moab::ParallelComm::tag_iface_entities(), moab::ScdInterface::tag_shared_vertices(), moab::ParallelMergeMesh::TagSharedElements(), test_spectral_hex(), test_spectral_quad(), TestErrorHandling_4(), moab::DualTool::traverse_hyperplane(), moab::unite(), moab::ParallelComm::unpack_sets(), moab::ReadNCDF::update(), moab::BoundBox::update(), moab::NestedRefine::update_global_ahf_1D(), moab::NestedRefine::update_global_ahf_1D_sub(), moab::NestedRefine::update_global_ahf_2D(), moab::NestedRefine::update_global_ahf_2D_sub(), moab::NestedRefine::update_global_ahf_3D(), moab::NestedRefine::update_local_ahf(), moab::ParallelComm::update_remote_data(), moab::ReorderTool::update_set_contents(), moab::NestedRefine::update_special_tags(), moab::Intx2MeshOnSphere::update_tracer_data(), moab::NestedRefine::update_tracking_verts(), upper_bound(), moab::RayIntersectSets::visit(), moab::TreeNodePrinter::visit(), moab::WriteHDF5::write_adjacencies(), MetisPartitioner::write_aggregationtag_partition(), moab::WriteVtk::write_bit_tag(), moab::WriteCCMIO::write_cells_and_faces(), moab::WriteNCDF::write_elementblocks(), moab::WriteHDF5::write_elems(), moab::WriteVtk::write_elems(), moab::WriteDamsel::write_entities(), moab::WriteCCMIO::write_external_faces(), moab::WriteNCDF::write_file(), moab::WriteAns::write_file(), moab::WriteGmsh::write_file(), moab::WriteSLAC::write_file(), moab::WriteSmf::write_file(), moab::WriteTemplate::write_file(), moab::WriteDamsel::write_file(), moab::Core::write_file(), moab::WriteNCDF::write_global_node_order_map(), moab::WriteSLAC::write_matsets(), moab::WriteTemplate::write_matsets(), moab::WriteHDF5::write_nodes(), moab::WriteCCMIO::write_nodes(), moab::WriteVtk::write_nodes(), moab::ScdNCWriteHelper::write_nonset_variables(), MetisPartitioner::write_partition(), ZoltanPartitioner::write_partition(), moab::WriteNCDF::write_poly_faces(), moab::WriteCCMIO::write_problem_description(), moab::WriteHDF5::write_set_data(), moab::WriteHDF5::write_sets(), moab::WriteHDF5::write_sparse_ids(), moab::WriteVtk::write_tag(), moab::WriteHDF5::write_tag_values(), moab::WriteHDF5::write_var_len_data(), moab::WriteHDF5::write_var_len_indices(), moab::WriteDamsel::write_vertices(), moab::IODebugTrack::~IODebugTrack(), and moab::Tqdcfr::~Tqdcfr().
void moab::Range::clear | ( | ) |
clears the contents of the list
Definition at line 197 of file Range.cpp.
References free_pair(), mHead, moab::Range::PairNode::mNext, and moab::Range::PairNode::mPrev.
Referenced by MetisPartitioner::assemble_graph(), ZoltanPartitioner::assemble_graph(), MetisPartitioner::assemble_taggedents_graph(), moab::ReadUtil::assign_ids(), moab::Skinner::classify_2d_boundary(), moab::TempestRemapper::clear(), moab::OrientedBoxTreeTool::closest_to_location(), moab::MeshTopoUtil::construct_aentities(), moab::DualTool::construct_hp_parent_child(), moab::GeomTopoTool::construct_vertex_ranges(), moab::TempestRemapper::convert_mesh_to_tempest_private(), moab::TempestRemapper::convert_overlap_mesh_sorted_by_source(), moab::WriteHDF5::count_set_size(), moab::ReadGmsh::create_geometric_topology(), moab::ParallelComm::create_iface_pc_links(), moab::ReadABAQUS::create_instance_of_part(), moab::AEntityFactory::create_vert_elem_adjacencies(), moab::Core::create_vertices(), moab::Core::delete_entities(), do_test_mode(), dot_get_sets(), moab::WriteHDF5Parallel::exchange_file_ids(), moab::ParallelComm::filter_pstatus(), moab::AdaptiveKDTree::find_close_triangle(), moab::MergeMesh::find_merged_to(), moab::ReadHDF5::find_sets_containing(), moab::DualTool::fs_get_quads(), moab::DualTool::fsr_get_fourth_quad(), moab::WriteHDF5Parallel::gather_interface_meshes(), moab::WriteVtk::gather_mesh(), moab::WriteHDF5::gather_mesh_info(), moab::ReadUtil::gather_related_ents(), get_adjacent_elems(), moab::MeshTopoUtil::get_bridge_adjacencies(), moab::GeomTopoTool::get_ct_children_by_dimension(), moab::DualTool::get_graphics_points(), moab::Coupler::get_matching_entities(), get_max_volume(), moab::DualTool::get_opposite_verts(), moab::ParallelComm::get_sent_ents(), moab::ReadABAQUS::get_set_elements(), moab::ReadABAQUS::get_set_nodes(), moab::WriteCCMIO::get_sets(), moab::ParallelComm::get_shared_entities(), moab::SharedSetData::get_shared_sets(), moab::WriteHDF5::get_sparse_tagged_entities(), moab::WriteHDF5::get_tag_data_length(), moab::WriteHDF5::get_write_entities(), moab::FBEngine::getAdjacentEntities(), moab::FBEngine::getEntities(), moab::TempestRemapper::GetOverlapAugmentedEntities(), iMOAB_UpdateMeshInfo(), ZoltanPartitioner::include_closure(), moab::SmoothFace::init_gradient(), moab::WriteHDF5::initialize_mesh(), moab::ReorderTool::int_order_from_sets_and_adj(), moab::AdaptiveKDTree::intersect_children_with_elems(), moab::Intx2Mesh::intersect_meshes(), moab::Intx2Mesh::intersect_meshes_kdtree(), moab::ReadCGNS::load_file(), moab::ReadGmsh::load_file(), moab::ReadHDF5::load_file_impl(), moab::ReadSms::load_file_impl(), moab::WriteGMV::local_write_mesh(), main(), make_tris_from_quads(), moab::GeomQueryTool::measure_area(), moab::GeomQueryTool::measure_volume(), moab::MergeMesh::merge_higher_dimensions(), moab::AdaptiveKDTree::merge_leaf(), moab::BSPTree::merge_leaf(), moab::DualTool::next_loop_vertex(), operator=(), moab::ParallelComm::pack_entities(), ZoltanPartitioner::partition_owned_cells(), moab::GeomQueryTool::point_in_volume_slow(), moab::ParallelMergeMesh::PopulateMySkinEnts(), moab::Core::print(), moab::AdaptiveKDTree::ray_intersect_triangles(), moab::OrientedBoxTreeTool::ray_intersect_triangles(), moab::ReadABAQUS::read_element_set(), DeformMeshRemap::read_file(), McnpData::read_mcnpfile(), moab::ReadABAQUS::read_node_set(), moab::ReadHDF5VarLen::read_offsets(), moab::ReadHDF5::read_sparse_tag_indices(), remove_entities_from_sets(), moab::ParallelComm::resolve_shared_sets(), moab::GeomTopoTool::restore_topology_from_adjacency(), moab::FBEngine::separate(), moab::ReadHDF5Dataset::set_all_file_ids(), moab::OrientedBoxTreeTool::sphere_intersect_triangles(), moab::split_box(), moab::MeshTopoUtil::split_entities_manifold(), moab::FBEngine::split_surface(), moab::MeshTopoUtil::star_next_entity(), tag_depth(), moab::ParallelComm::tag_iface_entities(), moab::DualTool::traverse_hyperplane(), moab::NestedRefine::update_global_ahf_1D_sub(), moab::NestedRefine::update_special_tags(), MetisPartitioner::write_aggregationtag_partition(), moab::WriteCCMIO::write_cells_and_faces(), moab::WriteNCDF::write_file(), moab::WriteSLAC::write_file(), moab::WriteTemplate::write_file(), moab::WriteHDF5::write_finished(), MetisPartitioner::write_partition(), ZoltanPartitioner::write_partition(), and ~Range().
|
inline |
Definition at line 952 of file Range.hpp.
References get_memory_use(), and size().
Referenced by moab::TempestRemapper::convert_mesh_to_tempest_private(), moab::ParallelComm::exchange_ghost_cells(), moab::ParallelComm::exchange_owned_mesh(), and moab::ParallelComm::send_entities().
|
inline |
Definition at line 748 of file Range.hpp.
References moab::Range::PairNode::mNext.
Referenced by moab::RangeSetIterator::build_pair_vec(), moab::BitTag::clear_data(), moab::convert_to_ranged_ids(), moab::ReadHDF5::delete_non_side_elements(), moab::BitTag::get_data(), index(), moab::MeshSet::insert_entity_ranges(), moab::ReadHDF5::insert_in_id_map(), moab::DebugOutput::list_ints_real(), moab::DebugOutput::list_range_real(), moab::merge_ranged_ids(), num_of_dimension(), num_of_type(), operator-=(), moab::operator==(), moab::PACK_RANGE(), moab::ParallelComm::pack_range_map(), psize(), moab::WriteHDF5::range_to_blocked_list(), moab::WriteHDF5::range_to_id_list(), moab::ReadHDF5VarLen::read_offsets(), moab::ReadHDF5::read_set_data(), moab::BitTag::remove_data(), moab::MeshSet::remove_entity_ranges(), moab::BitTag::set_data(), str_rep(), and moab::IODebugTrack::~IODebugTrack().
|
inline |
Definition at line 752 of file Range.hpp.
Referenced by moab::RangeSetIterator::build_pair_vec(), moab::BitTag::clear_data(), moab::convert_to_ranged_ids(), moab::VarLenDenseTag::find_entities_with_value(), moab::DenseTag::find_entities_with_value(), moab::BitTag::get_data(), index(), moab::MeshSet::insert_entity_ranges(), moab::ReadHDF5::insert_in_id_map(), moab::DebugOutput::list_ints_real(), moab::DebugOutput::list_range_real(), moab::merge_ranged_ids(), num_of_dimension(), num_of_type(), moab::operator==(), moab::PACK_RANGE(), moab::ParallelComm::pack_range_map(), psize(), moab::WriteHDF5::range_to_blocked_list(), moab::WriteHDF5::range_to_id_list(), moab::ReadHDF5VarLen::read_offsets(), moab::ReadHDF5::read_set_data(), moab::BitTag::remove_data(), moab::MeshSet::remove_entity_ranges(), moab::BitTag::set_data(), str_rep(), and moab::IODebugTrack::~IODebugTrack().
bool moab::Range::contains | ( | const Range & | other | ) | const |
Check if this range is a non-strict superset of some other range.
Definition at line 1005 of file Range.cpp.
References empty(), mHead, and moab::Range::PairNode::mNext.
Referenced by moab::WriteHDF5::check_dense_format_tag(), and moab::Skinner::find_skin_scd().
|
protected |
if dead_node is not mHead, remove it from the list and free it's memory.
Definition at line 471 of file Range.cpp.
References free_pair(), mHead, moab::Range::PairNode::mNext, and moab::Range::PairNode::mPrev.
Referenced by erase(), operator-=(), pop_back(), and pop_front().
|
inline |
return whether empty or not always use "if(!Ranges::empty())" instead of "if(Ranges::size())"
Definition at line 885 of file Range.hpp.
References mHead, and moab::Range::PairNode::mNext.
Referenced by moab::SpatialLocator::add_elems(), all_of_dimension(), all_of_type(), MetisPartitioner::assemble_graph(), ZoltanPartitioner::assemble_graph(), MetisPartitioner::assemble_taggedents_graph(), MetisPartitioner::assemble_taggedsets_graph(), moab::OrientedBoxTreeTool::build_tree(), moab::WriteHDF5::check_dense_format_tag(), moab::DualTool::check_dual_equiv_edges(), moab::NCHelperGCRM::check_existing_mesh(), moab::NCHelperHOMME::check_existing_mesh(), moab::NCHelperMPAS::check_existing_mesh(), moab::ParallelComm::check_global_ids(), moab::HalfFacetRep::check_mixed_entity_type(), moab::GeomTopoTool::check_model(), moab::ParallelComm::check_my_shared_handles(), moab::Skinner::classify_2d_boundary(), moab::MeshTag::clear_data(), moab::OrientedBoxTreeTool::closest_to_location(), moab::NCWriteGCRM::collect_mesh_info(), moab::ScdNCWriteHelper::collect_mesh_info(), moab::NCWriteHOMME::collect_mesh_info(), moab::NCWriteMPAS::collect_mesh_info(), moab::NCWriteGCRM::collect_variable_data(), moab::NCWriteMPAS::collect_variable_data(), moab::MeshTopoUtil::common_entity(), moab::ParCommGraph::compute_partition(), moab::TempestRemapper::ComputeOverlapMesh(), moab::DualTool::construct_dual_cells(), moab::DualTool::construct_dual_edges(), moab::DualTool::construct_dual_faces(), moab::DualTool::construct_dual_vertices(), moab::NestedRefine::construct_hm_2D(), moab::DualTool::construct_hp_parent_child(), moab::GeomTopoTool::construct_obb_tree(), contains(), moab::Tqdcfr::convert_nodesets_sidesets(), moab::SpectralMeshTool::convert_to_coarse(), moab::Coupler::Coupler(), moab::ReadIDEAS::create_elements(), moab::ReadGmsh::create_geometric_topology(), moab::ParallelComm::create_iface_pc_links(), create_mesh_no_holes(), moab::ReadParallel::create_partition_sets(), moab::ReadGmsh::create_sets(), moab::SpatialLocator::create_tree(), moab::DataCoupler::DataCoupler(), moab::Core::delete_entities(), moab::ReadHDF5::delete_non_side_elements(), moab::ReadParallel::delete_nonlocal_entities(), moab::ParallelComm::exchange_ghost_cells(), moab::ParallelComm::exchange_owned_mesh(), moab::ParallelComm::exchange_tags(), moab::DualTool::face_shrink(), moab::ParallelComm::filter_pstatus(), moab::Tree::find_all_trees(), moab::ScdInterface::find_boxes(), moab::AdaptiveKDTree::find_close_triangle(), moab::ParallelComm::find_existing_entity(), moab::Skinner::find_geometric_skin(), moab::MergeMesh::find_merged_to(), DeformMeshRemap::find_other_sets(), moab::ReadHDF5::find_sets_containing(), moab::Skinner::find_skin(), moab::Skinner::find_skin_noadj(), moab::DualTool::foc_delete_dual(), moab::DualTool::foc_get_ents(), moab::WriteCCMIO::gather_matset_info(), moab::WriteVtk::gather_mesh(), moab::WriteNCDF::gather_mesh_information(), moab::WriteUtil::gather_nodes_from_elements(), moab::ReadUtil::gather_related_ents(), moab::GeomTopoTool::geometrize_surface_set(), moab::get_adjacencies_intersection(), moab::MeshTopoUtil::get_average_position(), moab::MeshTopoUtil::get_bridge_adjacencies(), moab::MeshTag::get_data(), moab::WriteUtil::get_element_connect(), moab::Core::get_entities_by_type_and_tag(), moab::ReadUtil::get_gather_set(), moab::ParallelComm::get_ghosted_entities(), moab::ParallelData::get_interface_sets(), moab::Tqdcfr::get_mesh_entities(), moab::ParallelComm::get_remote_handles(), moab::ParallelComm::get_sent_ents(), moab::WriteHDF5::get_sparse_tagged_entities(), moab::ReadHDF5::get_subset_ids(), moab::ParallelComm::get_tag_send_list(), moab::Core::get_vertices(), moab::IntxUtils::global_gnomonic_projection(), iMOAB_CreateElements(), iMOAB_CreateVertices(), iMOAB_DeregisterApplication(), iMOAB_GetBlockInfo(), iMOAB_GetElementID(), iMOAB_GetElementOwnership(), iMOAB_UpdateMeshInfo(), ZoltanPartitioner::include_closure(), moab::NestedRefine::initialize(), moab::HiReconstruction::initialize(), moab::Coupler::initialize_spectral_elements(), moab::Coupler::initialize_tree(), moab::ReadHDF5::insert_in_id_map(), moab::Intx2Mesh::intersect_meshes(), moab::Intx2Mesh::intersect_meshes_kdtree(), moab::DualTool::is_blind(), moab::OrientedBoxTreeTool::join_trees(), moab::DebugOutput::list_ints_real(), moab::DebugOutput::list_range_real(), moab::ReadParallel::load_file(), moab::ReadHDF5::load_file_partial(), moab::SpatialLocator::locate_points(), main(), moab::WriteDamsel::map_sparse_tags(), moab::DualTool::next_loop_vertex(), moab::DualTool::order_chord(), moab::ParallelComm::pack_entities(), moab::ParallelComm::pack_sets(), ZoltanPartitioner::partition_mesh_and_geometry(), ZoltanPartitioner::partition_owned_cells(), perform_laplacian_smoothing(), perform_lloyd_relaxation(), moab::LloydSmoother::perform_smooth(), moab::ParallelMergeMesh::PopulateMySkinEnts(), moab::Core::print(), moab::TreeNodePrinter::print_contents(), process_partition_file(), moab::ReadHDF5::read_elems(), DeformMeshRemap::read_file(), moab::ReadHDF5::read_nodes(), moab::Tqdcfr::read_nodes(), moab::ReadHDF5VarLen::read_offsets(), moab::ReadHDF5::read_set_ids_recursive(), moab::ReadHDF5::read_sets(), moab::ReadHDF5::read_sparse_tag_indices(), moab::ParCommGraph::receive_mesh(), moab::FBEngine::redistribute_boundary_edges_to_faces(), moab::ParallelComm::reduce_tags(), moab::MeshTag::remove_data(), remove_entities_from_sets(), moab::ReorderTool::reorder_tag_data(), moab::ParallelComm::resolve_shared_ents(), moab::ParallelComm::resolve_shared_sets(), moab::GeomTopoTool::restore_topology_from_adjacency(), moab::DualTool::rev_face_shrink(), sanity_check(), moab::ScdBox::ScdBox(), moab::ParallelComm::send_entities(), moab::WriteHDF5::serial_create_file(), moab::MeshTag::set_data(), moab::Core::side_element(), moab::SpatialLocator::SpatialLocator(), moab::OrientedBoxTreeTool::sphere_intersect_triangles(), moab::MeshTopoUtil::star_entities_nonmanifold(), moab::MeshTopoUtil::star_next_entity(), str_rep(), moab::NestedRefine::subdivide_cells(), moab::NestedRefine::subdivide_tets(), tag_depth(), test_spectral_hex(), moab::ParallelComm::unpack_sets(), moab::ReorderTool::update_set_contents(), moab::RayIntersectSets::visit(), MetisPartitioner::write_aggregationtag_partition(), moab::WriteCCMIO::write_cells_and_faces(), moab::WriteNCDF::write_file(), moab::WriteGmsh::write_file(), moab::WriteSTL::write_file(), moab::WriteDamsel::write_file(), moab::WriteHDF5::write_file_impl(), moab::NCWriteGCRM::write_nonset_variables(), moab::NCWriteMPAS::write_nonset_variables(), MetisPartitioner::write_partition(), ZoltanPartitioner::write_partition(), moab::WriteNCDF::write_poly_faces(), moab::WriteCCMIO::write_problem_description(), moab::WriteVtk::write_tags(), and moab::IODebugTrack::~IODebugTrack().
|
inline |
return the ending const iterator for this range
Definition at line 872 of file Range.hpp.
References mHead.
Referenced by moab::SmoothFace::area(), moab::IntxAreaUtils::area_on_sphere(), moab::WriteSTL::ascii_write_triangles(), MetisPartitioner::assemble_graph(), ZoltanPartitioner::assemble_graph(), MetisPartitioner::assemble_taggedents_graph(), MetisPartitioner::assemble_taggedsets_graph(), moab::ScdInterface::assign_global_ids(), moab::WriteUtil::assign_ids(), moab::WriteSTL::binary_write_triangles(), moab::FBEngine::boundary_mesh_edges_on_face(), moab::FBEngine::boundary_nodes_on_face(), moab::box_from_axes(), moab::MeshGeneration::BrickInstance(), SphereDecomp::build_hexes(), moab::FBEngine::chain_able_edge(), moab::FBEngine::chain_two_edges(), moab::Core::check_adjacencies(), moab::DualTool::check_dual_adjs(), moab::DualTool::check_dual_equiv_edges(), moab::AEntityFactory::check_equiv_entities(), moab::NCHelperMPAS::check_existing_mesh(), moab::GeomTopoTool::check_model(), moab::ParallelComm::check_my_shared_handles(), moab::ParallelComm::check_sent_ents(), moab::Skinner::classify_2d_boundary(), moab::Intx2Mesh::clean(), moab::Core::clear_meshset(), moab::OrientedBoxTreeTool::closest_to_location(), moab::closest_to_triangles(), moab::WriteHDF5Parallel::communicate_shared_set_data(), moab::WriteHDF5Parallel::communicate_shared_set_ids(), moab::HiReconstruction::compute_average_vertex_normals_surf(), moab::HiReconstruction::compute_average_vertex_tangents_curve(), moab::SmoothFace::compute_control_points_on_edges(), compute_dual_mesh(), moab::SmoothFace::compute_internal_control_points_on_facets(), moab::FBEngine::compute_intersection_points(), SphereDecomp::compute_nodes(), moab::ParCommGraph::compute_partition(), moab::SmoothFace::compute_tangents_for_each_edge(), moab::TempestOnlineMap::ComputeAdjacencyRelations(), moab::TempestRemapper::ComputeOverlapMesh(), moab::DualTool::construct_dual(), moab::DualTool::construct_dual_cells(), moab::DualTool::construct_dual_edges(), moab::DualTool::construct_dual_faces(), moab::DualTool::construct_dual_vertices(), moab::BVHTree::construct_element_vec(), moab::DualTool::construct_hp_parent_child(), moab::GeomTopoTool::construct_obb_tree(), moab::GeomTopoTool::construct_obb_trees(), moab::GeomTopoTool::construct_vertex_ranges(), moab::TempestRemapper::convert_mesh_to_tempest_private(), moab::Tqdcfr::convert_nodesets_sidesets(), moab::TempestRemapper::convert_overlap_mesh_sorted_by_source(), moab::HigherOrderFactory::convert_sequence(), moab::SpectralMeshTool::convert_to_coarse(), moab::WriteHDF5::count_adjacencies(), moab::WriteHDF5::count_set_size(), moab::AEntityFactory::create_explicit_adjs(), moab::NCHelperGCRM::create_gather_set_vertices(), moab::NCHelperMPAS::create_gather_set_vertices(), moab::ParallelComm::create_iface_pc_links(), moab::ReadABAQUS::create_instance_of_part(), moab::NCHelperESMF::create_local_cells(), moab::NCHelperMPAS::create_local_cells(), moab::NCHelperGCRM::create_local_edges(), moab::NCHelperMPAS::create_local_edges(), moab::NCHelperESMF::create_local_vertices(), moab::NCHelperGCRM::create_local_vertices(), moab::NCHelperMPAS::create_local_vertices(), moab::ScdNCHelper::create_mesh(), moab::NCHelperDomain::create_mesh(), moab::NCHelperHOMME::create_mesh(), moab::NCHelperScrip::create_mesh(), moab::NCHelperESMF::create_padded_local_cells(), moab::NCHelperGCRM::create_padded_local_cells(), moab::NCHelperMPAS::create_padded_local_cells(), moab::ScdNCHelper::create_quad_coordinate_tag(), moab::ReadGmsh::create_sets(), moab::SpectralMeshTool::create_spectral_elems(), moab::AEntityFactory::create_vert_elem_adjacencies(), moab::Intx2Mesh::createTags(), moab::IntxUtils::deep_copy_set_with_quads(), moab::GeomTopoTool::delete_all_obb_trees(), moab::DualTool::delete_dual_entities(), moab::ParallelComm::delete_entities(), moab::ReadHDF5::delete_non_side_elements(), moab::ReadParallel::delete_nonlocal_entities(), moab::GeomTopoTool::delete_obb_tree(), moab::HalfFacetRep::determine_border_vertices(), moab::HalfFacetRep::determine_incident_halfedges(), moab::HalfFacetRep::determine_incident_halffaces(), moab::HalfFacetRep::determine_incident_halfverts(), moab::HalfFacetRep::determine_sibling_halfedges(), moab::HalfFacetRep::determine_sibling_halffaces(), moab::HalfFacetRep::determine_sibling_halfverts(), moab::Intx2Mesh::DetermineOrderedNeighbors(), do_test_mode(), dot_children(), dot_contained(), dot_nodes(), dot_write_id_nodes(), moab::SmoothFace::DumpModelControlPoints(), moab::GeomTopoTool::duplicate_model(), moab::IntxUtils::enforce_convexity(), equal_range(), erase(), moab::WriteHDF5Parallel::exchange_file_ids(), moab::DualTool::face_shrink(), moab::ParallelComm::filter_pstatus(), moab::Intx2Mesh::filterByMask(), find(), moab::ScdInterface::find_boxes(), moab::SparseTag::find_entities_with_value(), moab::VarLenSparseTag::find_entities_with_value(), moab::Skinner::find_geometric_skin(), moab::Skinner::find_inferred_edges(), moab::Skinner::find_skin_noadj(), moab::Skinner::find_skin_vertices_1D(), moab::Skinner::find_skin_vertices_2D(), moab::Skinner::find_skin_vertices_3D(), moab::HalfFacetRep::find_total_edges_2d(), moab::HalfFacetRep::find_total_edges_faces_3d(), moab::FBEngine::find_vertex_set_for_node(), moab::GeomQueryTool::find_volume_slow(), moab::Intx2Mesh::FindMaxEdgesInSet(), moab::IntxUtils::fix_degenerate_quads(), moab::DualTool::foc_delete_dual(), moab::ParCommGraph::form_tuples_to_migrate_mesh(), moab::ParallelComm::gather_data(), moab::WriteSLAC::gather_interior_exterior(), moab::WriteVtk::gather_mesh(), moab::WriteHDF5::gather_mesh_info(), moab::WriteNCDF::gather_mesh_information(), moab::WriteSLAC::gather_mesh_information(), moab::WriteTemplate::gather_mesh_information(), moab::WriteUtil::gather_nodes_from_elements(), moab::ReadUtil::gather_related_ents(), gather_set_stats(), moab::GeomTopoTool::generate_implicit_complement(), moab::TempestRemapper::GenerateMeshMetadata(), moab::GeomTopoTool::geometrize_surface_set(), moab::Core::get_adjacencies(), moab::get_adjacencies_intersection(), moab::MeshTopoUtil::get_average_position(), moab::Core::get_connectivity(), moab::Core::get_connectivity_by_type(), moab::GeomTopoTool::get_ct_children_by_dimension(), moab::AEntityFactory::get_down_adjacency_elements_poly(), moab::DualTool::get_dual_entities(), moab::WriteUtil::get_element_connect(), moab::Core::get_entities_by_handle(), moab::ParallelComm::get_ghosted_entities(), moab::Coupler::get_gl_points_on_elements(), moab::DualTool::get_graphics_points(), moab::ParallelComm::get_iface_entities(), moab::ParallelComm::get_interface_procs(), moab::ParallelComm::get_interface_sets(), moab::ParallelData::get_interface_sets(), moab::ParallelComm::get_local_handles(), moab::MeshTopoUtil::get_manifold(), get_max_volume(), moab::AEntityFactory::get_memory_use(), moab::WriteCCMIO::get_neuset_elems(), moab::WriteSLAC::get_neuset_elems(), moab::WriteTemplate::get_neuset_elems(), moab::HiReconstruction::get_normals_surf(), moab::ParallelComm::get_part_entities(), moab::ReadHDF5::get_partition(), moab::ParallelComm::get_proc_nvecs(), moab::ParallelComm::get_pstatus_entities(), moab::ParallelComm::get_remote_handles(), moab::ReorderTool::get_reordered_handles(), moab::ParallelComm::get_sent_ents(), moab::ReadABAQUS::get_set_by_name(), moab::WriteCCMIO::get_sets(), moab::WriteNCDF::get_sideset_elems(), get_signed_volume(), moab::WriteHDF5::get_tag_data_length(), moab::HiReconstruction::get_tangents_curve(), moab::MeshSetSequence::get_type(), moab::AEntityFactory::get_up_adjacency_elements(), moab::WriteSLAC::get_valid_sides(), moab::WriteTemplate::get_valid_sides(), moab::WriteNCDF::get_valid_sides(), moab::Core::get_vertex_coordinates(), moab::Core::get_vertices(), moab::FBEngine::getAdjacentEntities(), moab::FBEngine::getEntType(), moab::IntxUtils::global_gnomonic_projection(), hcFilter(), iMOAB_GetElementOwnership(), iMOAB_GetMeshInfo(), iMOAB_GetPointerToSurfaceBC(), iMOAB_GetPointerToVertexBC(), iMOAB_GetVertexOwnership(), iMOAB_GetVisibleElementsInfo(), iMOAB_SetDoubleTagStorageWithGid(), ZoltanPartitioner::include_closure(), moab::ReadHDF5Dataset::init(), moab::HalfFacetRep::init_curve(), moab::SmoothFace::init_gradient(), moab::HalfFacetRep::init_surface(), moab::HalfFacetRep::init_volume(), moab::WriteNCDF::initialize_exodus_file(), moab::FBEngine::initializeSmoothing(), insert(), moab::ReorderTool::int_order_from_sets_and_adj(), moab::intersect(), moab::AdaptiveKDTree::intersect_children_with_elems(), moab::Intx2Mesh::intersect_meshes(), moab::Intx2Mesh::intersect_meshes_kdtree(), moab::DualTool::is_blind(), moab::GeomTopoTool::is_owned_set(), moab::FBEngine::isEntAdj(), moab::OrientedBoxTreeTool::join_trees(), laplacianFilter(), moab::RayIntersectSets::leaf(), moab::ParallelComm::list_entities(), moab::Core::list_entities(), moab::ReadABAQUS::load_file(), moab::WriteGMV::local_write_mesh(), moab::Coupler::locate_points(), lower_bound(), main(), make_tris_from_quads(), moab::WriteDamsel::map_sparse_tags(), moab::GeomQueryTool::measure_area(), moab::GeomQueryTool::measure_volume(), merge(), moab::AEntityFactory::merge_adjust_adjacencies(), merge_duplicate_vertices(), moab::MergeMesh::merge_higher_dimensions(), merge_input_surfs(), moab::MergeMesh::merge_using_integer_tag(), moab::Coupler::nat_param(), obbstat_write(), obbvis_create(), moab::HiReconstruction::obtain_nring_ngbvs(), operator-=(), moab::DualTool::order_chord(), moab::ParallelComm::pack_entity_seq(), moab::ParallelComm::pack_sets(), ZoltanPartitioner::partition_owned_cells(), perform_laplacian_smoothing(), perform_lloyd_relaxation(), moab::LloydSmoother::perform_smooth(), moab::GeomQueryTool::point_in_volume_slow(), moab::IntxAreaUtils::positive_orientation(), moab::AdaptiveKDTree::print(), moab::Core::print(), moab::TreeNodePrinter::print_contents(), moab::WriteHDF5Parallel::print_set_sharing_data(), moab::HalfFacetRep::print_tags(), moab::print_type_sets(), moab::ReadDamsel::process_entity_tags(), quads_to_tris(), moab::AdaptiveKDTree::ray_intersect_triangles(), moab::OrientedBoxTreeTool::ray_intersect_triangles(), moab::ReadCCMIO::read_cells(), moab::ReadHDF5VarLen::read_data(), DeformMeshRemap::read_file(), moab::Tqdcfr::read_nodes(), moab::ReadNCDF::read_nodesets(), moab::ReadHDF5VarLen::read_offsets(), moab::ScdNCHelper::read_scd_variables_to_nonset_allocate(), moab::ReadHDF5::read_set_data(), moab::ReadNCDF::read_sidesets(), moab::ReadHDF5::read_sparse_tag_indices(), moab::ReadHDF5::read_tag_values_partial(), moab::NCHelperMPAS::read_ucd_variables_to_nonset(), moab::NCHelperGCRM::read_ucd_variables_to_nonset_allocate(), moab::NCHelperHOMME::read_ucd_variables_to_nonset_allocate(), moab::NCHelperMPAS::read_ucd_variables_to_nonset_allocate(), moab::ParCommGraph::receive_mesh(), moab::ParCommGraph::receive_tag_values(), moab::HiReconstruction::reconstruct3D_curve_geom(), moab::HiReconstruction::reconstruct3D_surf_geom(), moab::MeshSetSequence::recursive_get_sets(), remove_entities_from_sets(), remove_from_vector(), moab::IntxUtils::remove_padded_vertices(), moab::WriteHDF5Parallel::remove_remote_entities(), moab::WriteHDF5Parallel::remove_remote_sets(), moab::ReorderTool::reorder_tag_data(), moab::HalfFacetRep::resize_hf_maps(), moab::ParallelComm::resolve_shared_ents(), moab::ParallelComm::resolve_shared_sets(), moab::GeomTopoTool::restore_obb_index(), moab::GeomTopoTool::restore_topology_from_adjacency(), moab::GeomTopoTool::restore_topology_from_geometric_inclusion(), moab::DualTool::rev_atomic_pillow(), moab::DualTool::rev_face_shrink(), moab::IntxUtils::ScaleToRadius(), moab::ParCommGraph::send_tag_values(), moab::FBEngine::separate(), moab::GeomTopoTool::separate_by_dimension(), moab::Core::set_coords(), moab::ReadHDF5Dataset::set_file_ids(), moab::ParallelComm::settle_intersection_points(), moab::FBEngine::smooth_new_intx_points(), moab::OrientedBoxTreeTool::sphere_intersect_triangles(), moab::AdaptiveKDTree::sphere_intersect_triangles(), moab::FBEngine::split_bedge_at_new_mesh_node(), moab::FBEngine::split_boundary(), moab::FBEngine::split_edge_at_mesh_node(), moab::MeshTopoUtil::split_entities_manifold(), moab::MeshTopoUtil::split_entity_nonmanifold(), moab::FBEngine::split_internal_edge(), moab::DualTool::split_pair_nonmanifold(), moab::FBEngine::split_quads(), moab::FBEngine::split_surface_with_direction(), moab::MeshTopoUtil::star_entities_nonmanifold(), moab::MeshTopoUtil::star_next_entity(), subset_by_dimension(), moab::ParallelComm::tag_iface_entities(), moab::ScdInterface::tag_shared_vertices(), moab::ParallelComm::tag_shared_verts(), moab::ParallelMergeMesh::TagSharedElements(), test_spectral_hex(), test_spectral_quad(), TestErrorHandling_4(), moab::DualTool::traverse_hyperplane(), moab::unite(), moab::ParallelComm::unpack_sets(), moab::ReadNCDF::update(), moab::BoundBox::update(), moab::ParallelComm::update_remote_data(), moab::ReorderTool::update_set_contents(), moab::NestedRefine::update_special_tags(), moab::Intx2MeshOnSphere::update_tracer_data(), upper_bound(), moab::vector_remove_range(), moab::FBEngine::weave_lateral_face_from_edges(), moab::WriteHDF5::write_adjacencies(), MetisPartitioner::write_aggregationtag_partition(), moab::WriteVtk::write_bit_tag(), moab::WriteNCDF::write_elementblocks(), moab::WriteVtk::write_elems(), moab::WriteDamsel::write_entities(), moab::WriteCCMIO::write_external_faces(), moab::WriteNCDF::write_file(), moab::WriteAns::write_file(), moab::WriteGmsh::write_file(), moab::WriteSLAC::write_file(), moab::WriteSmf::write_file(), moab::WriteTemplate::write_file(), moab::WriteDamsel::write_file(), moab::Core::write_file(), moab::WriteNCDF::write_global_node_order_map(), moab::WriteSLAC::write_matsets(), moab::WriteTemplate::write_matsets(), moab::WriteCCMIO::write_nodes(), moab::WriteVtk::write_nodes(), moab::ScdNCWriteHelper::write_nonset_variables(), MetisPartitioner::write_partition(), ZoltanPartitioner::write_partition(), moab::WriteNCDF::write_poly_faces(), moab::WriteHDF5::write_set_data(), moab::WriteHDF5::write_sets(), moab::WriteVtk::write_tag(), moab::WriteHDF5::write_var_len_data(), moab::WriteDamsel::write_vertices(), and moab::Tqdcfr::~Tqdcfr().
std::pair< Range::const_iterator, Range::const_iterator > moab::Range::equal_range | ( | EntityType | type | ) | const |
Definition at line 867 of file Range.cpp.
References begin(), moab::CREATE_HANDLE(), end(), and lower_bound().
Referenced by moab::ParallelComm::add_verts(), moab::SparseTag::find_entities_with_value(), moab::VarLenSparseTag::find_entities_with_value(), moab::ReadUtil::gather_related_ents(), moab::BitTag::get_entities_with_bits(), moab::WriteHDF5::initialize_mesh(), and subset_by_type().
|
inline |
Range::iterator moab::Range::erase | ( | iterator | iter | ) |
remove an item from this list and return an iterator to the next item
erases an item from this list and returns an iterator to the next item
Definition at line 364 of file Range.cpp.
References alloc_pair(), end(), free_pair(), moab::Range::PairNode::mNext, moab::Range::const_iterator::mNode, moab::Range::PairNode::mPrev, and moab::Range::const_iterator::mValue.
Referenced by MetisPartitioner::assemble_taggedsets_graph(), moab::ParallelComm::check_clean_iface(), moab::DualTool::check_dual_equiv_edges(), moab::ParallelComm::check_my_shared_handles(), moab::Skinner::classify_2d_boundary(), moab::TempestRemapper::ComputeOverlapMesh(), moab::ReadCCMIO::create_cell_from_faces(), moab::ReadHDF5::delete_non_side_elements(), moab::GeomTopoTool::delete_obb_tree(), moab::ParallelComm::destroy_part(), erase(), moab::Skinner::find_skin_noadj(), moab::Skinner::find_skin_vertices_3D(), moab::DualTool::foc_delete_dual(), moab::DualTool::foc_get_addl_ents(), moab::DualTool::foc_get_ents(), moab::DualTool::fsr_get_fourth_quad(), moab::WriteVtk::gather_mesh(), moab::GeomTopoTool::geometrize_surface_set(), moab::get_adjacencies_intersection(), moab::MeshTopoUtil::get_bridge_adjacencies(), moab::ParallelComm::get_interface_sets(), moab::WriteCCMIO::get_neuset_elems(), moab::WriteSLAC::get_neuset_elems(), moab::WriteTemplate::get_neuset_elems(), moab::ReorderTool::get_new_handles(), moab::ReadHDF5::get_partition(), moab::WriteNCDF::get_sideset_elems(), moab::Core::get_vertices(), hcFilter(), moab::Intx2Mesh::intersect_meshes(), laplacianFilter(), merge_duplicate_vertices(), moab::MergeMesh::merge_higher_dimensions(), moab::DualTool::next_loop_vertex(), operator-=(), moab::DualTool::order_chord(), moab::ReadHDF5::read_sparse_tag(), moab::ReorderTool::reorder_tag_data(), moab::ParallelComm::resolve_shared_ents(), moab::ParallelComm::resolve_shared_sets(), moab::DualTool::rev_face_shrink(), moab::Core::serial_load_file(), moab::MeshTopoUtil::star_entities(), moab::MeshTopoUtil::star_entities_nonmanifold(), moab::MeshTopoUtil::star_next_entity(), moab::ParallelMergeMesh::TagSharedElements(), and moab::ReadNCDF::update().
Range::iterator moab::Range::erase | ( | iterator | iter1, |
iterator | iter2 | ||
) |
remove a range of items from the list
Definition at line 411 of file Range.cpp.
References alloc_pair(), delete_pair_node(), mHead, moab::Range::PairNode::mNext, moab::Range::const_iterator::mNode, moab::Range::PairNode::mPrev, and moab::Range::const_iterator::mValue.
Range::const_iterator moab::Range::find | ( | EntityHandle | val | ) | const |
find an item int the list and return an iterator at that value
finds a value in the list. this method is preferred over other algorithms because it can be found faster this way.
Definition at line 510 of file Range.cpp.
References end(), mHead, and moab::Range::PairNode::mNext.
Referenced by moab::Core::check_adjacencies(), moab::DualTool::check_dual_adjs(), moab::AEntityFactory::check_equiv_entities(), moab::Skinner::classify_2d_boundary(), moab::FBEngine::compute_intersection_points(), moab::GeomTopoTool::duplicate_model(), erase(), moab::Skinner::find_skin_noadj(), moab::ReadUtil::gather_related_ents(), moab::GeomTopoTool::geometrize_surface_set(), moab::FBEngine::getAdjacentEntities(), moab::FBEngine::getEntType(), ZoltanPartitioner::include_closure(), moab::Intx2Mesh::intersect_meshes(), moab::GeomTopoTool::is_owned_set(), moab::FBEngine::isEntAdj(), merge_duplicate_vertices(), moab::HiReconstruction::obtain_nring_ngbvs(), moab::WriteHDF5Parallel::remove_remote_entities(), moab::WriteHDF5Parallel::remove_remote_sets(), moab::FBEngine::separate(), moab::ParallelComm::settle_intersection_points(), moab::FBEngine::split_boundary(), moab::DualTool::split_pair_nonmanifold(), moab::DualTool::traverse_hyperplane(), moab::vector_remove_range(), moab::FBEngine::weave_lateral_face_from_edges(), and moab::WriteSLAC::write_matsets().
|
inline |
get first entity in range
Definition at line 907 of file Range.hpp.
References mHead, and moab::Range::PairNode::mNext.
Referenced by all_of_dimension(), all_of_type(), moab::ReadMCNP5::average_with_existing_tally(), moab::MeshTag::clear_data(), moab::OrientedBoxTreeTool::closest_to_location(), moab::ReadIDEAS::create_elements(), moab::ReadHDF5::delete_non_side_elements(), moab::GeomTopoTool::entity_by_id(), moab::WriteNCDF::gather_mesh_information(), moab::WriteUtil::gather_nodes_from_elements(), moab::MeshTag::get_data(), moab::ParallelComm::get_part_handle(), moab::ReadHDF5::get_partition(), moab::ReadHDF5::insert_in_id_map(), moab::Intx2Mesh::intersect_meshes(), pop_front(), moab::ReadHDF5::read_elems(), moab::ReadHDF5VarLen::read_offsets(), moab::MeshTag::remove_data(), moab::GeomTopoTool::resize_rootSets(), moab::GeomTopoTool::restore_topology_from_adjacency(), moab::MeshTag::set_data(), and moab::OrientedBoxTreeTool::sphere_intersect_triangles().
unsigned long moab::Range::get_memory_use | ( | ) | const |
Definition at line 997 of file Range.cpp.
References mHead, and moab::Range::PairNode::mNext.
Referenced by compactness().
|
inline |
Definition at line 936 of file Range.hpp.
References begin(), const_pair_begin(), const_pair_end(), and rbegin().
Referenced by MetisPartitioner::assemble_graph(), moab::NestedRefine::construct_hm_1D(), moab::NestedRefine::construct_hm_2D(), moab::TempestRemapper::convert_mesh_to_tempest_private(), moab::NCHelperESMF::create_local_cells(), moab::NCHelperMPAS::create_local_cells(), moab::NCHelperGCRM::create_local_edges(), moab::NCHelperMPAS::create_local_edges(), moab::NCHelperESMF::create_padded_local_cells(), moab::NCHelperGCRM::create_padded_local_cells(), moab::NCHelperMPAS::create_padded_local_cells(), moab::ParallelComm::create_part(), moab::ParallelComm::delete_entities(), moab::ParallelComm::destroy_part(), moab::HalfFacetRep::determine_incident_halfedges(), moab::HalfFacetRep::determine_sibling_halfedges(), moab::HalfFacetRep::determine_sibling_halffaces(), moab::HalfFacetRep::determine_sibling_halfverts(), moab::HalfFacetRep::find_total_edges_2d(), moab::HalfFacetRep::find_total_edges_faces_3d(), moab::Intx2MeshInPlane::findNodes(), moab::Intx2MeshOnSphere::findNodes(), moab::IntxRllCssphere::findNodes(), moab::HiReconstruction::get_fittings_data(), moab::HiReconstruction::get_normals_surf(), moab::HiReconstruction::get_tangents_curve(), hcFilter(), moab::HiReconstruction::hiproj_walf_around_vertex(), moab::HiReconstruction::hiproj_walf_in_element(), iMOAB_GetBlockElementConnectivities(), iMOAB_GetElementConnectivity(), iMOAB_GetElementID(), iMOAB_GetNeighborElements(), iMOAB_GetPointerToSurfaceBC(), iMOAB_GetPointerToVertexBC(), iMOAB_GetVisibleElementsInfo(), laplacianFilter(), main(), moab::HalfFacetRep::mark_halfedges(), moab::NestedRefine::parent_to_child(), moab::HiReconstruction::polyfit3d_walf_curve_vertex(), moab::HiReconstruction::polyfit3d_walf_surf_vertex(), moab::NCHelperMPAS::read_ucd_variables_to_nonset(), moab::HiReconstruction::reconstruct3D_curve_geom(), moab::HiReconstruction::reconstruct3D_surf_geom(), moab::NestedRefine::update_global_ahf_3D(), moab::Intx2MeshOnSphere::update_tracer_data(), moab::WriteNCDF::write_elementblocks(), and moab::WriteVtk::write_elems().
|
inline |
|
inline |
|
inline |
Range::iterator moab::Range::insert | ( | Range::iterator | prev, |
EntityHandle | first, | ||
EntityHandle | last | ||
) |
insert a range of items into this list and return the iterator for the first inserted item
Definition at line 295 of file Range.cpp.
References alloc_pair(), end(), moab::GeomUtil::first(), free_pair(), mHead, moab::Range::PairNode::mNext, moab::Range::const_iterator::mNode, and moab::Range::PairNode::mPrev.
Range::iterator moab::Range::insert | ( | Range::iterator | hint, |
EntityHandle | val | ||
) |
inserts a single value into this range
Definition at line 229 of file Range.cpp.
References alloc_pair(), end(), free_pair(), mHead, moab::Range::PairNode::mNext, moab::Range::const_iterator::mNode, and moab::Range::PairNode::mPrev.
Referenced by moab::ReadSms::add_entities(), moab::GeomTopoTool::add_geo_set(), moab::SharedSetData::append_local_handles(), moab::ParallelComm::assign_global_ids(), moab::ReadUtil::assign_ids(), moab::DualTool::atomic_pillow(), moab::FBEngine::BreakTriangle2(), moab::MeshGeneration::BrickInstance(), moab::ReadRTT::build_moab(), moab::ParallelComm::check_clean_iface(), moab::DualTool::check_dual_equiv_edges(), moab::ParallelComm::check_my_shared_handles(), moab::ParallelComm::check_sent_ents(), moab::Skinner::classify_2d_boundary(), moab::MeshTopoUtil::common_entity(), moab::WriteHDF5Parallel::communicate_shared_set_ids(), moab::TempestOnlineMap::ComputeAdjacencyRelations(), moab::TempestRemapper::ComputeOverlapMesh(), moab::DualTool::construct_dual_cells(), moab::DualTool::construct_dual_edges(), moab::DualTool::construct_dual_faces(), moab::DualTool::construct_dual_vertices(), moab::DualTool::construct_hp_parent_child(), moab::GeomTopoTool::construct_obb_tree(), moab::GeomTopoTool::construct_obb_trees(), moab::TempestRemapper::ConstructCoveringSet(), moab::DamselUtil::container_to_range(), moab::Tqdcfr::convert_nodesets_sidesets(), moab::ReadHDF5::convert_range_to_handle(), moab::convert_to_ranged_ids(), moab::BVHTree::convert_tree(), moab::copy_set_contents(), moab::copy_sorted_file_ids(), moab::ReadGmsh::create_elements(), moab::NCHelperMPAS::create_gather_set_cells(), moab::ReadABAQUS::create_instance_of_part(), moab::ParallelComm::create_interface_sets(), moab::NCHelperESMF::create_local_cells(), moab::NCHelperMPAS::create_local_cells(), moab::ScdNCHelper::create_mesh(), moab::NCHelperDomain::create_mesh(), moab::NCHelperESMF::create_mesh(), moab::NCHelperGCRM::create_mesh(), moab::NCHelperHOMME::create_mesh(), moab::NCHelperMPAS::create_mesh(), moab::NCHelperScrip::create_mesh(), moab::FBEngine::create_new_gedge(), moab::NCHelperESMF::create_padded_local_cells(), moab::NCHelperGCRM::create_padded_local_cells(), moab::NCHelperMPAS::create_padded_local_cells(), moab::ParallelComm::create_part(), moab::ScdInterface::create_scd_sequence(), moab::ReadGmsh::create_sets(), moab::ReadTemplate::create_sets(), moab::SpectralMeshTool::create_spectral_elems(), moab::ReadOBJ::create_tri_faces(), moab::Core::create_vertices(), moab::ReadIDEAS::create_vertices(), moab::DualTool::delete_dual_entities(), moab::Core::delete_entities(), moab::ParallelComm::delete_entities(), moab::ReadHDF5::delete_non_side_elements(), moab::ReadParallel::delete_nonlocal_entities(), moab::GeomTopoTool::delete_obb_tree(), moab::FBEngine::divide_triangle(), moab::GeomTopoTool::duplicate_model(), moab::Core::estimated_memory_use(), moab::DualTool::face_shrink(), moab::ParallelComm::filter_pstatus(), moab::Intx2Mesh::filterByMask(), moab::ScdInterface::find_boxes(), moab::Skinner::find_geometric_skin(), moab::Skinner::find_inferred_edges(), moab::find_map_values(), moab::ReadHDF5::find_sets_containing(), moab::Skinner::find_skin_noadj(), moab::Skinner::find_skin_vertices_1D(), moab::Skinner::find_skin_vertices_2D(), moab::Skinner::find_skin_vertices_3D(), moab::find_tag_values(), moab::DualTool::foc_delete_dual(), moab::DualTool::foc_get_ents(), moab::ParCommGraph::form_mesh_from_tuples(), moab::DualTool::fs_get_quad_loops(), moab::DualTool::fsr_get_fourth_quad(), moab::WriteSLAC::gather_interior_exterior(), moab::WriteVtk::gather_mesh(), moab::WriteHDF5::gather_mesh_info(), moab::WriteUtil::gather_nodes_from_elements(), moab::ReadUtil::gather_related_ents(), moab::GeomTopoTool::geometrize_surface_set(), moab::get_adjacencies_union(), moab::MeshTopoUtil::get_bridge_adjacencies(), moab::DualTool::get_cell_points(), moab::Core::get_connectivity(), moab::ReadDamsel::get_contents(), moab::GeomTopoTool::get_ct_children_by_dimension(), moab::AEntityFactory::get_down_adjacency_elements_poly(), moab::DualTool::get_dual_entities(), moab::ReadABAQUS::get_elements_by_id(), moab::TypeSequenceManager::get_entities(), moab::MeshSet::get_entities_by_dimension(), moab::MeshSet::get_entities_by_type(), moab::DualTool::get_graphics_points(), moab::MeshTopoUtil::get_manifold(), moab::ReadABAQUS::get_nodes_by_id(), moab::MeshSet::get_non_set_entities(), moab::ParallelComm::get_pstatus_entities(), moab::ReadHDF5::get_tagged_entities(), moab::TempestRemapper::GetOverlapAugmentedEntities(), moab::IntxUtils::global_gnomonic_projection(), iMOAB_DeregisterApplication(), insert(), insert_list(), moab::AdaptiveKDTree::intersect_children_with_elems(), moab::Intx2Mesh::intersect_meshes(), moab::ReadParallel::load_file(), moab::ReadHDF5::load_file_impl(), main(), make_tris_from_quads(), moab::MergeMesh::merge_entities(), moab::MergeMesh::merge_higher_dimensions(), moab::merge_ranged_ids(), moab::DualTool::next_loop_vertex(), moab::HiReconstruction::obtain_nring_ngbvs(), moab::range_inserter::operator=(), ZoltanPartitioner::partition_owned_cells(), perform_laplacian_smoothing(), perform_lloyd_relaxation(), moab::LloydSmoother::perform_smooth(), moab::AdaptiveKDTree::print(), moab::print_type_sets(), moab::ReadDamsel::process_ent_info(), process_partition_file(), moab::ReadCCMIO::read_cells(), moab::ReadHDF5::read_dense_tag(), moab::ReadTetGen::read_elem_file(), moab::ReadNCDF::read_elements(), moab::ReadTemplate::read_elements(), moab::ReadHDF5::read_elems(), McnpData::read_mcnpfile(), moab::ReadHDF5::read_node_adj_elems(), moab::ReadTetGen::read_node_file(), moab::ReadNCDF::read_nodes(), moab::ReadHDF5VarLen::read_offsets(), moab::ReadHDF5::read_set_data(), moab::ReadHDF5::read_sparse_tag_indices(), moab::ReadHDF5::read_tag_values_partial(), moab::ReadTemplate::read_vertices(), moab::MeshSetSequence::recursive_get_sets(), remove_entities_from_sets(), moab::IntxUtils::remove_padded_vertices(), moab::WriteHDF5Parallel::remove_remote_entities(), moab::WriteHDF5Parallel::remove_remote_sets(), moab::ParallelComm::resolve_shared_sets(), moab::DualTool::rev_atomic_pillow(), moab::DualTool::rev_face_shrink(), moab::BitPage::search(), moab::GeomTopoTool::separate_by_dimension(), moab::ReadHDF5Dataset::set_all_file_ids(), moab::Core::side_element(), moab::Skinner::skin_box(), moab::split_box(), moab::MeshTopoUtil::split_entities_manifold(), moab::MeshTopoUtil::split_entity_nonmanifold(), moab::FBEngine::split_internal_edge(), moab::ParCommGraph::split_owned_range(), moab::DualTool::split_pair_nonmanifold(), moab::FBEngine::split_surface_with_direction(), moab::MeshTopoUtil::star_next_entity(), moab::ReadHDF5::store_file_ids(), moab::ReadHDF5::store_sets_file_ids(), subset_by_dimension(), subset_by_type(), tag_depth(), moab::ParallelComm::tag_iface_entities(), moab::ScdInterface::tag_shared_vertices(), moab::DualTool::traverse_hyperplane(), moab::unite(), moab::UNPACK_RANGE(), moab::ParallelComm::unpack_sets(), moab::ReadNCDF::update(), moab::ReorderTool::update_set_contents(), moab::ReadVtk::vtk_read_polydata(), moab::ReadVtk::vtk_read_rectilinear_grid(), moab::ReadVtk::vtk_read_structured_grid(), moab::ReadVtk::vtk_read_structured_points(), moab::ReadVtk::vtk_read_unstructured_grid(), MetisPartitioner::write_aggregationtag_partition(), moab::WriteCCMIO::write_cells_and_faces(), moab::WriteGmsh::write_file(), MetisPartitioner::write_partition(), ZoltanPartitioner::write_partition(), moab::WriteHDF5::write_set_data(), and moab::IODebugTrack::~IODebugTrack().
void moab::Range::insert | ( | Range::const_iterator | begini, |
Range::const_iterator | endi | ||
) |
merges another Range with this one
Definition at line 523 of file Range.cpp.
References insert(), moab::Range::PairNode::mNext, and moab::Range::const_iterator::mNode.
|
inline |
Range::iterator moab::Range::insert_list | ( | Iterator | begin_iter, |
Iterator | end_iter | ||
) |
iterator moab::Range::insert_list | ( | T | begin_iter, |
T | end_iter | ||
) |
Referenced by moab::GeomTopoTool::delete_obb_tree().
|
static |
return an iterator to the first value >= val
Definition at line 809 of file Range.cpp.
References moab::GeomUtil::first(), moab::Range::PairNode::mNext, and moab::Range::const_iterator::mNode.
Referenced by MetisPartitioner::assemble_taggedsets_graph(), moab::OrientedBox::compute_from_vertices(), moab::OrientedBox::covariance_data_from_tris(), moab::ReadHDF5::delete_non_side_elements(), equal_range(), moab::VarLenDenseTag::find_entities_with_value(), moab::DenseTag::find_entities_with_value(), moab::WriteVtk::gather_mesh(), moab::ReadHDF5::get_partition(), moab::ParallelComm::get_shared_entities(), moab::intersect(), moab::AdaptiveKDTree::intersect_children_with_elems(), moab::ReadHDF5::load_file_partial(), moab::WriteGMV::local_write_mesh(), lower_bound(), moab::ReadHDF5::read_tag_values_partial(), moab::WriteHDF5Parallel::remove_remote_entities(), moab::ParallelComm::resolve_shared_ents(), moab::ParallelComm::resolve_shared_sets(), moab::OrientedBoxTreeTool::sphere_intersect_triangles(), subset_by_dimension(), moab::ParallelMergeMesh::TagSharedElements(), moab::BoundBox::update(), and upper_bound().
|
inline |
Range::const_iterator moab::Range::lower_bound | ( | EntityType | type | ) | const |
Definition at line 839 of file Range.cpp.
References begin(), moab::CREATE_HANDLE(), end(), and lower_bound().
Range::const_iterator moab::Range::lower_bound | ( | EntityType | type, |
const_iterator | first | ||
) | const |
Definition at line 845 of file Range.cpp.
References moab::CREATE_HANDLE(), end(), moab::GeomUtil::first(), and lower_bound().
|
inline |
merges this Range with another range
Definition at line 346 of file Range.hpp.
References begin(), and end().
Referenced by moab::DualTool::atomic_pillow(), moab::MeshGeneration::BrickInstance(), moab::AdaptiveKDTree::build_tree(), moab::AxisBox::calculate(), moab::ParallelComm::check_my_shared_handles(), moab::WriteHDF5Parallel::communicate_shared_set_data(), compute_dual_mesh(), moab::TempestRemapper::ComputeOverlapMesh(), moab::DualTool::construct_dual_faces(), moab::NCHelperHOMME::create_mesh(), DeformMeshRemap::deform_master(), dot_get_sets(), moab::WriteHDF5Parallel::exchange_file_ids(), moab::ParallelComm::exchange_owned_mesh(), DeformMeshRemap::execute(), moab::Skinner::find_skin(), moab::Skinner::find_skin_scd(), moab::Skinner::find_skin_vertices_3D(), moab::DualTool::foc_delete_dual(), moab::DualTool::foc_get_addl_ents(), moab::WriteHDF5Parallel::gather_interface_meshes(), moab::WriteVtk::gather_mesh(), moab::WriteHDF5::gather_mesh_info(), moab::ReadUtil::gather_related_ents(), moab::MeshTopoUtil::get_bridge_adjacencies(), moab::Core::get_connectivity(), moab::Core::get_entities_by_type_and_tag(), moab::ParallelComm::get_ghosted_entities(), moab::MeshTopoUtil::get_manifold(), moab::ParallelComm::get_part_entities(), moab::ParallelComm::get_sent_ents(), moab::ReadABAQUS::get_set_elements(), moab::ReadABAQUS::get_set_nodes(), moab::ParallelComm::get_shared_entities(), moab::WriteHDF5::get_sparse_tagged_entities(), moab::WriteHDF5::get_tag_data_length(), moab::WriteSTL::get_triangles(), moab::Core::get_vertices(), moab::WriteHDF5::get_write_entities(), moab::FBEngine::getEntities(), moab::TempestRemapper::GetOverlapAugmentedEntities(), iMOAB_CreateElements(), iMOAB_CreateVertices(), moab::ReorderTool::int_order_from_sets_and_adj(), moab::intersect(), moab::Intx2Mesh::intersect_meshes(), moab::Intx2Mesh::intersect_meshes_kdtree(), moab::ReadHDF5::load_file_partial(), main(), merge_duplicate_vertices(), obbvis_create(), moab::AdaptiveKDTree::print(), moab::ReadABAQUS::read_element_set(), DeformMeshRemap::read_file(), moab::ReadABAQUS::read_node_set(), moab::Tqdcfr::read_nodes(), moab::ReadHDF5::read_set_ids_recursive(), moab::ReadHDF5::read_sparse_tag(), moab::ReadHDF5::read_tag_values_partial(), moab::ReadCCMIO::read_vertices(), moab::ParCommGraph::receive_mesh(), moab::WriteHDF5Parallel::remove_remote_entities(), moab::GeomTopoTool::resize_rootSets(), moab::ParallelComm::resolve_shared_ents(), moab::ParallelComm::resolve_shared_sets(), moab::ParallelComm::send_entities(), moab::ParCommGraph::send_mesh_parts(), moab::Skinner::skin_box(), moab::MeshTopoUtil::split_entities_manifold(), moab::ParallelMergeMesh::TagSharedElements(), moab::WriteCCMIO::write_cells_and_faces(), moab::WriteVtk::write_elems(), moab::WriteGmsh::write_file(), moab::WriteHDF5::write_sparse_ids(), moab::WriteHDF5::write_tag_values(), moab::WriteVtk::write_tags(), and moab::WriteHDF5::write_var_len_indices().
|
inline |
unsigned moab::Range::num_of_dimension | ( | int | dim | ) | const |
Definition at line 911 of file Range.cpp.
References const_pair_begin(), const_pair_end(), moab::CREATE_HANDLE(), dim, moab::CN::Dimension(), moab::GeomUtil::first(), MB_END_ID, moab::TYPE_FROM_HANDLE(), and moab::CN::TypeDimensionMap.
Referenced by main().
unsigned moab::Range::num_of_type | ( | EntityType | type | ) | const |
Definition at line 890 of file Range.cpp.
References const_pair_begin(), const_pair_end(), moab::ID_FROM_HANDLE(), MB_END_ID, and moab::TYPE_FROM_HANDLE().
just like subtract, but as an operator
Definition at line 730 of file Range.cpp.
References alloc_pair(), begin(), const_pair_begin(), delete_pair_node(), end(), erase(), moab::Range::PairNode::mNext, moab::Range::PairNode::mPrev, moab::Range::pair_iterator::node(), and pair_begin().
operator=
Definition at line 210 of file Range.cpp.
References alloc_pair(), clear(), mHead, moab::Range::PairNode::mNext, and moab::Range::PairNode::mPrev.
|
inline |
|
inline |
Definition at line 739 of file Range.hpp.
References moab::Range::PairNode::mNext.
Referenced by moab::ReadUtil::assign_ids(), moab::NCHelperGCRM::create_local_edges(), moab::NCHelperMPAS::create_local_edges(), moab::NCHelperESMF::create_local_vertices(), moab::NCHelperGCRM::create_local_vertices(), moab::NCHelperMPAS::create_local_vertices(), moab::NCHelperESMF::create_mesh(), moab::NCHelperGCRM::create_mesh(), moab::NCHelperMPAS::create_mesh(), moab::NCHelperScrip::create_mesh(), moab::DataCoupler::DataCoupler(), moab::UcdNCWriteHelper::jik_to_kji_stride(), moab::UcdNCHelper::kji_to_jik_stride(), operator-=(), moab::NCHelperGCRM::read_ucd_variables_to_nonset(), moab::NCHelperHOMME::read_ucd_variables_to_nonset(), moab::NCHelperMPAS::read_ucd_variables_to_nonset(), moab::NCWriteGCRM::write_nonset_variables(), moab::NCWriteHOMME::write_nonset_variables(), and moab::NCWriteMPAS::write_nonset_variables().
|
inline |
Definition at line 756 of file Range.hpp.
References moab::Range::PairNode::mNext.
|
inline |
Definition at line 743 of file Range.hpp.
Referenced by moab::ReadUtil::assign_ids(), moab::NCHelperGCRM::create_local_edges(), moab::NCHelperMPAS::create_local_edges(), moab::NCHelperESMF::create_local_vertices(), moab::NCHelperGCRM::create_local_vertices(), moab::NCHelperMPAS::create_local_vertices(), moab::NCHelperESMF::create_mesh(), moab::NCHelperGCRM::create_mesh(), moab::NCHelperMPAS::create_mesh(), moab::NCHelperScrip::create_mesh(), moab::DataCoupler::DataCoupler(), moab::UcdNCWriteHelper::jik_to_kji_stride(), moab::UcdNCHelper::kji_to_jik_stride(), moab::NCHelperGCRM::read_ucd_variables_to_nonset(), moab::NCHelperHOMME::read_ucd_variables_to_nonset(), moab::NCHelperMPAS::read_ucd_variables_to_nonset(), moab::NCWriteGCRM::write_nonset_variables(), moab::NCWriteHOMME::write_nonset_variables(), and moab::NCWriteMPAS::write_nonset_variables().
|
inline |
EntityHandle moab::Range::pop_back | ( | ) |
remove last entity from range
Definition at line 494 of file Range.cpp.
References back(), delete_pair_node(), mHead, and moab::Range::PairNode::mPrev.
Referenced by MetisPartitioner::write_aggregationtag_partition(), MetisPartitioner::write_partition(), and ZoltanPartitioner::write_partition().
EntityHandle moab::Range::pop_front | ( | ) |
remove first entity from range
Definition at line 482 of file Range.cpp.
References delete_pair_node(), front(), mHead, and moab::Range::PairNode::mNext.
Referenced by moab::GeomTopoTool::insert_in_tree().
void moab::Range::print | ( | const char * | indent_prefix = NULL | ) | const |
Definition at line 618 of file Range.cpp.
Referenced by moab::MeshGeneration::BrickInstance(), moab::ParallelComm::list_entities(), main(), moab::operator<<(), and moab::ParallelComm::pack_tag().
void moab::Range::print | ( | std::ostream & | s, |
const char * | indent_prefix = NULL |
||
) | const |
|
inline |
return the number of range pairs in the list
Definition at line 979 of file Range.hpp.
References const_pair_begin(), and const_pair_end().
Referenced by moab::RangeSetIterator::build_pair_vec(), moab::SpectralMeshTool::convert_to_coarse(), moab::convert_to_ranged_ids(), moab::NCHelperGCRM::create_local_edges(), moab::NCHelperMPAS::create_local_edges(), moab::NCHelperESMF::create_local_vertices(), moab::NCHelperGCRM::create_local_vertices(), moab::NCHelperMPAS::create_local_vertices(), moab::NCHelperESMF::create_mesh(), moab::NCHelperGCRM::create_mesh(), moab::NCHelperMPAS::create_mesh(), moab::NCHelperScrip::create_mesh(), moab::ParallelComm::gather_data(), moab::ReadHDF5::load_file_partial(), moab::merge_ranged_ids(), moab::PACK_RANGE(), moab::ReadDamsel::process_ent_info(), moab::RANGE_SIZE(), moab::ReadHDF5::read_elems(), moab::ScdNCHelper::read_scd_variables_to_nonset_allocate(), moab::NCHelperGCRM::read_ucd_variables_to_nonset(), moab::NCHelperHOMME::read_ucd_variables_to_nonset(), moab::NCHelperMPAS::read_ucd_variables_to_nonset(), moab::NCHelperGCRM::read_ucd_variables_to_nonset_allocate(), moab::NCHelperHOMME::read_ucd_variables_to_nonset_allocate(), moab::NCHelperMPAS::read_ucd_variables_to_nonset_allocate(), moab::NCWriteGCRM::write_nonset_variables(), moab::NCWriteHOMME::write_nonset_variables(), and moab::NCWriteMPAS::write_nonset_variables().
|
inline |
return the beginning const reverse iterator of this range
Definition at line 866 of file Range.hpp.
References mHead, and moab::Range::PairNode::mPrev.
Referenced by moab::SpatialLocator::add_elems(), MetisPartitioner::assemble_taggedsets_graph(), moab::DualTool::atomic_pillow(), moab::ParallelComm::check_clean_iface(), moab::WriteHDF5Parallel::communicate_shared_set_ids(), SphereDecomp::compute_nodes(), moab::ParCommGraph::compute_partition(), moab::DualTool::construct_dual_cells(), moab::DualTool::construct_dual_edges(), moab::DualTool::construct_dual_faces(), moab::DualTool::construct_dual_vertices(), moab::ParallelComm::create_iface_pc_links(), moab::SpatialLocator::create_tree(), moab::DataCoupler::DataCoupler(), moab::Core::delete_entities(), moab::DualTool::delete_whole_dual(), DeformMeshRemap::find_other_sets(), moab::DualTool::foc_get_ents(), index(), moab::OrientedBoxTreeTool::join_trees(), moab::WriteGMV::local_write_mesh(), moab::SpatialLocator::locate_points(), moab::DualTool::next_loop_vertex(), moab::GeomTopoTool::other_entity(), moab::ParallelComm::pack_entities(), ZoltanPartitioner::partition_owned_cells(), moab::LloydSmoother::perform_smooth(), moab::ReadTemplate::read_elements(), DeformMeshRemap::read_file(), moab::Tqdcfr::read_nodes(), moab::ParallelComm::resolve_shared_ents(), moab::DualTool::rev_atomic_pillow(), moab::ParallelComm::set_pstatus_entities(), moab::SpatialLocator::SpatialLocator(), moab::MeshTopoUtil::split_entities_manifold(), tag_depth(), moab::ParallelMergeMesh::TagSharedElements(), moab::ParallelComm::unpack_sets(), moab::WriteCCMIO::write_external_faces(), and moab::WriteGmsh::write_file().
|
inline |
return the ending const reverse iterator for this range
Definition at line 878 of file Range.hpp.
References mHead.
Referenced by moab::ParallelComm::check_clean_iface(), moab::WriteHDF5Parallel::communicate_shared_set_ids(), moab::Core::delete_entities(), moab::DualTool::delete_whole_dual(), moab::OrientedBoxTreeTool::join_trees(), moab::DualTool::split_pair_nonmanifold(), tag_depth(), and moab::WriteGmsh::write_file().
void moab::Range::sanity_check | ( | ) | const |
check for internal consistency
Definition at line 554 of file Range.cpp.
References empty(), mHead, moab::Range::PairNode::mNext, and moab::Range::PairNode::mPrev.
size_t moab::Range::size | ( | ) | const |
return the number of values this Ranges represents
returns the number of values this list represents
Definition at line 77 of file Range.cpp.
References mHead, and moab::Range::PairNode::mNext.
Referenced by moab::TempestOnlineMap::ApplyWeights(), moab::IntxAreaUtils::area_on_sphere(), MetisPartitioner::assemble_graph(), ZoltanPartitioner::assemble_graph(), MetisPartitioner::assemble_taggedsets_graph(), moab::DualTool::atomic_pillow(), moab::ParallelComm::augment_default_sets_with_ghosts(), moab::ReadMCNP5::average_with_existing_tally(), moab::AdaptiveKDTree::best_subdivision_plane(), moab::AdaptiveKDTree::best_subdivision_snap_plane(), moab::AdaptiveKDTree::best_vertex_median_plane(), moab::AdaptiveKDTree::best_vertex_sample_plane(), moab::WriteSTL::binary_write_triangles(), moab::FBEngine::boundary_nodes_on_face(), moab::MeshGeneration::BrickInstance(), moab::OrientedBoxTreeTool::build_tree(), moab::AxisBox::calculate(), moab::WriteHDF5::check_dense_format_tag(), moab::DualTool::check_dual_adjs(), moab::DualTool::check_dual_equiv_edges(), moab::AEntityFactory::check_equiv_entities(), moab::NCHelperGCRM::check_existing_mesh(), moab::NCHelperHOMME::check_existing_mesh(), moab::NCHelperMPAS::check_existing_mesh(), moab::HalfFacetRep::check_mixed_entity_type(), moab::GeomTopoTool::check_model(), moab::ParallelComm::check_sent_ents(), moab::Skinner::classify_2d_boundary(), moab::FBEngine::clean(), moab::ParallelComm::clean_shared_tags(), moab::OrientedBoxTreeTool::closest_to_location(), moab::NCWriteGCRM::collect_mesh_info(), moab::NCWriteHOMME::collect_mesh_info(), moab::NCWriteMPAS::collect_mesh_info(), moab::NCWriteGCRM::collect_variable_data(), moab::NCWriteHOMME::collect_variable_data(), moab::NCWriteMPAS::collect_variable_data(), moab::ParallelComm::collective_sync_partition(), moab::WriteHDF5Parallel::communicate_shared_set_data(), moab::WriteHDF5Parallel::communicate_shared_set_ids(), compactness(), compute_dual_mesh(), moab::ParCommGraph::compute_partition(), moab::TempestRemapper::ComputeGlobalLocalMaps(), moab::TempestRemapper::ComputeOverlapMesh(), moab::DualTool::construct_dual(), moab::NestedRefine::construct_hm_1D(), moab::NestedRefine::construct_hm_2D(), moab::DualTool::construct_hp_parent_child(), moab::DualTool::construct_new_hyperplane(), moab::TempestRemapper::ConstructCoveringSet(), moab::TempestRemapper::convert_mesh_to_tempest_private(), moab::Tqdcfr::convert_nodesets_sidesets(), moab::TempestRemapper::convert_overlap_mesh_sorted_by_source(), moab::SpectralMeshTool::convert_to_coarse(), moab::TempestRemapper::ConvertMeshToTempest(), moab::NestedRefine::copy_vertices_from_prev_level(), moab::HalfFacetRep::count_subentities(), moab::ReadCCMIO::create_cell_from_faces(), moab::ReadMCNP5::create_elements(), moab::ReadIDEAS::create_elements(), moab::NCHelperMPAS::create_gather_set_cells(), moab::ParallelComm::create_iface_pc_links(), moab::ReadABAQUS::create_instance_of_part(), moab::NCHelperGCRM::create_local_edges(), moab::NCHelperMPAS::create_local_edges(), moab::NCHelperESMF::create_local_vertices(), moab::NCHelperGCRM::create_local_vertices(), moab::NCHelperMPAS::create_local_vertices(), moab::NCHelperESMF::create_mesh(), moab::NCHelperGCRM::create_mesh(), moab::NCHelperHOMME::create_mesh(), moab::NCHelperMPAS::create_mesh(), moab::NCHelperScrip::create_mesh(), moab::WriteHDF5Parallel::create_meshset_tables(), moab::FBEngine::create_new_gedge(), moab::WriteHDF5Parallel::create_node_table(), moab::ReadParallel::create_partition_sets(), moab::ScdNCHelper::create_quad_coordinate_tag(), moab::WriteHDF5Parallel::create_tag_tables(), moab::Coupler::create_tuples(), moab::ReadMCNP5::create_vertices(), moab::Intx2Mesh::createTags(), moab::IntxUtils::deep_copy_set_with_quads(), DeformMeshRemap::deform_master(), moab::ParallelComm::delete_entities(), moab::ReadHDF5::delete_non_side_elements(), moab::ReadParallel::delete_nonlocal_entities(), moab::GeomTopoTool::delete_obb_tree(), moab::HalfFacetRep::determine_incident_halfedges(), moab::HalfFacetRep::determine_sibling_halfedges(), moab::HalfFacetRep::determine_sibling_halffaces(), moab::HalfFacetRep::determine_sibling_halfverts(), do_test_mode(), moab::ParCommGraph::dump_comm_information(), moab::GeomTopoTool::duplicate_model(), moab::NestedRefine::estimate_hm_storage(), moab::WriteHDF5Parallel::exchange_file_ids(), moab::ParallelComm::exchange_ghost_cells(), moab::ParallelComm::exchange_owned_mesh(), moab::ParallelComm::exchange_tags(), DeformMeshRemap::execute(), moab::ParallelComm::filter_pstatus(), moab::Intx2Mesh::filterByMask(), moab::Tree::find_all_trees(), moab::ScdInterface::find_boxes(), moab::Skinner::find_geometric_skin(), moab::MergeMesh::find_merged_to(), moab::Skinner::find_skin_scd(), moab::HalfFacetRep::find_total_edges_2d(), moab::HalfFacetRep::find_total_edges_faces_3d(), moab::DualTool::foc_get_ents(), moab::ParCommGraph::form_mesh_from_tuples(), moab::ParCommGraph::form_tuples_to_migrate_mesh(), moab::DualTool::fs_get_quad_loops(), moab::DualTool::fs_get_quads(), moab::DualTool::fsr_get_fourth_quad(), moab::ParallelComm::gather_data(), moab::WriteSLAC::gather_interior_exterior(), moab::WriteNCDF::gather_mesh_information(), moab::WriteSLAC::gather_mesh_information(), moab::WriteTemplate::gather_mesh_information(), moab::ReadUtil::gather_related_ents(), moab::GeomTopoTool::geometrize_surface_set(), moab::ScdBox::get_adj_edge_or_face(), moab::MeshTopoUtil::get_average_position(), moab::DualTool::get_cell_points(), moab::Core::get_connectivity(), moab::Core::get_connectivity_by_type(), moab::AEntityFactory::get_down_adjacency_elements_poly(), moab::DualTool::get_dual_entities(), moab::ReadABAQUS::get_elements_by_id(), moab::Core::get_entities_by_handle(), moab::WriteCCMIO::get_gids(), moab::Coupler::get_gl_points_on_elements(), moab::DualTool::get_graphics_points(), get_imesh_mesh(), moab::ParallelComm::get_interface_procs(), moab::ParallelData::get_interface_sets(), moab::ParallelComm::get_local_handles(), moab::MeshTopoUtil::get_manifold(), moab::Coupler::get_matching_entities(), moab::Tqdcfr::get_mesh_entities(), moab::ReorderTool::get_new_handles(), moab::ReadABAQUS::get_nodes_by_id(), moab::WriteHDF5::get_num_sparse_tagged_entities(), moab::Core::get_number_entities_by_type_and_tag(), moab::DualTool::get_opposite_verts(), moab::ReadHDF5::get_partition(), moab::ParallelComm::get_pstatus_entities(), moab::ParallelComm::get_remote_handles(), moab::ReorderTool::get_reordered_handles(), moab::ReadABAQUS::get_set_elements_by_name(), moab::ReadABAQUS::get_set_nodes(), moab::WriteHDF5::get_tag_data_length(), get_vartag_data(), moab::Core::get_vertex_coordinates(), moab::NestedRefine::get_vertex_duplicates(), moab::TempestRemapper::GetIMasks(), moab::FBEngine::getNumOfType(), moab::TempestRemapper::GetOverlapAugmentedEntities(), moab::FBEngine::getPntRayIntsct(), hcFilter(), iMOAB_GetBlockID(), iMOAB_GetBlockInfo(), iMOAB_GetDoubleTagStorage(), iMOAB_GetElementConnectivity(), iMOAB_GetElementID(), iMOAB_GetElementOwnership(), iMOAB_GetIntTagStorage(), iMOAB_GetMeshInfo(), iMOAB_GetNeighborElements(), iMOAB_GetPointerToSurfaceBC(), iMOAB_GetPointerToVertexBC(), iMOAB_GetVertexID(), iMOAB_GetVertexOwnership(), iMOAB_GetVisibleVerticesCoordinates(), iMOAB_SetDoubleTagStorage(), iMOAB_SetDoubleTagStorageWithGid(), iMOAB_SetIntTagStorage(), iMOAB_UpdateMeshInfo(), ZoltanPartitioner::include_closure(), moab::SmoothFace::init_gradient(), moab::HalfFacetRep::initialize(), moab::LloydSmoother::initialize(), moab::HiReconstruction::initialize(), moab::WriteNCDF::initialize_exodus_file(), moab::WriteSLAC::initialize_file(), moab::FBEngine::initializeSmoothing(), moab::GeomTopoTool::insert_in_tree(), moab::DataCoupler::interpolate(), moab::AdaptiveKDTree::intersect_children_with_elems(), moab::Intx2Mesh::intersect_meshes(), moab::Intx2Mesh::intersect_meshes_kdtree(), moab::SmoothCurve::is_periodic(), laplacianFilter(), TriStats::leaf(), moab::ReadParallel::load_file(), moab::ReadVtk::load_file(), moab::ReadSms::load_file_impl(), moab::ReadHDF5::load_file_partial(), moab::ReadMCNP5::load_one_file(), moab::WriteGMV::local_write_mesh(), moab::SpatialLocator::locate_points(), moab::Coupler::locate_points(), main(), moab::WriteDamsel::map_sparse_tags(), moab::MergeMesh::merge_higher_dimensions(), merge_input_surfs(), moab::MergeMesh::merge_using_integer_tag(), moab::MeshSetSequence::num_dimension(), moab::MeshSetSequence::num_entities(), moab::GeomTopoTool::num_ents_of_dim(), moab::DenseTag::num_tagged_entities(), moab::MeshSetSequence::num_type(), obbstat_write(), obbvis_create(), moab::HiReconstruction::obtain_nring_ngbvs(), moab::DualTool::order_chord(), moab::GeomTopoTool::other_entity(), moab::ParallelComm::pack_entities(), moab::ParallelComm::pack_entity_seq(), moab::ParallelComm::pack_sets(), moab::ParallelComm::pack_tag(), moab::ParallelComm::packed_tag_size(), ZoltanPartitioner::partition_inferred_mesh(), ZoltanPartitioner::partition_mesh_and_geometry(), ZoltanPartitioner::partition_owned_cells(), perform_laplacian_smoothing(), perform_lloyd_relaxation(), moab::LloydSmoother::perform_smooth(), moab::HiReconstruction::polyfit3d_walf_curve_vertex(), moab::HiReconstruction::polyfit3d_walf_surf_vertex(), moab::AdaptiveKDTree::print(), moab::TreeNodePrinter::print_contents(), moab::print_type_sets(), moab::ReadDamsel::process_ent_info(), moab::ReadDamsel::process_entity_tags(), process_partition_file(), moab::putElementField(), moab::putSpectralElementField(), moab::putVertexField(), quads_to_tris(), moab::ReadHDF5VarLen::read_data(), moab::ReadABAQUS::read_element_list(), moab::ReadHDF5::read_elems(), DeformMeshRemap::read_file(), McnpData::read_mcnpfile(), moab::ReadABAQUS::read_node_list(), moab::ReadHDF5::read_nodes(), moab::ReadHDF5VarLen::read_offsets(), moab::ScdNCHelper::read_scd_variables_to_nonset_allocate(), moab::ReadHDF5::read_set_data(), moab::ReadHDF5::read_sets(), moab::ReadHDF5::read_sparse_tag_indices(), moab::ReadHDF5::read_tag_values_partial(), moab::NCHelperGCRM::read_ucd_variables_to_nonset_allocate(), moab::NCHelperHOMME::read_ucd_variables_to_nonset_allocate(), moab::NCHelperMPAS::read_ucd_variables_to_nonset_allocate(), moab::ParCommGraph::receive_mesh(), moab::ParCommGraph::receive_tag_values(), moab::ReorderTool::reorder_tag_data(), replace_surface(), moab::ParallelComm::resolve_shared_ents(), moab::GeomTopoTool::restore_topology_from_adjacency(), moab::GeomTopoTool::restore_topology_from_geometric_inclusion(), moab::DualTool::rev_atomic_pillow(), moab::ParallelComm::send_entities(), moab::ParCommGraph::send_tag_values(), moab::FBEngine::separate(), moab::GeomTopoTool::separate_by_dimension(), moab::WriteHDF5::serial_create_file(), moab::ParallelComm::set_pstatus_entities(), moab::ParCommGraph::set_split_ranges(), moab::ParallelComm::settle_intersection_points(), moab::FBEngine::smooth_new_intx_points(), moab::FBEngine::split_bedge_at_new_mesh_node(), moab::FBEngine::split_edge_at_mesh_node(), moab::FBEngine::split_internal_edge(), moab::FBEngine::split_surface_with_direction(), moab::MeshTopoUtil::star_next_entity(), moab::NestedRefine::subdivide_cells(), moab::NestedRefine::subdivide_tets(), tag_depth(), moab::Core::tag_get_by_ptr(), moab::ParallelComm::tag_iface_entities(), moab::Core::tag_set_by_ptr(), moab::ScdInterface::tag_shared_vertices(), test_spectral_hex(), test_spectral_quad(), moab::DualTool::traverse_hyperplane(), moab::ParallelComm::unpack_sets(), moab::ReadNCDF::update(), moab::NestedRefine::update_global_ahf_1D(), moab::NestedRefine::update_global_ahf_1D_sub(), moab::NestedRefine::update_global_ahf_2D(), moab::NestedRefine::update_global_ahf_2D_sub(), moab::NestedRefine::update_global_ahf_3D(), moab::ReorderTool::update_set_contents(), moab::Intx2MeshOnSphere::update_tracer_data(), moab::RayIntersectSets::visit(), moab::TreeNodePrinter::visit(), moab::ReadVtk::vtk_read_unstructured_grid(), MetisPartitioner::write_aggregationtag_partition(), moab::WriteCCMIO::write_cells_and_faces(), moab::WriteHDF5::write_dense_tag(), moab::WriteNCDF::write_elementblocks(), moab::WriteHDF5::write_elems(), moab::WriteVtk::write_elems(), moab::WriteDamsel::write_entities(), moab::WriteCCMIO::write_external_faces(), moab::WriteGmsh::write_file(), moab::WriteSmf::write_file(), moab::Core::write_file(), moab::WriteHDF5::write_file_impl(), moab::WriteSLAC::write_matsets(), moab::WriteHDF5::write_nodes(), moab::WriteCCMIO::write_nodes(), moab::WriteVtk::write_nodes(), moab::NCWriteGCRM::write_nonset_variables(), moab::ScdNCWriteHelper::write_nonset_variables(), moab::NCWriteHOMME::write_nonset_variables(), moab::NCWriteMPAS::write_nonset_variables(), MetisPartitioner::write_partition(), ZoltanPartitioner::write_partition(), moab::WriteNCDF::write_poly_faces(), moab::WriteHDF5::write_sets(), moab::WriteHDF5::write_sparse_ids(), moab::WriteHDF5::write_sparse_tag(), moab::WriteHDF5::write_tag_values(), DeformMeshRemap::write_to_coords(), moab::WriteHDF5::write_var_len_indices(), moab::WriteHDF5::write_var_len_tag(), and moab::WriteDamsel::write_vertices().
const std::string moab::Range::str_rep | ( | const char * | indent_prefix = NULL | ) | const |
for debugging
Definition at line 579 of file Range.cpp.
References const_pair_begin(), const_pair_end(), empty(), moab::CN::EntityTypeName(), moab::ID_FROM_HANDLE(), and moab::TYPE_FROM_HANDLE().
Referenced by print().
Range moab::Range::subset_by_dimension | ( | int | dim | ) | const |
return a subset of this range, by dimension
return a subset of this range, by type
Definition at line 966 of file Range.cpp.
References begin(), moab::CREATE_HANDLE(), end(), moab::GeomUtil::first(), insert(), lower_bound(), and moab::CN::TypeDimensionMap.
Referenced by moab::ParCommGraph::compute_partition(), moab::SpectralMeshTool::convert_to_coarse(), DeformMeshRemap::find_other_sets(), moab::WriteNCDF::gather_mesh_information(), moab::ParallelComm::get_ghosted_entities(), get_max_volume(), main(), and DeformMeshRemap::read_file().
Range moab::Range::subset_by_type | ( | EntityType | t | ) | const |
return a subset of this range, by type
Definition at line 957 of file Range.cpp.
References equal_range(), and insert().
Referenced by moab::ParallelComm::add_verts(), moab::HalfFacetRep::check_mixed_entity_type(), moab::ReadParallel::delete_nonlocal_entities(), moab::ParallelComm::exchange_owned_meshs(), iMOAB_DeregisterApplication(), moab::Coupler::locate_points(), moab::Tqdcfr::read_nodes(), moab::ParCommGraph::receive_mesh(), and moab::WriteVtk::write_elems().
void moab::Range::swap | ( | Range & | range | ) |
swap the contents of this range with another one
swap the contents of this range with another one THIS FUNCTION MUST NOT BE INLINED, THAT WILL ELIMINATE RANGE_EMPTY AND THIS_EMPTY BY SUBSTITUTION AND THE FUNCTION WON'T WORK RIGHT!
Definition at line 937 of file Range.cpp.
References mHead, moab::Range::PairNode::mNext, and moab::Range::PairNode::mPrev.
Referenced by moab::AdaptiveKDTree::best_subdivision_plane(), moab::AdaptiveKDTree::best_subdivision_snap_plane(), moab::AdaptiveKDTree::best_vertex_median_plane(), moab::AdaptiveKDTree::best_vertex_sample_plane(), moab::OrientedBoxTreeTool::build_tree(), moab::ReadParallel::delete_nonlocal_entities(), dot_nodes(), moab::MeshTopoUtil::equivalent_entities(), moab::ParallelComm::filter_pstatus(), moab::Skinner::find_skin(), moab::Core::get_entities_by_type_and_tag(), moab::ReadHDF5::get_subset_ids(), moab::Core::get_vertices(), moab::DualTool::next_loop_vertex(), moab::WriteHDF5Parallel::remove_remote_entities(), moab::WriteHDF5Parallel::remove_remote_sets(), moab::ParallelComm::resolve_shared_sets(), MetisPartitioner::write_aggregationtag_partition(), MetisPartitioner::write_partition(), and ZoltanPartitioner::write_partition().
|
static |
Definition at line 832 of file Range.cpp.
References moab::GeomUtil::first(), and lower_bound().
Referenced by MetisPartitioner::assemble_taggedsets_graph(), moab::ParallelComm::check_my_shared_handles(), moab::OrientedBox::compute_from_vertices(), moab::Skinner::find_skin_vertices_3D(), moab::ReadHDF5::get_partition(), moab::ParallelComm::get_shared_entities(), moab::Core::get_vertices(), moab::ParallelComm::resolve_shared_ents(), and moab::ParallelMergeMesh::TagSharedElements().
|
inline |
Range::const_iterator moab::Range::upper_bound | ( | EntityType | type | ) | const |
Definition at line 852 of file Range.cpp.
References begin(), moab::CREATE_HANDLE(), end(), and lower_bound().
Range::const_iterator moab::Range::upper_bound | ( | EntityType | type, |
const_iterator | first | ||
) | const |
Definition at line 859 of file Range.cpp.
References moab::CREATE_HANDLE(), end(), moab::GeomUtil::first(), and lower_bound().
|
protected |
the head of the list that contains pairs that represent the ranges this list is sorted and unique at all times
Definition at line 392 of file Range.hpp.
Referenced by back(), begin(), clear(), contains(), delete_pair_node(), empty(), end(), erase(), find(), front(), get_memory_use(), insert(), operator=(), pop_back(), pop_front(), Range(), rbegin(), rend(), sanity_check(), size(), and swap().