DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_STKMeshNodalShapeFunction.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_STKMeshNodalShapeFunction.hpp"
42 #include "DTK_DBC.hpp"
43 #include "DTK_STKMeshHelpers.hpp"
44 
45 #include <Shards_CellTopology.hpp>
46 
47 namespace DataTransferKit
48 {
49 //---------------------------------------------------------------------------//
50 // Constructor.
52  const Teuchos::RCP<stk::mesh::BulkData> &bulk_data )
53  : d_bulk_data( bulk_data )
54 { /* ... */
55 }
56 
57 //---------------------------------------------------------------------------//
58 // Given an entity, get the ids of the degrees of freedom in the vector space
59 // supporting its shape function.
61  const Entity &entity, Teuchos::Array<SupportId> &support_ids ) const
62 {
63  // Extract the stk entity.
64  const stk::mesh::Entity &stk_entity =
66  const stk::mesh::EntityRank entity_rank =
67  d_bulk_data->entity_rank( stk_entity );
68 
69  // If the entity is a node, return the id of the node as the support id.
70  if ( stk::topology::NODE_RANK == entity_rank )
71  {
72  support_ids.assign( 1, d_bulk_data->identifier( stk_entity ) );
73  }
74 
75  // Otherwise get the ids of the nodes supporting the entity.
76  else
77  {
78  const stk::mesh::Entity *begin = d_bulk_data->begin_nodes( stk_entity );
79  const stk::mesh::Entity *end = d_bulk_data->end_nodes( stk_entity );
80 
81  // Extract the node ids as the support ids.
82  int num_nodes = std::distance( begin, end );
83  support_ids.resize( num_nodes );
84  for ( int n = 0; n < num_nodes; ++n )
85  {
86  support_ids[n] = d_bulk_data->identifier( begin[n] );
87  }
88  }
89 }
90 
91 //---------------------------------------------------------------------------//
92 // Given an entity and a reference point, evaluate the shape function of the
93 // entity at that point.
95  const Entity &entity,
96  const Teuchos::ArrayView<const double> &reference_point,
97  Teuchos::Array<double> &values ) const
98 {
99  const stk::mesh::Entity &stk_entity =
101 
102  shards::CellTopology entity_topo = stk::mesh::get_cell_topology(
103  d_bulk_data->bucket( stk_entity ).topology() );
104 
105  d_intrepid_shape.evaluateValue( entity_topo, reference_point, values );
106 }
107 
108 //---------------------------------------------------------------------------//
109 // Given an entity and a reference point, evaluate the gradient of the shape
110 // function of the entity at that point.
112  const Entity &entity,
113  const Teuchos::ArrayView<const double> &reference_point,
114  Teuchos::Array<Teuchos::Array<double>> &gradients ) const
115 {
116  const stk::mesh::Entity &stk_entity =
118 
119  shards::CellTopology entity_topo = stk::mesh::get_cell_topology(
120  d_bulk_data->bucket( stk_entity ).topology() );
121 
122  d_intrepid_shape.evaluateGradient( entity_topo, reference_point,
123  gradients );
124 }
125 
126 //---------------------------------------------------------------------------//
127 
128 } // end namespace DataTransferKit
129 
130 //---------------------------------------------------------------------------//
131 // end DTK_STKMeshNodalShapeFunction.cpp
132 //---------------------------------------------------------------------------//
Geometric entity interface definition.
Definition: DTK_Entity.hpp:61
void evaluateGradient(const Entity &entity, const Teuchos::ArrayView< const double > &reference_point, Teuchos::Array< Teuchos::Array< double >> &gradients) const override
Given an entity and a reference point, evaluate the gradient of the shape function of the entity at t...
void evaluateValue(const Entity &entity, const Teuchos::ArrayView< const double > &reference_point, Teuchos::Array< double > &values) const override
Given an entity and a reference point, evaluate the shape function of the entity at that point...
static const stk::mesh::Entity & extractEntity(const Entity dtk_entity)
Given a DTK entity, extract the STK entity.
Assertions and Design-by-Contract for error handling.
STKMeshNodalShapeFunction(const Teuchos::RCP< stk::mesh::BulkData > &bulk_data)
Constructor.
void entitySupportIds(const Entity &entity, Teuchos::Array< SupportId > &support_ids) const override
Given an entity, get the ids of the support locations.
void evaluateGradient(const shards::CellTopology &topology, const Teuchos::ArrayView< const double > &reference_point, Teuchos::Array< Teuchos::Array< double >> &gradients) const
Given an topology and a reference point, evaluate the gradient of the shape function of the topology ...
void evaluateValue(const shards::CellTopology &topology, const Teuchos::ArrayView< const double > &reference_point, Teuchos::Array< double > &values) const
Given an topology and a reference point, evaluate the shape function of the topology at that point...
DTK_BasicEntitySet.cpp.