DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_STKMeshHelpers.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 "DTK_STKMeshHelpers.hpp"
42 #include "DTK_DBC.hpp"
43 #include "DTK_STKMeshEntityExtraData.hpp"
44 
45 #include <stk_mesh/base/CoordinateSystems.hpp>
46 #include <stk_mesh/base/Field.hpp>
47 #include <stk_mesh/base/FieldBase.hpp>
48 #include <stk_mesh/base/MetaData.hpp>
49 #include <stk_topology/topology.hpp>
50 
51 namespace DataTransferKit
52 {
53 //---------------------------------------------------------------------------//
54 // Given a DTK entity, extract the STK entity.
55 const stk::mesh::Entity &
57 {
58  return Teuchos::rcp_dynamic_cast<STKMeshEntityExtraData>(
59  dtk_entity.extraData() )
60  ->d_stk_entity;
61 }
62 
63 //---------------------------------------------------------------------------//
64 // Given a topological dimension, get the STK entity rank.
65 stk::mesh::EntityRank
67  const int space_dim )
68 {
69  stk::mesh::EntityRank stk_rank = stk::topology::INVALID_RANK;
70 
71  switch ( space_dim )
72  {
73  case 3:
74  switch ( topo_dim )
75  {
76  case 0:
77  stk_rank = stk::topology::NODE_RANK;
78  break;
79  case 1:
80  stk_rank = stk::topology::EDGE_RANK;
81  break;
82  case 2:
83  stk_rank = stk::topology::FACE_RANK;
84  break;
85  case 3:
86  stk_rank = stk::topology::ELEM_RANK;
87  break;
88  default:
89  DTK_CHECK( 0 == topo_dim || 1 == topo_dim || 2 == topo_dim ||
90  3 == topo_dim );
91  break;
92  }
93  break;
94 
95  case 2:
96  switch ( topo_dim )
97  {
98  case 0:
99  stk_rank = stk::topology::NODE_RANK;
100  break;
101  case 1:
102  stk_rank = stk::topology::EDGE_RANK;
103  break;
104  case 2:
105  stk_rank = stk::topology::ELEM_RANK;
106  break;
107  default:
108  DTK_CHECK( 0 == topo_dim || 1 == topo_dim || 2 == topo_dim );
109  break;
110  }
111  break;
112 
113  default:
114  stk_rank = stk::topology::NODE_RANK;
115  }
116 
117  return stk_rank;
118 }
119 
120 //---------------------------------------------------------------------------//
121 // Given a STK entity stk_rank, get the topological dimension.
123  const stk::mesh::EntityRank stk_rank, const int space_dim )
124 {
125  int topo_dim = 0;
126 
127  switch ( space_dim )
128  {
129  case 3:
130  switch ( stk_rank )
131  {
132  case stk::topology::NODE_RANK:
133  topo_dim = 0;
134  break;
135  case stk::topology::EDGE_RANK:
136  topo_dim = 1;
137  break;
138  case stk::topology::FACE_RANK:
139  topo_dim = 2;
140  break;
141  case stk::topology::ELEM_RANK:
142  topo_dim = 3;
143  break;
144  default:
145  DTK_CHECK( stk::topology::NODE_RANK == stk_rank ||
146  stk::topology::EDGE_RANK == stk_rank ||
147  stk::topology::FACE_RANK == stk_rank ||
148  stk::topology::ELEM_RANK == stk_rank );
149  break;
150  }
151  break;
152 
153  case 2:
154  switch ( stk_rank )
155  {
156  case stk::topology::NODE_RANK:
157  topo_dim = 0;
158  break;
159  case stk::topology::EDGE_RANK:
160  topo_dim = 1;
161  break;
162  case stk::topology::ELEM_RANK:
163  topo_dim = 2;
164  break;
165  default:
166  DTK_CHECK( stk::topology::NODE_RANK == stk_rank ||
167  stk::topology::EDGE_RANK == stk_rank ||
168  stk::topology::ELEM_RANK == stk_rank );
169  break;
170  }
171  break;
172 
173  default:
174  DTK_CHECK( 3 == space_dim || 2 == space_dim );
175  break;
176  }
177 
178  return topo_dim;
179 }
180 
181 //---------------------------------------------------------------------------//
182 // Given a DTK entity, return the corresponding STK entity key.
183 stk::mesh::EntityKey STKMeshHelpers::getKeyFromEntity( const Entity dtk_entity )
184 {
185  return stk::mesh::EntityKey(
187  dtk_entity.physicalDimension() ),
188  dtk_entity.id() );
189 }
190 
191 //---------------------------------------------------------------------------//
192 // Given a STK entity, return its shards topology.
193 shards::CellTopology
194 STKMeshHelpers::getShardsTopology( const stk::mesh::Entity stk_entity,
195  const stk::mesh::BulkData &bulk_data )
196 {
197  return stk::mesh::get_cell_topology(
198  bulk_data.bucket( stk_entity ).topology() );
199 }
200 
201 //---------------------------------------------------------------------------//
202 // Given a set of STK entities, return the coordinates of its nodes in a field
203 // container ordered by canonical node order (N,D).
204 Intrepid::FieldContainer<double> STKMeshHelpers::getEntityNodeCoordinates(
205  const Teuchos::Array<stk::mesh::Entity> &stk_entities,
206  const stk::mesh::BulkData &bulk_data )
207 {
208  int space_dim = bulk_data.mesh_meta_data().spatial_dimension();
209  switch ( space_dim )
210  {
211  case 3:
212  return STKMeshHelpers::extractEntityNodeCoordinates<
213  stk::mesh::Cartesian3d>( stk_entities, bulk_data, space_dim );
214  break;
215  case 2:
216  return STKMeshHelpers::extractEntityNodeCoordinates<
217  stk::mesh::Cartesian2d>( stk_entities, bulk_data, space_dim );
218  break;
219  default:
220  DTK_CHECK( 2 == space_dim || 3 == space_dim );
221  break;
222  }
223  return Intrepid::FieldContainer<double>();
224 }
225 
226 //---------------------------------------------------------------------------//
227 
228 } // end namespace DataTransferKit
229 
230 //---------------------------------------------------------------------------//
231 // end DTK_STKMeshHelpers.cpp
232 //---------------------------------------------------------------------------//
static Intrepid::FieldContainer< double > getEntityNodeCoordinates(const Teuchos::Array< stk::mesh::Entity > &stk_entities, const stk::mesh::BulkData &bulk_data)
Given a STK entity, return the coordinates of its nodes in a field container ordered by canonical nod...
Geometric entity interface definition.
Definition: DTK_Entity.hpp:61
static const stk::mesh::Entity & extractEntity(const Entity dtk_entity)
Given a DTK entity, extract the STK entity.
static int getTopologicalDimensionFromRank(const stk::mesh::EntityRank stk_rank, const int space_dim)
Given a STK entity rank, get the topological dimension.
Assertions and Design-by-Contract for error handling.
Teuchos::RCP< EntityExtraData > extraData() const
Get the extra data on the entity. This is a convenient helper for implementing the other interfaces...
Definition: DTK_Entity.cpp:138
static stk::mesh::EntityKey getKeyFromEntity(const Entity dtk_entity)
Given a DTK entity, return the corresponding STK entity key.
static shards::CellTopology getShardsTopology(const stk::mesh::Entity stk_entity, const stk::mesh::BulkData &bulk_data)
Given a STK entity, return its shards topology.
int physicalDimension() const
Return the physical dimension of the entity.
Definition: DTK_Entity.cpp:106
int topologicalDimension() const
Return the topological dimension of the entity.
Definition: DTK_Entity.cpp:98
A base class for setting extra data with entities.
static stk::mesh::EntityRank getRankFromTopologicalDimension(const int topo_dim, const int space_dim)
Given a topological dimension, get the STK entity rank.
DTK_BasicEntitySet.cpp.
EntityId id() const
Client interface.
Definition: DTK_Entity.cpp:82