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