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