Intrepid2
Intrepid2_HGRAD_LINE_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
16#ifndef __INTREPID2_HGRAD_LINE_CN_FEM_HPP__
17#define __INTREPID2_HGRAD_LINE_CN_FEM_HPP__
18
19#include "Intrepid2_Basis.hpp"
21
23#include "Teuchos_LAPACK.hpp"
24
25namespace Intrepid2 {
26
43 namespace Impl {
44
49 public:
50 typedef struct Line<2> cell_topology_type;
54 template<EOperator opType>
55 struct Serial {
56 template<typename outputValueViewType,
57 typename inputPointViewType,
58 typename workViewType,
59 typename vinvViewType>
60 KOKKOS_INLINE_FUNCTION
61 static void
62 getValues( outputValueViewType outputValues,
63 const inputPointViewType inputPoints,
64 workViewType work,
65 const vinvViewType vinv,
66 const ordinal_type operatorDn = 0 );
67 };
68
69 template<typename DeviceType, ordinal_type numPtsPerEval,
70 typename outputValueValueType, class ...outputValueProperties,
71 typename inputPointValueType, class ...inputPointProperties,
72 typename vinvValueType, class ...vinvProperties>
73 static void
74 getValues( const typename DeviceType::execution_space& space,
75 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
76 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
77 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
78 const EOperator operatorType );
79
83 template<typename outputValueViewType,
84 typename inputPointViewType,
85 typename vinvViewType,
86 typename workViewType,
87 EOperator opType,
88 ordinal_type numPtsEval>
89 struct Functor {
90 outputValueViewType _outputValues;
91 const inputPointViewType _inputPoints;
92 const vinvViewType _vinv;
93 workViewType _work;
94 const ordinal_type _opDn;
95
96 KOKKOS_INLINE_FUNCTION
97 Functor( outputValueViewType outputValues_,
98 inputPointViewType inputPoints_,
99 vinvViewType vinv_,
100 workViewType work_,
101 const ordinal_type opDn_ = 0 )
102 : _outputValues(outputValues_), _inputPoints(inputPoints_),
103 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
104
105 KOKKOS_INLINE_FUNCTION
106 void operator()(const size_type iter) const {
107 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
108 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
109
110 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
111 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
112
113 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
114
115 auto vcprop = Kokkos::common_view_alloc_prop(_work);
116 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
117
118 switch (opType) {
119 case OPERATOR_VALUE : {
120 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
121 Serial<opType>::getValues( output, input, work, _vinv );
122 break;
123 }
124 case OPERATOR_Dn : {
125 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
126 Serial<opType>::getValues( output, input, work, _vinv, _opDn );
127 break;
128 }
129 default: {
130 INTREPID2_TEST_FOR_ABORT( true,
131 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::Functor) operator is not supported");
132
133 }
134 }
135 }
136 };
137 };
138 }
139
140 template<typename DeviceType = void,
141 typename outputValueType = double,
142 typename pointValueType = double>
144 : public Basis<DeviceType,outputValueType,pointValueType> {
145 public:
147 using typename BasisBase::ExecutionSpace;
148
150
154
155 using typename BasisBase::OutputViewType;
156 using typename BasisBase::PointViewType ;
157 using typename BasisBase::ScalarViewType;
158
159 private:
160
163 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType> vinv_;
164 EPointType pointType_;
165 public:
168 Basis_HGRAD_LINE_Cn_FEM(const ordinal_type order,
169 const EPointType pointType = POINTTYPE_EQUISPACED);
170
172
173 virtual
174 void
176 OutputViewType outputValues,
177 const PointViewType inputPoints,
178 const EOperator operatorType = OPERATOR_VALUE ) const override {
179#ifdef HAVE_INTREPID2_DEBUG
181 inputPoints,
182 operatorType,
183 this->getBaseCellTopology(),
184 this->getCardinality() );
185#endif
186 constexpr ordinal_type numPtsPerEval = 1;
187 Impl::Basis_HGRAD_LINE_Cn_FEM::
188 getValues<DeviceType,numPtsPerEval>(space,
189 outputValues,
190 inputPoints,
191 this->vinv_,
192 operatorType);
193 }
194
195 virtual void
196 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
197 ordinal_type& perThreadSpaceSize,
198 const PointViewType inputPointsconst,
199 const EOperator operatorType = OPERATOR_VALUE) const override;
200
201 KOKKOS_INLINE_FUNCTION
202 virtual void
203 getValues(
204 OutputViewType outputValues,
205 const PointViewType inputPoints,
206 const EOperator operatorType,
207 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
208 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
209 const ordinal_type subcellDim = -1,
210 const ordinal_type subcellOrdinal = -1) const override;
211
212 virtual
213 void
214 getDofCoords( ScalarViewType dofCoords ) const override {
215#ifdef HAVE_INTREPID2_DEBUG
216 // Verify rank of output array.
217 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
218 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
219 // Verify 0th dimension of output array.
220 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
221 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
222 // Verify 1st dimension of output array.
223 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
224 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
225#endif
226 Kokkos::deep_copy(dofCoords, this->dofCoords_);
227 }
228
229 virtual
230 void
231 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
232#ifdef HAVE_INTREPID2_DEBUG
233 // Verify rank of output array.
234 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
235 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
236 // Verify 0th dimension of output array.
237 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
238 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
239#endif
240 Kokkos::deep_copy(dofCoeffs, 1.0);
241 }
242
243 virtual
244 const char*
245 getName() const override {
246 return "Intrepid2_HGRAD_LINE_Cn_FEM";
247 }
248
249 void
250 getVandermondeInverse( ScalarViewType vinv ) const {
251 // has to be same rank and dimensions
252 Kokkos::deep_copy(vinv, this->vinv_);
253 }
254
255 Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
256 getVandermondeInverse() const {
257 return vinv_;
258 }
259
260 ordinal_type
261 getWorkSizePerPoint(const EOperator operatorType) const {
262 return getPnCardinality<1>(this->basisDegree_);
263 }
264
269 virtual HostBasisPtr<outputValueType,pointValueType>
270 getHostBasis() const override {
271 auto hostBasis = Teuchos::rcp(new HostBasis(this->basisDegree_, pointType_));
272
273 return hostBasis;
274 }
275 };
276
277}// namespace Intrepid2
278
280
281#endif
Header file for the abstract base class Intrepid2::Basis.
void getValues_HGRAD_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 HGRAD-conforming FEM basis....
Definition file for FEM basis functions of degree n for H(grad) functions on LINE.
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM_JACOBI class.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates,...
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1,...
virtual const char * getName() const override
Returns basis name.
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
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.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
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::DynRankView< typename ScalarViewType::value_type, 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 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...
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.
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.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
See Intrepid2::Basis_HGRAD_LINE_Cn_FEM.
small utility functions