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