Intrepid2
Intrepid2_HCURL_TRI_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_TRI_I1_FEM_HPP
17 #define INTREPID2_HCURL_TRI_I1_FEM_HPP
18 
19 #include "Intrepid2_Basis.hpp"
21 
22 namespace Intrepid2 {
23 
64  namespace Impl {
65 
70  public:
71  typedef struct Triangle<3> cell_topology_type;
75  template<EOperator opType>
76  struct Serial {
77  template<typename OutputViewType,
78  typename inputViewType>
79  KOKKOS_INLINE_FUNCTION
80  static void
81  getValues( OutputViewType output,
82  const inputViewType input );
83 
84  };
85 
86  template<typename DeviceType,
87  typename outputValueValueType, class ...outputValueProperties,
88  typename inputPointValueType, class ...inputPointProperties>
89  static void
90  getValues( const typename DeviceType::execution_space& space,
91  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
92  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
93  const EOperator operatorType);
94 
98  template<typename outputValueViewType,
99  typename inputPointViewType,
100  EOperator opType>
101  struct Functor {
102  outputValueViewType _outputValues;
103  const inputPointViewType _inputPoints;
104 
105  KOKKOS_INLINE_FUNCTION
106  Functor( outputValueViewType outputValues_,
107  inputPointViewType inputPoints_ )
108  : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
109 
110  KOKKOS_INLINE_FUNCTION
111  void operator()(const ordinal_type pt) const {
112  switch (opType) {
113  case OPERATOR_VALUE : {
114  auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
115  const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
116  Serial<opType>::getValues( output, input );
117  break;
118  }
119  case OPERATOR_CURL : {
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_DIV: {
126  INTREPID2_TEST_FOR_ABORT( (opType == OPERATOR_DIV),
127  ">>> ERROR (Basis_HCURL_TRI_I1_FEM): DIV is invalid operator for HCURL Basis Functions");
128  break;
129  }
130  case OPERATOR_GRAD: {
131  INTREPID2_TEST_FOR_ABORT( (opType == OPERATOR_GRAD),
132  ">>> ERROR (Basis_HCURL_TRI_I1_FEM): GRAD is invalid operator for HCURL Basis Functions");
133  break;
134  }
135  default: {
136  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
137  opType != OPERATOR_CURL,
138  ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_I1_FEM::Serial::getValues) operator is not supported");
139  }
140  }
141  }
142  };
143 
144  };
145  }
146 
147  template< typename DeviceType = void,
148  typename outputValueType = double,
149  typename pointValueType = double >
150  class Basis_HCURL_TRI_I1_FEM : public Basis<DeviceType, outputValueType, pointValueType> {
151  public:
153  using typename BasisBase::ExecutionSpace;
154 
155  using typename BasisBase::OrdinalTypeArray1DHost;
156  using typename BasisBase::OrdinalTypeArray2DHost;
157  using typename BasisBase::OrdinalTypeArray3DHost;
158 
159  using typename BasisBase::OutputViewType;
160  using typename BasisBase::PointViewType ;
161  using typename BasisBase::ScalarViewType;
162 
166 
167  using BasisBase::getValues;
168 
169  virtual
170  void
171  getValues( const ExecutionSpace& space,
172  OutputViewType outputValues,
173  const PointViewType inputPoints,
174  const EOperator operatorType = OPERATOR_VALUE ) const override {
175 #ifdef HAVE_INTREPID2_DEBUG
176  // Verify arguments
177  Intrepid2::getValues_HCURL_Args(outputValues,
178  inputPoints,
179  operatorType,
180  this->getBaseCellTopology(),
181  this->getCardinality() );
182 #endif
183  Impl::Basis_HCURL_TRI_I1_FEM::
184  getValues<DeviceType>(space,
185  outputValues,
186  inputPoints,
187  operatorType);
188  }
189 
190  virtual
191  void
192  getDofCoords( ScalarViewType dofCoords ) const override {
193 #ifdef HAVE_INTREPID2_DEBUG
194  // Verify rank of output array.
195  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
196  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
197  // Verify 0th dimension of output array.
198  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
199  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
200  // Verify 1st dimension of output array.
201  INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
202  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
203 #endif
204  Kokkos::deep_copy(dofCoords, this->dofCoords_);
205  }
206 
207  virtual
208  void
209  getDofCoeffs( ScalarViewType dofCoeffs ) const override {
210 #ifdef HAVE_INTREPID2_DEBUG
211  // Verify rank of output array.
212  INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 2, std::invalid_argument,
213  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
214  // Verify 0th dimension of output array.
215  INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
216  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
217  // Verify 1st dimension of output array.
218  INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
219  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
220 #endif
221  Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
222  }
223 
224 
225  virtual
226  const char*
227  getName() const override {
228  return "Intrepid2_HCURL_TRI_I1_FEM";
229  }
230 
231  virtual
232  bool
233  requireOrientation() const override {
234  return true;
235  }
236 
246  BasisPtr<DeviceType,outputValueType,pointValueType>
247  getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
248  if(subCellDim == 1)
249  return Teuchos::rcp( new
250  Basis_HVOL_C0_FEM<DeviceType,outputValueType,pointValueType>(shards::getCellTopologyData<shards::Line<2> >()));
251 
252  INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
253  }
254 
255  BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
256  getHostBasis() const override{
258  }
259  };
260 
261 }// namespace Intrepid2
262 
264 
265 #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.
virtual void getValues(const ExecutionSpace &, OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Definition file for default FEM basis functions of degree 1 for H(curl) functions on Triangle cells...
See Intrepid2::Basis_HCURL_TRI_I1_FEM.
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::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.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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< 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 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. ...
virtual const char * getName() const override
Returns basis name.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
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...
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< 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.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the abstract base class Intrepid2::Basis.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.