Intrepid2
Intrepid2_HDIV_HEX_In_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_HDIV_HEX_IN_FEM_HPP__
50#define __INTREPID2_HDIV_HEX_IN_FEM_HPP__
51
52#include "Intrepid2_Basis.hpp"
54
55namespace Intrepid2 {
56
57 namespace Impl {
58
63 public:
64 typedef struct Hexahedron<8> cell_topology_type;
68 template<EOperator opType>
69 struct Serial {
70 template<typename outputValueViewType,
71 typename inputPointViewType,
72 typename workViewType,
73 typename vinvViewType>
74 KOKKOS_INLINE_FUNCTION
75 static void
76 getValues( outputValueViewType outputValues,
77 const inputPointViewType inputPoints,
78 workViewType work,
79 const vinvViewType vinvLine,
80 const vinvViewType vinvBubble );
81
82 KOKKOS_INLINE_FUNCTION
83 static ordinal_type
84 getWorkSizePerPoint(ordinal_type order) {
85 return 2*getPnCardinality<1>(order)+2*getPnCardinality<1>(order-1);
86 }
87 };
88
89 template<typename DeviceType, ordinal_type numPtsPerEval,
90 typename outputValueValueType, class ...outputValueProperties,
91 typename inputPointValueType, class ...inputPointProperties,
92 typename vinvValueType, class ...vinvProperties>
93 static void
94 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
95 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
96 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvLine,
97 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvBubble,
98 const EOperator operatorType );
99
103 template<typename outputValueViewType,
104 typename inputPointViewType,
105 typename vinvViewType,
106 typename workViewType,
107 EOperator opType,
108 ordinal_type numPtsEval>
109 struct Functor {
110 outputValueViewType _outputValues;
111 const inputPointViewType _inputPoints;
112 const vinvViewType _vinvLine;
113 const vinvViewType _vinvBubble;
114 workViewType _work;
115
116 KOKKOS_INLINE_FUNCTION
117 Functor( outputValueViewType outputValues_,
118 inputPointViewType inputPoints_,
119 vinvViewType vinvLine_,
120 vinvViewType vinvBubble_,
121 workViewType work_)
122 : _outputValues(outputValues_), _inputPoints(inputPoints_),
123 _vinvLine(vinvLine_), _vinvBubble(vinvBubble_), _work(work_) {}
124
125 KOKKOS_INLINE_FUNCTION
126 void operator()(const size_type iter) const {
127 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
128 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
129
130 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
131 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
132
133 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
134
135 auto vcprop = Kokkos::common_view_alloc_prop(_work);
136 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
137
138 switch (opType) {
139 case OPERATOR_VALUE : {
140 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
141 Serial<opType>::getValues( output, input, work, _vinvLine, _vinvBubble );
142 break;
143 }
144 case OPERATOR_DIV : {
145 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
146 Serial<opType>::getValues( output, input, work, _vinvLine, _vinvBubble );
147 break;
148 }
149 default: {
150 INTREPID2_TEST_FOR_ABORT( true,
151 ">>> ERROR: (Intrepid2::Basis_HDIV_HEX_In_FEM::Functor) operator is not supported." );
152
153 }
154 }
155 }
156 };
157 };
158 }
159
167 template<typename DeviceType = void,
168 typename outputValueType = double,
169 typename pointValueType = double>
171 : public Basis<DeviceType,outputValueType,pointValueType> {
172 public:
176
179 Basis_HDIV_HEX_In_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 = 1;
201 Impl::Basis_HDIV_HEX_In_FEM::
202 getValues<DeviceType,numPtsPerEval>( outputValues,
203 inputPoints,
204 this->vinvLine_,
205 this->vinvBubble_,
206 operatorType );
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_HDIV_HEX_In_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_HDIV_HEX_In_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_HDIV_HEX_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
222#endif
223 Kokkos::deep_copy(dofCoords, this->dofCoords_);
224 }
225
226
227 virtual
228 void
229 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
230#ifdef HAVE_INTREPID2_DEBUG
231 // Verify rank of output array.
232 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
233 ">>> ERROR: (Intrepid2::Basis_HDIV_HEX_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
234 // Verify 0th dimension of output array.
235 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
236 ">>> ERROR: (Intrepid2::Basis_HDIV_HEX_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
237 // Verify 1st dimension of output array.
238 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
239 ">>> ERROR: (Intrepid2::Basis_HDIV_HEX_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
240#endif
241 Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
242 }
243
244 virtual
245 const char*
246 getName() const override {
247 return "Intrepid2_HDIV_HEX_In_FEM";
248 }
249
250 virtual
251 bool
252 requireOrientation() const override {
253 return true;
254 }
255
266 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
267
268 if(subCellDim == 2) {
269 return Teuchos::rcp(new
271 (this->basisDegree_-1,POINTTYPE_GAUSS));
272 }
273 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
274 }
275
277 getHostBasis() const override{
279 }
280 private:
281
283 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType> vinvLine_, vinvBubble_;
284 EPointType pointType_;
285 };
286
287}// namespace Intrepid2
288
289
290
292
293#endif
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
void getValues_HDIV_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 HDIV-conforming FEM basis....
Definition file for FEM basis functions of degree n for H(div) functions on HEX cells.
Header file for the Intrepid2::Basis_HGRAD_QUAD_Cn_FEM class.
Implementation of the default H(div)-compatible FEM basis on Hexahedron cell.
virtual bool requireOrientation() const override
True if orientation is required.
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...
virtual const char * getName() const override
Returns basis name.
Basis_HDIV_HEX_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Kokkos::DynRankView< typename ScalarViewType::value_type, DeviceType > vinvLine_
inverse of Generalized Vandermonde matrix (isotropic order)
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
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...
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
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.
void 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::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
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_HDIV_HEX_In_FEM.
small utility functions