Intrepid2
Intrepid2_HGRAD_HEX_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_HEX_C2_FEM_HPP__
50#define __INTREPID2_HGRAD_HEX_C2_FEM_HPP__
51
52#include "Intrepid2_Basis.hpp"
54
55namespace Intrepid2 {
56
146
147 namespace Impl {
148
152 template<bool serendipity>
154 public:
155 typedef struct Hexahedron<serendipity ? 20 : 27> cell_topology_type;
159 template<EOperator opType>
160 struct Serial {
161 template<typename OutputViewType,
162 typename inputViewType>
163 KOKKOS_INLINE_FUNCTION
164 static void
165 getValues( OutputViewType output,
166 const inputViewType input );
167
168 };
169
170 template<typename DeviceType,
171 typename outputValueValueType, class ...outputValueProperties,
172 typename inputPointValueType, class ...inputPointProperties>
173 static void
174 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
175 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
176 const EOperator operatorType);
177
181 template<typename outputValueViewType,
182 typename inputPointViewType,
183 EOperator opType>
184 struct Functor {
185 outputValueViewType _outputValues;
186 const inputPointViewType _inputPoints;
187
188 KOKKOS_INLINE_FUNCTION
189 Functor( outputValueViewType outputValues_,
190 inputPointViewType inputPoints_ )
191 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
192
193 KOKKOS_INLINE_FUNCTION
194 void operator()(const ordinal_type pt) const {
195 switch (opType) {
196 case OPERATOR_VALUE : {
197 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
198 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
199 Serial<opType>::getValues( output, input );
200 break;
201 }
202 case OPERATOR_GRAD :
203 case OPERATOR_D1 :
204 case OPERATOR_D2 :
205 case OPERATOR_D3 :
206 case OPERATOR_D4 :
207 case OPERATOR_MAX : {
208 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
209 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
210 Serial<opType>::getValues( output, input );
211 break;
212 }
213 default: {
214 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
215 opType != OPERATOR_GRAD &&
216 opType != OPERATOR_D1 &&
217 opType != OPERATOR_D2 &&
218 opType != OPERATOR_D3 &&
219 opType != OPERATOR_D4 &&
220 opType != OPERATOR_MAX,
221 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::Serial::getValues) operator is not supported");
222 }
223 }
224 }
225 };
226 };
227 }
228
229 template<bool serendipity,
230 typename DeviceType,
231 typename outputValueType,
232 typename pointValueType>
233 class Basis_HGRAD_HEX_DEG2_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
234 public:
238
242
246
247 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
248
249 virtual
250 void
251 getValues( OutputViewType outputValues,
252 const PointViewType inputPoints,
253 const EOperator operatorType = OPERATOR_VALUE ) const override {
254#ifdef HAVE_INTREPID2_DEBUG
255 // Verify arguments
257 inputPoints,
258 operatorType,
259 this->getBaseCellTopology(),
260 this->getCardinality() );
261#endif
262 if constexpr (serendipity)
264 getValues<DeviceType>( outputValues,
265 inputPoints,
266 operatorType );
267 else
269 getValues<DeviceType>( outputValues,
270 inputPoints,
271 operatorType );
272 }
273
274 virtual
275 void
276 getDofCoords( ScalarViewType dofCoords ) const override {
277#ifdef HAVE_INTREPID2_DEBUG
278 // Verify rank of output array.
279 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
280 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getDofCoords) rank = 2 required for dofCoords array");
281 // Verify 0th dimension of output array.
282 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
283 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
284 // Verify 1st dimension of output array.
285 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
286 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
287#endif
288 Kokkos::deep_copy(dofCoords, this->dofCoords_);
289 }
290
291 virtual
292 void
293 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
294#ifdef HAVE_INTREPID2_DEBUG
295 // Verify rank of output array.
296 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
297 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
298 // Verify 0th dimension of output array.
299 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
300 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_DEG2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
301#endif
302 Kokkos::deep_copy(dofCoeffs, 1.0);
303 }
304
305 virtual
306 const char*
307 getName() const override {
308 if constexpr (serendipity)
309 return "Intrepid2_HGRAD_HEX_I2_Serendipity_FEM";
310 else
311 return "Intrepid2_HGRAD_HEX_C2_FEM";
312 }
313
323 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
324 if(subCellDim == 1) {
325 return Teuchos::rcp(new
327 } else if(subCellDim == 2) {
328 return Teuchos::rcp(new
330 }
331 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
332 }
333
338 };
339
340 template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
341 using Basis_HGRAD_HEX_C2_FEM = Basis_HGRAD_HEX_DEG2_FEM<false, DeviceType, outputValueType, pointValueType>;
342
343 template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
344 using Basis_HGRAD_HEX_I2_Serendipity_FEM = Basis_HGRAD_HEX_DEG2_FEM<true, DeviceType, outputValueType, pointValueType>;
345
346}// namespace Intrepid2
347
349
350#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 HEX cells.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C2_FEM class.
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
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.
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.
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.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HGRAD_HEX_DEG2_FEM.