1 #ifndef _COMPADRE_OPERATORS_HPP_
2 #define _COMPADRE_OPERATORS_HPP_
6 #define make_sampling_functional(input, output, targets, nontrivial, transform) SamplingFunctional(input, output, targets, nontrivial, transform, __COUNTER__)
117 KOKKOS_INLINE_FUNCTION
119 const bool use_target_site_weights_,
const bool nontrivial_nullspace_,
120 const int transform_type_,
const int id_) :
125 KOKKOS_INLINE_FUNCTION
130 KOKKOS_INLINE_FUNCTION
138 constexpr SamplingFunctional
Divergence-free vector polynomial basis.
Point evaluation of the curl of a curl of a vector (results in a vector)
constexpr SamplingFunctional FaceNormalIntegralSample
For integrating polynomial dotted with normal over an edge.
Standard GMLS problem type.
constexpr SamplingFunctional ScalarFaceAverageSample
For polynomial integrated on faces.
int output_rank
Rank of sampling functional output for each SamplingFunctional.
Point evaluation of the curl of a vector (results in a vector)
Point evaluation of the laplacian of each component of a vector.
Point evaluation of the partial with respect to y of a scalar.
constexpr SamplingFunctional ManifoldVectorPointSample
Point evaluations of the entire vector source function (but on a manifold, so it includes a transform...
Point evaluation of a vector (reconstructs entire vector at once, requiring a ReconstructionSpace hav...
a 2D manifold in 3D in local coordinates would have 2 components for a vector
Neumann Gradient Scalar Type.
SamplingTransformType
Describes the SamplingFunction relationship to targets, neighbors.
CoordinatesType
Coordinate type for input and output format of vector data on manifold problems.
KOKKOS_INLINE_FUNCTION constexpr bool operator==(const SamplingFunctional sf) const
constexpr SamplingFunctional VaryingManifoldVectorPointSample
For integrating polynomial dotted with normal over an edge.
Each target applies a different data transform, but the same to each neighbor.
Each target applies a different transform for each neighbor.
bool use_target_site_weights
Whether or not the SamplingTensor acts on the target site as well as the neighbors.
constexpr SamplingFunctional StaggeredEdgeAnalyticGradientIntegralSample
Analytical integral of a gradient source vector is just a difference of the scalar source at neighbor...
Scalar polynomial basis centered at the target site and scaled by sum of basis powers e...
Solve GMLS problem on a manifold (will use QR or SVD to solve the resultant GMLS problem dependent on...
Should be the total count of all available target functionals.
Point evaluation of Gaussian curvature.
constexpr int TargetOutputTensorRank[]
Rank of target functional output for each TargetOperation Rank of target functional input for each Ta...
constexpr SamplingFunctional FaceTangentIntegralSample
For integrating polynomial dotted with tangent over an edge.
QR factorization performed on P*sqrt(w) matrix.
ProblemType
Problem type, that optionally can handle manifolds.
constexpr SamplingFunctional FaceNormalPointSample
For polynomial dotted with normal on edge.
Point evaluation of a scalar.
int input_rank
Rank of sampling functional input for each SamplingFunctional.
constexpr SamplingFunctional StaggeredEdgeIntegralSample
Samples consist of the result of integrals of a vector dotted with the tangent along edges between ne...
Scalar basis reused as many times as there are components in the vector resulting in a much cheaper p...
SVD factorization performed on P*sqrt(w) matrix.
ReconstructionSpace
Space in which to reconstruct polynomial.
KOKKOS_INLINE_FUNCTION constexpr SamplingFunctional(const int input_rank_, const int output_rank_, const bool use_target_site_weights_, const bool nontrivial_nullspace_, const int transform_type_, const int id_)
Each neighbor for each target all apply the same transform.
WeightingFunctionType
Available weighting kernel function types.
KOKKOS_INLINE_FUNCTION constexpr bool operator!=(const SamplingFunctional sf) const
Point evaluation of the laplacian of a scalar (could be on a manifold or not)
Point evaluation of the chained staggered Laplacian acting on VectorTaylorPolynomial basis + Staggere...
Point evaluation of the divergence of a vector (results in a scalar)
Point evaluation of the gradient of a scalar.
constexpr SamplingFunctional FaceTangentPointSample
For polynomial dotted with tangent.
No action performed on data before GMLS target operation.
DenseSolverType
Dense solver type.
Point evaluation of the partial with respect to z of a scalar.
Point evaluation of the partial with respect to x of a scalar.
LU factorization performed on P^T*W*P matrix.
int transform_type
Describes the SamplingFunction relationship to targets, neighbors.
constexpr int ActualReconstructionSpaceRank[]
Number of actual components in the ReconstructionSpace.
constexpr SamplingFunctional PointSample
Available sampling functionals.
bool nontrivial_nullspace
Whether the SamplingFunctional + ReconstructionSpace results in a nontrivial nullspace requiring SVD...
a 2D manifold in 3D in ambient coordinates would have 3 components for a vector
Point evaluation of the gradient of a vector (results in a matrix, NOT CURRENTLY IMPLEMENTED) ...
#define make_sampling_functional(input, output, targets, nontrivial, transform)
Vector polynomial basis having # of components _dimensions, or (_dimensions-1) in the case of manifol...
constexpr SamplingFunctional VectorPointSample
Point evaluations of the entire vector source function.
TargetOperation
Available target functionals.
ConstraintType
Constraint type.
Average of values in a face of a cell using quadrature 2D in 3D problem, 1D in 2D problem...