Intrepid2
Intrepid2_HCURL_TRI_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_HCURL_TRI_IN_FEM_HPP__
50#define __INTREPID2_HCURL_TRI_IN_FEM_HPP__
51
52#include "Intrepid2_Basis.hpp"
55
57#include "Teuchos_LAPACK.hpp"
58
59namespace Intrepid2 {
60
81#define CardinalityHCurlTri(order) (order*(order+2))
82
83namespace Impl {
84
89public:
90 typedef struct Triangle<3> cell_topology_type;
91
95 template<EOperator opType>
96 struct Serial {
97 template<typename outputValueViewType,
98 typename inputPointViewType,
99 typename workViewType,
100 typename vinvViewType>
101 KOKKOS_INLINE_FUNCTION
102 static void
103 getValues( outputValueViewType outputValues,
104 const inputPointViewType inputPoints,
105 workViewType work,
106 const vinvViewType vinv );
107
108 KOKKOS_INLINE_FUNCTION
109 static ordinal_type
110 getWorkSizePerPoint(ordinal_type order) {
111 auto cardinality = CardinalityHCurlTri(order);
112 switch (opType) {
113 case OPERATOR_GRAD:
114 case OPERATOR_CURL:
115 case OPERATOR_D1:
116 return 5*cardinality;
117 default:
118 return getDkCardinality<opType,2>()*cardinality;
119 }
120 }
121 };
122
123 template<typename DeviceType, ordinal_type numPtsPerEval,
124 typename outputValueValueType, class ...outputValueProperties,
125 typename inputPointValueType, class ...inputPointProperties,
126 typename vinvValueType, class ...vinvProperties>
127 static void
128 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
129 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
130 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
131 const EOperator operatorType);
132
136 template<typename outputValueViewType,
137 typename inputPointViewType,
138 typename vinvViewType,
139 typename workViewType,
140 EOperator opType,
141 ordinal_type numPtsEval>
142 struct Functor {
143 outputValueViewType _outputValues;
144 const inputPointViewType _inputPoints;
145 const vinvViewType _coeffs;
146 workViewType _work;
147
148 KOKKOS_INLINE_FUNCTION
149 Functor( outputValueViewType outputValues_,
150 inputPointViewType inputPoints_,
151 vinvViewType coeffs_,
152 workViewType work_)
153 : _outputValues(outputValues_), _inputPoints(inputPoints_),
154 _coeffs(coeffs_), _work(work_) {}
155
156 KOKKOS_INLINE_FUNCTION
157 void operator()(const size_type iter) const {
158 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
159 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
160
161 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
162 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
163
164
165 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
166
167 auto vcprop = Kokkos::common_view_alloc_prop(_work);
168 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
169
170 switch (opType) {
171 case OPERATOR_VALUE : {
172 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
173 Serial<opType>::getValues( output, input, work, _coeffs );
174 break;
175 }
176 case OPERATOR_CURL: {
177 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
178 Serial<opType>::getValues( output, input, work, _coeffs );
179 break;
180 }
181 default: {
182 INTREPID2_TEST_FOR_ABORT( true,
183 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::Functor) operator is not supported");
184
185 }
186 }
187 }
188 };
189};
190}
191
192template<typename DeviceType = void,
193 typename outputValueType = double,
194 typename pointValueType = double>
196 : public Basis<DeviceType,outputValueType,pointValueType> {
197 public:
201
204 Basis_HCURL_TRI_In_FEM(const ordinal_type order,
205 const EPointType pointType = POINTTYPE_EQUISPACED);
206
208
212
214
215 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
216
217 virtual
218 void
219 getValues( OutputViewType outputValues,
220 const PointViewType inputPoints,
221 const EOperator operatorType = OPERATOR_VALUE) const override {
222#ifdef HAVE_INTREPID2_DEBUG
224 inputPoints,
225 operatorType,
226 this->getBaseCellTopology(),
227 this->getCardinality() );
228#endif
229 constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
230 Impl::Basis_HCURL_TRI_In_FEM::
231 getValues<DeviceType,numPtsPerEval>( outputValues,
232 inputPoints,
233 this->coeffs_,
234 operatorType);
235 }
236
237 virtual
238 void
239 getDofCoords( ScalarViewType dofCoords ) const override {
240#ifdef HAVE_INTREPID2_DEBUG
241 // Verify rank of output array.
242 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
243 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
244 // Verify 0th dimension of output array.
245 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
246 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
247 // Verify 1st dimension of output array.
248 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
249 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
250#endif
251 Kokkos::deep_copy(dofCoords, this->dofCoords_);
252 }
253
254 virtual
255 void
256 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
257#ifdef HAVE_INTREPID2_DEBUG
258 // Verify rank of output array.
259 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
260 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
261 // Verify 0th dimension of output array.
262 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
263 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
264 // Verify 1st dimension of output array.
265 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
266 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
267#endif
268 Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
269 }
270
271 void
272 getExpansionCoeffs( ScalarViewType coeffs ) const {
273 // has to be same rank and dimensions
274 Kokkos::deep_copy(coeffs, this->coeffs_);
275 }
276
277 virtual
278 const char*
279 getName() const override {
280 return "Intrepid2_HCURL_TRI_In_FEM";
281 }
282
283 virtual
284 bool
285 requireOrientation() const override {
286 return true;
287 }
288
299 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
300 if(subCellDim == 1) {
301 return Teuchos::rcp(new
303 ( this->basisDegree_ - 1, pointType_) );
304 }
305 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
306 }
307
309 getHostBasis() const override{
311 }
312 private:
313
316 Kokkos::DynRankView<scalarType,DeviceType> coeffs_;
317
319 EPointType pointType_;
320
321};
322
323}// namespace Intrepid2
324
326
327#endif
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
void getValues_HCURL_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 HCURL-conforming FEM basis....
Definition file for FEM basis functions of degree n for H(curl) functions on TRI.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates,...
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Tr...
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...
Kokkos::DynRankView< scalarType, DeviceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
virtual bool requireOrientation() const override
True if orientation is required.
virtual const char * getName() const override
Returns basis name.
EPointType pointType_
type of lattice used for creating the DoF coordinates
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 getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the 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...
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,...
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::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.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
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_HCURL_TRI_In_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
small utility functions