#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 |
default constructor
Definition at line 846 of file Range.hpp.
847 {
848 // set the head node to point to itself
849 mHead.mNext = mHead.mPrev = &mHead;
850 mHead.first = mHead.second = 0;
851 }
References mHead, moab::Range::PairNode::mNext, and moab::Range::PairNode::mPrev.
moab::Range::Range | ( | const Range & | copy | ) |
copy constructor
Definition at line 179 of file Range.cpp.
180 {
181 // set the head node to point to itself
182 mHead.mNext = mHead.mPrev = &mHead;
183 mHead.first = mHead.second = 0;
184
185 const PairNode* copy_node = copy.mHead.mNext;
186 PairNode* new_node = &mHead;
187 for( ; copy_node != &( copy.mHead ); copy_node = copy_node->mNext )
188 {
189 PairNode* tmp_node = alloc_pair( new_node->mNext, new_node, copy_node->first, copy_node->second );
190 new_node->mNext->mPrev = tmp_node;
191 new_node->mNext = tmp_node;
192 new_node = tmp_node;
193 }
194 }
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.
173 { 174 mHead.mNext = mHead.mPrev = alloc_pair( &mHead, &mHead, val1, val2 ); 175 mHead.first = mHead.second = 0; 176 }
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.
885 {
886 return empty() || ( CN::Dimension( TYPE_FROM_HANDLE( front() ) ) == dimension &&
887 CN::Dimension( TYPE_FROM_HANDLE( back() ) ) == dimension );
888 }
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.
880 {
881 return empty() || ( TYPE_FROM_HANDLE( front() ) == type && TYPE_FROM_HANDLE( back() ) == type );
882 }
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(), iMOAB_ReceiveMesh(), 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.
913 {
914 return mHead.mPrev->second;
915 }
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.
861 {
862 return const_iterator( mHead.mNext, mHead.mNext->first );
863 }
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::IntxUtils::max_diagonal(), 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.
198 {
199 PairNode* tmp_node = mHead.mNext;
200 while( tmp_node != &mHead )
201 {
202 PairNode* to_delete = tmp_node;
203 tmp_node = tmp_node->mNext;
204 free_pair( to_delete );
205 }
206 mHead.mNext = &mHead;
207 mHead.mPrev = &mHead;
208 }
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_MergeVertices(), 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.
953 {
954 unsigned int num_ents = size();
955 return ( num_ents ? ( (double)get_memory_use() / (double)( num_ents * sizeof( EntityHandle ) ) ) : -1 );
956 }
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.
749 {
750 return const_pair_iterator( mHead.mNext );
751 }
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.
753 {
754 return const_pair_iterator( &mHead );
755 }
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.
1006 {
1007 if( othr.empty() ) return true;
1008 if( empty() ) return false;
1009
1010 // neither range is empty, so both have valid pair nodes
1011 // other than dummy mHead
1012 const PairNode* this_node = mHead.mNext;
1013 const PairNode* othr_node = othr.mHead.mNext;
1014 for( ;; )
1015 {
1016 // Loop while the node in this list is entirely before
1017 // the node in the other list.
1018 while( this_node->second < othr_node->first )
1019 {
1020 this_node = this_node->mNext;
1021 if( this_node == &mHead ) return false;
1022 }
1023 // If other node is not entirely contained in this node
1024 // then other list is not contained in this list
1025 if( this_node->first > othr_node->first ) break;
1026 // Loop while other node is entirely contained in this node.
1027 while( othr_node->second <= this_node->second )
1028 {
1029 othr_node = othr_node->mNext;
1030 if( othr_node == &othr.mHead ) return true;
1031 }
1032 // If other node overlapped end of this node
1033 if( othr_node->first <= this_node->second ) break;
1034 }
1035
1036 // should be unreachable
1037 return false;
1038 }
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.
472 {
473 if( node != &mHead )
474 { // pop_front() and pop_back() rely on this check
475 node->mPrev->mNext = node->mNext;
476 node->mNext->mPrev = node->mPrev;
477 free_pair( node );
478 }
479 }
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.
886 {
887 return ( mHead.mNext == &mHead );
888 }
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::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_MergeVertices(), iMOAB_ResolveSharedEntities(), 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.
873 {
874 return const_iterator( &mHead, mHead.first );
875 }
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::IntxUtils::max_diagonal(), 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.
868 {
869 std::pair< Range::const_iterator, Range::const_iterator > result;
870 int err;
871 EntityHandle handle = CREATE_HANDLE( type, 0, err );
872 result.first = err ? end() : lower_bound( begin(), end(), handle );
873 // if (type+1) overflows, err will be true and we return end().
874 handle = CREATE_HANDLE( type + 1, 0, err );
875 result.second = err ? end() : lower_bound( result.first, end(), handle );
876 return result;
877 }
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.
365 {
366 // one of a few things could happen
367 // 1. shrink a range
368 // 2. split a range
369 // 3. remove a range
370
371 if( iter == end() ) return end();
372
373 // the iterator most likely to be returned
374 iterator new_iter = iter;
375 ++new_iter;
376
377 PairNode* kter = iter.mNode;
378
379 // just remove the range
380 if( kter->first == kter->second )
381 {
382 kter->mNext->mPrev = kter->mPrev;
383 kter->mPrev->mNext = kter->mNext;
384 free_pair( kter );
385 return new_iter;
386 }
387 // shrink it
388 else if( kter->first == iter.mValue )
389 {
390 kter->first++;
391 return new_iter;
392 }
393 // shrink it the other way
394 else if( kter->second == iter.mValue )
395 {
396 kter->second--;
397 return new_iter;
398 }
399 // split the range
400 else
401 {
402 PairNode* new_node = alloc_pair( iter.mNode->mNext, iter.mNode, iter.mValue + 1, kter->second );
403 new_node->mPrev->mNext = new_node->mNext->mPrev = new_node;
404 iter.mNode->second = iter.mValue - 1;
405 new_iter = const_iterator( new_node, new_node->first );
406 return new_iter;
407 }
408 }
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.
412 {
413 iterator result;
414
415 if( iter1.mNode == iter2.mNode )
416 {
417 if( iter2.mValue <= iter1.mValue )
418 {
419 // empty range OK, otherwise invalid input
420 return iter2;
421 }
422
423 // If both iterators reference the same pair node, then
424 // we're either removing values from the front of the
425 // node or splitting the node. We can never be removing
426 // the last value in the node in this case because iter2
427 // points *after* the last entry to be removed.
428
429 PairNode* node = iter1.mNode;
430 if( iter1.mValue == node->first )
431 {
432 node->first = iter2.mValue;
433 result = iter2;
434 }
435 else
436 {
437 PairNode* new_node = alloc_pair( node->mNext, node, iter2.mValue, node->second );
438 new_node->mNext->mPrev = new_node;
439 new_node->mPrev->mNext = new_node;
440 node->second = iter1.mValue - 1;
441 result = iterator( new_node, new_node->first );
442 }
443 }
444 else
445 {
446 if( iter1.mNode == &mHead ) return iter1;
447 PairNode* dn = iter1.mNode;
448 if( iter1.mValue > dn->first )
449 {
450 dn->second = iter1.mValue - 1;
451 dn = dn->mNext;
452 }
453 if( iter2.mNode != &mHead ) iter2.mNode->first = iter2.mValue;
454 while( dn != iter2.mNode )
455 {
456 PairNode* dead = dn;
457 dn = dn->mNext;
458
459 dead->mPrev->mNext = dead->mNext;
460 dead->mNext->mPrev = dead->mPrev;
461 // dead->mPrev = dead->mNext = 0;
462 delete_pair_node( dead );
463 }
464
465 result = iter2;
466 }
467
468 return result;
469 }
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.
511 {
512 // iterator through the list
513 PairNode* iter = mHead.mNext;
514 for( ; iter != &mHead && ( val > iter->second ); iter = iter->mNext )
515 ;
516 return ( ( iter->second >= val ) && ( iter->first <= val ) ) ? const_iterator( iter, val ) : end();
517 }
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.
908 {
909 return mHead.mNext->first;
910 }
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.
998 {
999 unsigned long result = 0;
1000 for( const PairNode* n = mHead.mNext; n != &mHead; n = n->mNext )
1001 result += sizeof( PairNode );
1002 return result;
1003 }
References mHead, and moab::Range::PairNode::mNext.
Referenced by compactness().
|
inline |
Definition at line 936 of file Range.hpp.
937 {
938 if( handle < *begin() || handle > *rbegin() ) return -1;
939
940 unsigned int i = 0;
941 Range::const_pair_iterator pit = const_pair_begin();
942 while( handle > ( *pit ).second && pit != const_pair_end() )
943 {
944 i += ( *pit ).second - ( *pit ).first + 1;
945 ++pit;
946 }
947 if( handle < ( *pit ).first || pit == const_pair_end() ) return -1;
948
949 return i + handle - ( *pit ).first;
950 }
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.
296 {
297 if( val1 == 0 || val1 > val2 ) return end();
298
299 // Empty
300 if( mHead.mNext == &mHead )
301 {
302 assert( prev == end() );
303 PairNode* new_node = alloc_pair( &mHead, &mHead, val1, val2 );
304 mHead.mNext = mHead.mPrev = new_node;
305 return iterator( mHead.mNext, val1 );
306 }
307
308 PairNode* iter = prev.mNode;
309 // If iterator is at end, set it to last.
310 // Thus if the hint was to append, we start searching
311 // at the end of the list.
312 if( iter == &mHead ) iter = mHead.mPrev;
313 // if hint (prev) is past insert position, reset it to the beginning.
314 if( iter != &mHead && iter->first > val2 + 1 ) iter = mHead.mNext;
315
316 // If hint is bogus then search backwards
317 while( iter != mHead.mNext && iter->mPrev->second >= val1 - 1 )
318 iter = iter->mPrev;
319
320 // Input range is before beginning?
321 if( iter->mPrev == &mHead && val2 < iter->first - 1 )
322 {
323 PairNode* new_node = alloc_pair( iter, &mHead, val1, val2 );
324 mHead.mNext = iter->mPrev = new_node;
325 return iterator( mHead.mNext, val1 );
326 }
327
328 // Find first intersecting list entry, or the next entry
329 // if none intersects.
330 while( iter != &mHead && iter->second + 1 < val1 )
331 iter = iter->mNext;
332
333 // Need to insert new pair (don't intersect any existing pair)?
334 if( iter == &mHead || iter->first - 1 > val2 )
335 {
336 PairNode* new_node = alloc_pair( iter, iter->mPrev, val1, val2 );
337 iter->mPrev = iter->mPrev->mNext = new_node;
338 return iterator( iter->mPrev, val1 );
339 }
340
341 // Make the first intersecting pair the union of itself with [val1,val2]
342 if( iter->first > val1 ) iter->first = val1;
343 if( iter->second >= val2 ) return iterator( iter, val1 );
344 iter->second = val2;
345
346 // Merge any remaining pairs that intersect [val1,val2]
347 while( iter->mNext != &mHead && iter->mNext->first <= val2 + 1 )
348 {
349 PairNode* dead = iter->mNext;
350 iter->mNext = dead->mNext;
351 dead->mNext->mPrev = iter;
352
353 if( dead->second > val2 ) iter->second = dead->second;
354 free_pair( dead );
355 }
356
357 return iterator( iter, val1 );
358 }
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.
230 {
231 // don't allow zero-valued handles in Range
232 if( val == 0 ) return end();
233
234 // if this is empty, just add it and return an iterator to it
235 if( &mHead == mHead.mNext )
236 {
237 mHead.mNext = mHead.mPrev = alloc_pair( &mHead, &mHead, val, val );
238 return iterator( mHead.mNext, val );
239 }
240
241 // find the location in the list where we can safely insert
242 // new items and keep it ordered
243 PairNode* hter = hint.mNode;
244 PairNode* jter = hter->first <= val ? hter : mHead.mNext;
245 for( ; ( jter != &mHead ) && ( jter->second < val ); jter = jter->mNext )
246 ;
247 PairNode* iter = jter;
248 jter = jter->mPrev;
249
250 // if this val is already in the list
251 if( ( iter->first <= val && iter->second >= val ) && ( iter != &mHead ) )
252 {
253 // return an iterator pointing to this location
254 return iterator( iter, val );
255 }
256
257 // one of a few things can happen at this point
258 // 1. this range needs to be backwardly extended
259 // 2. the previous range needs to be forwardly extended
260 // 3. a new range needs to be added
261
262 // extend this range back a bit
263 else if( ( iter->first == ( val + 1 ) ) && ( iter != &mHead ) )
264 {
265 iter->first = val;
266 // see if we need to merge two ranges
267 if( ( iter != mHead.mNext ) && ( jter->second == ( val - 1 ) ) )
268 {
269 jter->second = iter->second;
270 iter->mPrev->mNext = iter->mNext;
271 iter->mNext->mPrev = iter->mPrev;
272 free_pair( iter );
273 return iterator( jter, val );
274 }
275 else
276 {
277 return iterator( iter, val );
278 }
279 }
280 // extend the previous range forward a bit
281 else if( ( jter->second == ( val - 1 ) ) && ( iter != mHead.mNext ) )
282 {
283 jter->second = val;
284 return iterator( jter, val );
285 }
286 // make a new range
287 else
288 {
289 PairNode* new_node = alloc_pair( iter, iter->mPrev, val, val );
290 iter->mPrev = new_node->mPrev->mNext = new_node;
291 return iterator( new_node, val );
292 }
293 }
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.
524 {
525 if( begini == endi ) return;
526
527 PairNode* node = begini.mNode;
528 if( endi.mNode == node )
529 {
530 insert( *begini, ( *endi ) - 1 );
531 return;
532 }
533
534 Range::iterator hint = insert( *begini, node->second );
535 node = node->mNext;
536 while( node != endi.mNode )
537 {
538 hint = insert( hint, node->first, node->second );
539 node = node->mNext;
540 }
541
542 if( *endi > node->first )
543 {
544 if( *endi <= node->second )
545 insert( hint, node->first, *(endi)-1 );
546 else
547 insert( hint, node->first, node->second );
548 }
549 }
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 | ||
) |
Definition at line 959 of file Range.hpp.
960 {
961 size_t n = std::distance( begin_iter, end_iter );
962 EntityHandle* sorted = new EntityHandle[n];
963 std::copy( begin_iter, end_iter, sorted );
964 std::sort( sorted, sorted + n );
965 iterator hint = begin();
966 size_t i = 0;
967 while( i < n )
968 {
969 size_t j = i + 1;
970 while( j < n && sorted[j] == 1 + sorted[j - 1] )
971 ++j;
972 hint = insert( hint, sorted[i], sorted[i] + ( ( j - i ) - 1 ) );
973 i = j;
974 }
975 delete[] sorted;
976 return hint;
977 }
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.
810 {
811 // Find the first pair whose end is >= val
812 PairNode* iter;
813 for( iter = first.mNode; iter != last.mNode; iter = iter->mNext )
814 {
815 if( iter->second >= val )
816 {
817 // This is the correct pair. Either 'val' is in the range, or
818 // the range starts before 'val' and iter->first IS the lower_bound.
819 if( iter->first > val ) return const_iterator( iter, iter->first );
820 return const_iterator( iter, val );
821 }
822 }
823
824 if( iter->first >= val )
825 return const_iterator( iter, iter->first );
826 else if( *last > val )
827 return const_iterator( iter, val );
828 else
829 return last;
830 }
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.
840 {
841 int err;
842 EntityHandle handle = CREATE_HANDLE( type, 0, err );
843 return err ? end() : lower_bound( begin(), end(), handle );
844 }
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.
846 {
847 int err;
848 EntityHandle handle = CREATE_HANDLE( type, 0, err );
849 return err ? end() : lower_bound( first, end(), handle );
850 }
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.
347 {
348 insert( range.begin(), range.end() );
349 }
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::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.
912 {
913 const_pair_iterator iter = const_pair_begin();
914 while( iter != const_pair_end() && CN::Dimension( TYPE_FROM_HANDLE( ( *iter ).second ) ) < dim )
915 ++iter;
916
917 int junk;
918 unsigned count = 0;
919 for( ; iter != const_pair_end(); ++iter )
920 {
921 int start_dim = CN::Dimension( TYPE_FROM_HANDLE( ( *iter ).first ) );
922 int end_dim = CN::Dimension( TYPE_FROM_HANDLE( ( *iter ).second ) );
923 if( start_dim > dim ) break;
924
925 EntityHandle sh = start_dim < dim ? CREATE_HANDLE( CN::TypeDimensionMap[dim].first, 1, junk ) : ( *iter ).first;
926 EntityHandle eh =
927 end_dim > dim ? CREATE_HANDLE( CN::TypeDimensionMap[dim].second, MB_END_ID, junk ) : ( *iter ).second;
928 count += eh - sh + 1;
929 }
930
931 return count;
932 }
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.
891 {
892 const_pair_iterator iter = const_pair_begin();
893 while( iter != const_pair_end() && TYPE_FROM_HANDLE( ( *iter ).second ) < type )
894 ++iter;
895
896 unsigned count = 0;
897 for( ; iter != const_pair_end(); ++iter )
898 {
899 EntityType start_type = TYPE_FROM_HANDLE( ( *iter ).first );
900 EntityType end_type = TYPE_FROM_HANDLE( ( *iter ).second );
901 if( start_type > type ) break;
902
903 EntityID sid = start_type < type ? 1 : ID_FROM_HANDLE( ( *iter ).first );
904 EntityID eid = end_type > type ? MB_END_ID : ID_FROM_HANDLE( ( *iter ).second );
905 count += eid - sid + 1;
906 }
907
908 return count;
909 }
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.
731 {
732 const bool braindead = false;
733
734 if( braindead )
735 {
736 // brain-dead implementation right now
737 Range res( *this );
738 for( Range::const_iterator rit = range2.begin(); rit != range2.end(); ++rit )
739 res.erase( *rit );
740
741 return *this;
742 }
743 else
744 {
745 Range::pair_iterator r_it0 = this->pair_begin();
746 Range::const_pair_iterator r_it1 = range2.const_pair_begin();
747
748 // terminate the while loop when at least one "start" iterator is at the
749 // end of the list
750 while( r_it0 != this->end() && r_it1 != range2.end() )
751 {
752 // case a: pair wholly within subtracted pair
753 if( r_it0->first >= r_it1->first && r_it0->second <= r_it1->second )
754 {
755 Range::PairNode* rtmp = r_it0.node();
756 ++r_it0;
757 this->delete_pair_node( rtmp );
758 }
759 // case b: pair overlaps upper part of subtracted pair
760 else if( r_it0->first <= r_it1->second && r_it0->first >= r_it1->first )
761 {
762 r_it0->first = r_it1->second + 1;
763 ++r_it1;
764 }
765 // case c: pair overlaps lower part of subtracted pair
766 else if( r_it0->second >= r_it1->first && r_it0->second <= r_it1->second )
767 {
768 r_it0->second = r_it1->first - 1;
769 ++r_it0;
770 }
771 // case d: pair completely surrounds subtracted pair
772 else if( r_it0->first < r_it1->first && r_it0->second > r_it1->second )
773 {
774 Range::PairNode* new_node =
775 alloc_pair( r_it0.node(), r_it0.node()->mPrev, r_it0->first, r_it1->first - 1 );
776 new_node->mPrev->mNext = new_node->mNext->mPrev = new_node;
777 r_it0.node()->first = r_it1->second + 1;
778 ++r_it1;
779 }
780 else
781 {
782 while( r_it0->second < r_it1->first && r_it0 != this->end() )
783 ++r_it0;
784 if( r_it0 == this->end() ) break;
785 while( r_it1->second < r_it0->first && r_it1 != range2.end() )
786 ++r_it1;
787 }
788 }
789 return *this;
790 }
791 }
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.
211 {
212 clear();
213 const PairNode* copy_node = copy.mHead.mNext;
214 PairNode* new_node = &mHead;
215 for( ; copy_node != &( copy.mHead ); copy_node = copy_node->mNext )
216 {
217 PairNode* tmp_node = alloc_pair( new_node->mNext, new_node, copy_node->first, copy_node->second );
218 new_node->mNext->mPrev = tmp_node;
219 new_node->mNext = tmp_node;
220 new_node = tmp_node;
221 }
222 return *this;
223 }
References alloc_pair(), clear(), mHead, moab::Range::PairNode::mNext, and moab::Range::PairNode::mPrev.
|
inline |
|
inline |
Definition at line 739 of file Range.hpp.
740 {
741 return pair_iterator( mHead.mNext );
742 }
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.
757 {
758 return const_pair_iterator( mHead.mNext );
759 }
References moab::Range::PairNode::mNext.
|
inline |
Definition at line 743 of file Range.hpp.
744 {
745 return pair_iterator( &mHead );
746 }
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.
495 {
496 EntityHandle retval = back();
497 if( mHead.mPrev->first == mHead.mPrev->second ) // need to remove pair from range
498 delete_pair_node( mHead.mPrev );
499 else
500 --( mHead.mPrev->second ); // otherwise just adjust end value of pair
501
502 return retval;
503 }
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.
483 {
484 EntityHandle retval = front();
485 if( mHead.mNext->first == mHead.mNext->second ) // need to remove pair from range
486 delete_pair_node( mHead.mNext );
487 else
488 ++( mHead.mNext->first ); // otherwise just adjust start value of pair
489
490 return retval;
491 }
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.
619 {
620 print( std::cout, indent_prefix );
621 }
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.
980 {
981 size_t i = 0;
982 Range::const_pair_iterator pit;
983 for( pit = const_pair_begin(), i = 0; pit != const_pair_end(); ++pit, i++ )
984 ;
985
986 return i;
987 }
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.
867 {
868 return const_reverse_iterator( mHead.mPrev, mHead.mPrev->second );
869 }
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.
879 {
880 return const_reverse_iterator( &mHead, mHead.second );
881 }
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.
555 {
556 if( empty() ) return;
557
558 const PairNode* node = mHead.mNext;
559 std::vector< const PairNode* > seen_before;
560 bool stop_it = false;
561
562 for( ; stop_it == false; node = node->mNext )
563 {
564 // have we seen this node before?
565 assert( std::find( seen_before.begin(), seen_before.end(), node ) == seen_before.end() );
566 seen_before.push_back( node );
567
568 // is the connection correct?
569 assert( node->mNext->mPrev == node );
570
571 // are the values right?
572 assert( node->first <= node->second );
573 if( node != &mHead && node->mPrev != &mHead ) assert( node->mPrev->second < node->first );
574
575 if( node == &mHead ) stop_it = true;
576 }
577 }
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.
78 {
79 // go through each pair and add up the number of values
80 // we have.
81 size_t sz = 0;
82 for( PairNode* iter = mHead.mNext; iter != &mHead; iter = iter->mNext )
83 {
84 sz += ( ( iter->second - iter->first ) + 1 );
85 }
86 return sz;
87 }
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_ComputeCommGraph(), 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_ReceiveMesh(), iMOAB_ResolveSharedEntities(), iMOAB_SendMesh(), 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.
580 {
581 std::stringstream str_stream;
582 std::string indent_prefix_str;
583 if( NULL != indent_prefix )
584 {
585 indent_prefix_str += indent_prefix;
586 }
587
588 if( empty() )
589 {
590 str_stream << indent_prefix_str << "\tempty" << std::endl;
591 return str_stream.str().c_str();
592 }
593
594 for( const_pair_iterator i = const_pair_begin(); i != const_pair_end(); ++i )
595 {
596 EntityType t1 = TYPE_FROM_HANDLE( i->first );
597 EntityType t2 = TYPE_FROM_HANDLE( i->second );
598
599 str_stream << indent_prefix_str << "\t" << CN::EntityTypeName( t1 ) << " " << ID_FROM_HANDLE( i->first );
600 if( i->first != i->second )
601 {
602 str_stream << " - ";
603 if( t1 != t2 ) str_stream << CN::EntityTypeName( t2 ) << " ";
604 str_stream << ID_FROM_HANDLE( i->second );
605 }
606 str_stream << std::endl;
607 }
608
609 return str_stream.str();
610 }
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.
967 {
968 EntityHandle handle1 = CREATE_HANDLE( CN::TypeDimensionMap[d].first, 0 );
969 iterator st = lower_bound( begin(), end(), handle1 );
970
971 iterator en;
972 if( d < 4 )
973 { // dimension 4 is MBENTITYSET
974 EntityHandle handle2 = CREATE_HANDLE( CN::TypeDimensionMap[d + 1].first, 0 );
975 en = lower_bound( st, end(), handle2 );
976 }
977 else
978 {
979 en = end();
980 }
981
982 Range result;
983 result.insert( st, en );
984 return result;
985 }
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.
958 {
959 Range result;
960 std::pair< const_iterator, const_iterator > iters = equal_range( t );
961 result.insert( iters.first, iters.second );
962 return result;
963 }
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(), iMOAB_ReceiveMesh(), 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.
938 {
939 // update next/prev nodes of head of both ranges
940 bool range_empty = ( range.mHead.mNext == &( range.mHead ) );
941 bool this_empty = ( mHead.mNext == &mHead );
942
943 range.mHead.mNext->mPrev = ( range_empty ? &( range.mHead ) : &mHead );
944 range.mHead.mPrev->mNext = ( range_empty ? &( range.mHead ) : &mHead );
945 mHead.mNext->mPrev = ( this_empty ? &mHead : &( range.mHead ) );
946 mHead.mPrev->mNext = ( this_empty ? &mHead : &( range.mHead ) );
947
948 // switch data in head nodes of both ranges
949 PairNode *range_next = range.mHead.mNext, *range_prev = range.mHead.mPrev;
950 range.mHead.mNext = ( this_empty ? &( range.mHead ) : mHead.mNext );
951 range.mHead.mPrev = ( this_empty ? &( range.mHead ) : mHead.mPrev );
952 mHead.mNext = ( range_empty ? &mHead : range_next );
953 mHead.mPrev = ( range_empty ? &mHead : range_prev );
954 }
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.
833 {
834 Range::const_iterator result = lower_bound( first, last, val );
835 if( result != last && *result == val ) ++result;
836 return result;
837 }
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.
853 {
854 // if (type+1) overflows, err will be true and we return end().
855 int err;
856 EntityHandle handle = CREATE_HANDLE( type + 1, 0, err );
857 return err ? end() : lower_bound( begin(), end(), handle );
858 }
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.
860 {
861 // if (type+1) overflows, err will be true and we return end().
862 int err;
863 EntityHandle handle = CREATE_HANDLE( type + 1, 0, err );
864 return err ? end() : lower_bound( first, end(), handle );
865 }
References moab::CREATE_HANDLE(), end(), moab::GeomUtil::first(), and lower_bound().
intersect two ranges, placing the results in the return range
Definition at line 626 of file Range.cpp.
627 {
628 Range::const_pair_iterator r_it[2] = { range1.const_pair_begin(), range2.const_pair_begin() };
629 EntityHandle low_it, high_it;
630
631 Range lhs;
632 Range::iterator hint = lhs.begin();
633
634 // terminate the while loop when at least one "start" iterator is at the
635 // end of the list
636 while( r_it[0] != range1.end() && r_it[1] != range2.end() )
637 {
638
639 if( r_it[0]->second < r_it[1]->first )
640 // 1st subrange completely below 2nd subrange
641 ++r_it[0];
642 else if( r_it[1]->second < r_it[0]->first )
643 // 2nd subrange completely below 1st subrange
644 ++r_it[1];
645
646 else
647 {
648 // else ranges overlap; first find greater start and lesser end
649 low_it = MAX( r_it[0]->first, r_it[1]->first );
650 high_it = MIN( r_it[0]->second, r_it[1]->second );
651
652 // insert into result
653 hint = lhs.insert( hint, low_it, high_it );
654
655 // now find bounds of this insertion and increment corresponding iterator
656 if( high_it == r_it[0]->second ) ++r_it[0];
657 if( high_it == r_it[1]->second ) ++r_it[1];
658 }
659 }
660
661 return lhs;
662 }
subtract range2 from this, placing the results in the return range
Definition at line 664 of file Range.cpp.
665 {
666 const bool braindead = false;
667
668 if( braindead )
669 {
670 // brain-dead implementation right now
671 Range res( range1 );
672 for( Range::const_iterator rit = range2.begin(); rit != range2.end(); ++rit )
673 res.erase( *rit );
674
675 return res;
676 }
677 else
678 {
679 Range lhs( range1 );
680
681 Range::pair_iterator r_it0 = lhs.pair_begin();
682 Range::const_pair_iterator r_it1 = range2.const_pair_begin();
683
684 // terminate the while loop when at least one "start" iterator is at the
685 // end of the list
686 while( r_it0 != lhs.end() && r_it1 != range2.end() )
687 {
688 // case a: pair wholly within subtracted pair
689 if( r_it0->first >= r_it1->first && r_it0->second <= r_it1->second )
690 {
691 Range::PairNode* rtmp = r_it0.node();
692 ++r_it0;
693 lhs.delete_pair_node( rtmp );
694 }
695 // case b: pair overlaps upper part of subtracted pair
696 else if( r_it0->first <= r_it1->second && r_it0->first >= r_it1->first )
697 {
698 r_it0->first = r_it1->second + 1;
699 ++r_it1;
700 }
701 // case c: pair overlaps lower part of subtracted pair
702 else if( r_it0->second >= r_it1->first && r_it0->second <= r_it1->second )
703 {
704 r_it0->second = r_it1->first - 1;
705 ++r_it0;
706 }
707 // case d: pair completely surrounds subtracted pair
708 else if( r_it0->first < r_it1->first && r_it0->second > r_it1->second )
709 {
710 Range::PairNode* new_node =
711 alloc_pair( r_it0.node(), r_it0.node()->mPrev, r_it0->first, r_it1->first - 1 );
712 new_node->mPrev->mNext = new_node->mNext->mPrev = new_node;
713 r_it0.node()->first = r_it1->second + 1;
714 ++r_it1;
715 }
716 else
717 {
718 while( r_it0->second < r_it1->first && r_it0 != lhs.end() )
719 ++r_it0;
720 if( r_it0 == lhs.end() ) break;
721 while( r_it1->second < r_it0->first && r_it1 != range2.end() )
722 ++r_it1;
723 }
724 }
725
726 return lhs;
727 }
728 }
|
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().