Intrepid2
Intrepid2_HDIV_TET_In_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_HDIV_TET_IN_FEM_DEF_HPP__
50 #define __INTREPID2_HDIV_TET_IN_FEM_DEF_HPP__
51 
54 
55 namespace Intrepid2 {
56 
57 // -------------------------------------------------------------------------------------
58 
59 namespace Impl {
60 
61 template<EOperator opType>
62 template<typename OutputViewType,
63 typename inputViewType,
64 typename workViewType,
65 typename vinvViewType>
66 KOKKOS_INLINE_FUNCTION
67 void
68 Basis_HDIV_TET_In_FEM::Serial<opType>::
69 getValues( OutputViewType output,
70  const inputViewType input,
71  workViewType work,
72  const vinvViewType coeffs ) {
73 
74  constexpr ordinal_type spaceDim = 3;
75  const ordinal_type
76  cardPn = coeffs.extent(0)/spaceDim,
77  card = coeffs.extent(1),
78  npts = input.extent(0);
79 
80  // compute order
81  ordinal_type order = 0;
82  for (ordinal_type p=0;p<=Parameters::MaxOrder;++p) {
83  if (card == CardinalityHDivTet(p)) {
84  order = p;
85  break;
86  }
87  }
88 
89  typedef typename Kokkos::DynRankView<typename workViewType::value_type, typename workViewType::memory_space> viewType;
90  auto vcprop = Kokkos::common_view_alloc_prop(work);
91  auto ptr = work.data();
92 
93  switch (opType) {
94  case OPERATOR_VALUE: {
95  const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts);
96  workViewType dummyView;
97 
98  Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::
99  Serial<opType>::getValues(phis, input, dummyView, order);
100 
101  for (ordinal_type i=0;i<card;++i)
102  for (ordinal_type j=0;j<npts;++j)
103  for (ordinal_type d=0;d<spaceDim;++d) {
104  output.access(i,j,d) = 0.0;
105  for (ordinal_type k=0;k<cardPn;++k)
106  output.access(i,j,d) += coeffs(k+d*cardPn,i) * phis.access(k,j);
107  }
108  break;
109  }
110  case OPERATOR_DIV: {
111  const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim);
112  ptr += card*npts*spaceDim*get_dimension_scalar(work);
113  const viewType workView(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim+1);
114 
115  Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::
116  Serial<OPERATOR_GRAD>::getValues(phis, input, workView, order);
117 
118  for (ordinal_type i=0;i<card;++i)
119  for (ordinal_type j=0;j<npts;++j) {
120  output.access(i,j) = 0.0;
121  for (ordinal_type k=0; k<cardPn; ++k)
122  for (ordinal_type d=0; d<spaceDim; ++d)
123  output.access(i,j) += coeffs(k+d*cardPn,i)*phis.access(k,j,d);
124  }
125  break;
126  }
127  default: {
128  INTREPID2_TEST_FOR_ABORT( true,
129  ">>> ERROR (Basis_HDIV_TET_In_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_HDIV_TET_In_FEM::
140 getValues( /* */ Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
141  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
142  const Kokkos::DynRankView<vinvValueType, vinvProperties...> coeffs,
143  const EOperator operatorType) {
144  typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
145  typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
146  typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
147  typedef typename ExecSpace<typename inputPointViewType::execution_space,typename DT::execution_space>::ExecSpaceType ExecSpaceType;
148 
149  // loopSize corresponds to cardinality
150  const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
151  const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
152  const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
153  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
154 
155  typedef typename inputPointViewType::value_type inputPointType;
156 
157  const ordinal_type cardinality = outputValues.extent(0);
158  const ordinal_type spaceDim = 3;
159 
160  auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
161  typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
162 
163  switch (operatorType) {
164  case OPERATOR_VALUE: {
165  workViewType work(Kokkos::view_alloc("Basis_HDIV_TET_In_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0));
166  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
167  OPERATOR_VALUE,numPtsPerEval> FunctorType;
168  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, coeffs, work) );
169  break;
170  }
171  case OPERATOR_DIV: {
172  workViewType work(Kokkos::view_alloc("Basis_HDIV_TET_In_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
173  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
174  OPERATOR_DIV,numPtsPerEval> FunctorType;
175  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, coeffs, work) );
176  break;
177  }
178  default: {
179  INTREPID2_TEST_FOR_EXCEPTION( true , std::invalid_argument,
180  ">>> ERROR (Basis_HDIV_TET_In_FEM): Operator type not implemented" );
181  }
182  }
183 }
184 }
185 
186 // -------------------------------------------------------------------------------------
187 template<typename DT, typename OT, typename PT>
189 Basis_HDIV_TET_In_FEM( const ordinal_type order,
190  const EPointType pointType ) {
191 
192  constexpr ordinal_type spaceDim = 3;
193  this->basisCardinality_ = CardinalityHDivTet(order);
194  this->basisDegree_ = order; // small n
195  this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Tetrahedron<4> >() );
196  this->basisType_ = BASIS_FEM_LAGRANGIAN;
197  this->basisCoordinates_ = COORDINATES_CARTESIAN;
198  this->functionSpace_ = FUNCTION_SPACE_HDIV;
199  pointType_ = pointType;
200 
201  const ordinal_type card = this->basisCardinality_;
202 
203  const ordinal_type cardPn = Intrepid2::getPnCardinality<spaceDim>(order); // dim of (P_{n}) -- smaller space
204  const ordinal_type cardPnm1 = Intrepid2::getPnCardinality<spaceDim>(order-1); // dim of (P_{n-1}) -- smaller space
205  const ordinal_type cardPnm2 = Intrepid2::getPnCardinality<spaceDim>(order-2); // dim of (P_{n-2}) -- smaller space
206  const ordinal_type cardVecPn = spaceDim*cardPn; // dim of (P_{n})^3 -- larger space
207  const ordinal_type cardVecPnm1 = spaceDim*cardPnm1; // dim of (P_{n-1})^3 -- smaller space
208  const ordinal_type dim_PkH = cardPnm1 - cardPnm2;
209 
210  // 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.
211  // 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.)
212  INTREPID2_TEST_FOR_EXCEPTION( order > Parameters::MaxOrder, std::invalid_argument, "polynomial order exceeds the max supported by this class");
213 
214  // Basis-dependent initializations
215  constexpr ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
216  constexpr ordinal_type maxCard = CardinalityHDivTet(Parameters::MaxOrder);
217  ordinal_type tags[maxCard][tagSize];
218 
219  // points are computed in the host and will be copied
220  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
221  dofCoords("Hdiv::Tet::In::dofCoords", card, spaceDim);
222 
223  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
224  dofCoeffs("Hdiv::Tet::In::dofCoeffs", card, spaceDim);
225 
226  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
227  coeffs("Hdiv::Tet::In::coeffs", cardVecPn, card);
228 
229  // first, need to project the basis for RT space onto the
230  // orthogonal basis of degree n
231  // get coefficients of PkHx
232 
233  const ordinal_type lwork = card*card;
234  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
235  V1("Hdiv::Tet::In::V1", cardVecPn, card);
236 
237  // basis for the space is
238  // { (phi_i,0,0) }_{i=0}^{cardPnm1-1} ,
239  // { (0,phi_i,0) }_{i=0}^{cardPnm1-1} ,
240  // { (0,0,phi_i) }_{i=0}^{cardPnm1-1} ,
241  // { (x,y) . phi_i}_{i=cardPnm2}^{cardPnm1-1}
242  // columns of V1 are expansion of this basis in terms of the orthogonal basis
243  // for P_{n}^3
244 
245  // these two loops get the first two sets of basis functions
246  for (ordinal_type i=0;i<cardPnm1;i++) {
247  for (ordinal_type k=0; k<3;k++) {
248  V1(k*cardPn+i,k*cardPnm1+i) = 1.0;
249  }
250  }
251 
252  // now I need to integrate { (x,y,z) phi } against the big basis
253  // first, get a cubature rule.
255  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> cubPoints("Hdiv::Tet::In::cubPoints", myCub.getNumPoints() , spaceDim );
256  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> cubWeights("Hdiv::Tet::In::cubWeights", myCub.getNumPoints() );
257  myCub.getCubature( cubPoints , cubWeights );
258 
259  // tabulate the scalar orthonormal basis at cubature points
260  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> phisAtCubPoints("Hdiv::Tet::In::phisAtCubPoints", cardPn , myCub.getNumPoints() );
261  Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(phisAtCubPoints, cubPoints, order, OPERATOR_VALUE);
262 
263  // now do the integration
264  for (ordinal_type i=0;i<dim_PkH;i++) {
265  for (ordinal_type j=0;j<cardPn;j++) { // int (x,y,z) phi_i \cdot (phi_j,0,0)
266  V1(j,cardVecPnm1+i) = 0.0;
267  for (ordinal_type d=0; d< spaceDim; ++d)
268  for (ordinal_type k=0;k<myCub.getNumPoints();k++) {
269  V1(j+d*cardPn,cardVecPnm1+i) +=
270  cubWeights(k) * cubPoints(k,d)
271  * phisAtCubPoints(cardPnm2+i,k)
272  * phisAtCubPoints(j,k);
273  }
274  }
275  }
276 
277  // next, apply the RT nodes (rows) to the basis for (P_n)^3 (columns)
278  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
279  V2("Hdiv::Tet::In::V2", card ,cardVecPn);
280 
281  const ordinal_type numFaces = this->basisCellTopology_.getFaceCount();
282 
283  shards::CellTopology faceTop(shards::getCellTopologyData<shards::Triangle<3> >() );
284 
285  const int numPtsPerFace = PointTools::getLatticeSize( faceTop ,
286  order+2 ,
287  1 );
288 
289  // get the points on the tetrahedron face
290  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> triPts("Hdiv::Tet::In::triPts", numPtsPerFace , 2 );
291 
292  // construct lattice
293  const ordinal_type offset = 1;
294  PointTools::getLattice( triPts,
295  faceTop,
296  order+2,
297  offset,
298  pointType );
299 
300  // holds the image of the tet points
301  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> facePts("Hdiv::Tet::In::facePts", numPtsPerFace , spaceDim );
302  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> phisAtFacePoints("Hdiv::Tet::In::phisAtFacePoints", cardPn , numPtsPerFace );
303  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> faceNormal("Hcurl::Tet::In::faceNormal", spaceDim );
304 
305  // loop over faces
306  for (ordinal_type face=0;face<numFaces;face++) { // loop over faces
307 
308  // these are normal scaled by the appropriate face areas.
310  face ,
311  this->basisCellTopology_ );
312 
314  triPts ,
315  2 ,
316  face ,
317  this->basisCellTopology_ );
318 
319  // get phi values at face points
320  Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(phisAtFacePoints, facePts, order, OPERATOR_VALUE);
321 
322  // loop over points (rows of V2)
323  for (ordinal_type j=0;j<numPtsPerFace;j++) {
324 
325  const ordinal_type i_card = numPtsPerFace*face+j;
326 
327  // loop over orthonormal basis functions (columns of V2)
328  for (ordinal_type k=0;k<cardPn;k++) {
329  // loop over space dimension
330  for (ordinal_type l=0; l<spaceDim; l++)
331  V2(i_card,k+l*cardPn) = faceNormal(l) * phisAtFacePoints(k,j);
332  }
333 
334  //save dof coordinates and coefficients
335  for(ordinal_type l=0; l<spaceDim; ++l) {
336  dofCoords(i_card,l) = facePts(j,l);
337  dofCoeffs(i_card,l) = faceNormal(l);
338  }
339 
340  tags[i_card][0] = 2; // face dof
341  tags[i_card][1] = face; // face id
342  tags[i_card][2] = j; // local dof id
343  tags[i_card][3] = numPtsPerFace; // total vert dof
344 
345  }
346  }
347 
348  // remaining nodes point values of each vector component on interior
349  // points of a lattice of degree+2
350  // This way, RT0 --> degree = 1 and internal lattice has no points
351  // RT1 --> degree = 2, and internal lattice has one point (inside of quartic)
352  const ordinal_type numPtsPerCell = PointTools::getLatticeSize( this->basisCellTopology_ ,
353  order + 2 ,
354  1 );
355 
356  if (numPtsPerCell > 0) {
357  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
358  internalPoints( "Hdiv::Tet::In::internalPoints", numPtsPerCell , spaceDim );
359  PointTools::getLattice( internalPoints ,
360  this->basisCellTopology_ ,
361  order + 2 ,
362  1 ,
363  pointType );
364 
365  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
366  phisAtInternalPoints("Hdiv::Tet::In::phisAtInternalPoints", cardPn , numPtsPerCell );
367  Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>( phisAtInternalPoints , internalPoints , order, OPERATOR_VALUE );
368 
369  // copy values into right positions of V2
370  for (ordinal_type j=0;j<numPtsPerCell;j++) {
371 
372  const ordinal_type i_card = numFaces*numPtsPerFace+spaceDim*j;
373 
374  for (ordinal_type k=0;k<cardPn;k++) {
375  for (ordinal_type l=0;l<spaceDim;l++) {
376  V2(i_card+l,l*cardPn+k) = phisAtInternalPoints(k,j);
377  }
378  }
379 
380  //save dof coordinates and coefficients
381  for(ordinal_type d=0; d<spaceDim; ++d) {
382  for(ordinal_type l=0; l<spaceDim; ++l) {
383  dofCoords(i_card+d,l) = internalPoints(j,l);
384  dofCoeffs(i_card+d,l) = (l==d);
385  }
386 
387  tags[i_card+d][0] = spaceDim; // elem dof
388  tags[i_card+d][1] = 0; // elem id
389  tags[i_card+d][2] = spaceDim*j+d; // local dof id
390  tags[i_card+d][3] = spaceDim*numPtsPerCell; // total vert dof
391  }
392  }
393  }
394 
395  // form Vandermonde matrix. Actually, this is the transpose of the VDM,
396  // so we transpose on copy below.
397  Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
398  vmat("Hdiv::Tet::In::vmat", card, card),
399  work("Hdiv::Tet::In::work", lwork),
400  ipiv("Hdiv::Tet::In::ipiv", card);
401 
402  //vmat' = V2*V1;
403  for(ordinal_type i=0; i< card; ++i) {
404  for(ordinal_type j=0; j< card; ++j) {
405  scalarType s=0;
406  for(ordinal_type k=0; k< cardVecPn; ++k)
407  s += V2(i,k)*V1(k,j);
408  vmat(i,j) = s;
409  }
410  }
411 
412  ordinal_type info = 0;
413  Teuchos::LAPACK<ordinal_type,scalarType> lapack;
414 
415  lapack.GETRF(card, card,
416  vmat.data(), vmat.stride_1(),
417  (ordinal_type*)ipiv.data(),
418  &info);
419 
420  INTREPID2_TEST_FOR_EXCEPTION( info != 0,
421  std::runtime_error ,
422  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM) lapack.GETRF returns nonzero info." );
423 
424  lapack.GETRI(card,
425  vmat.data(), vmat.stride_1(),
426  (ordinal_type*)ipiv.data(),
427  work.data(), lwork,
428  &info);
429 
430  INTREPID2_TEST_FOR_EXCEPTION( info != 0,
431  std::runtime_error ,
432  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM) lapack.GETRI returns nonzero info." );
433 
434  for (ordinal_type i=0;i<cardVecPn;++i)
435  for (ordinal_type j=0;j<card;++j){
436  scalarType s=0;
437  for(ordinal_type k=0; k< card; ++k)
438  s += V1(i,k)*vmat(k,j);
439  coeffs(i,j) = s;
440  }
441 
442  this->coeffs_ = Kokkos::create_mirror_view(typename DT::memory_space(), coeffs);
443  Kokkos::deep_copy(this->coeffs_ , coeffs);
444 
445  this->dofCoords_ = Kokkos::create_mirror_view(typename DT::memory_space(), dofCoords);
446  Kokkos::deep_copy(this->dofCoords_, dofCoords);
447 
448  this->dofCoeffs_ = Kokkos::create_mirror_view(typename DT::memory_space(), dofCoeffs);
449  Kokkos::deep_copy(this->dofCoeffs_, dofCoeffs);
450 
451 
452  // set tags
453  {
454  // Basis-dependent initializations
455  const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
456  const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
457  const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
458 
459  OrdinalTypeArray1DHost tagView(&tags[0][0], card*tagSize);
460 
461  // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
462  // tags are constructed on host
463  this->setOrdinalTagData(this->tagToOrdinal_,
464  this->ordinalToTag_,
465  tagView,
466  this->basisCardinality_,
467  tagSize,
468  posScDim,
469  posScOrd,
470  posDfOrd);
471  }
472 }
473 } // namespace Intrepid2
474 #endif
virtual void getCubature(PointViewType cubPoints, weightViewType cubWeights) const override
Returns cubature points and weights (return arrays must be pre-sized/pre-allocated).
virtual ordinal_type getNumPoints() const override
Returns the number of cubature points.
static void getReferenceSideNormal(Kokkos::DynRankView< refSideNormalValueType, refSideNormalProperties... > refSideNormal, const ordinal_type sideOrd, const shards::CellTopology parentCell)
Computes constant normal vectors to sides of 2D or 3D reference cells.
Header file for the Intrepid2::CubatureDirectTetDefault class.
static void mapToReferenceSubcell(Kokkos::DynRankView< refSubcellPointValueType, refSubcellPointProperties... > refSubcellPoints, const Kokkos::DynRankView< paramPointValueType, paramPointProperties... > paramPoints, const ordinal_type subcellDim, const ordinal_type subcellOrd, const shards::CellTopology parentCell)
Computes parameterization maps of 1- and 2-subcells of reference cells.
EPointType
Enumeration of types of point distributions in Intrepid.
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...
Defines direct integration rules on a tetrahedron.
Basis_HDIV_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.
static ordinal_type getLatticeSize(const shards::CellTopology cellType, const ordinal_type order, const ordinal_type offset=0)
Computes the number of points in a lattice of a given order on a simplex (currently disabled for othe...