Intrepid2
Intrepid2_HGRAD_WEDGE_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
49#ifndef __INTREPID2_HGRAD_WEDGE_C2_FEM_HPP__
50#define __INTREPID2_HGRAD_WEDGE_C2_FEM_HPP__
51
52#include "Intrepid2_Basis.hpp"
53
54namespace Intrepid2 {
55
120 namespace Impl {
121
126 public:
127 typedef struct Wedge<18> cell_topology_type;
131 template<EOperator opType>
132 struct Serial {
133 template<typename OutputViewType,
134 typename inputViewType>
135 KOKKOS_INLINE_FUNCTION
136 static void
137 getValues( OutputViewType output,
138 const inputViewType input );
139
140 };
141
142 template<typename DeviceType,
143 typename outputValueValueType, class ...outputValueProperties,
144 typename inputPointValueType, class ...inputPointProperties>
145 static void
146 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
147 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
148 const EOperator operatorType);
149
153 template<typename outputValueViewType,
154 typename inputPointViewType,
155 EOperator opType>
156 struct Functor {
157 outputValueViewType _outputValues;
158 const inputPointViewType _inputPoints;
159
160 KOKKOS_INLINE_FUNCTION
161 Functor( outputValueViewType outputValues_,
162 inputPointViewType inputPoints_ )
163 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
164 KOKKOS_INLINE_FUNCTION
165 void operator()(const ordinal_type pt) const {
166 switch (opType) {
167 case OPERATOR_VALUE : {
168 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
169 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
170 Serial<opType>::getValues( output, input );
171 break;
172 }
173 case OPERATOR_GRAD :
174 case OPERATOR_D2 :
175 case OPERATOR_D3 :
176 case OPERATOR_D4 :
177 case OPERATOR_MAX : {
178 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
179 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
180 Serial<opType>::getValues( output, input );
181 break;
182 }
183 default: {
184 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
185 opType != OPERATOR_GRAD &&
186 opType != OPERATOR_D2 &&
187 opType != OPERATOR_MAX,
188 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C2_FEM::Serial::getValues) operator is not supported");
189 }
190 }
191 }
192 };
193 };
194 }
195
196 template<typename DeviceType = void,
197 typename outputValueType = double,
198 typename pointValueType = double>
199 class Basis_HGRAD_WEDGE_C2_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
200 public:
207
211
212 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
213
214 virtual
215 void
216 getValues( OutputViewType outputValues,
217 const PointViewType inputPoints,
218 const EOperator operatorType = OPERATOR_VALUE ) const override {
219#ifdef HAVE_INTREPID2_DEBUG
220 // Verify arguments
222 inputPoints,
223 operatorType,
224 this->getBaseCellTopology(),
225 this->getCardinality() );
226#endif
227 Impl::Basis_HGRAD_WEDGE_C2_FEM::
228 getValues<DeviceType>( outputValues,
229 inputPoints,
230 operatorType );
231 }
232
233 virtual
234 void
235 getDofCoords( ScalarViewType dofCoords ) const override {
236#ifdef HAVE_INTREPID2_DEBUG
237 // Verify rank of output array.
238 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
239 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C2_FEM::getDofCoords) rank = 2 required for dofCoords array");
240 // Verify 0th dimension of output array.
241 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
242 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
243 // Verify 1st dimension of output array.
244 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
245 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
246#endif
247 Kokkos::deep_copy(dofCoords, this->dofCoords_);
248 }
249
250 virtual
251 void
252 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
253#ifdef HAVE_INTREPID2_DEBUG
254 // Verify rank of output array.
255 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
256 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
257 // Verify 0th dimension of output array.
258 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
259 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
260#endif
261 Kokkos::deep_copy(dofCoeffs, 1.0);
262 }
263
264 virtual
265 const char*
266 getName() const override {
267 return "Intrepid2_HGRAD_WEDGE_C2_FEM";
268 }
269
271 getHostBasis() const override{
273 }
274 };
275
276}// namespace Intrepid2
277
279
280#endif
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
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....
Definition file for FEM basis functions of degree 2 for H(grad) functions on WEDGE cells.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Wedge cell.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the 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...
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...
virtual const char * getName() const override
Returns basis name.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
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.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Device DeviceType
(Kokkos) Device type on which Basis is templated. Does not necessarily return true for Kokkos::is_dev...
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HGRAD_WEDGE_C2_FEM.