Compadre
1.5.9
|
Namespaces | |
BernsteinPolynomialBasis | |
Definition of scalar Bernstein polynomial basis. | |
DivergenceFreePolynomialBasis | |
Definition of the divergence-free polynomial basis. | |
GMLS_LinearAlgebra | |
ScalarTaylorPolynomialBasis | |
Definition of scalar Taylor polynomial basis. | |
Classes | |
struct | SubviewND |
Creates 1D subviews of data from a 2D view, generally constructed with CreateNDSliceOnDeviceView. More... | |
class | Evaluator |
Lightweight Evaluator Helper This class is a lightweight wrapper for extracting and applying all relevant data from a GMLS class in order to transform data into a form that can be acted on by the GMLS operator, apply the action of the GMLS operator, and then transform data again (only if on a manifold) More... | |
struct | GMLSBasisData |
struct | GMLSSolutionData |
struct | ApplyTargets |
Functor to apply target evaluation to polynomial coefficients to store in _alphas. More... | |
struct | EvaluateStandardTargets |
Functor to evaluate targets operations on the basis. More... | |
struct | ComputePrestencilWeights |
Functor to calculate prestencil weights to apply to data to transform into a format expected by a GMLS stencil. More... | |
struct | AssembleStandardPsqrtW |
Functor to assemble the P*sqrt(weights) matrix and construct sqrt(weights)*Identity. More... | |
struct | ComputeCoarseTangentPlane |
Functor to create a coarse tangent approximation from a given neighborhood of points. More... | |
struct | AssembleCurvaturePsqrtW |
Functor to assemble the P*sqrt(weights) matrix and construct sqrt(weights)*Identity for curvature. More... | |
struct | GetAccurateTangentDirections |
Functor to evaluate curvature targets and construct accurate tangent direction approximation for manifolds. More... | |
struct | FixTangentDirectionOrdering |
Functor to determine if tangent directions need reordered, and to reorder them if needed We require that the normal is consistent with a right hand rule on the tangent vectors. More... | |
struct | ApplyCurvatureTargets |
Functor to evaluate curvature targets and apply to coefficients of curvature reconstruction. More... | |
struct | AssembleManifoldPsqrtW |
Functor to assemble the P*sqrt(weights) matrix and construct sqrt(weights)*Identity. More... | |
struct | EvaluateManifoldTargets |
Functor to evaluate targets on a manifold. More... | |
class | GMLS |
Generalized Moving Least Squares (GMLS) More... | |
class | KokkosParser |
Class handling Kokkos command line arguments and returning parameters. More... | |
struct | XYZ |
struct | NeighborLists |
NeighborLists assists in accessing entries of compressed row neighborhood lists. More... | |
struct | SamplingFunctional |
class | ParallelManager |
Parallel Manager. More... | |
class | RadiusResultSet |
Custom RadiusResultSet for nanoflann that uses pre-allocated space for indices and radii instead of using std::vec for std::pairs. More... | |
class | PointCloudSearch |
PointCloudSearch generates neighbor lists and window sizes for each target site. More... | |
struct | PointConnections |
Combines NeighborLists with the PointClouds from which it was derived Assumed that memory_space is the same as device, but it can be set to host, if desired. More... | |
class | Quadrature |
Quadrature. More... | |
struct | SolutionSet |
All vairables and functionality related to the layout and storage of GMLS solutions (alpha values) More... | |
struct | Extract |
Functions | |
template<typename SolutionData > | |
KOKKOS_INLINE_FUNCTION void | applyTargetsToCoefficients (const SolutionData &data, const member_type &teamMember, scratch_matrix_right_type Q, scratch_matrix_right_type P_target_row) |
For applying the evaluations from a target functional to the polynomial coefficients. More... | |
template<typename BasisData > | |
KOKKOS_INLINE_FUNCTION void | calcPij (const BasisData &data, const member_type &teamMember, double *delta, double *thread_workspace, const int target_index, int neighbor_index, const double alpha, const int dimension, const int poly_order, bool specific_order_only=false, const scratch_matrix_right_type *V=NULL, const ReconstructionSpace reconstruction_space=ReconstructionSpace::ScalarTaylorPolynomial, const SamplingFunctional polynomial_sampling_functional=PointSample, const int evaluation_site_local_index=0) |
Evaluates the polynomial basis under a particular sampling function. Generally used to fill a row of P. More... | |
template<typename BasisData > | |
KOKKOS_INLINE_FUNCTION void | calcGradientPij (const BasisData &data, const member_type &teamMember, double *delta, double *thread_workspace, const int target_index, int neighbor_index, const double alpha, const int partial_direction, const int dimension, const int poly_order, bool specific_order_only, const scratch_matrix_right_type *V, const ReconstructionSpace reconstruction_space, const SamplingFunctional polynomial_sampling_functional, const int evaluation_site_local_index=0) |
Evaluates the gradient of a polynomial basis under the Dirac Delta (pointwise) sampling function. More... | |
template<typename BasisData > | |
KOKKOS_INLINE_FUNCTION void | calcHessianPij (const BasisData &data, const member_type &teamMember, double *delta, double *thread_workspace, const int target_index, int neighbor_index, const double alpha, const int partial_direction_1, const int partial_direction_2, const int dimension, const int poly_order, bool specific_order_only, const scratch_matrix_right_type *V, const ReconstructionSpace reconstruction_space, const SamplingFunctional polynomial_sampling_functional, const int evaluation_site_local_index=0) |
Evaluates the Hessian of a polynomial basis under the Dirac Delta (pointwise) sampling function. More... | |
template<typename BasisData > | |
KOKKOS_INLINE_FUNCTION void | createWeightsAndP (const BasisData &data, const member_type &teamMember, scratch_vector_type delta, scratch_vector_type thread_workspace, scratch_matrix_right_type P, scratch_vector_type w, const int dimension, int polynomial_order, bool weight_p=false, scratch_matrix_right_type *V=NULL, const ReconstructionSpace reconstruction_space=ReconstructionSpace::ScalarTaylorPolynomial, const SamplingFunctional polynomial_sampling_functional=PointSample) |
Fills the _P matrix with either P or P*sqrt(w) More... | |
template<typename BasisData > | |
KOKKOS_INLINE_FUNCTION void | createWeightsAndPForCurvature (const BasisData &data, const member_type &teamMember, scratch_vector_type delta, scratch_vector_type thread_workspace, scratch_matrix_right_type P, scratch_vector_type w, const int dimension, bool only_specific_order, scratch_matrix_right_type *V=NULL) |
Fills the _P matrix with P*sqrt(w) for use in solving for curvature. More... | |
template<typename T , typename T2 > | |
struct SubviewND< T, T2, enable_if_t<(T::rank< 2)> >{T _data_in;T2 _data_original_view;bool _scalar_as_vector_if_needed;SubviewND(T data_in, T2 data_original_view, bool scalar_as_vector_if_needed){_data_in=data_in;_data_original_view=data_original_view;_scalar_as_vector_if_needed=scalar_as_vector_if_needed;}auto get1DView(const int column_num) -> | decltype (Kokkos::subview(_data_in, Kokkos::ALL)) |
Creates 1D subviews of data from a 1D view, generally constructed with CreateNDSliceOnDeviceView. More... | |
auto | get2DView (const int column_num, const int block_size) -> decltype(Kokkos::subview(_data_in, Kokkos::ALL)) |
T2 | copyToAndReturnOriginalView () |
template<typename T > | |
auto | CreateNDSliceOnDeviceView (T sampling_input_data_host_or_device, bool scalar_as_vector_if_needed) -> SubviewND< decltype(Kokkos::create_mirror_view(device_memory_space(), sampling_input_data_host_or_device)), T > |
Copies data_in to the device, and then allows for access to 1D columns of data on device. More... | |
KOKKOS_INLINE_FUNCTION double | MetricFactor (const scratch_vector_type a_, const double h, const double u1, const double u2) |
Metric factor (det(G)) at any point in the local chart. More... | |
KOKKOS_INLINE_FUNCTION double | GaussianCurvature (const scratch_vector_type a_, const double h, const double u1, const double u2) |
Gaussian curvature K at any point in the local chart. More... | |
KOKKOS_INLINE_FUNCTION double | SurfaceCurlOfScalar (const scratch_vector_type a_, const double h, const double u1, const double u2, int x_pow, int y_pow, const int component) |
Surface curl at any point in the local chart. More... | |
KOKKOS_INLINE_FUNCTION XYZ | operator+ (const XYZ &vecA, const XYZ &vecB) |
KOKKOS_INLINE_FUNCTION XYZ | operator- (const XYZ &vecA, const XYZ &vecB) |
KOKKOS_INLINE_FUNCTION XYZ | operator* (const XYZ &vecA, const XYZ &vecB) |
KOKKOS_INLINE_FUNCTION XYZ | operator+ (const XYZ &vecA, const scalar_type &constant) |
KOKKOS_INLINE_FUNCTION XYZ | operator+ (const scalar_type &constant, const XYZ &vecA) |
KOKKOS_INLINE_FUNCTION XYZ | operator- (const XYZ &vecA, const scalar_type &constant) |
KOKKOS_INLINE_FUNCTION XYZ | operator- (const scalar_type &constant, const XYZ &vecA) |
KOKKOS_INLINE_FUNCTION XYZ | operator* (const XYZ &vecA, const scalar_type &constant) |
KOKKOS_INLINE_FUNCTION XYZ | operator* (const scalar_type &constant, const XYZ &vecA) |
KOKKOS_INLINE_FUNCTION XYZ | operator/ (const XYZ &vecA, const scalar_type &constant) |
std::ostream & | operator<< (std::ostream &os, const XYZ &vec) |
KOKKOS_INLINE_FUNCTION int | pown (int n, unsigned p) |
n^p (n^0 returns 1, regardless of n) More... | |
KOKKOS_INLINE_FUNCTION int | getAdditionalAlphaSizeFromConstraint (DenseSolverType dense_solver_type, ConstraintType constraint_type) |
KOKKOS_INLINE_FUNCTION int | getAdditionalCoeffSizeFromConstraintAndSpace (DenseSolverType dense_solver_type, ConstraintType constraint_type, ReconstructionSpace reconstruction_space, const int dimension) |
KOKKOS_INLINE_FUNCTION void | getRHSDims (DenseSolverType dense_solver_type, ConstraintType constraint_type, ReconstructionSpace reconstruction_space, const int dimension, const int M, const int N, int &RHS_row, int &RHS_col) |
KOKKOS_INLINE_FUNCTION void | getPDims (DenseSolverType dense_solver_type, ConstraintType constraint_type, ReconstructionSpace reconstruction_space, const int dimension, const int M, const int N, int &out_row, int &out_col) |
KOKKOS_INLINE_FUNCTION int | getTargetOutputIndex (const int operation_num, const int output_component_axis_1, const int output_component_axis_2, const int dimensions) |
Helper function for finding alpha coefficients. More... | |
KOKKOS_INLINE_FUNCTION int | getSamplingOutputIndex (const SamplingFunctional sf, const int output_component_axis_1, const int output_component_axis_2) |
Helper function for finding alpha coefficients. More... | |
KOKKOS_INLINE_FUNCTION int | getInputRankOfSampling (SamplingFunctional sro) |
Input rank for sampling operation. More... | |
KOKKOS_INLINE_FUNCTION int | getOutputDimensionOfSampling (SamplingFunctional sro, const int local_dimensions) |
Dimensions ^ output rank for sampling operation (always in local chart if on a manifold, never ambient space) More... | |
KOKKOS_INLINE_FUNCTION int | getInputDimensionOfSampling (SamplingFunctional sro, const int global_dimensions) |
Dimensions ^ output rank for sampling operation (always in ambient space, never local chart on a manifold) More... | |
KOKKOS_INLINE_FUNCTION int | calculateBasisMultiplier (const ReconstructionSpace rs, const int local_dimensions) |
Calculate basis_multiplier. More... | |
KOKKOS_INLINE_FUNCTION int | calculateSamplingMultiplier (const ReconstructionSpace rs, const SamplingFunctional sro, const int local_dimensions) |
Calculate sampling_multiplier. More... | |
KOKKOS_INLINE_FUNCTION int | getOutputDimensionOfOperation (TargetOperation lro, const int local_dimensions) |
Dimensions ^ output rank for target operation. More... | |
KOKKOS_INLINE_FUNCTION int | getInputDimensionOfOperation (TargetOperation lro, SamplingFunctional sro, const int local_dimensions) |
Dimensions ^ input rank for target operation (always in local chart if on a manifold, never ambient space) More... | |
template<typename view_type > | |
NeighborLists< view_type > | CreateNeighborLists (view_type number_of_neighbors_list) |
CreateNeighborLists allows for the construction of an object of type NeighborLists with template deduction. More... | |
template<typename view_type > | |
NeighborLists< view_type > | CreateNeighborLists (view_type neighbor_lists, view_type number_of_neighbors_list) |
CreateNeighborLists allows for the construction of an object of type NeighborLists with template deduction. More... | |
template<typename view_type > | |
NeighborLists< view_type > | CreateNeighborLists (view_type neighbor_lists, view_type number_of_neighbors_list, view_type neighbor_lists_row_offsets) |
CreateNeighborLists allows for the construction of an object of type NeighborLists with template deduction. More... | |
template<typename view_type_2d , typename view_type_1d = Kokkos::View<int*, typename view_type_2d::memory_space, typename view_type_2d::memory_traits>> | |
NeighborLists< view_type_1d > | Convert2DToCompressedRowNeighborLists (view_type_2d neighbor_lists) |
Converts 2D neighbor lists to compressed row neighbor lists. More... | |
KOKKOS_INLINE_FUNCTION int | getTargetOutputTensorRank (const int &index) |
Rank of target functional output for each TargetOperation Rank of target functional input for each TargetOperation is based on the output rank of the SamplingFunctional used on the polynomial basis. More... | |
KOKKOS_INLINE_FUNCTION int | getActualReconstructionSpaceRank (const int &index) |
Number of actual components in the ReconstructionSpace. More... | |
template<typename view_type > | |
PointCloudSearch< view_type > | CreatePointCloudSearch (view_type src_view, const local_index_type dimensions=-1, const local_index_type max_leaf=-1) |
CreatePointCloudSearch allows for the construction of an object of type PointCloudSearch with template deduction. More... | |
template<typename TargetData > | |
KOKKOS_INLINE_FUNCTION void | computeTargetFunctionals (const TargetData &data, const member_type &teamMember, scratch_vector_type delta, scratch_vector_type thread_workspace, scratch_matrix_right_type P_target_row) |
Evaluates a polynomial basis with a target functional applied to each member of the basis. More... | |
template<typename TargetData > | |
KOKKOS_INLINE_FUNCTION void | computeCurvatureFunctionals (const TargetData &data, const member_type &teamMember, scratch_vector_type delta, scratch_vector_type thread_workspace, scratch_matrix_right_type P_target_row, const scratch_matrix_right_type *V, const local_index_type local_neighbor_index=-1) |
Evaluates a polynomial basis for the curvature with a gradient target functional applied. More... | |
template<typename TargetData > | |
KOKKOS_INLINE_FUNCTION void | computeTargetFunctionalsOnManifold (const TargetData &data, const member_type &teamMember, scratch_vector_type delta, scratch_vector_type thread_workspace, scratch_matrix_right_type P_target_row, scratch_matrix_right_type V, scratch_vector_type curvature_coefficients) |
Evaluates a polynomial basis with a target functional applied, using information from the manifold curvature. More... | |
KOKKOS_INLINE_FUNCTION void | getMidpointFromCellVertices (const member_type &teamMember, scratch_vector_type midpoint_storage, scratch_matrix_right_type cell_coordinates, const int cell_num, const int dim=3) |
template<typename view_type_1 , typename view_type_2 > | |
KOKKOS_INLINE_FUNCTION double | getAreaFromVectors (const member_type &teamMember, view_type_1 v1, view_type_2 v2) |
template<typename output_memory_space , typename view_type_input_data , typename output_array_layout = typename view_type_input_data::array_layout, typename index_type = int> | |
Kokkos::View< int *, output_array_layout, output_memory_space > | filterViewByID (view_type_input_data input_data_host_or_device, index_type filtered_value) |
KOKKOS_INLINE_FUNCTION void | evaluateConstraints (scratch_matrix_right_type M, scratch_matrix_right_type PsqrtW, const ConstraintType constraint_type, const ReconstructionSpace reconstruction_space, const int NP, const double cutoff_p, const int dimension, const int num_neighbors=0, scratch_matrix_right_type *T=NULL) |
Variables | |
constexpr SamplingFunctional | PointSample = make_sampling_functional(0,0,false,false,(int)Identity) |
Available sampling functionals. More... | |
constexpr SamplingFunctional | VectorPointSample = make_sampling_functional(1,1,false,false,(int)Identity) |
Point evaluations of the entire vector source function. More... | |
constexpr SamplingFunctional | ManifoldVectorPointSample = make_sampling_functional(1,1,false,false,(int)DifferentEachTarget) |
Point evaluations of the entire vector source function (but on a manifold, so it includes a transform into local coordinates) More... | |
constexpr SamplingFunctional | StaggeredEdgeAnalyticGradientIntegralSample = make_sampling_functional(0,0,true,true,(int)SameForAll) |
Analytical integral of a gradient source vector is just a difference of the scalar source at neighbor and target. More... | |
constexpr SamplingFunctional | StaggeredEdgeIntegralSample = make_sampling_functional(1,0,true,true,(int)DifferentEachNeighbor) |
Samples consist of the result of integrals of a vector dotted with the tangent along edges between neighbor and target. More... | |
constexpr SamplingFunctional | VaryingManifoldVectorPointSample = make_sampling_functional(1,1,false,false,(int)DifferentEachNeighbor) |
For integrating polynomial dotted with normal over an edge. More... | |
constexpr SamplingFunctional | FaceNormalIntegralSample = make_sampling_functional(1,0,false,false,(int)Identity) |
For integrating polynomial dotted with normal over an edge. More... | |
constexpr SamplingFunctional | FaceNormalAverageSample = make_sampling_functional(1,0,false,false,(int)Identity) |
For polynomial dotted with normal on edge. More... | |
constexpr SamplingFunctional | EdgeTangentIntegralSample = make_sampling_functional(1,0,false,false,(int)Identity) |
For integrating polynomial dotted with tangent over an edge. More... | |
constexpr SamplingFunctional | EdgeTangentAverageSample = make_sampling_functional(1,0,false,false,(int)Identity) |
For polynomial dotted with tangent. More... | |
constexpr SamplingFunctional | CellAverageSample = make_sampling_functional(0,0,false,false,(int)DifferentEachNeighbor) |
For polynomial integrated on cells. More... | |
constexpr SamplingFunctional | CellIntegralSample = make_sampling_functional(0,0,false,false,(int)DifferentEachNeighbor) |
For polynomial integrated on cells. More... | |
enum Compadre::ConstraintType : int |
Constraint type.
Enumerator | |
---|---|
NO_CONSTRAINT |
No constraint. |
NEUMANN_GRAD_SCALAR |
Neumann Gradient Scalar Type. |
Definition at line 234 of file Compadre_Operators.hpp.
enum Compadre::CoordinatesType : int |
Coordinate type for input and output format of vector data on manifold problems.
Anything without a manifold is always Ambient.
Enumerator | |
---|---|
Ambient |
a 2D manifold in 3D in ambient coordinates would have 3 components for a vector |
Local |
a 2D manifold in 3D in local coordinates would have 2 components for a vector |
Definition at line 253 of file Compadre_Operators.hpp.
enum Compadre::DenseSolverType : int |
Dense solver type.
Enumerator | |
---|---|
QR |
QR+Pivoting factorization performed on P*sqrt(w) matrix. |
LU |
LU factorization performed on P^T*W*P matrix. |
Definition at line 217 of file Compadre_Operators.hpp.
enum Compadre::ProblemType : int |
Problem type, that optionally can handle manifolds.
Enumerator | |
---|---|
STANDARD |
Standard GMLS problem type. |
MANIFOLD |
Solve GMLS problem on a manifold (will use QR or SVD to solve the resultant GMLS problem dependent on SamplingNontrivialNullspace. |
Definition at line 225 of file Compadre_Operators.hpp.
enum Compadre::QuadratureType : int |
Enumerator | |
---|---|
INVALID | |
LINE | |
TRI | |
QUAD | |
TET | |
HEX |
Definition at line 19 of file Compadre_Quadrature.hpp.
enum Compadre::ReconstructionSpace : int |
Space in which to reconstruct polynomial.
Definition at line 103 of file Compadre_Operators.hpp.
enum Compadre::SamplingTransformType : int |
Describes the SamplingFunction relationship to targets, neighbors.
Enumerator | |
---|---|
Identity |
No action performed on data before GMLS target operation. |
SameForAll |
Each neighbor for each target all apply the same transform. |
DifferentEachTarget |
Each target applies a different data transform, but the same to each neighbor. |
DifferentEachNeighbor |
Each target applies a different transform for each neighbor. |
Definition at line 133 of file Compadre_Operators.hpp.
enum Compadre::TargetOperation : int |
Available target functionals.
Definition at line 19 of file Compadre_Operators.hpp.
enum Compadre::WeightingFunctionType : int |
Available weighting kernel function types.
Enumerator | |
---|---|
Power | |
Gaussian | |
CubicSpline | |
CardinalCubicBSpline | |
Cosine | |
Sigmoid |
Definition at line 242 of file Compadre_Operators.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::applyTargetsToCoefficients | ( | const SolutionData & | data, |
const member_type & | teamMember, | ||
scratch_matrix_right_type | Q, | ||
scratch_matrix_right_type | P_target_row | ||
) |
For applying the evaluations from a target functional to the polynomial coefficients.
data | [out/in] - GMLSSolutionData struct (stores solution in data._d_ss._alphas) |
teamMember | [in] - Kokkos::TeamPolicy member type (created by parallel_for) |
Q | [in] - 2D Kokkos View containing the polynomial coefficients |
P_target_row | [in] - 1D Kokkos View where the evaluation of the polynomial basis is stored |
Definition at line 23 of file Compadre_ApplyTargetEvaluations.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::calcGradientPij | ( | const BasisData & | data, |
const member_type & | teamMember, | ||
double * | delta, | ||
double * | thread_workspace, | ||
const int | target_index, | ||
int | neighbor_index, | ||
const double | alpha, | ||
const int | partial_direction, | ||
const int | dimension, | ||
const int | poly_order, | ||
bool | specific_order_only, | ||
const scratch_matrix_right_type * | V, | ||
const ReconstructionSpace | reconstruction_space, | ||
const SamplingFunctional | polynomial_sampling_functional, | ||
const int | evaluation_site_local_index = 0 |
||
) |
Evaluates the gradient of a polynomial basis under the Dirac Delta (pointwise) sampling function.
data | [in] - GMLSBasisData struct |
teamMember | [in] - Kokkos::TeamPolicy member type (created by parallel_for) |
delta | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large is the _basis_multipler*the dimension of the polynomial basis. |
thread_workspace | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _poly_order*the spatial dimension of the polynomial basis. |
target_index | [in] - target number |
neighbor_index | [in] - index of neighbor for this target with respect to local numbering [0,...,number of neighbors for target] |
alpha | [in] - double to determine convex combination of target and neighbor site at which to evaluate polynomials. (1-alpha)*neighbor + alpha*target |
partial_direction | [in] - direction that partial is taken with respect to, e.g. 0 is x direction, 1 is y direction |
dimension | [in] - spatial dimension of basis to evaluate. e.g. dimension two basis of order one is 1, x, y, whereas for dimension 3 it is 1, x, y, z |
poly_order | [in] - polynomial basis degree |
specific_order_only | [in] - boolean for only evaluating one degree of polynomial when true |
V | [in] - orthonormal basis matrix size _dimensions * _dimensions whose first _dimensions-1 columns are an approximation of the tangent plane |
reconstruction_space | [in] - space of polynomial that a sampling functional is to evaluate |
sampling_strategy | [in] - sampling functional specification |
evaluation_site_local_index | [in] - local index for evaluation sites (0 is target site) |
Definition at line 681 of file Compadre_Basis.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::calcHessianPij | ( | const BasisData & | data, |
const member_type & | teamMember, | ||
double * | delta, | ||
double * | thread_workspace, | ||
const int | target_index, | ||
int | neighbor_index, | ||
const double | alpha, | ||
const int | partial_direction_1, | ||
const int | partial_direction_2, | ||
const int | dimension, | ||
const int | poly_order, | ||
bool | specific_order_only, | ||
const scratch_matrix_right_type * | V, | ||
const ReconstructionSpace | reconstruction_space, | ||
const SamplingFunctional | polynomial_sampling_functional, | ||
const int | evaluation_site_local_index = 0 |
||
) |
Evaluates the Hessian of a polynomial basis under the Dirac Delta (pointwise) sampling function.
data | [in] - GMLSBasisData struct |
teamMember | [in] - Kokkos::TeamPolicy member type (created by parallel_for) |
delta | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large is the _basis_multipler*the dimension of the polynomial basis. |
thread_workspace | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _poly_order*the spatial dimension of the polynomial basis. |
target_index | [in] - target number |
neighbor_index | [in] - index of neighbor for this target with respect to local numbering [0,...,number of neighbors for target] |
alpha | [in] - double to determine convex combination of target and neighbor site at which to evaluate polynomials. (1-alpha)*neighbor + alpha*target |
partial_direction_1 | [in] - first direction that partial is taken with respect to, e.g. 0 is x direction, 1 is y direction |
partial_direction_2 | [in] - second direction that partial is taken with respect to, e.g. 0 is x direction, 1 is y direction |
dimension | [in] - spatial dimension of basis to evaluate. e.g. dimension two basis of order one is 1, x, y, whereas for dimension 3 it is 1, x, y, z |
poly_order | [in] - polynomial basis degree |
specific_order_only | [in] - boolean for only evaluating one degree of polynomial when true |
V | [in] - orthonormal basis matrix size _dimensions * _dimensions whose first _dimensions-1 columns are an approximation of the tangent plane |
reconstruction_space | [in] - space of polynomial that a sampling functional is to evaluate |
sampling_strategy | [in] - sampling functional specification |
evaluation_site_local_index | [in] - local index for evaluation sites (0 is target site) |
Definition at line 769 of file Compadre_Basis.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::calcPij | ( | const BasisData & | data, |
const member_type & | teamMember, | ||
double * | delta, | ||
double * | thread_workspace, | ||
const int | target_index, | ||
int | neighbor_index, | ||
const double | alpha, | ||
const int | dimension, | ||
const int | poly_order, | ||
bool | specific_order_only = false , |
||
const scratch_matrix_right_type * | V = NULL , |
||
const ReconstructionSpace | reconstruction_space = ReconstructionSpace::ScalarTaylorPolynomial , |
||
const SamplingFunctional | polynomial_sampling_functional = PointSample , |
||
const int | evaluation_site_local_index = 0 |
||
) |
Evaluates the polynomial basis under a particular sampling function. Generally used to fill a row of P.
data | [in] - GMLSBasisData struct |
teamMember | [in] - Kokkos::TeamPolicy member type (created by parallel_for) |
delta | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _basis_multipler*the dimension of the polynomial basis. |
thread_workspace | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _poly_order*the spatial dimension of the polynomial basis. |
target_index | [in] - target number |
neighbor_index | [in] - index of neighbor for this target with respect to local numbering [0,...,number of neighbors for target] |
alpha | [in] - double to determine convex combination of target and neighbor site at which to evaluate polynomials. (1-alpha)*neighbor + alpha*target |
dimension | [in] - spatial dimension of basis to evaluate. e.g. dimension two basis of order one is 1, x, y, whereas for dimension 3 it is 1, x, y, z |
poly_order | [in] - polynomial basis degree |
specific_order_only | [in] - boolean for only evaluating one degree of polynomial when true |
V | [in] - orthonormal basis matrix size _dimensions * _dimensions whose first _dimensions-1 columns are an approximation of the tangent plane |
reconstruction_space | [in] - space of polynomial that a sampling functional is to evaluate |
sampling_strategy | [in] - sampling functional specification |
evaluation_site_local_index | [in] - local index for evaluation sites (0 is target site) |
Definition at line 34 of file Compadre_Basis.hpp.
KOKKOS_INLINE_FUNCTION int Compadre::calculateBasisMultiplier | ( | const ReconstructionSpace | rs, |
const int | local_dimensions | ||
) |
Calculate basis_multiplier.
Definition at line 268 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION int Compadre::calculateSamplingMultiplier | ( | const ReconstructionSpace | rs, |
const SamplingFunctional | sro, | ||
const int | local_dimensions | ||
) |
Calculate sampling_multiplier.
Definition at line 276 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::computeCurvatureFunctionals | ( | const TargetData & | data, |
const member_type & | teamMember, | ||
scratch_vector_type | delta, | ||
scratch_vector_type | thread_workspace, | ||
scratch_matrix_right_type | P_target_row, | ||
const scratch_matrix_right_type * | V, | ||
const local_index_type | local_neighbor_index = -1 |
||
) |
Evaluates a polynomial basis for the curvature with a gradient target functional applied.
data._operations is used by this function which is set through a modifier function
data | [in] - GMLSBasisData struct |
teamMember | [in] - Kokkos::TeamPolicy member type (created by parallel_for) |
delta | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large is the _basis_multipler*the dimension of the polynomial basis. |
thread_workspace | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _curvature_poly_order*the spatial dimension of the polynomial basis. |
P_target_row | [out] - 1D Kokkos View where the evaluation of the polynomial basis is stored |
V | [in] - orthonormal basis matrix size _dimensions * _dimensions whose first _dimensions-1 columns are an approximation of the tangent plane |
Definition at line 1055 of file Compadre_Targets.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::computeTargetFunctionals | ( | const TargetData & | data, |
const member_type & | teamMember, | ||
scratch_vector_type | delta, | ||
scratch_vector_type | thread_workspace, | ||
scratch_matrix_right_type | P_target_row | ||
) |
Evaluates a polynomial basis with a target functional applied to each member of the basis.
data | [in] - GMLSBasisData struct |
teamMember | [in] - Kokkos::TeamPolicy member type (created by parallel_for) |
delta | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large is the _basis_multipler*the dimension of the polynomial basis. |
thread_workspace | [in/out] - scratch space that is allocated so that each team has its own copy. Must be at least as large is the _poly_order*_global_dimensions. |
P_target_row | [out] - 1D Kokkos View where the evaluation of the polynomial basis is stored |
Definition at line 26 of file Compadre_Targets.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::computeTargetFunctionalsOnManifold | ( | const TargetData & | data, |
const member_type & | teamMember, | ||
scratch_vector_type | delta, | ||
scratch_vector_type | thread_workspace, | ||
scratch_matrix_right_type | P_target_row, | ||
scratch_matrix_right_type | V, | ||
scratch_vector_type | curvature_coefficients | ||
) |
Evaluates a polynomial basis with a target functional applied, using information from the manifold curvature.
data._operations is used by this function which is set through a modifier function
data | [in] - GMLSBasisData struct |
teamMember | [in] - Kokkos::TeamPolicy member type (created by parallel_for) |
delta | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large is the _basis_multipler*the dimension of the polynomial basis. |
thread_workspace | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _curvature_poly_order*the spatial dimension of the polynomial basis. |
P_target_row | [out] - 1D Kokkos View where the evaluation of the polynomial basis is stored |
V | [in] - orthonormal basis matrix size _dimensions * _dimensions whose first _dimensions-1 columns are an approximation of the tangent plane |
curvature_coefficients | [in] - polynomial coefficients for curvature |
Definition at line 1124 of file Compadre_Targets.hpp.
NeighborLists<view_type_1d> Compadre::Convert2DToCompressedRowNeighborLists | ( | view_type_2d | neighbor_lists | ) |
Converts 2D neighbor lists to compressed row neighbor lists.
Definition at line 332 of file Compadre_NeighborLists.hpp.
T2 Compadre::copyToAndReturnOriginalView | ( | ) |
Definition at line 95 of file Compadre_Evaluator.hpp.
auto Compadre::CreateNDSliceOnDeviceView | ( | T | sampling_input_data_host_or_device, |
bool | scalar_as_vector_if_needed | ||
) | -> SubviewND<decltype(Kokkos::create_mirror_view( device_memory_space(), sampling_input_data_host_or_device)), T> |
Copies data_in to the device, and then allows for access to 1D columns of data on device.
Handles either 2D or 1D views as input, and they can be on the host or the device.
Definition at line 106 of file Compadre_Evaluator.hpp.
NeighborLists<view_type> Compadre::CreateNeighborLists | ( | view_type | number_of_neighbors_list | ) |
CreateNeighborLists allows for the construction of an object of type NeighborLists with template deduction.
Definition at line 314 of file Compadre_NeighborLists.hpp.
NeighborLists<view_type> Compadre::CreateNeighborLists | ( | view_type | neighbor_lists, |
view_type | number_of_neighbors_list | ||
) |
CreateNeighborLists allows for the construction of an object of type NeighborLists with template deduction.
Definition at line 320 of file Compadre_NeighborLists.hpp.
NeighborLists<view_type> Compadre::CreateNeighborLists | ( | view_type | neighbor_lists, |
view_type | number_of_neighbors_list, | ||
view_type | neighbor_lists_row_offsets | ||
) |
CreateNeighborLists allows for the construction of an object of type NeighborLists with template deduction.
Definition at line 326 of file Compadre_NeighborLists.hpp.
PointCloudSearch<view_type> Compadre::CreatePointCloudSearch | ( | view_type | src_view, |
const local_index_type | dimensions = -1 , |
||
const local_index_type | max_leaf = -1 |
||
) |
CreatePointCloudSearch allows for the construction of an object of type PointCloudSearch with template deduction.
Definition at line 796 of file Compadre_PointCloudSearch.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::createWeightsAndP | ( | const BasisData & | data, |
const member_type & | teamMember, | ||
scratch_vector_type | delta, | ||
scratch_vector_type | thread_workspace, | ||
scratch_matrix_right_type | P, | ||
scratch_vector_type | w, | ||
const int | dimension, | ||
int | polynomial_order, | ||
bool | weight_p = false , |
||
scratch_matrix_right_type * | V = NULL , |
||
const ReconstructionSpace | reconstruction_space = ReconstructionSpace::ScalarTaylorPolynomial , |
||
const SamplingFunctional | polynomial_sampling_functional = PointSample |
||
) |
Fills the _P matrix with either P or P*sqrt(w)
data | [in] - GMLSBasisData struct |
teamMember | [in] - Kokkos::TeamPolicy member type (created by parallel_for) |
delta | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large is the _basis_multipler*the dimension of the polynomial basis. |
thread_workspace | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _poly_order*the spatial dimension of the polynomial basis. |
P | [out] - 2D Kokkos View which will contain evaluation of sampling functional on polynomial basis for each neighbor the target has (stored column major) |
w | [out] - 1D Kokkos View which will contain weighting kernel values for the target with each neighbor if weight_p = true |
dimension | [in] - spatial dimension of basis to evaluate. e.g. dimension two basis of order one is 1, x, y, whereas for dimension 3 it is 1, x, y, z |
polynomial_order | [in] - polynomial basis degree |
weight_p | [in] - boolean whether to fill w with kernel weights |
V | [in] - orthonormal basis matrix size _dimensions * _dimensions whose first _dimensions-1 columns are an approximation of the tangent plane |
reconstruction_space | [in] - space of polynomial that a sampling functional is to evaluate |
sampling_strategy | [in] - sampling functional specification |
Definition at line 850 of file Compadre_Basis.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::createWeightsAndPForCurvature | ( | const BasisData & | data, |
const member_type & | teamMember, | ||
scratch_vector_type | delta, | ||
scratch_vector_type | thread_workspace, | ||
scratch_matrix_right_type | P, | ||
scratch_vector_type | w, | ||
const int | dimension, | ||
bool | only_specific_order, | ||
scratch_matrix_right_type * | V = NULL |
||
) |
Fills the _P matrix with P*sqrt(w) for use in solving for curvature.
Uses _curvature_poly_order as the polynomial order of the basis \param data [in] - GMLSBasisData struct \param teamMember [in] - Kokkos::TeamPolicy member type (created by parallel_for) \param delta [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large is the
s_multipler*the dimension of the polynomial basis.
thread_workspace | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _order*the spatial dimension of the polynomial basis. |
P | [out] - 2D Kokkos View which will contain evaluation of sampling functional on polynomial basis for each neighbor the target has (stored column major) |
w | [out] - 1D Kokkos View which will contain weighting kernel values for the target with each neighbor if weight_p = true |
dimension | [in] - spatial dimension of basis to evaluate. e.g. dimension two basis of order one is 1, x, y, whereas for dimension 3 it is 1, x, y, z |
only_specific_order | [in] - boolean for only evaluating one degree of polynomial when true |
V | [in] - orthonormal basis matrix size _dimensions * _dimensions whose first _dimensions-1 columns are an approximation of the tangent plane |
Definition at line 947 of file Compadre_Basis.hpp.
struct SubviewND< T, T2, enable_if_t<(T::rank< 2)> >{T _data_in;T2 _data_original_view;bool _scalar_as_vector_if_needed;SubviewND(T data_in, T2 data_original_view, bool scalar_as_vector_if_needed){_data_in=data_in;_data_original_view=data_original_view;_scalar_as_vector_if_needed=scalar_as_vector_if_needed;}auto get1DView(const int column_num) -> Compadre::decltype | ( | Kokkos:: | subview_data_in, Kokkos::ALL | ) |
Creates 1D subviews of data from a 1D view, generally constructed with CreateNDSliceOnDeviceView.
Definition at line 80 of file Compadre_Evaluator.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::evaluateConstraints | ( | scratch_matrix_right_type | M, |
scratch_matrix_right_type | PsqrtW, | ||
const ConstraintType | constraint_type, | ||
const ReconstructionSpace | reconstruction_space, | ||
const int | NP, | ||
const double | cutoff_p, | ||
const int | dimension, | ||
const int | num_neighbors = 0 , |
||
scratch_matrix_right_type * | T = NULL |
||
) |
Definition at line 17 of file Compadre_CreateConstraints.hpp.
Kokkos::View<int*, output_array_layout, output_memory_space> Compadre::filterViewByID | ( | view_type_input_data | input_data_host_or_device, |
index_type | filtered_value | ||
) |
Definition at line 55 of file Compadre_Utilities.hpp.
KOKKOS_INLINE_FUNCTION double Compadre::GaussianCurvature | ( | const scratch_vector_type | a_, |
const double | h, | ||
const double | u1, | ||
const double | u2 | ||
) |
Gaussian curvature K at any point in the local chart.
Definition at line 58 of file Compadre_Manifold_Functions.hpp.
auto Compadre::get2DView | ( | const int | column_num, |
const int | block_size | ||
) | -> decltype(Kokkos::subview(_data_in, Kokkos::ALL)) |
Definition at line 90 of file Compadre_Evaluator.hpp.
KOKKOS_INLINE_FUNCTION int Compadre::getActualReconstructionSpaceRank | ( | const int & | index | ) |
Number of actual components in the ReconstructionSpace.
< ScalarTaylorPolynomial
< VectorTaylorPolynomial
< VectorOfScalarClonesTaylorPolynomial
< DivergenceFreeVectorTaylorPolynomial
< BernsteinPolynomial
Definition at line 121 of file Compadre_Operators.hpp.
KOKKOS_INLINE_FUNCTION int Compadre::getAdditionalAlphaSizeFromConstraint | ( | DenseSolverType | dense_solver_type, |
ConstraintType | constraint_type | ||
) |
Definition at line 169 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION int Compadre::getAdditionalCoeffSizeFromConstraintAndSpace | ( | DenseSolverType | dense_solver_type, |
ConstraintType | constraint_type, | ||
ReconstructionSpace | reconstruction_space, | ||
const int | dimension | ||
) |
Definition at line 179 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION double Compadre::getAreaFromVectors | ( | const member_type & | teamMember, |
view_type_1 | v1, | ||
view_type_2 | v2 | ||
) |
Definition at line 32 of file Compadre_Utilities.hpp.
KOKKOS_INLINE_FUNCTION int Compadre::getInputDimensionOfOperation | ( | TargetOperation | lro, |
SamplingFunctional | sro, | ||
const int | local_dimensions | ||
) |
Dimensions ^ input rank for target operation (always in local chart if on a manifold, never ambient space)
Definition at line 298 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION int Compadre::getInputDimensionOfSampling | ( | SamplingFunctional | sro, |
const int | global_dimensions | ||
) |
Dimensions ^ output rank for sampling operation (always in ambient space, never local chart on a manifold)
Definition at line 262 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION int Compadre::getInputRankOfSampling | ( | SamplingFunctional | sro | ) |
Input rank for sampling operation.
Definition at line 248 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::getMidpointFromCellVertices | ( | const member_type & | teamMember, |
scratch_vector_type | midpoint_storage, | ||
scratch_matrix_right_type | cell_coordinates, | ||
const int | cell_num, | ||
const int | dim = 3 |
||
) |
Definition at line 18 of file Compadre_Utilities.hpp.
KOKKOS_INLINE_FUNCTION int Compadre::getOutputDimensionOfOperation | ( | TargetOperation | lro, |
const int | local_dimensions | ||
) |
Dimensions ^ output rank for target operation.
Definition at line 292 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION int Compadre::getOutputDimensionOfSampling | ( | SamplingFunctional | sro, |
const int | local_dimensions | ||
) |
Dimensions ^ output rank for sampling operation (always in local chart if on a manifold, never ambient space)
Definition at line 255 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::getPDims | ( | DenseSolverType | dense_solver_type, |
ConstraintType | constraint_type, | ||
ReconstructionSpace | reconstruction_space, | ||
const int | dimension, | ||
const int | M, | ||
const int | N, | ||
int & | out_row, | ||
int & | out_col | ||
) |
Definition at line 210 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::getRHSDims | ( | DenseSolverType | dense_solver_type, |
ConstraintType | constraint_type, | ||
ReconstructionSpace | reconstruction_space, | ||
const int | dimension, | ||
const int | M, | ||
const int | N, | ||
int & | RHS_row, | ||
int & | RHS_col | ||
) |
Definition at line 190 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION int Compadre::getSamplingOutputIndex | ( | const SamplingFunctional | sf, |
const int | output_component_axis_1, | ||
const int | output_component_axis_2 | ||
) |
Helper function for finding alpha coefficients.
Definition at line 241 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION int Compadre::getTargetOutputIndex | ( | const int | operation_num, |
const int | output_component_axis_1, | ||
const int | output_component_axis_2, | ||
const int | dimensions | ||
) |
Helper function for finding alpha coefficients.
Definition at line 234 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION int Compadre::getTargetOutputTensorRank | ( | const int & | index | ) |
Rank of target functional output for each TargetOperation Rank of target functional input for each TargetOperation is based on the output rank of the SamplingFunctional used on the polynomial basis.
< PointEvaluation
< VectorPointEvaluation
< LaplacianOfScalarPointEvaluation
< VectorLaplacianPointEvaluation
< GradientOfScalarPointEvaluation
< GradientOfVectorPointEvaluation
< DivergenceOfVectorPointEvaluation
< CurlOfVectorPointEvaluation
< CurlCurlOfVectorPointEvaluation
< PartialXOfScalarPointEvaluation
< PartialYOfScalarPointEvaluation
< PartialZOfScalarPointEvaluation
< ChainedStaggeredLaplacianOfScalarPointEvaluation
< GaussianCurvaturePointEvaluation
< CellAverageEvaluation
< CellIntegralEvaluation
< FaceNormalAverageEvaluation
< FaceNormalIntegralEvaluation
< EdgeTangentAverageEvaluation
< EdgeTangentIntegralEvaluation
Definition at line 76 of file Compadre_Operators.hpp.
KOKKOS_INLINE_FUNCTION double Compadre::MetricFactor | ( | const scratch_vector_type | a_, |
const double | h, | ||
const double | u1, | ||
const double | u2 | ||
) |
Metric factor (det(G)) at any point in the local chart.
Definition at line 14 of file Compadre_Manifold_Functions.hpp.
KOKKOS_INLINE_FUNCTION XYZ Compadre::operator* | ( | const XYZ & | vecA, |
const XYZ & | vecB | ||
) |
Definition at line 113 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION XYZ Compadre::operator* | ( | const XYZ & | vecA, |
const scalar_type & | constant | ||
) |
Definition at line 133 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION XYZ Compadre::operator* | ( | const scalar_type & | constant, |
const XYZ & | vecA | ||
) |
Definition at line 137 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION XYZ Compadre::operator+ | ( | const XYZ & | vecA, |
const XYZ & | vecB | ||
) |
Definition at line 105 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION XYZ Compadre::operator+ | ( | const XYZ & | vecA, |
const scalar_type & | constant | ||
) |
Definition at line 117 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION XYZ Compadre::operator+ | ( | const scalar_type & | constant, |
const XYZ & | vecA | ||
) |
Definition at line 121 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION XYZ Compadre::operator- | ( | const XYZ & | vecA, |
const XYZ & | vecB | ||
) |
Definition at line 109 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION XYZ Compadre::operator- | ( | const XYZ & | vecA, |
const scalar_type & | constant | ||
) |
Definition at line 125 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION XYZ Compadre::operator- | ( | const scalar_type & | constant, |
const XYZ & | vecA | ||
) |
Definition at line 129 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION XYZ Compadre::operator/ | ( | const XYZ & | vecA, |
const scalar_type & | constant | ||
) |
Definition at line 141 of file Compadre_Misc.hpp.
|
inline |
Definition at line 144 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION int Compadre::pown | ( | int | n, |
unsigned | p | ||
) |
n^p (n^0 returns 1, regardless of n)
Definition at line 149 of file Compadre_Misc.hpp.
KOKKOS_INLINE_FUNCTION double Compadre::SurfaceCurlOfScalar | ( | const scratch_vector_type | a_, |
const double | h, | ||
const double | u1, | ||
const double | u2, | ||
int | x_pow, | ||
int | y_pow, | ||
const int | component | ||
) |
Surface curl at any point in the local chart.
Definition at line 104 of file Compadre_Manifold_Functions.hpp.
constexpr SamplingFunctional Compadre::CellAverageSample = make_sampling_functional(0,0,false,false,(int)DifferentEachNeighbor) |
For polynomial integrated on cells.
Definition at line 211 of file Compadre_Operators.hpp.
constexpr SamplingFunctional Compadre::CellIntegralSample = make_sampling_functional(0,0,false,false,(int)DifferentEachNeighbor) |
For polynomial integrated on cells.
Definition at line 214 of file Compadre_Operators.hpp.
constexpr SamplingFunctional Compadre::EdgeTangentAverageSample = make_sampling_functional(1,0,false,false,(int)Identity) |
For polynomial dotted with tangent.
Definition at line 208 of file Compadre_Operators.hpp.
constexpr SamplingFunctional Compadre::EdgeTangentIntegralSample = make_sampling_functional(1,0,false,false,(int)Identity) |
For integrating polynomial dotted with tangent over an edge.
Definition at line 205 of file Compadre_Operators.hpp.
constexpr SamplingFunctional Compadre::FaceNormalAverageSample = make_sampling_functional(1,0,false,false,(int)Identity) |
For polynomial dotted with normal on edge.
Definition at line 202 of file Compadre_Operators.hpp.
constexpr SamplingFunctional Compadre::FaceNormalIntegralSample = make_sampling_functional(1,0,false,false,(int)Identity) |
For integrating polynomial dotted with normal over an edge.
Definition at line 199 of file Compadre_Operators.hpp.
constexpr SamplingFunctional Compadre::ManifoldVectorPointSample = make_sampling_functional(1,1,false,false,(int)DifferentEachTarget) |
Point evaluations of the entire vector source function (but on a manifold, so it includes a transform into local coordinates)
Definition at line 187 of file Compadre_Operators.hpp.
constexpr SamplingFunctional Compadre::PointSample = make_sampling_functional(0,0,false,false,(int)Identity) |
Available sampling functionals.
Point evaluations of the scalar source function
Definition at line 180 of file Compadre_Operators.hpp.
constexpr SamplingFunctional Compadre::StaggeredEdgeAnalyticGradientIntegralSample = make_sampling_functional(0,0,true,true,(int)SameForAll) |
Analytical integral of a gradient source vector is just a difference of the scalar source at neighbor and target.
Definition at line 190 of file Compadre_Operators.hpp.
constexpr SamplingFunctional Compadre::StaggeredEdgeIntegralSample = make_sampling_functional(1,0,true,true,(int)DifferentEachNeighbor) |
Samples consist of the result of integrals of a vector dotted with the tangent along edges between neighbor and target.
Definition at line 193 of file Compadre_Operators.hpp.
constexpr SamplingFunctional Compadre::VaryingManifoldVectorPointSample = make_sampling_functional(1,1,false,false,(int)DifferentEachNeighbor) |
For integrating polynomial dotted with normal over an edge.
Definition at line 196 of file Compadre_Operators.hpp.
constexpr SamplingFunctional Compadre::VectorPointSample = make_sampling_functional(1,1,false,false,(int)Identity) |
Point evaluations of the entire vector source function.
Definition at line 183 of file Compadre_Operators.hpp.