Intrepid2
Intrepid2_HGRAD_PYR_I2_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_PYR_I2_SERENDIPITY_FEM_HPP__
17 #define __INTREPID2_HGRAD_PYR_I2_SERENDIPITY_FEM_HPP__
18 
19 #include "Intrepid2_Basis.hpp"
22 
23 namespace Intrepid2 {
24 
70  namespace Impl {
71 
76  public:
77  typedef struct Pyramid<5> cell_topology_type;
81  template<EOperator opType>
82  struct Serial {
83  template<typename OutputViewType,
84  typename inputViewType>
85  KOKKOS_INLINE_FUNCTION
86  static void
87  getValues( OutputViewType output,
88  const inputViewType input );
89 
90  };
91 
92  template<typename DeviceType,
93  typename outputValueValueType, class ...outputValueProperties,
94  typename inputPointValueType, class ...inputPointProperties>
95  static void
96  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
97  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
98  const EOperator operatorType);
99 
103  template<typename outputValueViewType,
104  typename inputPointViewType,
105  EOperator opType>
106  struct Functor {
107  outputValueViewType _outputValues;
108  const inputPointViewType _inputPoints;
109 
110  KOKKOS_INLINE_FUNCTION
111  Functor( outputValueViewType outputValues_,
112  inputPointViewType inputPoints_ )
113  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
114 
115  KOKKOS_INLINE_FUNCTION
116  void operator()(const ordinal_type pt) const {
117  switch (opType) {
118  case OPERATOR_VALUE : {
119  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
120  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
121  Serial<opType>::getValues( output, input );
122  break;
123  }
124  case OPERATOR_GRAD :
125  case OPERATOR_CURL :
126  case OPERATOR_D2 : {
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_CURL &&
136  opType != OPERATOR_D2,
137  ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_I2_FEM::Serial::getValues) operator is not supported");
138  }
139  }
140  }
141  };
142 
143  };
144  }
145  template<typename DeviceType = void,
146  typename outputValueType = double,
147  typename pointValueType = double>
148  class Basis_HGRAD_PYR_I2_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
149  public:
150  using OrdinalTypeArray1DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray1DHost;
151  using OrdinalTypeArray2DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray2DHost;
152  using OrdinalTypeArray3DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray3DHost;
153 
157 
161 
163 
164  virtual
165  void
166  getValues( OutputViewType outputValues,
167  const PointViewType inputPoints,
168  const EOperator operatorType = OPERATOR_VALUE ) const override {
169 #ifdef HAVE_INTREPID2_DEBUG
170  // Verify arguments
171  Intrepid2::getValues_HGRAD_Args(outputValues,
172  inputPoints,
173  operatorType,
174  this->getBaseCellTopology(),
175  this->getCardinality() );
176 #endif
177  Impl::Basis_HGRAD_PYR_I2_FEM::
178  getValues<DeviceType>( outputValues,
179  inputPoints,
180  operatorType );
181  }
182 
183  virtual void
184  getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
185  ordinal_type& perThreadSpaceSize,
186  const PointViewType inputPointsconst,
187  const EOperator operatorType = OPERATOR_VALUE) const override;
188 
189  KOKKOS_INLINE_FUNCTION
190  virtual void
191  getValues(
192  OutputViewType outputValues,
193  const PointViewType inputPoints,
194  const EOperator operatorType,
195  const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
196  const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
197  const ordinal_type subcellDim = -1,
198  const ordinal_type subcellOrdinal = -1) const override;
199 
200  virtual
201  void
202  getDofCoords( ScalarViewType dofCoords ) const override {
203 #ifdef HAVE_INTREPID2_DEBUG
204  // Verify rank of output array.
205  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
206  ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_I2_FEM::getDofCoords) rank = 2 required for dofCoords array");
207  // Verify 0th dimension of output array.
208  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
209  ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_I2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
210  // Verify 1st dimension of output array.
211  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
212  ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_I2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
213 #endif
214  Kokkos::deep_copy(dofCoords, this->dofCoords_);
215  }
216 
217  virtual
218  void
219  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
220 #ifdef HAVE_INTREPID2_DEBUG
221  // Verify rank of output array.
222  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
223  ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_I2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
224  // Verify 0th dimension of output array.
225  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
226  ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_I2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
227 #endif
228  Kokkos::deep_copy(dofCoeffs, 1.0);
229  }
230 
231  virtual
232  const char*
233  getName() const override {
234  return "Intrepid2_HGRAD_PYR_I2_FEM";
235  }
236 
245  BasisPtr<DeviceType,outputValueType,pointValueType>
246  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
247  if(subCellDim == 1) {
248  return Teuchos::rcp(new
250  } else if(subCellDim == 2) {
251  if(subCellOrd == 0) //pyramid base
252  return Teuchos::rcp(new
254  else
255  return Teuchos::rcp(new
257  }
258  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
259  }
260 
261  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
262  getHostBasis() const override{
264  }
265 
266  };
267 }// namespace Intrepid2
268 
270 
271 #endif
ordinal_type getCardinality() const
Returns cardinality of the basis.
Implementation of an H(grad)-compatible FEM basis of degree 2 on a Pyramid cell.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C2_FEM class.
Definition file for FEM basis functions of degree 1 for H(grad) functions on PYR cells.
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...
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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.
See Intrepid2::Basis_HGRAD_PYR_I2_FEM.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Line cell.
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...
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell...
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Header file for the abstract base class Intrepid2::Basis.