Intrepid2
Intrepid2_HGRAD_TET_C1_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_C1_FEM_HPP__
17 #define __INTREPID2_HGRAD_TET_C1_FEM_HPP__
18 
19 #include "Intrepid2_Basis.hpp"
21 
22 namespace Intrepid2 {
23 
49  namespace Impl {
50 
55  public:
56  typedef struct Tetrahedron<4> cell_topology_type;
60  template<EOperator opType>
61  struct Serial {
62  template<typename OutputViewType,
63  typename inputViewType>
64  KOKKOS_INLINE_FUNCTION
65  static void
66  getValues( OutputViewType output,
67  const inputViewType input );
68 
69  };
70 
71  template<typename DeviceType,
72  typename outputValueValueType, class ...outputValueProperties,
73  typename inputPointValueType, class ...inputPointProperties>
74  static void
75  getValues( const typename DeviceType::execution_space& space,
76  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
77  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
78  const EOperator operatorType);
79 
83  template<typename outputValueViewType,
84  typename inputPointViewType,
85  EOperator opType>
86  struct Functor {
87  outputValueViewType _outputValues;
88  const inputPointViewType _inputPoints;
89 
90  KOKKOS_INLINE_FUNCTION
91  Functor( outputValueViewType outputValues_,
92  inputPointViewType inputPoints_ )
93  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
94 
95  KOKKOS_INLINE_FUNCTION
96  void operator()(const ordinal_type pt) const {
97  switch (opType) {
98  case OPERATOR_VALUE : {
99  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
100  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
101  Serial<opType>::getValues( output, input );
102  break;
103  }
104  case OPERATOR_GRAD :
105  case OPERATOR_MAX : {
106  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
107  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
108  Serial<opType>::getValues( output, input );
109  break;
110  }
111  default: {
112  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
113  opType != OPERATOR_GRAD &&
114  opType != OPERATOR_MAX,
115  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C1_FEM::Serial::getValues) operator is not supported");
116  }
117  }
118  }
119  };
120 
121 
122  };
123  }
124 
125  template<typename DeviceType = void,
126  typename outputValueType = double,
127  typename pointValueType = double>
128  class Basis_HGRAD_TET_C1_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
129  public:
131  using typename BasisBase::ExecutionSpace;
132  using typename BasisBase::OrdinalTypeArray1DHost;
133  using typename BasisBase::OrdinalTypeArray2DHost;
134  using typename BasisBase::OrdinalTypeArray3DHost;
135 
139 
140  using typename BasisBase::OutputViewType;
141  using typename BasisBase::PointViewType;
142  using typename BasisBase::ScalarViewType;
143 
144  using BasisBase::getValues;
145 
146  virtual
147  void
148  getValues( const ExecutionSpace& space,
149  OutputViewType outputValues,
150  const PointViewType inputPoints,
151  const EOperator operatorType = OPERATOR_VALUE ) const override {
152 #ifdef HAVE_INTREPID2_DEBUG
153  // Verify arguments
154  Intrepid2::getValues_HGRAD_Args(outputValues,
155  inputPoints,
156  operatorType,
157  this->getBaseCellTopology(),
158  this->getCardinality() );
159 #endif
160  Impl::Basis_HGRAD_TET_C1_FEM::
161  getValues<DeviceType>(space,
162  outputValues,
163  inputPoints,
164  operatorType);
165  }
166 
167  virtual void
168  getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
169  ordinal_type& perThreadSpaceSize,
170  const PointViewType inputPointsconst,
171  const EOperator operatorType = OPERATOR_VALUE) const override;
172 
173  KOKKOS_INLINE_FUNCTION
174  virtual void
175  getValues(
176  OutputViewType outputValues,
177  const PointViewType inputPoints,
178  const EOperator operatorType,
179  const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
180  const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
181  const ordinal_type subcellDim = -1,
182  const ordinal_type subcellOrdinal = -1) const override;
183 
184  virtual
185  void
186  getDofCoords( ScalarViewType dofCoords ) const override {
187 #ifdef HAVE_INTREPID2_DEBUG
188  // Verify rank of output array.
189  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
190  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
191  // Verify 0th dimension of output array.
192  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
193  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
194  // Verify 1st dimension of output array.
195  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
196  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
197 #endif
198  Kokkos::deep_copy(dofCoords, this->dofCoords_);
199  }
200 
201  virtual
202  void
203  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
204 #ifdef HAVE_INTREPID2_DEBUG
205  // Verify rank of output array.
206  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
207  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
208  // Verify 0th dimension of output array.
209  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
210  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
211 #endif
212  Kokkos::deep_copy(dofCoeffs, 1.0);
213  }
214 
215  virtual
216  const char*
217  getName() const override {
218  return "Intrepid2_HGRAD_TET_C1_FEM";
219  }
220 
229  BasisPtr<DeviceType, outputValueType, pointValueType>
230  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
231  if(subCellDim == 1) {
232  return Teuchos::rcp(new
234  } else if(subCellDim == 1) {
235  return Teuchos::rcp(new
237  }
238  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
239  }
240 
241  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
242  getHostBasis() const override{
244  }
245  };
246 }// namespace Intrepid2
247 
249 
250 #endif
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...
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Header file for the Intrepid2::Basis_HGRAD_TRI_C1_FEM class.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
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...
Definition file for FEM basis functions of degree 1 for H(grad) functions on TET cells.
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...
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Triangle 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 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.
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.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
See Intrepid2::Basis_HGRAD_TET_C1_FEM.
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
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 Tetrahedron cell...
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference 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.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.