Intrepid2
Intrepid2_HGRAD_WEDGE_C2_FEM.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Intrepid2 Package
4//
5// Copyright 2007 NTESS and the Intrepid2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
16#ifndef __INTREPID2_HGRAD_WEDGE_C2_FEM_HPP__
17#define __INTREPID2_HGRAD_WEDGE_C2_FEM_HPP__
18
19#include "Intrepid2_Basis.hpp"
22
23namespace Intrepid2 {
24
97 namespace Impl {
98
102 template<bool serendipity>
104 public:
105 typedef struct Wedge<serendipity ? 15 : 18> cell_topology_type;
109 template<EOperator opType>
110 struct Serial {
111 template<typename OutputViewType,
112 typename inputViewType>
113 KOKKOS_INLINE_FUNCTION
114 static void
115 getValues( OutputViewType output,
116 const inputViewType input );
117
118 };
119
120 template<typename DeviceType,
121 typename outputValueValueType, class ...outputValueProperties,
122 typename inputPointValueType, class ...inputPointProperties>
123 static void
124 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
125 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
126 const EOperator operatorType);
127
131 template<typename outputValueViewType,
132 typename inputPointViewType,
133 EOperator opType>
134 struct Functor {
135 outputValueViewType _outputValues;
136 const inputPointViewType _inputPoints;
137
138 KOKKOS_INLINE_FUNCTION
139 Functor( outputValueViewType outputValues_,
140 inputPointViewType inputPoints_ )
141 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
142 KOKKOS_INLINE_FUNCTION
143 void operator()(const ordinal_type pt) const {
144 switch (opType) {
145 case OPERATOR_VALUE : {
146 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
147 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
148 Serial<opType>::getValues( output, input );
149 break;
150 }
151 case OPERATOR_GRAD :
152 case OPERATOR_D2 :
153 case OPERATOR_D3 :
154 case OPERATOR_D4 :
155 case OPERATOR_MAX : {
156 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
157 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
158 Serial<opType>::getValues( output, input );
159 break;
160 }
161 default: {
162 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
163 opType != OPERATOR_GRAD &&
164 opType != OPERATOR_D2 &&
165 opType != OPERATOR_MAX,
166 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::Serial::getValues) operator is not supported");
167 }
168 }
169 }
170 };
171 };
172 }
173
174 template<bool serendipity,
175 typename DeviceType = void,
176 typename outputValueType = double,
177 typename pointValueType = double>
178 class Basis_HGRAD_WEDGE_DEG2_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
179 public:
186
190
191 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
192
193 virtual
194 void
195 getValues( OutputViewType outputValues,
196 const PointViewType inputPoints,
197 const EOperator operatorType = OPERATOR_VALUE ) const override {
198#ifdef HAVE_INTREPID2_DEBUG
199 // Verify arguments
201 inputPoints,
202 operatorType,
203 this->getBaseCellTopology(),
204 this->getCardinality() );
205#endif
206 if constexpr (serendipity)
208 getValues<DeviceType>( outputValues,
209 inputPoints,
210 operatorType );
211 else
213 getValues<DeviceType>( outputValues,
214 inputPoints,
215 operatorType );;
216 }
217
218 virtual void
219 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
220 ordinal_type& perThreadSpaceSize,
221 const PointViewType inputPointsconst,
222 const EOperator operatorType = OPERATOR_VALUE) const override;
223
224 KOKKOS_INLINE_FUNCTION
225 virtual void
226 getValues(
227 OutputViewType outputValues,
228 const PointViewType inputPoints,
229 const EOperator operatorType,
230 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
231 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
232 const ordinal_type subcellDim = -1,
233 const ordinal_type subcellOrdinal = -1) const override;
234
235 virtual
236 void
237 getDofCoords( ScalarViewType dofCoords ) const override {
238#ifdef HAVE_INTREPID2_DEBUG
239 // Verify rank of output array.
240 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
241 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getDofCoords) rank = 2 required for dofCoords array");
242 // Verify 0th dimension of output array.
243 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
244 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
245 // Verify 1st dimension of output array.
246 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
247 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
248#endif
249 Kokkos::deep_copy(dofCoords, this->dofCoords_);
250 }
251
252 virtual
253 void
254 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
255#ifdef HAVE_INTREPID2_DEBUG
256 // Verify rank of output array.
257 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
258 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
259 // Verify 0th dimension of output array.
260 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
261 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_DEG2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
262#endif
263 Kokkos::deep_copy(dofCoeffs, 1.0);
264 }
265
266 virtual
267 const char*
268 getName() const override {
269 return serendipity ? "Intrepid2_HGRAD_WEDGE_I2_FEM" : "Intrepid2_HGRAD_WEDGE_C2_FEM";
270 }
271
281 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
282 if(subCellDim == 1) {
283 return Teuchos::rcp(new
285 } else if(subCellDim == 2) {
286 if(subCellOrd < 3) //lateral faces
288 else
290 }
291 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
292 }
293
298 };
299
300 template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
301 using Basis_HGRAD_WEDGE_C2_FEM = Basis_HGRAD_WEDGE_DEG2_FEM<false, DeviceType, outputValueType, pointValueType>;
302
303 template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
304 using Basis_HGRAD_WEDGE_I2_FEM = Basis_HGRAD_WEDGE_DEG2_FEM<true, DeviceType, outputValueType, pointValueType>;
305
306}// namespace Intrepid2
307
309
310#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....
Header file for the Intrepid2::Basis_HGRAD_QUAD_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TRI_C2_FEM class.
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 Line cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Triangle cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Wedge cell.
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...
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual const char * getName() const override
Returns basis name.
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_DEG2_FEM.