Intrepid2
Intrepid2_HGRAD_TRI_Cn_FEMDef.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Intrepid2 Package
5 // Copyright (2007) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Kyungjoo Kim (kyukim@sandia.gov), or
38 // Mauro Perego (mperego@sandia.gov)
39 //
40 // ************************************************************************
41 // @HEADER
42 
49 #ifndef __INTREPID2_HGRAD_TRI_CN_FEM_DEF_HPP__
50 #define __INTREPID2_HGRAD_TRI_CN_FEM_DEF_HPP__
51 
53 
54 namespace Intrepid2 {
55 
56 // -------------------------------------------------------------------------------------
57 namespace Impl {
58 
59 template<EOperator opType>
60 template<typename OutputViewType,
61 typename inputViewType,
62 typename workViewType,
63 typename vinvViewType>
64 KOKKOS_INLINE_FUNCTION
65 void
66 Basis_HGRAD_TRI_Cn_FEM::Serial<opType>::
67 getValues( OutputViewType output,
68  const inputViewType input,
69  workViewType work,
70  const vinvViewType vinv ) {
71 
72  constexpr ordinal_type spaceDim = 2;
73  const ordinal_type
74  card = vinv.extent(0),
75  npts = input.extent(0);
76 
77  // compute order
78  ordinal_type order = 0;
79  for (ordinal_type p=0;p<=Parameters::MaxOrder;++p) {
80  if (card == Intrepid2::getPnCardinality<spaceDim>(p) ) {
81  order = p;
82  break;
83  }
84  }
85 
86  typedef typename Kokkos::DynRankView<typename workViewType::value_type, typename workViewType::memory_space> viewType;
87  auto vcprop = Kokkos::common_view_alloc_prop(work);
88  auto ptr = work.data();
89 
90  switch (opType) {
91  case OPERATOR_VALUE: {
92  const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts);
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  output.access(i,j) = 0.0;
101  for (ordinal_type k=0;k<card;++k)
102  output.access(i,j) += vinv(k,i)*phis.access(k,j);
103  }
104  break;
105  }
106  case OPERATOR_GRAD:
107  case OPERATOR_D1: {
108  const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim);
109  ptr += card*npts*spaceDim*get_dimension_scalar(work);
110  const viewType workView(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim+1);
111  Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
112  Serial<opType>::getValues(phis, input, workView, order);
113 
114  for (ordinal_type i=0;i<card;++i)
115  for (ordinal_type j=0;j<npts;++j)
116  for (ordinal_type k=0;k<spaceDim;++k) {
117  output.access(i,j,k) = 0.0;
118  for (ordinal_type l=0;l<card;++l)
119  output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
120  }
121  break;
122  }
123  case OPERATOR_D2:
124  case OPERATOR_D3:
125  case OPERATOR_D4:
126  case OPERATOR_D5:
127  case OPERATOR_D6:
128  case OPERATOR_D7:
129  case OPERATOR_D8:
130  case OPERATOR_D9:
131  case OPERATOR_D10: {
132  const ordinal_type dkcard = getDkCardinality<opType,spaceDim>(); //(orDn + 1);
133  const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, dkcard);
134  viewType dummyView;
135 
136  Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
137  Serial<opType>::getValues(phis, input, dummyView, order);
138 
139  for (ordinal_type i=0;i<card;++i)
140  for (ordinal_type j=0;j<npts;++j)
141  for (ordinal_type k=0;k<dkcard;++k) {
142  output.access(i,j,k) = 0.0;
143  for (ordinal_type l=0;l<card;++l)
144  output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
145  }
146  break;
147  }
148  case OPERATOR_CURL: { // only works in 2d. first component is -d/dy, second is d/dx
149  const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim);
150  ptr += card*npts*spaceDim*get_dimension_scalar(work);
151  const viewType workView(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim+1);
152 
153 
154  Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
155  Serial<OPERATOR_D1>::getValues(phis, input, workView, order);
156 
157  for (ordinal_type i=0;i<card;++i)
158  for (ordinal_type j=0;j<npts;++j) {
159  output.access(i,j,0) = 0.0;
160  for (ordinal_type l=0;l<card;++l)
161  output.access(i,j,0) += vinv(l,i)*phis.access(l,j,1);
162  output.access(i,j,1) = 0.0;
163  for (ordinal_type l=0;l<card;++l)
164  output.access(i,j,1) -= vinv(l,i)*phis.access(l,j,0);
165  }
166  break;
167  }
168  default: {
169  INTREPID2_TEST_FOR_ABORT( true,
170  ">>> ERROR (Basis_HGRAD_TRI_Cn_FEM): Operator type not implemented");
171  }
172  }
173 }
174 
175 template<typename DT, ordinal_type numPtsPerEval,
176 typename outputValueValueType, class ...outputValueProperties,
177 typename inputPointValueType, class ...inputPointProperties,
178 typename vinvValueType, class ...vinvProperties>
179 void
180 Basis_HGRAD_TRI_Cn_FEM::
181 getValues(
182  const typename DT::execution_space& space,
183  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
184  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
185  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
186  const EOperator operatorType) {
187  typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
188  typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
189  typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
190  typedef typename ExecSpace<typename inputPointViewType::execution_space,typename DT::execution_space>::ExecSpaceType ExecSpaceType;
191 
192  // loopSize corresponds to cardinality
193  const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
194  const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
195  const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
196  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(space, 0, loopSize);
197 
198  typedef typename inputPointViewType::value_type inputPointType;
199 
200  const ordinal_type cardinality = outputValues.extent(0);
201  const ordinal_type spaceDim = 2;
202 
203  auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
204  typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
205 
206  switch (operatorType) {
207  case OPERATOR_VALUE: {
208  workViewType work(Kokkos::view_alloc(space, "Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0));
209  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
210  OPERATOR_VALUE,numPtsPerEval> FunctorType;
211  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
212  break;
213  }
214  case OPERATOR_GRAD:
215  case OPERATOR_D1: {
216  workViewType work(Kokkos::view_alloc(space, "Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
217  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
218  OPERATOR_D1,numPtsPerEval> FunctorType;
219  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
220  break;
221  }
222  case OPERATOR_CURL: {
223  workViewType work(Kokkos::view_alloc(space, "Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
224  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
225  OPERATOR_CURL,numPtsPerEval> FunctorType;
226  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
227  break;
228  }
229  case OPERATOR_D2: {
230  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
231  OPERATOR_D2,numPtsPerEval> FunctorType;
232  workViewType work(Kokkos::view_alloc(space, "Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality*outputValues.extent(2), inputPoints.extent(0));
233  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
234  break;
235  }
236  /* case OPERATOR_D3: {
237  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType
238  OPERATOR_D3,numPtsPerEval> FunctorType;
239  workViewType work(Kokkos::view_alloc("Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0), outputValues.extent(2));
240  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
241  break;
242  }*/
243  default: {
244  INTREPID2_TEST_FOR_EXCEPTION( true , std::invalid_argument,
245  ">>> ERROR (Basis_HGRAD_TRI_Cn_FEM): Operator type not implemented" );
246  }
247  }
248 }
249 }
250 
251 // -------------------------------------------------------------------------------------
252 template<typename DT, typename OT, typename PT>
254 Basis_HGRAD_TRI_Cn_FEM( const ordinal_type order,
255  const EPointType pointType ) {
256  constexpr ordinal_type spaceDim = 2;
257 
258  this->basisCardinality_ = Intrepid2::getPnCardinality<spaceDim>(order); // bigN
259  this->basisDegree_ = order; // small n
260  this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Triangle<3> >() );
261  this->basisType_ = BASIS_FEM_LAGRANGIAN;
262  this->basisCoordinates_ = COORDINATES_CARTESIAN;
263  this->functionSpace_ = FUNCTION_SPACE_HGRAD;
264 
265  pointType_ = (pointType == POINTTYPE_DEFAULT) ? POINTTYPE_EQUISPACED : pointType;
266  const ordinal_type card = this->basisCardinality_;
267 
268  // points are computed in the host and will be copied
269  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
270  dofCoords("Hgrad::Tri::Cn::dofCoords", card, spaceDim);
271 
272  // construct lattice
273  const ordinal_type offset = 0;
274  PointTools::getLattice( dofCoords,
275  this->basisCellTopology_,
276  order, offset,
277  pointType_ );
278 
279  this->dofCoords_ = Kokkos::create_mirror_view(typename DT::memory_space(), dofCoords);
280  Kokkos::deep_copy(this->dofCoords_, dofCoords);
281 
282  // form Vandermonde matrix. Actually, this is the transpose of the VDM,
283  // so we transpose on copy below.
284  const ordinal_type lwork = card*card;
285  Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
286  vmat("Hgrad::Tri::Cn::vmat", card, card),
287  work("Hgrad::Tri::Cn::work", lwork),
288  ipiv("Hgrad::Tri::Cn::ipiv", card);
289 
290  Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(typename Kokkos::HostSpace::execution_space{},
291  vmat,
292  dofCoords,
293  order,
294  OPERATOR_VALUE);
295 
296  ordinal_type info = 0;
297  Teuchos::LAPACK<ordinal_type,scalarType> lapack;
298 
299  lapack.GETRF(card, card,
300  vmat.data(), vmat.stride_1(),
301  (ordinal_type*)ipiv.data(),
302  &info);
303 
304  INTREPID2_TEST_FOR_EXCEPTION( info != 0,
305  std::runtime_error ,
306  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM) lapack.GETRF returns nonzero info." );
307 
308  lapack.GETRI(card,
309  vmat.data(), vmat.stride_1(),
310  (ordinal_type*)ipiv.data(),
311  work.data(), lwork,
312  &info);
313 
314  INTREPID2_TEST_FOR_EXCEPTION( info != 0,
315  std::runtime_error ,
316  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM) lapack.GETRI returns nonzero info." );
317 
318  // create host mirror
319  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
320  vinv("Hgrad::Line::Cn::vinv", card, card);
321 
322  for (ordinal_type i=0;i<card;++i)
323  for (ordinal_type j=0;j<card;++j)
324  vinv(i,j) = vmat(j,i);
325 
326  this->vinv_ = Kokkos::create_mirror_view(typename DT::memory_space(), vinv);
327  Kokkos::deep_copy(this->vinv_ , vinv);
328 
329  // initialize tags
330  {
331  // Basis-dependent initializations
332  constexpr ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
333  const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
334  const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
335  const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
336 
337  // 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.
338  // 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.)
339  INTREPID2_TEST_FOR_EXCEPTION( order > Parameters::MaxOrder, std::invalid_argument, "polynomial order exceeds the max supported by this class");
340  constexpr ordinal_type maxCard = Intrepid2::getPnCardinality<spaceDim, Parameters::MaxOrder>();
341  ordinal_type tags[maxCard][tagSize];
342 
343  const ordinal_type
344  numEdgeDof = Intrepid2::getPnCardinality<1>(order-2),
345  numElemDof = (order > 2 ? Intrepid2::getPnCardinality<2>(order-3) : 0);
346 
347  scalarType xi0, xi1, xi2;
348  const scalarType eps = threshold();
349 
350  ordinal_type edgeId[3] = {}, elemId = 0;
351  for (ordinal_type i=0;i<card;++i) {
352 
353  // compute barycentric coordinates
354  const auto x = dofCoords(i,0);
355  const auto y = dofCoords(i,1);
356  xi0 = 1.0 - x - y;
357  xi1= x;
358  xi2= y;
359 
360  // vertex
361  if ((1.0 - xi0) < eps) { // vert 0
362  tags[i][0] = 0; // vertex dof
363  tags[i][1] = 0; // vertex id
364  tags[i][2] = 0; // local dof id
365  tags[i][3] = 1; // total vert dof
366  }
367  else if ((1.0 - xi1) < eps) { // vert 1
368  tags[i][0] = 0; // vertex dof
369  tags[i][1] = 1; // vertex id
370  tags[i][2] = 0; // local dof id
371  tags[i][3] = 1; // total vert dof
372  }
373  else if ((1.0 - xi2) < eps) { // vert 2
374  tags[i][0] = 0; // vertex dof
375  tags[i][1] = 2; // vertex id
376  tags[i][2] = 0; // local dof id
377  tags[i][3] = 1; // total vert dof
378  }
379  else if (xi2 < eps) { // edge 0
380  tags[i][0] = 1; // edge dof
381  tags[i][1] = 0; // edge id
382  tags[i][2] = edgeId[0]++; // local dof id
383  tags[i][3] = numEdgeDof; // total vert dof
384  }
385  else if (xi0 < eps) { // edge 1
386  tags[i][0] = 1; // edge dof
387  tags[i][1] = 1; // edge id
388  tags[i][2] = edgeId[1]++; // local dof id
389  tags[i][3] = numEdgeDof; // total vert dof
390  }
391  else if (xi1 < eps) { // edge 2
392  tags[i][0] = 1; // edge dof
393  tags[i][1] = 2; // edge id
394  tags[i][2] = edgeId[2]++; // local dof id
395  tags[i][3] = numEdgeDof; // total vert dof
396  }
397  else { // elem
398  tags[i][0] = 2; // intr dof
399  tags[i][1] = 0; // intr id
400  tags[i][2] = elemId++; // local dof id
401  tags[i][3] = numElemDof; // total vert dof
402  }
403  }
404 
405  OrdinalTypeArray1DHost tagView(&tags[0][0], card*tagSize);
406 
407  // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
408  // tags are constructed on host
409  this->setOrdinalTagData(this->tagToOrdinal_,
410  this->ordinalToTag_,
411  tagView,
412  this->basisCardinality_,
413  tagSize,
414  posScDim,
415  posScOrd,
416  posDfOrd);
417  }
418 }
419 } // namespace Intrepid2
420 #endif
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
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.
Basis_HGRAD_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.