Intrepid2
Intrepid2_HGRAD_TET_COMP12_FEM.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Intrepid2 Package
4 //
5 // Copyright 2007 NTESS and the Intrepid2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
16 #ifndef __INTREPID2_HGRAD_TET_COMP12_FEM_HPP__
17 #define __INTREPID2_HGRAD_TET_COMP12_FEM_HPP__
18 
19 #include "Intrepid2_Basis.hpp"
20 
21 namespace Intrepid2 {
22 
67  namespace Impl {
68 
73  public:
74  typedef struct Tetrahedron<4> cell_topology_type;
75 
79  template<typename pointValueType>
80  KOKKOS_INLINE_FUNCTION
81  static ordinal_type
82  getLocalSubTet( const pointValueType x,
83  const pointValueType y,
84  const pointValueType z );
85 
89  template<EOperator opType>
90  struct Serial {
91  template<typename outputValueViewType,
92  typename inputPointViewType>
93  KOKKOS_INLINE_FUNCTION
94  static void
95  getValues( outputValueViewType outputValues,
96  const inputPointViewType inputPoints );
97 
98  };
99 
100  template<typename DeviceType,
101  typename outputValueValueType, class ...outputValueProperties,
102  typename inputPointValueType, class ...inputPointProperties>
103  static void
104  getValues( const typename DeviceType::execution_space& space,
105  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
106  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
107  const EOperator operatorType );
108 
112  template<typename outputValueViewType,
113  typename inputPointViewType,
114  EOperator opType>
115  struct Functor {
116  outputValueViewType _outputValues;
117  const inputPointViewType _inputPoints;
118 
119  KOKKOS_INLINE_FUNCTION
120  Functor( outputValueViewType outputValues_,
121  inputPointViewType inputPoints_ )
122  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
123 
124  KOKKOS_INLINE_FUNCTION
125  void operator()(const ordinal_type pt) const {
126  switch (opType) {
127  case OPERATOR_VALUE : {
128  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
129  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
130  Serial<opType>::getValues( output, input );
131  break;
132  }
133  case OPERATOR_GRAD :
134  case OPERATOR_MAX : {
135  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
136  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
137  Serial<opType>::getValues( output, input );
138  break;
139  }
140  default: {
141  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
142  opType != OPERATOR_GRAD &&
143  opType != OPERATOR_MAX,
144  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_COMP12_FEM::Functor::operator() operator is not supported");
145  }
146  }
147  }
148  };
149  };
150  }
151 
152  template<typename DeviceType = void,
153  typename outputValueType = double,
154  typename pointValueType = double>
155  class Basis_HGRAD_TET_COMP12_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
156  public:
158  using typename BasisBase::ExecutionSpace;
159 
160  using typename BasisBase::OrdinalTypeArray1DHost;
161  using typename BasisBase::OrdinalTypeArray2DHost;
162  using typename BasisBase::OrdinalTypeArray3DHost;
163 
164  using typename BasisBase::OutputViewType;
165  using typename BasisBase::PointViewType ;
166  using typename BasisBase::ScalarViewType;
167 
171 
172  using BasisBase::getValues;
173 
186  virtual
187  void
188  getValues( const ExecutionSpace& space,
189  OutputViewType outputValues,
190  const PointViewType inputPoints,
191  const EOperator operatorType = OPERATOR_VALUE ) const override {
192 #ifdef HAVE_INTREPID2_DEBUG
193  Intrepid2::getValues_HGRAD_Args(outputValues,
194  inputPoints,
195  operatorType,
196  this->getBaseCellTopology(),
197  this->getCardinality() );
198 #endif
199  Impl::Basis_HGRAD_TET_COMP12_FEM::
200  getValues<DeviceType>(space,
201  outputValues,
202  inputPoints,
203  operatorType);
204  }
205 
206  virtual void
207  getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
208  ordinal_type& perThreadSpaceSize,
209  const PointViewType inputPointsconst,
210  const EOperator operatorType = OPERATOR_VALUE) const override;
211 
212  KOKKOS_INLINE_FUNCTION
213  virtual void
214  getValues(
215  OutputViewType outputValues,
216  const PointViewType inputPoints,
217  const EOperator operatorType,
218  const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
219  const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
220  const ordinal_type subcellDim = -1,
221  const ordinal_type subcellOrdinal = -1) const override;
222 
223  virtual
224  void
225  getDofCoords( ScalarViewType dofCoords ) const override {
226 #ifdef HAVE_INTREPID2_DEBUG
227  // Verify rank of output array.
228  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
229  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_COMP12_FEM::getDofCoords) rank = 2 required for dofCoords array");
230  // Verify 0th dimension of output array.
231  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
232  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_COMP12_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
233  // Verify 1st dimension of output array.
234  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
235  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_COMP12_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
236 #endif
237  Kokkos::deep_copy(dofCoords, this->dofCoords_);
238  }
239 
240  virtual
241  const char*
242  getName() const override {
243  return "Intrepid2_HGRAD_TET_COMP12_FEM";
244  }
245 
246  };
247 }// namespace Intrepid2
248 
250 
251 #endif
See Intrepid2::Basis_HGRAD_TET_COMP12_FEM.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
ordinal_type getCardinality() const
Returns cardinality of the basis.
virtual void getValues(const ExecutionSpace &space, OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
FEM basis evaluation on a reference Tetrahedron 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...
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
virtual KOKKOS_INLINE_FUNCTION void getValues(OutputViewType, const PointViewType, const EOperator, const typename Kokkos::TeamPolicy< ExecutionSpace >::member_type &teamMember, const typename ExecutionSpace::scratch_memory_space &scratchStorage, const ordinal_type subcellDim=-1, const ordinal_type subcellOrdinal=-1) const
Team-level evaluation of basis functions on a reference cell.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
static KOKKOS_INLINE_FUNCTION ordinal_type getLocalSubTet(const pointValueType x, const pointValueType y, const pointValueType z)
See Intrepid2::Basis_HGRAD_TET_COMP12_FEM.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPointsconst, const EOperator operatorType=OPERATOR_VALUE) const override
Return the size of the scratch space, in bytes, needed for using the team-level implementation of get...
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Tetrahedron cell...
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Definition file for the composite H(grad)-compatible FEM basis of degree 1 on Tetrahedron cell with 1...
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.