55 template<
class Scalar,
class ArrayType>
57 const shards::CellTopology& cellType ,
60 const EPointType pointType )
63 case POINTTYPE_EQUISPACED:
64 getEquispacedLattice<Scalar,ArrayType>( cellType , pts , order , offset );
66 case POINTTYPE_WARPBLEND:
68 getWarpBlendLattice<Scalar,ArrayType>( cellType , pts , order , offset );
71 TEUCHOS_TEST_FOR_EXCEPTION(
true ,
72 std::invalid_argument ,
73 "PointTools::getLattice: invalid EPointType" );
78 template<
class Scalar,
class ArrayType>
83 Scalar *z =
new Scalar[order+1];
84 Scalar *w =
new Scalar[order+1];
87 for (
int i=0;i<order+1;i++) {
96 template<
class Scalar,
class ArrayType>
103 switch (cellType.getKey()) {
104 case shards::Tetrahedron<4>::key:
105 case shards::Tetrahedron<8>::key:
106 case shards::Tetrahedron<10>::key:
107 TEUCHOS_TEST_FOR_EXCEPTION( ( points.dimension(0) !=
getLatticeSize( cellType , order , offset ) )
108 || points.dimension(1) != 3 ,
109 std::invalid_argument ,
110 ">>> ERROR(PointTools::getEquispacedLattice): points argument is ill-sized." );
111 getEquispacedLatticeTetrahedron<Scalar,ArrayType>( points , order , offset );
113 case shards::Triangle<3>::key:
114 case shards::Triangle<4>::key:
115 case shards::Triangle<6>::key:
116 TEUCHOS_TEST_FOR_EXCEPTION( ( points.dimension(0) !=
getLatticeSize( cellType , order , offset ) )
117 || points.dimension(1) != 2 ,
118 std::invalid_argument ,
119 ">>> ERROR(PointTools::getEquispacedLattice): points argument is ill-sized." );
120 getEquispacedLatticeTriangle<Scalar,ArrayType>( points , order , offset );
122 case shards::Line<2>::key:
123 case shards::Line<3>::key:
124 TEUCHOS_TEST_FOR_EXCEPTION( ( points.dimension(0) !=
getLatticeSize( cellType , order , offset ) )
125 || points.dimension(1) != 1 ,
126 std::invalid_argument ,
127 ">>> ERROR(PointTools::getEquispacedLattice): points argument is ill-sized." );
128 getEquispacedLatticeLine<Scalar,ArrayType>( points , order , offset );
131 TEUCHOS_TEST_FOR_EXCEPTION(
true , std::invalid_argument ,
132 ">>> ERROR (Intrepid::PointTools::getEquispacedLattice): Illegal cell type" );
137 template<
class Scalar,
class ArrayType>
145 switch (cellType.getKey()) {
146 case shards::Tetrahedron<4>::key:
147 case shards::Tetrahedron<8>::key:
148 case shards::Tetrahedron<10>::key:
149 TEUCHOS_TEST_FOR_EXCEPTION( ( points.dimension(0) !=
getLatticeSize( cellType , order , offset ) )
150 || points.dimension(1) != 3 ,
151 std::invalid_argument ,
152 ">>> ERROR(PointTools::getWarpBlendLattice): points argument is ill-sized." );
154 getWarpBlendLatticeTetrahedron<Scalar,ArrayType>( points , order , offset );
156 case shards::Triangle<3>::key:
157 case shards::Triangle<4>::key:
158 case shards::Triangle<6>::key:
159 TEUCHOS_TEST_FOR_EXCEPTION( ( points.dimension(0) !=
getLatticeSize( cellType , order , offset ) )
160 || points.dimension(1) != 2 ,
161 std::invalid_argument ,
162 ">>> ERROR(PointTools::getWarpBlendLattice): points argument is ill-sized." );
164 getWarpBlendLatticeTriangle<Scalar,ArrayType>( points , order , offset );
166 case shards::Line<2>::key:
167 case shards::Line<3>::key:
168 TEUCHOS_TEST_FOR_EXCEPTION( ( points.dimension(0) !=
getLatticeSize( cellType , order , offset ) )
169 || points.dimension(1) != 1 ,
170 std::invalid_argument ,
171 ">>> ERROR(PointTools::getWarpBlendLattice): points argument is ill-sized." );
173 getWarpBlendLatticeLine<Scalar,ArrayType>( points , order , offset );
176 TEUCHOS_TEST_FOR_EXCEPTION(
true , std::invalid_argument ,
177 ">>> ERROR (Intrepid::PointTools::getWarpBlendLattice): Illegal cell type" );
182 template<
class Scalar,
class ArrayType>
187 TEUCHOS_TEST_FOR_EXCEPTION( order < 0 ,
188 std::invalid_argument ,
189 ">>> ERROR (Intrepid::PointTools::getEquispacedLatticeLine): order must be positive" );
194 const Scalar h = 2.0 / order;
196 for (
int i=offset;i<=order-offset;i++) {
197 points(i-offset,0) = -1.0 + h * (Scalar) i;
204 template<
class Scalar,
class ArrayType>
209 TEUCHOS_TEST_FOR_EXCEPTION( order <= 0 ,
210 std::invalid_argument ,
211 ">>> ERROR (Intrepid::PointTools::getEquispacedLatticeLine): order must be positive" );
213 const Scalar h = 1.0 / order;
216 for (
int i=offset;i<=order-offset;i++) {
217 for (
int j=offset;j<=order-i-offset;j++) {
218 points(cur,0) = (Scalar)0.0 + (Scalar) j * h ;
219 points(cur,1) = (Scalar)0.0 + (Scalar) i * h;
227 template<
class Scalar,
class ArrayType>
232 TEUCHOS_TEST_FOR_EXCEPTION( (order <= 0) ,
233 std::invalid_argument ,
234 ">>> ERROR (Intrepid::PointTools::getEquispacedLatticeTetrahedron): order must be positive" );
236 const Scalar h = 1.0 / order;
239 for (
int i=offset;i<=order-offset;i++) {
240 for (
int j=offset;j<=order-i-offset;j++) {
241 for (
int k=offset;k<=order-i-j-offset;k++) {
242 points(cur,0) = (Scalar) k * h;
243 points(cur,1) = (Scalar) j * h;
244 points(cur,2) = (Scalar) i * h;
253 template<
class Scalar,
class ArrayType>
258 Scalar *z =
new Scalar[order+1];
259 Scalar *w =
new Scalar[order+1];
265 for (
int i=offset;i<=order-offset;i++) {
266 points(i-offset,0) = z[i];
275 template<
class Scalar,
class ArrayType>
277 const ArrayType &xnodes ,
278 const ArrayType &xout ,
281 TEUCHOS_TEST_FOR_EXCEPTION( ( warp.dimension(0) != xout.dimension(0) ) ,
282 std::invalid_argument ,
283 ">>> ERROR (PointTools::warpFactor): xout and warp must be same size." );
291 PointTools::getEquispacedLatticeLine<Scalar,ArrayType>( xeq , order , 0 );
292 xeq.resize( order + 1 );
294 TEUCHOS_TEST_FOR_EXCEPTION( ( xeq.dimension(0) != xnodes.dimension(0) ) ,
295 std::invalid_argument ,
296 ">>> ERROR (PointTools::warpFactor): xeq and xnodes must be same size." );
298 for (
int i=0;i<=order;i++) {
300 for (
int k=0;k<d.dimension(0);k++) {
301 d(k) = xnodes(i) - xeq(i);
304 for (
int j=1;j<order;j++) {
306 for (
int k=0;k<d.dimension(0);k++) {
307 d(k) = d(k) * ( (xout(k)-xeq(j)) / (xeq(i)-xeq(j)) );
314 for (
int k=0;k<d.dimension(0);k++) {
315 d(k) = -d(k) / (xeq(i) - xeq(0));
320 for (
int k=0;k<d.dimension(0);k++) {
321 d(k) = d(k) / (xeq(i) - xeq(order));
325 for (
int k=0;k<d.dimension(0);k++) {
335 template<
class Scalar,
class ArrayType>
344 PointTools::getWarpBlendLatticeLine<Scalar,FieldContainer<Scalar> >( gaussX , order , 0 );
346 gaussX.resize(gaussX.dimension(0));
348 Scalar alpopt[] = {0.0000,0.0000,1.4152,0.1001,0.2751,0.9800,1.0999,
349 1.2832,1.3648, 1.4773, 1.4959, 1.5743, 1.5770, 1.6223, 1.6258};
353 if (order >= 1 && order < 16) {
354 alpha = alpopt[order-1];
361 int N = (p+1)*(p+2)/2;
371 for (
int n=1;n<=p+1;n++) {
372 for (
int m=1;m<=p+2-n;m++) {
373 L1(sk) = (n-1) / (Scalar)p;
374 L3(sk) = (m-1) / (Scalar)p;
375 L2(sk) = 1.0 - L1(sk) - L3(sk);
380 for (
int n=0;n<N;n++) {
381 X(n) = -L2(n) + L3(n);
382 Y(n) = (-L2(n) - L3(n) + 2*L1(n))/1.7320508075688772;
390 for (
int n=0;n<N;n++) {
391 blend1(n) = 4.0 * L2(n) * L3(n);
392 blend2(n) = 4.0 * L1(n) * L3(n);
393 blend3(n) = 4.0 * L1(n) * L2(n);
401 for (
int k=0;k<N;k++) {
402 L3mL2(k) = L3(k)-L2(k);
403 L1mL3(k) = L1(k)-L3(k);
404 L2mL1(k) = L2(k)-L1(k);
411 warpFactor<Scalar,FieldContainer<Scalar> >( order , gaussX , L3mL2 , warpfactor1 );
412 warpFactor<Scalar,FieldContainer<Scalar> >( order , gaussX , L1mL3 , warpfactor2 );
413 warpFactor<Scalar,FieldContainer<Scalar> >( order , gaussX , L2mL1 , warpfactor3 );
419 for (
int k=0;k<N;k++) {
420 warp1(k) = blend1(k) * warpfactor1(k) *
421 ( 1.0 + alpha * alpha * L1(k) * L1(k) );
422 warp2(k) = blend2(k) * warpfactor2(k) *
423 ( 1.0 + alpha * alpha * L2(k) * L2(k) );
424 warp3(k) = blend3(k) * warpfactor3(k) *
425 ( 1.0 + alpha * alpha * L3(k) * L3(k) );
428 for (
int k=0;k<N;k++) {
429 X(k) += 1.0 * warp1(k) + cos( 2.0 * M_PI / 3.0 ) * warp2(k) + cos(4*M_PI/3.0) * warp3(k);
430 Y(k) += 0.0 * warp1(k) + sin( 2.0 * M_PI / 3.0 ) * warp2(k) + sin( 4*M_PI/3.0) * warp3(k);
435 for (
int k=0;k<N;k++) {
443 warburtonVerts(0,0,0) = -1.0;
444 warburtonVerts(0,0,1) = -1.0/sqrt(3.0);
445 warburtonVerts(0,1,0) = 1.0;
446 warburtonVerts(0,1,1) = -1.0/sqrt(3.0);
447 warburtonVerts(0,2,0) = 0.0;
448 warburtonVerts(0,2,1) = 2.0/sqrt(3.0);
455 shards::getCellTopologyData< shards::Triangle<3> >(),
461 for (
int i=0;i<=order;i++) {
462 for (
int j=0;j<=order-i;j++) {
463 if ( (i >= offset) && (i <= order-offset) &&
464 (j >= offset) && (j <= order-i-offset) ) {
465 points(offcur,0) = refPts(noffcur,0);
466 points(offcur,1) = refPts(noffcur,1);
477 template<
class Scalar,
class ArrayType>
486 evalshift<Scalar,ArrayType>(order,pval,L2,L3,L4,dxy);
490 template<
class Scalar,
class ArrayType>
493 const ArrayType &L1 ,
494 const ArrayType &L2 ,
495 const ArrayType &L3 ,
500 PointTools::getWarpBlendLatticeLine<Scalar,FieldContainer<Scalar> >( gaussX , order , 0 );
501 gaussX.resize(order+1);
502 const int N = L1.dimension(0);
505 for (
int k=0;k<=order;k++) {
514 for (
int i=0;i<N;i++) {
515 blend1(i) = L2(i) * L3(i);
516 blend2(i) = L1(i) * L3(i);
517 blend3(i) = L1(i) * L2(i);
530 for (
int k=0;k<N;k++) {
531 L3mL2(k) = L3(k)-L2(k);
532 L1mL3(k) = L1(k)-L3(k);
533 L2mL1(k) = L2(k)-L1(k);
536 evalwarp<Scalar,FieldContainer<Scalar> >( warpfactor1 , order , gaussX , L3mL2 );
537 evalwarp<Scalar,FieldContainer<Scalar> >( warpfactor2 , order , gaussX , L1mL3 );
538 evalwarp<Scalar,FieldContainer<Scalar> >( warpfactor3 , order , gaussX , L2mL1 );
540 for (
int k=0;k<N;k++) {
541 warpfactor1(k) *= 4.0;
542 warpfactor2(k) *= 4.0;
543 warpfactor3(k) *= 4.0;
550 for (
int k=0;k<N;k++) {
551 warp1(k) = blend1(k) * warpfactor1(k) *
552 ( 1.0 + pval * pval * L1(k) * L1(k) );
553 warp2(k) = blend2(k) * warpfactor2(k) *
554 ( 1.0 + pval * pval * L2(k) * L2(k) );
555 warp3(k) = blend3(k) * warpfactor3(k) *
556 ( 1.0 + pval * pval * L3(k) * L3(k) );
559 for (
int k=0;k<N;k++) {
560 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);
561 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);
569 template<
class Scalar,
class ArrayType>
572 const ArrayType &xnodes ,
573 const ArrayType &xout )
580 for (
int i=0;i<=order;i++) {
581 xeq(i) = -1.0 + 2.0 * ( order - i ) / order;
586 for (
int i=0;i<=order;i++) {
588 for (
int j=1;j<order;j++) {
590 for (
int k=0;k<d.dimension(0);k++) {
591 d(k) = d(k) * (xout(k)-xeq(j))/(xeq(i)-xeq(j));
596 for (
int k=0;k<d.dimension(0);k++) {
597 d(k) = -d(k)/(xeq(i)-xeq(0));
601 for (
int k=0;k<d.dimension(0);k++) {
602 d(k) = d(k)/(xeq(i)-xeq(order));
606 for (
int k=0;k<d.dimension(0);k++) {
615 template<
class Scalar,
class ArrayType>
620 Scalar alphastore[] = { 0,0,0,0.1002, 1.1332,1.5608,1.3413,1.2577,1.1603,
621 1.10153,0.6080,0.4523,0.8856,0.8717,0.9655};
625 alpha = alphastore[order-1];
631 const int N = (order+1)*(order+2)*(order+3)/6;
640 for (
int n=0;n<=order;n++) {
641 for (
int m=0;m<=order-n;m++) {
642 for (
int q=0;q<=order-n-m;q++) {
643 equipoints(sk,0) = -1.0 + (q * 2.0 ) / order;
644 equipoints(sk,1) = -1.0 + (m * 2.0 ) / order;
645 equipoints(sk,2) = -1.0 + (n * 2.0 ) / order;
657 for (
int i=0;i<N;i++) {
658 L1(i) = (1.0 + equipoints(i,2)) / 2.0;
659 L2(i) = (1.0 + equipoints(i,1)) / 2.0;
660 L3(i) = -(1.0 + equipoints(i,0) + equipoints(i,1) + equipoints(i,2)) / 2.0;
661 L4(i) = (1.0 + equipoints(i,0)) / 2.0;
667 warVerts(0,0) = -1.0;
668 warVerts(0,1) = -1.0/sqrt(3.0);
669 warVerts(0,2) = -1.0/sqrt(6.0);
671 warVerts(1,1) = -1.0/sqrt(3.0);
672 warVerts(1,2) = -1.0/sqrt(6.0);
674 warVerts(2,1) = 2.0 / sqrt(3.0);
675 warVerts(2,2) = -1.0/sqrt(6.0);
678 warVerts(3,2) = 3.0 / sqrt(6.0);
684 for (
int i=0;i<3;i++) {
685 t1(0,i) = warVerts(1,i) - warVerts(0,i);
686 t1(1,i) = warVerts(1,i) - warVerts(0,i);
687 t1(2,i) = warVerts(2,i) - warVerts(1,i);
688 t1(3,i) = warVerts(2,i) - warVerts(0,i);
689 t2(0,i) = warVerts(2,i) - 0.5 * ( warVerts(0,i) + warVerts(1,i) );
690 t2(1,i) = warVerts(3,i) - 0.5 * ( warVerts(0,i) + warVerts(1,i) );
691 t2(2,i) = warVerts(3,i) - 0.5 * ( warVerts(1,i) + warVerts(2,i) );
692 t2(3,i) = warVerts(3,i) - 0.5 * ( warVerts(0,i) + warVerts(2,i) );
696 for (
int n=0;n<4;n++) {
698 Scalar normt1n = 0.0;
699 Scalar normt2n = 0.0;
700 for (
int i=0;i<3;i++) {
701 normt1n += (t1(n,i) * t1(n,i));
702 normt2n += (t2(n,i) * t2(n,i));
704 normt1n = sqrt(normt1n);
705 normt2n = sqrt(normt2n);
707 for (
int i=0;i<3;i++) {
715 for (
int i=0;i<N;i++) {
716 for (
int j=0;j<3;j++) {
717 XYZ(i,j) = L3(i)*warVerts(0,j) + L4(i)*warVerts(1,j) + L2(i)*warVerts(2,j) + L1(i)*warVerts(3,j);
721 for (
int face=1;face<=4;face++) {
728 La = L1; Lb = L2; Lc = L3; Ld = L4;
break;
730 La = L2; Lb = L1; Lc = L3; Ld = L4;
break;
732 La = L3; Lb = L1; Lc = L4; Ld = L2;
break;
734 La = L4; Lb = L1; Lc = L3; Ld = L2;
break;
738 warpShiftFace3D<Scalar,ArrayType>(order,alpha,La,Lb,Lc,Ld,warp);
740 for (
int k=0;k<N;k++) {
741 blend(k) = Lb(k) * Lc(k) * Ld(k);
744 for (
int k=0;k<N;k++) {
745 denom(k) = (Lb(k) + 0.5 * La(k)) * (Lc(k) + 0.5*La(k)) * (Ld(k) + 0.5 * La(k));
748 for (
int k=0;k<N;k++) {
749 if (denom(k) > tol) {
750 blend(k) *= ( 1.0 + alpha * alpha * La(k) * La(k) ) / denom(k);
756 for (
int k=0;k<N;k++) {
757 for (
int j=0;j<3;j++) {
758 shift(k,j) = shift(k,j) + blend(k) * warp(k,0) * t1(face-1,j)
759 + blend(k) * warp(k,1) * t2(face-1,j);
763 for (
int k=0;k<N;k++) {
764 if (La(k) < tol && ( Lb(k) < tol || Lc(k) < tol || Ld(k) < tol )) {
765 for (
int j=0;j<3;j++) {
766 shift(k,j) = warp(k,0) * t1(face-1,j) + warp(k,1) * t2(face-1,j);
774 for (
int k=0;k<N;k++) {
775 for (
int j=0;j<3;j++) {
776 updatedPoints(k,j) = XYZ(k,j) + shift(k,j);
780 warVerts.
resize( 1 , 4 , 3 );
786 shards::getCellTopologyData<shards::Tetrahedron<4> >() ,
792 for (
int i=0;i<=order;i++) {
793 for (
int j=0;j<=order-i;j++) {
794 for (
int k=0;k<=order-i-j;k++) {
795 if ( (i >= offset) && (i <= order-offset) &&
796 (j >= offset) && (j <= order-i-offset) &&
797 (k >= offset) && (k <= order-i-j-offset) ) {
798 points(offcur,0) = refPts(noffcur,0);
799 points(offcur,1) = refPts(noffcur,1);
800 points(offcur,2) = refPts(noffcur,2);
static void zwglj(Scalar *z, Scalar *w, const int np, const Scalar alpha, const Scalar beta)
Gauss-Lobatto-Jacobi zeros and weights with end point at z=-1,1.
void resize(const int dim0)
Resizes FieldContainer to a rank-1 container with the specified dimension, initialized by 0...
void initialize(const Scalar value=0)
Initializes a field container by assigning value to all its elements.
static void zwgj(Scalar *z, Scalar *w, const int np, const Scalar alpha, const Scalar beta)
Gauss-Jacobi zeros and weights.