Intrepid2
Intrepid2_HDIV_TET_In_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_TET_IN_FEM_HPP__
17#define __INTREPID2_HDIV_TET_IN_FEM_HPP__
18
19#include "Intrepid2_Basis.hpp"
22
24#include "Teuchos_LAPACK.hpp"
25
26namespace Intrepid2 {
27
54#define CardinalityHDivTet(order) (order*(order+1)*(order+3)/2)
55
56namespace Impl {
57
62public:
63
67 template<EOperator opType>
68 struct Serial {
69 template<typename outputValueViewType,
70 typename inputPointViewType,
71 typename workViewType,
72 typename vinvViewType>
73 KOKKOS_INLINE_FUNCTION
74 static void
75 getValues( outputValueViewType outputValues,
76 const inputPointViewType inputPoints,
77 workViewType work,
78 const vinvViewType vinv );
79
80
81 KOKKOS_INLINE_FUNCTION
82 static ordinal_type
83 getWorkSizePerPoint(ordinal_type order) {
84 auto cardinality = CardinalityHDivTet(order);
85 switch (opType) {
86 case OPERATOR_DIV:
87 case OPERATOR_D1:
88 return 7*cardinality;
89 default:
90 return getDkCardinality<opType,3>()*cardinality;
91 }
92 }
93 };
94
95 template<typename DeviceType, ordinal_type numPtsPerEval,
96 typename outputValueValueType, class ...outputValueProperties,
97 typename inputPointValueType, class ...inputPointProperties,
98 typename vinvValueType, class ...vinvProperties>
99 static void
100 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
101 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
102 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
103 const EOperator operatorType);
104
108 template<typename outputValueViewType,
109 typename inputPointViewType,
110 typename vinvViewType,
111 typename workViewType,
112 EOperator opType,
113 ordinal_type numPtsEval>
114 struct Functor {
115 outputValueViewType _outputValues;
116 const inputPointViewType _inputPoints;
117 const vinvViewType _coeffs;
118 workViewType _work;
119
120 KOKKOS_INLINE_FUNCTION
121 Functor( outputValueViewType outputValues_,
122 inputPointViewType inputPoints_,
123 vinvViewType coeffs_,
124 workViewType work_)
125 : _outputValues(outputValues_), _inputPoints(inputPoints_),
126 _coeffs(coeffs_), _work(work_) {}
127
128 KOKKOS_INLINE_FUNCTION
129 void operator()(const size_type iter) const {
130 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
131 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
132
133 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
134 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
135
136 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
137
138 auto vcprop = Kokkos::common_view_alloc_prop(_work);
139 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
140
141 switch (opType) {
142 case OPERATOR_VALUE : {
143 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
144 Serial<opType>::getValues( output, input, work, _coeffs );
145 break;
146 }
147 case OPERATOR_DIV: {
148 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
149 Serial<opType>::getValues( output, input, work, _coeffs );
150 break;
151 }
152 default: {
153 INTREPID2_TEST_FOR_ABORT( true,
154 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::Functor) operator is not supported");
155
156 }
157 }
158 }
159 };
160};
161}
162
163template<typename DeviceType = void,
164 typename outputValueType = double,
165 typename pointValueType = double>
167 : public Basis<DeviceType,outputValueType,pointValueType> {
168 public:
170 using OrdinalTypeArray1DHost = typename BasisBase::OrdinalTypeArray1DHost;
171 using OrdinalTypeArray2DHost = typename BasisBase::OrdinalTypeArray2DHost;
172 using OrdinalTypeArray3DHost = typename BasisBase::OrdinalTypeArray3DHost;
173
176 Basis_HDIV_TET_In_FEM(const ordinal_type order,
177 const EPointType pointType = POINTTYPE_EQUISPACED);
178
179 using OutputViewType = typename BasisBase::OutputViewType;
180 using PointViewType = typename BasisBase::PointViewType;
181 using ScalarViewType = typename BasisBase::ScalarViewType;
182 using scalarType = typename BasisBase::scalarType;
184
185 virtual
186 void
187 getValues( /* */ OutputViewType outputValues,
188 const PointViewType inputPoints,
189 const EOperator operatorType = OPERATOR_VALUE) const override {
190#ifdef HAVE_INTREPID2_DEBUG
192 inputPoints,
193 operatorType,
194 this->getBaseCellTopology(),
195 this->getCardinality() );
196#endif
197 constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
198 Impl::Basis_HDIV_TET_In_FEM::
199 getValues<DeviceType,numPtsPerEval>( outputValues,
200 inputPoints,
201 this->coeffs_,
202 operatorType);
203 }
204
205 virtual void
206 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
207 ordinal_type& perThreadSpaceSize,
208 const PointViewType inputPointsconst,
209 const EOperator operatorType = OPERATOR_VALUE) const override;
210
211 KOKKOS_INLINE_FUNCTION
212 virtual void
213 getValues(
214 OutputViewType outputValues,
215 const PointViewType inputPoints,
216 const EOperator operatorType,
217 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
218 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
219 const ordinal_type subcellDim = -1,
220 const ordinal_type subcellOrdinal = -1) const override;
221
222 virtual
223 void
224 getDofCoords( ScalarViewType dofCoords ) const override {
225#ifdef HAVE_INTREPID2_DEBUG
226 // Verify rank of output array.
227 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
228 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
229 // Verify 0th dimension of output array.
230 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
231 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
232 // Verify 1st dimension of output array.
233 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
234 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
235#endif
236 Kokkos::deep_copy(dofCoords, this->dofCoords_);
237 }
238
239 virtual
240 void
241 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
242 #ifdef HAVE_INTREPID2_DEBUG
243 // Verify rank of output array.
244 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
245 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
246 // Verify 0th dimension of output array.
247 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
248 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
249 // Verify 1st dimension of output array.
250 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
251 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
252#endif
253 Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
254 }
255
256 void
257 getExpansionCoeffs( ScalarViewType coeffs ) const {
258 // has to be same rank and dimensions
259 Kokkos::deep_copy(coeffs, this->coeffs_);
260 }
261
262 virtual
263 const char*
264 getName() const override {
265 return "Intrepid2_HDIV_TET_In_FEM";
266 }
267
268 virtual
269 bool
270 requireOrientation() const override {
271 return true;
272 }
273
284 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
285
286 if(subCellDim == 2) {
287 return Teuchos::rcp(new
289 (this->basisDegree_-1, pointType_));
290 }
291 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
292 }
293
298 private:
299
301 Kokkos::DynRankView<scalarType,DeviceType> coeffs_;
302
304 EPointType pointType_;
305 };
306
307}// namespace Intrepid2
308
310
311#endif
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
void getValues_HDIV_Args(const outputValueViewType outputValues, const inputPointViewType inputPoints, const EOperator operatorType, const shards::CellTopology cellTopo, const ordinal_type basisCard)
Runtime check of the arguments for the getValues method in an HDIV-conforming FEM basis....
Definition file for FEM basis functions of degree n for H(grad) functions on TET cells.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
Header file for the Intrepid2::Basis_HVOL_TRI_Cn_FEM class.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates,...
Implementation of the default H(div)-compatible Raviart-Thomas basis of arbitrary degree on Tetrahedr...
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 bool requireOrientation() const override
True if orientation is required.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual const char * getName() const override
Returns basis name.
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...
EPointType pointType_
type of lattice used for creating the DoF coordinates
Kokkos::DynRankView< scalarType, DeviceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
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 Lagrange basis of arbitrary degree on Triangle cell.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
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::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
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.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HDIV_TET_In_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
small utility functions