Intrepid2
Intrepid2_HGRAD_TET_C2_FEM.hpp
Go to the documentation of this file.
1// @HEADER
2// ************************************************************************
3//
4// Intrepid2 Package
5// Copyright (2007) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact Kyungjoo Kim (kyukim@sandia.gov), or
38// Mauro Perego (mperego@sandia.gov)
39//
40// ************************************************************************
41// @HEADER
42
48
49#ifndef __INTREPID2_HGRAD_TET_C2_FEM_HPP__
50#define __INTREPID2_HGRAD_TET_C2_FEM_HPP__
51
52#include "Intrepid2_Basis.hpp"
53
54namespace Intrepid2 {
55
99
100 namespace Impl {
101
106 public:
107 typedef struct Tetrahedron<10> cell_topology_type;
111 template<EOperator opType>
112 struct Serial {
113 template<typename OutputViewType,
114 typename inputViewType>
115 KOKKOS_INLINE_FUNCTION
116 static void
117 getValues( OutputViewType output,
118 const inputViewType input );
119
120 };
121
122 template<typename DeviceType,
123 typename outputValueValueType, class ...outputValueProperties,
124 typename inputPointValueType, class ...inputPointProperties>
125 static void
126 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
127 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
128 const EOperator operatorType);
129
133 template<typename outputValueViewType,
134 typename inputPointViewType,
135 EOperator opType>
136 struct Functor {
137 outputValueViewType _outputValues;
138 const inputPointViewType _inputPoints;
139
140 KOKKOS_INLINE_FUNCTION
141 Functor( outputValueViewType outputValues_,
142 inputPointViewType inputPoints_ )
143 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
144
145 KOKKOS_INLINE_FUNCTION
146 void operator()(const ordinal_type pt) const {
147 switch (opType) {
148 case OPERATOR_VALUE : {
149 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
150 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
151 Serial<opType>::getValues( output, input );
152 break;
153 }
154 case OPERATOR_GRAD :
155 case OPERATOR_D1 :
156 case OPERATOR_D2 :
157 case OPERATOR_MAX : {
158 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
159 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
160 Serial<opType>::getValues( output, input );
161 break;
162 }
163 default: {
164 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
165 opType != OPERATOR_GRAD &&
166 opType != OPERATOR_MAX,
167 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::Serial::getValues) operator is not supported");
168 }
169 }
170 }
171 };
172 };
173 }
174
175 template<typename DeviceType = void,
176 typename outputValueType = double,
177 typename pointValueType = double>
178 class Basis_HGRAD_TET_C2_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
179 public:
183
187
191
192 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
193
194 virtual
195 void
196 getValues( OutputViewType outputValues,
197 const PointViewType inputPoints,
198 const EOperator operatorType = OPERATOR_VALUE ) const override {
199#ifdef HAVE_INTREPID2_DEBUG
200 // Verify arguments
202 inputPoints,
203 operatorType,
204 this->getBaseCellTopology(),
205 this->getCardinality() );
206#endif
207 Impl::Basis_HGRAD_TET_C2_FEM::
208 getValues<DeviceType>( outputValues,
209 inputPoints,
210 operatorType );
211 }
212
213 virtual
214 void
215 getDofCoords( ScalarViewType dofCoords ) const override {
216#ifdef HAVE_INTREPID2_DEBUG
217 // Verify rank of output array.
218 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
219 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getDofCoords) rank = 2 required for dofCoords array");
220 // Verify 0th dimension of output array.
221 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
222 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
223 // Verify 1st dimension of output array.
224 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
225 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
226#endif
227 Kokkos::deep_copy(dofCoords, this->dofCoords_);
228 }
229
230
231 virtual
232 void
233 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
234#ifdef HAVE_INTREPID2_DEBUG
235 // Verify rank of output array.
236 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
237 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
238 // Verify 0th dimension of output array.
239 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
240 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
241#endif
242 Kokkos::deep_copy(dofCoeffs, 1.0);
243 }
244
245 virtual
246 const char*
247 getName() const override {
248 return "Intrepid2_HGRAD_TET_C2_FEM";
249 }
250
255 };
256}// namespace Intrepid2
257
259
260#endif
Header file for the abstract base class Intrepid2::Basis.
void getValues_HGRAD_Args(const outputValueViewType outputValues, const inputPointViewType inputPoints, const EOperator operatorType, const shards::CellTopology cellTopo, const ordinal_type basisCard)
Runtime check of the arguments for the getValues method in an HGRAD-conforming FEM basis....
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
Definition file for FEM basis functions of degree 2 for H(grad) functions on TET cells.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const override
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual const char * getName() const override
Returns basis name.
BasisPtr< typename Kokkos::HostSpace::device_type, outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
shards::CellTopology getBaseCellTopology() const
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HGRAD_TET_C2_FEM.