49#ifndef __INTREPID2_HGRAD_HEX_C1_FEM_HPP__
50#define __INTREPID2_HGRAD_HEX_C1_FEM_HPP__
103 template<EOperator opType>
105 template<
typename OutputViewType,
106 typename inputViewType>
107 KOKKOS_INLINE_FUNCTION
109 getValues( OutputViewType output,
110 const inputViewType input );
114 template<
typename DeviceType,
115 typename outputValueValueType,
class ...outputValueProperties,
116 typename inputPointValueType,
class ...inputPointProperties>
118 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
119 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
120 const EOperator operatorType);
125 template<
typename outputValueViewType,
126 typename inputPointViewType,
129 outputValueViewType _outputValues;
130 const inputPointViewType _inputPoints;
132 KOKKOS_INLINE_FUNCTION
133 Functor( outputValueViewType outputValues_,
134 inputPointViewType inputPoints_ )
135 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
137 KOKKOS_INLINE_FUNCTION
138 void operator()(
const ordinal_type pt)
const {
140 case OPERATOR_VALUE : {
141 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
142 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
149 case OPERATOR_MAX : {
150 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
151 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
156 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
157 opType != OPERATOR_GRAD &&
158 opType != OPERATOR_D2 &&
159 opType != OPERATOR_D3 &&
160 opType != OPERATOR_MAX,
161 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_C1_FEM::Serial::getValues) operator is not supported");
169 template<
typename DeviceType = void,
170 typename outputValueType = double,
171 typename pointValueType =
double>
191 const PointViewType inputPoints,
192 const EOperator operatorType = OPERATOR_VALUE )
const override {
193#ifdef HAVE_INTREPID2_DEBUG
201 Impl::Basis_HGRAD_HEX_C1_FEM::
202 getValues<DeviceType>( outputValues,
210#ifdef HAVE_INTREPID2_DEBUG
212 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
213 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
215 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
216 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
218 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
219 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
221 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
227#ifdef HAVE_INTREPID2_DEBUG
229 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
230 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
232 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
233 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
235 Kokkos::deep_copy(dofCoeffs, 1.0);
241 return "Intrepid2_HGRAD_HEX_C1_FEM";
254 if(subCellDim == 1) {
255 return Teuchos::rcp(
new
257 }
else if(subCellDim == 2) {
258 return Teuchos::rcp(
new
261 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
void getValues_HGRAD_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 HGRAD-conforming FEM basis....
Definition file for FEM basis functions of degree 1 for H(grad) functions on HEX cells.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C1_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Hexahedron cell.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual const char * getName() const override
Returns basis name.
Basis_HGRAD_HEX_C1_FEM()
Constructor.
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...
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(grad)-compatible FEM basis of degree 1 on Line cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Quadrilateral 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 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::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
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_HGRAD_HEX_C1_FEM.
See Intrepid2::Basis_HGRAD_HEX_C1_FEM.
See Intrepid2::Basis_HGRAD_HEX_C1_FEM.
Hexahedron topology, 8 nodes.