Compadre  1.5.9
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Compadre_Basis.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Compadre: COMpatible PArticle Discretization and REmap Toolkit
4 //
5 // Copyright 2018 NTESS and the Compadre contributors.
6 // SPDX-License-Identifier: BSD-2-Clause
7 // *****************************************************************************
8 // @HEADER
9 #ifndef _COMPADRE_BASIS_HPP_
10 #define _COMPADRE_BASIS_HPP_
11 
12 #include "Compadre_GMLS.hpp"
13 
14 namespace Compadre {
15 
16 /*! \brief Evaluates the polynomial basis under a particular sampling function. Generally used to fill a row of P.
17  \param data [in] - GMLSBasisData struct
18  \param teamMember [in] - Kokkos::TeamPolicy member type (created by parallel_for)
19  \param delta [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _basis_multipler*the dimension of the polynomial basis.
20  \param thread_workspace [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _poly_order*the spatial dimension of the polynomial basis.
21  \param target_index [in] - target number
22  \param neighbor_index [in] - index of neighbor for this target with respect to local numbering [0,...,number of neighbors for target]
23  \param alpha [in] - double to determine convex combination of target and neighbor site at which to evaluate polynomials. (1-alpha)*neighbor + alpha*target
24  \param dimension [in] - spatial dimension of basis to evaluate. e.g. dimension two basis of order one is 1, x, y, whereas for dimension 3 it is 1, x, y, z
25  \param poly_order [in] - polynomial basis degree
26  \param specific_order_only [in] - boolean for only evaluating one degree of polynomial when true
27  \param V [in] - orthonormal basis matrix size _dimensions * _dimensions whose first _dimensions-1 columns are an approximation of the tangent plane
28  \param reconstruction_space [in] - space of polynomial that a sampling functional is to evaluate
29  \param sampling_strategy [in] - sampling functional specification
30  \param evaluation_site_local_index [in] - local index for evaluation sites (0 is target site)
31 */
32 template <typename BasisData>
33 KOKKOS_INLINE_FUNCTION
34 void calcPij(const BasisData& data, const member_type& teamMember, double* delta, double* thread_workspace, const int target_index, int neighbor_index, const double alpha, const int dimension, const int poly_order, bool specific_order_only = false, const scratch_matrix_right_type* V = NULL, const ReconstructionSpace reconstruction_space = ReconstructionSpace::ScalarTaylorPolynomial, const SamplingFunctional polynomial_sampling_functional = PointSample, const int evaluation_site_local_index = 0) {
35 /*
36  * This class is under two levels of hierarchical parallelism, so we
37  * do not put in any finer grain parallelism in this function
38  */
39  const int my_num_neighbors = data._pc._nla.getNumberOfNeighborsDevice(target_index);
40 
41  // store precalculated factorials for speedup
42  const double factorial[15] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600, 6227020800, 87178291200};
43 
44  int component = 0;
45  if (neighbor_index >= my_num_neighbors) {
46  component = neighbor_index / my_num_neighbors;
47  neighbor_index = neighbor_index % my_num_neighbors;
48  } else if (neighbor_index < 0) {
49  // -1 maps to 0 component
50  // -2 maps to 1 component
51  // -3 maps to 2 component
52  component = -(neighbor_index+1);
53  }
54 
55  XYZ relative_coord;
56  if (neighbor_index > -1) {
57  // Evaluate at neighbor site
58  for (int i=0; i<dimension; ++i) {
59  // calculates (alpha*target+(1-alpha)*neighbor)-1*target = (alpha-1)*target + (1-alpha)*neighbor
60  relative_coord[i] = (alpha-1)*data._pc.getTargetCoordinate(target_index, i, V);
61  relative_coord[i] += (1-alpha)*data._pc.getNeighborCoordinate(target_index, neighbor_index, i, V);
62  }
63  } else if (evaluation_site_local_index > 0) {
64  // Extra evaluation site
65  for (int i=0; i<dimension; ++i) {
66  // evaluation_site_local_index is for evaluation site, which includes target site
67  // the -1 converts to the local index for ADDITIONAL evaluation sites
68  relative_coord[i] = data._additional_pc.getNeighborCoordinate(target_index, evaluation_site_local_index-1, i, V);
69  relative_coord[i] -= data._pc.getTargetCoordinate(target_index, i, V);
70  }
71  } else {
72  // Evaluate at the target site
73  for (int i=0; i<3; ++i) relative_coord[i] = 0;
74  }
75 
76  // basis ActualReconstructionSpaceRank is 0 (evaluated like a scalar) and sampling functional is traditional
77  if ((polynomial_sampling_functional == PointSample ||
78  polynomial_sampling_functional == VectorPointSample ||
79  polynomial_sampling_functional == ManifoldVectorPointSample ||
80  polynomial_sampling_functional == VaryingManifoldVectorPointSample)&&
81  (reconstruction_space == ScalarTaylorPolynomial || reconstruction_space == VectorOfScalarClonesTaylorPolynomial)) {
82 
83  double cutoff_p = data._epsilons(target_index);
84  const int start_index = specific_order_only ? poly_order : 0; // only compute specified order if requested
85 
86  ScalarTaylorPolynomialBasis::evaluate(teamMember, delta, thread_workspace, dimension, poly_order, cutoff_p, relative_coord.x, relative_coord.y, relative_coord.z, start_index);
87 
88  // basis ActualReconstructionSpaceRank is 1 (is a true vector basis) and sampling functional is traditional
89  } else if ((polynomial_sampling_functional == VectorPointSample ||
90  polynomial_sampling_functional == ManifoldVectorPointSample ||
91  polynomial_sampling_functional == VaryingManifoldVectorPointSample) &&
92  (reconstruction_space == VectorTaylorPolynomial)) {
93 
94  const int dimension_offset = GMLS::getNP(data._poly_order, dimension, reconstruction_space);
95  double cutoff_p = data._epsilons(target_index);
96  const int start_index = specific_order_only ? poly_order : 0; // only compute specified order if requested
97 
98  for (int d=0; d<dimension; ++d) {
99  if (d==component) {
100  ScalarTaylorPolynomialBasis::evaluate(teamMember, delta+component*dimension_offset, thread_workspace, dimension, poly_order, cutoff_p, relative_coord.x, relative_coord.y, relative_coord.z, start_index);
101  } else {
102  for (int n=0; n<dimension_offset; ++n) {
103  *(delta+d*dimension_offset+n) = 0;
104  }
105  }
106  }
107  } else if ((polynomial_sampling_functional == VectorPointSample) &&
108  (reconstruction_space == DivergenceFreeVectorTaylorPolynomial)) {
109  // Divergence free vector polynomial basis
110  double cutoff_p = data._epsilons(target_index);
111 
112  DivergenceFreePolynomialBasis::evaluate(teamMember, delta, thread_workspace, dimension, poly_order, component, cutoff_p, relative_coord.x, relative_coord.y, relative_coord.z);
113  } else if (reconstruction_space == BernsteinPolynomial) {
114  // Bernstein vector polynomial basis
115  double cutoff_p = data._epsilons(target_index);
116 
117  BernsteinPolynomialBasis::evaluate(teamMember, delta, thread_workspace, dimension, poly_order, component, cutoff_p, relative_coord.x, relative_coord.y, relative_coord.z);
118 
119  } else if ((polynomial_sampling_functional == StaggeredEdgeAnalyticGradientIntegralSample) &&
120  (reconstruction_space == ScalarTaylorPolynomial)) {
121  double cutoff_p = data._epsilons(target_index);
122  const int start_index = specific_order_only ? poly_order : 0; // only compute specified order if requested
123  // basis is actually scalar with staggered sampling functional
124  ScalarTaylorPolynomialBasis::evaluate(teamMember, delta, thread_workspace, dimension, poly_order, cutoff_p, relative_coord.x, relative_coord.y, relative_coord.z, start_index, 0.0, -1.0);
125  relative_coord.x = 0;
126  relative_coord.y = 0;
127  relative_coord.z = 0;
128  ScalarTaylorPolynomialBasis::evaluate(teamMember, delta, thread_workspace, dimension, poly_order, cutoff_p, relative_coord.x, relative_coord.y, relative_coord.z, start_index, 1.0, 1.0);
129  } else if (polynomial_sampling_functional == StaggeredEdgeIntegralSample) {
130  Kokkos::single(Kokkos::PerThread(teamMember), [&] () {
131  if (data._problem_type == ProblemType::MANIFOLD) {
132  double cutoff_p = data._epsilons(target_index);
133  int alphax, alphay;
134  double alphaf;
135  const int start_index = specific_order_only ? poly_order : 0; // only compute specified order if requested
136 
137  for (int quadrature = 0; quadrature<data._qm.getNumberOfQuadraturePoints(); ++quadrature) {
138 
139  int i = 0;
140 
141  XYZ tangent_quadrature_coord_2d;
142  XYZ quadrature_coord_2d;
143  for (int j=0; j<dimension; ++j) {
144  // calculates (alpha*target+(1-alpha)*neighbor)-1*target = (alpha-1)*target + (1-alpha)*neighbor
145  quadrature_coord_2d[j] = (data._qm.getSite(quadrature,0)-1)*data._pc.getTargetCoordinate(target_index, j, V);
146  quadrature_coord_2d[j] += (1-data._qm.getSite(quadrature,0))*data._pc.getNeighborCoordinate(target_index, neighbor_index, j, V);
147  tangent_quadrature_coord_2d[j] = data._pc.getTargetCoordinate(target_index, j, V);
148  tangent_quadrature_coord_2d[j] -= data._pc.getNeighborCoordinate(target_index, neighbor_index, j, V);
149  }
150  for (int j=0; j<data._basis_multiplier; ++j) {
151  for (int n = start_index; n <= poly_order; n++){
152  for (alphay = 0; alphay <= n; alphay++){
153  alphax = n - alphay;
154  alphaf = factorial[alphax]*factorial[alphay];
155 
156  // local evaluation of vector [0,p] or [p,0]
157  double v0, v1;
158  v0 = (j==0) ? std::pow(quadrature_coord_2d.x/cutoff_p,alphax)
159  *std::pow(quadrature_coord_2d.y/cutoff_p,alphay)/alphaf : 0;
160  v1 = (j==0) ? 0 : std::pow(quadrature_coord_2d.x/cutoff_p,alphax)
161  *std::pow(quadrature_coord_2d.y/cutoff_p,alphay)/alphaf;
162 
163  double dot_product = tangent_quadrature_coord_2d[0]*v0 + tangent_quadrature_coord_2d[1]*v1;
164 
165  // multiply by quadrature weight
166  if (quadrature==0) {
167  *(delta+i) = dot_product * data._qm.getWeight(quadrature);
168  } else {
169  *(delta+i) += dot_product * data._qm.getWeight(quadrature);
170  }
171  i++;
172  }
173  }
174  }
175  }
176  } else {
177  // Calculate basis matrix for NON MANIFOLD problems
178  double cutoff_p = data._epsilons(target_index);
179  int alphax, alphay, alphaz;
180  double alphaf;
181  const int start_index = specific_order_only ? poly_order : 0; // only compute specified order if requested
182 
183  for (int quadrature = 0; quadrature<data._qm.getNumberOfQuadraturePoints(); ++quadrature) {
184 
185  int i = 0;
186 
187  XYZ quadrature_coord_3d;
188  XYZ tangent_quadrature_coord_3d;
189  for (int j=0; j<dimension; ++j) {
190  // calculates (alpha*target+(1-alpha)*neighbor)-1*target = (alpha-1)*target + (1-alpha)*neighbor
191  quadrature_coord_3d[j] = (data._qm.getSite(quadrature,0)-1)*data._pc.getTargetCoordinate(target_index, j, NULL);
192  quadrature_coord_3d[j] += (1-data._qm.getSite(quadrature,0))*data._pc.getNeighborCoordinate(target_index, neighbor_index, j, NULL);
193  tangent_quadrature_coord_3d[j] = data._pc.getTargetCoordinate(target_index, j, NULL);
194  tangent_quadrature_coord_3d[j] -= data._pc.getNeighborCoordinate(target_index, neighbor_index, j, NULL);
195  }
196  for (int j=0; j<data._basis_multiplier; ++j) {
197  for (int n = start_index; n <= poly_order; n++) {
198  if (dimension == 3) {
199  for (alphaz = 0; alphaz <= n; alphaz++){
200  int s = n - alphaz;
201  for (alphay = 0; alphay <= s; alphay++){
202  alphax = s - alphay;
203  alphaf = factorial[alphax]*factorial[alphay]*factorial[alphaz];
204 
205  // local evaluation of vector [p, 0, 0], [0, p, 0] or [0, 0, p]
206  double v0, v1, v2;
207  switch(j) {
208  case 1:
209  v0 = 0.0;
210  v1 = std::pow(quadrature_coord_3d.x/cutoff_p,alphax)
211  *std::pow(quadrature_coord_3d.y/cutoff_p,alphay)
212  *std::pow(quadrature_coord_3d.z/cutoff_p,alphaz)/alphaf;
213  v2 = 0.0;
214  break;
215  case 2:
216  v0 = 0.0;
217  v1 = 0.0;
218  v2 = std::pow(quadrature_coord_3d.x/cutoff_p,alphax)
219  *std::pow(quadrature_coord_3d.y/cutoff_p,alphay)
220  *std::pow(quadrature_coord_3d.z/cutoff_p,alphaz)/alphaf;
221  break;
222  default:
223  v0 = std::pow(quadrature_coord_3d.x/cutoff_p,alphax)
224  *std::pow(quadrature_coord_3d.y/cutoff_p,alphay)
225  *std::pow(quadrature_coord_3d.z/cutoff_p,alphaz)/alphaf;
226  v1 = 0.0;
227  v2 = 0.0;
228  break;
229  }
230 
231  double dot_product = tangent_quadrature_coord_3d[0]*v0 + tangent_quadrature_coord_3d[1]*v1 + tangent_quadrature_coord_3d[2]*v2;
232 
233  // multiply by quadrature weight
234  if (quadrature == 0) {
235  *(delta+i) = dot_product * data._qm.getWeight(quadrature);
236  } else {
237  *(delta+i) += dot_product * data._qm.getWeight(quadrature);
238  }
239  i++;
240  }
241  }
242  } else if (dimension == 2) {
243  for (alphay = 0; alphay <= n; alphay++){
244  alphax = n - alphay;
245  alphaf = factorial[alphax]*factorial[alphay];
246 
247  // local evaluation of vector [p, 0] or [0, p]
248  double v0, v1;
249  switch(j) {
250  case 1:
251  v0 = 0.0;
252  v1 = std::pow(quadrature_coord_3d.x/cutoff_p,alphax)
253  *std::pow(quadrature_coord_3d.y/cutoff_p,alphay)/alphaf;
254  break;
255  default:
256  v0 = std::pow(quadrature_coord_3d.x/cutoff_p,alphax)
257  *std::pow(quadrature_coord_3d.y/cutoff_p,alphay)/alphaf;
258  v1 = 0.0;
259  break;
260  }
261 
262  double dot_product = tangent_quadrature_coord_3d[0]*v0 + tangent_quadrature_coord_3d[1]*v1;
263 
264  // multiply by quadrature weight
265  if (quadrature == 0) {
266  *(delta+i) = dot_product * data._qm.getWeight(quadrature);
267  } else {
268  *(delta+i) += dot_product * data._qm.getWeight(quadrature);
269  }
270  i++;
271  }
272  }
273  }
274  }
275  }
276  } // NON MANIFOLD PROBLEMS
277  });
278  } else if ((polynomial_sampling_functional == FaceNormalIntegralSample ||
279  polynomial_sampling_functional == EdgeTangentIntegralSample ||
280  polynomial_sampling_functional == FaceNormalAverageSample ||
281  polynomial_sampling_functional == EdgeTangentAverageSample) &&
282  reconstruction_space == VectorTaylorPolynomial) {
283 
284  compadre_kernel_assert_debug(data._local_dimensions==2 &&
285  "FaceNormalIntegralSample, EdgeTangentIntegralSample, FaceNormalAverageSample, \
286  and EdgeTangentAverageSample only support 2d or 3d with 2d manifold");
287 
288  compadre_kernel_assert_debug(data._qm.getDimensionOfQuadraturePoints()==1 \
289  && "Only 1d quadrature may be specified for edge integrals");
290 
291  compadre_kernel_assert_debug(data._qm.getNumberOfQuadraturePoints()>=1 \
292  && "Quadrature points not generated");
293 
294  compadre_kernel_assert_debug(data._source_extra_data.extent(0)>0 && "Extra data used but not set.");
295 
296  compadre_kernel_assert_debug(!specific_order_only &&
297  "Sampling functional does not support specific_order_only");
298 
299  double cutoff_p = data._epsilons(target_index);
300  auto global_neighbor_index = data._pc.getNeighborIndex(target_index, neighbor_index);
301  int alphax, alphay;
302  double alphaf;
303 
304  /*
305  * requires quadrature points defined on an edge, not a target/source edge (spoke)
306  *
307  * data._source_extra_data will contain the endpoints (2 for 2D, 3 for 3D) and then the unit normals
308  * (e0_x, e0_y, e1_x, e1_y, n_x, n_y, t_x, t_y)
309  */
310 
311  int quadrature_point_loop = data._qm.getNumberOfQuadraturePoints();
312 
313  const int TWO = 2; // used because of # of vertices on an edge
314  double G_data[3*TWO]; // max(2,3)*TWO
315  double edge_coords[3*TWO];
316  for (int i=0; i<data._global_dimensions*TWO; ++i) G_data[i] = 0;
317  for (int i=0; i<data._global_dimensions*TWO; ++i) edge_coords[i] = 0;
318  // 2 is for # vertices on an edge
319  scratch_matrix_right_type G(G_data, data._global_dimensions, TWO);
320  scratch_matrix_right_type edge_coords_matrix(edge_coords, data._global_dimensions, TWO);
321 
322  // neighbor coordinate is assumed to be midpoint
323  // could be calculated, but is correct for sphere
324  // and also for non-manifold problems
325  // uses given midpoint, rather than computing the midpoint from vertices
326  double radius = 0.0;
327  // this midpoint should lie on the sphere, or this will be the wrong radius
328  for (int j=0; j<data._global_dimensions; ++j) {
329  edge_coords_matrix(j, 0) = data._source_extra_data(global_neighbor_index, j);
330  edge_coords_matrix(j, 1) = data._source_extra_data(global_neighbor_index, data._global_dimensions + j) - edge_coords_matrix(j, 0);
331  radius += edge_coords_matrix(j, 0)*edge_coords_matrix(j, 0);
332  }
333  radius = std::sqrt(radius);
334 
335  // edge_coords now has:
336  // (v0_x, v0_y, v1_x-v0_x, v1_y-v0_y)
337  auto E = edge_coords_matrix;
338 
339  // get arc length of edge on manifold
340  double theta = 0.0;
341  if (data._problem_type == ProblemType::MANIFOLD) {
342  XYZ a = {E(0,0), E(1,0), E(2,0)};
343  XYZ b = {E(0,1)+E(0,0), E(1,1)+E(1,0), E(2,1)+E(2,0)};
344  double a_dot_b = a[0]*b[0] + a[1]*b[1] + a[2]*b[2];
345  double norm_a_cross_b = getAreaFromVectors(teamMember, a, b);
346  theta = std::atan(norm_a_cross_b / a_dot_b);
347  }
348 
349  // loop
350  double entire_edge_length = 0.0;
351  for (int quadrature = 0; quadrature<quadrature_point_loop; ++quadrature) {
352 
353  double G_determinant = 1.0;
354  double scaled_transformed_qp[3] = {0,0,0};
355  double unscaled_transformed_qp[3] = {0,0,0};
356  for (int j=0; j<data._global_dimensions; ++j) {
357  unscaled_transformed_qp[j] += E(j,1)*data._qm.getSite(quadrature, 0);
358  // adds back on shift by endpoint
359  unscaled_transformed_qp[j] += E(j,0);
360  }
361 
362  // project onto the sphere
363  if (data._problem_type == ProblemType::MANIFOLD) {
364  // unscaled_transformed_qp now lives on cell edge, but if on manifold,
365  // not directly on the sphere, just near by
366 
367  // normalize to project back onto sphere
368  double transformed_qp_norm = 0;
369  for (int j=0; j<data._global_dimensions; ++j) {
370  transformed_qp_norm += unscaled_transformed_qp[j]*unscaled_transformed_qp[j];
371  }
372  transformed_qp_norm = std::sqrt(transformed_qp_norm);
373  // transformed_qp made unit length
374  for (int j=0; j<data._global_dimensions; ++j) {
375  scaled_transformed_qp[j] = unscaled_transformed_qp[j] * radius / transformed_qp_norm;
376  }
377 
378  G_determinant = radius * theta;
379  XYZ qp = XYZ(scaled_transformed_qp[0], scaled_transformed_qp[1], scaled_transformed_qp[2]);
380  for (int j=0; j<data._local_dimensions; ++j) {
381  relative_coord[j] = data._pc.convertGlobalToLocalCoordinate(qp,j,*V)
382  - data._pc.getTargetCoordinate(target_index,j,V);
383  // shift quadrature point by target site
384  }
385  relative_coord[2] = 0;
386  } else { // not on a manifold, but still integrated
387  XYZ endpoints_difference = {E(0,1), E(1,1), 0};
388  G_determinant = data._pc.EuclideanVectorLength(endpoints_difference, 2);
389  for (int j=0; j<data._local_dimensions; ++j) {
390  relative_coord[j] = unscaled_transformed_qp[j]
391  - data._pc.getTargetCoordinate(target_index,j,V);
392  // shift quadrature point by target site
393  }
394  relative_coord[2] = 0;
395  }
396 
397  // get normal or tangent direction (ambient)
398  XYZ direction;
399  if (polynomial_sampling_functional == FaceNormalIntegralSample
400  || polynomial_sampling_functional == FaceNormalAverageSample) {
401  for (int j=0; j<data._global_dimensions; ++j) {
402  // normal direction
403  direction[j] = data._source_extra_data(global_neighbor_index, 2*data._global_dimensions + j);
404  }
405  } else {
406  if (data._problem_type == ProblemType::MANIFOLD) {
407  // generate tangent from outward normal direction of the sphere and edge normal
408  XYZ k = {scaled_transformed_qp[0], scaled_transformed_qp[1], scaled_transformed_qp[2]};
409  double k_norm = std::sqrt(k[0]*k[0]+k[1]*k[1]+k[2]*k[2]);
410  k[0] = k[0]/k_norm; k[1] = k[1]/k_norm; k[2] = k[2]/k_norm;
411  XYZ n = {data._source_extra_data(global_neighbor_index, 2*data._global_dimensions + 0),
412  data._source_extra_data(global_neighbor_index, 2*data._global_dimensions + 1),
413  data._source_extra_data(global_neighbor_index, 2*data._global_dimensions + 2)};
414 
415  double norm_k_cross_n = getAreaFromVectors(teamMember, k, n);
416  direction[0] = (k[1]*n[2] - k[2]*n[1]) / norm_k_cross_n;
417  direction[1] = (k[2]*n[0] - k[0]*n[2]) / norm_k_cross_n;
418  direction[2] = (k[0]*n[1] - k[1]*n[0]) / norm_k_cross_n;
419  } else {
420  for (int j=0; j<data._global_dimensions; ++j) {
421  // tangent direction
422  direction[j] = data._source_extra_data(global_neighbor_index, 3*data._global_dimensions + j);
423  }
424  }
425  }
426 
427  // convert direction to local chart (for manifolds)
428  XYZ local_direction;
429  if (data._problem_type == ProblemType::MANIFOLD) {
430  for (int j=0; j<data._basis_multiplier; ++j) {
431  // Project ambient normal direction onto local chart basis as a local direction.
432  // Using V alone to provide vectors only gives tangent vectors at
433  // the target site. This could result in accuracy < 3rd order.
434  local_direction[j] = data._pc.convertGlobalToLocalCoordinate(direction,j,*V);
435  }
436  }
437 
438  int i = 0;
439  for (int j=0; j<data._basis_multiplier; ++j) {
440  for (int n = 0; n <= poly_order; n++){
441  for (alphay = 0; alphay <= n; alphay++){
442  alphax = n - alphay;
443  alphaf = factorial[alphax]*factorial[alphay];
444 
445  // local evaluation of vector [0,p] or [p,0]
446  double v0, v1;
447  v0 = (j==0) ? std::pow(relative_coord.x/cutoff_p,alphax)
448  *std::pow(relative_coord.y/cutoff_p,alphay)/alphaf : 0;
449  v1 = (j==1) ? std::pow(relative_coord.x/cutoff_p,alphax)
450  *std::pow(relative_coord.y/cutoff_p,alphay)/alphaf : 0;
451 
452  // either n*v or t*v
453  double dot_product = 0.0;
454  if (data._problem_type == ProblemType::MANIFOLD) {
455  // alternate option for projection
456  dot_product = local_direction[0]*v0 + local_direction[1]*v1;
457  } else {
458  dot_product = direction[0]*v0 + direction[1]*v1;
459  }
460 
461  // multiply by quadrature weight
462  if (quadrature==0) {
463  *(delta+i) = dot_product * data._qm.getWeight(quadrature) * G_determinant;
464  } else {
465  *(delta+i) += dot_product * data._qm.getWeight(quadrature) * G_determinant;
466  }
467  i++;
468  }
469  }
470  }
471  entire_edge_length += G_determinant * data._qm.getWeight(quadrature);
472  }
473  if (polynomial_sampling_functional == FaceNormalAverageSample ||
474  polynomial_sampling_functional == EdgeTangentAverageSample) {
475  int k = 0;
476  for (int j=0; j<data._basis_multiplier; ++j) {
477  for (int n = 0; n <= poly_order; n++){
478  for (alphay = 0; alphay <= n; alphay++){
479  *(delta+k) /= entire_edge_length;
480  k++;
481  }
482  }
483  }
484  }
485  } else if (polynomial_sampling_functional == CellAverageSample ||
486  polynomial_sampling_functional == CellIntegralSample) {
487 
488  compadre_kernel_assert_debug(data._local_dimensions==2 &&
489  "CellAverageSample only supports 2d or 3d with 2d manifold");
490  auto global_neighbor_index = data._pc.getNeighborIndex(target_index, neighbor_index);
491  double cutoff_p = data._epsilons(target_index);
492  int alphax, alphay;
493  double alphaf;
494 
495  double G_data[3*3]; //data._global_dimensions*3
496  double triangle_coords[3*3]; //data._global_dimensions*3
497  for (int i=0; i<data._global_dimensions*3; ++i) G_data[i] = 0;
498  for (int i=0; i<data._global_dimensions*3; ++i) triangle_coords[i] = 0;
499  // 3 is for # vertices in sub-triangle
500  scratch_matrix_right_type G(G_data, data._global_dimensions, 3);
501  scratch_matrix_right_type triangle_coords_matrix(triangle_coords, data._global_dimensions, 3);
502 
503  // neighbor coordinate is assumed to be midpoint
504  // could be calculated, but is correct for sphere
505  // and also for non-manifold problems
506  // uses given midpoint, rather than computing the midpoint from vertices
507  double radius = 0.0;
508  // this midpoint should lie on the sphere, or this will be the wrong radius
509  for (int j=0; j<data._global_dimensions; ++j) {
510  // midpoint
511  triangle_coords_matrix(j, 0) = data._pc.getNeighborCoordinate(target_index, neighbor_index, j);
512  radius += triangle_coords_matrix(j, 0)*triangle_coords_matrix(j, 0);
513  }
514  radius = std::sqrt(radius);
515 
516  // NaN in entry (data._global_dimensions) is a convention for indicating fewer vertices
517  // for this cell and NaN is checked by entry!=entry
518  int num_vertices = 0;
519  for (int j=0; j<data._source_extra_data.extent_int(1); ++j) {
520  auto val = data._source_extra_data(global_neighbor_index, j);
521  if (val != val) {
522  break;
523  } else {
524  num_vertices++;
525  }
526  }
527  num_vertices = num_vertices / data._global_dimensions;
528  auto T = triangle_coords_matrix;
529 
530  // loop over each two vertices
531  // made for flat surfaces (either dim=2 or on a manifold)
532  double entire_cell_area = 0.0;
533  for (int v=0; v<num_vertices; ++v) {
534  int v1 = v;
535  int v2 = (v+1) % num_vertices;
536 
537  for (int j=0; j<data._global_dimensions; ++j) {
538  triangle_coords_matrix(j,1) = data._source_extra_data(global_neighbor_index,
539  v1*data._global_dimensions+j)
540  - triangle_coords_matrix(j,0);
541  triangle_coords_matrix(j,2) = data._source_extra_data(global_neighbor_index,
542  v2*data._global_dimensions+j)
543  - triangle_coords_matrix(j,0);
544  }
545 
546  // triangle_coords now has:
547  // (midpoint_x, midpoint_y, midpoint_z,
548  // v1_x-midpoint_x, v1_y-midpoint_y, v1_z-midpoint_z,
549  // v2_x-midpoint_x, v2_y-midpoint_y, v2_z-midpoint_z);
550  for (int quadrature = 0; quadrature<data._qm.getNumberOfQuadraturePoints(); ++quadrature) {
551  double unscaled_transformed_qp[3] = {0,0,0};
552  double scaled_transformed_qp[3] = {0,0,0};
553  for (int j=0; j<data._global_dimensions; ++j) {
554  for (int k=1; k<3; ++k) { // 3 is for # vertices in subtriangle
555  // uses vertex-midpoint as one direction
556  // and other vertex-midpoint as other direction
557  unscaled_transformed_qp[j] += T(j,k)*data._qm.getSite(quadrature, k-1);
558  }
559  // adds back on shift by midpoint
560  unscaled_transformed_qp[j] += T(j,0);
561  }
562 
563  // project onto the sphere
564  double G_determinant = 1.0;
565  if (data._problem_type == ProblemType::MANIFOLD) {
566  // unscaled_transformed_qp now lives on cell, but if on manifold,
567  // not directly on the sphere, just near by
568 
569  // normalize to project back onto sphere
570  double transformed_qp_norm = 0;
571  for (int j=0; j<data._global_dimensions; ++j) {
572  transformed_qp_norm += unscaled_transformed_qp[j]*unscaled_transformed_qp[j];
573  }
574  transformed_qp_norm = std::sqrt(transformed_qp_norm);
575  // transformed_qp made unit length
576  for (int j=0; j<data._global_dimensions; ++j) {
577  scaled_transformed_qp[j] = unscaled_transformed_qp[j] * radius / transformed_qp_norm;
578  }
579 
580 
581  // u_qp = midpoint + r_qp[1]*(v_1-midpoint) + r_qp[2]*(v_2-midpoint)
582  // s_qp = u_qp * radius / norm(u_qp) = radius * u_qp / norm(u_qp)
583  //
584  // so G(:,i) is \partial{s_qp}/ \partial{r_qp[i]}
585  // where r_qp is reference quadrature point (R^2 in 2D manifold in R^3)
586  //
587  // G(:,i) = radius * ( \partial{u_qp}/\partial{r_qp[i]} * (\sum_m u_qp[k]^2)^{-1/2}
588  // + u_qp * \partial{(\sum_m u_qp[k]^2)^{-1/2}}/\partial{r_qp[i]} )
589  //
590  // = radius * ( T(:,i)/norm(u_qp) + u_qp*(-1/2)*(\sum_m u_qp[k]^2)^{-3/2}
591  // *2*(\sum_k u_qp[k]*\partial{u_qp[k]}/\partial{r_qp[i]}) )
592  //
593  // = radius * ( T(:,i)/norm(u_qp) + u_qp*(-1/2)*(\sum_m u_qp[k]^2)^{-3/2}
594  // *2*(\sum_k u_qp[k]*T(k,i)) )
595  //
596  // NOTE: we do not multiply G by radius before determining area from vectors,
597  // so we multiply by radius**2 after calculation
598  double qp_norm_sq = transformed_qp_norm*transformed_qp_norm;
599  for (int j=0; j<data._global_dimensions; ++j) {
600  G(j,1) = T(j,1)/transformed_qp_norm;
601  G(j,2) = T(j,2)/transformed_qp_norm;
602  for (int k=0; k<data._global_dimensions; ++k) {
603  G(j,1) += unscaled_transformed_qp[j]*(-0.5)*std::pow(qp_norm_sq,-1.5)
604  *2*(unscaled_transformed_qp[k]*T(k,1));
605  G(j,2) += unscaled_transformed_qp[j]*(-0.5)*std::pow(qp_norm_sq,-1.5)
606  *2*(unscaled_transformed_qp[k]*T(k,2));
607  }
608  }
609  G_determinant = getAreaFromVectors(teamMember,
610  Kokkos::subview(G, Kokkos::ALL(), 1), Kokkos::subview(G, Kokkos::ALL(), 2));
611  G_determinant *= radius*radius;
612  XYZ qp = XYZ(scaled_transformed_qp[0], scaled_transformed_qp[1], scaled_transformed_qp[2]);
613  for (int j=0; j<data._local_dimensions; ++j) {
614  relative_coord[j] = data._pc.convertGlobalToLocalCoordinate(qp,j,*V)
615  - data._pc.getTargetCoordinate(target_index,j,V);
616  // shift quadrature point by target site
617  }
618  relative_coord[2] = 0;
619  } else {
620  G_determinant = getAreaFromVectors(teamMember,
621  Kokkos::subview(T, Kokkos::ALL(), 1), Kokkos::subview(T, Kokkos::ALL(), 2));
622  for (int j=0; j<data._local_dimensions; ++j) {
623  relative_coord[j] = unscaled_transformed_qp[j]
624  - data._pc.getTargetCoordinate(target_index,j,V);
625  // shift quadrature point by target site
626  }
627  relative_coord[2] = 0;
628  }
629 
630  int k = 0;
631  compadre_kernel_assert_debug(!specific_order_only &&
632  "CellAverageSample does not support specific_order_only");
633  for (int n = 0; n <= poly_order; n++){
634  for (alphay = 0; alphay <= n; alphay++){
635  alphax = n - alphay;
636  alphaf = factorial[alphax]*factorial[alphay];
637  double val_to_sum = G_determinant * (data._qm.getWeight(quadrature)
638  * std::pow(relative_coord.x/cutoff_p,alphax)
639  * std::pow(relative_coord.y/cutoff_p,alphay) / alphaf);
640  if (quadrature==0 && v==0) *(delta+k) = val_to_sum;
641  else *(delta+k) += val_to_sum;
642  k++;
643  }
644  }
645  entire_cell_area += G_determinant * data._qm.getWeight(quadrature);
646  }
647  }
648  if (polynomial_sampling_functional == CellAverageSample) {
649  int k = 0;
650  for (int n = 0; n <= poly_order; n++){
651  for (alphay = 0; alphay <= n; alphay++){
652  *(delta+k) /= entire_cell_area;
653  k++;
654  }
655  }
656  }
657  } else {
658  compadre_kernel_assert_release((false) && "Sampling and basis space combination not defined.");
659  }
660 }
661 
662 /*! \brief Evaluates the gradient of a polynomial basis under the Dirac Delta (pointwise) sampling function.
663  \param data [in] - GMLSBasisData struct
664  \param teamMember [in] - Kokkos::TeamPolicy member type (created by parallel_for)
665  \param delta [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large is the _basis_multipler*the dimension of the polynomial basis.
666  \param thread_workspace [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _poly_order*the spatial dimension of the polynomial basis.
667  \param target_index [in] - target number
668  \param neighbor_index [in] - index of neighbor for this target with respect to local numbering [0,...,number of neighbors for target]
669  \param alpha [in] - double to determine convex combination of target and neighbor site at which to evaluate polynomials. (1-alpha)*neighbor + alpha*target
670  \param partial_direction [in] - direction that partial is taken with respect to, e.g. 0 is x direction, 1 is y direction
671  \param dimension [in] - spatial dimension of basis to evaluate. e.g. dimension two basis of order one is 1, x, y, whereas for dimension 3 it is 1, x, y, z
672  \param poly_order [in] - polynomial basis degree
673  \param specific_order_only [in] - boolean for only evaluating one degree of polynomial when true
674  \param V [in] - orthonormal basis matrix size _dimensions * _dimensions whose first _dimensions-1 columns are an approximation of the tangent plane
675  \param reconstruction_space [in] - space of polynomial that a sampling functional is to evaluate
676  \param sampling_strategy [in] - sampling functional specification
677  \param evaluation_site_local_index [in] - local index for evaluation sites (0 is target site)
678 */
679 template <typename BasisData>
680 KOKKOS_INLINE_FUNCTION
681 void calcGradientPij(const BasisData& data, const member_type& teamMember, double* delta, double* thread_workspace, const int target_index, int neighbor_index, const double alpha, const int partial_direction, const int dimension, const int poly_order, bool specific_order_only, const scratch_matrix_right_type* V, const ReconstructionSpace reconstruction_space, const SamplingFunctional polynomial_sampling_functional, const int evaluation_site_local_index = 0) {
682 /*
683  * This class is under two levels of hierarchical parallelism, so we
684  * do not put in any finer grain parallelism in this function
685  */
686 
687  const int my_num_neighbors = data._pc._nla.getNumberOfNeighborsDevice(target_index);
688 
689  int component = 0;
690  if (neighbor_index >= my_num_neighbors) {
691  component = neighbor_index / my_num_neighbors;
692  neighbor_index = neighbor_index % my_num_neighbors;
693  } else if (neighbor_index < 0) {
694  // -1 maps to 0 component
695  // -2 maps to 1 component
696  // -3 maps to 2 component
697  component = -(neighbor_index+1);
698  }
699 
700  // alpha corresponds to the linear combination of target_index and neighbor_index coordinates
701  // coordinate to evaluate = alpha*(target_index's coordinate) + (1-alpha)*(neighbor_index's coordinate)
702 
703  // partial_direction - 0=x, 1=y, 2=z
704  XYZ relative_coord;
705  if (neighbor_index > -1) {
706  for (int i=0; i<dimension; ++i) {
707  // calculates (alpha*target+(1-alpha)*neighbor)-1*target = (alpha-1)*target + (1-alpha)*neighbor
708  relative_coord[i] = (alpha-1)*data._pc.getTargetCoordinate(target_index, i, V);
709  relative_coord[i] += (1-alpha)*data._pc.getNeighborCoordinate(target_index, neighbor_index, i, V);
710  }
711  } else if (evaluation_site_local_index > 0) {
712  for (int i=0; i<dimension; ++i) {
713  // evaluation_site_local_index is for evaluation site, which includes target site
714  // the -1 converts to the local index for ADDITIONAL evaluation sites
715  relative_coord[i] = data._additional_pc.getNeighborCoordinate(target_index, evaluation_site_local_index-1, i, V);
716  relative_coord[i] -= data._pc.getTargetCoordinate(target_index, i, V);
717  }
718  } else {
719  for (int i=0; i<3; ++i) relative_coord[i] = 0;
720  }
721 
722  double cutoff_p = data._epsilons(target_index);
723  const int start_index = specific_order_only ? poly_order : 0; // only compute specified order if requested
724 
725  if ((polynomial_sampling_functional == PointSample ||
726  polynomial_sampling_functional == VectorPointSample ||
727  polynomial_sampling_functional == ManifoldVectorPointSample ||
728  polynomial_sampling_functional == VaryingManifoldVectorPointSample) &&
729  (reconstruction_space == ScalarTaylorPolynomial || reconstruction_space == VectorOfScalarClonesTaylorPolynomial)) {
730 
731  ScalarTaylorPolynomialBasis::evaluatePartialDerivative(teamMember, delta, thread_workspace, dimension, poly_order, partial_direction, cutoff_p, relative_coord.x, relative_coord.y, relative_coord.z, start_index);
732 
733  } else if ((polynomial_sampling_functional == VectorPointSample) &&
734  (reconstruction_space == DivergenceFreeVectorTaylorPolynomial)) {
735 
736  // Divergence free vector polynomial basis
737  DivergenceFreePolynomialBasis::evaluatePartialDerivative(teamMember, delta, thread_workspace, dimension, poly_order, component, partial_direction, cutoff_p, relative_coord.x, relative_coord.y, relative_coord.z);
738 
739  } else if (reconstruction_space == BernsteinPolynomial) {
740 
741  // Bernstein vector polynomial basis
742  BernsteinPolynomialBasis::evaluatePartialDerivative(teamMember, delta, thread_workspace, dimension, poly_order, component, partial_direction, cutoff_p, relative_coord.x, relative_coord.y, relative_coord.z);
743 
744  } else {
745  compadre_kernel_assert_release((false) && "Sampling and basis space combination not defined.");
746  }
747 }
748 
749 /*! \brief Evaluates the Hessian of a polynomial basis under the Dirac Delta (pointwise) sampling function.
750  \param data [in] - GMLSBasisData struct
751  \param teamMember [in] - Kokkos::TeamPolicy member type (created by parallel_for)
752  \param delta [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large is the _basis_multipler*the dimension of the polynomial basis.
753  \param thread_workspace [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _poly_order*the spatial dimension of the polynomial basis.
754  \param target_index [in] - target number
755  \param neighbor_index [in] - index of neighbor for this target with respect to local numbering [0,...,number of neighbors for target]
756  \param alpha [in] - double to determine convex combination of target and neighbor site at which to evaluate polynomials. (1-alpha)*neighbor + alpha*target
757  \param partial_direction_1 [in] - first direction that partial is taken with respect to, e.g. 0 is x direction, 1 is y direction
758  \param partial_direction_2 [in] - second direction that partial is taken with respect to, e.g. 0 is x direction, 1 is y direction
759  \param dimension [in] - spatial dimension of basis to evaluate. e.g. dimension two basis of order one is 1, x, y, whereas for dimension 3 it is 1, x, y, z
760  \param poly_order [in] - polynomial basis degree
761  \param specific_order_only [in] - boolean for only evaluating one degree of polynomial when true
762  \param V [in] - orthonormal basis matrix size _dimensions * _dimensions whose first _dimensions-1 columns are an approximation of the tangent plane
763  \param reconstruction_space [in] - space of polynomial that a sampling functional is to evaluate
764  \param sampling_strategy [in] - sampling functional specification
765  \param evaluation_site_local_index [in] - local index for evaluation sites (0 is target site)
766 */
767 template <typename BasisData>
768 KOKKOS_INLINE_FUNCTION
769 void calcHessianPij(const BasisData& data, const member_type& teamMember, double* delta, double* thread_workspace, const int target_index, int neighbor_index, const double alpha, const int partial_direction_1, const int partial_direction_2, const int dimension, const int poly_order, bool specific_order_only, const scratch_matrix_right_type* V, const ReconstructionSpace reconstruction_space, const SamplingFunctional polynomial_sampling_functional, const int evaluation_site_local_index = 0) {
770 /*
771  * This class is under two levels of hierarchical parallelism, so we
772  * do not put in any finer grain parallelism in this function
773  */
774 
775  const int my_num_neighbors = data._pc._nla.getNumberOfNeighborsDevice(target_index);
776 
777  int component = 0;
778  if (neighbor_index >= my_num_neighbors) {
779  component = neighbor_index / my_num_neighbors;
780  neighbor_index = neighbor_index % my_num_neighbors;
781  } else if (neighbor_index < 0) {
782  // -1 maps to 0 component
783  // -2 maps to 1 component
784  // -3 maps to 2 component
785  component = -(neighbor_index+1);
786  }
787 
788  // alpha corresponds to the linear combination of target_index and neighbor_index coordinates
789  // coordinate to evaluate = alpha*(target_index's coordinate) + (1-alpha)*(neighbor_index's coordinate)
790 
791  // partial_direction - 0=x, 1=y, 2=z
792  XYZ relative_coord;
793  if (neighbor_index > -1) {
794  for (int i=0; i<dimension; ++i) {
795  // calculates (alpha*target+(1-alpha)*neighbor)-1*target = (alpha-1)*target + (1-alpha)*neighbor
796  relative_coord[i] = (alpha-1)*data._pc.getTargetCoordinate(target_index, i, V);
797  relative_coord[i] += (1-alpha)*data._pc.getNeighborCoordinate(target_index, neighbor_index, i, V);
798  }
799  } else if (evaluation_site_local_index > 0) {
800  for (int i=0; i<dimension; ++i) {
801  // evaluation_site_local_index is for evaluation site, which includes target site
802  // the -1 converts to the local index for ADDITIONAL evaluation sites
803  relative_coord[i] = data._additional_pc.getNeighborCoordinate(target_index, evaluation_site_local_index-1, i, V);
804  relative_coord[i] -= data._pc.getTargetCoordinate(target_index, i, V);
805  }
806  } else {
807  for (int i=0; i<3; ++i) relative_coord[i] = 0;
808  }
809 
810  double cutoff_p = data._epsilons(target_index);
811 
812  if ((polynomial_sampling_functional == PointSample ||
813  polynomial_sampling_functional == VectorPointSample ||
814  polynomial_sampling_functional == ManifoldVectorPointSample ||
815  polynomial_sampling_functional == VaryingManifoldVectorPointSample) &&
816  (reconstruction_space == ScalarTaylorPolynomial || reconstruction_space == VectorOfScalarClonesTaylorPolynomial)) {
817 
818  ScalarTaylorPolynomialBasis::evaluateSecondPartialDerivative(teamMember, delta, thread_workspace, dimension, poly_order, partial_direction_1, partial_direction_2, cutoff_p, relative_coord.x, relative_coord.y, relative_coord.z);
819 
820  } else if ((polynomial_sampling_functional == VectorPointSample) &&
821  (reconstruction_space == DivergenceFreeVectorTaylorPolynomial)) {
822 
823  DivergenceFreePolynomialBasis::evaluateSecondPartialDerivative(teamMember, delta, thread_workspace, dimension, poly_order, component, partial_direction_1, partial_direction_2, cutoff_p, relative_coord.x, relative_coord.y, relative_coord.z);
824 
825  } else if (reconstruction_space == BernsteinPolynomial) {
826 
827  BernsteinPolynomialBasis::evaluateSecondPartialDerivative(teamMember, delta, thread_workspace, dimension, poly_order, component, partial_direction_1, partial_direction_2, cutoff_p, relative_coord.x, relative_coord.y, relative_coord.z);
828 
829  } else {
830  compadre_kernel_assert_release((false) && "Sampling and basis space combination not defined.");
831  }
832 }
833 
834 /*! \brief Fills the _P matrix with either P or P*sqrt(w)
835  \param data [in] - GMLSBasisData struct
836  \param teamMember [in] - Kokkos::TeamPolicy member type (created by parallel_for)
837  \param delta [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large is the _basis_multipler*the dimension of the polynomial basis.
838  \param thread_workspace [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _poly_order*the spatial dimension of the polynomial basis.
839  \param P [out] - 2D Kokkos View which will contain evaluation of sampling functional on polynomial basis for each neighbor the target has (stored column major)
840  \param w [out] - 1D Kokkos View which will contain weighting kernel values for the target with each neighbor if weight_p = true
841  \param dimension [in] - spatial dimension of basis to evaluate. e.g. dimension two basis of order one is 1, x, y, whereas for dimension 3 it is 1, x, y, z
842  \param polynomial_order [in] - polynomial basis degree
843  \param weight_p [in] - boolean whether to fill w with kernel weights
844  \param V [in] - orthonormal basis matrix size _dimensions * _dimensions whose first _dimensions-1 columns are an approximation of the tangent plane
845  \param reconstruction_space [in] - space of polynomial that a sampling functional is to evaluate
846  \param sampling_strategy [in] - sampling functional specification
847 */
848 template <typename BasisData>
849 KOKKOS_INLINE_FUNCTION
850 void createWeightsAndP(const BasisData& data, const member_type& teamMember, scratch_vector_type delta, scratch_vector_type thread_workspace, scratch_matrix_right_type P, scratch_vector_type w, const int dimension, int polynomial_order, bool weight_p = false, scratch_matrix_right_type* V = NULL, const ReconstructionSpace reconstruction_space = ReconstructionSpace::ScalarTaylorPolynomial, const SamplingFunctional polynomial_sampling_functional = PointSample) {
851  /*
852  * Creates sqrt(W)*P
853  */
854  const int target_index = data._initial_index_for_batch + teamMember.league_rank();
855 // printf("specific order: %d\n", specific_order);
856 // {
857 // const int storage_size = (specific_order > 0) ? GMLS::getNP(specific_order, dimension)-GMLS::getNP(specific_order-1, dimension) : GMLS::getNP(data._poly_order, dimension);
858 // printf("storage size: %d\n", storage_size);
859 // }
860 // printf("weight_p: %d\n", weight_p);
861 
862  // not const b.c. of gcc 7.2 issue
863  int my_num_neighbors = data._pc._nla.getNumberOfNeighborsDevice(target_index);
864 
865  // storage_size needs to change based on the size of the basis
866  int storage_size = GMLS::getNP(polynomial_order, dimension, reconstruction_space);
867  storage_size *= data._basis_multiplier;
868  for (int j = 0; j < delta.extent_int(0); ++j) {
869  delta(j) = 0;
870  }
871  for (int j = 0; j < thread_workspace.extent_int(0); ++j) {
872  thread_workspace(j) = 0;
873  }
874  Kokkos::parallel_for(Kokkos::TeamThreadRange(teamMember,data._pc._nla.getNumberOfNeighborsDevice(target_index)),
875  [&] (const int i) {
876 
877  for (int d=0; d<data._sampling_multiplier; ++d) {
878  // in 2d case would use distance between SVD coordinates
879 
880  // ignores V when calculating weights from a point, i.e. uses actual point values
881  double r;
882 
883  // coefficient muliplied by relative distance (allows for mid-edge weighting if applicable)
884  double alpha_weight = 1;
885  if (data._polynomial_sampling_functional==StaggeredEdgeIntegralSample
886  || data._polynomial_sampling_functional==StaggeredEdgeAnalyticGradientIntegralSample) {
887  alpha_weight = 0.5;
888  }
889 
890  // get Euchlidean distance of scaled relative coordinate from the origin
891  if (V==NULL) {
892  r = data._pc.EuclideanVectorLength(data._pc.getRelativeCoord(target_index, i, dimension) * alpha_weight, dimension);
893  } else {
894  r = data._pc.EuclideanVectorLength(data._pc.getRelativeCoord(target_index, i, dimension, V) * alpha_weight, dimension);
895  }
896 
897  // generate weight vector from distances and window sizes
898  w(i+my_num_neighbors*d) = GMLS::Wab(r, data._epsilons(target_index), data._weighting_type, data._weighting_p, data._weighting_n);
899 
900  calcPij<BasisData>(data, teamMember, delta.data(), thread_workspace.data(), target_index, i + d*my_num_neighbors, 0 /*alpha*/, dimension, polynomial_order, false /*bool on only specific order*/, V, reconstruction_space, polynomial_sampling_functional);
901 
902  if (weight_p) {
903  for (int j = 0; j < storage_size; ++j) {
904  // no need to convert offsets to global indices because the sum will never be large
905  P(i+my_num_neighbors*d, j) = delta[j] * std::sqrt(w(i+my_num_neighbors*d));
906  compadre_kernel_assert_extreme_debug(delta[j]==delta[j] && "NaN in sqrt(W)*P matrix.");
907  }
908 
909  } else {
910  for (int j = 0; j < storage_size; ++j) {
911  // no need to convert offsets to global indices because the sum will never be large
912  P(i+my_num_neighbors*d, j) = delta[j];
913 
914  compadre_kernel_assert_extreme_debug(delta[j]==delta[j] && "NaN in P matrix.");
915  }
916  }
917  }
918  });
919  teamMember.team_barrier();
920 // Kokkos::single(Kokkos::PerTeam(teamMember), [=] () {
921 // for (int k=0; k<data._pc._nla.getNumberOfNeighborsDevice(target_index); k++) {
922 // for (int l=0; l<_NP; l++) {
923 // printf("%i %i %0.16f\n", k, l, P(k,l) );// << " " << l << " " << R(k,l) << std::endl;
924 // }
925 // }
926 // });
927 }
928 
929 /*! \brief Fills the _P matrix with P*sqrt(w) for use in solving for curvature
930 
931  Uses _curvature_poly_order as the polynomial order of the basis
932 
933  \param data [in] - GMLSBasisData struct
934  \param teamMember [in] - Kokkos::TeamPolicy member type (created by parallel_for)
935  \param delta [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large is the
936 s_multipler*the dimension of the polynomial basis.
937  \param thread_workspace [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the
938 _order*the spatial dimension of the polynomial basis.
939  \param P [out] - 2D Kokkos View which will contain evaluation of sampling functional on polynomial basis for each neighbor the target has (stored column major)
940  \param w [out] - 1D Kokkos View which will contain weighting kernel values for the target with each neighbor if weight_p = true
941  \param dimension [in] - spatial dimension of basis to evaluate. e.g. dimension two basis of order one is 1, x, y, whereas for dimension 3 it is 1, x, y, z
942  \param only_specific_order [in] - boolean for only evaluating one degree of polynomial when true
943  \param V [in] - orthonormal basis matrix size _dimensions * _dimensions whose first _dimensions-1 columns are an approximation of the tangent plane
944 */
945 template <typename BasisData>
946 KOKKOS_INLINE_FUNCTION
947 void createWeightsAndPForCurvature(const BasisData& data, const member_type& teamMember, scratch_vector_type delta, scratch_vector_type thread_workspace, scratch_matrix_right_type P, scratch_vector_type w, const int dimension, bool only_specific_order, scratch_matrix_right_type* V = NULL) {
948 /*
949  * This function has two purposes
950  * 1.) Used to calculate specifically for 1st order polynomials, from which we can reconstruct a tangent plane
951  * 2.) Used to calculate a polynomial of data._curvature_poly_order, which we use to calculate curvature of the manifold
952  */
953 
954  const int target_index = data._initial_index_for_batch + teamMember.league_rank();
955  int storage_size = only_specific_order ? GMLS::getNP(1, dimension)-GMLS::getNP(0, dimension) : GMLS::getNP(data._curvature_poly_order, dimension);
956  for (int j = 0; j < delta.extent_int(0); ++j) {
957  delta(j) = 0;
958  }
959  for (int j = 0; j < thread_workspace.extent_int(0); ++j) {
960  thread_workspace(j) = 0;
961  }
962  Kokkos::parallel_for(Kokkos::TeamThreadRange(teamMember,data._pc._nla.getNumberOfNeighborsDevice(target_index)),
963  [&] (const int i) {
964 
965  // ignores V when calculating weights from a point, i.e. uses actual point values
966  double r;
967 
968  // get Euclidean distance of scaled relative coordinate from the origin
969  if (V==NULL) {
970  r = data._pc.EuclideanVectorLength(data._pc.getRelativeCoord(target_index, i, dimension), dimension);
971  } else {
972  r = data._pc.EuclideanVectorLength(data._pc.getRelativeCoord(target_index, i, dimension, V), dimension);
973  }
974 
975  // generate weight vector from distances and window sizes
976  if (only_specific_order) {
977  w(i) = GMLS::Wab(r, data._epsilons(target_index), data._curvature_weighting_type, data._curvature_weighting_p, data._curvature_weighting_n);
978  calcPij<BasisData>(data, teamMember, delta.data(), thread_workspace.data(), target_index, i, 0 /*alpha*/, dimension, 1, true /*bool on only specific order*/);
979  } else {
980  w(i) = GMLS::Wab(r, data._epsilons(target_index), data._curvature_weighting_type, data._curvature_weighting_p, data._curvature_weighting_n);
981  calcPij<BasisData>(data, teamMember, delta.data(), thread_workspace.data(), target_index, i, 0 /*alpha*/, dimension, data._curvature_poly_order, false /*bool on only specific order*/, V);
982  }
983 
984  for (int j = 0; j < storage_size; ++j) {
985  P(i, j) = delta[j] * std::sqrt(w(i));
986  }
987 
988  });
989  teamMember.team_barrier();
990 }
991 
992 } // Compadre
993 #endif
constexpr SamplingFunctional FaceNormalIntegralSample
For integrating polynomial dotted with normal over an edge.
constexpr SamplingFunctional EdgeTangentAverageSample
For polynomial dotted with tangent.
KOKKOS_INLINE_FUNCTION void evaluateSecondPartialDerivative(const member_type &teamMember, double *delta, double *workspace, const int dimension, const int max_degree, const int component, const int partial_direction_1, const int partial_direction_2, const double h, const double x, const double y, const double z, const int starting_order=0, const double weight_of_original_value=0.0, const double weight_of_new_value=1.0)
Evaluates the second partial derivatives of the Bernstein polynomial basis delta[j] = weight_of_origi...
ReconstructionSpace
Space in which to reconstruct polynomial.
Scalar basis reused as many times as there are components in the vector resulting in a much cheaper p...
#define compadre_kernel_assert_release(condition)
compadre_kernel_assert_release is similar to compadre_assert_release, but is a call on the device...
constexpr SamplingFunctional ManifoldVectorPointSample
Point evaluations of the entire vector source function (but on a manifold, so it includes a transform...
KOKKOS_INLINE_FUNCTION void evaluate(const member_type &teamMember, double *delta, double *workspace, const int dimension, const int max_degree, const int component, const double h, const double x, const double y, const double z, const int starting_order=0, const double weight_of_original_value=0.0, const double weight_of_new_value=1.0)
Evaluates the divergence-free polynomial basis delta[j] = weight_of_original_value * delta[j] + weigh...
team_policy::member_type member_type
#define compadre_kernel_assert_extreme_debug(condition)
constexpr SamplingFunctional VaryingManifoldVectorPointSample
For integrating polynomial dotted with normal over an edge.
KOKKOS_INLINE_FUNCTION double getAreaFromVectors(const member_type &teamMember, view_type_1 v1, view_type_2 v2)
KOKKOS_INLINE_FUNCTION void calcHessianPij(const BasisData &data, const member_type &teamMember, double *delta, double *thread_workspace, const int target_index, int neighbor_index, const double alpha, const int partial_direction_1, const int partial_direction_2, const int dimension, const int poly_order, bool specific_order_only, const scratch_matrix_right_type *V, const ReconstructionSpace reconstruction_space, const SamplingFunctional polynomial_sampling_functional, const int evaluation_site_local_index=0)
Evaluates the Hessian of a polynomial basis under the Dirac Delta (pointwise) sampling function...
KOKKOS_INLINE_FUNCTION void evaluate(const member_type &teamMember, double *delta, double *workspace, const int dimension, const int max_degree, const double h, const double x, const double y, const double z, const int starting_order=0, const double weight_of_original_value=0.0, const double weight_of_new_value=1.0)
Evaluates the scalar Taylor polynomial basis delta[j] = weight_of_original_value * delta[j] + weight_...
Vector polynomial basis having # of components _dimensions, or (_dimensions-1) in the case of manifol...
constexpr SamplingFunctional StaggeredEdgeAnalyticGradientIntegralSample
Analytical integral of a gradient source vector is just a difference of the scalar source at neighbor...
KOKKOS_INLINE_FUNCTION void calcGradientPij(const BasisData &data, const member_type &teamMember, double *delta, double *thread_workspace, const int target_index, int neighbor_index, const double alpha, const int partial_direction, const int dimension, const int poly_order, bool specific_order_only, const scratch_matrix_right_type *V, const ReconstructionSpace reconstruction_space, const SamplingFunctional polynomial_sampling_functional, const int evaluation_site_local_index=0)
Evaluates the gradient of a polynomial basis under the Dirac Delta (pointwise) sampling function...
Solve GMLS problem on a manifold (will use QR or SVD to solve the resultant GMLS problem dependent on...
static KOKKOS_INLINE_FUNCTION int getNP(const int m, const int dimension=3, const ReconstructionSpace r_space=ReconstructionSpace::ScalarTaylorPolynomial)
Returns size of the basis for a given polynomial order and dimension General to dimension 1...
constexpr SamplingFunctional StaggeredEdgeIntegralSample
Samples consist of the result of integrals of a vector dotted with the tangent along edges between ne...
Kokkos::View< double **, layout_right, Kokkos::MemoryTraits< Kokkos::Unmanaged > > scratch_matrix_right_type
Kokkos::View< double *, Kokkos::MemoryTraits< Kokkos::Unmanaged > > scratch_vector_type
Scalar polynomial basis centered at the target site and scaled by sum of basis powers e...
constexpr SamplingFunctional FaceNormalAverageSample
For polynomial dotted with normal on edge.
import subprocess import os import re import math import sys import argparse d d d
constexpr SamplingFunctional EdgeTangentIntegralSample
For integrating polynomial dotted with tangent over an edge.
KOKKOS_INLINE_FUNCTION void calcPij(const BasisData &data, const member_type &teamMember, double *delta, double *thread_workspace, const int target_index, int neighbor_index, const double alpha, const int dimension, const int poly_order, bool specific_order_only=false, const scratch_matrix_right_type *V=NULL, const ReconstructionSpace reconstruction_space=ReconstructionSpace::ScalarTaylorPolynomial, const SamplingFunctional polynomial_sampling_functional=PointSample, const int evaluation_site_local_index=0)
Evaluates the polynomial basis under a particular sampling function. Generally used to fill a row of ...
constexpr SamplingFunctional CellAverageSample
For polynomial integrated on cells.
KOKKOS_INLINE_FUNCTION void evaluatePartialDerivative(const member_type &teamMember, double *delta, double *workspace, const int dimension, const int max_degree, const int component, const int partial_direction, const double h, const double x, const double y, const double z, const int starting_order=0, const double weight_of_original_value=0.0, const double weight_of_new_value=1.0)
Evaluates the first partial derivatives of the Bernstein polynomial basis delta[j] = weight_of_origin...
KOKKOS_INLINE_FUNCTION void createWeightsAndP(const BasisData &data, const member_type &teamMember, scratch_vector_type delta, scratch_vector_type thread_workspace, scratch_matrix_right_type P, scratch_vector_type w, const int dimension, int polynomial_order, bool weight_p=false, scratch_matrix_right_type *V=NULL, const ReconstructionSpace reconstruction_space=ReconstructionSpace::ScalarTaylorPolynomial, const SamplingFunctional polynomial_sampling_functional=PointSample)
Fills the _P matrix with either P or P*sqrt(w)
Divergence-free vector polynomial basis.
KOKKOS_INLINE_FUNCTION void createWeightsAndPForCurvature(const BasisData &data, const member_type &teamMember, scratch_vector_type delta, scratch_vector_type thread_workspace, scratch_matrix_right_type P, scratch_vector_type w, const int dimension, bool only_specific_order, scratch_matrix_right_type *V=NULL)
Fills the _P matrix with P*sqrt(w) for use in solving for curvature.
KOKKOS_INLINE_FUNCTION void evaluate(const member_type &teamMember, double *delta, double *workspace, const int dimension, const int max_degree, const int component, const double h, const double x, const double y, const double z, const int starting_order=0, const double weight_of_original_value=0.0, const double weight_of_new_value=1.0)
Evaluates the Bernstein polynomial basis delta[j] = weight_of_original_value * delta[j] + weight_of_n...
constexpr SamplingFunctional PointSample
Available sampling functionals.
constexpr SamplingFunctional CellIntegralSample
For polynomial integrated on cells.
constexpr SamplingFunctional VectorPointSample
Point evaluations of the entire vector source function.
#define compadre_kernel_assert_debug(condition)
Bernstein polynomial basis.