Intrepid2
Intrepid2_HGRAD_HEX_Cn_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_HEX_CN_FEM_HPP__
50#define __INTREPID2_HGRAD_HEX_CN_FEM_HPP__
51
52#include "Intrepid2_Basis.hpp"
54
55namespace Intrepid2 {
56
57 namespace Impl {
62 public:
63 typedef struct Hexahedron<8> cell_topology_type;
67 template<EOperator opType>
68 struct Serial {
69 template<typename outputValueViewType,
70 typename inputPointViewType,
71 typename workViewType,
72 typename vinvViewType>
73 KOKKOS_INLINE_FUNCTION
74 static void
75 getValues( outputValueViewType outputValues,
76 const inputPointViewType inputPoints,
77 workViewType work,
78 const vinvViewType vinv,
79 const ordinal_type operatorDn = 0 );
80 };
81
82 template<typename DeviceType, ordinal_type numPtsPerEval,
83 typename outputValueValueType, class ...outputValueProperties,
84 typename inputPointValueType, class ...inputPointProperties,
85 typename vinvValueType, class ...vinvProperties>
86 static void
87 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
88 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
89 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
90 const EOperator operatorType );
91
95 template<typename outputValueViewType,
96 typename inputPointViewType,
97 typename vinvViewType,
98 typename workViewType,
99 EOperator opType,
100 ordinal_type numPtsEval>
101 struct Functor {
102 outputValueViewType _outputValues;
103 const inputPointViewType _inputPoints;
104 const vinvViewType _vinv;
105 workViewType _work;
106 const ordinal_type _opDn;
107
108 KOKKOS_INLINE_FUNCTION
109 Functor( outputValueViewType outputValues_,
110 inputPointViewType inputPoints_,
111 vinvViewType vinv_,
112 workViewType work_,
113 const ordinal_type opDn_ = 0 )
114 : _outputValues(outputValues_), _inputPoints(inputPoints_),
115 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
116
117 KOKKOS_INLINE_FUNCTION
118 void operator()(const size_type iter) const {
119 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
120 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
121
122 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
123 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
124
125 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
126
127 auto vcprop = Kokkos::common_view_alloc_prop(_work);
128 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
129
130 switch (opType) {
131 case OPERATOR_VALUE : {
132 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
133 Serial<opType>::getValues( output, input, work, _vinv );
134 break;
135 }
136 case OPERATOR_CURL :
137 case OPERATOR_Dn : {
138 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
139 Serial<opType>::getValues( output, input, work, _vinv, _opDn );
140 break;
141 }
142 default: {
143 INTREPID2_TEST_FOR_ABORT( true,
144 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::Functor) operator is not supported");
145
146 }
147 }
148 }
149 };
150 };
151 }
152
164 template<typename DeviceType = void,
165 typename outputValueType = double,
166 typename pointValueType = double>
168 : public Basis<DeviceType,outputValueType,pointValueType> {
169 public:
173
177
178 private:
180 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType> vinv_;
181
183 EPointType pointType_;
184
185 public:
186
189 Basis_HGRAD_HEX_Cn_FEM(const ordinal_type order,
190 const EPointType pointType = POINTTYPE_EQUISPACED);
191
192 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
193
194 virtual
195 void
196 getValues( OutputViewType outputValues,
197 const PointViewType inputPoints,
198 const EOperator operatorType = OPERATOR_VALUE ) const override {
199#ifdef HAVE_INTREPID2_DEBUG
201 inputPoints,
202 operatorType,
203 this->getBaseCellTopology(),
204 this->getCardinality() );
205#endif
206 constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
207 Impl::Basis_HGRAD_HEX_Cn_FEM::
208 getValues<DeviceType,numPtsPerEval>( outputValues,
209 inputPoints,
210 this->vinv_,
211 operatorType );
212 }
213
214
215 virtual
216 void
217 getDofCoords( ScalarViewType dofCoords ) const override {
218#ifdef HAVE_INTREPID2_DEBUG
219 // Verify rank of output array.
220 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
221 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
222 // Verify 0th dimension of output array.
223 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
224 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
225 // Verify 1st dimension of output array.
226 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
227 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
228#endif
229 Kokkos::deep_copy(dofCoords, this->dofCoords_);
230 }
231
232 virtual
233 void
234 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
235#ifdef HAVE_INTREPID2_DEBUG
236 // Verify rank of output array.
237 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
238 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
239 // Verify 0th dimension of output array.
240 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
241 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
242#endif
243 Kokkos::deep_copy(dofCoeffs, 1.0);
244 }
245
246 virtual
247 const char*
248 getName() const override {
249 return "Intrepid2_HGRAD_HEX_Cn_FEM";
250 }
251
252 virtual
253 bool
254 requireOrientation() const override {
255 return (this->basisDegree_ > 2);
256 }
257
258 Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
259 getVandermondeInverse() {
260 return vinv_;
261 }
262
263 ordinal_type
264 getWorkSizePerPoint(const EOperator operatorType) {
265 return 4*getPnCardinality<1>(this->basisDegree_);
266 }
267
276 BasisPtr<DeviceType,outputValueType,pointValueType>
277 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
278 if(subCellDim == 1) {
279 return Teuchos::rcp(new
281 (this->basisDegree_, pointType_));
282 } else if(subCellDim == 2) {
283 return Teuchos::rcp(new
285 (this->basisDegree_, pointType_));
286 }
287 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
288 }
289
291 getHostBasis() const override{
293 }
294 };
295
296}// namespace Intrepid2
297
299
300#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 basis function of degree n for H(grad) functions on HEX cells.
Header file for the Intrepid2::Basis_HGRAD_QUAD_Cn_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron cell.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a 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...
Kokkos::DynRankView< typename ScalarViewType::value_type, DeviceType > vinv_
inverse of Generalized Vandermonde matrix (isotropic order)
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...
EPointType pointType_
type of lattice used for creating the DoF coordinates
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.
virtual const char * getName() const override
Returns basis name.
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1,...
Implementation of the default H(grad)-compatible FEM basis of degree n on Quadrilateral cell Implemen...
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 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.
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_HEX_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
small utility functions