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