Compadre  1.2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Compadre_Operators.hpp
Go to the documentation of this file.
1 #ifndef _COMPADRE_OPERATORS_HPP_
2 #define _COMPADRE_OPERATORS_HPP_
3 
4 #include "Compadre_Typedefs.hpp"
5 
6 #define make_sampling_functional(input, output, targets, nontrivial, transform) SamplingFunctional(input, output, targets, nontrivial, transform, __COUNTER__)
7 
8 namespace Compadre {
9 
10  //! Available target functionals
12  //! Point evaluation of a scalar
14  //! Point evaluation of a vector (reconstructs entire vector at once, requiring a
15  //! ReconstructionSpace having a sufficient number of components in the basis)
17  //! Point evaluation of the laplacian of a scalar (could be on a manifold or not)
19  //! Point evaluation of the laplacian of each component of a vector
21  //! Point evaluation of the gradient of a scalar
23  //! Point evaluation of the gradient of a vector (results in a matrix, NOT CURRENTLY IMPLEMENTED)
25  //! Point evaluation of the divergence of a vector (results in a scalar)
27  //! Point evaluation of the curl of a vector (results in a vector)
29  //! Point evaluation of the curl of a curl of a vector (results in a vector)
31  //! Point evaluation of the partial with respect to x of a scalar
33  //! Point evaluation of the partial with respect to y of a scalar
35  //! Point evaluation of the partial with respect to z of a scalar
37  //! Point evaluation of the chained staggered Laplacian acting on VectorTaylorPolynomial
38  //! basis + StaggeredEdgeIntegralSample sampling functional
40  //! Point evaluation of Gaussian curvature
42  //! Average of values in a face of a cell using quadrature
43  //! 2D in 3D problem, 1D in 2D problem
45  //! Should be the total count of all available target functionals
46  COUNT=15,
47  };
48 
49  //! Rank of target functional output for each TargetOperation
50  //! Rank of target functional input for each TargetOperation is based on the output
51  //! rank of the SamplingFunctional used on the polynomial basis
52  constexpr int TargetOutputTensorRank[] {
53  0, ///< PointEvaluation
54  1, ///< VectorPointEvaluation
55  0, ///< LaplacianOfScalarPointEvaluation
56  1, ///< VectorLaplacianPointEvaluation
57  1, ///< GradientOfScalarPointEvaluation
58  2, ///< GradientOfVectorPointEvaluation
59  0, ///< DivergenceOfVectorPointEvaluation
60  1, ///< CurlOfVectorPointEvaluation
61  1, ///< CurlCurlOfVectorPointEvaluation
62  0, ///< PartialXOfScalarPointEvaluation
63  0, ///< PartialYOfScalarPointEvaluation
64  0, ///< PartialZOfScalarPointEvaluation
65  0, ///< ChainedStaggeredLaplacianOfScalarPointEvaluation
66  0, ///< GaussianCurvaturePointEvaluation
67  0, ///< ScalarFaceAverageEvaluation
68  };
69 
70  //! Space in which to reconstruct polynomial
72  //! Scalar polynomial basis centered at the target site and scaled by sum of basis powers
73  //! e.g. \f$(x-x_t)^2*(y-y_t)*(z-z_t)^3/factorial(2+1+3)\f$ would be a member of 3rd order in 3D, where
74  //! \f$(x_t,y_t,z_t)\f$ is the coordinate of the target site in 3D coordinates.
76  //! Vector polynomial basis having # of components _dimensions, or (_dimensions-1) in the case of manifolds)
78  //! Scalar basis reused as many times as there are components in the vector
79  //! resulting in a much cheaper polynomial reconstruction
81  //! Divergence-free vector polynomial basis
83  };
84 
85  //! Number of actual components in the ReconstructionSpace
86  constexpr int ActualReconstructionSpaceRank[] = {
87  0, ///< ScalarTaylorPolynomial
88  1, ///< VectorTaylorPolynomial
89  0, ///< VectorOfScalarClonesTaylorPolynomial
90  0, ///< DivergenceFreeVectorTaylorPolynomial
91  };
92 
93  //! Describes the SamplingFunction relationship to targets, neighbors
95  Identity, ///< No action performed on data before GMLS target operation
96  SameForAll, ///< Each neighbor for each target all apply the same transform
97  DifferentEachTarget, ///< Each target applies a different data transform, but the same to each neighbor
98  DifferentEachNeighbor, ///< Each target applies a different transform for each neighbor
99  };
100 
102  //! for uniqueness
103  size_t id;
104  //! Rank of sampling functional input for each SamplingFunctional
106  //! Rank of sampling functional output for each SamplingFunctional
108  //! Whether or not the SamplingTensor acts on the target site as well as the neighbors.
109  //! This makes sense only in staggered schemes, when each target site is also a source site
111  //! Whether the SamplingFunctional + ReconstructionSpace results in a nontrivial nullspace requiring SVD
113  //! Describes the SamplingFunction relationship to targets, neighbors
115 
116  //SamplingFunctional(std::string name, const int input_rank_, const int output_rank_,
117  KOKKOS_INLINE_FUNCTION
118  constexpr SamplingFunctional(const int input_rank_, const int output_rank_,
119  const bool use_target_site_weights_, const bool nontrivial_nullspace_,
120  const int transform_type_, const int id_) :
121  id(id_), input_rank(input_rank_), output_rank(output_rank_),
122  use_target_site_weights(use_target_site_weights_), nontrivial_nullspace(nontrivial_nullspace_),
123  transform_type(transform_type_) {}
124 
125  KOKKOS_INLINE_FUNCTION
126  constexpr bool operator == (const SamplingFunctional sf) const {
127  return id == sf.id;
128  }
129 
130  KOKKOS_INLINE_FUNCTION
131  constexpr bool operator != (const SamplingFunctional sf) const {
132  return id != sf.id;
133  }
134 
135  };
136 
137  //! Available sampling functionals
138  constexpr SamplingFunctional
139 
140  //! Point evaluations of the scalar source function
142 
143  //! Point evaluations of the entire vector source function
145 
146  //! Point evaluations of the entire vector source function
147  //! (but on a manifold, so it includes a transform into local coordinates)
149 
150  //! Analytical integral of a gradient source vector is just a difference of the scalar source at neighbor and target
152 
153  //! Samples consist of the result of integrals of a vector dotted with the tangent along edges between neighbor and target
155 
156  //! For integrating polynomial dotted with normal over an edge
158 
159  //! For integrating polynomial dotted with normal over an edge
161 
162  //! For polynomial dotted with normal on edge
164 
165  //! For integrating polynomial dotted with tangent over an edge
167 
168  //! For polynomial dotted with tangent
170 
171  //! For polynomial integrated on faces
173 
174  //! Dense solver type
176  //! QR factorization performed on P*sqrt(w) matrix
177  QR,
178  //! SVD factorization performed on P*sqrt(w) matrix
179  SVD,
180  //! LU factorization performed on P^T*W*P matrix
181  LU,
182  };
183 
184  //! Problem type, that optionally can handle manifolds
185  enum ProblemType {
186  //! Standard GMLS problem type
188  //! Solve GMLS problem on a manifold (will use QR or SVD to solve the resultant GMLS
189  //! problem dependent on SamplingNontrivialNullspace
191  };
192 
193  //! Constraint type
195  //! No constraint
197  //! Neumann Gradient Scalar Type
199  };
200 
201  //! Available weighting kernel function types
206  };
207 
208  //! Coordinate type for input and output format of vector data on manifold problems.
209  //! Anything without a manifold is always Ambient.
211  Ambient, ///< a 2D manifold in 3D in ambient coordinates would have 3 components for a vector
212  Local, ///< a 2D manifold in 3D in local coordinates would have 2 components for a vector
213  };
214 
215  /*
216  \page pageGMLSOperators GMLS Target Operations, Spaces, and Sampling Functionals
217  \section sectionGMLSTarget Target Operations
218 
219  \section sectionGMLSSpace Reconstruction Spaces
220 
221  \section sectionGMLSSampling Sampling Functionals
222  */
223 
224 
225 } // Compadre
226 
227 #endif
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...