Intrepid2
Intrepid2_PointToolsDef.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 
48 #ifndef __INTREPID2_POINTTOOLS_DEF_HPP__
49 #define __INTREPID2_POINTTOOLS_DEF_HPP__
50 
51 #ifdef _MSC_VER
52 #include "winmath.h"
53 #endif
54 
55 
56 namespace Intrepid2 {
57 
58 // -----------------------------------------------------------------------------------------
59 // Front interface
60 // -----------------------------------------------------------------------------------------
61 
62 inline
63 ordinal_type
65 getLatticeSize( const shards::CellTopology cellType,
66  const ordinal_type order,
67  const ordinal_type offset ) {
68 #ifdef HAVE_INTREPID2_DEBUG
69  INTREPID2_TEST_FOR_EXCEPTION( order < 0 || offset < 0,
70  std::invalid_argument ,
71  ">>> ERROR (PointTools::getLatticeSize): order and offset must be positive values." );
72 #endif
73  ordinal_type r_val = 0;
74  switch (cellType.getBaseKey()) {
75  case shards::Tetrahedron<>::key: {
76  const auto effectiveOrder = order - 4 * offset;
77  r_val = (effectiveOrder < 0 ? 0 :(effectiveOrder+1)*(effectiveOrder+2)*(effectiveOrder+3)/6);
78  break;
79  }
80  case shards::Triangle<>::key: {
81  const auto effectiveOrder = order - 3 * offset;
82  r_val = (effectiveOrder < 0 ? 0 : (effectiveOrder+1)*(effectiveOrder+2)/2);
83  break;
84  }
85  case shards::Line<>::key: {
86  const auto effectiveOrder = order - 2 * offset;
87  r_val = (effectiveOrder < 0 ? 0 : (effectiveOrder+1));
88  break;
89  }
90  default: {
91  INTREPID2_TEST_FOR_EXCEPTION( true , std::invalid_argument ,
92  ">>> ERROR (Intrepid2::PointTools::getLatticeSize): the specified cell type is not supported." );
93  }
94  }
95  return r_val;
96 }
97 
98 template<typename pointValueType, class ...pointProperties>
99 void
101 getLattice( Kokkos::DynRankView<pointValueType,pointProperties...> points,
102  const shards::CellTopology cell,
103  const ordinal_type order,
104  const ordinal_type offset,
105  const EPointType pointType ) {
106 #ifdef HAVE_INTREPID2_DEBUG
107  INTREPID2_TEST_FOR_EXCEPTION( points.rank() != 2,
108  std::invalid_argument ,
109  ">>> ERROR (PointTools::getLattice): points rank must be 2." );
110  INTREPID2_TEST_FOR_EXCEPTION( order < 0 || offset < 0,
111  std::invalid_argument ,
112  ">>> ERROR (PointTools::getLattice): order and offset must be positive values." );
113 
114  const size_type latticeSize = getLatticeSize( cell, order, offset );
115  const size_type spaceDim = cell.getDimension();
116 
117  INTREPID2_TEST_FOR_EXCEPTION( points.extent(0) != latticeSize ||
118  points.extent(1) != spaceDim,
119  std::invalid_argument ,
120  ">>> ERROR (PointTools::getLattice): dimension does not match to lattice size." );
121 #endif
122 
123  // const auto latticeSize = getLatticeSize( cell, order, offset );
124  // const auto spaceDim = cell.getDimension();
125 
126  // // the interface assumes that the input array follows the cell definition
127  // // so, let's match all dimensions according to the cell specification
128  // typedef Kokkos::pair<ordinal_type,ordinal_type> range_type;
129  // auto pts = Kokkos::subview( points,
130  // range_type(0, latticeSize),
131  // range_type(0, spaceDim) );
132  switch (pointType) {
133  case POINTTYPE_EQUISPACED: getEquispacedLattice( points, cell, order, offset ); break;
134  case POINTTYPE_WARPBLEND: getWarpBlendLattice ( points, cell, order, offset ); break;
135  default: {
136  INTREPID2_TEST_FOR_EXCEPTION( true ,
137  std::invalid_argument ,
138  ">>> ERROR (PointTools::getLattice): invalid EPointType." );
139  }
140  }
141 }
142 
143 template<typename pointValueType, class ...pointProperties>
144 void
146 getGaussPoints( Kokkos::DynRankView<pointValueType,pointProperties...> points,
147  const ordinal_type order ) {
148 #ifdef HAVE_INTREPID2_DEBUG
149  INTREPID2_TEST_FOR_EXCEPTION( points.rank() != 2,
150  std::invalid_argument ,
151  ">>> ERROR (PointTools::getGaussPoints): points rank must be 1." );
152  INTREPID2_TEST_FOR_EXCEPTION( order < 0,
153  std::invalid_argument ,
154  ">>> ERROR (PointTools::getGaussPoints): order must be positive value." );
155 #endif
156  const ordinal_type np = order + 1;
157  const double alpha = 0.0, beta = 0.0;
158 
159  // until view and dynrankview inter-operatible, we use views in a consistent way
160  Kokkos::View<pointValueType*,Kokkos::HostSpace>
161  zHost("PointTools::getGaussPoints::z", np),
162  wHost("PointTools::getGaussPoints::w", np);
163 
164  // sequential means that the code is decorated with KOKKOS_INLINE_FUNCTION
165  // and it does not invoke parallel for inside (cheap operation), which means
166  // that gpu memory is not accessible unless this is called inside of functor.
167  Polylib::Serial::Cubature<POLYTYPE_GAUSS>::getValues(zHost, wHost, np, alpha, beta);
168 
169  typedef Kokkos::pair<ordinal_type,ordinal_type> range_type;
170  auto pts = Kokkos::subview( points, range_type(0,np), 0 );
171  // should be fixed after view and dynrankview are inter-operatible
172  auto z = Kokkos::DynRankView<pointValueType,Kokkos::HostSpace>(zHost.data(), np);
173  Kokkos::deep_copy(pts, z);
174 }
175 
176 // -----------------------------------------------------------------------------------------
177 // Internal implementation
178 // -----------------------------------------------------------------------------------------
179 
180 template<typename pointValueType, class ...pointProperties>
181 void PointTools::
182 getEquispacedLattice( Kokkos::DynRankView<pointValueType,pointProperties...> points,
183  const shards::CellTopology cell,
184  const ordinal_type order,
185  const ordinal_type offset ) {
186  switch (cell.getBaseKey()) {
187  case shards::Tetrahedron<>::key: getEquispacedLatticeTetrahedron( points, order, offset ); break;
188  case shards::Triangle<>::key: getEquispacedLatticeTriangle ( points, order, offset ); break;
189  case shards::Line<>::key: getEquispacedLatticeLine ( points, order, offset ); break;
190  default: {
191  INTREPID2_TEST_FOR_EXCEPTION( true , std::invalid_argument ,
192  ">>> ERROR (Intrepid2::PointTools::getEquispacedLattice): the specified cell type is not supported." );
193  }
194  }
195 
196 }
197 
198 template<typename pointValueType, class ...pointProperties>
199 void PointTools::
200 getWarpBlendLattice( Kokkos::DynRankView<pointValueType,pointProperties...> points,
201  const shards::CellTopology cell,
202  const ordinal_type order,
203  const ordinal_type offset ) {
204  switch (cell.getBaseKey()) {
205  case shards::Tetrahedron<>::key: getWarpBlendLatticeTetrahedron( points, order, offset ); break;
206  case shards::Triangle<>::key: getWarpBlendLatticeTriangle ( points, order, offset ); break;
207  case shards::Line<>::key: getWarpBlendLatticeLine ( points, order, offset ); break;
208  default: {
209  INTREPID2_TEST_FOR_EXCEPTION( true , std::invalid_argument ,
210  ">>> ERROR (Intrepid2::PointTools::getWarpBlendLattice): the specified cell type is not supported." );
211  }
212  }
213 }
214 
215 // -----------------------------------------------------------------------------------------
216 
217 template<typename pointValueType, class ...pointProperties>
218 void
220 getEquispacedLatticeLine( Kokkos::DynRankView<pointValueType,pointProperties...> points,
221  const ordinal_type order,
222  const ordinal_type offset ) {
223  auto pointsHost = Kokkos::create_mirror_view(Kokkos::HostSpace::memory_space(), points);
224 
225  if (order == 0)
226  pointsHost(0,0) = 0.0;
227  else {
228  const pointValueType h = 2.0 / order;
229  const ordinal_type ibeg = offset, iend = order-offset+1;
230  for (ordinal_type i=ibeg;i<iend;++i)
231  pointsHost(i-ibeg, 0) = -1.0 + h * (pointValueType) i;
232  }
233 
234  Kokkos::deep_copy(points, pointsHost);
235 }
236 
237 template<typename pointValueType, class ...pointProperties>
238 void
240 getWarpBlendLatticeLine( Kokkos::DynRankView<pointValueType,pointProperties...> points,
241  const ordinal_type order,
242  const ordinal_type offset ) {
243  // order is order of polynomial degree. The Gauss-Lobatto points are accurate
244  // up to degree 2*i-1
245  const ordinal_type np = order + 1;
246  const double alpha = 0.0, beta = 0.0;
247  const ordinal_type s = np - 2*offset;
248 
249  if (s > 0) {
250  // until view and dynrankview inter-operatible, we use views in a consistent way
251  Kokkos::View<pointValueType*,Kokkos::HostSpace>
252  zHost("PointTools::getGaussPoints::z", np),
253  wHost("PointTools::getGaussPoints::w", np);
254 
255  // sequential means that the code is decorated with KOKKOS_INLINE_FUNCTION
256  // and it does not invoke parallel for inside (cheap operation), which means
257  // that gpu memory is not accessible unless this is called inside of functor.
259 
260  typedef Kokkos::pair<ordinal_type,ordinal_type> range_type;
261  auto pts = Kokkos::subview( points, range_type(0, s), 0 );
262 
263  // this should be fixed after view and dynrankview is interoperatable
264  auto z = Kokkos::DynRankView<pointValueType,Kokkos::HostSpace>(zHost.data() + offset, np-offset);
265 
266  Kokkos::deep_copy(pts, z);
267  }
268 }
269 
270 template<typename pointValueType, class ...pointProperties>
271 void
273 getEquispacedLatticeTriangle( Kokkos::DynRankView<pointValueType,pointProperties...> points,
274  const ordinal_type order,
275  const ordinal_type offset ) {
276  TEUCHOS_TEST_FOR_EXCEPTION( order <= 0 ,
277  std::invalid_argument ,
278  ">>> ERROR (Intrepid2::PointTools::getEquispacedLatticeTriangle): order must be positive" );
279 
280  auto pointsHost = Kokkos::create_mirror_view(Kokkos::HostSpace::memory_space(), points);
281 
282  const pointValueType h = 1.0 / order;
283  ordinal_type cur = 0;
284 
285  for (ordinal_type i=offset;i<=order-offset;i++) {
286  for (ordinal_type j=offset;j<=order-i-offset;j++) {
287  pointsHost(cur,0) = j * h ;
288  pointsHost(cur,1) = i * h;
289  cur++;
290  }
291  }
292  Kokkos::deep_copy(points, pointsHost);
293 }
294 
295 template<typename pointValueType, class ...pointProperties>
296 void
298 getEquispacedLatticeTetrahedron( Kokkos::DynRankView<pointValueType,pointProperties...> points,
299  const ordinal_type order ,
300  const ordinal_type offset ) {
301  TEUCHOS_TEST_FOR_EXCEPTION( (order <= 0) ,
302  std::invalid_argument ,
303  ">>> ERROR (Intrepid2::PointTools::getEquispacedLatticeTetrahedron): order must be positive" );
304 
305  auto pointsHost = Kokkos::create_mirror_view(Kokkos::HostSpace::memory_space(), points);
306 
307  const pointValueType h = 1.0 / order;
308  ordinal_type cur = 0;
309 
310  for (ordinal_type i=offset;i<=order-offset;i++) {
311  for (ordinal_type j=offset;j<=order-i-offset;j++) {
312  for (ordinal_type k=offset;k<=order-i-j-offset;k++) {
313  pointsHost(cur,0) = k * h;
314  pointsHost(cur,1) = j * h;
315  pointsHost(cur,2) = i * h;
316  cur++;
317  }
318  }
319  }
320  Kokkos::deep_copy(points, pointsHost);
321 }
322 
323 
324 template<typename pointValueType, class ...pointProperties>
325 void
327 warpFactor( Kokkos::DynRankView<pointValueType,pointProperties...> warp,
328  const ordinal_type order,
329  const Kokkos::DynRankView<pointValueType,pointProperties...> xnodes ,
330  const Kokkos::DynRankView<pointValueType,pointProperties...> xout
331  )
332  {
333  TEUCHOS_TEST_FOR_EXCEPTION( ( warp.extent(0) != xout.extent(0) ) ,
334  std::invalid_argument ,
335  ">>> ERROR (PointTools::warpFactor): xout and warp must be same size." );
336 
337  Kokkos::deep_copy(warp, pointValueType(0.0));
338 
339  ordinal_type xout_dim0 = xout.extent(0);
340  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> d("d", xout_dim0 );
341 
342  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> xeq_("xeq", order + 1 ,1);
343  PointTools::getEquispacedLatticeLine( xeq_ , order , 0 );
344  const auto xeq = Kokkos::subview(xeq_, Kokkos::ALL(),0);
345 
346  TEUCHOS_TEST_FOR_EXCEPTION( ( xeq.extent(0) != xnodes.extent(0) ) ,
347  std::invalid_argument ,
348  ">>> ERROR (PointTools::warpFactor): xeq and xnodes must be same size." );
349 
350 
351  for (ordinal_type i=0;i<=order;i++) {
352 
353  Kokkos::deep_copy(d, xnodes(i,0) - xeq(i));
354 
355  for (ordinal_type j=1;j<order;j++) {
356  if (i != j) {
357  for (ordinal_type k=0;k<xout_dim0;k++) {
358  d(k) = d(k) * ( (xout(k)-xeq(j)) / (xeq(i)-xeq(j)) );
359  }
360  }
361  }
362 
363  // deflate end roots
364  if ( i != 0 ) {
365  for (ordinal_type k=0;k<xout_dim0;k++) {
366  d(k) = -d(k) / (xeq(i) - xeq(0));
367  }
368  }
369 
370  if (i != order ) {
371  for (ordinal_type k=0;k<xout_dim0;k++) {
372  d(k) = d(k) / (xeq(i) - xeq(order));
373  }
374  }
375 
376  for (ordinal_type k=0;k<xout_dim0;k++) {
377  warp(k) += d(k);
378  }
379 
380  }
381  }
382 
383 template<typename pointValueType, class ...pointProperties>
384 void
386 getWarpBlendLatticeTriangle( Kokkos::DynRankView<pointValueType,pointProperties...> points,
387  const ordinal_type order ,
388  const ordinal_type offset)
389  {
390  /* get Gauss-Lobatto points */
391 
392  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> gaussX("gaussX", order + 1 , 1 );
393 
394  PointTools::getWarpBlendLatticeLine( gaussX , order , 0 );
395  //auto gaussX = Kokkos::subdynrankview(gaussX_, Kokkos::ALL(),0);
396 
397  // gaussX.resize(gaussX.extent(0));
398 
399  pointValueType alpopt[] = {0.0000,0.0000,1.4152,0.1001,0.2751,0.9800,1.0999,
400  1.2832,1.3648, 1.4773, 1.4959, 1.5743, 1.5770, 1.6223, 1.6258};
401 
402  pointValueType alpha;
403 
404  if (order >= 1 && order < 16) {
405  alpha = alpopt[order-1];
406  }
407  else {
408  alpha = 5.0 / 3.0;
409  }
410 
411  const ordinal_type p = order; /* switch to Warburton's notation */
412  ordinal_type N = (p+1)*(p+2)/2;
413 
414  /* equidistributed nodes on equilateral triangle */
415  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> L1("L1", N );
416  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> L2("L2", N );
417  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> L3("L3", N );
418  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> X("X", N);
419  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> Y("Y", N);
420 
421  ordinal_type sk = 0;
422  for (ordinal_type n=1;n<=p+1;n++) {
423  for (ordinal_type m=1;m<=p+2-n;m++) {
424  L1(sk) = (n-1) / (pointValueType)p;
425  L3(sk) = (m-1) / (pointValueType)p;
426  L2(sk) = 1.0 - L1(sk) - L3(sk);
427  sk++;
428  }
429  }
430 
431  for (ordinal_type n=0;n<N;n++) {
432  X(n) = -L2(n) + L3(n);
433  Y(n) = (-L2(n) - L3(n) + 2*L1(n))/1.7320508075688772;
434  }
435 
436  /* get blending function for each node at each edge */
437  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> blend1("blend1", N);
438  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> blend2("blend2", N);
439  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> blend3("blend3", N);
440 
441  for (ordinal_type n=0;n<N;n++) {
442  blend1(n) = 4.0 * L2(n) * L3(n);
443  blend2(n) = 4.0 * L1(n) * L3(n);
444  blend3(n) = 4.0 * L1(n) * L2(n);
445  }
446 
447  /* get difference of each barycentric coordinate */
448  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> L3mL2("L3mL2", N);
449  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> L1mL3("L1mL3", N);
450  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> L2mL1("L2mL1", N);
451 
452  for (ordinal_type k=0;k<N;k++) {
453  L3mL2(k) = L3(k)-L2(k);
454  L1mL3(k) = L1(k)-L3(k);
455  L2mL1(k) = L2(k)-L1(k);
456  }
457 
458  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> warpfactor1("warpfactor1", N);
459  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> warpfactor2("warpfactor2", N);
460  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> warpfactor3("warpfactor3", N);
461 
462  warpFactor( warpfactor1, order , gaussX , L3mL2 );
463  warpFactor( warpfactor2, order , gaussX , L1mL3 );
464  warpFactor( warpfactor3, order , gaussX , L2mL1 );
465 
466  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> warp1("warp1", N);
467  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> warp2("warp2", N);
468  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> warp3("warp3", N);
469 
470  for (ordinal_type k=0;k<N;k++) {
471  warp1(k) = blend1(k) * warpfactor1(k) *
472  ( 1.0 + alpha * alpha * L1(k) * L1(k) );
473  warp2(k) = blend2(k) * warpfactor2(k) *
474  ( 1.0 + alpha * alpha * L2(k) * L2(k) );
475  warp3(k) = blend3(k) * warpfactor3(k) *
476  ( 1.0 + alpha * alpha * L3(k) * L3(k) );
477  }
478 
479  for (ordinal_type k=0;k<N;k++) {
480  X(k) += 1.0 * warp1(k) + cos( 2.0 * M_PI / 3.0 ) * warp2(k) + cos(4*M_PI/3.0) * warp3(k);
481  Y(k) += 0.0 * warp1(k) + sin( 2.0 * M_PI / 3.0 ) * warp2(k) + sin( 4*M_PI/3.0) * warp3(k);
482  }
483 
484  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> warXY("warXY", 1, N,2);
485 
486  for (ordinal_type k=0;k<N;k++) {
487  warXY(0, k,0) = X(k);
488  warXY(0, k,1) = Y(k);
489  }
490 
491 
492  // finally, convert the warp-blend points to the correct triangle
493  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> warburtonVerts("warburtonVerts", 1,3,2);
494  warburtonVerts(0,0,0) = -1.0;
495  warburtonVerts(0,0,1) = -1.0/std::sqrt(3.0);
496  warburtonVerts(0,1,0) = 1.0;
497  warburtonVerts(0,1,1) = -1.0/std::sqrt(3.0);
498  warburtonVerts(0,2,0) = 0.0;
499  warburtonVerts(0,2,1) = 2.0/std::sqrt(3.0);
500 
501  Kokkos::DynRankView<pointValueType, Kokkos::DefaultHostExecutionSpace> refPts("refPts", 1, N,2);
502 
504  warXY ,
505  warburtonVerts ,
506  shards::getCellTopologyData< shards::Triangle<3> >()
507  );
508 
509 
510  auto pointsHost = Kokkos::create_mirror_view(Kokkos::HostSpace::memory_space(), points);
511  // now write from refPts into points, taking care of offset
512  ordinal_type noffcur = 0; // index into refPts
513  ordinal_type offcur = 0; // index ordinal_type points
514  for (ordinal_type i=0;i<=order;i++) {
515  for (ordinal_type j=0;j<=order-i;j++) {
516  if ( (i >= offset) && (i <= order-offset) &&
517  (j >= offset) && (j <= order-i-offset) ) {
518  pointsHost(offcur,0) = refPts(0, noffcur,0);
519  pointsHost(offcur,1) = refPts(0, noffcur,1);
520  offcur++;
521  }
522  noffcur++;
523  }
524  }
525 
526  Kokkos::deep_copy(points, pointsHost);
527 
528  }
529 
530 
531 template<typename pointValueType, class ...pointProperties>
532 void
534 warpShiftFace3D( Kokkos::DynRankView<pointValueType,pointProperties...> dxy,
535  const ordinal_type order ,
536  const pointValueType pval ,
537  const Kokkos::DynRankView<pointValueType,pointProperties...> /* L1 */,
538  const Kokkos::DynRankView<pointValueType,pointProperties...> L2,
539  const Kokkos::DynRankView<pointValueType,pointProperties...> L3,
540  const Kokkos::DynRankView<pointValueType,pointProperties...> L4
541  )
542  {
543  evalshift(dxy,order,pval,L2,L3,L4);
544  return;
545  }
546 
547 template<typename pointValueType, class ...pointProperties>
548 void
550 evalshift( Kokkos::DynRankView<pointValueType,pointProperties...> dxy,
551  const ordinal_type order ,
552  const pointValueType pval ,
553  const Kokkos::DynRankView<pointValueType,pointProperties...> L1 ,
554  const Kokkos::DynRankView<pointValueType,pointProperties...> L2 ,
555  const Kokkos::DynRankView<pointValueType,pointProperties...> L3
556  )
557  {
558  // get Gauss-Lobatto-nodes
559  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> gaussX("gaussX",order+1,1);
560  PointTools::getWarpBlendLatticeLine( gaussX , order , 0 );
561  //gaussX.resize(order+1);
562  const ordinal_type N = L1.extent(0);
563 
564  // Warburton code reverses them
565  for (ordinal_type k=0;k<=order;k++) {
566  gaussX(k,0) *= -1.0;
567  }
568 
569  // blending function at each node for each edge
570  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> blend1("blend1",N);
571  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> blend2("blend2",N);
572  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> blend3("blend3",N);
573 
574  for (ordinal_type i=0;i<N;i++) {
575  blend1(i) = L2(i) * L3(i);
576  blend2(i) = L1(i) * L3(i);
577  blend3(i) = L1(i) * L2(i);
578  }
579 
580  // amount of warp for each node for each edge
581  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> warpfactor1("warpfactor1",N);
582  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> warpfactor2("warpfactor2",N);
583  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> warpfactor3("warpfactor3",N);
584 
585  // differences of barycentric coordinates
586  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> L3mL2("L3mL2",N);
587  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> L1mL3("L1mL3",N);
588  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> L2mL1("L2mL1",N);
589 
590  for (ordinal_type k=0;k<N;k++) {
591  L3mL2(k) = L3(k)-L2(k);
592  L1mL3(k) = L1(k)-L3(k);
593  L2mL1(k) = L2(k)-L1(k);
594  }
595 
596  evalwarp( warpfactor1 , order , gaussX , L3mL2 );
597  evalwarp( warpfactor2 , order , gaussX , L1mL3 );
598  evalwarp( warpfactor3 , order , gaussX , L2mL1 );
599 
600  for (ordinal_type k=0;k<N;k++) {
601  warpfactor1(k) *= 4.0;
602  warpfactor2(k) *= 4.0;
603  warpfactor3(k) *= 4.0;
604  }
605 
606  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> warp1("warp1",N);
607  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> warp2("warp2",N);
608  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> warp3("warp3",N);
609 
610  for (ordinal_type k=0;k<N;k++) {
611  warp1(k) = blend1(k) * warpfactor1(k) *
612  ( 1.0 + pval * pval * L1(k) * L1(k) );
613  warp2(k) = blend2(k) * warpfactor2(k) *
614  ( 1.0 + pval * pval * L2(k) * L2(k) );
615  warp3(k) = blend3(k) * warpfactor3(k) *
616  ( 1.0 + pval * pval * L3(k) * L3(k) );
617  }
618 
619  for (ordinal_type k=0;k<N;k++) {
620  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);
621  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);
622  }
623  }
624 
625  /* one-d edge warping function */
626 template<typename pointValueType, class ...pointProperties>
627 void
629 evalwarp(Kokkos::DynRankView<pointValueType,pointProperties...> warp ,
630  const ordinal_type order ,
631  const Kokkos::DynRankView<pointValueType,pointProperties...> xnodes,
632  const Kokkos::DynRankView<pointValueType,pointProperties...> xout )
633  {
634  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> xeq("xeq",order+1);
635 
636  ordinal_type xout_dim0 = xout.extent(0);
637  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> d("d",xout_dim0);
638 
639  //Kokkos::deep_copy(d, 0.0);
640 
641  for (ordinal_type i=0;i<=order;i++) {
642  xeq(i) = -1.0 + 2.0 * ( order - i ) / order;
643  }
644 
645  for (ordinal_type i=0;i<=order;i++) {
646  Kokkos::deep_copy(d, xnodes(i,0) - xeq(i));
647  for (ordinal_type j=1;j<order;j++) {
648  if (i!=j) {
649  for (ordinal_type k=0;k<xout_dim0;k++) {
650  d(k) = d(k) * (xout(k)-xeq(j))/(xeq(i)-xeq(j));
651  }
652  }
653  }
654  if (i!=0) {
655  for (ordinal_type k=0;k<xout_dim0;k++) {
656  d(k) = -d(k)/(xeq(i)-xeq(0));
657  }
658  }
659  if (i!=order) {
660  for (ordinal_type k=0;k<xout_dim0;k++) {
661  d(k) = d(k)/(xeq(i)-xeq(order));
662  }
663  }
664 
665  for (ordinal_type k=0;k<xout_dim0;k++) {
666  warp(k) += d(k);
667  }
668  }
669  }
670 
671 
672 template<typename pointValueType, class ...pointProperties>
673 void
675 getWarpBlendLatticeTetrahedron(Kokkos::DynRankView<pointValueType,pointProperties...> points,
676  const ordinal_type order ,
677  const ordinal_type offset )
678  {
679  pointValueType alphastore[] = { 0,0,0,0.1002, 1.1332,1.5608,1.3413,1.2577,1.1603,
680  1.10153,0.6080,0.4523,0.8856,0.8717,0.9655};
681  pointValueType alpha;
682 
683  if (order <= 15) {
684  alpha = alphastore[std::max(order-1,0)];
685  }
686  else {
687  alpha = 1.0;
688  }
689 
690  const ordinal_type N = (order+1)*(order+2)*(order+3)/6;
691  pointValueType tol = 1.e-10;
692 
693  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> shift("shift",N,3);
694  Kokkos::deep_copy(shift,0.0);
695 
696  /* create 3d equidistributed nodes on Warburton tet */
697  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> equipoints("equipoints", N,3);
698  ordinal_type sk = 0;
699  for (ordinal_type n=0;n<=order;n++) {
700  for (ordinal_type m=0;m<=order-n;m++) {
701  for (ordinal_type q=0;q<=order-n-m;q++) {
702  equipoints(sk,0) = -1.0 + (q * 2.0 ) / order;
703  equipoints(sk,1) = -1.0 + (m * 2.0 ) / order;
704  equipoints(sk,2) = -1.0 + (n * 2.0 ) / order;
705  sk++;
706  }
707  }
708  }
709 
710 
711  /* construct barycentric coordinates for equispaced lattice */
712  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> L1("L1",N);
713  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> L2("L2",N);
714  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> L3("L3",N);
715  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> L4("L4",N);
716  for (ordinal_type i=0;i<N;i++) {
717  L1(i) = (1.0 + equipoints(i,2)) / 2.0;
718  L2(i) = (1.0 + equipoints(i,1)) / 2.0;
719  L3(i) = -(1.0 + equipoints(i,0) + equipoints(i,1) + equipoints(i,2)) / 2.0;
720  L4(i) = (1.0 + equipoints(i,0)) / 2.0;
721  }
722 
723 
724  /* vertices of Warburton tet */
725  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> warVerts_("warVerts",1,4,3);
726  auto warVerts = Kokkos::subview(warVerts_,0,Kokkos::ALL(),Kokkos::ALL());
727  warVerts(0,0) = -1.0;
728  warVerts(0,1) = -1.0/sqrt(3.0);
729  warVerts(0,2) = -1.0/sqrt(6.0);
730  warVerts(1,0) = 1.0;
731  warVerts(1,1) = -1.0/sqrt(3.0);
732  warVerts(1,2) = -1.0/sqrt(6.0);
733  warVerts(2,0) = 0.0;
734  warVerts(2,1) = 2.0 / sqrt(3.0);
735  warVerts(2,2) = -1.0/sqrt(6.0);
736  warVerts(3,0) = 0.0;
737  warVerts(3,1) = 0.0;
738  warVerts(3,2) = 3.0 / sqrt(6.0);
739 
740 
741  /* tangents to faces */
742  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> t1("t1",4,3);
743  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> t2("t2",4,3);
744  for (ordinal_type i=0;i<3;i++) {
745  t1(0,i) = warVerts(1,i) - warVerts(0,i);
746  t1(1,i) = warVerts(1,i) - warVerts(0,i);
747  t1(2,i) = warVerts(2,i) - warVerts(1,i);
748  t1(3,i) = warVerts(2,i) - warVerts(0,i);
749  t2(0,i) = warVerts(2,i) - 0.5 * ( warVerts(0,i) + warVerts(1,i) );
750  t2(1,i) = warVerts(3,i) - 0.5 * ( warVerts(0,i) + warVerts(1,i) );
751  t2(2,i) = warVerts(3,i) - 0.5 * ( warVerts(1,i) + warVerts(2,i) );
752  t2(3,i) = warVerts(3,i) - 0.5 * ( warVerts(0,i) + warVerts(2,i) );
753  }
754 
755  /* normalize tangents */
756  for (ordinal_type n=0;n<4;n++) {
757  /* Compute norm of t1(n) and t2(n) */
758  pointValueType normt1n = 0.0;
759  pointValueType normt2n = 0.0;
760  for (ordinal_type i=0;i<3;i++) {
761  normt1n += (t1(n,i) * t1(n,i));
762  normt2n += (t2(n,i) * t2(n,i));
763  }
764  normt1n = sqrt(normt1n);
765  normt2n = sqrt(normt2n);
766  /* normalize each tangent now */
767  for (ordinal_type i=0;i<3;i++) {
768  t1(n,i) /= normt1n;
769  t2(n,i) /= normt2n;
770  }
771  }
772 
773  /* undeformed coordinates */
774  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> XYZ("XYZ",N,3);
775  for (ordinal_type i=0;i<N;i++) {
776  for (ordinal_type j=0;j<3;j++) {
777  XYZ(i,j) = L3(i)*warVerts(0,j) + L4(i)*warVerts(1,j) + L2(i)*warVerts(2,j) + L1(i)*warVerts(3,j);
778  }
779  }
780 
781  for (ordinal_type face=1;face<=4;face++) {
782  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> La, Lb, Lc, Ld;
783  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> warp("warp",N,2);
784  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> blend("blend",N);
785  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> denom("denom",N);
786  switch (face) {
787  case 1:
788  La = L1; Lb = L2; Lc = L3; Ld = L4; break;
789  case 2:
790  La = L2; Lb = L1; Lc = L3; Ld = L4; break;
791  case 3:
792  La = L3; Lb = L1; Lc = L4; Ld = L2; break;
793  case 4:
794  La = L4; Lb = L1; Lc = L3; Ld = L2; break;
795  }
796 
797  /* get warp tangential to face */
798  warpShiftFace3D(warp,order,alpha,La,Lb,Lc,Ld);
799 
800  for (ordinal_type k=0;k<N;k++) {
801  blend(k) = Lb(k) * Lc(k) * Ld(k);
802  }
803 
804  for (ordinal_type k=0;k<N;k++) {
805  denom(k) = (Lb(k) + 0.5 * La(k)) * (Lc(k) + 0.5*La(k)) * (Ld(k) + 0.5 * La(k));
806  }
807 
808  for (ordinal_type k=0;k<N;k++) {
809  if (denom(k) > tol) {
810  blend(k) *= ( 1.0 + alpha * alpha * La(k) * La(k) ) / denom(k);
811  }
812  }
813 
814 
815  // compute warp and blend
816  for (ordinal_type k=0;k<N;k++) {
817  for (ordinal_type j=0;j<3;j++) {
818  shift(k,j) = shift(k,j) + blend(k) * warp(k,0) * t1(face-1,j)
819  + blend(k) * warp(k,1) * t2(face-1,j);
820  }
821  }
822 
823  for (ordinal_type k=0;k<N;k++) {
824  if (La(k) < tol && ( Lb(k) < tol || Lc(k) < tol || Ld(k) < tol )) {
825  for (ordinal_type j=0;j<3;j++) {
826  shift(k,j) = warp(k,0) * t1(face-1,j) + warp(k,1) * t2(face-1,j);
827  }
828  }
829  }
830 
831  }
832 
833  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> updatedPoints("updatedPoints",1,N,3);
834  for (ordinal_type k=0;k<N;k++) {
835  for (ordinal_type j=0;j<3;j++) {
836  updatedPoints(0,k,j) = XYZ(k,j) + shift(k,j);
837  }
838  }
839 
840  //warVerts.resize( 1 , 4 , 3 );
841 
842  // now we convert to Pavel's reference triangle!
843  Kokkos::DynRankView<pointValueType,Kokkos::DefaultHostExecutionSpace> refPts("refPts",1,N,3);
845  warVerts_ ,
846  shards::getCellTopologyData<shards::Tetrahedron<4> >()
847  );
848 
849  auto pointsHost = Kokkos::create_mirror_view(Kokkos::HostSpace::memory_space(), points);
850  // now write from refPts into points, taking offset into account
851  ordinal_type noffcur = 0;
852  ordinal_type offcur = 0;
853  for (ordinal_type i=0;i<=order;i++) {
854  for (ordinal_type j=0;j<=order-i;j++) {
855  for (ordinal_type k=0;k<=order-i-j;k++) {
856  if ( (i >= offset) && (i <= order-offset) &&
857  (j >= offset) && (j <= order-i-offset) &&
858  (k >= offset) && (k <= order-i-j-offset) ) {
859  pointsHost(offcur,0) = refPts(0,noffcur,0);
860  pointsHost(offcur,1) = refPts(0,noffcur,1);
861  pointsHost(offcur,2) = refPts(0,noffcur,2);
862  offcur++;
863  }
864  noffcur++;
865  }
866  }
867  }
868 
869  Kokkos::deep_copy(points, pointsHost);
870  }
871 
872 
873 } // face Intrepid
874 #endif
static void evalwarp(Kokkos::DynRankView< pointValueType, pointProperties...> warp, const ordinal_type order, const Kokkos::DynRankView< pointValueType, pointProperties...> xnodes, const Kokkos::DynRankView< pointValueType, pointProperties...> xout)
Used internally to compute the warp on edges of a triangle in warp-blend points.
static void evalshift(Kokkos::DynRankView< pointValueType, pointProperties...> dxy, const ordinal_type order, const pointValueType pval, const Kokkos::DynRankView< pointValueType, pointProperties...> L1, const Kokkos::DynRankView< pointValueType, pointProperties...> L2, const Kokkos::DynRankView< pointValueType, pointProperties...> L3)
Used internally to evaluate the point shift for warp-blend points on faces of tets.
static void getWarpBlendLattice(Kokkos::DynRankView< pointValueType, pointProperties...> points, const shards::CellTopology cellType, const ordinal_type order, const ordinal_type offset=0)
Computes a warped lattice (ala Warburton&#39;s warp-blend points of a given order on a reference simplex ...
static void warpFactor(Kokkos::DynRankView< pointValueType, pointProperties...> warp, const ordinal_type order, const Kokkos::DynRankView< pointValueType, pointProperties...> xnodes, const Kokkos::DynRankView< pointValueType, pointProperties...> xout)
interpolates Warburton&#39;s warp function on the line
Gauss-Jacobi/Gauss-Radau-Jacobi/Gauss-Lobatto zeros and weights.
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 (currently disabled for other ce...
static void getEquispacedLatticeTetrahedron(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0)
Computes an equispaced lattice of a given order on the reference tetrahedron. The output array is (P...
static void getWarpBlendLatticeLine(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0)
Returns the Gauss-Lobatto points of a given order on the reference line [-1,1]. The output array is (...
static void mapToReferenceFrame(Kokkos::DynRankView< refPointValueType, refPointProperties...> refPoints, const Kokkos::DynRankView< physPointValueType, physPointProperties...> physPoints, const Kokkos::DynRankView< worksetCellValueType, worksetCellProperties...> worksetCell, const shards::CellTopology cellTopo)
Computes , the inverse of the reference-to-physical frame map using a default initial guess...
static void getWarpBlendLatticeTriangle(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0)
Returns Warburton&#39;s warp-blend points of a given order on the reference triangle. The output array is...
static void getEquispacedLatticeLine(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0)
Computes an equispaced lattice of a given order on the reference line [-1,1]. The output array is (P...
static void getEquispacedLattice(Kokkos::DynRankView< pointValueType, pointProperties...> points, const shards::CellTopology cellType, const ordinal_type order, const ordinal_type offset=0)
Computes an equispaced lattice of a given order on a reference simplex (currently disabled for other ...
static void getGaussPoints(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order)
static void warpShiftFace3D(Kokkos::DynRankView< pointValueType, pointProperties...> dxy, const ordinal_type order, const pointValueType pval, const Kokkos::DynRankView< pointValueType, pointProperties...> L1, const Kokkos::DynRankView< pointValueType, pointProperties...> L2, const Kokkos::DynRankView< pointValueType, pointProperties...> L3, const Kokkos::DynRankView< pointValueType, pointProperties...> L4)
This is used internally to compute the tetrahedral warp-blend points one each face.
static void getWarpBlendLatticeTetrahedron(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0)
Returns Warburton&#39;s warp-blend points of a given order on the reference tetrahedron. The output array is (P,3), where.
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 getEquispacedLatticeTriangle(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0)
Computes an equispaced lattice of a given order on the reference triangle. The output array is (P...