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;
204 const ordinal_type card = this->basisCardinality_;
206 const ordinal_type cardPn = Intrepid2::getPnCardinality<spaceDim>(order);
207 const ordinal_type cardPnm1 = Intrepid2::getPnCardinality<spaceDim>(order-1);
208 const ordinal_type cardPnm2 = Intrepid2::getPnCardinality<spaceDim>(order-2);
209 const ordinal_type cardVecPn = spaceDim*cardPn;
210 const ordinal_type cardVecPnm1 = spaceDim*cardPnm1;
211 const ordinal_type cardPnm1H = cardPnm1-cardPnm2;
216 constexpr ordinal_type tagSize = 4;
218 ordinal_type tags[maxCard][tagSize];
221 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
222 dofCoords(
"Hcurl::Tet::In::dofCoords", card, spaceDim);
224 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
225 coeffs(
"Hcurl::Tet::In::coeffs", cardVecPn, card);
227 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
228 dofCoeffs(
"Hcurl::Tet::In::dofCoeffs", card, spaceDim);
234 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
235 V1(
"Hcurl::Tet::In::V1", cardVecPn, cardVecPnm1 + spaceDim*cardPnm1H);
239 for (ordinal_type i=0;i<cardPnm1;i++)
240 for (ordinal_type d=0;d<spaceDim;d++)
241 V1(i+d*cardPn,i+d*cardPnm1) = 1.0;
247 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> cubPoints(
"Hcurl::Tet::In::cubPoints", myCub.
getNumPoints() , spaceDim );
248 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> cubWeights(
"Hcurl::Tet::In::cubWeights", myCub.
getNumPoints() );
252 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> phisAtCubPoints(
"Hcurl::Tet::In::phisAtCubPoints", cardPn , myCub.
getNumPoints() );
253 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(phisAtCubPoints, cubPoints, order, OPERATOR_VALUE);
256 for (ordinal_type i=0;i<cardPn;i++) {
257 for (ordinal_type j=0;j<cardPnm1H;j++) {
258 for (ordinal_type d=0; d< spaceDim; ++d) {
259 scalarType integral(0);
261 integral += cubWeights(k) * cubPoints(k,d)
262 * phisAtCubPoints(cardPnm2+j,k)
263 * phisAtCubPoints(i,k);
264 ordinal_type d1 = (d+1) % spaceDim, d2 = (d+2) % spaceDim;
265 V1(i+d2*cardPn,cardVecPnm1+d1*cardPnm1H + j) = -integral;
266 V1(i+d1*cardPn,cardVecPnm1+d2*cardPnm1H + j) = integral;
276 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
277 S(
"Hcurl::Tet::In::S", cardVecPn,1),
278 U(
"Hcurl::Tet::In::U", cardVecPn, cardVecPn),
279 Vt(
"Hcurl::Tet::In::Vt", cardVecPn, cardVecPn),
280 work(
"Hcurl::Tet::In::work", 5*cardVecPn,1),
281 rWork(
"Hcurl::Tet::In::rW", 1,1);
285 ordinal_type info = 0;
286 Teuchos::LAPACK<ordinal_type,scalarType> lapack;
306 #ifdef HAVE_INTREPID2_DEBUG
307 ordinal_type num_nonzero_sv = 0;
308 for (
int i=0;i<cardVecPn;i++)
309 num_nonzero_sv += (S(i,0) > tolerence());
311 INTREPID2_TEST_FOR_EXCEPTION( num_nonzero_sv != card, std::invalid_argument,
312 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM( order, pointType), Matrix V1 should have rank equal to the cardinality of HCURL space");
316 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
317 V2(
"Hcurl::Tet::In::V2", card ,cardVecPn);
319 const ordinal_type numEdges = this->basisCellTopology_.getEdgeCount();
320 const ordinal_type numFaces = this->basisCellTopology_.getFaceCount();
325 shards::CellTopology edgeTop(shards::getCellTopologyData<shards::Line<2> >() );
326 shards::CellTopology faceTop(shards::getCellTopologyData<shards::Triangle<3> >() );
340 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> linePts(
"Hcurl::Tet::In::linePts", numPtsPerEdge , 1 );
341 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> triPts(
"Hcurl::Tet::In::triPts", numPtsPerFace , 2 );
344 const ordinal_type offset = 1;
359 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> edgePts(
"Hcurl::Tet::In::edgePts", numPtsPerEdge , spaceDim );
360 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> facePts(
"Hcurl::Tet::In::facePts", numPtsPerFace , spaceDim );
361 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> phisAtEdgePoints(
"Hcurl::Tet::In::phisAtEdgePoints", cardPn , numPtsPerEdge );
362 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> phisAtFacePoints(
"Hcurl::Tet::In::phisAtFacePoints", cardPn , numPtsPerFace);
364 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> edgeTan(
"Hcurl::Tet::In::edgeTan", spaceDim );
367 for (ordinal_type i=0;i<numEdges;i++) {
370 this->basisCellTopology_ );
372 const scalarType refEdgeMeasure = 2.0;
373 for (ordinal_type j=0;j<spaceDim;j++)
374 edgeTan(j) *= refEdgeMeasure;
380 this->basisCellTopology_ );
382 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(phisAtEdgePoints , edgePts, order, OPERATOR_VALUE);
385 for (ordinal_type j=0;j<numPtsPerEdge;j++) {
387 const ordinal_type i_card = numPtsPerEdge*i+j;
390 for (ordinal_type k=0;k<cardPn;k++)
391 for (ordinal_type d=0;d<spaceDim;d++)
392 V2(i_card,k+d*cardPn) = edgeTan(d) * phisAtEdgePoints(k,j);
395 for(ordinal_type k=0; k<spaceDim; ++k) {
396 dofCoords(i_card,k) = edgePts(j,k);
397 dofCoeffs(i_card,k) = edgeTan(k);
403 tags[i_card][3] = numPtsPerEdge;
408 if(numPtsPerFace >0) {
409 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> faceTan1(
"Hcurl::Tet::In::edgeTan", spaceDim );
410 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace> faceTan2(
"Hcurl::Tet::In::edgeTan", spaceDim );
412 for (ordinal_type i=0;i<numFaces;i++) {
416 this->basisCellTopology_ );
422 this->basisCellTopology_ );
424 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(phisAtFacePoints , facePts, order, OPERATOR_VALUE);
427 for (ordinal_type j=0;j<numPtsPerFace;j++) {
429 const ordinal_type i_card = numEdges*numPtsPerEdge+2*numPtsPerFace*i+2*j;
430 const ordinal_type i_card_p1 = i_card+1;
433 for (ordinal_type k=0;k<cardPn;k++)
434 for (ordinal_type d=0;d<spaceDim;d++) {
435 V2(i_card,k+d*cardPn) = faceTan1(d) * phisAtFacePoints(k,j);
436 V2(i_card_p1,k+d*cardPn) = faceTan2(d) * phisAtFacePoints(k,j);
440 for(ordinal_type k=0; k<spaceDim; ++k) {
441 dofCoords(i_card,k) = facePts(j,k);
442 dofCoords(i_card_p1,k) = facePts(j,k);
443 dofCoeffs(i_card,k) = faceTan1(k);
444 dofCoeffs(i_card_p1,k) = faceTan2(k);
449 tags[i_card][2] = 2*j;
450 tags[i_card][3] = 2*numPtsPerFace;
452 tags[i_card_p1][0] = 2;
453 tags[i_card_p1][1] = i;
454 tags[i_card_p1][2] = 2*j+1;
455 tags[i_card_p1][3] = 2*numPtsPerFace;
463 if (numPtsPerCell > 0) {
464 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
465 cellPoints(
"Hcurl::Tet::In::cellPoints", numPtsPerCell , spaceDim );
467 this->basisCellTopology_ ,
472 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
473 phisAtCellPoints(
"Hcurl::Tet::In::phisAtCellPoints", cardPn , numPtsPerCell );
474 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>( phisAtCellPoints , cellPoints , order, OPERATOR_VALUE );
477 for (ordinal_type j=0;j<numPtsPerCell;j++) {
479 const ordinal_type i_card = numEdges*numPtsPerEdge+2*numFaces*numPtsPerFace+spaceDim*j;
481 for (ordinal_type k=0;k<cardPn;k++)
482 for (ordinal_type d=0;d<spaceDim;d++)
483 V2(i_card+d,d*cardPn+k) = phisAtCellPoints(k,j);
487 for(ordinal_type d=0; d<spaceDim; ++d) {
488 for(ordinal_type dim=0; dim<spaceDim; ++dim) {
489 dofCoords(i_card+d,dim) = cellPoints(j,dim);
490 dofCoeffs(i_card+d,dim) = (d==dim);
493 tags[i_card+d][0] = spaceDim;
494 tags[i_card+d][1] = 0;
495 tags[i_card+d][2] = spaceDim*j+d;
496 tags[i_card+d][3] = spaceDim*numPtsPerCell;
503 const ordinal_type lwork = card*card;
504 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
505 vmat(
"Hcurl::Tet::In::vmat", card, card),
506 work1(
"Hcurl::Tet::In::work", lwork),
507 ipiv(
"Hcurl::Tet::In::ipiv", card);
510 for(ordinal_type i=0; i< card; ++i) {
511 for(ordinal_type j=0; j< card; ++j) {
513 for(ordinal_type k=0; k< cardVecPn; ++k)
521 lapack.GETRF(card, card,
522 vmat.data(), vmat.stride_1(),
523 (ordinal_type*)ipiv.data(),
526 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
528 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM) lapack.GETRF returns nonzero info." );
531 vmat.data(), vmat.stride_1(),
532 (ordinal_type*)ipiv.data(),
536 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
538 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM) lapack.GETRI returns nonzero info." );
540 for (ordinal_type i=0;i<cardVecPn;++i) {
541 for (ordinal_type j=0;j<card;++j){
543 for(ordinal_type k=0; k< card; ++k)
544 s += U(i,k)*vmat(k,j);
549 this->coeffs_ = Kokkos::create_mirror_view(
typename SpT::memory_space(), coeffs);
550 Kokkos::deep_copy(this->coeffs_ , coeffs);
552 this->dofCoords_ = Kokkos::create_mirror_view(
typename SpT::memory_space(), dofCoords);
553 Kokkos::deep_copy(this->dofCoords_, dofCoords);
555 this->dofCoeffs_ = Kokkos::create_mirror_view(
typename SpT::memory_space(), dofCoeffs);
556 Kokkos::deep_copy(this->dofCoeffs_, dofCoeffs);
562 const ordinal_type posScDim = 0;
563 const ordinal_type posScOrd = 1;
564 const ordinal_type posDfOrd = 2;
570 this->setOrdinalTagData(this->tagToOrdinal_,
573 this->basisCardinality_,
Basis_HCURL_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Kokkos::View< ordinal_type *,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_1d_host
View type for 1d host array.
Header file for the Intrepid2::CubatureDirectTetDefault class.
virtual void getCubature(pointViewType cubPoints, weightViewType cubWeights) const
Returns cubature points and weights (return arrays must be pre-sized/pre-allocated).
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.