Intrepid2
Intrepid2_HGRAD_HEX_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_HEX_C2_FEM_HPP__
17#define __INTREPID2_HGRAD_HEX_C2_FEM_HPP__
18
19#include "Intrepid2_Basis.hpp"
21
22namespace Intrepid2 {
23
114 namespace Impl {
115
119 template<bool serendipity>
121 public:
122 typedef struct Hexahedron<serendipity ? 20 : 27> cell_topology_type;
126 template<EOperator opType>
127 struct Serial {
128 template<typename OutputViewType,
129 typename inputViewType>
130 KOKKOS_INLINE_FUNCTION
131 static void
132 getValues( OutputViewType output,
133 const inputViewType input );
134
135 };
136
137 template<typename DeviceType,
138 typename outputValueValueType, class ...outputValueProperties,
139 typename inputPointValueType, class ...inputPointProperties>
140 static void
141 getValues( const typename DeviceType::execution_space& space,
142 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
143 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
144 const EOperator operatorType);
145
149 template<typename outputValueViewType,
150 typename inputPointViewType,
151 EOperator opType>
152 struct Functor {
153 outputValueViewType _outputValues;
154 const inputPointViewType _inputPoints;
155
156 KOKKOS_INLINE_FUNCTION
157 Functor( outputValueViewType outputValues_,
158 inputPointViewType inputPoints_ )
159 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
160
161 KOKKOS_INLINE_FUNCTION
162 void operator()(const ordinal_type pt) const {
163 switch (opType) {
164 case OPERATOR_VALUE : {
165 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
166 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
167 Serial<opType>::getValues( output, input );
168 break;
169 }
170 case OPERATOR_GRAD :
171 case OPERATOR_D1 :
172 case OPERATOR_D2 :
173 case OPERATOR_D3 :
174 case OPERATOR_D4 :
175 case OPERATOR_MAX : {
176 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
177 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
178 Serial<opType>::getValues( output, input );
179 break;
180 }
181 default: {
182 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
183 opType != OPERATOR_GRAD &&
184 opType != OPERATOR_D1 &&
185 opType != OPERATOR_D2 &&
186 opType != OPERATOR_D3 &&
187 opType != OPERATOR_D4 &&
188 opType != OPERATOR_MAX,
189 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::Serial::getValues) operator is not supported");
190 }
191 }
192 }
193 };
194 };
195 }
196
197 template<bool serendipity,
198 typename DeviceType,
199 typename outputValueType,
200 typename pointValueType>
201 class Basis_HGRAD_HEX_DEG2_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
202 public:
204 using typename BasisBase::ExecutionSpace;
205
209
210 using typename BasisBase::OutputViewType;
211 using typename BasisBase::PointViewType ;
212 using typename BasisBase::ScalarViewType;
213
217
219
220 virtual
221 void
223 OutputViewType outputValues,
224 const PointViewType inputPoints,
225 const EOperator operatorType = OPERATOR_VALUE ) const override {
226#ifdef HAVE_INTREPID2_DEBUG
227 // Verify arguments
229 inputPoints,
230 operatorType,
231 this->getBaseCellTopology(),
232 this->getCardinality());
233#endif
234 if constexpr (serendipity)
236 getValues<DeviceType>(space,
237 outputValues,
238 inputPoints,
239 operatorType);
240 else
243 outputValues,
244 inputPoints,
245 operatorType);
246 }
247
248 virtual void
249 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
250 ordinal_type& perThreadSpaceSize,
251 const PointViewType inputPointsconst,
252 const EOperator operatorType = OPERATOR_VALUE) const override;
253
254 KOKKOS_INLINE_FUNCTION
255 virtual void
256 getValues(
257 OutputViewType outputValues,
258 const PointViewType inputPoints,
259 const EOperator operatorType,
260 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
261 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
262 const ordinal_type subcellDim = -1,
263 const ordinal_type subcellOrdinal = -1) const override;
264
265 virtual
266 void
267 getDofCoords( ScalarViewType dofCoords ) const override {
268#ifdef HAVE_INTREPID2_DEBUG
269 // Verify rank of output array.
270 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
271 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getDofCoords) rank = 2 required for dofCoords array");
272 // Verify 0th dimension of output array.
273 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
274 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
275 // Verify 1st dimension of output array.
276 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
277 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
278#endif
279 Kokkos::deep_copy(dofCoords, this->dofCoords_);
280 }
281
282 virtual
283 void
284 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
285#ifdef HAVE_INTREPID2_DEBUG
286 // Verify rank of output array.
287 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
288 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
289 // Verify 0th dimension of output array.
290 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
291 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
292#endif
293 Kokkos::deep_copy(dofCoeffs, 1.0);
294 }
295
296 virtual
297 const char*
298 getName() const override {
299 return serendipity ? "Intrepid2_HGRAD_HEX_I2_FEM" : "Intrepid2_HGRAD_HEX_C2_FEM";
300 }
301
311 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
312 if(subCellDim == 1) {
313 return Teuchos::rcp(new
315 } else if(subCellDim == 2) {
316 return Teuchos::rcp(new
318 }
319 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
320 }
321
326 };
327
328 template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
329 using Basis_HGRAD_HEX_C2_FEM = Basis_HGRAD_HEX_DEG2_FEM<false, DeviceType, outputValueType, pointValueType>;
330
331 template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
332 using Basis_HGRAD_HEX_I2_FEM = Basis_HGRAD_HEX_DEG2_FEM<true, DeviceType, outputValueType, pointValueType>;
333
334}// namespace Intrepid2
335
337
338#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 HEX cells.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C2_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron cell.
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPointsconst, const EOperator operatorType=OPERATOR_VALUE) const override
Return the size of the scratch space, in bytes, needed for using the team-level implementation of get...
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 void getValues(const ExecutionSpace &space, OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
virtual const char * getName() const override
Returns basis name.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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.
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.
virtual KOKKOS_INLINE_FUNCTION void getValues(OutputViewType, const PointViewType, const EOperator, const typename Kokkos::TeamPolicy< ExecutionSpace >::member_type &teamMember, const typename ExecutionSpace::scratch_memory_space &scratchStorage, const ordinal_type subcellDim=-1, const ordinal_type subcellOrdinal=-1) const
Team-level evaluation of basis functions on a reference cell.
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.
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.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
See Intrepid2::Basis_HGRAD_HEX_DEG2_FEM.