Intrepid2
Intrepid2_HVOL_QUAD_Cn_FEMDef.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_QUAD_CN_FEM_DEF_HPP__
16#define __INTREPID2_HVOL_QUAD_CN_FEM_DEF_HPP__
17
18namespace Intrepid2 {
19
20 // -------------------------------------------------------------------------------------
21 namespace Impl {
22
23 template<EOperator OpType>
24 template<typename OutputViewType,
25 typename InputViewType,
26 typename WorkViewType,
27 typename VinvViewType>
28 KOKKOS_INLINE_FUNCTION
29 void
30 Basis_HVOL_QUAD_Cn_FEM::Serial<OpType>::
31 getValues( OutputViewType output,
32 const InputViewType input,
33 WorkViewType work,
34 const VinvViewType vinv,
35 const ordinal_type operatorDn ) {
36 ordinal_type opDn = operatorDn;
37
38 const ordinal_type cardLine = vinv.extent(0);
39 const ordinal_type npts = input.extent(0);
40
41 typedef Kokkos::pair<ordinal_type,ordinal_type> range_type;
42 const auto input_x = Kokkos::subview(input, Kokkos::ALL(), range_type(0,1));
43 const auto input_y = Kokkos::subview(input, Kokkos::ALL(), range_type(1,2));
44
45 const ordinal_type dim_s = get_dimension_scalar(input);
46 auto ptr0 = work.data();
47 auto ptr1 = work.data()+cardLine*npts*dim_s;
48 auto ptr2 = work.data()+2*cardLine*npts*dim_s;
49
50 typedef typename Kokkos::DynRankView<typename InputViewType::value_type, typename WorkViewType::memory_space> ViewType;
51 auto vcprop = Kokkos::common_view_alloc_prop(input);
52
53 switch (OpType) {
54 case OPERATOR_VALUE: {
55 ViewType work_line(Kokkos::view_wrap(ptr0, vcprop), cardLine, npts);
56 ViewType output_x(Kokkos::view_wrap(ptr1, vcprop), cardLine, npts);
57 ViewType output_y(Kokkos::view_wrap(ptr2, vcprop), cardLine, npts);
58
59 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
60 getValues(output_x, input_x, work_line, vinv);
61
62 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
63 getValues(output_y, input_y, work_line, vinv);
64
65 // tensor product
66 ordinal_type idx = 0;
67 for (ordinal_type j=0;j<cardLine;++j) // y
68 for (ordinal_type i=0;i<cardLine;++i,++idx) // x
69 for (ordinal_type k=0;k<npts;++k)
70 output.access(idx,k) = output_x.access(i,k)*output_y.access(j,k);
71 break;
72 }
73 case OPERATOR_GRAD:
74 case OPERATOR_D1:
75 case OPERATOR_D2:
76 case OPERATOR_D3:
77 case OPERATOR_D4:
78 case OPERATOR_D5:
79 case OPERATOR_D6:
80 case OPERATOR_D7:
81 case OPERATOR_D8:
82 case OPERATOR_D9:
83 case OPERATOR_D10:
84 opDn = getOperatorOrder(OpType);
85 case OPERATOR_Dn: {
86 const auto dkcard = opDn + 1;
87 for (auto l=0;l<dkcard;++l) {
88 ViewType work_line(Kokkos::view_wrap(ptr0, vcprop), cardLine, npts);
89
90 ViewType output_x, output_y;
91
92 const auto mult_x = opDn - l;
93 const auto mult_y = l;
94
95 if (mult_x) {
96 output_x = ViewType(Kokkos::view_wrap(ptr1, vcprop), cardLine, npts, 1);
97 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_Dn>::
98 getValues(output_x, input_x, work_line, vinv, mult_x);
99 } else {
100 output_x = ViewType(Kokkos::view_wrap(ptr1, vcprop), cardLine, npts);
101 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
102 getValues(output_x, input_x, work_line, vinv);
103 }
104
105 if (mult_y) {
106 output_y = ViewType(Kokkos::view_wrap(ptr2, vcprop), cardLine, npts, 1);
107 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_Dn>::
108 getValues(output_y, input_y, work_line, vinv, mult_y);
109 } else {
110 output_y = ViewType(Kokkos::view_wrap(ptr2, vcprop), cardLine, npts);
111 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
112 getValues(output_y, input_y, work_line, vinv);
113 }
114
115 // tensor product (extra dimension of ouput x and y are ignored)
116 ordinal_type idx = 0;
117 for (ordinal_type j=0;j<cardLine;++j) // y
118 for (ordinal_type i=0;i<cardLine;++i,++idx) // x
119 for (ordinal_type k=0;k<npts;++k)
120 output.access(idx,k,l) = output_x.access(i,k,0)*output_y.access(j,k,0);
121 }
122 break;
123 }
124 default: {
125 INTREPID2_TEST_FOR_ABORT( true,
126 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::Serial::getValues) operator is not supported" );
127 }
128 }
129 }
130
131 template<typename DT, ordinal_type numPtsPerEval,
132 typename outputValueValueType, class ...outputValueProperties,
133 typename inputPointValueType, class ...inputPointProperties,
134 typename vinvValueType, class ...vinvProperties>
135 void
136 Basis_HVOL_QUAD_Cn_FEM::
137 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
138 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
139 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
140 const EOperator operatorType ) {
141 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
142 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
143 typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
144 typedef typename ExecSpace<typename inputPointViewType::execution_space,typename DT::execution_space>::ExecSpaceType ExecSpaceType;
145
146 // loopSize corresponds to cardinality
147 const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
148 const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
149 const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
150 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
151
152 typedef typename inputPointViewType::value_type inputPointType;
153
154 const ordinal_type cardinality = outputValues.extent(0);
155 const ordinal_type cardLine = std::sqrt(cardinality);
156 const ordinal_type workSize = 3*cardLine;
157
158 auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
159 typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
160 workViewType work(Kokkos::view_alloc("Basis_HVOL_QUAD_Cn_FEM::getValues::work", vcprop), workSize, inputPoints.extent(0));
161
162 switch (operatorType) {
163 case OPERATOR_VALUE: {
164 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
165 OPERATOR_VALUE,numPtsPerEval> FunctorType;
166 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
167 break;
168 }
169 case OPERATOR_GRAD:
170 case OPERATOR_D1:
171 case OPERATOR_D2:
172 case OPERATOR_D3:
173 case OPERATOR_D4:
174 case OPERATOR_D5:
175 case OPERATOR_D6:
176 case OPERATOR_D7:
177 case OPERATOR_D8:
178 case OPERATOR_D9:
179 case OPERATOR_D10: {
180 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
181 OPERATOR_Dn,numPtsPerEval> FunctorType;
182 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work,
183 getOperatorOrder(operatorType)) );
184 break;
185 }
186 default: {
187 INTREPID2_TEST_FOR_EXCEPTION( true , std::invalid_argument,
188 ">>> ERROR (Basis_HVOL_QUAD_Cn_FEM): Operator type not implemented" );
189 // break;commented out because exception
190 }
191 }
192 }
193 }
194
195 // -------------------------------------------------------------------------------------
196 template<typename DT, typename OT, typename PT>
198 Basis_HVOL_QUAD_Cn_FEM( const ordinal_type order,
199 const EPointType pointType ) {
200 // INTREPID2_TEST_FOR_EXCEPTION( !(pointType == POINTTYPE_EQUISPACED ||
201 // pointType == POINTTYPE_WARPBLEND), std::invalid_argument,
202 // ">>> ERROR (Basis_HVOL_QUAD_Cn_FEM): pointType must be either equispaced or warpblend." );
203
204 // this should be in host
205 Basis_HVOL_LINE_Cn_FEM<DT,OT,PT> lineBasis( order, pointType );
206 const auto cardLine = lineBasis.getCardinality();
207
208 this->pointType_ = pointType;
209 this->vinv_ = Kokkos::DynRankView<typename ScalarViewType::value_type,DT>("HVOL::Quad::Cn::vinv", cardLine, cardLine);
210 lineBasis.getVandermondeInverse(this->vinv_);
211
212 const ordinal_type spaceDim = 2;
213 this->basisCardinality_ = cardLine*cardLine;
214 this->basisDegree_ = order;
215 this->basisCellTopologyKey_ = shards::Quadrilateral<4>::key;
216 this->basisType_ = BASIS_FEM_LAGRANGIAN;
217 this->basisCoordinates_ = COORDINATES_CARTESIAN;
218 this->functionSpace_ = FUNCTION_SPACE_HVOL;
219
220 // initialize tags
221 {
222 // Basis-dependent initializations
223 const ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
224 const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
225 const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
226 const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
227
228 // An array with local DoF tags assigned to the basis functions, in the order of their local enumeration
229 constexpr ordinal_type maxCardLine = Parameters::MaxOrder + 1;
230 ordinal_type tags[maxCardLine*maxCardLine][4];
231
232 {
233 ordinal_type idx = 0;
234 for (ordinal_type j=0;j<cardLine;++j) { // y
235 const auto tag_y = lineBasis.getDofTag(j);
236 for (ordinal_type i=0;i<cardLine;++i,++idx) { // x
237 const auto tag_x = lineBasis.getDofTag(i);
238
239 // interior
240 tags[idx][0] = 2; // interior dof
241 tags[idx][1] = 0;
242 tags[idx][2] = tag_x(2) + tag_x(3)*tag_y(2); // local dof id
243 tags[idx][3] = tag_x(3)*tag_y(3); // total number of dofs in this vertex
244 }
245 }
246 }
247
248 OrdinalTypeArray1DHost tagView(&tags[0][0], this->basisCardinality_*4);
249
250 // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
251 // tags are constructed on host
252 this->setOrdinalTagData(this->tagToOrdinal_,
253 this->ordinalToTag_,
254 tagView,
255 this->basisCardinality_,
256 tagSize,
257 posScDim,
258 posScOrd,
259 posDfOrd);
260 }
261
262 // dofCoords on host and create its mirror view to device
263 Kokkos::DynRankView<typename ScalarViewType::value_type,typename DT::execution_space::array_layout,Kokkos::HostSpace>
264 dofCoordsHost("dofCoordsHost", this->basisCardinality_, spaceDim);
265
266 Kokkos::DynRankView<typename ScalarViewType::value_type,DT>
267 dofCoordsLine("dofCoordsLine", cardLine, 1);
268
269 lineBasis.getDofCoords(dofCoordsLine);
270 auto dofCoordsLineHost = Kokkos::create_mirror_view(dofCoordsLine);
271 Kokkos::deep_copy(dofCoordsLineHost, dofCoordsLine);
272 {
273 ordinal_type idx = 0;
274 for (ordinal_type j=0;j<cardLine;++j) { // y
275 for (ordinal_type i=0;i<cardLine;++i,++idx) { // x
276 dofCoordsHost(idx,0) = dofCoordsLineHost(i,0);
277 dofCoordsHost(idx,1) = dofCoordsLineHost(j,0);
278 }
279 }
280 }
281
282 this->dofCoords_ = Kokkos::create_mirror_view(typename DT::memory_space(), dofCoordsHost);
283 Kokkos::deep_copy(this->dofCoords_, dofCoordsHost);
284 }
285
286 template<typename DT, typename OT, typename PT>
287 void
289 ordinal_type& perTeamSpaceSize,
290 ordinal_type& perThreadSpaceSize,
291 const PointViewType inputPoints,
292 const EOperator operatorType) const {
293 perTeamSpaceSize = 0;
294 perThreadSpaceSize = 3*this->vinv_.extent(0)*get_dimension_scalar(inputPoints)*sizeof(typename BasisBase::scalarType);
295 }
296
297 template<typename DT, typename OT, typename PT>
298 KOKKOS_INLINE_FUNCTION
299 void
301 OutputViewType outputValues,
302 const PointViewType inputPoints,
303 const EOperator operatorType,
304 const typename Kokkos::TeamPolicy<typename DT::execution_space>::member_type& team_member,
305 const typename DT::execution_space::scratch_memory_space & scratchStorage,
306 const ordinal_type subcellDim,
307 const ordinal_type subcellOrdinal) const {
308
309 INTREPID2_TEST_FOR_ABORT( !((subcellDim == -1) && (subcellOrdinal == -1)),
310 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::getValues), The capability of selecting subsets of basis functions has not been implemented yet.");
311
312 const int numPoints = inputPoints.extent(0);
313 using ScalarType = typename ScalarTraits<typename PointViewType::value_type>::scalar_type;
314 using WorkViewType = Kokkos::DynRankView< ScalarType,typename DT::execution_space::scratch_memory_space,Kokkos::MemoryTraits<Kokkos::Unmanaged> >;
315 auto sizePerPoint = 3*this->vinv_.extent(0)*get_dimension_scalar(inputPoints);
316 WorkViewType workView(scratchStorage, sizePerPoint*team_member.team_size());
317 using range_type = Kokkos::pair<ordinal_type,ordinal_type>;
318 switch(operatorType) {
319 case OPERATOR_VALUE:
320 Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=, &vinv_ = this->vinv_, basisDegree_ = this->basisDegree_] (ordinal_type& pt) {
321 auto output = Kokkos::subview( outputValues, Kokkos::ALL(), range_type (pt,pt+1), Kokkos::ALL() );
322 const auto input = Kokkos::subview( inputPoints, range_type(pt, pt+1), Kokkos::ALL() );
323 WorkViewType work(workView.data() + sizePerPoint*team_member.team_rank(), sizePerPoint);
324 Impl::Basis_HVOL_QUAD_Cn_FEM::Serial<OPERATOR_VALUE>::getValues( output, input, work, vinv_, basisDegree_);
325 });
326 break;
327 default: {
328 INTREPID2_TEST_FOR_ABORT( true,
329 ">>> ERROR (Basis_HVOL_QUAD_Cn_FEM): getValues not implemented for this operator");
330 }
331 }
332 }
333
334} // namespace Intrepid2
335
336#endif
KOKKOS_INLINE_FUNCTION ordinal_type getOperatorOrder(const EOperator operatorType)
Returns order of an operator.
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,...
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Implementation of the default HVOL-compatible FEM basis of degree n on Quadrilateral cell Implements ...
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...
Basis_HVOL_QUAD_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
const OrdinalTypeArrayStride1DHost getDofTag(const ordinal_type dofOrd) const
DoF ordinal to DoF tag lookup.
ordinal_type getCardinality() const
Returns cardinality of the basis.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.