Intrepid2
Intrepid2_HGRAD_WEDGE_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_WEDGE_C1_FEM_HPP__
17 #define __INTREPID2_HGRAD_WEDGE_C1_FEM_HPP__
18 
19 #include "Intrepid2_Basis.hpp"
22 
23 namespace Intrepid2 {
24 
55  namespace Impl {
56 
61  public:
62  typedef struct Wedge<6> cell_topology_type;
66  template<EOperator opType>
67  struct Serial {
68  template<typename OutputViewType,
69  typename inputViewType>
70  KOKKOS_INLINE_FUNCTION
71  static void
72  getValues( OutputViewType output,
73  const inputViewType input );
74 
75  };
76 
77  template<typename DeviceType,
78  typename outputValueValueType, class ...outputValueProperties,
79  typename inputPointValueType, class ...inputPointProperties>
80  static void
81  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
82  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
83  const EOperator operatorType);
84 
88  template<typename outputValueViewType,
89  typename inputPointViewType,
90  EOperator opType>
91  struct Functor {
92  outputValueViewType _outputValues;
93  const inputPointViewType _inputPoints;
94 
95  KOKKOS_INLINE_FUNCTION
96  Functor( outputValueViewType outputValues_,
97  inputPointViewType inputPoints_ )
98  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
99 
100  KOKKOS_INLINE_FUNCTION
101  void operator()(const ordinal_type pt) const {
102  switch (opType) {
103  case OPERATOR_VALUE : {
104  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
105  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
106  Serial<opType>::getValues( output, input );
107  break;
108  }
109  case OPERATOR_GRAD :
110  case OPERATOR_D2 :
111  case OPERATOR_MAX : {
112  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
113  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
114  Serial<opType>::getValues( output, input );
115  break;
116  }
117  default: {
118  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
119  opType != OPERATOR_GRAD &&
120  opType != OPERATOR_D2 &&
121  opType != OPERATOR_MAX,
122  ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::Serial::getValues) operator is not supported");
123  }
124  }
125  }
126  };
127 
128 
129  };
130  }
131  template<typename DeviceType = void,
132  typename outputValueType = double,
133  typename pointValueType = double>
134  class Basis_HGRAD_WEDGE_C1_FEM : public Basis<DeviceType,outputValueType,pointValueType> {
135  public:
136  using OrdinalTypeArray1DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray1DHost;
137  using OrdinalTypeArray2DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray2DHost;
138  using OrdinalTypeArray3DHost = typename Basis<DeviceType,outputValueType,pointValueType>::OrdinalTypeArray3DHost;
139 
143 
147 
149 
150  virtual
151  void
152  getValues( 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_WEDGE_C1_FEM::
164  getValues<DeviceType>( outputValues,
165  inputPoints,
166  operatorType );
167  }
168 
169  virtual void
170  getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
171  ordinal_type& perThreadSpaceSize,
172  const PointViewType inputPointsconst,
173  const EOperator operatorType = OPERATOR_VALUE) const override;
174 
175  KOKKOS_INLINE_FUNCTION
176  virtual void
177  getValues(
178  OutputViewType outputValues,
179  const PointViewType inputPoints,
180  const EOperator operatorType,
181  const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
182  const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
183  const ordinal_type subcellDim = -1,
184  const ordinal_type subcellOrdinal = -1) const override;
185 
186  virtual
187  void
188  getDofCoords( ScalarViewType dofCoords ) const override {
189 #ifdef HAVE_INTREPID2_DEBUG
190  // Verify rank of output array.
191  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
192  ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
193  // Verify 0th dimension of output array.
194  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
195  ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
196  // Verify 1st dimension of output array.
197  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
198  ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
199 #endif
200  Kokkos::deep_copy(dofCoords, this->dofCoords_);
201  }
202 
203  virtual
204  void
205  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
206 #ifdef HAVE_INTREPID2_DEBUG
207  // Verify rank of output array.
208  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
209  ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
210  // Verify 0th dimension of output array.
211  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
212  ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
213 #endif
214  Kokkos::deep_copy(dofCoeffs, 1.0);
215  }
216 
217  virtual
218  const char*
219  getName() const override {
220  return "Intrepid2_HGRAD_WEDGE_C1_FEM";
221  }
222 
231  BasisPtr<DeviceType,outputValueType,pointValueType>
232  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
233  if(subCellDim == 1) {
234  return Teuchos::rcp(new
236  } else if(subCellDim == 2) {
237  if(subCellOrd < 3) //lateral faces
239  else
241  }
242  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
243  }
244 
245  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
246  getHostBasis() const override{
248  }
249  };
250 }// namespace Intrepid2
251 
253 
254 #endif
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Quadrilateral cell...
Header file for the Intrepid2::Basis_HGRAD_TRI_C1_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Wedge cell.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C1_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Triangle cell...
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< 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...
virtual const char * getName() const override
Returns basis name.
Definition file for FEM basis functions of degree 1 for H(grad) functions on WEDGE cells...
See Intrepid2::Basis_HGRAD_WEDGE_C1_FEM.
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.
Header file for the abstract base class Intrepid2::Basis.