Intrepid2
Intrepid2_HCURL_WEDGE_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_HCURL_WEDGE_I1_FEM_HPP__
17 #define __INTREPID2_HCURL_WEDGE_I1_FEM_HPP__
18 
19 #include "Intrepid2_Basis.hpp"
22 
23 namespace Intrepid2 {
24 
78  namespace Impl {
79 
84  public:
85  typedef struct Wedge<6> cell_topology_type;
86 
90  template<EOperator opType>
91  struct Serial {
92  template<typename OutputViewType,
93  typename inputViewType>
94  KOKKOS_INLINE_FUNCTION
95  static void
96  getValues( OutputViewType output,
97  const inputViewType input );
98 
99  };
100 
101  template<typename DeviceType,
102  typename outputValueValueType, class ...outputValueProperties,
103  typename inputPointValueType, class ...inputPointProperties>
104  static void
105  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
106  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
107  const EOperator operatorType);
108 
112  template<typename outputValueViewType,
113  typename inputPointViewType,
114  EOperator opType>
115  struct Functor {
116  outputValueViewType _outputValues;
117  const inputPointViewType _inputPoints;
118 
119  KOKKOS_INLINE_FUNCTION
120  Functor( outputValueViewType outputValues_,
121  inputPointViewType inputPoints_ )
122  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
123 
124  KOKKOS_INLINE_FUNCTION
125  void operator()(const ordinal_type pt) const {
126  switch (opType) {
127  case OPERATOR_VALUE : {
128  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
129  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
130  Serial<opType>::getValues( output, input );
131  break;
132  }
133  case OPERATOR_CURL : {
134  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
135  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
136  Serial<opType>::getValues( output, input );
137  break;
138  }
139  default: {
140  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
141  opType != OPERATOR_CURL,
142  ">>> ERROR: (Intrepid2::Basis_HCURL_WEDGE_I1_FEM::Serial::getValues) operator is not supported");
143  }
144  }
145  }
146  };
147  };
148  }
149 
150  template<typename DeviceType = void,
151  typename outputValueType = double,
152  typename pointValueType = double>
153  class Basis_HCURL_WEDGE_I1_FEM : public Basis<DeviceType, outputValueType, pointValueType> {
154  public:
155  using OrdinalTypeArray1DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray1DHost;
156  using OrdinalTypeArray2DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray2DHost;
157  using OrdinalTypeArray3DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray3DHost;
158 
162 
166 
168 
169  virtual
170  void
171  getValues( OutputViewType outputValues,
172  const PointViewType inputPoints,
173  const EOperator operatorType = OPERATOR_VALUE ) const override {
174 #ifdef HAVE_INTREPID2_DEBUG
175  // Verify arguments
176  Intrepid2::getValues_HCURL_Args(outputValues,
177  inputPoints,
178  operatorType,
179  this->getBaseCellTopology(),
180  this->getCardinality() );
181 #endif
182  Impl::Basis_HCURL_WEDGE_I1_FEM::
183  getValues<DeviceType>( outputValues,
184  inputPoints,
185  operatorType );
186  }
187 
188  virtual void
189  getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
190  ordinal_type& perThreadSpaceSize,
191  const PointViewType inputPointsconst,
192  const EOperator operatorType = OPERATOR_VALUE) const override;
193 
194  KOKKOS_INLINE_FUNCTION
195  virtual void
196  getValues(
197  OutputViewType outputValues,
198  const PointViewType inputPoints,
199  const EOperator operatorType,
200  const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
201  const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
202  const ordinal_type subcellDim = -1,
203  const ordinal_type subcellOrdinal = -1) const override;
204 
205  virtual
206  void
207  getDofCoords( ScalarViewType dofCoords ) const override {
208 #ifdef HAVE_INTREPID2_DEBUG
209  // Verify rank of output array.
210  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
211  ">>> ERROR: (Intrepid2::Basis_HCURL_WEDGE_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
212  // Verify 0th dimension of output array.
213  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
214  ">>> ERROR: (Intrepid2::Basis_HCURL_WEDGE_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
215  // Verify 1st dimension of output array.
216  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
217  ">>> ERROR: (Intrepid2::Basis_HCURL_WEDGE_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
218 #endif
219  Kokkos::deep_copy(dofCoords, this->dofCoords_);
220  }
221 
222  virtual
223  void
224  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
225  #ifdef HAVE_INTREPID2_DEBUG
226  // Verify rank of output array.
227  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
228  ">>> ERROR: (Intrepid2::Basis_HCURL_WEDGE_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
229  // Verify 0th dimension of output array.
230  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
231  ">>> ERROR: (Intrepid2::Basis_HCURL_WEDGE_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
232  // Verify 1st dimension of output array.
233  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
234  ">>> ERROR: (Intrepid2::Basis_HCURL_WEDGE_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
235  #endif
236  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
237  }
238 
239  virtual
240  const char*
241  getName() const override {
242  return "Intrepid2_HCURL_WEDGE_I1_FEM";
243  }
244 
245  virtual
246  bool
247  requireOrientation() const override {
248  return true;
249  }
250 
260  BasisPtr<DeviceType,outputValueType,pointValueType>
261  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
262 
263  if(subCellDim == 1)
264  return Teuchos::rcp( new
265  Basis_HVOL_C0_FEM<DeviceType,outputValueType,pointValueType>(shards::getCellTopologyData<shards::Line<2> >()));
266  else if(subCellDim == 2) {
267  if((subCellOrd >=0) && (subCellOrd<3))
268  return Teuchos::rcp(new
270  if((subCellOrd==3)||(subCellOrd==4))
271  return Teuchos::rcp(new
273  }
274 
275  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
276  }
277 
278  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
279  getHostBasis() const override{
281  }
282  };
283 }// namespace Intrepid2
284 
286 
287 #endif
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Wedge cell.
ordinal_type getCardinality() const
Returns cardinality of the basis.
See Intrepid2::Basis_HCURL_WEDGE_I1_FEM.
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 bool requireOrientation() const override
True if orientation is required.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Quadrilateral cell...
virtual const char * getName() const override
Returns basis name.
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 HVOL-compatible FEM contstant basis on triangle, quadrilateral, hexahedron and tetrahedron cells.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Triangle cell...
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
Header file for the Intrepid2::Basis_HCURL_QUAD_I1_FEM class.
Definition file for FEM basis functions of degree 1 for H(curl) functions on WEDGE cells...
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Header file for the Intrepid2::Basis_HCURL_TRI_I1_FEM class.
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.