15 #ifndef __INTREPID2_POINTTOOLS_DEF_HPP__
16 #define __INTREPID2_POINTTOOLS_DEF_HPP__
18 #if defined(_MSC_VER) || defined(_WIN32) && defined(__ICL)
20 #ifndef _USE_MATH_DEFINES
21 #define _USE_MATH_DEFINES
36 const ordinal_type order,
37 const ordinal_type offset ) {
38 #ifdef HAVE_INTREPID2_DEBUG
39 INTREPID2_TEST_FOR_EXCEPTION( order < 0 || offset < 0,
40 std::invalid_argument ,
41 ">>> ERROR (PointTools::getLatticeSize): order and offset must be positive values." );
43 ordinal_type r_val = 0;
44 switch (cellType.getBaseKey()) {
45 case shards::Tetrahedron<>::key: {
46 const auto effectiveOrder = order - 4 * offset;
47 r_val = (effectiveOrder < 0 ? 0 :(effectiveOrder+1)*(effectiveOrder+2)*(effectiveOrder+3)/6);
50 case shards::Triangle<>::key: {
51 const auto effectiveOrder = order - 3 * offset;
52 r_val = (effectiveOrder < 0 ? 0 : (effectiveOrder+1)*(effectiveOrder+2)/2);
55 case shards::Line<>::key: {
56 const auto effectiveOrder = order - 2 * offset;
57 r_val = (effectiveOrder < 0 ? 0 : (effectiveOrder+1));
60 case shards::Quadrilateral<>::key: {
61 const auto effectiveOrder = order - 2 * offset;
62 r_val = std::pow(effectiveOrder < 0 ? 0 : (effectiveOrder+1),2);
65 case shards::Hexahedron<>::key: {
66 const auto effectiveOrder = order - 2 * offset;
67 r_val = std::pow(effectiveOrder < 0 ? 0 : (effectiveOrder+1),3);
70 case shards::Pyramid<>::key: {
71 const auto effectiveOrder = order - 2 * offset;
72 r_val = (effectiveOrder < 0 ? 0 : (effectiveOrder+1)*(effectiveOrder+2)*(2*effectiveOrder+3)/6);
76 INTREPID2_TEST_FOR_EXCEPTION(
true , std::invalid_argument ,
77 ">>> ERROR (Intrepid2::PointTools::getLatticeSize): the specified cell type is not supported." );
83 template<
typename pointValueType,
class ...pointProperties>
86 getLattice( Kokkos::DynRankView<pointValueType,pointProperties...> points,
87 const shards::CellTopology cell,
88 const ordinal_type order,
89 const ordinal_type offset,
90 const EPointType pointType ) {
91 #ifdef HAVE_INTREPID2_DEBUG
92 INTREPID2_TEST_FOR_EXCEPTION( points.rank() != 2,
93 std::invalid_argument ,
94 ">>> ERROR (PointTools::getLattice): points rank must be 2." );
95 INTREPID2_TEST_FOR_EXCEPTION( order < 0 || offset < 0,
96 std::invalid_argument ,
97 ">>> ERROR (PointTools::getLattice): order and offset must be positive values." );
99 const size_type latticeSize =
getLatticeSize( cell, order, offset );
100 const size_type spaceDim = cell.getDimension();
102 INTREPID2_TEST_FOR_EXCEPTION( points.extent(0) != latticeSize ||
103 points.extent(1) != spaceDim,
104 std::invalid_argument ,
105 ">>> ERROR (PointTools::getLattice): dimension does not match to lattice size." );
108 switch (cell.getBaseKey()) {
110 case shards::Pyramid<>::key:
getLatticePyramid ( points, order, offset, pointType );
break;
111 case shards::Triangle<>::key:
getLatticeTriangle ( points, order, offset, pointType );
break;
112 case shards::Line<>::key:
getLatticeLine ( points, order, offset, pointType );
break;
113 case shards::Quadrilateral<>::key: {
114 auto hostPoints = Kokkos::create_mirror_view(points);
115 shards::CellTopology line(shards::getCellTopologyData<shards::Line<2> >());
116 const ordinal_type numPoints =
getLatticeSize( line, order, offset );
117 auto linePoints = getMatchingViewWithLabel(hostPoints,
"linePoints", numPoints, 1);
120 for (ordinal_type j=0; j<numPoints; ++j) {
121 for (ordinal_type i=0; i<numPoints; ++i, ++idx) {
122 hostPoints(idx,0) = linePoints(i,0);
123 hostPoints(idx,1) = linePoints(j,0);
126 Kokkos::deep_copy(points,hostPoints);
129 case shards::Hexahedron<>::key: {
130 auto hostPoints = Kokkos::create_mirror_view(points);
131 shards::CellTopology line(shards::getCellTopologyData<shards::Line<2> >());
132 const ordinal_type numPoints =
getLatticeSize( line, order, offset );
133 auto linePoints = getMatchingViewWithLabel(hostPoints,
"linePoints", numPoints, 1);
136 for (ordinal_type k=0; k<numPoints; ++k) {
137 for (ordinal_type j=0; j<numPoints; ++j) {
138 for (ordinal_type i=0; i<numPoints; ++i, ++idx) {
139 hostPoints(idx,0) = linePoints(i,0);
140 hostPoints(idx,1) = linePoints(j,0);
141 hostPoints(idx,2) = linePoints(k,0);
145 Kokkos::deep_copy(points,hostPoints);
149 INTREPID2_TEST_FOR_EXCEPTION(
true , std::invalid_argument ,
150 ">>> ERROR (Intrepid2::PointTools::getLattice): the specified cell type is not supported." );
155 template<
typename pointValueType,
class ...pointProperties>
159 const ordinal_type order ) {
160 #ifdef HAVE_INTREPID2_DEBUG
161 INTREPID2_TEST_FOR_EXCEPTION( points.rank() != 2,
162 std::invalid_argument ,
163 ">>> ERROR (PointTools::getGaussPoints): points rank must be 1." );
164 INTREPID2_TEST_FOR_EXCEPTION( order < 0,
165 std::invalid_argument ,
166 ">>> ERROR (PointTools::getGaussPoints): order must be positive value." );
168 const ordinal_type np = order + 1;
169 const double alpha = 0.0, beta = 0.0;
172 Kokkos::View<pointValueType*,Kokkos::HostSpace>
173 zHost(
"PointTools::getGaussPoints::z", np),
174 wHost(
"PointTools::getGaussPoints::w", np);
181 typedef Kokkos::pair<ordinal_type,ordinal_type> range_type;
182 auto pts = Kokkos::subview( points, range_type(0,np), 0 );
184 auto z = Kokkos::DynRankView<pointValueType,Kokkos::HostSpace>(zHost.data(), np);
185 Kokkos::deep_copy(pts, z);
192 template<
typename pointValueType,
class ...pointProperties>
196 const ordinal_type order,
197 const ordinal_type offset,
198 const EPointType pointType ) {
203 INTREPID2_TEST_FOR_EXCEPTION(
true ,
204 std::invalid_argument ,
205 ">>> ERROR (PointTools::getLattice): invalid EPointType." );
210 template<
typename pointValueType,
class ...pointProperties>
214 const ordinal_type order,
215 const ordinal_type offset,
216 const EPointType pointType ) {
221 INTREPID2_TEST_FOR_EXCEPTION(
true ,
222 std::invalid_argument ,
223 ">>> ERROR (PointTools::getLattice): invalid EPointType." );
228 template<
typename pointValueType,
class ...pointProperties>
232 const ordinal_type order,
233 const ordinal_type offset,
234 const EPointType pointType ) {
239 INTREPID2_TEST_FOR_EXCEPTION(
true ,
240 std::invalid_argument ,
241 ">>> ERROR (PointTools::getLattice): invalid EPointType." );
246 template<
typename pointValueType,
class ...pointProperties>
250 const ordinal_type order,
251 const ordinal_type offset,
252 const EPointType pointType )
257 INTREPID2_TEST_FOR_EXCEPTION(
true ,
258 std::invalid_argument ,
259 ">>> ERROR (PointTools::getLattice): invalid EPointType." );
266 template<
typename pointValueType,
class ...pointProperties>
270 const ordinal_type order,
271 const ordinal_type offset ) {
272 auto pointsHost = Kokkos::create_mirror_view(points);
275 pointsHost(0,0) = 0.0;
277 const pointValueType h = 2.0 / order;
278 const ordinal_type ibeg = offset, iend = order-offset+1;
279 for (ordinal_type i=ibeg;i<iend;++i)
280 pointsHost(i-ibeg, 0) = -1.0 + h * (pointValueType) i;
283 Kokkos::deep_copy(points, pointsHost);
286 template<
typename pointValueType,
class ...pointProperties>
290 const ordinal_type order,
291 const ordinal_type offset ) {
294 const ordinal_type np = order + 1;
295 const double alpha = 0.0, beta = 0.0;
296 const ordinal_type s = np - 2*offset;
300 Kokkos::View<pointValueType*,Kokkos::HostSpace>
301 zHost(
"PointTools::getGaussPoints::z", np),
302 wHost(
"PointTools::getGaussPoints::w", np);
309 typedef Kokkos::pair<ordinal_type,ordinal_type> range_type;
310 auto pts = Kokkos::subview( points, range_type(0, s), 0 );
313 auto z = Kokkos::DynRankView<pointValueType,Kokkos::HostSpace>(zHost.data() + offset, np-offset);
315 const auto common_range = range_type(0, std::min(pts.extent(0), z.extent(0)));
316 Kokkos::deep_copy(Kokkos::subview(pts, common_range),
317 Kokkos::subview(z, common_range));
321 template<
typename pointValueType,
class ...pointProperties>
325 const ordinal_type order,
326 const ordinal_type offset ) {
327 TEUCHOS_TEST_FOR_EXCEPTION( order <= 0 ,
328 std::invalid_argument ,
329 ">>> ERROR (Intrepid2::PointTools::getEquispacedLatticeTriangle): order must be positive" );
331 auto pointsHost = Kokkos::create_mirror_view(points);
333 const pointValueType h = 1.0 / order;
334 ordinal_type cur = 0;
336 for (ordinal_type i=offset;i<=order-offset;i++) {
337 for (ordinal_type j=offset;j<=order-i-offset;j++) {
338 pointsHost(cur,0) = j * h ;
339 pointsHost(cur,1) = i * h;
343 Kokkos::deep_copy(points, pointsHost);
346 template<
typename pointValueType,
class ...pointProperties>
350 const ordinal_type order ,
351 const ordinal_type offset ) {
352 TEUCHOS_TEST_FOR_EXCEPTION( (order <= 0) ,
353 std::invalid_argument ,
354 ">>> ERROR (Intrepid2::PointTools::getEquispacedLatticeTetrahedron): order must be positive" );
356 auto pointsHost = Kokkos::create_mirror_view(points);
358 const pointValueType h = 1.0 / order;
359 ordinal_type cur = 0;
361 for (ordinal_type i=offset;i<=order-offset;i++) {
362 for (ordinal_type j=offset;j<=order-i-offset;j++) {
363 for (ordinal_type k=offset;k<=order-i-j-offset;k++) {
364 pointsHost(cur,0) = k * h;
365 pointsHost(cur,1) = j * h;
366 pointsHost(cur,2) = i * h;
371 Kokkos::deep_copy(points, pointsHost);
374 template<
typename pointValueType,
class ...pointProperties>
378 const ordinal_type order ,
379 const ordinal_type offset ) {
380 TEUCHOS_TEST_FOR_EXCEPTION( (order <= 0) ,
381 std::invalid_argument ,
382 ">>> ERROR (Intrepid2::PointTools::getEquispacedLatticePyramid): order must be positive" );
384 auto pointsHost = Kokkos::create_mirror_view(points);
386 const pointValueType h = 1.0 / order;
387 ordinal_type cur = 0;
389 for (ordinal_type i=offset;i<=order-offset;i++) {
390 pointValueType z = i * h;
391 for (ordinal_type j=offset;j<=order-i-offset;j++) {
392 for (ordinal_type k=offset;k<=order-i-offset;k++) {
393 pointsHost(cur,0) = (2 * k * h - 1.) * (1-z);
394 pointsHost(cur,1) = (2 * j * h - 1.) * (1-z);
395 pointsHost(cur,2) = z;
400 Kokkos::deep_copy(points, pointsHost);
403 template<
typename pointValueType,
class ...pointProperties>
406 warpFactor( Kokkos::DynRankView<pointValueType,pointProperties...> warp,
407 const ordinal_type order,
408 const Kokkos::DynRankView<pointValueType,pointProperties...> xnodes ,
409 const Kokkos::DynRankView<pointValueType,pointProperties...> xout
412 TEUCHOS_TEST_FOR_EXCEPTION( ( warp.extent(0) != xout.extent(0) ) ,
413 std::invalid_argument ,
414 ">>> ERROR (PointTools::warpFactor): xout and warp must be same size." );
416 Kokkos::deep_copy(warp, pointValueType(0.0));
418 ordinal_type xout_dim0 = xout.extent(0);
419 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> d(
"d", xout_dim0 );
421 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> xeq_(
"xeq", order + 1 ,1);
423 const auto xeq = Kokkos::subview(xeq_, Kokkos::ALL(),0);
425 TEUCHOS_TEST_FOR_EXCEPTION( ( xeq.extent(0) != xnodes.extent(0) ) ,
426 std::invalid_argument ,
427 ">>> ERROR (PointTools::warpFactor): xeq and xnodes must be same size." );
430 for (ordinal_type i=0;i<=order;i++) {
432 Kokkos::deep_copy(d, xnodes(i,0) - xeq(i));
434 for (ordinal_type j=1;j<order;j++) {
436 for (ordinal_type k=0;k<xout_dim0;k++) {
437 d(k) = d(k) * ( (xout(k)-xeq(j)) / (xeq(i)-xeq(j)) );
444 for (ordinal_type k=0;k<xout_dim0;k++) {
445 d(k) = -d(k) / (xeq(i) - xeq(0));
450 for (ordinal_type k=0;k<xout_dim0;k++) {
451 d(k) = d(k) / (xeq(i) - xeq(order));
455 for (ordinal_type k=0;k<xout_dim0;k++) {
462 template<
typename pointValueType,
class ...pointProperties>
466 const ordinal_type order ,
467 const ordinal_type offset)
471 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> gaussX(
"gaussX", order + 1 , 1 );
478 pointValueType alpopt[] = {0.0000,0.0000,1.4152,0.1001,0.2751,0.9800,1.0999,
479 1.2832,1.3648, 1.4773, 1.4959, 1.5743, 1.5770, 1.6223, 1.6258};
481 pointValueType alpha;
483 if (order >= 1 && order < 16) {
484 alpha = alpopt[order-1];
490 const ordinal_type p = order;
491 ordinal_type N = (p+1)*(p+2)/2;
494 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> L1(
"L1", N );
495 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> L2(
"L2", N );
496 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> L3(
"L3", N );
497 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> X(
"X", N);
498 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> Y(
"Y", N);
501 for (ordinal_type n=1;n<=p+1;n++) {
502 for (ordinal_type m=1;m<=p+2-n;m++) {
503 L1(sk) = (n-1) / (pointValueType)p;
504 L3(sk) = (m-1) / (pointValueType)p;
505 L2(sk) = 1.0 - L1(sk) - L3(sk);
510 for (ordinal_type n=0;n<N;n++) {
511 X(n) = -L2(n) + L3(n);
512 Y(n) = (-L2(n) - L3(n) + 2*L1(n))/1.7320508075688772;
516 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> blend1(
"blend1", N);
517 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> blend2(
"blend2", N);
518 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> blend3(
"blend3", N);
520 for (ordinal_type n=0;n<N;n++) {
521 blend1(n) = 4.0 * L2(n) * L3(n);
522 blend2(n) = 4.0 * L1(n) * L3(n);
523 blend3(n) = 4.0 * L1(n) * L2(n);
527 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> L3mL2(
"L3mL2", N);
528 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> L1mL3(
"L1mL3", N);
529 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> L2mL1(
"L2mL1", N);
531 for (ordinal_type k=0;k<N;k++) {
532 L3mL2(k) = L3(k)-L2(k);
533 L1mL3(k) = L1(k)-L3(k);
534 L2mL1(k) = L2(k)-L1(k);
537 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> warpfactor1(
"warpfactor1", N);
538 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> warpfactor2(
"warpfactor2", N);
539 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> warpfactor3(
"warpfactor3", N);
541 warpFactor( warpfactor1, order , gaussX , L3mL2 );
542 warpFactor( warpfactor2, order , gaussX , L1mL3 );
543 warpFactor( warpfactor3, order , gaussX , L2mL1 );
545 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> warp1(
"warp1", N);
546 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> warp2(
"warp2", N);
547 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> warp3(
"warp3", N);
549 for (ordinal_type k=0;k<N;k++) {
550 warp1(k) = blend1(k) * warpfactor1(k) *
551 ( 1.0 + alpha * alpha * L1(k) * L1(k) );
552 warp2(k) = blend2(k) * warpfactor2(k) *
553 ( 1.0 + alpha * alpha * L2(k) * L2(k) );
554 warp3(k) = blend3(k) * warpfactor3(k) *
555 ( 1.0 + alpha * alpha * L3(k) * L3(k) );
558 for (ordinal_type k=0;k<N;k++) {
559 X(k) += 1.0 * warp1(k) + cos( 2.0 * M_PI / 3.0 ) * warp2(k) + cos(4*M_PI/3.0) * warp3(k);
560 Y(k) += 0.0 * warp1(k) + sin( 2.0 * M_PI / 3.0 ) * warp2(k) + sin( 4*M_PI/3.0) * warp3(k);
563 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> warXY(
"warXY", 1, N,2);
565 for (ordinal_type k=0;k<N;k++) {
566 warXY(0, k,0) = X(k);
567 warXY(0, k,1) = Y(k);
572 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> warburtonVerts(
"warburtonVerts", 1,3,2);
573 warburtonVerts(0,0,0) = -1.0;
574 warburtonVerts(0,0,1) = -1.0/std::sqrt(3.0);
575 warburtonVerts(0,1,0) = 1.0;
576 warburtonVerts(0,1,1) = -1.0/std::sqrt(3.0);
577 warburtonVerts(0,2,0) = 0.0;
578 warburtonVerts(0,2,1) = 2.0/std::sqrt(3.0);
580 Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> refPts(
"refPts", 1, N,2);
585 shards::getCellTopologyData< shards::Triangle<3> >()
589 auto pointsHost = Kokkos::create_mirror_view(points);
591 ordinal_type noffcur = 0;
592 ordinal_type offcur = 0;
593 for (ordinal_type i=0;i<=order;i++) {
594 for (ordinal_type j=0;j<=order-i;j++) {
595 if ( (i >= offset) && (i <= order-offset) &&
596 (j >= offset) && (j <= order-i-offset) ) {
597 pointsHost(offcur,0) = refPts(0, noffcur,0);
598 pointsHost(offcur,1) = refPts(0, noffcur,1);
605 Kokkos::deep_copy(points, pointsHost);
610 template<
typename pointValueType,
class ...pointProperties>
614 const ordinal_type order ,
615 const pointValueType pval ,
616 const Kokkos::DynRankView<pointValueType,pointProperties...> ,
617 const Kokkos::DynRankView<pointValueType,pointProperties...> L2,
618 const Kokkos::DynRankView<pointValueType,pointProperties...> L3,
619 const Kokkos::DynRankView<pointValueType,pointProperties...> L4
626 template<
typename pointValueType,
class ...pointProperties>
629 evalshift( Kokkos::DynRankView<pointValueType,pointProperties...> dxy,
630 const ordinal_type order ,
631 const pointValueType pval ,
632 const Kokkos::DynRankView<pointValueType,pointProperties...> L1 ,
633 const Kokkos::DynRankView<pointValueType,pointProperties...> L2 ,
634 const Kokkos::DynRankView<pointValueType,pointProperties...> L3
638 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> gaussX(
"gaussX",order+1,1);
641 const ordinal_type N = L1.extent(0);
644 for (ordinal_type k=0;k<=order;k++) {
649 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> blend1(
"blend1",N);
650 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> blend2(
"blend2",N);
651 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> blend3(
"blend3",N);
653 for (ordinal_type i=0;i<N;i++) {
654 blend1(i) = L2(i) * L3(i);
655 blend2(i) = L1(i) * L3(i);
656 blend3(i) = L1(i) * L2(i);
660 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> warpfactor1(
"warpfactor1",N);
661 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> warpfactor2(
"warpfactor2",N);
662 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> warpfactor3(
"warpfactor3",N);
665 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> L3mL2(
"L3mL2",N);
666 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> L1mL3(
"L1mL3",N);
667 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> L2mL1(
"L2mL1",N);
669 for (ordinal_type k=0;k<N;k++) {
670 L3mL2(k) = L3(k)-L2(k);
671 L1mL3(k) = L1(k)-L3(k);
672 L2mL1(k) = L2(k)-L1(k);
675 evalwarp( warpfactor1 , order , gaussX , L3mL2 );
676 evalwarp( warpfactor2 , order , gaussX , L1mL3 );
677 evalwarp( warpfactor3 , order , gaussX , L2mL1 );
679 for (ordinal_type k=0;k<N;k++) {
680 warpfactor1(k) *= 4.0;
681 warpfactor2(k) *= 4.0;
682 warpfactor3(k) *= 4.0;
685 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> warp1(
"warp1",N);
686 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> warp2(
"warp2",N);
687 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> warp3(
"warp3",N);
689 for (ordinal_type k=0;k<N;k++) {
690 warp1(k) = blend1(k) * warpfactor1(k) *
691 ( 1.0 + pval * pval * L1(k) * L1(k) );
692 warp2(k) = blend2(k) * warpfactor2(k) *
693 ( 1.0 + pval * pval * L2(k) * L2(k) );
694 warp3(k) = blend3(k) * warpfactor3(k) *
695 ( 1.0 + pval * pval * L3(k) * L3(k) );
698 for (ordinal_type k=0;k<N;k++) {
699 dxy(k,0) = 1.0 * warp1(k) + cos( 2.0 * M_PI / 3.0 ) * warp2(k) + cos( 4.0*M_PI/3.0 ) * warp3(k);
700 dxy(k,1) = 0.0 * warp1(k) + sin( 2.0 * M_PI / 3.0 ) * warp2(k) + sin( 4.0*M_PI/3.0 ) * warp3(k);
705 template<
typename pointValueType,
class ...pointProperties>
708 evalwarp(Kokkos::DynRankView<pointValueType,pointProperties...> warp ,
709 const ordinal_type order ,
710 const Kokkos::DynRankView<pointValueType,pointProperties...> xnodes,
711 const Kokkos::DynRankView<pointValueType,pointProperties...> xout )
713 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> xeq(
"xeq",order+1);
715 ordinal_type xout_dim0 = xout.extent(0);
716 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> d(
"d",xout_dim0);
720 for (ordinal_type i=0;i<=order;i++) {
721 xeq(i) = -1.0 + 2.0 * ( order - i ) / order;
724 for (ordinal_type i=0;i<=order;i++) {
725 Kokkos::deep_copy(d, xnodes(i,0) - xeq(i));
726 for (ordinal_type j=1;j<order;j++) {
728 for (ordinal_type k=0;k<xout_dim0;k++) {
729 d(k) = d(k) * (xout(k)-xeq(j))/(xeq(i)-xeq(j));
734 for (ordinal_type k=0;k<xout_dim0;k++) {
735 d(k) = -d(k)/(xeq(i)-xeq(0));
739 for (ordinal_type k=0;k<xout_dim0;k++) {
740 d(k) = d(k)/(xeq(i)-xeq(order));
744 for (ordinal_type k=0;k<xout_dim0;k++) {
751 template<
typename pointValueType,
class ...pointProperties>
755 const ordinal_type order ,
756 const ordinal_type offset )
758 pointValueType alphastore[] = { 0,0,0,0.1002, 1.1332,1.5608,1.3413,1.2577,1.1603,
759 1.10153,0.6080,0.4523,0.8856,0.8717,0.9655};
760 pointValueType alpha;
763 alpha = alphastore[std::max(order-1,0)];
769 const ordinal_type N = (order+1)*(order+2)*(order+3)/6;
770 pointValueType tol = 1.e-10;
772 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> shift(
"shift",N,3);
773 Kokkos::deep_copy(shift,0.0);
776 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> equipoints(
"equipoints", N,3);
778 for (ordinal_type n=0;n<=order;n++) {
779 for (ordinal_type m=0;m<=order-n;m++) {
780 for (ordinal_type q=0;q<=order-n-m;q++) {
781 equipoints(sk,0) = -1.0 + (q * 2.0 ) / order;
782 equipoints(sk,1) = -1.0 + (m * 2.0 ) / order;
783 equipoints(sk,2) = -1.0 + (n * 2.0 ) / order;
791 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> L1(
"L1",N);
792 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> L2(
"L2",N);
793 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> L3(
"L3",N);
794 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> L4(
"L4",N);
795 for (ordinal_type i=0;i<N;i++) {
796 L1(i) = (1.0 + equipoints(i,2)) / 2.0;
797 L2(i) = (1.0 + equipoints(i,1)) / 2.0;
798 L3(i) = -(1.0 + equipoints(i,0) + equipoints(i,1) + equipoints(i,2)) / 2.0;
799 L4(i) = (1.0 + equipoints(i,0)) / 2.0;
804 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> warVerts_(
"warVerts",1,4,3);
805 auto warVerts = Kokkos::subview(warVerts_,0,Kokkos::ALL(),Kokkos::ALL());
806 warVerts(0,0) = -1.0;
807 warVerts(0,1) = -1.0/sqrt(3.0);
808 warVerts(0,2) = -1.0/sqrt(6.0);
810 warVerts(1,1) = -1.0/sqrt(3.0);
811 warVerts(1,2) = -1.0/sqrt(6.0);
813 warVerts(2,1) = 2.0 / sqrt(3.0);
814 warVerts(2,2) = -1.0/sqrt(6.0);
817 warVerts(3,2) = 3.0 / sqrt(6.0);
821 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> t1(
"t1",4,3);
822 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> t2(
"t2",4,3);
823 for (ordinal_type i=0;i<3;i++) {
824 t1(0,i) = warVerts(1,i) - warVerts(0,i);
825 t1(1,i) = warVerts(1,i) - warVerts(0,i);
826 t1(2,i) = warVerts(2,i) - warVerts(1,i);
827 t1(3,i) = warVerts(2,i) - warVerts(0,i);
828 t2(0,i) = warVerts(2,i) - 0.5 * ( warVerts(0,i) + warVerts(1,i) );
829 t2(1,i) = warVerts(3,i) - 0.5 * ( warVerts(0,i) + warVerts(1,i) );
830 t2(2,i) = warVerts(3,i) - 0.5 * ( warVerts(1,i) + warVerts(2,i) );
831 t2(3,i) = warVerts(3,i) - 0.5 * ( warVerts(0,i) + warVerts(2,i) );
835 for (ordinal_type n=0;n<4;n++) {
837 pointValueType normt1n = 0.0;
838 pointValueType normt2n = 0.0;
839 for (ordinal_type i=0;i<3;i++) {
840 normt1n += (t1(n,i) * t1(n,i));
841 normt2n += (t2(n,i) * t2(n,i));
843 normt1n = sqrt(normt1n);
844 normt2n = sqrt(normt2n);
846 for (ordinal_type i=0;i<3;i++) {
853 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> XYZ(
"XYZ",N,3);
854 for (ordinal_type i=0;i<N;i++) {
855 for (ordinal_type j=0;j<3;j++) {
856 XYZ(i,j) = L3(i)*warVerts(0,j) + L4(i)*warVerts(1,j) + L2(i)*warVerts(2,j) + L1(i)*warVerts(3,j);
860 for (ordinal_type face=1;face<=4;face++) {
861 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> La, Lb, Lc, Ld;
862 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> warp(
"warp",N,2);
863 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> blend(
"blend",N);
864 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> denom(
"denom",N);
867 La = L1; Lb = L2; Lc = L3; Ld = L4;
break;
869 La = L2; Lb = L1; Lc = L3; Ld = L4;
break;
871 La = L3; Lb = L1; Lc = L4; Ld = L2;
break;
873 La = L4; Lb = L1; Lc = L3; Ld = L2;
break;
879 for (ordinal_type k=0;k<N;k++) {
880 blend(k) = Lb(k) * Lc(k) * Ld(k);
883 for (ordinal_type k=0;k<N;k++) {
884 denom(k) = (Lb(k) + 0.5 * La(k)) * (Lc(k) + 0.5*La(k)) * (Ld(k) + 0.5 * La(k));
887 for (ordinal_type k=0;k<N;k++) {
888 if (denom(k) > tol) {
889 blend(k) *= ( 1.0 + alpha * alpha * La(k) * La(k) ) / denom(k);
895 for (ordinal_type k=0;k<N;k++) {
896 for (ordinal_type j=0;j<3;j++) {
897 shift(k,j) = shift(k,j) + blend(k) * warp(k,0) * t1(face-1,j)
898 + blend(k) * warp(k,1) * t2(face-1,j);
902 for (ordinal_type k=0;k<N;k++) {
903 if (La(k) < tol && ( Lb(k) < tol || Lc(k) < tol || Ld(k) < tol )) {
904 for (ordinal_type j=0;j<3;j++) {
905 shift(k,j) = warp(k,0) * t1(face-1,j) + warp(k,1) * t2(face-1,j);
912 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> updatedPoints(
"updatedPoints",1,N,3);
913 for (ordinal_type k=0;k<N;k++) {
914 for (ordinal_type j=0;j<3;j++) {
915 updatedPoints(0,k,j) = XYZ(k,j) + shift(k,j);
922 Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> refPts(
"refPts",1,N,3);
925 shards::getCellTopologyData<shards::Tetrahedron<4> >()
928 auto pointsHost = Kokkos::create_mirror_view(points);
930 ordinal_type noffcur = 0;
931 ordinal_type offcur = 0;
932 for (ordinal_type i=0;i<=order;i++) {
933 for (ordinal_type j=0;j<=order-i;j++) {
934 for (ordinal_type k=0;k<=order-i-j;k++) {
935 if ( (i >= offset) && (i <= order-offset) &&
936 (j >= offset) && (j <= order-i-offset) &&
937 (k >= offset) && (k <= order-i-j-offset) ) {
938 pointsHost(offcur,0) = refPts(0,noffcur,0);
939 pointsHost(offcur,1) = refPts(0,noffcur,1);
940 pointsHost(offcur,2) = refPts(0,noffcur,2);
948 Kokkos::deep_copy(points, pointsHost);
Gauss-Jacobi/Gauss-Radau-Jacobi/Gauss-Lobatto zeros and weights.