Intrepid2
Intrepid2_HGRAD_QUAD_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_QUAD_C2_FEM_HPP__
17 #define __INTREPID2_HGRAD_QUAD_C2_FEM_HPP__
18 
19 #include "Intrepid2_Basis.hpp"
21 
22 namespace Intrepid2 {
23 
69  namespace Impl {
70 
74  template<bool serendipity>
76  public:
77  typedef struct Quadrilateral<serendipity ? 8 : 9> 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( const typename DeviceType::execution_space& space,
97  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
98  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
99  const EOperator operatorType);
100 
104  template<typename outputValueViewType,
105  typename inputPointViewType,
106  EOperator opType>
107  struct Functor {
108  outputValueViewType _outputValues;
109  const inputPointViewType _inputPoints;
110 
111  KOKKOS_INLINE_FUNCTION
112  Functor( outputValueViewType outputValues_,
113  inputPointViewType inputPoints_ )
114  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
115 
116  KOKKOS_INLINE_FUNCTION
117  void operator()(const ordinal_type pt) const {
118  switch (opType) {
119  case OPERATOR_VALUE : {
120  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
121  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
122  Serial<opType>::getValues( output, input );
123  break;
124  }
125  case OPERATOR_GRAD :
126  case OPERATOR_CURL :
127  case OPERATOR_D1 :
128  case OPERATOR_D2 :
129  case OPERATOR_D3 :
130  case OPERATOR_D4 :
131  case OPERATOR_MAX : {
132  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
133  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
134  Serial<opType>::getValues( output, input );
135  break;
136  }
137  default: {
138  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
139  opType != OPERATOR_GRAD &&
140  opType != OPERATOR_CURL &&
141  opType != OPERATOR_D1 &&
142  opType != OPERATOR_D2 &&
143  opType != OPERATOR_D3 &&
144  opType != OPERATOR_D4 &&
145  opType != OPERATOR_MAX,
146  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::Serial::getValues) operator is not supported");
147  }
148  }
149  }
150  };
151  };
152  }
153 
154  template<bool serendipity,
155  typename DeviceType,
156  typename outputValueType,
157  typename pointValueType>
158  class Basis_HGRAD_QUAD_DEG2_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
159  public:
161  using typename BasisBase::ExecutionSpace;
162  using typename BasisBase::OrdinalTypeArray1DHost;
163  using typename BasisBase::OrdinalTypeArray2DHost;
164  using typename BasisBase::OrdinalTypeArray3DHost;
165 
169 
170  using typename BasisBase::OutputViewType;
171  using typename BasisBase::PointViewType;
172  using typename BasisBase::ScalarViewType;
173 
174  using BasisBase::getValues;
175 
176  virtual
177  void
178  getValues( const ExecutionSpace& space,
179  OutputViewType outputValues,
180  const PointViewType inputPoints,
181  const EOperator operatorType = OPERATOR_VALUE ) const override {
182 #ifdef HAVE_INTREPID2_DEBUG
183  // Verify arguments
184  Intrepid2::getValues_HGRAD_Args(outputValues,
185  inputPoints,
186  operatorType,
187  this->getBaseCellTopology(),
188  this->getCardinality() );
189 #endif
191  template getValues<DeviceType>(space,
192  outputValues,
193  inputPoints,
194  operatorType);
195  }
196 
197  virtual void
198  getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
199  ordinal_type& perThreadSpaceSize,
200  const PointViewType inputPointsconst,
201  const EOperator operatorType = OPERATOR_VALUE) const override;
202 
203  KOKKOS_INLINE_FUNCTION
204  virtual void
205  getValues(
206  OutputViewType outputValues,
207  const PointViewType inputPoints,
208  const EOperator operatorType,
209  const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
210  const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
211  const ordinal_type subcellDim = -1,
212  const ordinal_type subcellOrdinal = -1) const override;
213 
214  virtual
215  void
216  getDofCoords( ScalarViewType dofCoords ) const override {
217 #ifdef HAVE_INTREPID2_DEBUG
218  // Verify rank of output array.
219  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
220  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::getDofCoords) rank = 2 required for dofCoords array");
221  // Verify 0th dimension of output array.
222  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
223  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
224  // Verify 1st dimension of output array.
225  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
226  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
227 #endif
228  Kokkos::deep_copy(dofCoords, this->dofCoords_);
229  }
230 
231  virtual
232  void
233  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
234 #ifdef HAVE_INTREPID2_DEBUG
235  // Verify rank of output array.
236  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
237  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
238  // Verify 0th dimension of output array.
239  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
240  ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
241 #endif
242  Kokkos::deep_copy(dofCoeffs, 1.0);
243  }
244 
245  virtual
246  const char*
247  getName() const override {
248  return serendipity ? "Intrepid2_HGRAD_QUAD_I2_FEM" : "Intrepid2_HGRAD_QUAD_C2_FEM";
249  }
250 
259  BasisPtr<DeviceType, outputValueType, pointValueType>
260  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
261  if(subCellDim == 1) {
262  return Teuchos::rcp(new
264  }
265  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
266  }
267 
268  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
269  getHostBasis() const override{
271  }
272  };
273 
274  template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
275  using Basis_HGRAD_QUAD_C2_FEM = Basis_HGRAD_QUAD_DEG2_FEM<false, DeviceType, outputValueType, pointValueType>;
276 
277  template<typename DeviceType = void, typename outputValueType = double, typename pointValueType = double>
278  using Basis_HGRAD_QUAD_I2_FEM = Basis_HGRAD_QUAD_DEG2_FEM<true, DeviceType, outputValueType, pointValueType>;
279 
280 }// namespace Intrepid2
281 
283 
284 #endif
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.
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...
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, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Header file for the Intrepid2::Basis_HGRAD_LINE_C2_FEM class.
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< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
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::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.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell...
See Intrepid2::Basis_HGRAD_QUAD_DEG2_FEM.
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::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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 2 for H(grad) functions on QUAD cells...
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference 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.
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.