Intrepid2
Intrepid2_HGRAD_TRI_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 
16 #ifndef __INTREPID2_HGRAD_TRI_CN_FEM_DEF_HPP__
17 #define __INTREPID2_HGRAD_TRI_CN_FEM_DEF_HPP__
18 
20 
21 namespace Intrepid2 {
22 
23 // -------------------------------------------------------------------------------------
24 namespace Impl {
25 
26 template<EOperator OpType>
27 template<typename OutputViewType,
28 typename InputViewType,
29 typename WorkViewType,
30 typename VinvViewType>
31 KOKKOS_INLINE_FUNCTION
32 void
33 Basis_HGRAD_TRI_Cn_FEM::Serial<OpType>::
34 getValues( OutputViewType output,
35  const InputViewType input,
36  WorkViewType work,
37  const VinvViewType vinv,
38  const ordinal_type order ) {
39 
40  constexpr ordinal_type spaceDim = 2;
41  const ordinal_type
42  card = vinv.extent(0),
43  npts = input.extent(0);
44 
45  typedef typename Kokkos::DynRankView<typename InputViewType::value_type, typename WorkViewType::memory_space> ViewType;
46  auto vcprop = Kokkos::common_view_alloc_prop(input);
47  auto ptr = work.data();
48 
49  switch (OpType) {
50  case OPERATOR_VALUE: {
51  const ViewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts);
52  ViewType dummyView;
53 
54  Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
55  Serial<OpType>::getValues(phis, input, dummyView, order);
56 
57  for (ordinal_type i=0;i<card;++i)
58  for (ordinal_type j=0;j<npts;++j) {
59  output.access(i,j) = 0.0;
60  for (ordinal_type k=0;k<card;++k)
61  output.access(i,j) += vinv(k,i)*phis.access(k,j);
62  }
63  break;
64  }
65  case OPERATOR_GRAD:
66  case OPERATOR_D1: {
67  const ViewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim);
68  ptr += card*npts*spaceDim*get_dimension_scalar(input);
69  const ViewType workView(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim+1);
70  Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
71  Serial<OpType>::getValues(phis, input, workView, order);
72 
73  for (ordinal_type i=0;i<card;++i)
74  for (ordinal_type j=0;j<npts;++j)
75  for (ordinal_type k=0;k<spaceDim;++k) {
76  output.access(i,j,k) = 0.0;
77  for (ordinal_type l=0;l<card;++l)
78  output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
79  }
80  break;
81  }
82  case OPERATOR_D2:
83  case OPERATOR_D3:
84  case OPERATOR_D4:
85  case OPERATOR_D5:
86  case OPERATOR_D6:
87  case OPERATOR_D7:
88  case OPERATOR_D8:
89  case OPERATOR_D9:
90  case OPERATOR_D10: {
91  const ordinal_type dkcard = getDkCardinality<OpType,spaceDim>(); //(orDn + 1);
92  const ViewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, dkcard);
93  ViewType dummyView;
94 
95  Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
96  Serial<OpType>::getValues(phis, input, dummyView, order);
97 
98  for (ordinal_type i=0;i<card;++i)
99  for (ordinal_type j=0;j<npts;++j)
100  for (ordinal_type k=0;k<dkcard;++k) {
101  output.access(i,j,k) = 0.0;
102  for (ordinal_type l=0;l<card;++l)
103  output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
104  }
105  break;
106  }
107  case OPERATOR_CURL: { // only works in 2d. first component is -d/dy, second is d/dx
108  const ViewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim);
109  ptr += card*npts*spaceDim*get_dimension_scalar(input);
110  const ViewType workView(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim+1);
111 
112 
113  Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
114  Serial<OPERATOR_D1>::getValues(phis, input, workView, order);
115 
116  for (ordinal_type i=0;i<card;++i)
117  for (ordinal_type j=0;j<npts;++j) {
118  output.access(i,j,0) = 0.0;
119  for (ordinal_type l=0;l<card;++l)
120  output.access(i,j,0) += vinv(l,i)*phis.access(l,j,1);
121  output.access(i,j,1) = 0.0;
122  for (ordinal_type l=0;l<card;++l)
123  output.access(i,j,1) -= vinv(l,i)*phis.access(l,j,0);
124  }
125  break;
126  }
127  default: {
128  INTREPID2_TEST_FOR_ABORT( true,
129  ">>> ERROR (Basis_HGRAD_TRI_Cn_FEM): Operator type not implemented");
130  }
131  }
132 }
133 
134 template<typename DT, ordinal_type numPtsPerEval,
135 typename outputValueValueType, class ...outputValueProperties,
136 typename inputPointValueType, class ...inputPointProperties,
137 typename vinvValueType, class ...vinvProperties>
138 void
139 Basis_HGRAD_TRI_Cn_FEM::
140 getValues(
141  const typename DT::execution_space& space,
142  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
143  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
144  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
145  const ordinal_type order,
146  const EOperator operatorType) {
147  typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
148  typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
149  typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
150  typedef typename ExecSpace<typename inputPointViewType::execution_space,typename DT::execution_space>::ExecSpaceType ExecSpaceType;
151 
152  // loopSize corresponds to cardinality
153  const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
154  const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
155  const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
156  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(space, 0, loopSize);
157 
158  typedef typename inputPointViewType::value_type inputPointType;
159 
160  const ordinal_type cardinality = outputValues.extent(0);
161  const ordinal_type spaceDim = 2;
162 
163  auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
164  typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
165 
166  switch (operatorType) {
167  case OPERATOR_VALUE: {
168  workViewType work(Kokkos::view_alloc(space, "Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0));
169  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
170  OPERATOR_VALUE,numPtsPerEval> FunctorType;
171  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work, order) );
172  break;
173  }
174  case OPERATOR_GRAD:
175  case OPERATOR_D1: {
176  workViewType work(Kokkos::view_alloc(space, "Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
177  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
178  OPERATOR_D1,numPtsPerEval> FunctorType;
179  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work, order) );
180  break;
181  }
182  case OPERATOR_CURL: {
183  workViewType work(Kokkos::view_alloc(space, "Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
184  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
185  OPERATOR_CURL,numPtsPerEval> FunctorType;
186  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work, order) );
187  break;
188  }
189  case OPERATOR_D2: {
190  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
191  OPERATOR_D2,numPtsPerEval> FunctorType;
192  workViewType work(Kokkos::view_alloc(space, "Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality*outputValues.extent(2), inputPoints.extent(0));
193  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work, order) );
194  break;
195  }
196  default: {
197  INTREPID2_TEST_FOR_EXCEPTION( true , std::invalid_argument,
198  ">>> ERROR (Basis_HGRAD_TRI_Cn_FEM): Operator type not implemented" );
199  }
200  }
201 }
202 }
203 
204 // -------------------------------------------------------------------------------------
205 template<typename DT, typename OT, typename PT>
207 Basis_HGRAD_TRI_Cn_FEM( const ordinal_type order,
208  const EPointType pointType ) {
209  constexpr ordinal_type spaceDim = 2;
210 
211  this->basisCardinality_ = Intrepid2::getPnCardinality<spaceDim>(order); // bigN
212  this->basisDegree_ = order; // small n
213  this->basisCellTopologyKey_ = shards::Triangle<3>::key;
214  this->basisType_ = BASIS_FEM_LAGRANGIAN;
215  this->basisCoordinates_ = COORDINATES_CARTESIAN;
216  this->functionSpace_ = FUNCTION_SPACE_HGRAD;
217 
218  pointType_ = (pointType == POINTTYPE_DEFAULT) ? POINTTYPE_EQUISPACED : pointType;
219  const ordinal_type card = this->basisCardinality_;
220 
221  // points are computed in the host and will be copied
222  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
223  dofCoords("Hgrad::Tri::Cn::dofCoords", card, spaceDim);
224 
225  // construct lattice
226  const shards::CellTopology cellTopo(shards::getCellTopologyData<shards::Triangle<3>>());
227  const ordinal_type offset = 0;
228  PointTools::getLattice( dofCoords,
229  cellTopo,
230  order, offset,
231  this->pointType_ );
232 
233  this->dofCoords_ = Kokkos::create_mirror_view(typename DT::memory_space(), dofCoords);
234  Kokkos::deep_copy(this->dofCoords_, dofCoords);
235 
236  // form Vandermonde matrix. Actually, this is the transpose of the VDM,
237  // so we transpose on copy below.
238  const ordinal_type lwork = card*card;
239  Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
240  vmat("Hgrad::Tri::Cn::vmat", card, card),
241  work("Hgrad::Tri::Cn::work", lwork),
242  ipiv("Hgrad::Tri::Cn::ipiv", card);
243 
244  Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(typename Kokkos::HostSpace::execution_space{},
245  vmat,
246  dofCoords,
247  order,
248  OPERATOR_VALUE);
249 
250  ordinal_type info = 0;
251  Teuchos::LAPACK<ordinal_type,scalarType> lapack;
252 
253  lapack.GETRF(card, card,
254  vmat.data(), vmat.stride_1(),
255  (ordinal_type*)ipiv.data(),
256  &info);
257 
258  INTREPID2_TEST_FOR_EXCEPTION( info != 0,
259  std::runtime_error ,
260  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM) lapack.GETRF returns nonzero info." );
261 
262  lapack.GETRI(card,
263  vmat.data(), vmat.stride_1(),
264  (ordinal_type*)ipiv.data(),
265  work.data(), lwork,
266  &info);
267 
268  INTREPID2_TEST_FOR_EXCEPTION( info != 0,
269  std::runtime_error ,
270  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM) lapack.GETRI returns nonzero info." );
271 
272  // create host mirror
273  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
274  vinv("Hgrad::Line::Cn::vinv", card, card);
275 
276  for (ordinal_type i=0;i<card;++i)
277  for (ordinal_type j=0;j<card;++j)
278  vinv(i,j) = vmat(j,i);
279 
280  this->vinv_ = Kokkos::create_mirror_view(typename DT::memory_space(), vinv);
281  Kokkos::deep_copy(this->vinv_ , vinv);
282 
283  // initialize tags
284  {
285  // Basis-dependent initializations
286  constexpr ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
287  const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
288  const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
289  const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
290 
291  // Note: the only reason why equispaced can't support higher order than Parameters::MaxOrder appears to be the fact that the tags below get stored into a fixed-length array.
292  // TODO: relax the maximum order requirement by setting up tags in a different container, perhaps directly into an OrdinalTypeArray1DHost (tagView, below). (As of this writing (1/25/22), looks like other nodal bases do this in a similar way -- those should be fixed at the same time; maybe search for Parameters::MaxOrder.)
293  INTREPID2_TEST_FOR_EXCEPTION( order > Parameters::MaxOrder, std::invalid_argument, "polynomial order exceeds the max supported by this class");
294  constexpr ordinal_type maxCard = Intrepid2::getPnCardinality<spaceDim, Parameters::MaxOrder>();
295  ordinal_type tags[maxCard][tagSize];
296 
297  const ordinal_type
298  numEdgeDof = Intrepid2::getPnCardinality<1>(order-2),
299  numElemDof = (order > 2 ? Intrepid2::getPnCardinality<2>(order-3) : 0);
300 
301  scalarType xi0, xi1, xi2;
302  const scalarType eps = threshold();
303 
304  ordinal_type edgeId[3] = {}, elemId = 0;
305  for (ordinal_type i=0;i<card;++i) {
306 
307  // compute barycentric coordinates
308  const auto x = dofCoords(i,0);
309  const auto y = dofCoords(i,1);
310  xi0 = 1.0 - x - y;
311  xi1= x;
312  xi2= y;
313 
314  // vertex
315  if ((1.0 - xi0) < eps) { // vert 0
316  tags[i][0] = 0; // vertex dof
317  tags[i][1] = 0; // vertex id
318  tags[i][2] = 0; // local dof id
319  tags[i][3] = 1; // total vert dof
320  }
321  else if ((1.0 - xi1) < eps) { // vert 1
322  tags[i][0] = 0; // vertex dof
323  tags[i][1] = 1; // vertex id
324  tags[i][2] = 0; // local dof id
325  tags[i][3] = 1; // total vert dof
326  }
327  else if ((1.0 - xi2) < eps) { // vert 2
328  tags[i][0] = 0; // vertex dof
329  tags[i][1] = 2; // vertex id
330  tags[i][2] = 0; // local dof id
331  tags[i][3] = 1; // total vert dof
332  }
333  else if (xi2 < eps) { // edge 0
334  tags[i][0] = 1; // edge dof
335  tags[i][1] = 0; // edge id
336  tags[i][2] = edgeId[0]++; // local dof id
337  tags[i][3] = numEdgeDof; // total vert dof
338  }
339  else if (xi0 < eps) { // edge 1
340  tags[i][0] = 1; // edge dof
341  tags[i][1] = 1; // edge id
342  tags[i][2] = edgeId[1]++; // local dof id
343  tags[i][3] = numEdgeDof; // total vert dof
344  }
345  else if (xi1 < eps) { // edge 2
346  tags[i][0] = 1; // edge dof
347  tags[i][1] = 2; // edge id
348  tags[i][2] = edgeId[2]++; // local dof id
349  tags[i][3] = numEdgeDof; // total vert dof
350  }
351  else { // elem
352  tags[i][0] = 2; // intr dof
353  tags[i][1] = 0; // intr id
354  tags[i][2] = elemId++; // local dof id
355  tags[i][3] = numElemDof; // total vert dof
356  }
357  }
358 
359  OrdinalTypeArray1DHost tagView(&tags[0][0], card*tagSize);
360 
361  // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
362  // tags are constructed on host
363  this->setOrdinalTagData(this->tagToOrdinal_,
364  this->ordinalToTag_,
365  tagView,
366  this->basisCardinality_,
367  tagSize,
368  posScDim,
369  posScOrd,
370  posDfOrd);
371  }
372 }
373 
374  template<typename DT, typename OT, typename PT>
375  void
377  ordinal_type& perTeamSpaceSize,
378  ordinal_type& perThreadSpaceSize,
379  const PointViewType inputPoints,
380  const EOperator operatorType) const {
381  perTeamSpaceSize = 0;
382  perThreadSpaceSize = getWorkSizePerPoint(operatorType)*get_dimension_scalar(inputPoints)*sizeof(typename BasisBase::scalarType);
383  }
384 
385  template<typename DT, typename OT, typename PT>
386  KOKKOS_INLINE_FUNCTION
387  void
389  OutputViewType outputValues,
390  const PointViewType inputPoints,
391  const EOperator operatorType,
392  const typename Kokkos::TeamPolicy<typename DT::execution_space>::member_type& team_member,
393  const typename DT::execution_space::scratch_memory_space & scratchStorage,
394  const ordinal_type subcellDim,
395  const ordinal_type subcellOrdinal) const {
396 
397  INTREPID2_TEST_FOR_ABORT( !((subcellDim == -1) && (subcellOrdinal == -1)),
398  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getValues), The capability of selecting subsets of basis functions has not been implemented yet.");
399 
400  const int numPoints = inputPoints.extent(0);
401  using ScalarType = typename ScalarTraits<typename PointViewType::value_type>::scalar_type;
402  using WorkViewType = Kokkos::DynRankView< ScalarType,typename DT::execution_space::scratch_memory_space,Kokkos::MemoryTraits<Kokkos::Unmanaged> >;
403  constexpr ordinal_type spaceDim = 2;
404  auto sizePerPoint = (operatorType==OPERATOR_VALUE) ?
405  this->vinv_.extent(0)*get_dimension_scalar(inputPoints) :
406  (2*spaceDim+1)*this->vinv_.extent(0)*get_dimension_scalar(inputPoints);
407  WorkViewType workView(scratchStorage, sizePerPoint*team_member.team_size());
408  using range_type = Kokkos::pair<ordinal_type,ordinal_type>;
409  switch(operatorType) {
410  case OPERATOR_VALUE:
411  Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=] (ordinal_type& pt) {
412  auto output = Kokkos::subview( outputValues, Kokkos::ALL(), range_type (pt,pt+1), Kokkos::ALL() );
413  const auto input = Kokkos::subview( inputPoints, range_type(pt, pt+1), Kokkos::ALL() );
414  WorkViewType work(workView.data() + sizePerPoint*team_member.team_rank(), sizePerPoint);
415  Impl::Basis_HGRAD_TRI_Cn_FEM::Serial<OPERATOR_VALUE>::getValues( output, input, work, this->vinv_, this->basisDegree_);
416  });
417  break;
418  case OPERATOR_GRAD:
419  Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=] (ordinal_type& pt) {
420  auto output = Kokkos::subview( outputValues, Kokkos::ALL(), range_type(pt,pt+1), Kokkos::ALL() );
421  const auto input = Kokkos::subview( inputPoints, range_type(pt,pt+1), Kokkos::ALL() );
422  WorkViewType work(workView.data() + sizePerPoint*team_member.team_rank(), sizePerPoint);
423  Impl::Basis_HGRAD_TRI_Cn_FEM::Serial<OPERATOR_GRAD>::getValues( output, input, work, this->vinv_, this->basisDegree_);
424  });
425  break;
426  case OPERATOR_CURL:
427  Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=] (ordinal_type& pt) {
428  auto output = Kokkos::subview( outputValues, Kokkos::ALL(), range_type(pt,pt+1), Kokkos::ALL() );
429  const auto input = Kokkos::subview( inputPoints, range_type(pt,pt+1), Kokkos::ALL() );
430  WorkViewType work(workView.data() + sizePerPoint*team_member.team_rank(), sizePerPoint);
431  Impl::Basis_HGRAD_TRI_Cn_FEM::Serial<OPERATOR_CURL>::getValues( output, input, work, this->vinv_, this->basisDegree_);
432  });
433  break;
434  default: {
435  INTREPID2_TEST_FOR_ABORT( true,
436  ">>> ERROR (Basis_HGRAD_TRI_Cn_FEM): getValues not implemented for this operator");
437  }
438  }
439  }
440 
441 } // namespace Intrepid2
442 
443 #endif
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
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.
static void getLattice(Kokkos::DynRankView< pointValueType, pointProperties...> points, const shards::CellTopology cellType, const ordinal_type order, const ordinal_type offset=0, const EPointType pointType=POINTTYPE_EQUISPACED)
Computes a lattice of points of a given order on a reference simplex, quadrilateral or hexahedron (cu...
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
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...
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM work is a rank 1 view having the same value_type of inputPoints...
Basis_HGRAD_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.