DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_STKMeshEntitySet.cpp
1 //---------------------------------------------------------------------------//
2 /*
3  Copyright (c) 2012, Stuart R. Slattery
4  All rights reserved.
5 
6  Redistribution and use in source and binary forms, with or without
7  modification, are permitted provided that the following conditions are
8  met:
9 
10  *: Redistributions of source code must retain the above copyright
11  notice, this list of conditions and the following disclaimer.
12 
13  *: Redistributions in binary form must reproduce the above copyright
14  notice, this list of conditions and the following disclaimer in the
15  documentation and/or other materials provided with the distribution.
16 
17  *: Neither the name of the University of Wisconsin - Madison nor the
18  names of its contributors may be used to endorse or promote products
19  derived from this software without specific prior written permission.
20 
21  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
22  "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
23  LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
24  A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
25  HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
26  SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
27  LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
28  DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
29  THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
30  (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
31  OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
32 */
33 //---------------------------------------------------------------------------//
39 //---------------------------------------------------------------------------//
40 
41 #include <vector>
42 
43 #include "DTK_DBC.hpp"
44 #include "DTK_STKMeshEntity.hpp"
45 #include "DTK_STKMeshEntityIterator.hpp"
46 #include "DTK_STKMeshEntityIteratorRange.hpp"
47 #include "DTK_STKMeshEntitySet.hpp"
48 #include "DTK_STKMeshHelpers.hpp"
49 
50 #include <stk_mesh/base/GetEntities.hpp>
51 #include <stk_mesh/base/MetaData.hpp>
52 #include <stk_topology/topology.hpp>
53 
54 #include <Teuchos_DefaultMpiComm.hpp>
55 
56 namespace DataTransferKit
57 {
58 //---------------------------------------------------------------------------//
59 // Constructor.
61  const Teuchos::RCP<stk::mesh::BulkData> &bulk_data )
62  : d_bulk_data( bulk_data )
63 {
64  DTK_REQUIRE( Teuchos::nonnull( d_bulk_data ) );
65 }
66 
67 //---------------------------------------------------------------------------//
68 // Get the parallel communicator for the entity set.
69 Teuchos::RCP<const Teuchos::Comm<int>> STKMeshEntitySet::communicator() const
70 {
71  return Teuchos::rcp( new Teuchos::MpiComm<int>( d_bulk_data->parallel() ) );
72 }
73 
74 //---------------------------------------------------------------------------//
75 // Return the largest physical dimension of the entities in the set.
77 {
78  return d_bulk_data->mesh_meta_data().spatial_dimension();
79 }
80 
81 //---------------------------------------------------------------------------//
82 // Given an EntityId, get the entity.
83 void STKMeshEntitySet::getEntity( const EntityId entity_id,
84  const int topological_dimension,
85  Entity &entity ) const
86 {
87  stk::mesh::Entity stk_entity = d_bulk_data->get_entity(
90  entity_id );
91  entity = STKMeshEntity( stk_entity, d_bulk_data.ptr() );
92 }
93 
94 //---------------------------------------------------------------------------//
95 // Get an iterator over a subset of the entity set that satisfies the given
96 // predicate.
98 STKMeshEntitySet::entityIterator( const int topological_dimension,
99  const PredicateFunction &predicate ) const
100 {
101  stk::mesh::EntityRank rank =
103  physicalDimension() );
104  Teuchos::RCP<STKMeshEntityIteratorRange> iterator_range =
105  Teuchos::rcp( new STKMeshEntityIteratorRange() );
106  stk::mesh::get_entities( *d_bulk_data, rank,
107  iterator_range->d_stk_entities );
108  return STKMeshEntityIterator( iterator_range, d_bulk_data.ptr(),
109  predicate );
110 }
111 
112 //---------------------------------------------------------------------------//
113 // Given an entity, get the entities of the given type that are adjacent to
114 // it.
116  const Entity &entity, const int adjacent_dimension,
117  Teuchos::Array<Entity> &adjacent_entities ) const
118 {
119  const stk::mesh::Entity &stk_entity =
121  stk::mesh::EntityRank rank =
123  physicalDimension() );
124  const stk::mesh::Entity *begin = d_bulk_data->begin( stk_entity, rank );
125  const stk::mesh::Entity *end = d_bulk_data->end( stk_entity, rank );
126  Teuchos::Array<stk::mesh::Entity> stk_adjacencies( begin, end );
127  adjacent_entities.resize( stk_adjacencies.size() );
128  Teuchos::Array<Entity>::iterator entity_it;
129  Teuchos::Array<stk::mesh::Entity>::iterator stk_it;
130  for ( entity_it = adjacent_entities.begin(),
131  stk_it = stk_adjacencies.begin();
132  entity_it != adjacent_entities.end(); ++entity_it, ++stk_it )
133  {
134  *entity_it = STKMeshEntity( *stk_it, d_bulk_data.ptr() );
135  }
136 }
137 
138 //---------------------------------------------------------------------------//
139 
140 } // end namespace DataTransferKit
141 
142 //---------------------------------------------------------------------------//
143 // end DTK_STKMeshEntitySet.cpp
144 //---------------------------------------------------------------------------//
Geometric entity interface definition.
Definition: DTK_Entity.hpp:61
STKMesh entity interface definition.
A container for entity vectors that can be reference-counted.
Entity iterator interface.
static const stk::mesh::Entity & extractEntity(const Entity dtk_entity)
Given a DTK entity, extract the STK entity.
STK mesh entity iterator implementation.
void getEntity(const EntityId entity_id, const int topological_dimension, Entity &entity) const override
Entity access functions.
STKMeshEntitySet(const Teuchos::RCP< stk::mesh::BulkData > &bulk_data)
Constructor.
EntityIterator entityIterator(const int topological_dimension, const PredicateFunction &predicate) const override
Get a iterator of the given entity type that satisfy the given predicate.
Assertions and Design-by-Contract for error handling.
void getAdjacentEntities(const Entity &entity, const int adjacent_dimension, Teuchos::Array< Entity > &adjacent_entities) const override
Given an entity, get the entities of the given topological dimension that are adjacent to it...
std::function< bool(Entity)> PredicateFunction
Predicate function typedef.
Definition: DTK_Types.hpp:64
unsigned long int EntityId
Entity id type.
Definition: DTK_Types.hpp:50
static stk::mesh::EntityRank getRankFromTopologicalDimension(const int topo_dim, const int space_dim)
Given a topological dimension, get the STK entity rank.
int physicalDimension() const override
Geometric data functions.
Teuchos::RCP< const Teuchos::Comm< int > > communicator() const override
Parallel functions.
DTK_BasicEntitySet.cpp.