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 DT, ordinal_type numPtsPerEval,
140 typename outputValueValueType,
class ...outputValueProperties,
141 typename inputPointValueType,
class ...inputPointProperties,
142 typename vinvValueType,
class ...vinvProperties>
144 Basis_HCURL_TET_In_FEM::
146 const typename DT::execution_space& space,
147 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
148 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
149 const Kokkos::DynRankView<vinvValueType, vinvProperties...> coeffs,
150 const EOperator operatorType) {
151 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
152 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
153 typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
154 typedef typename ExecSpace<typename inputPointViewType::execution_space,typename DT::execution_space>::ExecSpaceType ExecSpaceType;
157 const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
158 const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
159 const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
160 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(space, 0, loopSize);
162 typedef typename inputPointViewType::value_type inputPointType;
164 const ordinal_type cardinality = outputValues.extent(0);
165 const ordinal_type spaceDim = 3;
167 auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
168 typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
170 switch (operatorType) {
171 case OPERATOR_VALUE: {
172 workViewType work(Kokkos::view_alloc(space,
"Basis_HCURL_TET_In_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0));
173 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
174 OPERATOR_VALUE,numPtsPerEval> FunctorType;
175 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, coeffs, work) );
178 case OPERATOR_CURL: {
179 workViewType work(Kokkos::view_alloc(space,
"Basis_HCURL_TET_In_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
180 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
181 OPERATOR_CURL,numPtsPerEval> FunctorType;
182 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, coeffs, work) );
186 INTREPID2_TEST_FOR_EXCEPTION(
true , std::invalid_argument,
187 ">>> ERROR (Basis_HCURL_TET_In_FEM): Operator type not implemented" );
194 template<
typename DT,
typename OT,
typename PT>
197 const EPointType pointType ) {
199 constexpr ordinal_type spaceDim = 3;
200 this->basisCardinality_ = CardinalityHCurlTet(order);
201 this->basisDegree_ = order;
202 this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Tetrahedron<4> >() );
203 this->basisType_ = BASIS_FEM_LAGRANGIAN;
204 this->basisCoordinates_ = COORDINATES_CARTESIAN;
205 this->functionSpace_ = FUNCTION_SPACE_HCURL;
206 pointType_ = pointType;
207 const ordinal_type card = this->basisCardinality_;
209 const ordinal_type cardPn = Intrepid2::getPnCardinality<spaceDim>(order);
210 const ordinal_type cardPnm1 = Intrepid2::getPnCardinality<spaceDim>(order-1);
211 const ordinal_type cardPnm2 = Intrepid2::getPnCardinality<spaceDim>(order-2);
212 const ordinal_type cardVecPn = spaceDim*cardPn;
213 const ordinal_type cardVecPnm1 = spaceDim*cardPnm1;
214 const ordinal_type cardPnm1H = cardPnm1-cardPnm2;
218 INTREPID2_TEST_FOR_EXCEPTION( order >
Parameters::MaxOrder, std::invalid_argument,
"polynomial order exceeds the max supported by this class");
221 constexpr ordinal_type tagSize = 4;
223 ordinal_type tags[maxCard][tagSize];
226 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
227 dofCoords(
"Hcurl::Tet::In::dofCoords", card, spaceDim);
229 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
230 coeffs(
"Hcurl::Tet::In::coeffs", cardVecPn, card);
232 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
233 dofCoeffs(
"Hcurl::Tet::In::dofCoeffs", card, spaceDim);
239 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
240 V1(
"Hcurl::Tet::In::V1", cardVecPn, cardVecPnm1 + spaceDim*cardPnm1H);
244 for (ordinal_type i=0;i<cardPnm1;i++)
245 for (ordinal_type d=0;d<spaceDim;d++)
246 V1(i+d*cardPn,i+d*cardPnm1) = 1.0;
252 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> cubPoints(
"Hcurl::Tet::In::cubPoints", myCub.
getNumPoints() , spaceDim );
253 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> cubWeights(
"Hcurl::Tet::In::cubWeights", myCub.
getNumPoints() );
257 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> phisAtCubPoints(
"Hcurl::Tet::In::phisAtCubPoints", cardPn , myCub.
getNumPoints() );
258 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(
typename Kokkos::HostSpace::execution_space{},
265 for (ordinal_type i=0;i<cardPn;i++) {
266 for (ordinal_type j=0;j<cardPnm1H;j++) {
267 for (ordinal_type d=0; d< spaceDim; ++d) {
268 scalarType integral(0);
270 integral += cubWeights(k) * cubPoints(k,d)
271 * phisAtCubPoints(cardPnm2+j,k)
272 * phisAtCubPoints(i,k);
273 ordinal_type d1 = (d+1) % spaceDim, d2 = (d+2) % spaceDim;
274 V1(i+d2*cardPn,cardVecPnm1+d1*cardPnm1H + j) = -integral;
275 V1(i+d1*cardPn,cardVecPnm1+d2*cardPnm1H + j) = integral;
285 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
286 S(
"Hcurl::Tet::In::S", cardVecPn,1),
287 U(
"Hcurl::Tet::In::U", cardVecPn, cardVecPn),
288 Vt(
"Hcurl::Tet::In::Vt", cardVecPn, cardVecPn),
289 work(
"Hcurl::Tet::In::work", 5*cardVecPn,1),
290 rWork(
"Hcurl::Tet::In::rW", 1,1);
294 ordinal_type info = 0;
295 Teuchos::LAPACK<ordinal_type,scalarType> lapack;
315 #ifdef HAVE_INTREPID2_DEBUG
316 ordinal_type num_nonzero_sv = 0;
317 for (
int i=0;i<cardVecPn;i++)
318 num_nonzero_sv += (S(i,0) > tolerence());
320 INTREPID2_TEST_FOR_EXCEPTION( num_nonzero_sv != card, std::invalid_argument,
321 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM( order, pointType), Matrix V1 should have rank equal to the cardinality of HCURL space");
325 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
326 V2(
"Hcurl::Tet::In::V2", card ,cardVecPn);
328 const ordinal_type numEdges = this->basisCellTopology_.getEdgeCount();
329 const ordinal_type numFaces = this->basisCellTopology_.getFaceCount();
334 shards::CellTopology edgeTop(shards::getCellTopologyData<shards::Line<2> >() );
335 shards::CellTopology faceTop(shards::getCellTopologyData<shards::Triangle<3> >() );
349 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> linePts(
"Hcurl::Tet::In::linePts", numPtsPerEdge , 1 );
350 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> triPts(
"Hcurl::Tet::In::triPts", numPtsPerFace , 2 );
353 const ordinal_type offset = 1;
368 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> edgePts(
"Hcurl::Tet::In::edgePts", numPtsPerEdge , spaceDim );
369 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> facePts(
"Hcurl::Tet::In::facePts", numPtsPerFace , spaceDim );
370 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> phisAtEdgePoints(
"Hcurl::Tet::In::phisAtEdgePoints", cardPn , numPtsPerEdge );
371 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> phisAtFacePoints(
"Hcurl::Tet::In::phisAtFacePoints", cardPn , numPtsPerFace);
373 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> edgeTan(
"Hcurl::Tet::In::edgeTan", spaceDim );
376 for (ordinal_type i=0;i<numEdges;i++) {
379 this->basisCellTopology_ );
385 this->basisCellTopology_);
387 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(
typename Kokkos::HostSpace::execution_space{},
394 for (ordinal_type j=0;j<numPtsPerEdge;j++) {
396 const ordinal_type i_card = numPtsPerEdge*i+j;
399 for (ordinal_type k=0;k<cardPn;k++)
400 for (ordinal_type d=0;d<spaceDim;d++)
401 V2(i_card,k+d*cardPn) = edgeTan(d) * phisAtEdgePoints(k,j);
404 for(ordinal_type k=0; k<spaceDim; ++k) {
405 dofCoords(i_card,k) = edgePts(j,k);
406 dofCoeffs(i_card,k) = edgeTan(k);
412 tags[i_card][3] = numPtsPerEdge;
417 if(numPtsPerFace >0) {
418 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> faceTan1(
"Hcurl::Tet::In::edgeTan", spaceDim );
419 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> faceTan2(
"Hcurl::Tet::In::edgeTan", spaceDim );
421 for (ordinal_type i=0;i<numFaces;i++) {
425 this->basisCellTopology_ );
431 this->basisCellTopology_ );
433 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(
typename Kokkos::HostSpace::execution_space{},
440 for (ordinal_type j=0;j<numPtsPerFace;j++) {
442 const ordinal_type i_card = numEdges*numPtsPerEdge+2*numPtsPerFace*i+2*j;
443 const ordinal_type i_card_p1 = i_card+1;
446 for (ordinal_type k=0;k<cardPn;k++)
447 for (ordinal_type d=0;d<spaceDim;d++) {
448 V2(i_card,k+d*cardPn) = faceTan1(d) * phisAtFacePoints(k,j);
449 V2(i_card_p1,k+d*cardPn) = faceTan2(d) * phisAtFacePoints(k,j);
453 for(ordinal_type k=0; k<spaceDim; ++k) {
454 dofCoords(i_card,k) = facePts(j,k);
455 dofCoords(i_card_p1,k) = facePts(j,k);
456 dofCoeffs(i_card,k) = faceTan1(k);
457 dofCoeffs(i_card_p1,k) = faceTan2(k);
462 tags[i_card][2] = 2*j;
463 tags[i_card][3] = 2*numPtsPerFace;
465 tags[i_card_p1][0] = 2;
466 tags[i_card_p1][1] = i;
467 tags[i_card_p1][2] = 2*j+1;
468 tags[i_card_p1][3] = 2*numPtsPerFace;
476 if (numPtsPerCell > 0) {
477 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
478 cellPoints(
"Hcurl::Tet::In::cellPoints", numPtsPerCell , spaceDim );
480 this->basisCellTopology_ ,
485 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
486 phisAtCellPoints(
"Hcurl::Tet::In::phisAtCellPoints", cardPn , numPtsPerCell );
487 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(
typename Kokkos::HostSpace::execution_space{},
494 for (ordinal_type j=0;j<numPtsPerCell;j++) {
496 const ordinal_type i_card = numEdges*numPtsPerEdge+2*numFaces*numPtsPerFace+spaceDim*j;
498 for (ordinal_type k=0;k<cardPn;k++)
499 for (ordinal_type d=0;d<spaceDim;d++)
500 V2(i_card+d,d*cardPn+k) = phisAtCellPoints(k,j);
504 for(ordinal_type d=0; d<spaceDim; ++d) {
505 for(ordinal_type dim=0; dim<spaceDim; ++dim) {
506 dofCoords(i_card+d,dim) = cellPoints(j,dim);
507 dofCoeffs(i_card+d,dim) = (d==dim);
510 tags[i_card+d][0] = spaceDim;
511 tags[i_card+d][1] = 0;
512 tags[i_card+d][2] = spaceDim*j+d;
513 tags[i_card+d][3] = spaceDim*numPtsPerCell;
520 const ordinal_type lwork = card*card;
521 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
522 vmat(
"Hcurl::Tet::In::vmat", card, card),
523 work1(
"Hcurl::Tet::In::work", lwork),
524 ipiv(
"Hcurl::Tet::In::ipiv", card);
527 for(ordinal_type i=0; i< card; ++i) {
528 for(ordinal_type j=0; j< card; ++j) {
530 for(ordinal_type k=0; k< cardVecPn; ++k)
538 lapack.GETRF(card, card,
539 vmat.data(), vmat.stride_1(),
540 (ordinal_type*)ipiv.data(),
543 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
545 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM) lapack.GETRF returns nonzero info." );
548 vmat.data(), vmat.stride_1(),
549 (ordinal_type*)ipiv.data(),
553 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
555 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM) lapack.GETRI returns nonzero info." );
557 for (ordinal_type i=0;i<cardVecPn;++i) {
558 for (ordinal_type j=0;j<card;++j){
560 for(ordinal_type k=0; k< card; ++k)
561 s += U(i,k)*vmat(k,j);
566 this->coeffs_ = Kokkos::create_mirror_view(
typename DT::memory_space(), coeffs);
567 Kokkos::deep_copy(this->coeffs_ , coeffs);
569 this->dofCoords_ = Kokkos::create_mirror_view(
typename DT::memory_space(), dofCoords);
570 Kokkos::deep_copy(this->dofCoords_, dofCoords);
572 this->dofCoeffs_ = Kokkos::create_mirror_view(
typename DT::memory_space(), dofCoeffs);
573 Kokkos::deep_copy(this->dofCoeffs_, dofCoeffs);
579 const ordinal_type posScDim = 0;
580 const ordinal_type posScOrd = 1;
581 const ordinal_type posDfOrd = 2;
587 this->setOrdinalTagData(this->tagToOrdinal_,
590 this->basisCardinality_,
virtual void getCubature(PointViewType cubPoints, weightViewType cubWeights) const override
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.
virtual ordinal_type getNumPoints() const override
Returns the number of cubature points.
Header file for the Intrepid2::CubatureDirectTetDefault class.
Defines direct integration rules on a tetrahedron.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.