49 #ifndef __INTREPID2_HCURL_TRI_IN_FEM_HPP__
50 #define __INTREPID2_HCURL_TRI_IN_FEM_HPP__
56 #include "Teuchos_LAPACK.hpp"
80 #define CardinalityHCurlTri(order) (order*(order+2))
89 typedef struct Triangle<3> cell_topology_type;
94 template<EOperator opType>
96 template<
typename outputValueViewType,
97 typename inputPointViewType,
98 typename workViewType,
99 typename vinvViewType>
100 KOKKOS_INLINE_FUNCTION
102 getValues( outputValueViewType outputValues,
103 const inputPointViewType inputPoints,
105 const vinvViewType vinv );
107 KOKKOS_INLINE_FUNCTION
109 getWorkSizePerPoint(ordinal_type order) {
110 auto cardinality = CardinalityHCurlTri(order);
115 return 5*cardinality;
117 return getDkCardinality<opType,2>()*cardinality;
122 template<
typename ExecSpaceType, ordinal_type numPtsPerEval,
123 typename outputValueValueType,
class ...outputValueProperties,
124 typename inputPointValueType,
class ...inputPointProperties,
125 typename vinvValueType,
class ...vinvProperties>
127 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
128 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
129 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
130 const EOperator operatorType);
135 template<
typename outputValueViewType,
136 typename inputPointViewType,
137 typename vinvViewType,
138 typename workViewType,
140 ordinal_type numPtsEval>
142 outputValueViewType _outputValues;
143 const inputPointViewType _inputPoints;
144 const vinvViewType _coeffs;
147 KOKKOS_INLINE_FUNCTION
148 Functor( outputValueViewType outputValues_,
149 inputPointViewType inputPoints_,
150 vinvViewType coeffs_,
152 : _outputValues(outputValues_), _inputPoints(inputPoints_),
153 _coeffs(coeffs_), _work(work_) {}
155 KOKKOS_INLINE_FUNCTION
156 void operator()(
const size_type iter)
const {
160 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
161 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
164 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
166 auto vcprop = Kokkos::common_view_alloc_prop(_work);
167 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
170 case OPERATOR_VALUE : {
171 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
175 case OPERATOR_CURL: {
176 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
181 INTREPID2_TEST_FOR_ABORT(
true,
182 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::Functor) operator is not supported");
191 template<
typename ExecSpaceType = void,
192 typename outputValueType = double,
193 typename pointValueType =
double>
195 :
public Basis<ExecSpaceType,outputValueType,pointValueType> {
204 const EPointType pointType = POINTTYPE_EQUISPACED);
219 const EOperator operatorType = OPERATOR_VALUE)
const {
220 #ifdef HAVE_INTREPID2_DEBUG
221 Intrepid2::getValues_HCURL_Args(outputValues,
228 Impl::Basis_HCURL_TRI_In_FEM::
229 getValues<ExecSpaceType,numPtsPerEval>( outputValues,
238 #ifdef HAVE_INTREPID2_DEBUG
240 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
241 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
243 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
244 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
246 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
247 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
249 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
255 #ifdef HAVE_INTREPID2_DEBUG
257 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
258 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
260 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
261 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
263 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
264 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
266 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
270 getExpansionCoeffs( ScalarViewType coeffs )
const {
272 Kokkos::deep_copy(coeffs, this->
coeffs_);
278 return "Intrepid2_HCURL_TRI_In_FEM";
291 Kokkos::DynRankView<scalarType,ExecSpaceType>
coeffs_;
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Kokkos::View< ordinal_type **, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
virtual void getDofCoords(ScalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
ordinal_type getCardinality() const
Returns cardinality of the basis.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
See Intrepid2::Basis_HCURL_TRI_In_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
virtual const char * getName() const
Returns basis name.
Basis_HCURL_TRI_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
virtual bool requireOrientation() const
True if orientation is required.
Kokkos::DynRankView< scalarType, ExecSpaceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
See Intrepid2::Basis_HCURL_TRI_In_FEM.
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
See Intrepid2::Basis_HCURL_TRI_In_FEM.
Definition file for FEM basis functions of degree n for H(curl) functions on TRI. ...
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Tr...
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
Header file for the abstract base class Intrepid2::Basis.