49 #ifndef __INTREPID2_HCURL_TET_IN_FEM_DEF_HPP__
50 #define __INTREPID2_HCURL_TET_IN_FEM_DEF_HPP__
54 #include "Teuchos_SerialDenseMatrix.hpp"
62 template<EOperator opType>
63 template<
typename OutputViewType,
64 typename inputViewType,
65 typename workViewType,
66 typename vinvViewType>
67 KOKKOS_INLINE_FUNCTION
69 Basis_HCURL_TET_In_FEM::Serial<opType>::
70 getValues( OutputViewType output,
71 const inputViewType input,
73 const vinvViewType coeffs ) {
75 constexpr ordinal_type spaceDim = 3;
77 cardPn = coeffs.extent(0)/spaceDim,
78 card = coeffs.extent(1),
79 npts = input.extent(0);
82 ordinal_type order = 0;
84 if (card == CardinalityHCurlTet(p)) {
90 typedef typename Kokkos::DynRankView<typename workViewType::value_type, typename workViewType::memory_space> viewType;
91 auto vcprop = Kokkos::common_view_alloc_prop(work);
92 auto ptr = work.data();
95 case OPERATOR_VALUE: {
96 const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts);
97 workViewType dummyView;
99 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::
100 Serial<opType>::getValues(phis, input, dummyView, order);
102 for (ordinal_type i=0;i<card;++i)
103 for (ordinal_type j=0;j<npts;++j)
104 for (ordinal_type d=0;d<spaceDim;++d) {
105 output.access(i,j,d) = 0.0;
106 for (ordinal_type k=0;k<cardPn;++k)
107 output.access(i,j,d) += coeffs(k+d*cardPn,i) * phis(k,j);
111 case OPERATOR_CURL: {
112 const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim);
113 ptr += card*npts*spaceDim*get_dimension_scalar(work);
114 const viewType workView(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim+1);
116 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::
117 Serial<OPERATOR_GRAD>::getValues(phis, input, workView, order);
119 for (ordinal_type i=0;i<card;++i) {
120 for (ordinal_type j=0;j<npts;++j) {
121 for (ordinal_type d=0; d< spaceDim; ++d) {
122 output.access(i,j,d) = 0.0;
123 ordinal_type d1 = (d+1) % spaceDim, d2 = (d+2) % spaceDim;
124 for (ordinal_type k=0; k<cardPn; ++k)
125 output.access(i,j,d) += coeffs(k+d2*cardPn,i)*phis(k,j,d1)
126 -coeffs(k+d1*cardPn,i)*phis(k,j,d2);
133 INTREPID2_TEST_FOR_ABORT(
true,
134 ">>> ERROR (Basis_HCURL_TET_In_FEM): Operator type not implemented");
139 template<
typename SpT, ordinal_type numPtsPerEval,
140 typename outputValueValueType,
class ...outputValueProperties,
141 typename inputPointValueType,
class ...inputPointProperties,
142 typename vinvValueType,
class ...vinvProperties>
144 Basis_HCURL_TET_In_FEM::
145 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
146 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
147 const Kokkos::DynRankView<vinvValueType, vinvProperties...> coeffs,
148 const EOperator operatorType) {
149 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
150 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
151 typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
152 typedef typename ExecSpace<typename inputPointViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
155 const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
156 const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
157 const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
158 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
160 typedef typename inputPointViewType::value_type inputPointType;
162 const ordinal_type cardinality = outputValues.extent(0);
163 const ordinal_type spaceDim = 3;
165 auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
166 typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
168 switch (operatorType) {
169 case OPERATOR_VALUE: {
170 workViewType work(Kokkos::view_alloc(
"Basis_HCURL_TET_In_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0));
171 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
172 OPERATOR_VALUE,numPtsPerEval> FunctorType;
173 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, coeffs, work) );
176 case OPERATOR_CURL: {
177 workViewType work(Kokkos::view_alloc(
"Basis_HCURL_TET_In_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
178 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
179 OPERATOR_CURL,numPtsPerEval> FunctorType;
180 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, coeffs, work) );
184 INTREPID2_TEST_FOR_EXCEPTION(
true , std::invalid_argument,
185 ">>> ERROR (Basis_HCURL_TET_In_FEM): Operator type not implemented" );
192 template<
typename SpT,
typename OT,
typename PT>
195 const EPointType pointType ) {
197 constexpr ordinal_type spaceDim = 3;
198 this->basisCardinality_ = CardinalityHCurlTet(order);
199 this->basisDegree_ = order;
200 this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Tetrahedron<4> >() );
201 this->basisType_ = BASIS_FEM_FIAT;
202 this->basisCoordinates_ = COORDINATES_CARTESIAN;
203 this->functionSpace_ = FUNCTION_SPACE_HCURL;
205 const ordinal_type card = this->basisCardinality_;
207 const ordinal_type cardPn = Intrepid2::getPnCardinality<spaceDim>(order);
208 const ordinal_type cardPnm1 = Intrepid2::getPnCardinality<spaceDim>(order-1);
209 const ordinal_type cardPnm2 = Intrepid2::getPnCardinality<spaceDim>(order-2);
210 const ordinal_type cardVecPn = spaceDim*cardPn;
211 const ordinal_type cardVecPnm1 = spaceDim*cardPnm1;
212 const ordinal_type cardPnm1H = cardPnm1-cardPnm2;
217 constexpr ordinal_type tagSize = 4;
219 ordinal_type tags[maxCard][tagSize];
222 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
223 dofCoords(
"Hcurl::Tet::In::dofCoords", card, spaceDim);
225 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
226 coeffs(
"Hcurl::Tet::In::coeffs", cardVecPn, card);
228 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
229 dofCoeffs(
"Hcurl::Tet::In::dofCoeffs", card, spaceDim);
235 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
236 V1(
"Hcurl::Tet::In::V1", cardVecPn, cardVecPnm1 + spaceDim*cardPnm1H);
240 for (ordinal_type i=0;i<cardPnm1;i++)
241 for (ordinal_type d=0;d<spaceDim;d++)
242 V1(i+d*cardPn,i+d*cardPnm1) = 1.0;
248 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> cubPoints(
"Hcurl::Tet::In::cubPoints", myCub.
getNumPoints() , spaceDim );
249 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> cubWeights(
"Hcurl::Tet::In::cubWeights", myCub.
getNumPoints() );
253 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> phisAtCubPoints(
"Hcurl::Tet::In::phisAtCubPoints", cardPn , myCub.
getNumPoints() );
254 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(phisAtCubPoints, cubPoints, order, OPERATOR_VALUE);
257 for (ordinal_type i=0;i<cardPn;i++) {
258 for (ordinal_type j=0;j<cardPnm1H;j++) {
259 for (ordinal_type d=0; d< spaceDim; ++d) {
260 scalarType integral(0);
262 integral += cubWeights(k) * cubPoints(k,d)
263 * phisAtCubPoints(cardPnm2+j,k)
264 * phisAtCubPoints(i,k);
265 ordinal_type d1 = (d+1) % spaceDim, d2 = (d+2) % spaceDim;
266 V1(i+d2*cardPn,cardVecPnm1+d1*cardPnm1H + j) = -integral;
267 V1(i+d1*cardPn,cardVecPnm1+d2*cardPnm1H + j) = integral;
277 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
278 S(
"Hcurl::Tet::In::S", cardVecPn,1),
279 U(
"Hcurl::Tet::In::U", cardVecPn, cardVecPn),
280 Vt(
"Hcurl::Tet::In::Vt", cardVecPn, cardVecPn),
281 work(
"Hcurl::Tet::In::work", 5*cardVecPn,1),
282 rWork(
"Hcurl::Tet::In::rW", 1,1);
286 ordinal_type info = 0;
287 Teuchos::LAPACK<ordinal_type,scalarType> lapack;
307 #ifdef HAVE_INTREPID2_DEBUG
308 ordinal_type num_nonzero_sv = 0;
309 for (
int i=0;i<cardVecPn;i++)
310 num_nonzero_sv += (S(i,0) > tolerence());
312 INTREPID2_TEST_FOR_EXCEPTION( num_nonzero_sv != card, std::invalid_argument,
313 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM( order, pointType), Matrix V1 should have rank equal to the cardinality of HCURL space");
317 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
318 V2(
"Hcurl::Tet::In::V2", card ,cardVecPn);
320 const ordinal_type numEdges = this->basisCellTopology_.getEdgeCount();
321 const ordinal_type numFaces = this->basisCellTopology_.getFaceCount();
326 shards::CellTopology edgeTop(shards::getCellTopologyData<shards::Line<2> >() );
327 shards::CellTopology faceTop(shards::getCellTopologyData<shards::Triangle<3> >() );
341 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> linePts(
"Hcurl::Tet::In::linePts", numPtsPerEdge , 1 );
342 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> triPts(
"Hcurl::Tet::In::triPts", numPtsPerFace , 2 );
345 const ordinal_type offset = 1;
360 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> edgePts(
"Hcurl::Tet::In::edgePts", numPtsPerEdge , spaceDim );
361 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> facePts(
"Hcurl::Tet::In::facePts", numPtsPerFace , spaceDim );
362 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> phisAtEdgePoints(
"Hcurl::Tet::In::phisAtEdgePoints", cardPn , numPtsPerEdge );
363 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> phisAtFacePoints(
"Hcurl::Tet::In::phisAtFacePoints", cardPn , numPtsPerFace);
365 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> edgeTan(
"Hcurl::Tet::In::edgeTan", spaceDim );
368 for (ordinal_type i=0;i<numEdges;i++) {
371 this->basisCellTopology_ );
373 const scalarType refEdgeMeasure = 2.0;
374 for (ordinal_type j=0;j<spaceDim;j++)
375 edgeTan(j) *= refEdgeMeasure;
381 this->basisCellTopology_ );
383 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(phisAtEdgePoints , edgePts, order, OPERATOR_VALUE);
386 for (ordinal_type j=0;j<numPtsPerEdge;j++) {
388 const ordinal_type i_card = numPtsPerEdge*i+j;
391 for (ordinal_type k=0;k<cardPn;k++)
392 for (ordinal_type d=0;d<spaceDim;d++)
393 V2(i_card,k+d*cardPn) = edgeTan(d) * phisAtEdgePoints(k,j);
396 for(ordinal_type k=0; k<spaceDim; ++k) {
397 dofCoords(i_card,k) = edgePts(j,k);
398 dofCoeffs(i_card,k) = edgeTan(k);
404 tags[i_card][3] = numPtsPerEdge;
409 if(numPtsPerFace >0) {
410 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> faceTan1(
"Hcurl::Tet::In::edgeTan", spaceDim );
411 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> faceTan2(
"Hcurl::Tet::In::edgeTan", spaceDim );
413 for (ordinal_type i=0;i<numFaces;i++) {
417 this->basisCellTopology_ );
423 this->basisCellTopology_ );
425 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(phisAtFacePoints , facePts, order, OPERATOR_VALUE);
428 for (ordinal_type j=0;j<numPtsPerFace;j++) {
430 const ordinal_type i_card = numEdges*numPtsPerEdge+2*numPtsPerFace*i+2*j;
431 const ordinal_type i_card_p1 = i_card+1;
434 for (ordinal_type k=0;k<cardPn;k++)
435 for (ordinal_type d=0;d<spaceDim;d++) {
436 V2(i_card,k+d*cardPn) = faceTan1(d) * phisAtFacePoints(k,j);
437 V2(i_card_p1,k+d*cardPn) = faceTan2(d) * phisAtFacePoints(k,j);
441 for(ordinal_type k=0; k<spaceDim; ++k) {
442 dofCoords(i_card,k) = facePts(j,k);
443 dofCoords(i_card_p1,k) = facePts(j,k);
444 dofCoeffs(i_card,k) = faceTan1(k);
445 dofCoeffs(i_card_p1,k) = faceTan2(k);
450 tags[i_card][2] = 2*j;
451 tags[i_card][3] = 2*numPtsPerFace;
453 tags[i_card_p1][0] = 2;
454 tags[i_card_p1][1] = i;
455 tags[i_card_p1][2] = 2*j+1;
456 tags[i_card_p1][3] = 2*numPtsPerFace;
464 if (numPtsPerCell > 0) {
465 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
466 cellPoints(
"Hcurl::Tet::In::cellPoints", numPtsPerCell , spaceDim );
468 this->basisCellTopology_ ,
473 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
474 phisAtCellPoints(
"Hcurl::Tet::In::phisAtCellPoints", cardPn , numPtsPerCell );
475 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>( phisAtCellPoints , cellPoints , order, OPERATOR_VALUE );
478 for (ordinal_type j=0;j<numPtsPerCell;j++) {
480 const ordinal_type i_card = numEdges*numPtsPerEdge+2*numFaces*numPtsPerFace+spaceDim*j;
482 for (ordinal_type k=0;k<cardPn;k++)
483 for (ordinal_type d=0;d<spaceDim;d++)
484 V2(i_card+d,d*cardPn+k) = phisAtCellPoints(k,j);
488 for(ordinal_type d=0; d<spaceDim; ++d) {
489 for(ordinal_type dim=0; dim<spaceDim; ++dim) {
490 dofCoords(i_card+d,dim) = cellPoints(j,dim);
491 dofCoeffs(i_card+d,dim) = (d==dim);
494 tags[i_card+d][0] = spaceDim;
495 tags[i_card+d][1] = 0;
496 tags[i_card+d][2] = spaceDim*j+d;
497 tags[i_card+d][3] = spaceDim*numPtsPerCell;
504 const ordinal_type lwork = card*card;
505 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
506 vmat(
"Hcurl::Tet::In::vmat", card, card),
507 work1(
"Hcurl::Tet::In::work", lwork),
508 ipiv(
"Hcurl::Tet::In::ipiv", card);
511 for(ordinal_type i=0; i< card; ++i) {
512 for(ordinal_type j=0; j< card; ++j) {
514 for(ordinal_type k=0; k< cardVecPn; ++k)
522 lapack.GETRF(card, card,
523 vmat.data(), vmat.stride_1(),
524 (ordinal_type*)ipiv.data(),
527 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
529 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM) lapack.GETRF returns nonzero info." );
532 vmat.data(), vmat.stride_1(),
533 (ordinal_type*)ipiv.data(),
537 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
539 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM) lapack.GETRI returns nonzero info." );
541 for (ordinal_type i=0;i<cardVecPn;++i) {
542 for (ordinal_type j=0;j<card;++j){
544 for(ordinal_type k=0; k< card; ++k)
545 s += U(i,k)*vmat(k,j);
550 this->coeffs_ = Kokkos::create_mirror_view(
typename SpT::memory_space(), coeffs);
551 Kokkos::deep_copy(this->coeffs_ , coeffs);
553 this->dofCoords_ = Kokkos::create_mirror_view(
typename SpT::memory_space(), dofCoords);
554 Kokkos::deep_copy(this->dofCoords_, dofCoords);
556 this->dofCoeffs_ = Kokkos::create_mirror_view(
typename SpT::memory_space(), dofCoeffs);
557 Kokkos::deep_copy(this->dofCoeffs_, dofCoeffs);
563 const ordinal_type posScDim = 0;
564 const ordinal_type posScOrd = 1;
565 const ordinal_type posDfOrd = 2;
571 this->setOrdinalTagData(this->tagToOrdinal_,
574 this->basisCardinality_,
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
virtual void getCubature(PointViewType cubPoints, weightViewType cubWeights) const
Returns cubature points and weights (return arrays must be pre-sized/pre-allocated).
Basis_HCURL_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Header file for the Intrepid2::CubatureDirectTetDefault class.
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.