34void 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) {
39 const int my_num_neighbors = data._pc._nla.getNumberOfNeighborsDevice(target_index);
42 const double factorial[15] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600, 6227020800, 87178291200};
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) {
52 component = -(neighbor_index+1);
56 if (neighbor_index > -1) {
58 for (
int i=0; i<dimension; ++i) {
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);
63 }
else if (evaluation_site_local_index > 0) {
65 for (
int i=0; i<dimension; ++i) {
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);
73 for (
int i=0; i<3; ++i) relative_coord[i] = 0;
77 if ((polynomial_sampling_functional ==
PointSample ||
83 double cutoff_p = data._epsilons(target_index);
84 const int start_index = specific_order_only ? poly_order : 0;
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;
98 for (
int d=0; d<dimension; ++d) {
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);
102 for (
int n=0; n<dimension_offset; ++n) {
103 *(delta+d*dimension_offset+n) = 0;
110 double cutoff_p = data._epsilons(target_index);
115 double cutoff_p = data._epsilons(target_index);
121 double cutoff_p = data._epsilons(target_index);
122 const int start_index = specific_order_only ? poly_order : 0;
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);
130 Kokkos::single(Kokkos::PerThread(teamMember), [&] () {
132 double cutoff_p = data._epsilons(target_index);
135 const int start_index = specific_order_only ? poly_order : 0;
137 for (int quadrature = 0; quadrature<data._qm.getNumberOfQuadraturePoints(); ++quadrature) {
141 XYZ tangent_quadrature_coord_2d;
142 XYZ quadrature_coord_2d;
143 for (int j=0; j<dimension; ++j) {
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);
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++){
154 alphaf = factorial[alphax]*factorial[alphay];
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;
163 double dot_product = tangent_quadrature_coord_2d[0]*v0 + tangent_quadrature_coord_2d[1]*v1;
167 *(delta+i) = dot_product * data._qm.getWeight(quadrature);
169 *(delta+i) += dot_product * data._qm.getWeight(quadrature);
178 double cutoff_p = data._epsilons(target_index);
179 int alphax, alphay, alphaz;
181 const int start_index = specific_order_only ? poly_order : 0;
183 for (
int quadrature = 0; quadrature<data._qm.getNumberOfQuadraturePoints(); ++quadrature) {
187 XYZ quadrature_coord_3d;
188 XYZ tangent_quadrature_coord_3d;
189 for (
int j=0; j<dimension; ++j) {
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);
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++){
201 for (alphay = 0; alphay <= s; alphay++){
203 alphaf = factorial[alphax]*factorial[alphay]*factorial[alphaz];
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;
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;
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;
231 double dot_product = tangent_quadrature_coord_3d[0]*v0 + tangent_quadrature_coord_3d[1]*v1 + tangent_quadrature_coord_3d[2]*v2;
234 if (quadrature == 0) {
235 *(delta+i) = dot_product * data._qm.getWeight(quadrature);
237 *(delta+i) += dot_product * data._qm.getWeight(quadrature);
242 }
else if (dimension == 2) {
243 for (alphay = 0; alphay <= n; alphay++){
245 alphaf = factorial[alphax]*factorial[alphay];
252 v1 = std::pow(quadrature_coord_3d.
x/cutoff_p,alphax)
253 *std::pow(quadrature_coord_3d.
y/cutoff_p,alphay)/alphaf;
256 v0 = std::pow(quadrature_coord_3d.
x/cutoff_p,alphax)
257 *std::pow(quadrature_coord_3d.
y/cutoff_p,alphay)/alphaf;
262 double dot_product = tangent_quadrature_coord_3d[0]*v0 + tangent_quadrature_coord_3d[1]*v1;
265 if (quadrature == 0) {
266 *(delta+i) = dot_product * data._qm.getWeight(quadrature);
268 *(delta+i) += dot_product * data._qm.getWeight(quadrature);
285 "FaceNormalIntegralSample, EdgeTangentIntegralSample, FaceNormalAverageSample, \
286 and EdgeTangentAverageSample only support 2d or 3d with 2d manifold");
289 &&
"Only 1d quadrature may be specified for edge integrals");
292 &&
"Quadrature points not generated");
297 "Sampling functional does not support specific_order_only");
299 double cutoff_p = data._epsilons(target_index);
300 auto global_neighbor_index = data._pc.getNeighborIndex(target_index, neighbor_index);
311 int quadrature_point_loop = data._qm.getNumberOfQuadraturePoints();
314 double G_data[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;
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);
333 radius = std::sqrt(radius);
337 auto E = edge_coords_matrix;
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];
346 theta = std::atan(norm_a_cross_b / a_dot_b);
350 double entire_edge_length = 0.0;
351 for (
int quadrature = 0; quadrature<quadrature_point_loop; ++quadrature) {
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);
359 unscaled_transformed_qp[j] += E(j,0);
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];
372 transformed_qp_norm = std::sqrt(transformed_qp_norm);
374 for (
int j=0; j<data._global_dimensions; ++j) {
375 scaled_transformed_qp[j] = unscaled_transformed_qp[j] * radius / transformed_qp_norm;
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);
385 relative_coord[2] = 0;
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);
394 relative_coord[2] = 0;
401 for (
int j=0; j<data._global_dimensions; ++j) {
403 direction[j] = data._source_extra_data(global_neighbor_index, 2*data._global_dimensions + j);
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)};
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;
420 for (
int j=0; j<data._global_dimensions; ++j) {
422 direction[j] = data._source_extra_data(global_neighbor_index, 3*data._global_dimensions + j);
430 for (
int j=0; j<data._basis_multiplier; ++j) {
434 local_direction[j] = data._pc.convertGlobalToLocalCoordinate(direction,j,*V);
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++){
443 alphaf = factorial[alphax]*factorial[alphay];
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;
453 double dot_product = 0.0;
456 dot_product = local_direction[0]*v0 + local_direction[1]*v1;
458 dot_product = direction[0]*v0 + direction[1]*v1;
463 *(delta+i) = dot_product * data._qm.getWeight(quadrature) * G_determinant;
465 *(delta+i) += dot_product * data._qm.getWeight(quadrature) * G_determinant;
471 entire_edge_length += G_determinant * data._qm.getWeight(quadrature);
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;
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);
496 double triangle_coords[3*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;
509 for (
int j=0; j<data._global_dimensions; ++j) {
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);
514 radius = std::sqrt(radius);
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);
527 num_vertices = num_vertices / data._global_dimensions;
528 auto T = triangle_coords_matrix;
532 double entire_cell_area = 0.0;
533 for (
int v=0; v<num_vertices; ++v) {
535 int v2 = (v+1) % num_vertices;
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);
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) {
557 unscaled_transformed_qp[j] += T(j,k)*data._qm.getSite(quadrature, k-1);
560 unscaled_transformed_qp[j] += T(j,0);
564 double G_determinant = 1.0;
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];
574 transformed_qp_norm = std::sqrt(transformed_qp_norm);
576 for (
int j=0; j<data._global_dimensions; ++j) {
577 scaled_transformed_qp[j] = unscaled_transformed_qp[j] * radius / transformed_qp_norm;
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));
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);
618 relative_coord[2] = 0;
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);
627 relative_coord[2] = 0;
632 "CellAverageSample does not support specific_order_only");
633 for (
int n = 0; n <= poly_order; n++){
634 for (alphay = 0; alphay <= 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;
645 entire_cell_area += G_determinant * data._qm.getWeight(quadrature);
650 for (
int n = 0; n <= poly_order; n++){
651 for (alphay = 0; alphay <= n; alphay++){
652 *(delta+k) /= entire_cell_area;