1 #ifndef _COMPADRE_GMLS_BASIS_HPP_
2 #define _COMPADRE_GMLS_BASIS_HPP_
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 {
14 const int my_num_neighbors = this->
getNNeighbors(target_index);
17 const double factorial[15] = {1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600, 6227020800, 87178291200};
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) {
27 component = -(neighbor_index+1);
31 if (neighbor_index > -1) {
33 for (
int i=0; i<dimension; ++i) {
38 }
else if (additional_evaluation_local_index > 0) {
40 for (
int i=0; i<dimension; ++i) {
46 for (
int i=0; i<3; ++i) relative_coord[i] = 0;
50 if ((polynomial_sampling_functional ==
PointSample ||
56 double cutoff_p =
_epsilons(target_index);
57 const int start_index = specific_order_only ? poly_order : 0;
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;
71 for (
int d=0; d<dimension; ++d) {
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);
75 for (
int n=0; n<dimension_offset; ++n) {
76 *(delta+d*dimension_offset+n) = 0;
83 double cutoff_p =
_epsilons(target_index);
89 double cutoff_p =
_epsilons(target_index);
90 const int start_index = specific_order_only ? poly_order : 0;
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);
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);
98 Kokkos::single(Kokkos::PerThread(teamMember), [&] () {
100 double cutoff_p =
_epsilons(target_index);
103 const int start_index = specific_order_only ? poly_order : 0;
109 XYZ tangent_quadrature_coord_2d;
110 XYZ quadrature_coord_2d;
111 for (
int j=0; j<dimension; ++j) {
119 for (
int n = start_index; n <= poly_order; n++){
120 for (alphay = 0; alphay <= n; alphay++){
122 alphaf = factorial[alphax]*factorial[alphay];
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;
131 double dot_product = tangent_quadrature_coord_2d[0]*v0 + tangent_quadrature_coord_2d[1]*v1;
146 double cutoff_p =
_epsilons(target_index);
147 int alphax, alphay, alphaz;
149 const int start_index = specific_order_only ? poly_order : 0;
155 XYZ quadrature_coord_3d;
156 XYZ tangent_quadrature_coord_3d;
157 for (
int j=0; j<dimension; ++j) {
165 for (
int n = start_index; n <= poly_order; n++) {
166 if (dimension == 3) {
167 for (alphaz = 0; alphaz <= n; alphaz++){
169 for (alphay = 0; alphay <= s; alphay++){
171 alphaf = factorial[alphax]*factorial[alphay]*factorial[alphaz];
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;
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;
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;
199 double dot_product = tangent_quadrature_coord_3d[0]*v0 + tangent_quadrature_coord_3d[1]*v1 + tangent_quadrature_coord_3d[2]*v2;
202 if (quadrature == 0) {
210 }
else if (dimension == 2) {
211 for (alphay = 0; alphay <= n; alphay++){
213 alphaf = factorial[alphax]*factorial[alphay];
220 v1 = std::pow(quadrature_coord_3d.
x/cutoff_p,alphax)
221 *std::pow(quadrature_coord_3d.
y/cutoff_p,alphay)/alphaf;
224 v0 = std::pow(quadrature_coord_3d.
x/cutoff_p,alphax)
225 *std::pow(quadrature_coord_3d.
y/cutoff_p,alphay)/alphaf;
230 double dot_product = tangent_quadrature_coord_3d[0]*v0 + tangent_quadrature_coord_3d[1]*v1;
233 if (quadrature == 0) {
251 double cutoff_p =
_epsilons(target_index);
256 int neighbor_index_in_source =
getNeighborIndex(target_index, neighbor_index);
271 XYZ endpoints_difference;
272 for (
int j=0; j<dimension; ++j) {
279 const int start_index = specific_order_only ? poly_order : 0;
282 for (
int quadrature = 0; quadrature<quadrature_point_loop; ++quadrature) {
287 XYZ quadrature_coord_2d;
288 for (
int j=0; j<dimension; ++j) {
298 quadrature_coord_2d[j] = relative_coord[j];
314 for (
int n = start_index; n <= poly_order; n++){
315 for (alphay = 0; alphay <= n; alphay++){
317 alphaf = factorial[alphax]*factorial[alphay];
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;
327 double dot_product = direction_2d[0]*v0 + direction_2d[1]*v1;
334 *(delta+i) = dot_product *
_qm.
getWeight(quadrature) * magnitude;
337 *(delta+i) = dot_product;
341 *(delta+i) += dot_product *
_qm.
getWeight(quadrature) * magnitude;
349 auto global_neighbor_index =
getNeighborIndex(target_index, neighbor_index);
350 double cutoff_p =
_epsilons(target_index);
351 int alphax, alphay, alphaz;
356 double triangle_coords[3*3];
364 triangle_coords_matrix(j, 0) = midpoint(j);
368 double reference_cell_area = 0.5;
369 double entire_cell_area = 0.0;
370 auto T=triangle_coords_matrix;
372 for (
size_t v=0; v<num_vertices; ++v) {
374 int v2 = (v+1) % num_vertices;
380 Kokkos::subview(T, Kokkos::ALL(), 1), Kokkos::subview(T, Kokkos::ALL(), 2));
385 for (
size_t v=0; v<num_vertices; ++v) {
387 int v2 = (v+1) % num_vertices;
398 double transformed_qp[3] = {0,0,0};
400 for (
int k=1; k<3; ++k) {
401 transformed_qp[j] += T(j,k)*
_qm.
getSite(quadrature, k-1);
403 transformed_qp[j] += T(j,0);
408 Kokkos::subview(T, Kokkos::ALL(), 1), Kokkos::subview(T, Kokkos::ALL(), 2));
409 double scaling_factor = sub_cell_area / reference_cell_area;
412 XYZ qp =
XYZ(transformed_qp[0], transformed_qp[1], transformed_qp[2]);
413 for (
int j=0; j<2; ++j) {
415 relative_coord[2] = 0;
418 for (
int j=0; j<dimension; ++j) {
421 for (
int j=dimension; j<3; ++j) {
422 relative_coord[j] = 0.0;
427 const int start_index = specific_order_only ? poly_order : 0;
428 if (dimension == 3) {
429 for (
int n = start_index; n <= poly_order; n++){
430 for (alphaz = 0; alphaz <= n; alphaz++){
432 for (alphay = 0; alphay <= 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;
445 }
else if (dimension == 2) {
446 for (
int n = start_index; n <= poly_order; n++){
447 for (alphay = 0; alphay <= 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;
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 {
474 const int my_num_neighbors = this->
getNNeighbors(target_index);
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) {
484 component = -(neighbor_index+1);
492 if (neighbor_index > -1) {
493 for (
int i=0; i<dimension; ++i) {
498 }
else if (additional_evaluation_index > 0) {
499 for (
int i=0; i<dimension; ++i) {
504 for (
int i=0; i<3; ++i) relative_coord[i] = 0;
507 double cutoff_p =
_epsilons(target_index);
508 const int start_index = specific_order_only ? poly_order : 0;
510 if ((polynomial_sampling_functional ==
PointSample ||
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);
521 double cutoff_p =
_epsilons(target_index);
523 DivergenceFreePolynomialBasis::evaluatePartialDerivative(teamMember, delta, thread_workspace, dimension, poly_order, component, partial_direction, cutoff_p, relative_coord.
x, relative_coord.
y, relative_coord.
z);
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 {
537 const int my_num_neighbors = this->
getNNeighbors(target_index);
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) {
547 component = -(neighbor_index+1);
555 if (neighbor_index > -1) {
556 for (
int i=0; i<dimension; ++i) {
561 }
else if (additional_evaluation_index > 0) {
562 for (
int i=0; i<dimension; ++i) {
567 for (
int i=0; i<3; ++i) relative_coord[i] = 0;
570 double cutoff_p =
_epsilons(target_index);
572 if ((polynomial_sampling_functional ==
PointSample ||
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);
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);
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 {
603 const int my_num_neighbors = this->
getNNeighbors(target_index);
605 teamMember.team_barrier();
606 Kokkos::parallel_for(Kokkos::TeamThreadRange(teamMember,this->
getNNeighbors(target_index)),
616 double alpha_weight = 1;
632 this->
calcPij(teamMember, delta.data(), thread_workspace.data(), target_index, i + d*my_num_neighbors, 0 , dimension, polynomial_order,
false , V, reconstruction_space, polynomial_sampling_functional);
635 int storage_size = this->
getNP(polynomial_order, dimension, reconstruction_space);
639 for (
int j = 0; j < storage_size; ++j) {
641 P(i+my_num_neighbors*d, j) = delta[j] * std::sqrt(w(i+my_num_neighbors*d));
646 for (
int j = 0; j < storage_size; ++j) {
648 P(i+my_num_neighbors*d, j) = delta[j];
656 teamMember.team_barrier();
666 KOKKOS_INLINE_FUNCTION
676 teamMember.team_barrier();
677 Kokkos::parallel_for(Kokkos::TeamThreadRange(teamMember,this->
getNNeighbors(target_index)),
691 if (only_specific_order) {
693 this->
calcPij(teamMember, delta.data(), thread_workspace.data(), target_index, i, 0 , dimension, 1,
true );
701 for (
int j = 0; j < storage_size; ++j) {
702 P(i, j) = delta[j] * std::sqrt(w(i));
706 teamMember.team_barrier();
#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)