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 SpT, 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,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 = 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 SpT, 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_FIAT;
197  this->basisCoordinates_ = COORDINATES_CARTESIAN;
198  this->functionSpace_ = FUNCTION_SPACE_HDIV;
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})^3 -- larger space
206  const ordinal_type cardVecPnm1 = spaceDim*cardPnm1; // dim of (P_{n-1})^3 -- smaller space
207  const ordinal_type dim_PkH = cardPnm1 - cardPnm2;
208 
209 
210  // Basis-dependent initializations
211  constexpr ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
212  constexpr ordinal_type maxCard = CardinalityHDivTet(Parameters::MaxOrder);
213  ordinal_type tags[maxCard][tagSize];
214 
215  // points are computed in the host and will be copied
216  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
217  dofCoords("Hdiv::Tet::In::dofCoords", card, spaceDim);
218 
219  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
220  dofCoeffs("Hdiv::Tet::In::dofCoeffs", card, spaceDim);
221 
222  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
223  coeffs("Hdiv::Tet::In::coeffs", cardVecPn, card);
224 
225  // first, need to project the basis for RT space onto the
226  // orthogonal basis of degree n
227  // get coefficients of PkHx
228 
229  const ordinal_type lwork = card*card;
230  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
231  V1("Hdiv::Tet::In::V1", cardVecPn, card);
232 
233  // basis for the space is
234  // { (phi_i,0,0) }_{i=0}^{cardPnm1-1} ,
235  // { (0,phi_i,0) }_{i=0}^{cardPnm1-1} ,
236  // { (0,0,phi_i) }_{i=0}^{cardPnm1-1} ,
237  // { (x,y) . phi_i}_{i=cardPnm2}^{cardPnm1-1}
238  // columns of V1 are expansion of this basis in terms of the orthogonal basis
239  // for P_{n}^3
240 
241  // these two loops get the first two sets of basis functions
242  for (ordinal_type i=0;i<cardPnm1;i++) {
243  for (ordinal_type k=0; k<3;k++) {
244  V1(k*cardPn+i,k*cardPnm1+i) = 1.0;
245  }
246  }
247 
248  // now I need to integrate { (x,y,z) phi } against the big basis
249  // first, get a cubature rule.
251  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> cubPoints("Hdiv::Tet::In::cubPoints", myCub.getNumPoints() , spaceDim );
252  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> cubWeights("Hdiv::Tet::In::cubWeights", myCub.getNumPoints() );
253  myCub.getCubature( cubPoints , cubWeights );
254 
255  // tabulate the scalar orthonormal basis at cubature points
256  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> phisAtCubPoints("Hdiv::Tet::In::phisAtCubPoints", cardPn , myCub.getNumPoints() );
257  Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(phisAtCubPoints, cubPoints, order, OPERATOR_VALUE);
258 
259  // now do the integration
260  for (ordinal_type i=0;i<dim_PkH;i++) {
261  for (ordinal_type j=0;j<cardPn;j++) { // int (x,y,z) phi_i \cdot (phi_j,0,0)
262  V1(j,cardVecPnm1+i) = 0.0;
263  for (ordinal_type d=0; d< spaceDim; ++d)
264  for (ordinal_type k=0;k<myCub.getNumPoints();k++) {
265  V1(j+d*cardPn,cardVecPnm1+i) +=
266  cubWeights(k) * cubPoints(k,d)
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)^3 (columns)
274  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
275  V2("Hdiv::Tet::In::V2", card ,cardVecPn);
276 
277  const ordinal_type numFaces = this->basisCellTopology_.getFaceCount();
278 
279  shards::CellTopology faceTop(shards::getCellTopologyData<shards::Triangle<3> >() );
280 
281  const int numPtsPerFace = PointTools::getLatticeSize( faceTop ,
282  order+2 ,
283  1 );
284 
285  // get the points on the tetrahedron face
286  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> triPts("Hdiv::Tet::In::triPts", numPtsPerFace , 2 );
287 
288  // construct lattice
289  const ordinal_type offset = 1;
290  PointTools::getLattice( triPts,
291  faceTop,
292  order+2,
293  offset,
294  pointType );
295 
296  // holds the image of the tet points
297  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> facePts("Hdiv::Tet::In::facePts", numPtsPerFace , spaceDim );
298  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> phisAtFacePoints("Hdiv::Tet::In::phisAtFacePoints", cardPn , numPtsPerFace );
299  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> faceNormal("Hcurl::Tet::In::faceNormal", spaceDim );
300 
301  // loop over faces
302  for (ordinal_type face=0;face<numFaces;face++) { // loop over faces
303 
304  // these are normal scaled by the appropriate face areas.
306  face ,
307  this->basisCellTopology_ );
308 
309  // multiply to account for reference face areas, so that magnitude of faceNormal is equal to the face measure
310  const scalarType refTriangleMeasure = 0.5;
311  for (ordinal_type j=0;j<spaceDim;j++)
312  faceNormal(j) *= refTriangleMeasure;
313 
314 
316  triPts ,
317  2 ,
318  face ,
319  this->basisCellTopology_ );
320 
321  // get phi values at face points
322  Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(phisAtFacePoints, facePts, order, OPERATOR_VALUE);
323 
324  // loop over points (rows of V2)
325  for (ordinal_type j=0;j<numPtsPerFace;j++) {
326 
327  const ordinal_type i_card = numPtsPerFace*face+j;
328 
329  // loop over orthonormal basis functions (columns of V2)
330  for (ordinal_type k=0;k<cardPn;k++) {
331  // loop over space dimension
332  for (ordinal_type l=0; l<spaceDim; l++)
333  V2(i_card,k+l*cardPn) = faceNormal(l) * phisAtFacePoints(k,j);
334  }
335 
336  //save dof coordinates and coefficients
337  for(ordinal_type l=0; l<spaceDim; ++l) {
338  dofCoords(i_card,l) = facePts(j,l);
339  dofCoeffs(i_card,l) = faceNormal(l);
340  }
341 
342  tags[i_card][0] = 2; // face dof
343  tags[i_card][1] = face; // face id
344  tags[i_card][2] = j; // local dof id
345  tags[i_card][3] = numPtsPerFace; // total vert dof
346 
347  }
348  }
349 
350  // remaining nodes point values of each vector component on interior
351  // points of a lattice of degree+2
352  // This way, RT0 --> degree = 1 and internal lattice has no points
353  // RT1 --> degree = 2, and internal lattice has one point (inside of quartic)
354  const ordinal_type numPtsPerCell = PointTools::getLatticeSize( this->basisCellTopology_ ,
355  order + 2 ,
356  1 );
357 
358  if (numPtsPerCell > 0) {
359  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
360  internalPoints( "Hdiv::Tet::In::internalPoints", numPtsPerCell , spaceDim );
361  PointTools::getLattice( internalPoints ,
362  this->basisCellTopology_ ,
363  order + 2 ,
364  1 ,
365  pointType );
366 
367  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
368  phisAtInternalPoints("Hdiv::Tet::In::phisAtInternalPoints", cardPn , numPtsPerCell );
369  Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>( phisAtInternalPoints , internalPoints , order, OPERATOR_VALUE );
370 
371  // copy values into right positions of V2
372  for (ordinal_type j=0;j<numPtsPerCell;j++) {
373 
374  const ordinal_type i_card = numFaces*numPtsPerFace+spaceDim*j;
375 
376  for (ordinal_type k=0;k<cardPn;k++) {
377  for (ordinal_type l=0;l<spaceDim;l++) {
378  V2(i_card+l,l*cardPn+k) = phisAtInternalPoints(k,j);
379  }
380  }
381 
382  //save dof coordinates and coefficients
383  for(ordinal_type d=0; d<spaceDim; ++d) {
384  for(ordinal_type l=0; l<spaceDim; ++l) {
385  dofCoords(i_card+d,l) = internalPoints(j,l);
386  dofCoeffs(i_card+d,l) = (l==d);
387  }
388 
389  tags[i_card+d][0] = spaceDim; // elem dof
390  tags[i_card+d][1] = 0; // elem id
391  tags[i_card+d][2] = spaceDim*j+d; // local dof id
392  tags[i_card+d][3] = spaceDim*numPtsPerCell; // total vert dof
393  }
394  }
395  }
396 
397  // form Vandermonde matrix. Actually, this is the transpose of the VDM,
398  // so we transpose on copy below.
399  Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
400  vmat("Hdiv::Tet::In::vmat", card, card),
401  work("Hdiv::Tet::In::work", lwork),
402  ipiv("Hdiv::Tet::In::ipiv", card);
403 
404  //vmat' = V2*V1;
405  for(ordinal_type i=0; i< card; ++i) {
406  for(ordinal_type j=0; j< card; ++j) {
407  scalarType s=0;
408  for(ordinal_type k=0; k< cardVecPn; ++k)
409  s += V2(i,k)*V1(k,j);
410  vmat(i,j) = s;
411  }
412  }
413 
414  ordinal_type info = 0;
415  Teuchos::LAPACK<ordinal_type,scalarType> lapack;
416 
417  lapack.GETRF(card, card,
418  vmat.data(), vmat.stride_1(),
419  (ordinal_type*)ipiv.data(),
420  &info);
421 
422  INTREPID2_TEST_FOR_EXCEPTION( info != 0,
423  std::runtime_error ,
424  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM) lapack.GETRF returns nonzero info." );
425 
426  lapack.GETRI(card,
427  vmat.data(), vmat.stride_1(),
428  (ordinal_type*)ipiv.data(),
429  work.data(), lwork,
430  &info);
431 
432  INTREPID2_TEST_FOR_EXCEPTION( info != 0,
433  std::runtime_error ,
434  ">>> ERROR: (Intrepid2::Basis_HDIV_TET_In_FEM) lapack.GETRI returns nonzero info." );
435 
436  for (ordinal_type i=0;i<cardVecPn;++i)
437  for (ordinal_type j=0;j<card;++j){
438  scalarType s=0;
439  for(ordinal_type k=0; k< card; ++k)
440  s += V1(i,k)*vmat(k,j);
441  coeffs(i,j) = s;
442  }
443 
444  this->coeffs_ = Kokkos::create_mirror_view(typename SpT::memory_space(), coeffs);
445  Kokkos::deep_copy(this->coeffs_ , coeffs);
446 
447  this->dofCoords_ = Kokkos::create_mirror_view(typename SpT::memory_space(), dofCoords);
448  Kokkos::deep_copy(this->dofCoords_, dofCoords);
449 
450  this->dofCoeffs_ = Kokkos::create_mirror_view(typename SpT::memory_space(), dofCoeffs);
451  Kokkos::deep_copy(this->dofCoeffs_, dofCoeffs);
452 
453 
454  // set tags
455  {
456  // Basis-dependent initializations
457  const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
458  const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
459  const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
460 
461  OrdinalTypeArray1DHost tagView(&tags[0][0], card*tagSize);
462 
463  // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
464  // tags are constructed on host
465  this->setOrdinalTagData(this->tagToOrdinal_,
466  this->ordinalToTag_,
467  tagView,
468  this->basisCardinality_,
469  tagSize,
470  posScDim,
471  posScOrd,
472  posDfOrd);
473  }
474 }
475 } // namespace Intrepid2
476 #endif
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
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.
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).
Header file for the Intrepid2::CubatureDirectTetDefault class.
Basis_HDIV_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Defines direct integration rules on a tetrahedron.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
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...