16#ifndef __INTREPID2_HCURL_TRI_IN_FEM_HPP__
17#define __INTREPID2_HCURL_TRI_IN_FEM_HPP__
24#include "Teuchos_LAPACK.hpp"
48#define CardinalityHCurlTri(order) (order*(order+2))
62 template<EOperator opType>
64 template<
typename outputValueViewType,
65 typename inputPointViewType,
66 typename workViewType,
67 typename vinvViewType>
68 KOKKOS_INLINE_FUNCTION
70 getValues( outputValueViewType outputValues,
71 const inputPointViewType inputPoints,
73 const vinvViewType vinv );
75 KOKKOS_INLINE_FUNCTION
77 getWorkSizePerPoint(ordinal_type order) {
78 auto cardinality = CardinalityHCurlTri(order);
85 return getDkCardinality<opType,2>()*cardinality;
90 template<
typename DeviceType, ordinal_type numPtsPerEval,
91 typename outputValueValueType,
class ...outputValueProperties,
92 typename inputPointValueType,
class ...inputPointProperties,
93 typename vinvValueType,
class ...vinvProperties>
96 const typename DeviceType::execution_space& space,
97 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
98 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
99 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
100 const EOperator operatorType);
105 template<
typename outputValueViewType,
106 typename inputPointViewType,
107 typename vinvViewType,
108 typename workViewType,
110 ordinal_type numPtsEval>
112 outputValueViewType _outputValues;
113 const inputPointViewType _inputPoints;
114 const vinvViewType _coeffs;
117 KOKKOS_INLINE_FUNCTION
118 Functor( outputValueViewType outputValues_,
119 inputPointViewType inputPoints_,
120 vinvViewType coeffs_,
122 : _outputValues(outputValues_), _inputPoints(inputPoints_),
123 _coeffs(coeffs_), _work(work_) {}
125 KOKKOS_INLINE_FUNCTION
126 void operator()(
const size_type iter)
const {
130 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
131 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
134 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
136 auto vcprop = Kokkos::common_view_alloc_prop(_work);
137 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
140 case OPERATOR_VALUE : {
141 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
145 case OPERATOR_CURL: {
146 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
151 INTREPID2_TEST_FOR_ABORT(
true,
152 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::Functor) operator is not supported");
161template<
typename DeviceType = void,
162 typename outputValueType = double,
163 typename pointValueType =
double>
165 :
public Basis<DeviceType,outputValueType,pointValueType> {
176 const EPointType pointType = POINTTYPE_EQUISPACED);
194 const EOperator operatorType = OPERATOR_VALUE)
const override {
195#ifdef HAVE_INTREPID2_DEBUG
203 Impl::Basis_HCURL_TRI_In_FEM::
204 getValues<DeviceType,numPtsPerEval>(
214 ordinal_type& perThreadSpaceSize,
216 const EOperator operatorType = OPERATOR_VALUE)
const override;
218 KOKKOS_INLINE_FUNCTION
223 const EOperator operatorType,
224 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
225 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
226 const ordinal_type subcellDim = -1,
227 const ordinal_type subcellOrdinal = -1)
const override;
232#ifdef HAVE_INTREPID2_DEBUG
234 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
235 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
237 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
238 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
240 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
241 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
243 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
249#ifdef HAVE_INTREPID2_DEBUG
251 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 2, std::invalid_argument,
252 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
254 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
255 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
257 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
258 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
260 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
266 Kokkos::deep_copy(coeffs, this->
coeffs_);
272 return "Intrepid2_HCURL_TRI_In_FEM";
292 if(subCellDim == 1) {
293 return Teuchos::rcp(
new
297 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
308 Kokkos::DynRankView<scalarType,DeviceType>
coeffs_;
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.
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 void getValues(const ExecutionSpace &space, OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
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 getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPointsconst, const EOperator operatorType=OPERATOR_VALUE) const override
Return the size of the scratch space, in bytes, needed for using the team-level implementation of get...
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.
virtual KOKKOS_INLINE_FUNCTION void getValues(OutputViewType, const PointViewType, const EOperator, const typename Kokkos::TeamPolicy< ExecutionSpace >::member_type &teamMember, const typename ExecutionSpace::scratch_memory_space &scratchStorage, const ordinal_type subcellDim=-1, const ordinal_type subcellOrdinal=-1) const
Team-level evaluation of basis functions on a reference cell.
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.
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.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
See Intrepid2::Basis_HCURL_TRI_In_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
See Intrepid2::Basis_HCURL_TRI_In_FEM.
See Intrepid2::Basis_HCURL_TRI_In_FEM.
Triangle topology, 3 nodes.