Compadre  1.5.5
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends 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
11  enum TargetOperation : int {
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 values of a cell using quadrature
43  //! Supported on 2D faces in 3D problems (manifold) and 2D cells in 2D problems
45  //! Integral values over cell using quadrature
46  //! Supported on 2D faces in 3D problems (manifold) and 2D cells in 2D problems
48  //! Average value of vector dotted with normal direction
49  //! Supported on 1D edges in 3D problems (2D-manifold) and 1D edges on 2D cells
51  //! Integral value of vector dotted with normal direction
52  //! Supported on 1D edges in 3D problems (2D-manifold) and 1D edges on 2D cells
54  //! Average value of vector dotted with tangent directions
55  //! Supported on 1D edges in 3D problems (2D-manifold) and 1D edges on 2D cells
57  //! Integral value of vector dotted with tangent directions
58  //! Supported on 1D edges in 3D problems (2D-manifold) and 1D edges on 2D cells
60  //! Should be the total count of all available target functionals
61  COUNT=20,
62  };
63 
64  //! Rank of target functional output for each TargetOperation
65  //! Rank of target functional input for each TargetOperation is based on the output
66  //! rank of the SamplingFunctional used on the polynomial basis
67  KOKKOS_INLINE_FUNCTION
68  int getTargetOutputTensorRank(const int& index) {
69  constexpr int TargetOutputTensorRank[] {
70  0, ///< PointEvaluation
71  1, ///< VectorPointEvaluation
72  0, ///< LaplacianOfScalarPointEvaluation
73  1, ///< VectorLaplacianPointEvaluation
74  1, ///< GradientOfScalarPointEvaluation
75  2, ///< GradientOfVectorPointEvaluation
76  0, ///< DivergenceOfVectorPointEvaluation
77  1, ///< CurlOfVectorPointEvaluation
78  1, ///< CurlCurlOfVectorPointEvaluation
79  0, ///< PartialXOfScalarPointEvaluation
80  0, ///< PartialYOfScalarPointEvaluation
81  0, ///< PartialZOfScalarPointEvaluation
82  0, ///< ChainedStaggeredLaplacianOfScalarPointEvaluation
83  0, ///< GaussianCurvaturePointEvaluation
84  0, ///< CellAverageEvaluation
85  0, ///< CellIntegralEvaluation
86  0, ///< FaceNormalAverageEvaluation
87  0, ///< FaceNormalIntegralEvaluation
88  0, ///< EdgeTangentAverageEvaluation
89  0, ///< EdgeTangentIntegralEvaluation
90  };
91  return TargetOutputTensorRank[index];
92  }
93 
94  //! Space in which to reconstruct polynomial
95  enum ReconstructionSpace : int {
96  //! Scalar polynomial basis centered at the target site and scaled by sum of basis powers
97  //! 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
98  //! \f$(x_t,y_t,z_t)\f$ is the coordinate of the target site in 3D coordinates.
100  //! Vector polynomial basis having # of components _dimensions, or (_dimensions-1) in the case of manifolds)
102  //! Scalar basis reused as many times as there are components in the vector
103  //! resulting in a much cheaper polynomial reconstruction
105  //! Divergence-free vector polynomial basis
107  //! Bernstein polynomial basis
109  };
110 
111  //! Number of actual components in the ReconstructionSpace
112  KOKKOS_INLINE_FUNCTION
113  int getActualReconstructionSpaceRank(const int& index) {
114  constexpr int ActualReconstructionSpaceRank[] = {
115  0, ///< ScalarTaylorPolynomial
116  1, ///< VectorTaylorPolynomial
117  0, ///< VectorOfScalarClonesTaylorPolynomial
118  0, ///< DivergenceFreeVectorTaylorPolynomial
119  0, ///< BernsteinPolynomial
120  };
121  return ActualReconstructionSpaceRank[index];
122  }
123 
124  //! Describes the SamplingFunction relationship to targets, neighbors
126  Identity, ///< No action performed on data before GMLS target operation
127  SameForAll, ///< Each neighbor for each target all apply the same transform
128  DifferentEachTarget, ///< Each target applies a different data transform, but the same to each neighbor
129  DifferentEachNeighbor, ///< Each target applies a different transform for each neighbor
130  };
131 
133  //! for uniqueness
134  size_t id;
135  //! Rank of sampling functional input for each SamplingFunctional
137  //! Rank of sampling functional output for each SamplingFunctional
139  //! Whether or not the SamplingTensor acts on the target site as well as the neighbors.
140  //! This makes sense only in staggered schemes, when each target site is also a source site
142  //! Whether the SamplingFunctional + ReconstructionSpace results in a nontrivial nullspace
144  //! Describes the SamplingFunction relationship to targets, neighbors
146 
147  //SamplingFunctional(std::string name, const int input_rank_, const int output_rank_,
148  KOKKOS_INLINE_FUNCTION
149  constexpr SamplingFunctional(const int input_rank_, const int output_rank_,
150  const bool use_target_site_weights_, const bool nontrivial_nullspace_,
151  const int transform_type_, const int id_) :
152  id(id_), input_rank(input_rank_), output_rank(output_rank_),
153  use_target_site_weights(use_target_site_weights_), nontrivial_nullspace(nontrivial_nullspace_),
154  transform_type(transform_type_) {}
155 
156  KOKKOS_INLINE_FUNCTION
157  constexpr bool operator == (const SamplingFunctional sf) const {
158  return id == sf.id;
159  }
160 
161  KOKKOS_INLINE_FUNCTION
162  constexpr bool operator != (const SamplingFunctional sf) const {
163  return id != sf.id;
164  }
165 
166  };
167 
168  //! Available sampling functionals
169  constexpr SamplingFunctional
170 
171  //! Point evaluations of the scalar source function
173 
174  //! Point evaluations of the entire vector source function
176 
177  //! Point evaluations of the entire vector source function
178  //! (but on a manifold, so it includes a transform into local coordinates)
180 
181  //! Analytical integral of a gradient source vector is just a difference of the scalar source at neighbor and target
183 
184  //! Samples consist of the result of integrals of a vector dotted with the tangent along edges between neighbor and target
186 
187  //! For integrating polynomial dotted with normal over an edge
189 
190  //! For integrating polynomial dotted with normal over an edge
192 
193  //! For polynomial dotted with normal on edge
195 
196  //! For integrating polynomial dotted with tangent over an edge
198 
199  //! For polynomial dotted with tangent
201 
202  //! For polynomial integrated on cells
204 
205  //! For polynomial integrated on cells
207 
208  //! Dense solver type
209  enum DenseSolverType : int {
210  //! QR+Pivoting factorization performed on P*sqrt(w) matrix
211  QR,
212  //! LU factorization performed on P^T*W*P matrix
213  LU,
214  };
215 
216  //! Problem type, that optionally can handle manifolds
217  enum ProblemType : int {
218  //! Standard GMLS problem type
220  //! Solve GMLS problem on a manifold (will use QR or SVD to solve the resultant GMLS
221  //! problem dependent on SamplingNontrivialNullspace
223  };
224 
225  //! Constraint type
226  enum ConstraintType : int {
227  //! No constraint
229  //! Neumann Gradient Scalar Type
231  };
232 
233  //! Available weighting kernel function types
240  };
241 
242  //! Coordinate type for input and output format of vector data on manifold problems.
243  //! Anything without a manifold is always Ambient.
244  enum CoordinatesType : int {
245  Ambient, ///< a 2D manifold in 3D in ambient coordinates would have 3 components for a vector
246  Local, ///< a 2D manifold in 3D in local coordinates would have 2 components for a vector
247  };
248 
249  /*
250  \page pageGMLSOperators GMLS Target Operations, Spaces, and Sampling Functionals
251  \section sectionGMLSTarget Target Operations
252 
253  \section sectionGMLSSpace Reconstruction Spaces
254 
255  \section sectionGMLSSampling Sampling Functionals
256  */
257 
258 
259 } // Compadre
260 
261 #endif
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.