Intrepid2
Intrepid2_HDIV_QUAD_I1_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_HDIV_QUAD_I1_FEM_HPP__
17 #define __INTREPID2_HDIV_QUAD_I1_FEM_HPP__
18 
19 #include "Intrepid2_Basis.hpp"
21 
22 namespace Intrepid2 {
23 
72  namespace Impl {
73 
78  public:
79  typedef struct Quadrilateral<4> cell_topology_type;
83  template<EOperator opType>
84  struct Serial {
85  template<typename OutputViewType,
86  typename inputViewType>
87  KOKKOS_INLINE_FUNCTION
88  static void
89  getValues( OutputViewType output,
90  const inputViewType input );
91 
92  };
93 
94  template<typename DeviceType,
95  typename outputValueValueType, class ...outputValueProperties,
96  typename inputPointValueType, class ...inputPointProperties>
97  static void
98  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
99  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
100  const EOperator operatorType);
101 
105  template<typename outputValueViewType,
106  typename inputPointViewType,
107  EOperator opType>
108  struct Functor {
109  outputValueViewType _outputValues;
110  const inputPointViewType _inputPoints;
111 
112  KOKKOS_INLINE_FUNCTION
113  Functor( outputValueViewType outputValues_,
114  inputPointViewType inputPoints_ )
115  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
116 
117  KOKKOS_INLINE_FUNCTION
118  void operator()(const ordinal_type pt) const {
119  switch (opType) {
120  case OPERATOR_VALUE : {
121  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
122  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
123  Serial<opType>::getValues( output, input );
124  break;
125  }
126  case OPERATOR_DIV :{
127  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt);
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_DIV,
135  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::Serial::getValues) operator is not supported");
136  }
137  }
138  }
139  };
140  };
141  }
142 
143  template<typename DeviceType = void,
144  typename outputValueType = double,
145  typename pointValueType = double>
146  class Basis_HDIV_QUAD_I1_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
147  public:
149  using OrdinalTypeArray1DHost = typename BasisBase::OrdinalTypeArray1DHost;
150  using OrdinalTypeArray2DHost = typename BasisBase::OrdinalTypeArray2DHost;
151  using OrdinalTypeArray3DHost = typename BasisBase::OrdinalTypeArray3DHost;
152 
153 
157 
158  using OutputViewType = typename BasisBase::OutputViewType;
159  using PointViewType = typename BasisBase::PointViewType;
160  using ScalarViewType = typename BasisBase::ScalarViewType;
161 
162  using BasisBase::getValues;
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_HDIV_Args(outputValues,
172  inputPoints,
173  operatorType,
174  this->getBaseCellTopology(),
175  this->getCardinality() );
176 #endif
177  Impl::Basis_HDIV_QUAD_I1_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_HDIV_QUAD_I1_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->basisCardinality_, std::invalid_argument,
209  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_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_HDIV_QUAD_I1_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() != 2, std::invalid_argument,
223  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoeffs) rank = 2 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_HDIV_QUAD_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
227  // Verify 1st dimension of output array.
228  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
229  ">>> ERROR: (Intrepid2::Basis_HDIV_QUAD_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
230 #endif
231  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
232  }
233 
234  virtual
235  const char*
236  getName() const override {
237  return "Intrepid2_HDIV_QUAD_I1_FEM";
238  }
239 
240  virtual
241  bool
242  requireOrientation() const override {
243  return true;
244  }
245 
255  BasisPtr<DeviceType,outputValueType,pointValueType>
256  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
257  if(subCellDim == 1)
258  return Teuchos::rcp( new
259  Basis_HVOL_C0_FEM<DeviceType,outputValueType,pointValueType>(shards::getCellTopologyData<shards::Line<2> >()));
260 
261  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
262  }
263 
264  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
265  getHostBasis() const override{
267  }
268  };
269 
270 }// namespace Intrepid2
271 
273 
274 #endif
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
virtual bool requireOrientation() const override
True if orientation is required.
Implementation of the default H(div)-compatible FEM basis of degree 1 on Quadrilateral cell...
ordinal_type getCardinality() const
Returns cardinality of the basis.
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 const char * getName() const override
Returns basis name.
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::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Header file for the Intrepid2::Basis_HVOL_C0_FEM class.
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::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Definition file for FEM basis functions of degree 1 for H(div) functions on QUAD cells.
Implementation of the default HVOL-compatible FEM contstant basis on triangle, quadrilateral, hexahedron and tetrahedron cells.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
See Intrepid2::Basis_HDIV_QUAD_I1_FEM.
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.