DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_STKMeshNodalShapeFunction.hpp
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 #ifndef DTK_STKMESHNODALSHAPEFUNCTION
42 #define DTK_STKMESHNODALSHAPEFUNCTION
43 
44 #include "DTK_EntityShapeFunction.hpp"
45 #include "DTK_IntrepidShapeFunction.hpp"
46 #include "DTK_Types.hpp"
47 
48 #include <Teuchos_Array.hpp>
49 #include <Teuchos_RCP.hpp>
50 
51 #include <stk_mesh/base/BulkData.hpp>
52 
53 namespace DataTransferKit
54 {
55 //---------------------------------------------------------------------------//
66 //---------------------------------------------------------------------------//
68 {
69  public:
74  const Teuchos::RCP<stk::mesh::BulkData> &bulk_data );
75 
82  void
83  entitySupportIds( const Entity &entity,
84  Teuchos::Array<SupportId> &support_ids ) const override;
85 
95  void evaluateValue( const Entity &entity,
96  const Teuchos::ArrayView<const double> &reference_point,
97  Teuchos::Array<double> &values ) const override;
98 
111  void evaluateGradient(
112  const Entity &entity,
113  const Teuchos::ArrayView<const double> &reference_point,
114  Teuchos::Array<Teuchos::Array<double>> &gradients ) const override;
115 
116  private:
117  // Bulk data for the mesh over which the shape function is defined.
118  Teuchos::RCP<stk::mesh::BulkData> d_bulk_data;
119 
120  // Intrepid shape function.
121  IntrepidShapeFunction d_intrepid_shape;
122 };
123 
124 //---------------------------------------------------------------------------//
125 
126 } // end namespace DataTransferKit
127 
128 //---------------------------------------------------------------------------//
129 
130 #endif // end DTK_STKMESHNODALSHAPEFUNCTION
131 
132 //---------------------------------------------------------------------------//
133 // end DTK_STKMeshNodalShapeFunction.hpp
134 //---------------------------------------------------------------------------//
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...
STKMeshNodalShapeFunction(const Teuchos::RCP< stk::mesh::BulkData > &bulk_data)
Constructor.
Nodal shape function implementation for STK mesh.
void entitySupportIds(const Entity &entity, Teuchos::Array< SupportId > &support_ids) const override
Given an entity, get the ids of the support locations.
DTK_BasicEntitySet.cpp.