49#ifndef __INTREPID2_HDIV_TET_IN_FEM_HPP__
50#define __INTREPID2_HDIV_TET_IN_FEM_HPP__
57#include "Teuchos_LAPACK.hpp"
87#define CardinalityHDivTet(order) (order*(order+1)*(order+3)/2)
100 template<EOperator opType>
116 getWorkSizePerPoint(ordinal_type
order) {
133 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...>
outputValues,
134 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...>
inputPoints,
135 const Kokkos::DynRankView<vinvValueType, vinvProperties...>
vinv,
159 _coeffs(coeffs_), _work(
work_) {}
162 void operator()(
const size_type
iter)
const {
167 const auto input = Kokkos::subview( _inputPoints,
ptRange, Kokkos::ALL() );
169 typename workViewType::pointer_type
ptr = _work.data() + _work.extent(0)*
ptBegin*get_dimension_scalar(_work);
171 auto vcprop = Kokkos::common_view_alloc_prop(_work);
175 case OPERATOR_VALUE : {
176 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(),
ptRange, Kokkos::ALL() );
181 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(),
ptRange);
186 INTREPID2_TEST_FOR_ABORT(
true,
187 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::Functor) operator is not supported");
196template<
typename DeviceType = void,
197 typename outputValueType = double,
198 typename pointValueType =
double>
200 :
public Basis<DeviceType,outputValueType,pointValueType> {
209 const EPointType pointType = POINTTYPE_EQUISPACED);
223 const PointViewType inputPoints,
224 const EOperator operatorType = OPERATOR_VALUE)
const override {
225#ifdef HAVE_INTREPID2_DEBUG
233Impl::Basis_HDIV_TET_In_FEM::
234getValues<DeviceType,numPtsPerEval>( outputValues,
243#ifdef HAVE_INTREPID2_DEBUG
245 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
246 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
248 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
249 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
251 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
252 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
254 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
260#ifdef HAVE_INTREPID2_DEBUG
262 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
263 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
265 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
266 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
268 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
269 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
271 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
275 getExpansionCoeffs( ScalarViewType coeffs )
const {
277 Kokkos::deep_copy(coeffs, this->
coeffs_);
283 return "Intrepid2_HDIV_TET_In_FEM";
304 if(subCellDim == 2) {
305 return Teuchos::rcp(
new
309 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
320 Kokkos::DynRankView<scalarType,DeviceType>
coeffs_;
Header file for the abstract base class Intrepid2::Basis.
void getValues_HDIV_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 HDIV-conforming FEM basis....
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
Definition file for FEM basis functions of degree n for H(grad) functions on TET cells.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
Header file for the Intrepid2::Basis_HVOL_TRI_Cn_FEM class.
Implementation of the default H(div)-compatible Raviart-Thomas basis of arbitrary degree on Tetrahedr...
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...
virtual bool requireOrientation() const override
True if orientation is required.
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.
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...
Basis_HDIV_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
EPointType pointType_
type of lattice used for creating the DoF coordinates
Kokkos::DynRankView< scalarType, DeviceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Implementation of the default HVOL-compatible Lagrange basis of arbitrary degree on Triangle cell.
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_HDIV_TET_In_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
See Intrepid2::Basis_HDIV_TET_In_FEM.
See Intrepid2::Basis_HDIV_TET_In_FEM.