49 #ifndef __INTREPID2_HCURL_TET_IN_FEM_HPP__
50 #define __INTREPID2_HCURL_TET_IN_FEM_HPP__
56 #include "Teuchos_LAPACK.hpp"
92 #define CardinalityHCurlTet(order) (order*(order+2)*(order+3)/2)
105 template<EOperator opType>
107 template<
typename outputValueViewType,
108 typename inputPointViewType,
109 typename workViewType,
110 typename vinvViewType>
111 KOKKOS_INLINE_FUNCTION
113 getValues( outputValueViewType outputValues,
114 const inputPointViewType inputPoints,
116 const vinvViewType vinv );
119 KOKKOS_INLINE_FUNCTION
121 getWorkSizePerPoint(ordinal_type order) {
122 auto cardinality = CardinalityHCurlTet(order);
126 return 7*cardinality;
128 return getDkCardinality<opType,3>()*cardinality;
133 template<
typename ExecSpaceType, ordinal_type numPtsPerEval,
134 typename outputValueValueType,
class ...outputValueProperties,
135 typename inputPointValueType,
class ...inputPointProperties,
136 typename vinvValueType,
class ...vinvProperties>
138 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
139 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
140 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
141 const EOperator operatorType);
146 template<
typename outputValueViewType,
147 typename inputPointViewType,
148 typename vinvViewType,
149 typename workViewType,
151 ordinal_type numPtsEval>
153 outputValueViewType _outputValues;
154 const inputPointViewType _inputPoints;
155 const vinvViewType _coeffs;
158 KOKKOS_INLINE_FUNCTION
159 Functor( outputValueViewType outputValues_,
160 inputPointViewType inputPoints_,
161 vinvViewType coeffs_,
163 : _outputValues(outputValues_), _inputPoints(inputPoints_),
164 _coeffs(coeffs_), _work(work_) {}
166 KOKKOS_INLINE_FUNCTION
167 void operator()(
const size_type iter)
const {
171 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
172 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
174 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
176 auto vcprop = Kokkos::common_view_alloc_prop(_work);
177 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
180 case OPERATOR_VALUE : {
181 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
185 case OPERATOR_CURL: {
186 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
191 INTREPID2_TEST_FOR_ABORT(
true,
192 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::Functor) operator is not supported");
201 template<
typename ExecSpaceType = void,
202 typename outputValueType = double,
203 typename pointValueType =
double>
205 :
public Basis<ExecSpaceType,outputValueType,pointValueType> {
214 const EPointType pointType = POINTTYPE_EQUISPACED);
228 const EOperator operatorType = OPERATOR_VALUE)
const {
229 #ifdef HAVE_INTREPID2_DEBUG
230 Intrepid2::getValues_HCURL_Args(outputValues,
237 Impl::Basis_HCURL_TET_In_FEM::
238 getValues<ExecSpaceType,numPtsPerEval>( outputValues,
247 #ifdef HAVE_INTREPID2_DEBUG
249 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
250 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
252 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
253 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
255 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
256 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
258 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
264 #ifdef HAVE_INTREPID2_DEBUG
266 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
267 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
269 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
270 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
272 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
273 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
275 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
279 getExpansionCoeffs( scalarViewType coeffs )
const {
281 Kokkos::deep_copy(coeffs, this->
coeffs_);
287 return "Intrepid2_HCURL_TET_In_FEM";
300 Kokkos::DynRankView<scalarType,ExecSpaceType>
coeffs_;
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > scalarViewType
View type for scalars.
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Te...
virtual void getDofCoeffs(scalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< outputValueType, Kokkos::LayoutStride, ExecSpaceType > outputViewType
View type for basis value output.
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.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_3d_host
View type for 3d host array.
virtual const char * getName() const
Returns basis name.
Kokkos::DynRankView< scalarType, ExecSpaceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
Basis_HCURL_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
See Intrepid2::Basis_HCURL_TET_In_FEM.
Kokkos::View< ordinal_type *,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_1d_host
View type for 1d host array.
virtual void getDofCoords(scalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Kokkos::DynRankView< pointValueType, Kokkos::LayoutStride, ExecSpaceType > pointViewType
View type for input points.
virtual bool requireOrientation() const
True if orientation is required.
virtual void getValues(outputViewType outputValues, const pointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
Definition file for FEM basis functions of degree n for H(curl) functions on TET. ...
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
See Intrepid2::Basis_HCURL_TET_In_FEM.
See Intrepid2::Basis_HCURL_TET_In_FEM.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
Header file for the abstract base class Intrepid2::Basis.
Kokkos::View< ordinal_type **,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_2d_host
View type for 2d host array.