Intrepid2
Intrepid2_HGRAD_TET_C2_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_C2_FEM_HPP__
17 #define __INTREPID2_HGRAD_TET_C2_FEM_HPP__
18 
19 #include "Intrepid2_Basis.hpp"
21 
22 namespace Intrepid2 {
23 
68  namespace Impl {
69 
74  public:
75  typedef struct Tetrahedron<10> cell_topology_type;
79  template<EOperator opType>
80  struct Serial {
81  template<typename OutputViewType,
82  typename inputViewType>
83  KOKKOS_INLINE_FUNCTION
84  static void
85  getValues( OutputViewType output,
86  const inputViewType input );
87 
88  };
89 
90  template<typename DeviceType,
91  typename outputValueValueType, class ...outputValueProperties,
92  typename inputPointValueType, class ...inputPointProperties>
93  static void
94  getValues( const typename DeviceType::execution_space& space,
95  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
96  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
97  const EOperator operatorType);
98 
102  template<typename outputValueViewType,
103  typename inputPointViewType,
104  EOperator opType>
105  struct Functor {
106  outputValueViewType _outputValues;
107  const inputPointViewType _inputPoints;
108 
109  KOKKOS_INLINE_FUNCTION
110  Functor( outputValueViewType outputValues_,
111  inputPointViewType inputPoints_ )
112  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
113 
114  KOKKOS_INLINE_FUNCTION
115  void operator()(const ordinal_type pt) const {
116  switch (opType) {
117  case OPERATOR_VALUE : {
118  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
119  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
120  Serial<opType>::getValues( output, input );
121  break;
122  }
123  case OPERATOR_GRAD :
124  case OPERATOR_D1 :
125  case OPERATOR_D2 :
126  case OPERATOR_MAX : {
127  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
128  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
129  Serial<opType>::getValues( output, input );
130  break;
131  }
132  default: {
133  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
134  opType != OPERATOR_GRAD &&
135  opType != OPERATOR_MAX,
136  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::Serial::getValues) operator is not supported");
137  }
138  }
139  }
140  };
141  };
142  }
143 
144  template<typename DeviceType = void,
145  typename outputValueType = double,
146  typename pointValueType = double>
147  class Basis_HGRAD_TET_C2_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
148  public:
150  using typename BasisBase::ExecutionSpace;
151  using typename BasisBase::OrdinalTypeArray1DHost;
152  using typename BasisBase::OrdinalTypeArray2DHost;
153  using typename BasisBase::OrdinalTypeArray3DHost;
154 
158 
159  using typename BasisBase::OutputViewType;
160  using typename BasisBase::PointViewType;
161  using typename BasisBase::ScalarViewType;
162 
163  using BasisBase::getValues;
164 
165  virtual
166  void
167  getValues( const ExecutionSpace& space,
168  OutputViewType outputValues,
169  const PointViewType inputPoints,
170  const EOperator operatorType = OPERATOR_VALUE ) const override {
171 #ifdef HAVE_INTREPID2_DEBUG
172  // Verify arguments
173  Intrepid2::getValues_HGRAD_Args(outputValues,
174  inputPoints,
175  operatorType,
176  this->getBaseCellTopology(),
177  this->getCardinality() );
178 #endif
179  Impl::Basis_HGRAD_TET_C2_FEM::
180  getValues<DeviceType>(space,
181  outputValues,
182  inputPoints,
183  operatorType);
184  }
185 
186  virtual void
187  getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
188  ordinal_type& perThreadSpaceSize,
189  const PointViewType inputPointsconst,
190  const EOperator operatorType = OPERATOR_VALUE) const override;
191 
192  KOKKOS_INLINE_FUNCTION
193  virtual void
194  getValues(
195  OutputViewType outputValues,
196  const PointViewType inputPoints,
197  const EOperator operatorType,
198  const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
199  const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
200  const ordinal_type subcellDim = -1,
201  const ordinal_type subcellOrdinal = -1) const override;
202 
203  virtual
204  void
205  getDofCoords( ScalarViewType dofCoords ) const override {
206 #ifdef HAVE_INTREPID2_DEBUG
207  // Verify rank of output array.
208  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
209  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getDofCoords) rank = 2 required for dofCoords array");
210  // Verify 0th dimension of output array.
211  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
212  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
213  // Verify 1st dimension of output array.
214  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
215  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
216 #endif
217  Kokkos::deep_copy(dofCoords, this->dofCoords_);
218  }
219 
220 
221  virtual
222  void
223  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
224 #ifdef HAVE_INTREPID2_DEBUG
225  // Verify rank of output array.
226  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
227  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
228  // Verify 0th dimension of output array.
229  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
230  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
231 #endif
232  Kokkos::deep_copy(dofCoeffs, 1.0);
233  }
234 
235  virtual
236  const char*
237  getName() const override {
238  return "Intrepid2_HGRAD_TET_C2_FEM";
239  }
240 
249  BasisPtr<DeviceType,outputValueType,pointValueType>
250  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
251  if(subCellDim == 1) {
253  } else if(subCellDim == 2) {
255  }
256 
257  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
258  }
259 
260  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
261  getHostBasis() const override{
263  }
264  };
265 }// namespace Intrepid2
266 
268 
269 #endif
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...
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.
See Intrepid2::Basis_HGRAD_TET_C2_FEM.
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...
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
virtual const char * getName() const override
Returns basis name.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Triangle cell...
Header file for the Intrepid2::Basis_HGRAD_TRI_C2_FEM class.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Line cell.
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::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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.
Definition file for FEM basis functions of degree 2 for H(grad) functions on TET cells.
virtual void getValues(const ExecutionSpace &space, 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.
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...
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
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< 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.
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.
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...