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 
56 namespace Intrepid2 {
57 
58 // -------------------------------------------------------------------------------------
59 
60 namespace Impl {
61 
62 template<EOperator opType>
63 template<typename OutputViewType,
64 typename inputViewType,
65 typename workViewType,
66 typename vinvViewType>
67 KOKKOS_INLINE_FUNCTION
68 void
69 Basis_HCURL_TET_In_FEM::Serial<opType>::
70 getValues( 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 
139 template<typename DT, ordinal_type numPtsPerEval,
140 typename outputValueValueType, class ...outputValueProperties,
141 typename inputPointValueType, class ...inputPointProperties,
142 typename vinvValueType, class ...vinvProperties>
143 void
144 Basis_HCURL_TET_In_FEM::
145 getValues(
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;
155 
156  // loopSize corresponds to cardinality
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);
161 
162  typedef typename inputPointViewType::value_type inputPointType;
163 
164  const ordinal_type cardinality = outputValues.extent(0);
165  const ordinal_type spaceDim = 3;
166 
167  auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
168  typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
169 
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) );
176  break;
177  }
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) );
183  break;
184  }
185  default: {
186  INTREPID2_TEST_FOR_EXCEPTION( true , std::invalid_argument,
187  ">>> ERROR (Basis_HCURL_TET_In_FEM): Operator type not implemented" );
188  }
189  }
190 }
191 }
192 
193 // -------------------------------------------------------------------------------------
194 template<typename DT, typename OT, typename PT>
196 Basis_HCURL_TET_In_FEM( const ordinal_type order,
197  const EPointType pointType ) {
198 
199  constexpr ordinal_type spaceDim = 3;
200  this->basisCardinality_ = CardinalityHCurlTet(order);
201  this->basisDegree_ = order; // small n
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_;
208 
209  const ordinal_type cardPn = Intrepid2::getPnCardinality<spaceDim>(order); // dim of (P_{n}) -- smaller space
210  const ordinal_type cardPnm1 = Intrepid2::getPnCardinality<spaceDim>(order-1); // dim of (P_{n-1}) -- smaller space
211  const ordinal_type cardPnm2 = Intrepid2::getPnCardinality<spaceDim>(order-2); // dim of (P_{n-2}) -- smaller space
212  const ordinal_type cardVecPn = spaceDim*cardPn; // dim of (P_{n})^2 -- larger space
213  const ordinal_type cardVecPnm1 = spaceDim*cardPnm1; // dim of (P_{n-1})^2 -- smaller space
214  const ordinal_type cardPnm1H = cardPnm1-cardPnm2; //Homogeneous polynomial of order (n-1)
215 
216  // 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.
217  // 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.)
218  INTREPID2_TEST_FOR_EXCEPTION( order > Parameters::MaxOrder, std::invalid_argument, "polynomial order exceeds the max supported by this class");
219 
220  // Basis-dependent initializations
221  constexpr ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
222  constexpr ordinal_type maxCard = CardinalityHCurlTet(Parameters::MaxOrder);
223  ordinal_type tags[maxCard][tagSize];
224 
225  // points are computed in the host and will be copied
226  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
227  dofCoords("Hcurl::Tet::In::dofCoords", card, spaceDim);
228 
229  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
230  coeffs("Hcurl::Tet::In::coeffs", cardVecPn, card);
231 
232  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
233  dofCoeffs("Hcurl::Tet::In::dofCoeffs", card, spaceDim);
234 
235  // first, need to project the basis for RT space onto the
236  // orthogonal basis of degree n
237  // get coefficients of PkHx
238 
239  Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace> //use LayoutLeft for Lapack
240  V1("Hcurl::Tet::In::V1", cardVecPn, cardVecPnm1 + spaceDim*cardPnm1H);
241 
242 
243  // these two loops get the first three sets of basis functions
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;
247 
248 
249  // now I need to integrate { (x,y) \times phi } against the big basis
250  // first, get a cubature rule.
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() );
254  myCub.getCubature( cubPoints , cubWeights );
255 
256  // tabulate the scalar orthonormal basis at cubature points
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{},
259  phisAtCubPoints,
260  cubPoints,
261  order,
262  OPERATOR_VALUE);
263 
264  // 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)
265  for (ordinal_type i=0;i<cardPn;i++) {
266  for (ordinal_type j=0;j<cardPnm1H;j++) { // loop over homogeneous polynomials
267  for (ordinal_type d=0; d< spaceDim; ++d) {
268  scalarType integral(0);
269  for (ordinal_type k=0;k<myCub.getNumPoints();k++)
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;
276  }
277  }
278  }
279 
280 
281 
282 
283 
284  // now I need to set up an SVD to get a basis for the space
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);
291 
292 
293 
294  ordinal_type info = 0;
295  Teuchos::LAPACK<ordinal_type,scalarType> lapack;
296 
297 
298  lapack.GESVD( 'A',
299  'N',
300  V1.extent(0) ,
301  V1.extent(1) ,
302  V1.data() ,
303  V1.stride_1() ,
304  S.data() ,
305  U.data() ,
306  U.stride_1() ,
307  Vt.data() ,
308  Vt.stride_1() ,
309  work.data() ,
310  5*cardVecPn ,
311  rWork.data() ,
312  &info );
313 
314 
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());
319 
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");
322 #endif
323 
324  // next, apply the RT nodes (rows) to the basis for (P_n)^2 (columns)
325  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
326  V2("Hcurl::Tet::In::V2", card ,cardVecPn);
327 
328  const ordinal_type numEdges = this->basisCellTopology_.getEdgeCount();
329  const ordinal_type numFaces = this->basisCellTopology_.getFaceCount();
330 
331  // first numEdges * degree nodes are normals at each edge
332  // get the points on the line
333 
334  shards::CellTopology edgeTop(shards::getCellTopologyData<shards::Line<2> >() );
335  shards::CellTopology faceTop(shards::getCellTopologyData<shards::Triangle<3> >() );
336 
337  const int numPtsPerEdge = PointTools::getLatticeSize( edgeTop ,
338  order+1 ,
339  1 );
340 
341  const int numPtsPerFace = PointTools::getLatticeSize( faceTop ,
342  order+1 ,
343  1 );
344 
345  const int numPtsPerCell = PointTools::getLatticeSize( this->basisCellTopology_ ,
346  order+1 ,
347  1 );
348 
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 );
351 
352  // construct lattice
353  const ordinal_type offset = 1;
354 
355 
356 
357  PointTools::getLattice( linePts,
358  edgeTop,
359  order+1, offset,
360  pointType );
361 
362  PointTools::getLattice( triPts,
363  faceTop,
364  order+1, offset,
365  pointType );
366 
367  // holds the image of the line points
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);
372 
373  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace> edgeTan("Hcurl::Tet::In::edgeTan", spaceDim );
374 
375  // these are tangents scaled by the appropriate edge lengths.
376  for (ordinal_type i=0;i<numEdges;i++) { // loop over edges
378  i ,
379  this->basisCellTopology_ );
380 
382  linePts ,
383  1 ,
384  i ,
385  this->basisCellTopology_);
386 
387  Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(typename Kokkos::HostSpace::execution_space{},
388  phisAtEdgePoints,
389  edgePts,
390  order,
391  OPERATOR_VALUE);
392 
393  // loop over points (rows of V2)
394  for (ordinal_type j=0;j<numPtsPerEdge;j++) {
395 
396  const ordinal_type i_card = numPtsPerEdge*i+j;
397 
398  // loop over orthonormal basis functions (columns of V2)
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);
402 
403  //save dof coordinates and coefficients
404  for(ordinal_type k=0; k<spaceDim; ++k) {
405  dofCoords(i_card,k) = edgePts(j,k);
406  dofCoeffs(i_card,k) = edgeTan(k);
407  }
408 
409  tags[i_card][0] = 1; // edge dof
410  tags[i_card][1] = i; // edge id
411  tags[i_card][2] = j; // local dof id
412  tags[i_card][3] = numPtsPerEdge; // total vert dof
413 
414  }
415  }
416 
417  if(numPtsPerFace >0) {//handle faces if needed (order >1)
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 );
420 
421  for (ordinal_type i=0;i<numFaces;i++) { // loop over faces
423  faceTan2,
424  i ,
425  this->basisCellTopology_ );
426 
428  triPts ,
429  2 ,
430  i ,
431  this->basisCellTopology_ );
432 
433  Impl::Basis_HGRAD_TET_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(typename Kokkos::HostSpace::execution_space{},
434  phisAtFacePoints,
435  facePts,
436  order,
437  OPERATOR_VALUE);
438 
439  // loop over points (rows of V2)
440  for (ordinal_type j=0;j<numPtsPerFace;j++) {
441 
442  const ordinal_type i_card = numEdges*numPtsPerEdge+2*numPtsPerFace*i+2*j;
443  const ordinal_type i_card_p1 = i_card+1; // creating a temp otherwise nvcc gets confused
444 
445  // loop over orthonormal basis functions (columns of V2)
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);
450  }
451 
452  //save dof coordinates
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);
458  }
459 
460  tags[i_card][0] = 2; // face dof
461  tags[i_card][1] = i; // face id
462  tags[i_card][2] = 2*j; // local face id
463  tags[i_card][3] = 2*numPtsPerFace; // total face dof
464 
465  tags[i_card_p1][0] = 2; // face dof
466  tags[i_card_p1][1] = i; // face id
467  tags[i_card_p1][2] = 2*j+1; // local face id
468  tags[i_card_p1][3] = 2*numPtsPerFace; // total face dof
469 
470  }
471  }
472  }
473 
474 
475  // internal dof, if needed
476  if (numPtsPerCell > 0) {
477  Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
478  cellPoints( "Hcurl::Tet::In::cellPoints", numPtsPerCell , spaceDim );
479  PointTools::getLattice( cellPoints ,
480  this->basisCellTopology_ ,
481  order + 1 ,
482  1 ,
483  pointType );
484 
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{},
488  phisAtCellPoints,
489  cellPoints,
490  order,
491  OPERATOR_VALUE);
492 
493  // copy values into right positions of V2
494  for (ordinal_type j=0;j<numPtsPerCell;j++) {
495 
496  const ordinal_type i_card = numEdges*numPtsPerEdge+2*numFaces*numPtsPerFace+spaceDim*j;
497 
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);
501 
502 
503  //save dof coordinates
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);
508  }
509 
510  tags[i_card+d][0] = spaceDim; // elem dof
511  tags[i_card+d][1] = 0; // elem id
512  tags[i_card+d][2] = spaceDim*j+d; // local dof id
513  tags[i_card+d][3] = spaceDim*numPtsPerCell; // total vert dof
514  }
515  }
516  }
517 
518  // form Vandermonde matrix. Actually, this is the transpose of the VDM,
519  // so we transpose on copy below.
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);
525 
526  //vmat = V2*U;
527  for(ordinal_type i=0; i< card; ++i) {
528  for(ordinal_type j=0; j< card; ++j) {
529  scalarType s=0;
530  for(ordinal_type k=0; k< cardVecPn; ++k)
531  s += V2(i,k)*U(k,j);
532  vmat(i,j) = s;
533  }
534  }
535 
536  info = 0;
537 
538  lapack.GETRF(card, card,
539  vmat.data(), vmat.stride_1(),
540  (ordinal_type*)ipiv.data(),
541  &info);
542 
543  INTREPID2_TEST_FOR_EXCEPTION( info != 0,
544  std::runtime_error ,
545  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM) lapack.GETRF returns nonzero info." );
546 
547  lapack.GETRI(card,
548  vmat.data(), vmat.stride_1(),
549  (ordinal_type*)ipiv.data(),
550  work1.data(), lwork,
551  &info);
552 
553  INTREPID2_TEST_FOR_EXCEPTION( info != 0,
554  std::runtime_error ,
555  ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM) lapack.GETRI returns nonzero info." );
556 
557  for (ordinal_type i=0;i<cardVecPn;++i) {
558  for (ordinal_type j=0;j<card;++j){
559  scalarType s=0;
560  for(ordinal_type k=0; k< card; ++k)
561  s += U(i,k)*vmat(k,j);
562  coeffs(i,j) = s;
563  }
564  }
565 
566  this->coeffs_ = Kokkos::create_mirror_view(typename DT::memory_space(), coeffs);
567  Kokkos::deep_copy(this->coeffs_ , coeffs);
568 
569  this->dofCoords_ = Kokkos::create_mirror_view(typename DT::memory_space(), dofCoords);
570  Kokkos::deep_copy(this->dofCoords_, dofCoords);
571 
572  this->dofCoeffs_ = Kokkos::create_mirror_view(typename DT::memory_space(), dofCoeffs);
573  Kokkos::deep_copy(this->dofCoeffs_, dofCoeffs);
574 
575 
576  // set tags
577  {
578  // Basis-dependent initializations
579  const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
580  const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
581  const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
582 
583  OrdinalTypeArray1DHost tagView(&tags[0][0], card*tagSize);
584 
585  // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
586  // tags are constructed on host
587  this->setOrdinalTagData(this->tagToOrdinal_,
588  this->ordinalToTag_,
589  tagView,
590  this->basisCardinality_,
591  tagSize,
592  posScDim,
593  posScOrd,
594  posDfOrd);
595  }
596 }
597 } // namespace Intrepid2
598 #endif
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.
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...
virtual ordinal_type getNumPoints() const override
Returns the number of cubature points.
static void getReferenceEdgeTangent(RefEdgeTangentViewType refEdgeTangent, const ordinal_type edgeOrd, const shards::CellTopology parentCell)
Computes constant tangent vectors to edges of 2D or 3D reference cells.
Header file for the Intrepid2::CubatureDirectTetDefault class.
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.
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.
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 getReferenceFaceTangents(RefFaceTanViewType refFaceTanU, RefFaceTanViewType refFaceTanV, const ordinal_type faceOrd, const shards::CellTopology parentCell)
Computes pairs of constant tangent vectors to faces of a 3D reference cells.