DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_IntrepidShapeFunction.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_INTREPIDSHAPEFUNCTION
42 #define DTK_INTREPIDSHAPEFUNCTION
43 
44 #include <unordered_map>
45 
46 #include <Teuchos_Array.hpp>
47 #include <Teuchos_RCP.hpp>
48 
49 #include <Intrepid_Basis.hpp>
50 #include <Intrepid_FieldContainer.hpp>
51 
52 #include <Shards_CellTopology.hpp>
53 
54 namespace DataTransferKit
55 {
56 //---------------------------------------------------------------------------//
61 //---------------------------------------------------------------------------//
63 {
64  public:
74  void evaluateValue( const shards::CellTopology &topology,
75  const Teuchos::ArrayView<const double> &reference_point,
76  Teuchos::Array<double> &values ) const;
77 
90  void
91  evaluateGradient( const shards::CellTopology &topology,
92  const Teuchos::ArrayView<const double> &reference_point,
93  Teuchos::Array<Teuchos::Array<double>> &gradients ) const;
94 
95  private:
96  // Get the basis of a topology.
97  Teuchos::RCP<Intrepid::Basis<double, Intrepid::FieldContainer<double>>>
98  getIntrepidBasis( const shards::CellTopology &topology ) const;
99 
100  private:
101  // Map of already created shape functions.
102  mutable std::unordered_map<
103  unsigned,
104  Teuchos::RCP<Intrepid::Basis<double, Intrepid::FieldContainer<double>>>>
105  d_basis;
106 };
107 
108 //---------------------------------------------------------------------------//
109 
110 } // end namespace DataTransferKit
111 
112 //---------------------------------------------------------------------------//
113 
114 #endif // end DTK_INTREPIDSHAPEFUNCTION
115 
116 //---------------------------------------------------------------------------//
117 // end DTK_IntrepidShapeFunction.hpp
118 //---------------------------------------------------------------------------//
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.