Compadre  1.5.9
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Compadre_Operators.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Compadre: COMpatible PArticle Discretization and REmap Toolkit
4 //
5 // Copyright 2018 NTESS and the Compadre contributors.
6 // SPDX-License-Identifier: BSD-2-Clause
7 // *****************************************************************************
8 // @HEADER
9 #ifndef _COMPADRE_OPERATORS_HPP_
10 #define _COMPADRE_OPERATORS_HPP_
11 
12 #include "Compadre_Typedefs.hpp"
13 
14 #define make_sampling_functional(input, output, targets, nontrivial, transform) SamplingFunctional(input, output, targets, nontrivial, transform, __COUNTER__)
15 
16 namespace Compadre {
17 
18  //! Available target functionals
19  enum TargetOperation : int {
20  //! Point evaluation of a scalar
22  //! Point evaluation of a vector (reconstructs entire vector at once, requiring a
23  //! ReconstructionSpace having a sufficient number of components in the basis)
25  //! Point evaluation of the laplacian of a scalar (could be on a manifold or not)
27  //! Point evaluation of the laplacian of each component of a vector
29  //! Point evaluation of the gradient of a scalar
31  //! Point evaluation of the gradient of a vector (results in a matrix, NOT CURRENTLY IMPLEMENTED)
33  //! Point evaluation of the divergence of a vector (results in a scalar)
35  //! Point evaluation of the curl of a vector (results in a vector)
37  //! Point evaluation of the curl of a curl of a vector (results in a vector)
39  //! Point evaluation of the partial with respect to x of a scalar
41  //! Point evaluation of the partial with respect to y of a scalar
43  //! Point evaluation of the partial with respect to z of a scalar
45  //! Point evaluation of the chained staggered Laplacian acting on VectorTaylorPolynomial
46  //! basis + StaggeredEdgeIntegralSample sampling functional
48  //! Point evaluation of Gaussian curvature
50  //! Average values of a cell using quadrature
51  //! Supported on 2D faces in 3D problems (manifold) and 2D cells in 2D problems
53  //! Integral values over cell using quadrature
54  //! Supported on 2D faces in 3D problems (manifold) and 2D cells in 2D problems
56  //! Average value of vector dotted with normal direction
57  //! Supported on 1D edges in 3D problems (2D-manifold) and 1D edges on 2D cells
59  //! Integral value of vector dotted with normal direction
60  //! Supported on 1D edges in 3D problems (2D-manifold) and 1D edges on 2D cells
62  //! Average value of vector dotted with tangent directions
63  //! Supported on 1D edges in 3D problems (2D-manifold) and 1D edges on 2D cells
65  //! Integral value of vector dotted with tangent directions
66  //! Supported on 1D edges in 3D problems (2D-manifold) and 1D edges on 2D cells
68  //! Should be the total count of all available target functionals
69  COUNT=20,
70  };
71 
72  //! Rank of target functional output for each TargetOperation
73  //! Rank of target functional input for each TargetOperation is based on the output
74  //! rank of the SamplingFunctional used on the polynomial basis
75  KOKKOS_INLINE_FUNCTION
76  int getTargetOutputTensorRank(const int& index) {
77  constexpr int TargetOutputTensorRank[] {
78  0, ///< PointEvaluation
79  1, ///< VectorPointEvaluation
80  0, ///< LaplacianOfScalarPointEvaluation
81  1, ///< VectorLaplacianPointEvaluation
82  1, ///< GradientOfScalarPointEvaluation
83  2, ///< GradientOfVectorPointEvaluation
84  0, ///< DivergenceOfVectorPointEvaluation
85  1, ///< CurlOfVectorPointEvaluation
86  1, ///< CurlCurlOfVectorPointEvaluation
87  0, ///< PartialXOfScalarPointEvaluation
88  0, ///< PartialYOfScalarPointEvaluation
89  0, ///< PartialZOfScalarPointEvaluation
90  0, ///< ChainedStaggeredLaplacianOfScalarPointEvaluation
91  0, ///< GaussianCurvaturePointEvaluation
92  0, ///< CellAverageEvaluation
93  0, ///< CellIntegralEvaluation
94  0, ///< FaceNormalAverageEvaluation
95  0, ///< FaceNormalIntegralEvaluation
96  0, ///< EdgeTangentAverageEvaluation
97  0, ///< EdgeTangentIntegralEvaluation
98  };
99  return TargetOutputTensorRank[index];
100  }
101 
102  //! Space in which to reconstruct polynomial
103  enum ReconstructionSpace : int {
104  //! Scalar polynomial basis centered at the target site and scaled by sum of basis powers
105  //! 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
106  //! \f$(x_t,y_t,z_t)\f$ is the coordinate of the target site in 3D coordinates.
108  //! Vector polynomial basis having # of components _dimensions, or (_dimensions-1) in the case of manifolds)
110  //! Scalar basis reused as many times as there are components in the vector
111  //! resulting in a much cheaper polynomial reconstruction
113  //! Divergence-free vector polynomial basis
115  //! Bernstein polynomial basis
117  };
118 
119  //! Number of actual components in the ReconstructionSpace
120  KOKKOS_INLINE_FUNCTION
121  int getActualReconstructionSpaceRank(const int& index) {
122  constexpr int ActualReconstructionSpaceRank[] = {
123  0, ///< ScalarTaylorPolynomial
124  1, ///< VectorTaylorPolynomial
125  0, ///< VectorOfScalarClonesTaylorPolynomial
126  0, ///< DivergenceFreeVectorTaylorPolynomial
127  0, ///< BernsteinPolynomial
128  };
129  return ActualReconstructionSpaceRank[index];
130  }
131 
132  //! Describes the SamplingFunction relationship to targets, neighbors
134  Identity, ///< No action performed on data before GMLS target operation
135  SameForAll, ///< Each neighbor for each target all apply the same transform
136  DifferentEachTarget, ///< Each target applies a different data transform, but the same to each neighbor
137  DifferentEachNeighbor, ///< Each target applies a different transform for each neighbor
138  };
139 
141  //! for uniqueness
142  size_t id;
143  //! Rank of sampling functional input for each SamplingFunctional
145  //! Rank of sampling functional output for each SamplingFunctional
147  //! Whether or not the SamplingTensor acts on the target site as well as the neighbors.
148  //! This makes sense only in staggered schemes, when each target site is also a source site
150  //! Whether the SamplingFunctional + ReconstructionSpace results in a nontrivial nullspace
152  //! Describes the SamplingFunction relationship to targets, neighbors
154 
155  //SamplingFunctional(std::string name, const int input_rank_, const int output_rank_,
156  KOKKOS_INLINE_FUNCTION
157  constexpr SamplingFunctional(const int input_rank_, const int output_rank_,
158  const bool use_target_site_weights_, const bool nontrivial_nullspace_,
159  const int transform_type_, const int id_) :
160  id(id_), input_rank(input_rank_), output_rank(output_rank_),
161  use_target_site_weights(use_target_site_weights_), nontrivial_nullspace(nontrivial_nullspace_),
162  transform_type(transform_type_) {}
163 
164  KOKKOS_INLINE_FUNCTION
165  constexpr bool operator == (const SamplingFunctional sf) const {
166  return id == sf.id;
167  }
168 
169  KOKKOS_INLINE_FUNCTION
170  constexpr bool operator != (const SamplingFunctional sf) const {
171  return id != sf.id;
172  }
173 
174  };
175 
176  //! Available sampling functionals
177  constexpr SamplingFunctional
178 
179  //! Point evaluations of the scalar source function
181 
182  //! Point evaluations of the entire vector source function
184 
185  //! Point evaluations of the entire vector source function
186  //! (but on a manifold, so it includes a transform into local coordinates)
188 
189  //! Analytical integral of a gradient source vector is just a difference of the scalar source at neighbor and target
191 
192  //! Samples consist of the result of integrals of a vector dotted with the tangent along edges between neighbor and target
194 
195  //! For integrating polynomial dotted with normal over an edge
197 
198  //! For integrating polynomial dotted with normal over an edge
200 
201  //! For polynomial dotted with normal on edge
203 
204  //! For integrating polynomial dotted with tangent over an edge
206 
207  //! For polynomial dotted with tangent
209 
210  //! For polynomial integrated on cells
212 
213  //! For polynomial integrated on cells
215 
216  //! Dense solver type
217  enum DenseSolverType : int {
218  //! QR+Pivoting factorization performed on P*sqrt(w) matrix
219  QR,
220  //! LU factorization performed on P^T*W*P matrix
221  LU,
222  };
223 
224  //! Problem type, that optionally can handle manifolds
225  enum ProblemType : int {
226  //! Standard GMLS problem type
228  //! Solve GMLS problem on a manifold (will use QR or SVD to solve the resultant GMLS
229  //! problem dependent on SamplingNontrivialNullspace
231  };
232 
233  //! Constraint type
234  enum ConstraintType : int {
235  //! No constraint
237  //! Neumann Gradient Scalar Type
239  };
240 
241  //! Available weighting kernel function types
249  };
250 
251  //! Coordinate type for input and output format of vector data on manifold problems.
252  //! Anything without a manifold is always Ambient.
253  enum CoordinatesType : int {
254  Ambient, ///< a 2D manifold in 3D in ambient coordinates would have 3 components for a vector
255  Local, ///< a 2D manifold in 3D in local coordinates would have 2 components for a vector
256  };
257 
258  /*
259  \page pageGMLSOperators GMLS Target Operations, Spaces, and Sampling Functionals
260  \section sectionGMLSTarget Target Operations
261 
262  \section sectionGMLSSpace Reconstruction Spaces
263 
264  \section sectionGMLSSampling Sampling Functionals
265  */
266 
267 
268 } // Compadre
269 
270 #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.