Intrepid2
Intrepid2_HVOL_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
48#ifndef __INTREPID2_HVOL_HEX_CN_FEM_HPP__
49#define __INTREPID2_HVOL_HEX_CN_FEM_HPP__
50
51#include "Intrepid2_Basis.hpp"
53
54namespace Intrepid2 {
55
56 namespace Impl {
61 public:
62 typedef struct Hexahedron<8> cell_topology_type;
66 template<EOperator opType>
67 struct Serial {
68 template<typename outputValueViewType,
69 typename inputPointViewType,
70 typename workViewType,
71 typename vinvViewType>
73 static void
77 const vinvViewType vinv,
78 const ordinal_type operatorDn = 0 );
79
81 static ordinal_type
82 getWorkSizePerPoint(ordinal_type order) {return 4*getPnCardinality<1>(order); }
83 };
84
85 template<typename DeviceType, ordinal_type numPtsPerEval,
86 typename outputValueValueType, class ...outputValueProperties,
87 typename inputPointValueType, class ...inputPointProperties,
88 typename vinvValueType, class ...vinvProperties>
89 static void
90 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
91 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
92 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
93 const EOperator operatorType );
94
98 template<typename outputValueViewType,
99 typename inputPointViewType,
100 typename vinvViewType,
101 typename workViewType,
102 EOperator opType,
103 ordinal_type numPtsEval>
104 struct Functor {
105 outputValueViewType _outputValues;
106 const inputPointViewType _inputPoints;
107 const vinvViewType _vinv;
108 workViewType _work;
109 const ordinal_type _opDn;
110
113 inputPointViewType inputPoints_,
114 vinvViewType vinv_,
116 const ordinal_type opDn_ = 0 )
117 : _outputValues(outputValues_), _inputPoints(inputPoints_),
118 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
119
121 void operator()(const size_type iter) const {
122 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
123 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
124
125 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
126 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
127
128 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
129
130 auto vcprop = Kokkos::common_view_alloc_prop(_work);
131 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
132
133 switch (opType) {
134 case OPERATOR_VALUE : {
135 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
136 Serial<opType>::getValues( output, input, work, _vinv );
137 break;
138 }
139 case OPERATOR_Dn : {
140 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
141 Serial<opType>::getValues( output, input, work, _vinv, _opDn );
142 break;
143 }
144 default: {
145 INTREPID2_TEST_FOR_ABORT( true,
146 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::Functor) operator is not supported");
147
148 }
149 }
150 }
151 };
152 };
153 }
154
167 template<typename DeviceType = void,
168 typename outputValueType = double,
169 typename pointValueType = double>
171 : public Basis<DeviceType,outputValueType,pointValueType> {
172 public:
176
179 Basis_HVOL_HEX_Cn_FEM(const ordinal_type order,
180 const EPointType pointType = POINTTYPE_EQUISPACED);
181
185
186 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
187
188 virtual
189 void
190 getValues( OutputViewType outputValues,
191 const PointViewType inputPoints,
192 const EOperator operatorType = OPERATOR_VALUE ) const override {
193#ifdef HAVE_INTREPID2_DEBUG
195 inputPoints,
196 operatorType,
197 this->getBaseCellTopology(),
198 this->getCardinality() );
199#endif
200 constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
201 Impl::Basis_HVOL_HEX_Cn_FEM::
202 getValues<DeviceType,numPtsPerEval>( outputValues,
203 inputPoints,
204 this->vinv_,
205 operatorType );
206 }
207
208
209 virtual
210 void
211 getDofCoords( ScalarViewType dofCoords ) const override {
212#ifdef HAVE_INTREPID2_DEBUG
213 // Verify rank of output array.
214 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
215 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
216 // Verify 0th dimension of output array.
217 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
218 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
219 // Verify 1st dimension of output array.
220 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
221 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
222#endif
223 Kokkos::deep_copy(dofCoords, this->dofCoords_);
224 }
225
226 virtual
227 void
228 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
229#ifdef HAVE_INTREPID2_DEBUG
230 // Verify rank of output array.
231 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
232 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
233 // Verify 0th dimension of output array.
234 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
235 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
236#endif
237 Kokkos::deep_copy(dofCoeffs, 1.0);
238 }
239
240 virtual
241 const char*
242 getName() const override {
243 return "Intrepid2_HVOL_HEX_Cn_FEM";
244 }
245
246 virtual
247 bool
248 requireOrientation() const override {
249 return false;
250 }
251
256
257 private:
258
260 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType> vinv_;
261 EPointType pointType_;
262 };
263
264}// namespace Intrepid2
265
267
268#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.
Basis_HVOL_HEX_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
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)
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.
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_HVOL_HEX_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
small utility functions