Intrepid2
Intrepid2_HVOL_TET_Cn_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
15#ifndef __INTREPID2_HVOL_TET_CN_FEM_HPP__
16#define __INTREPID2_HVOL_TET_CN_FEM_HPP__
17
18#include "Intrepid2_Basis.hpp"
19
21#include "Teuchos_LAPACK.hpp"
22
23
24namespace Intrepid2 {
25
41 namespace Impl {
42
47 public:
48 typedef struct Tetrahedron<4> cell_topology_type;
52 template<EOperator opType>
53 struct Serial {
54 template<typename outputValueViewType,
55 typename inputPointViewType,
56 typename workViewType,
57 typename vinvViewType>
58 KOKKOS_INLINE_FUNCTION
59 static void
60 getValues( outputValueViewType outputValues,
61 const inputPointViewType inputPoints,
62 workViewType work,
63 const vinvViewType vinv );
64
65
66 KOKKOS_INLINE_FUNCTION
67 static ordinal_type
68 getWorkSizePerPoint(ordinal_type order) {
69 auto cardinality = getPnCardinality<3>(order);
70 switch (opType) {
71 case OPERATOR_GRAD:
72 case OPERATOR_CURL:
73 case OPERATOR_D1:
74 return 7*cardinality;
75 default:
76 return getDkCardinality<opType,3>()*cardinality;
77 }
78 }
79 };
80
81 template<typename DeviceType, ordinal_type numPtsPerEval,
82 typename outputValueValueType, class ...outputValueProperties,
83 typename inputPointValueType, class ...inputPointProperties,
84 typename vinvValueType, class ...vinvProperties>
85 static void
86 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
87 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
88 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
89 const EOperator operatorType);
90
94 template<typename outputValueViewType,
95 typename inputPointViewType,
96 typename vinvViewType,
97 typename workViewType,
98 EOperator opType,
99 ordinal_type numPtsEval>
100 struct Functor {
101 outputValueViewType _outputValues;
102 const inputPointViewType _inputPoints;
103 const vinvViewType _vinv;
104 workViewType _work;
105
106 KOKKOS_INLINE_FUNCTION
107 Functor( outputValueViewType outputValues_,
108 inputPointViewType inputPoints_,
109 vinvViewType vinv_,
110 workViewType work_)
111 : _outputValues(outputValues_), _inputPoints(inputPoints_),
112 _vinv(vinv_), _work(work_) {}
113
114 KOKKOS_INLINE_FUNCTION
115 void operator()(const size_type iter) const {
116 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
117 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
118
119 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
120 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
121
122 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
123
124 auto vcprop = Kokkos::common_view_alloc_prop(_work);
125 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
126
127 switch (opType) {
128 case OPERATOR_VALUE : {
129 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
130 Serial<opType>::getValues( output, input, work, _vinv );
131 break;
132 }
133 case OPERATOR_GRAD :
134 case OPERATOR_D1 :
135 case OPERATOR_D2 :
136 //case OPERATOR_D3 :
137 {
138 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
139 Serial<opType>::getValues( output, input, work, _vinv );
140 break;
141 }
142 default: {
143 INTREPID2_TEST_FOR_ABORT( true,
144 ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::Functor) operator is not supported");
145
146 }
147 }
148 }
149 };
150 };
151 }
152
153 template<typename DeviceType = void,
154 typename outputValueType = double,
155 typename pointValueType = double>
157 : public Basis<DeviceType,outputValueType,pointValueType> {
158 public:
160
164
165 using typename BasisBase::OutputViewType;
166 using typename BasisBase::PointViewType ;
167 using typename BasisBase::ScalarViewType;
168
171 Basis_HVOL_TET_Cn_FEM(const ordinal_type order,
172 const EPointType pointType = POINTTYPE_EQUISPACED);
173
174 using scalarType = typename BasisBase::scalarType;
176
177 virtual
178 void
179 getValues( OutputViewType outputValues,
180 const PointViewType inputPoints,
181 const EOperator operatorType = OPERATOR_VALUE) const override {
182#ifdef HAVE_INTREPID2_DEBUG
184 inputPoints,
185 operatorType,
186 this->getBaseCellTopology(),
187 this->getCardinality() );
188#endif
189 constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
190 Impl::Basis_HVOL_TET_Cn_FEM::
191 getValues<DeviceType,numPtsPerEval>( outputValues,
192 inputPoints,
193 this->vinv_,
194 operatorType);
195 }
196
197 virtual void
198 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
199 ordinal_type& perThreadSpaceSize,
200 const PointViewType inputPoints,
201 const EOperator operatorType = OPERATOR_VALUE) const override;
202
203 KOKKOS_INLINE_FUNCTION
204 virtual void
205 getValues(
206 OutputViewType outputValues,
207 const PointViewType inputPoints,
208 const EOperator operatorType,
209 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
210 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
211 const ordinal_type subcellDim = -1,
212 const ordinal_type subcellOrdinal = -1) const override;
213
214
215 virtual
216 void
217 getDofCoords( ScalarViewType dofCoords ) const override {
218#ifdef HAVE_INTREPID2_DEBUG
219 // Verify rank of output array.
220 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
221 ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
222 // Verify 0th dimension of output array.
223 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
224 ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
225 // Verify 1st dimension of output array.
226 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
227 ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
228#endif
229 Kokkos::deep_copy(dofCoords, this->dofCoords_);
230 }
231
232 virtual
233 void
234 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
235#ifdef HAVE_INTREPID2_DEBUG
236 // Verify rank of output array.
237 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
238 ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
239 // Verify 0th dimension of output array.
240 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
241 ">>> ERROR: (Intrepid2::Basis_HVOL_TET_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
242#endif
243 Kokkos::deep_copy(dofCoeffs, 1.0);
244 }
245
246 void
247 getVandermondeInverse( ScalarViewType vinv ) const {
248 // has to be same rank and dimensions
249 Kokkos::deep_copy(vinv, this->vinv_);
250 }
251
252 virtual
253 const char*
254 getName() const override {
255 return "Intrepid2_HVOL_TET_Cn_FEM";
256 }
257
258 virtual
259 bool
260 requireOrientation() const override {
261 return false;
262 }
263
268
269 private:
270
273 Kokkos::DynRankView<scalarType,DeviceType> vinv_;
274 EPointType pointType_;
275
276 };
277
278}// namespace Intrepid2
279
281
282#endif
Header file for the abstract base class Intrepid2::Basis.
BasisPtr< typename Kokkos::HostSpace::device_type, OutputType, PointType > HostBasisPtr
Pointer to a Basis whose device type is on the host (Kokkos::HostSpace::device_type),...
void getValues_HVOL_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 HVOL-conforming FEM basis....
Definition file for FEM basis functions of degree n for H(vol) functions on TET.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates,...
Implementation of the default HVOL-compatible Lagrange basis of arbitrary degree on Tetrahedron 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...
virtual bool requireOrientation() const override
True if orientation is required.
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
Kokkos::DynRankView< scalarType, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPoints, 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...
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference 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::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_HVOL_TET_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
small utility functions