Intrepid2
Intrepid2_HVOL_QUAD_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
47
48#ifndef __INTREPID2_HVOL_QUAD_CN_FEM_HPP__
49#define __INTREPID2_HVOL_QUAD_CN_FEM_HPP__
50
51#include "Intrepid2_Basis.hpp"
53
54namespace Intrepid2 {
55
56 namespace Impl {
57
62 public:
63 typedef struct Quadrilateral<4> 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 KOKKOS_INLINE_FUNCTION
82 static ordinal_type
83 getWorkSizePerPoint(ordinal_type order) {return 3*getPnCardinality<1>(order); }
84 };
85
86 template<typename DeviceType, ordinal_type numPtsPerEval,
87 typename outputValueValueType, class ...outputValueProperties,
88 typename inputPointValueType, class ...inputPointProperties,
89 typename vinvValueType, class ...vinvProperties>
90 static void
91 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
92 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
93 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
94 const EOperator operatorType);
95
99 template<typename outputValueViewType,
100 typename inputPointViewType,
101 typename vinvViewType,
102 typename workViewType,
103 EOperator opType,
104 ordinal_type numPtsEval>
105 struct Functor {
106 outputValueViewType _outputValues;
107 const inputPointViewType _inputPoints;
108 const vinvViewType _vinv;
109 workViewType _work;
110 const ordinal_type _opDn;
111
112 KOKKOS_INLINE_FUNCTION
113 Functor( outputValueViewType outputValues_,
114 inputPointViewType inputPoints_,
115 vinvViewType vinv_,
116 workViewType work_,
117 const ordinal_type opDn_ = 0 )
118 : _outputValues(outputValues_), _inputPoints(inputPoints_),
119 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
120
121 KOKKOS_INLINE_FUNCTION
122 void operator()(const size_type iter) const {
123 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
124 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
125
126 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
127 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
128
129 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
130
131 auto vcprop = Kokkos::common_view_alloc_prop(_work);
132 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
133
134 switch (opType) {
135 case OPERATOR_VALUE : {
136 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
137 Serial<opType>::getValues( output, input, work, _vinv );
138 break;
139 }
140 case OPERATOR_Dn : {
141 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
142 Serial<opType>::getValues( output, input, work, _vinv, _opDn );
143 break;
144 }
145 default: {
146 INTREPID2_TEST_FOR_ABORT( true,
147 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::Functor) operator is not supported");
148
149 }
150 }
151 }
152 };
153 };
154 }
155
162 template<typename DeviceType = void,
163 typename outputValueType = double,
164 typename pointValueType = double>
166 : public Basis<DeviceType,outputValueType,pointValueType> {
167 public:
171
174 Basis_HVOL_QUAD_Cn_FEM(const ordinal_type order,
175 const EPointType pointType = POINTTYPE_EQUISPACED);
176
180
181 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
182
183 virtual
184 void
185 getValues( OutputViewType outputValues,
186 const PointViewType inputPoints,
187 const EOperator operatorType = OPERATOR_VALUE ) const override {
188#ifdef HAVE_INTREPID2_DEBUG
190 inputPoints,
191 operatorType,
192 this->getBaseCellTopology(),
193 this->getCardinality() );
194#endif
195 constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
196 Impl::Basis_HVOL_QUAD_Cn_FEM::
197 getValues<DeviceType,numPtsPerEval>( outputValues,
198 inputPoints,
199 this->vinv_,
200 operatorType );
201 }
202
203 virtual
204 void
205 getDofCoords( ScalarViewType dofCoords ) const override {
206#ifdef HAVE_INTREPID2_DEBUG
207 // Verify rank of output array.
208 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
209 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
210 // Verify 0th dimension of output array.
211 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
212 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
213 // Verify 1st dimension of output array.
214 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
215 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
216#endif
217 Kokkos::deep_copy(dofCoords, this->dofCoords_);
218 }
219
220 virtual
221 void
222 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
223#ifdef HAVE_INTREPID2_DEBUG
224 // Verify rank of output array.
225 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
227 // Verify 0th dimension of output array.
228 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
229 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
230#endif
231 Kokkos::deep_copy(dofCoeffs, 1.0);
232 }
233
234 virtual
235 const char*
236 getName() const override {
237 return "Intrepid2_HVOL_QUAD_Cn_FEM";
238 }
239
240 virtual
241 bool
242 requireOrientation() const override {
243 return false;
244 }
245
250
251 private:
252
254 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType> vinv_;
255 EPointType pointType_;
256 };
257
258}// namespace Intrepid2
259
261
262#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),...
KOKKOS_INLINE_FUNCTION ordinal_type getPnCardinality(ordinal_type n)
Returns cardinality of Polynomials of order n (P^n).
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....
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Definition file for FEM basis functions of degree n for H(vol) functions on QUAD.
Kokkos::DynRankView< typename ScalarViewType::value_type, ExecutionSpace > vinv_
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...
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 bool requireOrientation() const override
True if orientation is required.
Basis_HVOL_QUAD_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference 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.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by 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
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HVOL_QUAD_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.