Compadre  1.2.0
Compadre_GMLS_Basis.hpp
Go to the documentation of this file.
1 #ifndef _COMPADRE_GMLS_BASIS_HPP_
2 #define _COMPADRE_GMLS_BASIS_HPP_
3 
4 #include "Compadre_GMLS.hpp"
5 
6 namespace Compadre {
7 
8 KOKKOS_INLINE_FUNCTION
9 void GMLS::calcPij(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, const scratch_matrix_right_type* V, const ReconstructionSpace reconstruction_space, const SamplingFunctional polynomial_sampling_functional, const int additional_evaluation_local_index) const {
10 /*
11  * This class is under two levels of hierarchical parallelism, so we
12  * do not put in any finer grain parallelism in this function
13  */
14  const int my_num_neighbors = this->getNNeighbors(target_index);
15 
16  // store precalculated factorials for speedup
17  const double factorial[15] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600, 6227020800, 87178291200};
18 
19  int component = 0;
20  if (neighbor_index >= my_num_neighbors) {
21  component = neighbor_index / my_num_neighbors;
22  neighbor_index = neighbor_index % my_num_neighbors;
23  } else if (neighbor_index < 0) {
24  // -1 maps to 0 component
25  // -2 maps to 1 component
26  // -3 maps to 2 component
27  component = -(neighbor_index+1);
28  }
29 
30  XYZ relative_coord;
31  if (neighbor_index > -1) {
32  // Evaluate at neighbor site
33  for (int i=0; i<dimension; ++i) {
34  // calculates (alpha*target+(1-alpha)*neighbor)-1*target = (alpha-1)*target + (1-alpha)*neighbor
35  relative_coord[i] = (alpha-1)*getTargetCoordinate(target_index, i, V);
36  relative_coord[i] += (1-alpha)*getNeighborCoordinate(target_index, neighbor_index, i, V);
37  }
38  } else if (additional_evaluation_local_index > 0) {
39  // Extra evaluation site
40  for (int i=0; i<dimension; ++i) {
41  relative_coord[i] = getTargetAuxiliaryCoordinate(target_index, additional_evaluation_local_index, i, V);
42  relative_coord[i] -= getTargetCoordinate(target_index, i, V);
43  }
44  } else {
45  // Evaluate at the target site
46  for (int i=0; i<3; ++i) relative_coord[i] = 0;
47  }
48 
49  // basis ActualReconstructionSpaceRank is 0 (evaluated like a scalar) and sampling functional is traditional
50  if ((polynomial_sampling_functional == PointSample ||
51  polynomial_sampling_functional == VectorPointSample ||
52  polynomial_sampling_functional == ManifoldVectorPointSample ||
53  polynomial_sampling_functional == VaryingManifoldVectorPointSample)&&
54  (reconstruction_space == ScalarTaylorPolynomial || reconstruction_space == VectorOfScalarClonesTaylorPolynomial)) {
55 
56  double cutoff_p = _epsilons(target_index);
57  const int start_index = specific_order_only ? poly_order : 0; // only compute specified order if requested
58 
59  ScalarTaylorPolynomialBasis::evaluate(teamMember, delta, thread_workspace, dimension, poly_order, cutoff_p, relative_coord.x, relative_coord.y, relative_coord.z, start_index);
60 
61  // basis ActualReconstructionSpaceRank is 1 (is a true vector basis) and sampling functional is traditional
62  } else if ((polynomial_sampling_functional == VectorPointSample ||
63  polynomial_sampling_functional == ManifoldVectorPointSample ||
64  polynomial_sampling_functional == VaryingManifoldVectorPointSample) &&
65  (reconstruction_space == VectorTaylorPolynomial)) {
66 
67  const int dimension_offset = this->getNP(_poly_order, dimension, reconstruction_space);
68  double cutoff_p = _epsilons(target_index);
69  const int start_index = specific_order_only ? poly_order : 0; // only compute specified order if requested
70 
71  for (int d=0; d<dimension; ++d) {
72  if (d==component) {
73  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);
74  } else {
75  for (int n=0; n<dimension_offset; ++n) {
76  *(delta+d*dimension_offset+n) = 0;
77  }
78  }
79  }
80  } else if ((polynomial_sampling_functional == VectorPointSample) &&
81  (reconstruction_space == DivergenceFreeVectorTaylorPolynomial)) {
82  // Divergence free vector polynomial basis
83  double cutoff_p = _epsilons(target_index);
84 
85  DivergenceFreePolynomialBasis::evaluate(teamMember, delta, thread_workspace, dimension, poly_order, component, cutoff_p, relative_coord.x, relative_coord.y, relative_coord.z);
86 
87  } else if ((polynomial_sampling_functional == StaggeredEdgeAnalyticGradientIntegralSample) &&
88  (reconstruction_space == ScalarTaylorPolynomial)) {
89  double cutoff_p = _epsilons(target_index);
90  const int start_index = specific_order_only ? poly_order : 0; // only compute specified order if requested
91  // basis is actually scalar with staggered sampling functional
92  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);
93  relative_coord.x = 0;
94  relative_coord.y = 0;
95  relative_coord.z = 0;
96  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);
97  } else if (polynomial_sampling_functional == StaggeredEdgeIntegralSample) {
98  Kokkos::single(Kokkos::PerThread(teamMember), [&] () {
100  double cutoff_p = _epsilons(target_index);
101  int alphax, alphay;
102  double alphaf;
103  const int start_index = specific_order_only ? poly_order : 0; // only compute specified order if requested
104 
105  for (int quadrature = 0; quadrature<_qm.getNumberOfQuadraturePoints(); ++quadrature) {
106 
107  int i = 0;
108 
109  XYZ tangent_quadrature_coord_2d;
110  XYZ quadrature_coord_2d;
111  for (int j=0; j<dimension; ++j) {
112  // calculates (alpha*target+(1-alpha)*neighbor)-1*target = (alpha-1)*target + (1-alpha)*neighbor
113  quadrature_coord_2d[j] = (_qm.getSite(quadrature,0)-1)*getTargetCoordinate(target_index, j, V);
114  quadrature_coord_2d[j] += (1-_qm.getSite(quadrature,0))*getNeighborCoordinate(target_index, neighbor_index, j, V);
115  tangent_quadrature_coord_2d[j] = getTargetCoordinate(target_index, j, V);
116  tangent_quadrature_coord_2d[j] += -getNeighborCoordinate(target_index, neighbor_index, j, V);
117  }
118  for (int j=0; j<_basis_multiplier; ++j) {
119  for (int n = start_index; n <= poly_order; n++){
120  for (alphay = 0; alphay <= n; alphay++){
121  alphax = n - alphay;
122  alphaf = factorial[alphax]*factorial[alphay];
123 
124  // local evaluation of vector [0,p] or [p,0]
125  double v0, v1;
126  v0 = (j==0) ? std::pow(quadrature_coord_2d.x/cutoff_p,alphax)
127  *std::pow(quadrature_coord_2d.y/cutoff_p,alphay)/alphaf : 0;
128  v1 = (j==0) ? 0 : std::pow(quadrature_coord_2d.x/cutoff_p,alphax)
129  *std::pow(quadrature_coord_2d.y/cutoff_p,alphay)/alphaf;
130 
131  double dot_product = tangent_quadrature_coord_2d[0]*v0 + tangent_quadrature_coord_2d[1]*v1;
132 
133  // multiply by quadrature weight
134  if (quadrature==0) {
135  *(delta+i) = dot_product * _qm.getWeight(quadrature);
136  } else {
137  *(delta+i) += dot_product * _qm.getWeight(quadrature);
138  }
139  i++;
140  }
141  }
142  }
143  }
144  } else {
145  // Calculate basis matrix for NON MANIFOLD problems
146  double cutoff_p = _epsilons(target_index);
147  int alphax, alphay, alphaz;
148  double alphaf;
149  const int start_index = specific_order_only ? poly_order : 0; // only compute specified order if requested
150 
151  for (int quadrature = 0; quadrature<_qm.getNumberOfQuadraturePoints(); ++quadrature) {
152 
153  int i = 0;
154 
155  XYZ quadrature_coord_3d;
156  XYZ tangent_quadrature_coord_3d;
157  for (int j=0; j<dimension; ++j) {
158  // calculates (alpha*target+(1-alpha)*neighbor)-1*target = (alpha-1)*target + (1-alpha)*neighbor
159  quadrature_coord_3d[j] = (_qm.getSite(quadrature,0)-1)*getTargetCoordinate(target_index, j, NULL);
160  quadrature_coord_3d[j] += (1-_qm.getSite(quadrature,0))*getNeighborCoordinate(target_index, neighbor_index, j, NULL);
161  tangent_quadrature_coord_3d[j] = getTargetCoordinate(target_index, j, NULL);
162  tangent_quadrature_coord_3d[j] += -getNeighborCoordinate(target_index, neighbor_index, j, NULL);
163  }
164  for (int j=0; j<_basis_multiplier; ++j) {
165  for (int n = start_index; n <= poly_order; n++) {
166  if (dimension == 3) {
167  for (alphaz = 0; alphaz <= n; alphaz++){
168  int s = n - alphaz;
169  for (alphay = 0; alphay <= s; alphay++){
170  alphax = s - alphay;
171  alphaf = factorial[alphax]*factorial[alphay]*factorial[alphaz];
172 
173  // local evaluation of vector [p, 0, 0], [0, p, 0] or [0, 0, p]
174  double v0, v1, v2;
175  switch(j) {
176  case 1:
177  v0 = 0.0;
178  v1 = std::pow(quadrature_coord_3d.x/cutoff_p,alphax)
179  *std::pow(quadrature_coord_3d.y/cutoff_p,alphay)
180  *std::pow(quadrature_coord_3d.z/cutoff_p,alphaz)/alphaf;
181  v2 = 0.0;
182  break;
183  case 2:
184  v0 = 0.0;
185  v1 = 0.0;
186  v2 = std::pow(quadrature_coord_3d.x/cutoff_p,alphax)
187  *std::pow(quadrature_coord_3d.y/cutoff_p,alphay)
188  *std::pow(quadrature_coord_3d.z/cutoff_p,alphaz)/alphaf;
189  break;
190  default:
191  v0 = std::pow(quadrature_coord_3d.x/cutoff_p,alphax)
192  *std::pow(quadrature_coord_3d.y/cutoff_p,alphay)
193  *std::pow(quadrature_coord_3d.z/cutoff_p,alphaz)/alphaf;
194  v1 = 0.0;
195  v2 = 0.0;
196  break;
197  }
198 
199  double dot_product = tangent_quadrature_coord_3d[0]*v0 + tangent_quadrature_coord_3d[1]*v1 + tangent_quadrature_coord_3d[2]*v2;
200 
201  // multiply by quadrature weight
202  if (quadrature == 0) {
203  *(delta+i) = dot_product * _qm.getWeight(quadrature);
204  } else {
205  *(delta+i) += dot_product * _qm.getWeight(quadrature);
206  }
207  i++;
208  }
209  }
210  } else if (dimension == 2) {
211  for (alphay = 0; alphay <= n; alphay++){
212  alphax = n - alphay;
213  alphaf = factorial[alphax]*factorial[alphay];
214 
215  // local evaluation of vector [p, 0] or [0, p]
216  double v0, v1;
217  switch(j) {
218  case 1:
219  v0 = 0.0;
220  v1 = std::pow(quadrature_coord_3d.x/cutoff_p,alphax)
221  *std::pow(quadrature_coord_3d.y/cutoff_p,alphay)/alphaf;
222  break;
223  default:
224  v0 = std::pow(quadrature_coord_3d.x/cutoff_p,alphax)
225  *std::pow(quadrature_coord_3d.y/cutoff_p,alphay)/alphaf;
226  v1 = 0.0;
227  break;
228  }
229 
230  double dot_product = tangent_quadrature_coord_3d[0]*v0 + tangent_quadrature_coord_3d[1]*v1;
231 
232  // multiply by quadrature weight
233  if (quadrature == 0) {
234  *(delta+i) = dot_product * _qm.getWeight(quadrature);
235  } else {
236  *(delta+i) += dot_product * _qm.getWeight(quadrature);
237  }
238  i++;
239  }
240  }
241  }
242  }
243  }
244  } // NON MANIFOLD PROBLEMS
245  });
246  } else if (polynomial_sampling_functional == FaceNormalIntegralSample ||
247  polynomial_sampling_functional == FaceTangentIntegralSample ||
248  polynomial_sampling_functional == FaceNormalPointSample ||
249  polynomial_sampling_functional == FaceTangentPointSample) {
250 
251  double cutoff_p = _epsilons(target_index);
252 
253  compadre_kernel_assert_debug(_dimensions==2 && "Only written for 2D");
254  compadre_kernel_assert_debug(_source_extra_data.extent(0)>0 && "Extra data used but not set.");
255 
256  int neighbor_index_in_source = getNeighborIndex(target_index, neighbor_index);
257 
258  /*
259  * requires quadrature points defined on an edge, not a target/source edge (spoke)
260  *
261  * _extra_data will contain the endpoints (2 for 2D, 3 for 3D) and then the unit normals
262  * (e0_x, e0_y, e1_x, e1_y, n_x, n_y, t_x, t_y)
263  */
264 
265  // if not integrating, set to 1
266  int quadrature_point_loop = (polynomial_sampling_functional == FaceNormalIntegralSample
267  || polynomial_sampling_functional == FaceTangentIntegralSample) ?
269 
270  // only used for integrated quantities
271  XYZ endpoints_difference;
272  for (int j=0; j<dimension; ++j) {
273  endpoints_difference[j] = _source_extra_data(neighbor_index_in_source, j) - _source_extra_data(neighbor_index_in_source, j+2);
274  }
275  double magnitude = EuclideanVectorLength(endpoints_difference, 2);
276 
277  int alphax, alphay;
278  double alphaf;
279  const int start_index = specific_order_only ? poly_order : 0; // only compute specified order if requested
280 
281  // loop
282  for (int quadrature = 0; quadrature<quadrature_point_loop; ++quadrature) {
283 
284  int i = 0;
285 
286  XYZ direction_2d;
287  XYZ quadrature_coord_2d;
288  for (int j=0; j<dimension; ++j) {
289 
290  if (polynomial_sampling_functional == FaceNormalIntegralSample
291  || polynomial_sampling_functional == FaceTangentIntegralSample) {
292  // quadrature coord site
293  quadrature_coord_2d[j] = _qm.getSite(quadrature,0)*_source_extra_data(neighbor_index_in_source, j);
294  quadrature_coord_2d[j] += (1-_qm.getSite(quadrature,0))*_source_extra_data(neighbor_index_in_source, j+2);
295  quadrature_coord_2d[j] -= getTargetCoordinate(target_index, j);
296  } else {
297  // traditional coord
298  quadrature_coord_2d[j] = relative_coord[j];
299  }
300 
301  // normal direction or tangent direction
302  if (polynomial_sampling_functional == FaceNormalIntegralSample
303  || polynomial_sampling_functional == FaceNormalPointSample) {
304  // normal direction
305  direction_2d[j] = _source_extra_data(neighbor_index_in_source, 4 + j);
306  } else {
307  // tangent direction
308  direction_2d[j] = _source_extra_data(neighbor_index_in_source, 6 + j);
309  }
310 
311  }
312 
313  for (int j=0; j<_basis_multiplier; ++j) {
314  for (int n = start_index; n <= poly_order; n++){
315  for (alphay = 0; alphay <= n; alphay++){
316  alphax = n - alphay;
317  alphaf = factorial[alphax]*factorial[alphay];
318 
319  // local evaluation of vector [0,p] or [p,0]
320  double v0, v1;
321  v0 = (j==0) ? std::pow(quadrature_coord_2d.x/cutoff_p,alphax)
322  *std::pow(quadrature_coord_2d.y/cutoff_p,alphay)/alphaf : 0;
323  v1 = (j==0) ? 0 : std::pow(quadrature_coord_2d.x/cutoff_p,alphax)
324  *std::pow(quadrature_coord_2d.y/cutoff_p,alphay)/alphaf;
325 
326  // either n*v or t*v
327  double dot_product = direction_2d[0]*v0 + direction_2d[1]*v1;
328 
329  // multiply by quadrature weight
330  if (quadrature==0) {
331  if (polynomial_sampling_functional == FaceNormalIntegralSample
332  || polynomial_sampling_functional == FaceTangentIntegralSample) {
333  // integral
334  *(delta+i) = dot_product * _qm.getWeight(quadrature) * magnitude;
335  } else {
336  // point
337  *(delta+i) = dot_product;
338  }
339  } else {
340  // non-integrated quantities never satisfy this condition
341  *(delta+i) += dot_product * _qm.getWeight(quadrature) * magnitude;
342  }
343  i++;
344  }
345  }
346  }
347  }
348  } else if (polynomial_sampling_functional == ScalarFaceAverageSample) {
349  auto global_neighbor_index = getNeighborIndex(target_index, neighbor_index);
350  double cutoff_p = _epsilons(target_index);
351  int alphax, alphay, alphaz;
352  double alphaf;
353 
354  // global dimension cannot be determined in a constexpr way, so we use a largest case scenario
355  // of dimensions 3 for _global_dimension
356  double triangle_coords[3/*_global_dimensions*/*3];
357  for (int i=0; i<_global_dimensions*3; ++i) triangle_coords[i] = 0;
358  // 3 is for # vertices in sub-triangle
359  scratch_matrix_right_type triangle_coords_matrix(triangle_coords, _global_dimensions, 3);
360 
361  scratch_vector_type midpoint(delta, _global_dimensions);
362  getMidpointFromCellVertices(teamMember, midpoint, _source_extra_data, global_neighbor_index, _global_dimensions /*dim*/);
363  for (int j=0; j<_global_dimensions; ++j) {
364  triangle_coords_matrix(j, 0) = midpoint(j);
365  }
366 
367  size_t num_vertices = _source_extra_data.extent(1) / _global_dimensions;
368  double reference_cell_area = 0.5;
369  double entire_cell_area = 0.0;
370  auto T=triangle_coords_matrix;
371 
372  for (size_t v=0; v<num_vertices; ++v) {
373  int v1 = v;
374  int v2 = (v+1) % num_vertices;
375  for (int j=0; j<_global_dimensions; ++j) {
376  triangle_coords_matrix(j,1) = _source_extra_data(global_neighbor_index, v1*_global_dimensions+j) - triangle_coords_matrix(j,0);
377  triangle_coords_matrix(j,2) = _source_extra_data(global_neighbor_index, v2*_global_dimensions+j) - triangle_coords_matrix(j,0);
378  }
379  entire_cell_area += 0.5 * getAreaFromVectors(teamMember,
380  Kokkos::subview(T, Kokkos::ALL(), 1), Kokkos::subview(T, Kokkos::ALL(), 2));
381  }
382 
383  // loop over each two vertices
384  // made for flat surfaces (either dim=2 or on a manifold)
385  for (size_t v=0; v<num_vertices; ++v) {
386  int v1 = v;
387  int v2 = (v+1) % num_vertices;
388 
389  for (int j=0; j<_global_dimensions; ++j) {
390  triangle_coords_matrix(j,1) = _source_extra_data(global_neighbor_index, v1*_global_dimensions+j) - triangle_coords_matrix(j,0);
391  triangle_coords_matrix(j,2) = _source_extra_data(global_neighbor_index, v2*_global_dimensions+j) - triangle_coords_matrix(j,0);
392  }
393  // triangle_coords now has:
394  // (midpoint_x, midpoint_y, midpoint_z,
395  // v1_x-midpoint_x, v1_y-midpoint_y, v1_z-midpoint_z,
396  // v2_x-midpoint_x, v2_y-midpoint_y, v2_z-midpoint_z);
397  for (int quadrature = 0; quadrature<_qm.getNumberOfQuadraturePoints(); ++quadrature) {
398  double transformed_qp[3] = {0,0,0};
399  for (int j=0; j<_global_dimensions; ++j) {
400  for (int k=1; k<3; ++k) { // 3 is for # vertices in subtriangle
401  transformed_qp[j] += T(j,k)*_qm.getSite(quadrature, k-1);
402  }
403  transformed_qp[j] += T(j,0);
404  }
405  // half the norm of the cross-product is the area of the triangle
406  // so scaling is area / reference area (0.5) = the norm of the cross-product
407  double sub_cell_area = 0.5 * getAreaFromVectors(teamMember,
408  Kokkos::subview(T, Kokkos::ALL(), 1), Kokkos::subview(T, Kokkos::ALL(), 2));
409  double scaling_factor = sub_cell_area / reference_cell_area;
410 
412  XYZ qp = XYZ(transformed_qp[0], transformed_qp[1], transformed_qp[2]);
413  for (int j=0; j<2; ++j) {
414  relative_coord[j] = convertGlobalToLocalCoordinate(qp,j,V) - getTargetCoordinate(target_index,j,V); // shift quadrature point by target site
415  relative_coord[2] = 0;
416  }
417  } else {
418  for (int j=0; j<dimension; ++j) {
419  relative_coord[j] = transformed_qp[j] - getTargetCoordinate(target_index,j,V); // shift quadrature point by target site
420  }
421  for (int j=dimension; j<3; ++j) {
422  relative_coord[j] = 0.0;
423  }
424  }
425 
426  int k = 0;
427  const int start_index = specific_order_only ? poly_order : 0; // only compute specified order if requested
428  if (dimension == 3) {
429  for (int n = start_index; n <= poly_order; n++){
430  for (alphaz = 0; alphaz <= n; alphaz++){
431  int s = n - alphaz;
432  for (alphay = 0; alphay <= s; alphay++){
433  alphax = s - alphay;
434  alphaf = factorial[alphax]*factorial[alphay]*factorial[alphaz];
435  double val_to_sum = (scaling_factor * _qm.getWeight(quadrature)
436  * std::pow(relative_coord.x/cutoff_p,alphax)
437  * std::pow(relative_coord.y/cutoff_p,alphay)
438  * std::pow(relative_coord.z/cutoff_p,alphaz)/alphaf) / entire_cell_area;
439  if (quadrature==0 && v==0) *(delta+k) = val_to_sum;
440  else *(delta+k) += val_to_sum;
441  k++;
442  }
443  }
444  }
445  } else if (dimension == 2) {
446  for (int n = start_index; n <= poly_order; n++){
447  for (alphay = 0; alphay <= n; alphay++){
448  alphax = n - alphay;
449  alphaf = factorial[alphax]*factorial[alphay];
450  double val_to_sum = (scaling_factor * _qm.getWeight(quadrature)
451  * std::pow(relative_coord.x/cutoff_p,alphax)
452  * std::pow(relative_coord.y/cutoff_p,alphay)/alphaf) / entire_cell_area;
453  if (quadrature==0 && v==0) *(delta+k) = val_to_sum;
454  else *(delta+k) += val_to_sum;
455  k++;
456  }
457  }
458  }
459  }
460  }
461  } else {
462  compadre_kernel_assert_release((false) && "Sampling and basis space combination not defined.");
463  }
464 }
465 
466 
467 KOKKOS_INLINE_FUNCTION
468 void GMLS::calcGradientPij(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 additional_evaluation_index) const {
469 /*
470  * This class is under two levels of hierarchical parallelism, so we
471  * do not put in any finer grain parallelism in this function
472  */
473 
474  const int my_num_neighbors = this->getNNeighbors(target_index);
475 
476  int component = 0;
477  if (neighbor_index >= my_num_neighbors) {
478  component = neighbor_index / my_num_neighbors;
479  neighbor_index = neighbor_index % my_num_neighbors;
480  } else if (neighbor_index < 0) {
481  // -1 maps to 0 component
482  // -2 maps to 1 component
483  // -3 maps to 2 component
484  component = -(neighbor_index+1);
485  }
486 
487  // alpha corresponds to the linear combination of target_index and neighbor_index coordinates
488  // coordinate to evaluate = alpha*(target_index's coordinate) + (1-alpha)*(neighbor_index's coordinate)
489 
490  // partial_direction - 0=x, 1=y, 2=z
491  XYZ relative_coord;
492  if (neighbor_index > -1) {
493  for (int i=0; i<dimension; ++i) {
494  // calculates (alpha*target+(1-alpha)*neighbor)-1*target = (alpha-1)*target + (1-alpha)*neighbor
495  relative_coord[i] = (alpha-1)*getTargetCoordinate(target_index, i, V);
496  relative_coord[i] += (1-alpha)*getNeighborCoordinate(target_index, neighbor_index, i, V);
497  }
498  } else if (additional_evaluation_index > 0) {
499  for (int i=0; i<dimension; ++i) {
500  relative_coord[i] = getTargetAuxiliaryCoordinate(target_index, additional_evaluation_index, i, V);
501  relative_coord[i] -= getTargetCoordinate(target_index, i, V);
502  }
503  } else {
504  for (int i=0; i<3; ++i) relative_coord[i] = 0;
505  }
506 
507  double cutoff_p = _epsilons(target_index);
508  const int start_index = specific_order_only ? poly_order : 0; // only compute specified order if requested
509 
510  if ((polynomial_sampling_functional == PointSample ||
511  polynomial_sampling_functional == VectorPointSample ||
512  polynomial_sampling_functional == ManifoldVectorPointSample ||
513  polynomial_sampling_functional == VaryingManifoldVectorPointSample) &&
514  (reconstruction_space == ScalarTaylorPolynomial || reconstruction_space == VectorOfScalarClonesTaylorPolynomial)) {
515 
516  ScalarTaylorPolynomialBasis::evaluatePartialDerivative(teamMember, delta, thread_workspace, dimension, poly_order, partial_direction, cutoff_p, relative_coord.x, relative_coord.y, relative_coord.z, start_index);
517 
518  } else if ((polynomial_sampling_functional == VectorPointSample) &&
519  (reconstruction_space == DivergenceFreeVectorTaylorPolynomial)) {
520  // Divergence free vector polynomial basis
521  double cutoff_p = _epsilons(target_index);
522 
523  DivergenceFreePolynomialBasis::evaluatePartialDerivative(teamMember, delta, thread_workspace, dimension, poly_order, component, partial_direction, cutoff_p, relative_coord.x, relative_coord.y, relative_coord.z);
524 
525  } else {
526  compadre_kernel_assert_release((false) && "Sampling and basis space combination not defined.");
527  }
528 }
529 
530 KOKKOS_INLINE_FUNCTION
531 void GMLS::calcHessianPij(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 additional_evaluation_index) const {
532 /*
533  * This class is under two levels of hierarchical parallelism, so we
534  * do not put in any finer grain parallelism in this function
535  */
536 
537  const int my_num_neighbors = this->getNNeighbors(target_index);
538 
539  int component = 0;
540  if (neighbor_index >= my_num_neighbors) {
541  component = neighbor_index / my_num_neighbors;
542  neighbor_index = neighbor_index % my_num_neighbors;
543  } else if (neighbor_index < 0) {
544  // -1 maps to 0 component
545  // -2 maps to 1 component
546  // -3 maps to 2 component
547  component = -(neighbor_index+1);
548  }
549 
550  // alpha corresponds to the linear combination of target_index and neighbor_index coordinates
551  // coordinate to evaluate = alpha*(target_index's coordinate) + (1-alpha)*(neighbor_index's coordinate)
552 
553  // partial_direction - 0=x, 1=y, 2=z
554  XYZ relative_coord;
555  if (neighbor_index > -1) {
556  for (int i=0; i<dimension; ++i) {
557  // calculates (alpha*target+(1-alpha)*neighbor)-1*target = (alpha-1)*target + (1-alpha)*neighbor
558  relative_coord[i] = (alpha-1)*getTargetCoordinate(target_index, i, V);
559  relative_coord[i] += (1-alpha)*getNeighborCoordinate(target_index, neighbor_index, i, V);
560  }
561  } else if (additional_evaluation_index > 0) {
562  for (int i=0; i<dimension; ++i) {
563  relative_coord[i] = getTargetAuxiliaryCoordinate(target_index, additional_evaluation_index, i, V);
564  relative_coord[i] -= getTargetCoordinate(target_index, i, V);
565  }
566  } else {
567  for (int i=0; i<3; ++i) relative_coord[i] = 0;
568  }
569 
570  double cutoff_p = _epsilons(target_index);
571 
572  if ((polynomial_sampling_functional == PointSample ||
573  polynomial_sampling_functional == VectorPointSample ||
574  polynomial_sampling_functional == ManifoldVectorPointSample ||
575  polynomial_sampling_functional == VaryingManifoldVectorPointSample) &&
576  (reconstruction_space == ScalarTaylorPolynomial || reconstruction_space == VectorOfScalarClonesTaylorPolynomial)) {
577 
578  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);
579 
580  } else if ((polynomial_sampling_functional == VectorPointSample) &&
581  (reconstruction_space == DivergenceFreeVectorTaylorPolynomial)) {
582 
583  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);
584 
585  } else {
586  compadre_kernel_assert_release((false) && "Sampling and basis space combination not defined.");
587  }
588 }
589 
590 
591 KOKKOS_INLINE_FUNCTION
592 void GMLS::createWeightsAndP(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, scratch_matrix_right_type* V, const ReconstructionSpace reconstruction_space, const SamplingFunctional polynomial_sampling_functional) const {
593  /*
594  * Creates sqrt(W)*P
595  */
596  const int target_index = _initial_index_for_batch + teamMember.league_rank();
597 // printf("specific order: %d\n", specific_order);
598 // {
599 // const int storage_size = (specific_order > 0) ? this->getNP(specific_order, dimension)-this->getNP(specific_order-1, dimension) : this->getNP(_poly_order, dimension);
600 // printf("storage size: %d\n", storage_size);
601 // }
602 // printf("weight_p: %d\n", weight_p);
603  const int my_num_neighbors = this->getNNeighbors(target_index);
604 
605  teamMember.team_barrier();
606  Kokkos::parallel_for(Kokkos::TeamThreadRange(teamMember,this->getNNeighbors(target_index)),
607  [=] (const int i) {
608 
609  for (int d=0; d<_sampling_multiplier; ++d) {
610  // in 2d case would use distance between SVD coordinates
611 
612  // ignores V when calculating weights from a point, i.e. uses actual point values
613  double r;
614 
615  // coefficient muliplied by relative distance (allows for mid-edge weighting if applicable)
616  double alpha_weight = 1;
619  alpha_weight = 0.5;
620  }
621 
622  // get Euchlidean distance of scaled relative coordinate from the origin
623  if (V==NULL) {
624  r = this->EuclideanVectorLength(this->getRelativeCoord(target_index, i, dimension) * alpha_weight, dimension);
625  } else {
626  r = this->EuclideanVectorLength(this->getRelativeCoord(target_index, i, dimension, V) * alpha_weight, dimension);
627  }
628 
629  // generate weight vector from distances and window sizes
630  w(i+my_num_neighbors*d) = this->Wab(r, _epsilons(target_index), _weighting_type, _weighting_power);
631 
632  this->calcPij(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);
633 
634  // storage_size needs to change based on the size of the basis
635  int storage_size = this->getNP(polynomial_order, dimension, reconstruction_space);
636  storage_size *= _basis_multiplier;
637 
638  if (weight_p) {
639  for (int j = 0; j < storage_size; ++j) {
640  // no need to convert offsets to global indices because the sum will never be large
641  P(i+my_num_neighbors*d, j) = delta[j] * std::sqrt(w(i+my_num_neighbors*d));
642  compadre_kernel_assert_extreme_debug(delta[j]==delta[j] && "NaN in sqrt(W)*P matrix.");
643  }
644 
645  } else {
646  for (int j = 0; j < storage_size; ++j) {
647  // no need to convert offsets to global indices because the sum will never be large
648  P(i+my_num_neighbors*d, j) = delta[j];
649 
650  compadre_kernel_assert_extreme_debug(delta[j]==delta[j] && "NaN in P matrix.");
651  }
652  }
653  }
654  });
655 
656  teamMember.team_barrier();
657 // Kokkos::single(Kokkos::PerTeam(teamMember), [=] () {
658 // for (int k=0; k<this->getNNeighbors(target_index); k++) {
659 // for (int l=0; l<_NP; l++) {
660 // printf("%i %i %0.16f\n", k, l, P(k,l) );// << " " << l << " " << R(k,l) << std::endl;
661 // }
662 // }
663 // });
664 }
665 
666 KOKKOS_INLINE_FUNCTION
667 void GMLS::createWeightsAndPForCurvature(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) const {
668 /*
669  * This function has two purposes
670  * 1.) Used to calculate specifically for 1st order polynomials, from which we can reconstruct a tangent plane
671  * 2.) Used to calculate a polynomial of _curvature_poly_order, which we use to calculate curvature of the manifold
672  */
673 
674  const int target_index = _initial_index_for_batch + teamMember.league_rank();
675 
676  teamMember.team_barrier();
677  Kokkos::parallel_for(Kokkos::TeamThreadRange(teamMember,this->getNNeighbors(target_index)),
678  [=] (const int i) {
679 
680  // ignores V when calculating weights from a point, i.e. uses actual point values
681  double r;
682 
683  // get Euclidean distance of scaled relative coordinate from the origin
684  if (V==NULL) {
685  r = this->EuclideanVectorLength(this->getRelativeCoord(target_index, i, dimension), dimension);
686  } else {
687  r = this->EuclideanVectorLength(this->getRelativeCoord(target_index, i, dimension, V), dimension);
688  }
689 
690  // generate weight vector from distances and window sizes
691  if (only_specific_order) {
692  w(i) = this->Wab(r, _epsilons(target_index), _curvature_weighting_type, _curvature_weighting_power);
693  this->calcPij(teamMember, delta.data(), thread_workspace.data(), target_index, i, 0 /*alpha*/, dimension, 1, true /*bool on only specific order*/);
694  } else {
695  w(i) = this->Wab(r, _epsilons(target_index), _curvature_weighting_type, _curvature_weighting_power);
696  this->calcPij(teamMember, delta.data(), thread_workspace.data(), target_index, i, 0 /*alpha*/, dimension, _curvature_poly_order, false /*bool on only specific order*/, V);
697  }
698 
699  int storage_size = only_specific_order ? this->getNP(1, dimension)-this->getNP(0, dimension) : this->getNP(_curvature_poly_order, dimension);
700 
701  for (int j = 0; j < storage_size; ++j) {
702  P(i, j) = delta[j] * std::sqrt(w(i));
703  }
704 
705  });
706  teamMember.team_barrier();
707 }
708 
709 } // Compadre
710 #endif
#define compadre_kernel_assert_debug(condition)
Kokkos::View< double *, Kokkos::MemoryTraits< Kokkos::Unmanaged > > scratch_vector_type
#define compadre_kernel_assert_release(condition)
compadre_kernel_assert_release is similar to compadre_assert_release, but is a call on the device,...
team_policy::member_type member_type
#define compadre_kernel_assert_extreme_debug(condition)
Kokkos::View< double **, layout_right, Kokkos::MemoryTraits< Kokkos::Unmanaged > > scratch_matrix_right_type
int _sampling_multiplier
actual dimension of the sampling functional e.g.
int _weighting_power
power to be used for weighting kernel
const SamplingFunctional _polynomial_sampling_functional
polynomial sampling functional used to construct P matrix, set at GMLS class instantiation
static KOKKOS_INLINE_FUNCTION double EuclideanVectorLength(const XYZ &delta_vector, const int dimension)
Returns Euclidean norm of a vector.
int _basis_multiplier
dimension of the reconstructed function e.g.
KOKKOS_INLINE_FUNCTION XYZ getRelativeCoord(const int target_index, const int neighbor_list_num, const int dimension, const scratch_matrix_right_type *V=NULL) const
Returns the relative coordinate as a vector between the target site and the neighbor site.
KOKKOS_INLINE_FUNCTION void createWeightsAndP(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 sampling_strategy=PointSample) const
Fills the _P matrix with either P or P*sqrt(w)
int _global_dimensions
spatial dimension of the points, set at class instantiation only
Kokkos::View< double * > _epsilons
h supports determined through neighbor search (device)
WeightingFunctionType _weighting_type
weighting kernel type for GMLS
ProblemType _problem_type
problem type for GMLS problem, can also be set to STANDARD for normal or MANIFOLD for manifold proble...
KOKKOS_INLINE_FUNCTION void createWeightsAndPForCurvature(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) const
Fills the _P matrix with P*sqrt(w) for use in solving for curvature.
Kokkos::View< double **, layout_right > _source_extra_data
Extra data available to basis functions (optional)
KOKKOS_INLINE_FUNCTION int getNNeighbors(const int target_index) const
Returns number of neighbors for a particular target.
int _initial_index_for_batch
initial index for current batch
int _curvature_weighting_power
power to be used for weighting kernel for curvature
KOKKOS_INLINE_FUNCTION void calcPij(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 sampling_strategy=PointSample, const int additional_evaluation_local_index=0) const
Evaluates the polynomial basis under a particular sampling function. Generally used to fill a row of ...
KOKKOS_INLINE_FUNCTION double convertGlobalToLocalCoordinate(const XYZ global_coord, const int dim, const scratch_matrix_right_type *V) const
Returns a component of the local coordinate after transformation from global to local under the ortho...
static KOKKOS_INLINE_FUNCTION double Wab(const double r, const double h, const WeightingFunctionType &weighting_type, const int power)
Evaluates the weighting kernel.
int _poly_order
order of basis for polynomial reconstruction
KOKKOS_INLINE_FUNCTION int getNeighborIndex(const int target_index, const int neighbor_list_num) const
Mapping from [0,number of neighbors for a target] to the row that contains the source coordinates for...
int _dimensions
dimension of the problem, set at class instantiation only
KOKKOS_INLINE_FUNCTION double getTargetAuxiliaryCoordinate(const int target_index, const int additional_list_num, const int dim, const scratch_matrix_right_type *V=NULL) const
(OPTIONAL) Returns one component of the additional evaluation coordinates.
int _curvature_poly_order
order of basis for curvature reconstruction
WeightingFunctionType _curvature_weighting_type
weighting kernel type for curvature problem
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....
KOKKOS_INLINE_FUNCTION void calcHessianPij(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 sampling_strategy, const int additional_evaluation_local_index=0) const
Evaluates the Hessian of a polynomial basis under the Dirac Delta (pointwise) sampling function.
KOKKOS_INLINE_FUNCTION void calcGradientPij(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 sampling_strategy, const int additional_evaluation_local_index=0) const
Evaluates the gradient of a polynomial basis under the Dirac Delta (pointwise) sampling function.
KOKKOS_INLINE_FUNCTION double getNeighborCoordinate(const int target_index, const int neighbor_list_num, const int dim, const scratch_matrix_right_type *V=NULL) const
Returns one component of the neighbor coordinate for a particular target.
KOKKOS_INLINE_FUNCTION double getTargetCoordinate(const int target_index, const int dim, const scratch_matrix_right_type *V=NULL) const
Returns one component of the target coordinate for a particular target.
Quadrature _qm
manages and calculates quadrature
KOKKOS_INLINE_FUNCTION double getWeight(const int index) const
KOKKOS_INLINE_FUNCTION int getNumberOfQuadraturePoints() const
KOKKOS_INLINE_FUNCTION double getSite(const int index, const int component) const
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...
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 divergence-free polynomial basis delta[j] = weight_of_...
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 divergence-free polynomial basis delta[j] = weight_of...
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_...
KOKKOS_INLINE_FUNCTION void evaluatePartialDerivative(const member_type &teamMember, double *delta, double *workspace, const int dimension, const int max_degree, 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 scalar Taylor polynomial basis delta[j] = weight_of_origin...
KOKKOS_INLINE_FUNCTION void evaluateSecondPartialDerivative(const member_type &teamMember, double *delta, double *workspace, const int dimension, const int max_degree, 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 scalar Taylor polynomial basis delta[j] = weight_of_origi...
@ MANIFOLD
Solve GMLS problem on a manifold (will use QR or SVD to solve the resultant GMLS problem dependent on...
constexpr SamplingFunctional VaryingManifoldVectorPointSample
For integrating polynomial dotted with normal over an edge.
constexpr SamplingFunctional FaceNormalIntegralSample
For integrating polynomial dotted with normal over an edge.
constexpr SamplingFunctional StaggeredEdgeIntegralSample
Samples consist of the result of integrals of a vector dotted with the tangent along edges between ne...
constexpr SamplingFunctional PointSample
Available sampling functionals.
constexpr SamplingFunctional StaggeredEdgeAnalyticGradientIntegralSample
Analytical integral of a gradient source vector is just a difference of the scalar source at neighbor...
constexpr SamplingFunctional FaceTangentIntegralSample
For integrating polynomial dotted with tangent over an edge.
constexpr SamplingFunctional FaceNormalPointSample
For polynomial dotted with normal on edge.
KOKKOS_INLINE_FUNCTION void getMidpointFromCellVertices(const member_type &teamMember, scratch_vector_type midpoint_storage, scratch_matrix_right_type cell_coordinates, const int cell_num, const int dim=3)
constexpr SamplingFunctional ManifoldVectorPointSample
Point evaluations of the entire vector source function (but on a manifold, so it includes a transform...
constexpr SamplingFunctional FaceTangentPointSample
For polynomial dotted with tangent.
constexpr SamplingFunctional VectorPointSample
Point evaluations of the entire vector source function.
constexpr SamplingFunctional ScalarFaceAverageSample
For polynomial integrated on faces.
ReconstructionSpace
Space in which to reconstruct polynomial.
@ DivergenceFreeVectorTaylorPolynomial
Divergence-free vector polynomial basis.
@ VectorTaylorPolynomial
Vector polynomial basis having # of components _dimensions, or (_dimensions-1) in the case of manifol...
@ ScalarTaylorPolynomial
Scalar polynomial basis centered at the target site and scaled by sum of basis powers e....
@ VectorOfScalarClonesTaylorPolynomial
Scalar basis reused as many times as there are components in the vector resulting in a much cheaper p...
KOKKOS_INLINE_FUNCTION double getAreaFromVectors(const member_type &teamMember, view_type_1 v1, view_type_2 v2)