Intrepid2
Intrepid2_HCURL_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_HCURL_TET_IN_FEM_DEF_HPP__
50#define __INTREPID2_HCURL_TET_IN_FEM_DEF_HPP__
51
54#include "Teuchos_SerialDenseMatrix.hpp"
55
56namespace Intrepid2 {
57
58// -------------------------------------------------------------------------------------
59
60namespace Impl {
61
62template<EOperator opType>
63template<typename OutputViewType,
64typename inputViewType,
65typename workViewType,
66typename vinvViewType>
67KOKKOS_INLINE_FUNCTION
68void
69Basis_HCURL_TET_In_FEM::Serial<opType>::
70getValues( OutputViewType output,
71 const inputViewType input,
72 workViewType work,
73 const vinvViewType coeffs ) {
74
75 constexpr ordinal_type spaceDim = 3;
76 const ordinal_type
77 cardPn = coeffs.extent(0)/spaceDim,
78 card = coeffs.extent(1),
79 npts = input.extent(0);
80
81 // compute order
82 ordinal_type order = 0;
83 for (ordinal_type p=0;p<=Parameters::MaxOrder;++p) {
84 if (card == CardinalityHCurlTet(p)) {
85 order = p;
86 break;
87 }
88 }
89
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();
93
94 switch (opType) {
95 case OPERATOR_VALUE: {
96 const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts);
97 workViewType dummyView;
98
99 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::
100 Serial<opType>::getValues(phis, input, dummyView, order);
101
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);
108 }
109 break;
110 }
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);
115
116 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::
117 Serial<OPERATOR_GRAD>::getValues(phis, input, workView, order);
118
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) //\sum_k (coeffs_k, coeffs_{k+cardPn}, coeffs_{k+2 cardPn}) \times phis_kj (cross product)
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);
127 }
128 }
129 }
130 break;
131 }
132 default: {
133 INTREPID2_TEST_FOR_ABORT( true,
134 ">>> ERROR (Basis_HCURL_TET_In_FEM): Operator type not implemented");
135 }
136 }
137}
138
139template<typename DT, ordinal_type numPtsPerEval,
140typename outputValueValueType, class ...outputValueProperties,
141typename inputPointValueType, class ...inputPointProperties,
142typename vinvValueType, class ...vinvProperties>
143void
144Basis_HCURL_TET_In_FEM::
145getValues( 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,typename DT::execution_space>::ExecSpaceType ExecSpaceType;
153
154 // loopSize corresponds to cardinality
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);
159
160 typedef typename inputPointViewType::value_type inputPointType;
161
162 const ordinal_type cardinality = outputValues.extent(0);
163 const ordinal_type spaceDim = 3;
164
165 auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
166 typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
167
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) );
174 break;
175 }
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) );
181 break;
182 }
183 default: {
184 INTREPID2_TEST_FOR_EXCEPTION( true , std::invalid_argument,
185 ">>> ERROR (Basis_HCURL_TET_In_FEM): Operator type not implemented" );
186 }
187 }
188}
189}
190
191// -------------------------------------------------------------------------------------
192template<typename DT, typename OT, typename PT>
194Basis_HCURL_TET_In_FEM( const ordinal_type order,
195 const EPointType pointType ) {
196
197 constexpr ordinal_type spaceDim = 3;
198 this->basisCardinality_ = CardinalityHCurlTet(order);
199 this->basisDegree_ = order; // small n
200 this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Tetrahedron<4> >() );
201 this->basisType_ = BASIS_FEM_LAGRANGIAN;
202 this->basisCoordinates_ = COORDINATES_CARTESIAN;
203 this->functionSpace_ = FUNCTION_SPACE_HCURL;
204 pointType_ = pointType;
205 const ordinal_type card = this->basisCardinality_;
206
207 const ordinal_type cardPn = Intrepid2::getPnCardinality<spaceDim>(order); // dim of (P_{n}) -- smaller space
208 const ordinal_type cardPnm1 = Intrepid2::getPnCardinality<spaceDim>(order-1); // dim of (P_{n-1}) -- smaller space
209 const ordinal_type cardPnm2 = Intrepid2::getPnCardinality<spaceDim>(order-2); // dim of (P_{n-2}) -- smaller space
210 const ordinal_type cardVecPn = spaceDim*cardPn; // dim of (P_{n})^2 -- larger space
211 const ordinal_type cardVecPnm1 = spaceDim*cardPnm1; // dim of (P_{n-1})^2 -- smaller space
212 const ordinal_type cardPnm1H = cardPnm1-cardPnm2; //Homogeneous polynomial of order (n-1)
213
214
215
216 // Basis-dependent initializations
217 constexpr ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
218 constexpr ordinal_type maxCard = CardinalityHCurlTet(Parameters::MaxOrder);
219 ordinal_type tags[maxCard][tagSize];
220
221 // points are computed in the host and will be copied
222 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
223 dofCoords("Hcurl::Tet::In::dofCoords", card, spaceDim);
224
225 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
226 coeffs("Hcurl::Tet::In::coeffs", cardVecPn, card);
227
228 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
229 dofCoeffs("Hcurl::Tet::In::dofCoeffs", card, spaceDim);
230
231 // first, need to project the basis for RT space onto the
232 // orthogonal basis of degree n
233 // get coefficients of PkHx
234
235 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace> //use LayoutLeft for Lapack
236 V1("Hcurl::Tet::In::V1", cardVecPn, cardVecPnm1 + spaceDim*cardPnm1H);
237
238
239 // these two loops get the first three sets of basis functions
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;
243
244
245 // now I need to integrate { (x,y) \times phi } against the big basis
246 // first, get a cubature rule.
248 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> cubPoints("Hcurl::Tet::In::cubPoints", myCub.getNumPoints() , spaceDim );
249 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> cubWeights("Hcurl::Tet::In::cubWeights", myCub.getNumPoints() );
250 myCub.getCubature( cubPoints , cubWeights );
251
252 // tabulate the scalar orthonormal basis at cubature points
253 Kokkos::DynRankView<scalarType,typename DT::execution_space::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);
255
256 // Integrate (x psi_j, y psi_j, z psi_j) \times (phi_i, phi_{i+cardPn}, phi_{i+2 cardPn}) cross product. psi are homogeneous polynomials of order (n-1)
257 for (ordinal_type i=0;i<cardPn;i++) {
258 for (ordinal_type j=0;j<cardPnm1H;j++) { // loop over homogeneous polynomials
259 for (ordinal_type d=0; d< spaceDim; ++d) {
260 scalarType integral(0);
261 for (ordinal_type k=0;k<myCub.getNumPoints();k++)
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;
268 }
269 }
270 }
271
272
273
274
275
276 // now I need to set up an SVD to get a basis for the space
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);
283
284
285
286 ordinal_type info = 0;
287 Teuchos::LAPACK<ordinal_type,scalarType> lapack;
288
289
290 lapack.GESVD( 'A',
291 'N',
292 V1.extent(0) ,
293 V1.extent(1) ,
294 V1.data() ,
295 V1.stride_1() ,
296 S.data() ,
297 U.data() ,
298 U.stride_1() ,
299 Vt.data() ,
300 Vt.stride_1() ,
301 work.data() ,
302 5*cardVecPn ,
303 rWork.data() ,
304 &info );
305
306
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());
311
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");
314#endif
315
316 // next, apply the RT nodes (rows) to the basis for (P_n)^2 (columns)
317 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
318 V2("Hcurl::Tet::In::V2", card ,cardVecPn);
319
320 const ordinal_type numEdges = this->basisCellTopology_.getEdgeCount();
321 const ordinal_type numFaces = this->basisCellTopology_.getFaceCount();
322
323 // first numEdges * degree nodes are normals at each edge
324 // get the points on the line
325
326 shards::CellTopology edgeTop(shards::getCellTopologyData<shards::Line<2> >() );
327 shards::CellTopology faceTop(shards::getCellTopologyData<shards::Triangle<3> >() );
328
329 const int numPtsPerEdge = PointTools::getLatticeSize( edgeTop ,
330 order+1 ,
331 1 );
332
333 const int numPtsPerFace = PointTools::getLatticeSize( faceTop ,
334 order+1 ,
335 1 );
336
337 const int numPtsPerCell = PointTools::getLatticeSize( this->basisCellTopology_ ,
338 order+1 ,
339 1 );
340
341 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> linePts("Hcurl::Tet::In::linePts", numPtsPerEdge , 1 );
342 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> triPts("Hcurl::Tet::In::triPts", numPtsPerFace , 2 );
343
344 // construct lattice
345 const ordinal_type offset = 1;
346
347
348
349 PointTools::getLattice( linePts,
350 edgeTop,
351 order+1, offset,
352 pointType );
353
355 faceTop,
356 order+1, offset,
357 pointType );
358
359 // holds the image of the line points
360 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> edgePts("Hcurl::Tet::In::edgePts", numPtsPerEdge , spaceDim );
361 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> facePts("Hcurl::Tet::In::facePts", numPtsPerFace , spaceDim );
362 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> phisAtEdgePoints("Hcurl::Tet::In::phisAtEdgePoints", cardPn , numPtsPerEdge );
363 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> phisAtFacePoints("Hcurl::Tet::In::phisAtFacePoints", cardPn , numPtsPerFace);
364
365 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> edgeTan("Hcurl::Tet::In::edgeTan", spaceDim );
366
367 // these are tangents scaled by the appropriate edge lengths.
368 for (ordinal_type i=0;i<numEdges;i++) { // loop over edges
370 i ,
371 this->basisCellTopology_ );
372
374 linePts ,
375 1 ,
376 i ,
377 this->basisCellTopology_);
378
379 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace,Parameters::MaxNumPtsPerBasisEval>(phisAtEdgePoints , edgePts, order, OPERATOR_VALUE);
380
381 // loop over points (rows of V2)
382 for (ordinal_type j=0;j<numPtsPerEdge;j++) {
383
384 const ordinal_type i_card = numPtsPerEdge*i+j;
385
386 // loop over orthonormal basis functions (columns of V2)
387 for (ordinal_type k=0;k<cardPn;k++)
388 for (ordinal_type d=0;d<spaceDim;d++)
389 V2(i_card,k+d*cardPn) = edgeTan(d) * phisAtEdgePoints(k,j);
390
391 //save dof coordinates and coefficients
392 for(ordinal_type k=0; k<spaceDim; ++k) {
393 dofCoords(i_card,k) = edgePts(j,k);
394 dofCoeffs(i_card,k) = edgeTan(k);
395 }
396
397 tags[i_card][0] = 1; // edge dof
398 tags[i_card][1] = i; // edge id
399 tags[i_card][2] = j; // local dof id
400 tags[i_card][3] = numPtsPerEdge; // total vert dof
401
402 }
403 }
404
405 if(numPtsPerFace >0) {//handle faces if needed (order >1)
406 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> faceTan1("Hcurl::Tet::In::edgeTan", spaceDim );
407 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> faceTan2("Hcurl::Tet::In::edgeTan", spaceDim );
408
409 for (ordinal_type i=0;i<numFaces;i++) { // loop over faces
411 faceTan2,
412 i ,
413 this->basisCellTopology_ );
414
416 triPts ,
417 2 ,
418 i ,
419 this->basisCellTopology_ );
420
421 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace,Parameters::MaxNumPtsPerBasisEval>(phisAtFacePoints , facePts, order, OPERATOR_VALUE);
422
423 // loop over points (rows of V2)
424 for (ordinal_type j=0;j<numPtsPerFace;j++) {
425
426 const ordinal_type i_card = numEdges*numPtsPerEdge+2*numPtsPerFace*i+2*j;
427 const ordinal_type i_card_p1 = i_card+1; // creating a temp otherwise nvcc gets confused
428
429 // loop over orthonormal basis functions (columns of V2)
430 for (ordinal_type k=0;k<cardPn;k++)
431 for (ordinal_type d=0;d<spaceDim;d++) {
432 V2(i_card,k+d*cardPn) = faceTan1(d) * phisAtFacePoints(k,j);
433 V2(i_card_p1,k+d*cardPn) = faceTan2(d) * phisAtFacePoints(k,j);
434 }
435
436 //save dof coordinates
437 for(ordinal_type k=0; k<spaceDim; ++k) {
438 dofCoords(i_card,k) = facePts(j,k);
439 dofCoords(i_card_p1,k) = facePts(j,k);
440 dofCoeffs(i_card,k) = faceTan1(k);
441 dofCoeffs(i_card_p1,k) = faceTan2(k);
442 }
443
444 tags[i_card][0] = 2; // face dof
445 tags[i_card][1] = i; // face id
446 tags[i_card][2] = 2*j; // local face id
447 tags[i_card][3] = 2*numPtsPerFace; // total face dof
448
449 tags[i_card_p1][0] = 2; // face dof
450 tags[i_card_p1][1] = i; // face id
451 tags[i_card_p1][2] = 2*j+1; // local face id
452 tags[i_card_p1][3] = 2*numPtsPerFace; // total face dof
453
454 }
455 }
456 }
457
458
459 // internal dof, if needed
460 if (numPtsPerCell > 0) {
461 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
462 cellPoints( "Hcurl::Tet::In::cellPoints", numPtsPerCell , spaceDim );
463 PointTools::getLattice( cellPoints ,
464 this->basisCellTopology_ ,
465 order + 1 ,
466 1 ,
467 pointType );
468
469 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
470 phisAtCellPoints("Hcurl::Tet::In::phisAtCellPoints", cardPn , numPtsPerCell );
471 Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>( phisAtCellPoints , cellPoints , order, OPERATOR_VALUE );
472
473 // copy values into right positions of V2
474 for (ordinal_type j=0;j<numPtsPerCell;j++) {
475
476 const ordinal_type i_card = numEdges*numPtsPerEdge+2*numFaces*numPtsPerFace+spaceDim*j;
477
478 for (ordinal_type k=0;k<cardPn;k++)
479 for (ordinal_type d=0;d<spaceDim;d++)
480 V2(i_card+d,d*cardPn+k) = phisAtCellPoints(k,j);
481
482
483 //save dof coordinates
484 for(ordinal_type d=0; d<spaceDim; ++d) {
485 for(ordinal_type dim=0; dim<spaceDim; ++dim) {
486 dofCoords(i_card+d,dim) = cellPoints(j,dim);
487 dofCoeffs(i_card+d,dim) = (d==dim);
488 }
489
490 tags[i_card+d][0] = spaceDim; // elem dof
491 tags[i_card+d][1] = 0; // elem id
492 tags[i_card+d][2] = spaceDim*j+d; // local dof id
493 tags[i_card+d][3] = spaceDim*numPtsPerCell; // total vert dof
494 }
495 }
496 }
497
498 // form Vandermonde matrix. Actually, this is the transpose of the VDM,
499 // so we transpose on copy below.
500 const ordinal_type lwork = card*card;
501 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
502 vmat("Hcurl::Tet::In::vmat", card, card),
503 work1("Hcurl::Tet::In::work", lwork),
504 ipiv("Hcurl::Tet::In::ipiv", card);
505
506 //vmat = V2*U;
507 for(ordinal_type i=0; i< card; ++i) {
508 for(ordinal_type j=0; j< card; ++j) {
509 scalarType s=0;
510 for(ordinal_type k=0; k< cardVecPn; ++k)
511 s += V2(i,k)*U(k,j);
512 vmat(i,j) = s;
513 }
514 }
515
516 info = 0;
517
518 lapack.GETRF(card, card,
519 vmat.data(), vmat.stride_1(),
520 (ordinal_type*)ipiv.data(),
521 &info);
522
523 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
524 std::runtime_error ,
525 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM) lapack.GETRF returns nonzero info." );
526
527 lapack.GETRI(card,
528 vmat.data(), vmat.stride_1(),
529 (ordinal_type*)ipiv.data(),
530 work1.data(), lwork,
531 &info);
532
533 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
534 std::runtime_error ,
535 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM) lapack.GETRI returns nonzero info." );
536
537 for (ordinal_type i=0;i<cardVecPn;++i) {
538 for (ordinal_type j=0;j<card;++j){
539 scalarType s=0;
540 for(ordinal_type k=0; k< card; ++k)
541 s += U(i,k)*vmat(k,j);
542 coeffs(i,j) = s;
543 }
544 }
545
546 this->coeffs_ = Kokkos::create_mirror_view(typename DT::memory_space(), coeffs);
547 Kokkos::deep_copy(this->coeffs_ , coeffs);
548
549 this->dofCoords_ = Kokkos::create_mirror_view(typename DT::memory_space(), dofCoords);
550 Kokkos::deep_copy(this->dofCoords_, dofCoords);
551
552 this->dofCoeffs_ = Kokkos::create_mirror_view(typename DT::memory_space(), dofCoeffs);
553 Kokkos::deep_copy(this->dofCoeffs_, dofCoeffs);
554
555
556 // set tags
557 {
558 // Basis-dependent initializations
559 const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
560 const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
561 const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
562
563 OrdinalTypeArray1DHost tagView(&tags[0][0], card*tagSize);
564
565 // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
566 // tags are constructed on host
567 this->setOrdinalTagData(this->tagToOrdinal_,
568 this->ordinalToTag_,
569 tagView,
570 this->basisCardinality_,
571 tagSize,
572 posScDim,
573 posScOrd,
574 posDfOrd);
575 }
576}
577} // namespace Intrepid2
578#endif
Header file for the Intrepid2::CubatureDirectTetDefault class.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
Basis_HCURL_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
static void getReferenceEdgeTangent(Kokkos::DynRankView< refEdgeTangentValueType, refEdgeTangentProperties... > refEdgeTangent, const ordinal_type edgeOrd, const shards::CellTopology parentCell)
Computes constant tangent vectors to edges of 2D or 3D reference cells.
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.
static void getReferenceFaceTangents(Kokkos::DynRankView< refFaceTanValueType, refFaceTanProperties... > refFaceTanU, Kokkos::DynRankView< refFaceTanValueType, refFaceTanProperties... > refFaceTanV, const ordinal_type faceOrd, const shards::CellTopology parentCell)
Computes pairs of constant tangent vectors to faces of a 3D reference cells.
Defines direct integration rules on a tetrahedron.
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...