9 #ifndef _COMPADRE_OPERATORS_HPP_
10 #define _COMPADRE_OPERATORS_HPP_
14 #define make_sampling_functional(input, output, targets, nontrivial, transform) SamplingFunctional(input, output, targets, nontrivial, transform, __COUNTER__)
75 KOKKOS_INLINE_FUNCTION
77 constexpr
int TargetOutputTensorRank[] {
99 return TargetOutputTensorRank[index];
120 KOKKOS_INLINE_FUNCTION
122 constexpr
int ActualReconstructionSpaceRank[] = {
129 return ActualReconstructionSpaceRank[index];
156 KOKKOS_INLINE_FUNCTION
158 const bool use_target_site_weights_,
const bool nontrivial_nullspace_,
159 const int transform_type_,
const int id_) :
164 KOKKOS_INLINE_FUNCTION
169 KOKKOS_INLINE_FUNCTION
177 constexpr SamplingFunctional
constexpr SamplingFunctional FaceNormalIntegralSample
For integrating polynomial dotted with normal over an edge.
Each neighbor for each target all apply the same transform.
constexpr SamplingFunctional EdgeTangentAverageSample
For polynomial dotted with tangent.
Average value of vector dotted with tangent directions Supported on 1D edges in 3D problems (2D-manif...
a 2D manifold in 3D in local coordinates would have 2 components for a vector
Point evaluation of Gaussian curvature.
ConstraintType
Constraint type.
int output_rank
Rank of sampling functional output for each SamplingFunctional.
Neumann Gradient Scalar Type.
Point evaluation of a scalar.
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...
constexpr SamplingFunctional ManifoldVectorPointSample
Point evaluations of the entire vector source function (but on a manifold, so it includes a transform...
Point evaluation of the gradient of a scalar.
SamplingTransformType
Describes the SamplingFunction relationship to targets, neighbors.
WeightingFunctionType
Available weighting kernel function types.
Point evaluation of the laplacian of a scalar (could be on a manifold or not)
Point evaluation of the partial with respect to y of a scalar.
KOKKOS_INLINE_FUNCTION constexpr bool operator==(const SamplingFunctional sf) const
Point evaluation of the partial with respect to z of a scalar.
Integral value of vector dotted with tangent directions Supported on 1D edges in 3D problems (2D-mani...
constexpr SamplingFunctional VaryingManifoldVectorPointSample
For integrating polynomial dotted with normal over an edge.
Point evaluation of the curl of a vector (results in a vector)
a 2D manifold in 3D in ambient coordinates would have 3 components for a vector
KOKKOS_INLINE_FUNCTION int getActualReconstructionSpaceRank(const int &index)
Number of actual components in the ReconstructionSpace.
DenseSolverType
Dense solver type.
bool use_target_site_weights
Whether or not the SamplingTensor acts on the target site as well as the neighbors.
Vector polynomial basis having # of components _dimensions, or (_dimensions-1) in the case of manifol...
Point evaluation of the curl of a curl of a vector (results in a vector)
constexpr SamplingFunctional StaggeredEdgeAnalyticGradientIntegralSample
Analytical integral of a gradient source vector is just a difference of the scalar source at neighbor...
Average value of vector dotted with normal direction Supported on 1D edges in 3D problems (2D-manifol...
Should be the total count of all available target functionals.
TargetOperation
Available target functionals.
Each target applies a different transform for each neighbor.
Solve GMLS problem on a manifold (will use QR or SVD to solve the resultant GMLS problem dependent on...
int input_rank
Rank of sampling functional input for each SamplingFunctional.
Point evaluation of the partial with respect to x of a scalar.
constexpr SamplingFunctional StaggeredEdgeIntegralSample
Samples consist of the result of integrals of a vector dotted with the tangent along edges between ne...
LU factorization performed on P^T*W*P matrix.
Standard GMLS problem type.
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_)
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.
Point evaluation of the divergence of a vector (results in a scalar)
QR+Pivoting factorization performed on P*sqrt(w) matrix.
Average values of a cell using quadrature Supported on 2D faces in 3D problems (manifold) and 2D cell...
KOKKOS_INLINE_FUNCTION constexpr bool operator!=(const SamplingFunctional sf) const
CoordinatesType
Coordinate type for input and output format of vector data on manifold problems.
constexpr SamplingFunctional EdgeTangentIntegralSample
For integrating polynomial dotted with tangent over an edge.
ProblemType
Problem type, that optionally can handle manifolds.
constexpr SamplingFunctional CellAverageSample
For polynomial integrated on cells.
No action performed on data before GMLS target operation.
Point evaluation of the laplacian of each component of a vector.
KOKKOS_INLINE_FUNCTION int getTargetOutputTensorRank(const int &index)
Rank of target functional output for each TargetOperation Rank of target functional input for each Ta...
Divergence-free vector polynomial basis.
int transform_type
Describes the SamplingFunction relationship to targets, neighbors.
constexpr SamplingFunctional PointSample
Available sampling functionals.
constexpr SamplingFunctional CellIntegralSample
For polynomial integrated on cells.
Each target applies a different data transform, but the same to each neighbor.
bool nontrivial_nullspace
Whether the SamplingFunctional + ReconstructionSpace results in a nontrivial nullspace.
Integral value of vector dotted with normal direction Supported on 1D edges in 3D problems (2D-manifo...
#define make_sampling_functional(input, output, targets, nontrivial, transform)
Integral values over cell using quadrature Supported on 2D faces in 3D problems (manifold) and 2D cel...
constexpr SamplingFunctional VectorPointSample
Point evaluations of the entire vector source function.
Point evaluation of the gradient of a vector (results in a matrix, NOT CURRENTLY IMPLEMENTED) ...
Point evaluation of the chained staggered Laplacian acting on VectorTaylorPolynomial basis + Staggere...
Point evaluation of a vector (reconstructs entire vector at once, requiring a ReconstructionSpace hav...
Bernstein polynomial basis.