Intrepid2
Intrepid2_HCURL_TRI_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_HCURL_TRI_IN_FEM_DEF_HPP__
50 #define __INTREPID2_HCURL_TRI_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_HCURL_TRI_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 = 2;
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 == CardinalityHCurlTri(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_TRI_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(k,j);
107  }
108  break;
109  }
110  case OPERATOR_CURL: {
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_TRI_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  output.access(i,j) += - coeffs(k,i)*phis(k,j,1) // - dy of x component
123  + coeffs(k+cardPn,i)*phis(k,j,0); // dx of y component
124  }
125  break;
126  }
127  default: {
128  INTREPID2_TEST_FOR_ABORT( true,
129  ">>> ERROR (Basis_HCURL_TRI_In_FEM): Operator type not implemented");
130  }
131  }
132  }
133 
134  template<typename SpT, ordinal_type numPtsPerEval,
135  typename outputValueValueType, class ...outputValueProperties,
136  typename inputPointValueType, class ...inputPointProperties,
137  typename vinvValueType, class ...vinvProperties>
138  void
139  Basis_HCURL_TRI_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,SpT>::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 = 2;
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_HCURL_TRI_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_CURL: {
172  workViewType work(Kokkos::view_alloc("Basis_HCURL_TRI_In_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
173  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
174  OPERATOR_CURL,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_HCURL_TRI_In_FEM): Operator type not implemented" );
181  }
182  }
183  }
184  }
185 
186  // -------------------------------------------------------------------------------------
187  template<typename SpT, typename OT, typename PT>
189  Basis_HCURL_TRI_In_FEM( const ordinal_type order,
190  const EPointType pointType ) {
191 
192  constexpr ordinal_type spaceDim = 2;
193  this->basisCardinality_ = CardinalityHCurlTri(order);
194  this->basisDegree_ = order; // small n
195  this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Triangle<3> >() );
196  this->basisType_ = BASIS_FEM_FIAT;
197  this->basisCoordinates_ = COORDINATES_CARTESIAN;
198  this->functionSpace_ = FUNCTION_SPACE_HCURL;
199 
200  const ordinal_type card = this->basisCardinality_;
201 
202  const ordinal_type cardPn = Intrepid2::getPnCardinality<spaceDim>(order); // dim of (P_{n}) -- smaller space
203  const ordinal_type cardPnm1 = Intrepid2::getPnCardinality<spaceDim>(order-1); // dim of (P_{n-1}) -- smaller space
204  const ordinal_type cardPnm2 = Intrepid2::getPnCardinality<spaceDim>(order-2); // dim of (P_{n-2}) -- smaller space
205  const ordinal_type cardVecPn = spaceDim*cardPn; // dim of (P_{n})^2 -- larger space
206  const ordinal_type cardVecPnm1 = spaceDim*cardPnm1; // dim of (P_{n-1})^2 -- smaller space
207 
208 
209  // Basis-dependent initializations
210  constexpr ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
211  constexpr ordinal_type maxCard = CardinalityHCurlTri(Parameters::MaxOrder);
212  ordinal_type tags[maxCard][tagSize];
213 
214  // points are computed in the host and will be copied
215  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
216  dofCoords("Hcurl::Tri::In::dofCoords", card, spaceDim);
217 
218  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
219  coeffs("Hcurl::Tri::In::coeffs", cardVecPn, card);
220 
221  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
222  dofCoeffs("Hcurl::Tri::In::dofCoeffs", card, spaceDim);
223 
224  // first, need to project the basis for RT space onto the
225  // orthogonal basis of degree n
226  // get coefficients of PkHx
227 
228  const ordinal_type lwork = card*card;
229  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
230  V1("Hcurl::Tri::In::V1", cardVecPn, card);
231 
232  // basis for the space is
233  // { (phi_i,0) }_{i=0}^{cardPnm1-1} ,
234  // { (0,phi_i) }_{i=0}^{cardPnm1-1} ,
235  // { (x,y) \times phi_i}_{i=cardPnm2}^{cardPnm1-1}
236  // { (x,y) \times phi = (y phi , -x \phi)
237  // columns of V1 are expansion of this basis in terms of the basis
238  // for P_{n}^2
239 
240  // these two loops get the first two sets of basis functions
241  for (ordinal_type i=0;i<cardPnm1;i++)
242  for (ordinal_type d=0;d<spaceDim;d++)
243  V1(d*cardPn+i,d*cardPnm1+i) = 1.0;
244 
245 
246  // now I need to integrate { (x,y) \times phi } against the big basis
247  // first, get a cubature rule.
249  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> cubPoints("Hcurl::Tri::In::cubPoints", myCub.getNumPoints() , spaceDim );
250  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> cubWeights("Hcurl::Tri::In::cubWeights", myCub.getNumPoints() );
251  myCub.getCubature( cubPoints , cubWeights );
252 
253  // tabulate the scalar orthonormal basis at cubature points
254  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> phisAtCubPoints("Hcurl::Tri::In::phisAtCubPoints", cardPn , myCub.getNumPoints() );
255  Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(phisAtCubPoints, cubPoints, order, OPERATOR_VALUE);
256 
257  // now do the integration
258  for (ordinal_type i=0;i<order;i++) {
259  for (ordinal_type j=0;j<cardPn;j++) { // int (x,y) phi_i \cdot (phi_j,phi_{j+cardPn})
260  for (ordinal_type k=0;k<myCub.getNumPoints();k++) {
261  V1(j,cardVecPnm1+i) -=
262  cubWeights(k) * cubPoints(k,1)
263  * phisAtCubPoints(cardPnm2+i,k)
264  * phisAtCubPoints(j,k);
265  V1(j+cardPn,cardVecPnm1+i) +=
266  cubWeights(k) * cubPoints(k,0)
267  * phisAtCubPoints(cardPnm2+i,k)
268  * phisAtCubPoints(j,k);
269  }
270  }
271  }
272 
273  // next, apply the RT nodes (rows) to the basis for (P_n)^2 (columns)
274  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
275  V2("Hcurl::Tri::In::V2", card ,cardVecPn);
276 
277  const ordinal_type numEdges = this->basisCellTopology_.getEdgeCount();
278 
279  shards::CellTopology edgeTop(shards::getCellTopologyData<shards::Line<2> >() );
280 
281  const int numPtsPerEdge = PointTools::getLatticeSize( edgeTop ,
282  order+1 ,
283  1 );
284 
285  // first numEdges * degree nodes are tangents at each edge
286  // get the points on the line
287  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> linePts("Hcurl::Tri::In::linePts", numPtsPerEdge , 1 );
288 
289  // construct lattice
290  const ordinal_type offset = 1;
291  PointTools::getLattice( linePts,
292  edgeTop,
293  order+1, offset,
294  pointType );
295 
296  // holds the image of the line points
297  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> edgePts("Hcurl::Tri::In::edgePts", numPtsPerEdge , spaceDim );
298  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> phisAtEdgePoints("Hcurl::Tri::In::phisAtEdgePoints", cardPn , numPtsPerEdge );
299  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> edgeTan("Hcurl::Tri::In::edgeTan", spaceDim );
300 
301  // these are tangents scaled by the appropriate edge lengths.
302  for (ordinal_type edge=0;edge<numEdges;edge++) { // loop over edges
304  edge ,
305  this->basisCellTopology_ );
306  /* multiply by measure of reference edge so that magnitude of the edgeTan is equal to the edge measure */
307  const scalarType refEdgeMeasure = 2.0;
308  for (ordinal_type j=0;j<spaceDim;j++)
309  edgeTan(j) *= refEdgeMeasure;
310 
312  linePts ,
313  1 ,
314  edge ,
315  this->basisCellTopology_ );
316 
317  Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(phisAtEdgePoints , edgePts, order, OPERATOR_VALUE);
318 
319  // loop over points (rows of V2)
320  for (ordinal_type j=0;j<numPtsPerEdge;j++) {
321 
322  const ordinal_type i_card = numPtsPerEdge*edge+j;
323 
324  // loop over orthonormal basis functions (columns of V2)
325  for (ordinal_type k=0;k<cardPn;k++) {
326  V2(i_card,k) = edgeTan(0) * phisAtEdgePoints(k,j);
327  V2(i_card,k+cardPn) = edgeTan(1) * phisAtEdgePoints(k,j);
328  }
329 
330 
331  //save dof coordinates
332  for(ordinal_type k=0; k<spaceDim; ++k) {
333  dofCoords(i_card,k) = edgePts(j,k);
334  dofCoeffs(i_card,k) = edgeTan(k);
335  }
336 
337  tags[i_card][0] = 1; // edge dof
338  tags[i_card][1] = edge; // edge id
339  tags[i_card][2] = j; // local dof id
340  tags[i_card][3] = numPtsPerEdge; // total edge dof
341 
342  }
343 
344 
345  }
346 
347  // remaining nodes are x- and y- components at internal points (this code is same as HDIV).
348  //These are evaluated at the interior of a lattice of degree + 1, For then
349  // the degree == 1 space corresponds classicaly to RT0 and so gets
350  // no internal nodes, and degree == 2 corresponds to RT1 and needs
351  // one internal node per vector component.
352  const ordinal_type numPtsPerCell = PointTools::getLatticeSize( this->basisCellTopology_ ,
353  order + 1 ,
354  1 );
355 
356  if (numPtsPerCell > 0) {
357  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
358  internalPoints( "Hcurl::Tri::In::internalPoints", numPtsPerCell , spaceDim );
359  PointTools::getLattice( internalPoints ,
360  this->basisCellTopology_ ,
361  order + 1 ,
362  1 ,
363  pointType );
364 
365  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
366  phisAtInternalPoints("Hcurl::Tri::In::phisAtInternalPoints", cardPn , numPtsPerCell );
367  Impl::Basis_HGRAD_TRI_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 = numEdges*order+spaceDim*j;
373 
374  for (ordinal_type k=0;k<cardPn;k++) {
375  // x component
376  V2(i_card,k) = phisAtInternalPoints(k,j);
377  // y component
378  V2(i_card+1,cardPn+k) = phisAtInternalPoints(k,j);
379  }
380 
381  //save dof coordinates
382  for(ordinal_type d=0; d<spaceDim; ++d) {
383  for(ordinal_type dim=0; dim<spaceDim; ++dim) {
384  dofCoords(i_card+d,dim) = internalPoints(j,dim);
385  dofCoeffs(i_card+d,dim) = (d==dim);
386  }
387 
388  tags[i_card+d][0] = spaceDim; // elem dof
389  tags[i_card+d][1] = 0; // elem id
390  tags[i_card+d][2] = spaceDim*j+d; // local dof id
391  tags[i_card+d][3] = spaceDim*numPtsPerCell; // total vert dof
392  }
393  }
394  }
395 
396  // form Vandermonde matrix. Actually, this is the transpose of the VDM,
397  // so we transpose on copy below.
398  Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
399  vmat("Hcurl::Tri::In::vmat", card, card),
400  work("Hcurl::Tri::In::work", lwork),
401  ipiv("Hcurl::Tri::In::ipiv", card);
402 
403  //vmat' = V2*V1;
404  for(ordinal_type i=0; i< card; ++i) {
405  for(ordinal_type j=0; j< card; ++j) {
406  scalarType s=0;
407  for(ordinal_type k=0; k< cardVecPn; ++k)
408  s += V2(i,k)*V1(k,j);
409  vmat(i,j) = s;
410  }
411  }
412 
413  ordinal_type info = 0;
414  Teuchos::LAPACK<ordinal_type,scalarType> lapack;
415 
416  lapack.GETRF(card, card,
417  vmat.data(), vmat.stride_1(),
418  (ordinal_type*)ipiv.data(),
419  &info);
420 
421  INTREPID2_TEST_FOR_EXCEPTION( info != 0,
422  std::runtime_error ,
423  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM) lapack.GETRF returns nonzero info." );
424 
425  lapack.GETRI(card,
426  vmat.data(), vmat.stride_1(),
427  (ordinal_type*)ipiv.data(),
428  work.data(), lwork,
429  &info);
430 
431  INTREPID2_TEST_FOR_EXCEPTION( info != 0,
432  std::runtime_error ,
433  ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM) lapack.GETRI returns nonzero info." );
434 
435  for (ordinal_type i=0;i<cardVecPn;++i)
436  for (ordinal_type j=0;j<card;++j){
437  scalarType s=0;
438  for(ordinal_type k=0; k< card; ++k)
439  s += V1(i,k)*vmat(k,j);
440  coeffs(i,j) = s;
441  }
442 
443  this->coeffs_ = Kokkos::create_mirror_view(typename SpT::memory_space(), coeffs);
444  Kokkos::deep_copy(this->coeffs_ , coeffs);
445 
446  this->dofCoords_ = Kokkos::create_mirror_view(typename SpT::memory_space(), dofCoords);
447  Kokkos::deep_copy(this->dofCoords_, dofCoords);
448 
449  this->dofCoeffs_ = Kokkos::create_mirror_view(typename SpT::memory_space(), dofCoeffs);
450  Kokkos::deep_copy(this->dofCoeffs_, dofCoeffs);
451 
452 
453  // set tags
454  {
455  // Basis-dependent initializations
456  const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
457  const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
458  const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
459 
460  OrdinalTypeArray1DHost tagView(&tags[0][0], card*tagSize);
461 
462  // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
463  // tags are constructed on host
464  this->setOrdinalTagData(this->tagToOrdinal_,
465  this->ordinalToTag_,
466  tagView,
467  this->basisCardinality_,
468  tagSize,
469  posScDim,
470  posScOrd,
471  posDfOrd);
472  }
473  }
474 } // namespace Intrepid2
475 #endif
static void getReferenceEdgeTangent(Kokkos::DynRankView< refEdgeTangentValueType, refEdgeTangentProperties...> refEdgeTangent, const ordinal_type edgeOrd, const shards::CellTopology parentCell)
Computes constant tangent vectors to edges of 2D or 3D reference cells.
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
Header file for the Intrepid2::CubatureDirectTriDefault class.
Defines direct integration rules on a triangle.
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 (currently disabled for other ce...
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.
virtual void getCubature(PointViewType cubPoints, weightViewType cubWeights) const
Returns cubature points and weights (return arrays must be pre-sized/pre-allocated).
Basis_HCURL_TRI_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.
virtual ordinal_type getNumPoints() const
Returns the number of cubature points.
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...