49 #ifndef __INTREPID2_HDIV_TET_I1_FEM_HPP__ 
   50 #define __INTREPID2_HDIV_TET_I1_FEM_HPP__ 
  113       template<EOperator opType>
 
  115         template<
typename OutputViewType,
 
  116                  typename inputViewType>
 
  117         KOKKOS_INLINE_FUNCTION
 
  119         getValues(       OutputViewType output,
 
  120                    const inputViewType input );
 
  124       template<
typename ExecSpaceType,
 
  125                typename outputValueValueType, 
class ...outputValueProperties,
 
  126                typename inputPointValueType,  
class ...inputPointProperties>
 
  128       getValues(       Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
 
  129                  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...>  inputPoints,
 
  130                  const EOperator operatorType);
 
  135       template<
typename outputValueViewType,
 
  136                typename inputPointViewType,
 
  139               outputValueViewType _outputValues;
 
  140         const inputPointViewType  _inputPoints;
 
  142         KOKKOS_INLINE_FUNCTION
 
  143         Functor(       outputValueViewType outputValues_,
 
  144                        inputPointViewType  inputPoints_ )
 
  145           : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
 
  147         KOKKOS_INLINE_FUNCTION
 
  148         void operator()(
const ordinal_type pt)
 const {
 
  150           case OPERATOR_VALUE : {
 
  151             auto       output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
 
  152             const auto input  = Kokkos::subview( _inputPoints,                 pt, Kokkos::ALL() );
 
  156           case OPERATOR_DIV : {
 
  157             auto       output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
 
  158             const auto input  = Kokkos::subview( _inputPoints,                 pt, Kokkos::ALL() );
 
  163             INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
 
  164                                       opType != OPERATOR_DIV,
 
  165                                       ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::Serial::getValues) operator is not supported");
 
  174   template<
typename ExecSpaceType = void,
 
  175            typename outputValueType = double,
 
  176            typename pointValueType = 
double>
 
  197                const EOperator operatorType = OPERATOR_VALUE )
 const {
 
  198 #ifdef HAVE_INTREPID2_DEBUG 
  200       Intrepid2::getValues_HDIV_Args(outputValues,
 
  206       Impl::Basis_HDIV_TET_I1_FEM::
 
  207         getValues<ExecSpaceType>( outputValues,
 
  215 #ifdef HAVE_INTREPID2_DEBUG 
  217       INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
 
  218                                     ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
 
  220       INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
basisCardinality_, std::invalid_argument,
 
  221                                     ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
 
  223       INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
basisCellTopology_.getDimension(), std::invalid_argument,
 
  224                                     ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
 
  226       Kokkos::deep_copy(dofCoords, this->
dofCoords_);
 
  232 #ifdef HAVE_INTREPID2_DEBUG 
  234     INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
 
  235         ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
 
  237     INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
 
  238         ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
 
  240     INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
 
  241         ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
 
  243     Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
 
  249       return "Intrepid2_HDIV_TET_I1_FEM";
 
Definition file for FEM basis functions of degree 1 for H(div) functions on TET cells. 
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. 
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
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. 
See Intrepid2::Basis_HDIV_TET_I1_FEM. 
virtual void getDofCoords(ScalarViewType dofCoords) const 
Returns spatial locations (coordinates) of degrees of freedom on the reference cell. 
ordinal_type getCardinality() const 
Returns cardinality of the basis. 
See Intrepid2::Basis_HDIV_TET_I1_FEM. 
Implementation of the default H(div)-compatible FEM basis of degree 1 on a Tetrahedron cell...
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output. 
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const 
Evaluation of a FEM basis on a reference cell. 
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array. 
virtual const char * getName() const 
Returns basis name. 
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points. 
virtual bool requireOrientation() const 
True if orientation is required. 
See Intrepid2::Basis_HDIV_TET_I1_FEM. 
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...
Basis_HDIV_TET_I1_FEM()
Constructor. 
Header file for the abstract base class Intrepid2::Basis. 
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const 
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package for definition of b...