48 #ifndef __INTREPID2_HVOL_C0_FEM_HPP__ 
   49 #define __INTREPID2_HVOL_C0_FEM_HPP__ 
   66       template<EOperator opType>
 
   68         template<
typename OutputViewType,
 
   69                  typename inputViewType>
 
   70         KOKKOS_INLINE_FUNCTION
 
   72         getValues(       OutputViewType output,
 
   73                    const inputViewType input );
 
   77       template<
typename ExecSpaceType,
 
   78                typename outputValueValueType, 
class ...outputValueProperties,
 
   79                typename inputPointValueType,  
class ...inputPointProperties>
 
   81       getValues(       Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
 
   82                  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...>  inputPoints,
 
   83                  const EOperator operatorType);
 
   88       template<
typename outputValueViewType,
 
   89                typename inputPointViewType,
 
   92               outputValueViewType _outputValues;
 
   93         const inputPointViewType  _inputPoints;
 
   95         KOKKOS_INLINE_FUNCTION
 
   96         Functor(       outputValueViewType outputValues_,
 
   97                        inputPointViewType  inputPoints_ )
 
   98           : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
 
  100         KOKKOS_INLINE_FUNCTION
 
  101         void operator()(
const ordinal_type pt)
 const {
 
  103           case OPERATOR_VALUE : {
 
  104             auto       output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
 
  105             const auto input  = Kokkos::subview( _inputPoints,                 pt, Kokkos::ALL() );
 
  109           case OPERATOR_MAX : {
 
  110             auto       output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
 
  111             const auto input  = Kokkos::subview( _inputPoints,                 pt, Kokkos::ALL() );
 
  116             INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
 
  117                                       opType != OPERATOR_MAX,
 
  118                                       ">>> ERROR: (Intrepid2::Basis_HVOL_C0_FEM::Serial::getValues) operator is not supported");
 
  129   template<
typename ExecSpaceType = void,
 
  130            typename outputValueType = double,
 
  131            typename pointValueType = 
double>
 
  156                const EOperator operatorType = OPERATOR_VALUE )
 const {
 
  157 #ifdef HAVE_INTREPID2_DEBUG 
  159       Intrepid2::getValues_HGRAD_Args(outputValues,
 
  165       Impl::Basis_HVOL_C0_FEM::
 
  166         getValues<ExecSpaceType>( outputValues,
 
  174 #ifdef HAVE_INTREPID2_DEBUG 
  176       INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
 
  177                                     ">>> ERROR: (Intrepid2::Basis_HVOL_C0_FEM::getDofCoords) rank = 2 required for dofCoords array");
 
  179       INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
 
  180                                     ">>> ERROR: (Intrepid2::Basis_HVOL_C0_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
 
  182       INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
 
  183                                     ">>> ERROR: (Intrepid2::Basis_HVOL_C0_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
 
  185       Kokkos::deep_copy(dofCoords, this->
dofCoords_);
 
  191 #ifdef HAVE_INTREPID2_DEBUG 
  193       INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
 
  194                                     ">>> ERROR: (Intrepid2::Basis_HVOL_C0_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
 
  196       INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
 
  197                                     ">>> ERROR: (Intrepid2::Basis_HVOL_C0_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
 
  199       Kokkos::deep_copy(dofCoeffs, 1.0);
 
  204       return "Intrepid2_HVOL_C0_FEM";
 
See Intrepid2::Basis_HVOL_C0_FEM. 
 
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. 
 
virtual void getDofCoords(ScalarViewType dofCoords) const 
Returns spatial locations (coordinates) of degrees of freedom on the reference cell. 
 
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. 
 
ordinal_type getCardinality() const 
Returns cardinality of the basis. 
 
See Intrepid2::Basis_HVOL_C0_FEM. 
 
See Intrepid2::Basis_HVOL_C0_FEM. 
 
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const 
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
 
virtual const char * getName() const 
Returns basis name. 
 
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output. 
 
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array. 
 
Implementation of the default HVOL-compatible FEM contstant basis on triangle, quadrilateral, hexahedron and tetrahedron cells. 
 
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points. 
 
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars. 
 
Basis_HVOL_C0_FEM()=delete
Constructor. 
 
Definition file FEM basis functions of degree 0 for H(vol) functions on all supported topologies...
 
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const 
Evaluation of a FEM basis on a reference cell. 
 
Header file for the abstract base class Intrepid2::Basis.