49#ifndef __INTREPID2_HCURL_TET_IN_FEM_HPP__
50#define __INTREPID2_HCURL_TET_IN_FEM_HPP__
57#include "Teuchos_LAPACK.hpp"
93#define CardinalityHCurlTet(order) (order*(order+2)*(order+3)/2)
106 template<EOperator opType>
122 getWorkSizePerPoint(ordinal_type
order) {
139 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...>
outputValues,
140 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...>
inputPoints,
141 const Kokkos::DynRankView<vinvValueType, vinvProperties...>
vinv,
165 _coeffs(coeffs_), _work(
work_) {}
168 void operator()(
const size_type
iter)
const {
173 const auto input = Kokkos::subview( _inputPoints,
ptRange, Kokkos::ALL() );
175 typename workViewType::pointer_type
ptr = _work.data() + _work.extent(0)*
ptBegin*get_dimension_scalar(_work);
177 auto vcprop = Kokkos::common_view_alloc_prop(_work);
181 case OPERATOR_VALUE : {
182 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(),
ptRange, Kokkos::ALL() );
186 case OPERATOR_CURL: {
187 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(),
ptRange, Kokkos::ALL() );
192 INTREPID2_TEST_FOR_ABORT(
true,
193 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::Functor) operator is not supported");
202template<
typename DeviceType = void,
203 typename outputValueType = double,
204 typename pointValueType =
double>
206 :
public Basis<DeviceType,outputValueType,pointValueType> {
215 const EPointType pointType = POINTTYPE_EQUISPACED);
229 const PointViewType inputPoints,
230 const EOperator operatorType = OPERATOR_VALUE)
const override {
231#ifdef HAVE_INTREPID2_DEBUG
239 Impl::Basis_HCURL_TET_In_FEM::
240 getValues<DeviceType,numPtsPerEval>( outputValues,
249#ifdef HAVE_INTREPID2_DEBUG
251 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
252 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
254 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
255 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
257 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
258 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
260 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
266#ifdef HAVE_INTREPID2_DEBUG
268 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
269 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
271 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
272 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
274 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
275 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
277 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
281 getExpansionCoeffs( ScalarViewType coeffs )
const {
283 Kokkos::deep_copy(coeffs, this->
coeffs_);
289 return "Intrepid2_HCURL_TET_In_FEM";
309 if(subCellDim == 1) {
310 return Teuchos::rcp(
new
313 }
else if(subCellDim == 2) {
314 return Teuchos::rcp(
new
318 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
330 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 TET.
Header file for the Intrepid2::Basis_HCURL_TRI_In_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Te...
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 getDofCoeffs(ScalarViewType dofCoeffs) const override
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual const char * getName() const override
Returns basis name.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Basis_HCURL_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
EPointType pointType_
type of lattice used for creating the DoF coordinates
virtual bool requireOrientation() const override
True if orientation is required.
Kokkos::DynRankView< scalarType, DeviceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
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...
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Tr...
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_TET_In_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
See Intrepid2::Basis_HCURL_TET_In_FEM.
See Intrepid2::Basis_HCURL_TET_In_FEM.
Tetrahedron topology, 4 nodes.