Intrepid2
Intrepid2_HVOL_HEX_Cn_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
15#ifndef __INTREPID2_HVOL_HEX_CN_FEM_HPP__
16#define __INTREPID2_HVOL_HEX_CN_FEM_HPP__
17
18#include "Intrepid2_Basis.hpp"
20
21namespace Intrepid2 {
22
23 namespace Impl {
28 public:
29 typedef struct Hexahedron<8> cell_topology_type;
33 template<EOperator opType>
34 struct Serial {
35 template<typename outputValueViewType,
36 typename inputPointViewType,
37 typename workViewType,
38 typename vinvViewType>
39 KOKKOS_INLINE_FUNCTION
40 static void
41 getValues( outputValueViewType outputValues,
42 const inputPointViewType inputPoints,
43 workViewType work,
44 const vinvViewType vinv,
45 const ordinal_type operatorDn = 0 );
46
47 KOKKOS_INLINE_FUNCTION
48 static ordinal_type
49 getWorkSizePerPoint(ordinal_type order) {return 4*getPnCardinality<1>(order); }
50 };
51
52 template<typename DeviceType, ordinal_type numPtsPerEval,
53 typename outputValueValueType, class ...outputValueProperties,
54 typename inputPointValueType, class ...inputPointProperties,
55 typename vinvValueType, class ...vinvProperties>
56 static void
57 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
58 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
59 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
60 const EOperator operatorType );
61
65 template<typename outputValueViewType,
66 typename inputPointViewType,
67 typename vinvViewType,
68 typename workViewType,
69 EOperator opType,
70 ordinal_type numPtsEval>
71 struct Functor {
72 outputValueViewType _outputValues;
73 const inputPointViewType _inputPoints;
74 const vinvViewType _vinv;
75 workViewType _work;
76 const ordinal_type _opDn;
77
78 KOKKOS_INLINE_FUNCTION
79 Functor( outputValueViewType outputValues_,
80 inputPointViewType inputPoints_,
81 vinvViewType vinv_,
82 workViewType work_,
83 const ordinal_type opDn_ = 0 )
84 : _outputValues(outputValues_), _inputPoints(inputPoints_),
85 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
86
87 KOKKOS_INLINE_FUNCTION
88 void operator()(const size_type iter) const {
89 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
90 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
91
92 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
93 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
94
95 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
96
97 auto vcprop = Kokkos::common_view_alloc_prop(_work);
98 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
99
100 switch (opType) {
101 case OPERATOR_VALUE : {
102 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
103 Serial<opType>::getValues( output, input, work, _vinv );
104 break;
105 }
106 case OPERATOR_Dn : {
107 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
108 Serial<opType>::getValues( output, input, work, _vinv, _opDn );
109 break;
110 }
111 default: {
112 INTREPID2_TEST_FOR_ABORT( true,
113 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::Functor) operator is not supported");
114
115 }
116 }
117 }
118 };
119 };
120 }
121
134 template<typename DeviceType = void,
135 typename outputValueType = double,
136 typename pointValueType = double>
138 : public Basis<DeviceType,outputValueType,pointValueType> {
139 public:
141
145
146 using typename BasisBase::OutputViewType;
147 using typename BasisBase::PointViewType ;
148 using typename BasisBase::ScalarViewType;
149
152 Basis_HVOL_HEX_Cn_FEM(const ordinal_type order,
153 const EPointType pointType = POINTTYPE_EQUISPACED);
154
156
157 virtual
158 void
159 getValues( OutputViewType outputValues,
160 const PointViewType inputPoints,
161 const EOperator operatorType = OPERATOR_VALUE ) const override {
162#ifdef HAVE_INTREPID2_DEBUG
164 inputPoints,
165 operatorType,
166 this->getBaseCellTopology(),
167 this->getCardinality() );
168#endif
169 constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
170 Impl::Basis_HVOL_HEX_Cn_FEM::
171 getValues<DeviceType,numPtsPerEval>( outputValues,
172 inputPoints,
173 this->vinv_,
174 operatorType );
175 }
176
177 virtual void
178 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
179 ordinal_type& perThreadSpaceSize,
180 const PointViewType inputPoints,
181 const EOperator operatorType = OPERATOR_VALUE) const override;
182
183 KOKKOS_INLINE_FUNCTION
184 virtual void
185 getValues(
186 OutputViewType outputValues,
187 const PointViewType inputPoints,
188 const EOperator operatorType,
189 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
190 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
191 const ordinal_type subcellDim = -1,
192 const ordinal_type subcellOrdinal = -1) const override;
193
194
195 virtual
196 void
197 getDofCoords( ScalarViewType dofCoords ) const override {
198#ifdef HAVE_INTREPID2_DEBUG
199 // Verify rank of output array.
200 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
201 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
202 // Verify 0th dimension of output array.
203 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
204 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
205 // Verify 1st dimension of output array.
206 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
207 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
208#endif
209 Kokkos::deep_copy(dofCoords, this->dofCoords_);
210 }
211
212 virtual
213 void
214 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
215#ifdef HAVE_INTREPID2_DEBUG
216 // Verify rank of output array.
217 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
218 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
219 // Verify 0th dimension of output array.
220 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
221 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
222#endif
223 Kokkos::deep_copy(dofCoeffs, 1.0);
224 }
225
226 virtual
227 const char*
228 getName() const override {
229 return "Intrepid2_HVOL_HEX_Cn_FEM";
230 }
231
232 virtual
233 bool
234 requireOrientation() const override {
235 return false;
236 }
237
242
243 private:
244
246 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType> vinv_;
247 EPointType pointType_;
248 };
249
250}// namespace Intrepid2
251
253
254#endif
Header file for the abstract base class Intrepid2::Basis.
BasisPtr< typename Kokkos::HostSpace::device_type, OutputType, PointType > HostBasisPtr
Pointer to a Basis whose device type is on the host (Kokkos::HostSpace::device_type),...
void getValues_HVOL_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 HVOL-conforming FEM basis....
Definition file for FEM basis functions of degree n for H(vol) functions on HEX cells.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Implementation of the default HVOL-compatible FEM basis of degree n on Hexahedron cell.
virtual const char * getName() const override
Returns basis name.
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
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...
Kokkos::DynRankView< typename ScalarViewType::value_type, DeviceType > vinv_
inverse of Generalized Vandermonde matrix (isotropic order)
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPoints, 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 bool requireOrientation() const override
True if orientation is required.
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 basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
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.
See Intrepid2::Basis_HVOL_HEX_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
small utility functions