DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_IntrepidShapeFunction.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_IntrepidShapeFunction.hpp"
42 #include "DTK_DBC.hpp"
44 
45 namespace DataTransferKit
46 {
47 //---------------------------------------------------------------------------//
48 // Given an topology and a reference point, evaluate the shape function of the
49 // topology at that point.
51  const shards::CellTopology &topology,
52  const Teuchos::ArrayView<const double> &reference_point,
53  Teuchos::Array<double> &values ) const
54 {
55  // Get the basis for the topology.
56  Teuchos::RCP<Intrepid::Basis<double, Intrepid::FieldContainer<double>>>
57  basis = getIntrepidBasis( topology );
58 
59  // Wrap the reference point.
60  Teuchos::Array<int> point_dims( 2 );
61  point_dims[0] = 1;
62  point_dims[1] = reference_point.size();
63  Intrepid::FieldContainer<double> point_container(
64  point_dims, const_cast<double *>( reference_point.getRawPtr() ) );
65 
66  // Wrap the evaluations.
67  values.resize( basis->getCardinality() );
68  Teuchos::Array<int> value_dims( 2 );
69  value_dims[0] = basis->getCardinality();
70  value_dims[1] = 1;
71  Intrepid::FieldContainer<double> value_container( value_dims,
72  values.getRawPtr() );
73 
74  // Evaluate the basis function.
75  basis->getValues( value_container, point_container,
76  Intrepid::OPERATOR_VALUE );
77 }
78 
79 //---------------------------------------------------------------------------//
80 // Given an topology and a reference point, evaluate the gradient of the shape
81 // function of the topology at that point.
83  const shards::CellTopology &topology,
84  const Teuchos::ArrayView<const double> &reference_point,
85  Teuchos::Array<Teuchos::Array<double>> &gradients ) const
86 {
87  // Get the basis for the topology.
88  Teuchos::RCP<Intrepid::Basis<double, Intrepid::FieldContainer<double>>>
89  basis = getIntrepidBasis( topology );
90 
91  // Wrap the reference point.
92  int space_dim = reference_point.size();
93  Teuchos::Array<int> point_dims( 2 );
94  point_dims[0] = 1;
95  point_dims[1] = space_dim;
96  Intrepid::FieldContainer<double> point_container(
97  point_dims, const_cast<double *>( reference_point.getRawPtr() ) );
98 
99  // Evaluate the basis function.
100  int cardinality = basis->getCardinality();
101  Intrepid::FieldContainer<double> grad_container( cardinality, 1,
102  space_dim );
103  basis->getValues( grad_container, point_container,
104  Intrepid::OPERATOR_GRAD );
105 
106  // Extract the evaluations.
107  gradients.resize( cardinality );
108  for ( int n = 0; n < cardinality; ++n )
109  {
110  gradients[n].resize( space_dim );
111  for ( int d = 0; d < space_dim; ++d )
112  {
113  gradients[n][d] = grad_container( n, 0, d );
114  }
115  }
116 }
117 
118 //---------------------------------------------------------------------------//
119 // Given a topology, get the intrepid basis function.
120 Teuchos::RCP<Intrepid::Basis<double, Intrepid::FieldContainer<double>>>
121 IntrepidShapeFunction::getIntrepidBasis(
122  const shards::CellTopology &topology ) const
123 {
124  // Either make a new basis for this topology or return an existing one.
125  unsigned basis_key = topology.getKey();
126  Teuchos::RCP<Intrepid::Basis<double, Intrepid::FieldContainer<double>>>
127  basis;
128  if ( d_basis.count( basis_key ) )
129  {
130  basis = d_basis.find( basis_key )->second;
131  }
132  else
133  {
134  basis = IntrepidBasisFactory::create( topology );
135  d_basis.emplace( basis_key, basis );
136  }
137  return basis;
138 }
139 
140 //---------------------------------------------------------------------------//
141 
142 } // end namespace DataTransferKit
143 
144 //---------------------------------------------------------------------------//
145 // end DTK_IntrepidShapeFunction.cpp
146 //---------------------------------------------------------------------------//
Assertions and Design-by-Contract for error handling.
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 ...
Factory for Intrepid basis functions.
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.