Compadre  1.5.7
Public Types | Private Attributes | Friends | List of all members
Compadre::GMLS Class Reference

Generalized Moving Least Squares (GMLS) More...

Detailed Description

Generalized Moving Least Squares (GMLS)

This class sets up a batch of GMLS problems from a given set of neighbor lists, target sites, and source sites. GMLS requires a target functional, reconstruction space, and sampling functional to be specified. For a given choice of reconstruction space and sampling functional, multiple targets can be generated with very little additional computation, which is why this class allows for multiple target functionals to be specified.

Definition at line 32 of file Compadre_GMLS.hpp.

#include <Compadre_GMLS.hpp>

Public Types

typedef PointConnections< Kokkos::View< double **, layout_right >, Kokkos::View< double **, layout_right >, NeighborLists< Kokkos::View< int * > > > point_connections_type
 
typedef NeighborLists< Kokkos::View< int * > > neighbor_lists_type
 
typedef Kokkos::View< double **, layout_rightcoordinates_type
 

Public Member Functions

Instantiation / Destruction
 GMLS (ReconstructionSpace reconstruction_space, const SamplingFunctional polynomial_sampling_strategy, const SamplingFunctional data_sampling_strategy, const int poly_order, const int dimensions, const DenseSolverType dense_solver_type, const ProblemType problem_type, const ConstraintType constraint_type, const int manifold_curvature_poly_order)
 Maximal constructor, not intended for users. More...
 
 GMLS (ReconstructionSpace reconstruction_space, const SamplingFunctional polynomial_sampling_strategy, const SamplingFunctional data_sampling_strategy, const int poly_order, const int dimensions=3, const std::string dense_solver_type=std::string("QR"), const std::string problem_type=std::string("STANDARD"), const std::string constraint_type=std::string("NO_CONSTRAINT"), const int manifold_curvature_poly_order=2)
 Maximal constructor, but with string arguments. More...
 
 GMLS (const int poly_order, const int dimensions=3, const std::string dense_solver_type=std::string("QR"), const std::string problem_type=std::string("STANDARD"), const std::string constraint_type=std::string("NO_CONSTRAINT"), const int manifold_curvature_poly_order=2)
 Constructor for the case when the data sampling functional does not match the polynomial sampling functional. More...
 
 GMLS (ReconstructionSpace reconstruction_space, SamplingFunctional dual_sampling_strategy, const int poly_order, const int dimensions=3, const std::string dense_solver_type=std::string("QR"), const std::string problem_type=std::string("STANDARD"), const std::string constraint_type=std::string("NO_CONSTRAINT"), const int manifold_curvature_poly_order=2)
 Constructor for the case when nonstandard sampling functionals or reconstruction spaces are to be used. More...
 
Accessors

Retrieve member variables through public member functions

host_managed_local_index_type getPolynomialCoefficientsDomainRangeSize () const
 Returns (size of the basis used in instance's polynomial reconstruction) x (data input dimension) More...
 
int getPolynomialCoefficientsSize () const
 Returns size of the basis used in instance's polynomial reconstruction. More...
 
host_managed_local_index_type getPolynomialCoefficientsMemorySize () const
 Returns 2D array size in memory on which coefficients are stored. More...
 
int getDimensions () const
 Dimension of the GMLS problem, set only at class instantiation. More...
 
int getGlobalDimensions () const
 Dimension of the GMLS problem's point data (spatial description of points in ambient space), set only at class instantiation. More...
 
int getLocalDimensions () const
 Local dimension of the GMLS problem (less than global dimension if on a manifold), set only at class instantiation. More...
 
DenseSolverType getDenseSolverType () const
 Get dense solver type. More...
 
ProblemType getProblemType () const
 Get problem type. More...
 
ConstraintType getConstraintType () const
 Get constraint type. More...
 
int getPolynomialOrder () const
 Get basis order used for reconstruction. More...
 
int getCurvaturePolynomialOrder () const
 Get basis order used for curvature reconstruction. More...
 
WeightingFunctionType getWeightingType () const
 Type for weighting kernel for GMLS problem. More...
 
WeightingFunctionType getManifoldWeightingType () const
 Type for weighting kernel for curvature. More...
 
int getWeightingParameter (const int index=0) const
 Get parameter for weighting kernel for GMLS problem. More...
 
int getManifoldWeightingParameter (const int index=0) const
 Get parameter for weighting kernel for curvature. More...
 
int getNumberOfQuadraturePoints () const
 Number of quadrature points. More...
 
int getOrderOfQuadraturePoints () const
 Order of quadrature points. More...
 
int getDimensionOfQuadraturePoints () const
 Dimensions of quadrature points. More...
 
std::string getQuadratureType () const
 Type of quadrature points. More...
 
neighbor_lists_typegetNeighborLists () const
 Get neighbor list accessor. More...
 
decltype(_pc) * getPointConnections ()
 Get a view (device) of all point connection info. More...
 
neighbor_lists_typegetAdditionalEvaluationIndices () const
 (OPTIONAL) Get additional evaluation sites neighbor list-like accessor More...
 
decltype(_additional_pc) * getAdditionalPointConnections ()
 (OPTIONAL) Get a view (device) of all additional evaluation point connection info More...
 
decltype(_epsilons) * getWindowSizes ()
 Get a view (device) of all window sizes. More...
 
decltype(_T) * getTangentDirections ()
 Get a view (device) of all tangent direction bundles. More...
 
decltype(_ref_N) * getReferenceNormalDirections ()
 Get a view (device) of all reference outward normal directions. More...
 
double getTangentBundle (const int target_index, const int direction, const int component) const
 Get component of tangent or normal directions for manifold problems. More...
 
double getReferenceNormalDirection (const int target_index, const int component) const
 Get component of tangent or normal directions for manifold problems. More...
 
decltype(_prestencil_weightsgetPrestencilWeights () const
 Get a view (device) of all rank 2 preprocessing tensors This is a rank 5 tensor that is able to provide data transformation into a form that GMLS is able to operate on. More...
 
decltype(_RHSgetFullPolynomialCoefficientsBasis () const
 Get a view (device) of all polynomial coefficients basis. More...
 
SamplingFunctional getPolynomialSamplingFunctional () const
 Get the polynomial sampling functional specified at instantiation. More...
 
SamplingFunctional getDataSamplingFunctional () const
 Get the data sampling functional specified at instantiation (often the same as the polynomial sampling functional) More...
 
ReconstructionSpace getReconstructionSpace () const
 Get the reconstruction space specified at instantiation. More...
 
double getPreStencilWeight (SamplingFunctional sro, const int target_index, const int neighbor_index, bool for_target, const int output_component=0, const int input_component=0) const
 Returns a stencil to transform data from its existing state into the input expected for some sampling functionals. More...
 
decltype(_h_ss) * getSolutionSetHost (bool alpha_validity_check=true)
 Get solution set on host. More...
 
decltype(_d_ss) * getSolutionSetDevice (bool alpha_validity_check=true)
 Get solution set on device. More...
 
bool containsValidAlphas () const
 Check if GMLS solution set contains valid alpha values (has generateAlphas been called) More...
 
const GMLSSolutionData extractSolutionData () const
 Get GMLS solution data. More...
 
const GMLSBasisData extractBasisData () const
 Get GMLS basis data. More...
 
Modifiers

Changed member variables through public member functions

void resetCoefficientData ()
 
template<typename view_type_1 , typename view_type_2 , typename view_type_3 , typename view_type_4 >
void setProblemData (view_type_1 neighbor_lists, view_type_2 source_coordinates, view_type_3 target_coordinates, view_type_4 epsilons)
 Sets basic problem data (neighbor lists, source coordinates, and target coordinates) More...
 
template<typename view_type_1 , typename view_type_2 , typename view_type_3 , typename view_type_4 >
void setProblemData (view_type_1 cr_neighbor_lists, view_type_1 number_of_neighbors_list, view_type_2 source_coordinates, view_type_3 target_coordinates, view_type_4 epsilons)
 Sets basic problem data (neighbor lists data, number of neighbors list, source coordinates, and target coordinates) More...
 
template<typename view_type_1 , typename view_type_2 >
void setAdditionalEvaluationSitesData (view_type_1 additional_evaluation_indices, view_type_2 additional_evaluation_coordinates)
 (OPTIONAL) Sets additional evaluation sites for each target site More...
 
template<typename view_type_1 , typename view_type_2 >
void setAdditionalEvaluationSitesData (view_type_1 cr_additional_evaluation_indices, view_type_1 number_of_additional_evaluation_indices, view_type_2 additional_evaluation_coordinates)
 (OPTIONAL) Sets additional evaluation sites for each target site More...
 
template<typename view_type >
std::enable_if< view_type::rank==1 &&std::is_same< neighbor_lists_type::internal_view_type, view_type >::value==1, void >::type setNeighborLists (view_type neighbor_lists, view_type number_of_neighbors_list)
 Sets neighbor list information from compressed row neighborhood lists data (if same view_type). More...
 
template<typename view_type >
std::enable_if< view_type::rank==1 &&std::is_same< neighbor_lists_type::internal_view_type, view_type >::value==0, void >::type setNeighborLists (view_type neighbor_lists, view_type number_of_neighbors_list)
 Sets neighbor list information from compressed row neighborhood lists data (if different view_type). More...
 
template<typename view_type >
std::enable_if< view_type::rank==2, void >::type setNeighborLists (view_type neighbor_lists)
 Sets neighbor list information. More...
 
template<typename view_type >
void setSourceSites (view_type source_coordinates)
 Sets source coordinate information. More...
 
template<typename view_type >
void setSourceSites (coordinates_type source_coordinates)
 Sets source coordinate information. More...
 
template<typename view_type >
void setTargetSites (view_type target_coordinates)
 Sets target coordinate information. Rows of this 2D-array should correspond to rows of the neighbor lists. More...
 
template<typename view_type >
void setTargetSites (coordinates_type target_coordinates)
 Sets target coordinate information. Rows of this 2D-array should correspond to rows of the neighbor lists. More...
 
template<typename view_type >
void setWindowSizes (view_type epsilons)
 Sets window sizes, also called the support of the kernel. More...
 
template<typename view_type >
void setWindowSizes (decltype(_epsilons) epsilons)
 Sets window sizes, also called the support of the kernel (device) More...
 
template<typename view_type >
void setTangentBundle (view_type tangent_directions)
 (OPTIONAL) Sets orthonormal tangent directions for reconstruction on a manifold. More...
 
template<typename view_type >
void setReferenceOutwardNormalDirection (view_type outward_normal_directions, bool use_to_orient_surface=true)
 (OPTIONAL) Sets outward normal direction. More...
 
template<typename view_type >
void setSourceExtraData (view_type extra_data)
 (OPTIONAL) Sets extra data to be used by sampling functionals in certain instances. More...
 
template<typename view_type >
void setSourceExtraData (decltype(_source_extra_data) extra_data)
 (OPTIONAL) Sets extra data to be used by sampling functionals in certain instances. More...
 
template<typename view_type >
void setTargetExtraData (view_type extra_data)
 (OPTIONAL) Sets extra data to be used by target operations in certain instances. More...
 
template<typename view_type >
void setTargetExtraData (decltype(_target_extra_data) extra_data)
 (OPTIONAL) Sets extra data to be used by target operations in certain instances. More...
 
void setDenseSolverType (const DenseSolverType dst)
 Set dense solver type type. More...
 
void setConstraintType (const ConstraintType ct)
 Set constraint type. More...
 
void setWeightingType (const std::string &wt)
 Type for weighting kernel for GMLS problem. More...
 
void setWeightingType (const WeightingFunctionType wt)
 Type for weighting kernel for GMLS problem. More...
 
void setCurvatureWeightingType (const std::string &wt)
 Type for weighting kernel for curvature. More...
 
void setCurvatureWeightingType (const WeightingFunctionType wt)
 Type for weighting kernel for curvature. More...
 
void setPolynomialOrder (const int poly_order)
 Sets basis order to be used when reconstructing any function. More...
 
void setCurvaturePolynomialOrder (const int curvature_poly_order)
 Sets basis order to be used when reconstruction curvature. More...
 
void setWeightingParameter (int wp, int index=0)
 Parameter for weighting kernel for GMLS problem index = 0 sets p paramater for weighting kernel index = 1 sets n paramater for weighting kernel. More...
 
void setCurvatureWeightingParameter (int wp, int index=0)
 Parameter for weighting kernel for curvature index = 0 sets p paramater for weighting kernel index = 1 sets n paramater for weighting kernel. More...
 
void setOrderOfQuadraturePoints (int order)
 Number quadrature points to use. More...
 
void setDimensionOfQuadraturePoints (int dim)
 Dimensions of quadrature points to use. More...
 
void setQuadratureType (std::string quadrature_type)
 Type of quadrature points. More...
 
void addTargets (TargetOperation lro)
 Adds a target to the vector of target functional to be applied to the reconstruction. More...
 
void addTargets (std::vector< TargetOperation > lro)
 Adds a vector of target functionals to the vector of target functionals already to be applied to the reconstruction. More...
 
void clearTargets ()
 Empties the vector of target functionals to apply to the reconstruction. More...
 
bool verifyPointConnections (bool assert_valid=false)
 Verify whether _pc is valid. More...
 
bool verifyAdditionalPointConnections (bool assert_valid=false)
 Verify whether _additional_pc is valid. More...
 
void generatePolynomialCoefficients (const int number_of_batches=1, const bool keep_coefficients=false, const bool clear_cache=true)
 Generates polynomial coefficients by setting up and solving least squares problems Sets up the batch of GMLS problems to be solved for. Provides alpha values that can later be contracted against data or degrees of freedom to form a global linear system. More...
 
void generateAlphas (const int number_of_batches=1, const bool keep_coefficients=false, const bool clear_cache=true)
 Meant to calculate target operations and apply the evaluations to the previously constructed polynomial coefficients. But now that is inside of generatePolynomialCoefficients because it must be to handle number_of_batches>1. Effectively, this just calls generatePolynomialCoefficients. More...
 

Static Public Member Functions

Public Utility
static KOKKOS_INLINE_FUNCTION int getNP (const int m, const int dimension=3, const ReconstructionSpace r_space=ReconstructionSpace::ScalarTaylorPolynomial)
 Returns size of the basis for a given polynomial order and dimension General to dimension 1..3 and polynomial order m The divfree options will return the divergence-free basis if true. More...
 
static KOKKOS_INLINE_FUNCTION int getNN (const int m, const int dimension=3, const ReconstructionSpace r_space=ReconstructionSpace::ScalarTaylorPolynomial)
 Returns number of neighbors needed for unisolvency for a given basis order and dimension. More...
 
static KOKKOS_INLINE_FUNCTION double Wab (const double r, const double h, const WeightingFunctionType &weighting_type, const int p, const int n)
 Evaluates the weighting kernel. More...
 

Private Member Functions

Private Modifiers

Private function because information lives on the device

template<typename view_type >
void setAuxiliaryEvaluationCoordinates (view_type evaluation_coordinates)
 (OPTIONAL) Sets additional points for evaluation of target operation on polynomial reconstruction. More...
 
template<typename view_type >
void setAuxiliaryEvaluationCoordinates (coordinates_type evaluation_coordinates)
 (OPTIONAL) Sets additional points for evaluation of target operation on polynomial reconstruction. More...
 
template<typename view_type >
std::enable_if< view_type::rank==1 &&std::is_same< neighbor_lists_type::internal_view_type, view_type >::value==1, void >::type setAuxiliaryEvaluationIndicesLists (view_type additional_evaluation_indices, view_type number_of_neighbors_list)
 (OPTIONAL) Sets the additional target evaluation indices list information from compressed row format (if same view_type) More...
 
template<typename view_type >
std::enable_if< view_type::rank==1 &&std::is_same< neighbor_lists_type::internal_view_type, view_type >::value==0, void >::type setAuxiliaryEvaluationIndicesLists (view_type additional_evaluation_indices, view_type number_of_neighbors_list)
 (OPTIONAL) Sets the additional target evaluation indices list information from compressed row format (if different view_type) More...
 
template<typename view_type >
std::enable_if< view_type::rank==2, void >::type setAuxiliaryEvaluationIndicesLists (view_type additional_evaluation_indices)
 (OPTIONAL) Sets the additional target evaluation indices list information. More...
 

Static Private Member Functions

Private Utility
static DenseSolverType parseSolverType (const std::string &dense_solver_type)
 Parses a string to determine solver type. More...
 
static ProblemType parseProblemType (const std::string &problem_type)
 Parses a string to determine problem type. More...
 
static ConstraintType parseConstraintType (const std::string &constraint_type)
 Parses a string to determine constraint type. More...
 

Private Attributes

Kokkos::View< double * > _w
 contains weights for all problems More...
 
Kokkos::View< double * > _P
 P*sqrt(w) matrix for all problems. More...
 
Kokkos::View< double * > _RHS
 sqrt(w)*Identity matrix for all problems, later holds polynomial coefficients for all problems More...
 
Kokkos::View< double * > _Z
 stores evaluations of targets applied to basis More...
 
Kokkos::View< double * > _T
 Rank 3 tensor for high order approximation of tangent vectors for all problems. More...
 
Kokkos::View< double * > _ref_N
 Rank 2 tensor for high order approximation of tangent vectors for all problems. More...
 
Kokkos::View< double * >::HostMirror _host_T
 tangent vectors information (host) More...
 
Kokkos::View< double * >::HostMirror _host_ref_N
 reference outward normal vectors information (host) More...
 
Kokkos::View< double * > _manifold_metric_tensor_inverse
 metric tensor inverse for all problems More...
 
Kokkos::View< double * > _manifold_curvature_coefficients
 curvature polynomial coefficients for all problems More...
 
Kokkos::View< double * > _manifold_curvature_gradient
 _dimension-1 gradient values for curvature for all problems More...
 
Kokkos::View< double **, layout_right_source_extra_data
 Extra data available to basis functions (optional) More...
 
Kokkos::View< double **, layout_right_target_extra_data
 Extra data available to target operations (optional) More...
 
point_connections_type _pc
 connections between points and neighbors More...
 
Kokkos::View< double * > _epsilons
 h supports determined through neighbor search (device) More...
 
Kokkos::View< double *****, layout_right_prestencil_weights
 generated weights for nontraditional samples required to transform data into expected sampling functional form (device). More...
 
Kokkos::View< const double *****, layout_right >::HostMirror _host_prestencil_weights
 generated weights for nontraditional samples required to transform data into expected sampling functional form (host) More...
 
point_connections_type _additional_pc
 (OPTIONAL) connections between additional points and neighbors More...
 
SolutionSet< host_memory_space_h_ss
 Solution Set (contains all alpha values from solution and alpha layout methods) More...
 
SolutionSet< device_memory_space_d_ss
 
int _poly_order
 order of basis for polynomial reconstruction More...
 
int _curvature_poly_order
 order of basis for curvature reconstruction More...
 
int _NP
 dimension of basis for polynomial reconstruction More...
 
int _global_dimensions
 spatial dimension of the points, set at class instantiation only More...
 
int _local_dimensions
 dimension of the problem, set at class instantiation only. For manifolds, generally _global_dimensions-1 More...
 
int _dimensions
 dimension of the problem, set at class instantiation only More...
 
ReconstructionSpace _reconstruction_space
 reconstruction space for GMLS problems, set at GMLS class instantiation More...
 
int _reconstruction_space_rank
 actual rank of reconstruction basis More...
 
DenseSolverType _dense_solver_type
 solver type for GMLS problem - can be QR, SVD or LU More...
 
ProblemType _problem_type
 problem type for GMLS problem, can also be set to STANDARD for normal or MANIFOLD for manifold problems
NOTE: can only be set at object instantiation More...
 
ConstraintType _constraint_type
 constraint type for GMLS problem More...
 
SamplingFunctional _polynomial_sampling_functional
 polynomial sampling functional used to construct P matrix, set at GMLS class instantiation
NOTE: can only be set at object instantiation More...
 
SamplingFunctional _data_sampling_functional
 generally the same as _polynomial_sampling_functional, but can differ if specified at More...
 
Kokkos::View< TargetOperation * > _curvature_support_operations
 vector containing target functionals to be applied for curvature More...
 
Kokkos::View< TargetOperation * > _operations
 vector containing target functionals to be applied for reconstruction problem (device) More...
 
Kokkos::View< TargetOperation * >::HostMirror _host_operations
 vector containing target functionals to be applied for reconstruction problem (host) More...
 
WeightingFunctionType _weighting_type
 weighting kernel type for GMLS More...
 
WeightingFunctionType _curvature_weighting_type
 weighting kernel type for curvature problem More...
 
int _weighting_p
 first parameter to be used for weighting kernel More...
 
int _weighting_n
 second parameter to be used for weighting kernel More...
 
int _curvature_weighting_p
 first parameter to be used for weighting kernel for curvature More...
 
int _curvature_weighting_n
 second parameter to be used for weighting kernel for curvature More...
 
int _basis_multiplier
 dimension of the reconstructed function e.g. More...
 
int _sampling_multiplier
 actual dimension of the sampling functional e.g. More...
 
int _data_sampling_multiplier
 effective dimension of the data sampling functional e.g. More...
 
bool _orthonormal_tangent_space_provided
 whether or not the orthonormal tangent directions were provided by the user. More...
 
bool _reference_outward_normal_direction_provided
 whether or not the reference outward normal directions were provided by the user. More...
 
bool _use_reference_outward_normal_direction_provided_to_orient_surface
 whether or not to use reference outward normal directions to orient the surface in a manifold problem. More...
 
bool _entire_batch_computed_at_once
 whether entire calculation was computed at once the alternative is that it was broken up over many smaller groups, in which case this is false, and so the _RHS matrix can not be stored or requested More...
 
bool _store_PTWP_inv_PTW
 whether polynomial coefficients were requested to be stored (in a state not yet applied to data) More...
 
int _initial_index_for_batch
 initial index for current batch More...
 
ParallelManager _pm
 determines scratch level spaces and is used to call kernels More...
 
int _order_of_quadrature_points
 order of exact polynomial integration for quadrature rule More...
 
int _dimension_of_quadrature_points
 dimension of quadrature rule More...
 
std::string _quadrature_type
 quadrature rule type More...
 
Quadrature _qm
 manages and calculates quadrature More...
 

Friends

class Evaluator
 

Member Typedef Documentation

◆ coordinates_type

typedef Kokkos::View<double**, layout_right> Compadre::GMLS::coordinates_type

Definition at line 45 of file Compadre_GMLS.hpp.

◆ neighbor_lists_type

typedef NeighborLists<Kokkos::View<int*> > Compadre::GMLS::neighbor_lists_type

Definition at line 43 of file Compadre_GMLS.hpp.

◆ point_connections_type

typedef PointConnections<Kokkos::View<double**, layout_right>, Kokkos::View<double**, layout_right>, NeighborLists<Kokkos::View<int*> > > Compadre::GMLS::point_connections_type

Definition at line 41 of file Compadre_GMLS.hpp.

Constructor & Destructor Documentation

◆ GMLS() [1/4]

Compadre::GMLS::GMLS ( ReconstructionSpace  reconstruction_space,
const SamplingFunctional  polynomial_sampling_strategy,
const SamplingFunctional  data_sampling_strategy,
const int  poly_order,
const int  dimensions,
const DenseSolverType  dense_solver_type,
const ProblemType  problem_type,
const ConstraintType  constraint_type,
const int  manifold_curvature_poly_order 
)
inline

Maximal constructor, not intended for users.

Definition at line 381 of file Compadre_GMLS.hpp.

◆ GMLS() [2/4]

Compadre::GMLS::GMLS ( ReconstructionSpace  reconstruction_space,
const SamplingFunctional  polynomial_sampling_strategy,
const SamplingFunctional  data_sampling_strategy,
const int  poly_order,
const int  dimensions = 3,
const std::string  dense_solver_type = std::string("QR"),
const std::string  problem_type = std::string("STANDARD"),
const std::string  constraint_type = std::string("NO_CONSTRAINT"),
const int  manifold_curvature_poly_order = 2 
)
inline

Maximal constructor, but with string arguments.

Definition at line 460 of file Compadre_GMLS.hpp.

◆ GMLS() [3/4]

Compadre::GMLS::GMLS ( const int  poly_order,
const int  dimensions = 3,
const std::string  dense_solver_type = std::string("QR"),
const std::string  problem_type = std::string("STANDARD"),
const std::string  constraint_type = std::string("NO_CONSTRAINT"),
const int  manifold_curvature_poly_order = 2 
)
inline

Constructor for the case when the data sampling functional does not match the polynomial sampling functional.

Only case anticipated is staggered Laplacian.

Definition at line 473 of file Compadre_GMLS.hpp.

◆ GMLS() [4/4]

Compadre::GMLS::GMLS ( ReconstructionSpace  reconstruction_space,
SamplingFunctional  dual_sampling_strategy,
const int  poly_order,
const int  dimensions = 3,
const std::string  dense_solver_type = std::string("QR"),
const std::string  problem_type = std::string("STANDARD"),
const std::string  constraint_type = std::string("NO_CONSTRAINT"),
const int  manifold_curvature_poly_order = 2 
)
inline

Constructor for the case when nonstandard sampling functionals or reconstruction spaces are to be used.

Reconstruction space and sampling strategy can only be set at instantiation.

Definition at line 483 of file Compadre_GMLS.hpp.

Member Function Documentation

◆ addTargets() [1/2]

void Compadre::GMLS::addTargets ( std::vector< TargetOperation lro)
inline

Adds a vector of target functionals to the vector of target functionals already to be applied to the reconstruction.

Definition at line 1272 of file Compadre_GMLS.hpp.

◆ addTargets() [2/2]

void Compadre::GMLS::addTargets ( TargetOperation  lro)
inline

Adds a target to the vector of target functional to be applied to the reconstruction.

Definition at line 1266 of file Compadre_GMLS.hpp.

◆ clearTargets()

void Compadre::GMLS::clearTargets ( )
inline

Empties the vector of target functionals to apply to the reconstruction.

Definition at line 1278 of file Compadre_GMLS.hpp.

◆ containsValidAlphas()

bool Compadre::GMLS::containsValidAlphas ( ) const
inline

Check if GMLS solution set contains valid alpha values (has generateAlphas been called)

Definition at line 808 of file Compadre_GMLS.hpp.

◆ extractBasisData()

const GMLSBasisData Compadre::GMLS::extractBasisData ( ) const

Get GMLS basis data.

Definition at line 521 of file Compadre_GMLS.cpp.

◆ extractSolutionData()

const GMLSSolutionData Compadre::GMLS::extractSolutionData ( ) const

Get GMLS solution data.

Definition at line 473 of file Compadre_GMLS.cpp.

◆ generateAlphas()

void Compadre::GMLS::generateAlphas ( const int  number_of_batches = 1,
const bool  keep_coefficients = false,
const bool  clear_cache = true 
)

Meant to calculate target operations and apply the evaluations to the previously constructed polynomial coefficients. But now that is inside of generatePolynomialCoefficients because it must be to handle number_of_batches>1. Effectively, this just calls generatePolynomialCoefficients.

Parameters
number_of_batches[in] - how many batches to break up the total workload into (for storage)
keep_coefficients[in] - whether to store (P^T W P)^-1 * P^T * W
clear_cache[in] - clear whatever temporary memory was used for calculations that
keep_coefficients hasn't indicated needs to be saved

clear_cache should be false to be used for debugging as it will leave all data structures used to compute the solution intact. If clear_cache is set true and keep_coefficients is true, it will allow the polynomial coefficients to persist after calculation.

Definition at line 467 of file Compadre_GMLS.cpp.

◆ generatePolynomialCoefficients()

void Compadre::GMLS::generatePolynomialCoefficients ( const int  number_of_batches = 1,
const bool  keep_coefficients = false,
const bool  clear_cache = true 
)

Generates polynomial coefficients by setting up and solving least squares problems Sets up the batch of GMLS problems to be solved for. Provides alpha values that can later be contracted against data or degrees of freedom to form a global linear system.

Parameters
number_of_batches[in] - how many batches to break up the total workload into (for storage)
keep_coefficients[in] - whether to store (P^T W P)^-1 * P^T * W
clear_cache[in] - clear whatever temporary memory was used for calculations that
keep_coefficients hasn't indicated needs to be saved

clear_cache should be false to be used for debugging as it will leave all data structures used to compute the solution intact. If clear_cache is set true and keep_coefficients is true, it will allow the polynomial coefficients to persist after calculation.

Definition at line 6 of file Compadre_GMLS.cpp.

◆ getAdditionalEvaluationIndices()

neighbor_lists_type* Compadre::GMLS::getAdditionalEvaluationIndices ( ) const
inline

(OPTIONAL) Get additional evaluation sites neighbor list-like accessor

Definition at line 690 of file Compadre_GMLS.hpp.

◆ getAdditionalPointConnections()

decltype(_additional_pc) * Compadre::GMLS::getAdditionalPointConnections ( )
inline

(OPTIONAL) Get a view (device) of all additional evaluation point connection info

Definition at line 695 of file Compadre_GMLS.hpp.

◆ getConstraintType()

ConstraintType Compadre::GMLS::getConstraintType ( ) const
inline

Get constraint type.

Definition at line 637 of file Compadre_GMLS.hpp.

◆ getCurvaturePolynomialOrder()

int Compadre::GMLS::getCurvaturePolynomialOrder ( ) const
inline

Get basis order used for curvature reconstruction.

Definition at line 643 of file Compadre_GMLS.hpp.

◆ getDataSamplingFunctional()

SamplingFunctional Compadre::GMLS::getDataSamplingFunctional ( ) const
inline

Get the data sampling functional specified at instantiation (often the same as the polynomial sampling functional)

Definition at line 759 of file Compadre_GMLS.hpp.

◆ getDenseSolverType()

DenseSolverType Compadre::GMLS::getDenseSolverType ( ) const
inline

Get dense solver type.

Definition at line 631 of file Compadre_GMLS.hpp.

◆ getDimensionOfQuadraturePoints()

int Compadre::GMLS::getDimensionOfQuadraturePoints ( ) const
inline

Dimensions of quadrature points.

Definition at line 676 of file Compadre_GMLS.hpp.

◆ getDimensions()

int Compadre::GMLS::getDimensions ( ) const
inline

Dimension of the GMLS problem, set only at class instantiation.

Definition at line 622 of file Compadre_GMLS.hpp.

◆ getFullPolynomialCoefficientsBasis()

decltype(_RHS) Compadre::GMLS::getFullPolynomialCoefficientsBasis ( ) const
inline

Get a view (device) of all polynomial coefficients basis.

Definition at line 743 of file Compadre_GMLS.hpp.

◆ getGlobalDimensions()

int Compadre::GMLS::getGlobalDimensions ( ) const
inline

Dimension of the GMLS problem's point data (spatial description of points in ambient space), set only at class instantiation.

Definition at line 625 of file Compadre_GMLS.hpp.

◆ getLocalDimensions()

int Compadre::GMLS::getLocalDimensions ( ) const
inline

Local dimension of the GMLS problem (less than global dimension if on a manifold), set only at class instantiation.

Definition at line 628 of file Compadre_GMLS.hpp.

◆ getManifoldWeightingParameter()

int Compadre::GMLS::getManifoldWeightingParameter ( const int  index = 0) const
inline

Get parameter for weighting kernel for curvature.

Definition at line 661 of file Compadre_GMLS.hpp.

◆ getManifoldWeightingType()

WeightingFunctionType Compadre::GMLS::getManifoldWeightingType ( ) const
inline

Type for weighting kernel for curvature.

Definition at line 649 of file Compadre_GMLS.hpp.

◆ getNeighborLists()

neighbor_lists_type* Compadre::GMLS::getNeighborLists ( ) const
inline

Get neighbor list accessor.

Definition at line 682 of file Compadre_GMLS.hpp.

◆ getNN()

static KOKKOS_INLINE_FUNCTION int Compadre::GMLS::getNN ( const int  m,
const int  dimension = 3,
const ReconstructionSpace  r_space = ReconstructionSpace::ScalarTaylorPolynomial 
)
inlinestatic

Returns number of neighbors needed for unisolvency for a given basis order and dimension.

Definition at line 516 of file Compadre_GMLS.hpp.

◆ getNP()

static KOKKOS_INLINE_FUNCTION int Compadre::GMLS::getNP ( const int  m,
const int  dimension = 3,
const ReconstructionSpace  r_space = ReconstructionSpace::ScalarTaylorPolynomial 
)
inlinestatic

Returns size of the basis for a given polynomial order and dimension General to dimension 1..3 and polynomial order m The divfree options will return the divergence-free basis if true.

Examples
GMLS Tutorial, and Manifold GMLS Tutorial.

Definition at line 504 of file Compadre_GMLS.hpp.

◆ getNumberOfQuadraturePoints()

int Compadre::GMLS::getNumberOfQuadraturePoints ( ) const
inline

Number of quadrature points.

Definition at line 670 of file Compadre_GMLS.hpp.

◆ getOrderOfQuadraturePoints()

int Compadre::GMLS::getOrderOfQuadraturePoints ( ) const
inline

Order of quadrature points.

Definition at line 673 of file Compadre_GMLS.hpp.

◆ getPointConnections()

decltype(_pc) * Compadre::GMLS::getPointConnections ( )
inline

Get a view (device) of all point connection info.

Definition at line 687 of file Compadre_GMLS.hpp.

◆ getPolynomialCoefficientsDomainRangeSize()

host_managed_local_index_type Compadre::GMLS::getPolynomialCoefficientsDomainRangeSize ( ) const
inline

Returns (size of the basis used in instance's polynomial reconstruction) x (data input dimension)

Definition at line 592 of file Compadre_GMLS.hpp.

◆ getPolynomialCoefficientsMemorySize()

host_managed_local_index_type Compadre::GMLS::getPolynomialCoefficientsMemorySize ( ) const
inline

Returns 2D array size in memory on which coefficients are stored.

Definition at line 606 of file Compadre_GMLS.hpp.

◆ getPolynomialCoefficientsSize()

int Compadre::GMLS::getPolynomialCoefficientsSize ( ) const
inline

Returns size of the basis used in instance's polynomial reconstruction.

Definition at line 600 of file Compadre_GMLS.hpp.

◆ getPolynomialOrder()

int Compadre::GMLS::getPolynomialOrder ( ) const
inline

Get basis order used for reconstruction.

Definition at line 640 of file Compadre_GMLS.hpp.

◆ getPolynomialSamplingFunctional()

SamplingFunctional Compadre::GMLS::getPolynomialSamplingFunctional ( ) const
inline

Get the polynomial sampling functional specified at instantiation.

Definition at line 756 of file Compadre_GMLS.hpp.

◆ getPreStencilWeight()

double Compadre::GMLS::getPreStencilWeight ( SamplingFunctional  sro,
const int  target_index,
const int  neighbor_index,
bool  for_target,
const int  output_component = 0,
const int  input_component = 0 
) const
inline

Returns a stencil to transform data from its existing state into the input expected for some sampling functionals.

Definition at line 766 of file Compadre_GMLS.hpp.

◆ getPrestencilWeights()

decltype(_prestencil_weights) Compadre::GMLS::getPrestencilWeights ( ) const
inline

Get a view (device) of all rank 2 preprocessing tensors This is a rank 5 tensor that is able to provide data transformation into a form that GMLS is able to operate on.

The ranks are as follows:

1 - Either size 2 if it operates on the target site and neighbor site (staggered schemes) or 1 if it operates only on the neighbor sites (almost every scheme)

2 - If the data transform varies with each target site (but could be the same for each neighbor of that target site), then this is the number of target sites

3 - If the data transform varies with each neighbor of each target site, then this is the number of neighbors for each respective target (max number of neighbors for all target sites is its uniform size)

4 - Data transform resulting in rank 1 data for the GMLS operator will have size _local_dimensions, otherwise 1

5 - Data transform taking in rank 1 data will have size _global_dimensions, otherwise 1

Definition at line 738 of file Compadre_GMLS.hpp.

◆ getProblemType()

ProblemType Compadre::GMLS::getProblemType ( ) const
inline

Get problem type.

Definition at line 634 of file Compadre_GMLS.hpp.

◆ getQuadratureType()

std::string Compadre::GMLS::getQuadratureType ( ) const
inline

Type of quadrature points.

Definition at line 679 of file Compadre_GMLS.hpp.

◆ getReconstructionSpace()

ReconstructionSpace Compadre::GMLS::getReconstructionSpace ( ) const
inline

Get the reconstruction space specified at instantiation.

Definition at line 762 of file Compadre_GMLS.hpp.

◆ getReferenceNormalDirection()

double Compadre::GMLS::getReferenceNormalDirection ( const int  target_index,
const int  component 
) const
inline

Get component of tangent or normal directions for manifold problems.

Definition at line 716 of file Compadre_GMLS.hpp.

◆ getReferenceNormalDirections()

decltype(_ref_N) * Compadre::GMLS::getReferenceNormalDirections ( )
inline

Get a view (device) of all reference outward normal directions.

Definition at line 704 of file Compadre_GMLS.hpp.

◆ getSolutionSetDevice()

decltype(_d_ss) * Compadre::GMLS::getSolutionSetDevice ( bool  alpha_validity_check = true)
inline

Get solution set on device.

Definition at line 801 of file Compadre_GMLS.hpp.

◆ getSolutionSetHost()

decltype(_h_ss) * Compadre::GMLS::getSolutionSetHost ( bool  alpha_validity_check = true)
inline

Get solution set on host.

Definition at line 789 of file Compadre_GMLS.hpp.

◆ getTangentBundle()

double Compadre::GMLS::getTangentBundle ( const int  target_index,
const int  direction,
const int  component 
) const
inline

Get component of tangent or normal directions for manifold problems.

Definition at line 707 of file Compadre_GMLS.hpp.

◆ getTangentDirections()

decltype(_T) * Compadre::GMLS::getTangentDirections ( )
inline

Get a view (device) of all tangent direction bundles.

Definition at line 701 of file Compadre_GMLS.hpp.

◆ getWeightingParameter()

int Compadre::GMLS::getWeightingParameter ( const int  index = 0) const
inline

Get parameter for weighting kernel for GMLS problem.

Definition at line 652 of file Compadre_GMLS.hpp.

◆ getWeightingType()

WeightingFunctionType Compadre::GMLS::getWeightingType ( ) const
inline

Type for weighting kernel for GMLS problem.

Definition at line 646 of file Compadre_GMLS.hpp.

◆ getWindowSizes()

decltype(_epsilons) * Compadre::GMLS::getWindowSizes ( )
inline

Get a view (device) of all window sizes.

Definition at line 698 of file Compadre_GMLS.hpp.

◆ parseConstraintType()

static ConstraintType Compadre::GMLS::parseConstraintType ( const std::string &  constraint_type)
inlinestaticprivate

Parses a string to determine constraint type.

Definition at line 359 of file Compadre_GMLS.hpp.

◆ parseProblemType()

static ProblemType Compadre::GMLS::parseProblemType ( const std::string &  problem_type)
inlinestaticprivate

Parses a string to determine problem type.

Definition at line 346 of file Compadre_GMLS.hpp.

◆ parseSolverType()

static DenseSolverType Compadre::GMLS::parseSolverType ( const std::string &  dense_solver_type)
inlinestaticprivate

Parses a string to determine solver type.

Definition at line 335 of file Compadre_GMLS.hpp.

◆ resetCoefficientData()

void Compadre::GMLS::resetCoefficientData ( )
inline

Definition at line 824 of file Compadre_GMLS.hpp.

◆ setAdditionalEvaluationSitesData() [1/2]

template<typename view_type_1 , typename view_type_2 >
void Compadre::GMLS::setAdditionalEvaluationSitesData ( view_type_1  additional_evaluation_indices,
view_type_2  additional_evaluation_coordinates 
)
inline

(OPTIONAL) Sets additional evaluation sites for each target site

Definition at line 860 of file Compadre_GMLS.hpp.

◆ setAdditionalEvaluationSitesData() [2/2]

template<typename view_type_1 , typename view_type_2 >
void Compadre::GMLS::setAdditionalEvaluationSitesData ( view_type_1  cr_additional_evaluation_indices,
view_type_1  number_of_additional_evaluation_indices,
view_type_2  additional_evaluation_coordinates 
)
inline

(OPTIONAL) Sets additional evaluation sites for each target site

Definition at line 869 of file Compadre_GMLS.hpp.

◆ setAuxiliaryEvaluationCoordinates() [1/2]

template<typename view_type >
void Compadre::GMLS::setAuxiliaryEvaluationCoordinates ( coordinates_type  evaluation_coordinates)
inlineprivate

(OPTIONAL) Sets additional points for evaluation of target operation on polynomial reconstruction.

If this is never called, then the target sites are the only locations where the target operations will be evaluated and applied to polynomial reconstructions. (device)

Definition at line 279 of file Compadre_GMLS.hpp.

◆ setAuxiliaryEvaluationCoordinates() [2/2]

template<typename view_type >
void Compadre::GMLS::setAuxiliaryEvaluationCoordinates ( view_type  evaluation_coordinates)
inlineprivate

(OPTIONAL) Sets additional points for evaluation of target operation on polynomial reconstruction.

If this is never called, then the target sites are the only locations where the target operations will be evaluated and applied to polynomial reconstructions.

Definition at line 249 of file Compadre_GMLS.hpp.

◆ setAuxiliaryEvaluationIndicesLists() [1/3]

template<typename view_type >
std::enable_if<view_type::rank==2, void>::type Compadre::GMLS::setAuxiliaryEvaluationIndicesLists ( view_type  additional_evaluation_indices)
inlineprivate

(OPTIONAL) Sets the additional target evaluation indices list information.

Should be # targets x maximum number of indices evaluation indices for any target + 1. first entry in every row should be the number of indices for the corresponding target.

Definition at line 318 of file Compadre_GMLS.hpp.

◆ setAuxiliaryEvaluationIndicesLists() [2/3]

template<typename view_type >
std::enable_if<view_type::rank==1&&std::is_same<neighbor_lists_type::internal_view_type,view_type>::value==1, void>::type Compadre::GMLS::setAuxiliaryEvaluationIndicesLists ( view_type  additional_evaluation_indices,
view_type  number_of_neighbors_list 
)
inlineprivate

(OPTIONAL) Sets the additional target evaluation indices list information from compressed row format (if same view_type)

Definition at line 288 of file Compadre_GMLS.hpp.

◆ setAuxiliaryEvaluationIndicesLists() [3/3]

template<typename view_type >
std::enable_if<view_type::rank==1&&std::is_same<neighbor_lists_type::internal_view_type,view_type>::value==0, void>::type Compadre::GMLS::setAuxiliaryEvaluationIndicesLists ( view_type  additional_evaluation_indices,
view_type  number_of_neighbors_list 
)
inlineprivate

(OPTIONAL) Sets the additional target evaluation indices list information from compressed row format (if different view_type)

Definition at line 300 of file Compadre_GMLS.hpp.

◆ setConstraintType()

void Compadre::GMLS::setConstraintType ( const ConstraintType  ct)
inline

Set constraint type.

Definition at line 1149 of file Compadre_GMLS.hpp.

◆ setCurvaturePolynomialOrder()

void Compadre::GMLS::setCurvaturePolynomialOrder ( const int  curvature_poly_order)
inline

Sets basis order to be used when reconstruction curvature.

Definition at line 1217 of file Compadre_GMLS.hpp.

◆ setCurvatureWeightingParameter()

void Compadre::GMLS::setCurvatureWeightingParameter ( int  wp,
int  index = 0 
)
inline

Parameter for weighting kernel for curvature index = 0 sets p paramater for weighting kernel index = 1 sets n paramater for weighting kernel.

Definition at line 1238 of file Compadre_GMLS.hpp.

◆ setCurvatureWeightingType() [1/2]

void Compadre::GMLS::setCurvatureWeightingType ( const std::string &  wt)
inline

Type for weighting kernel for curvature.

Definition at line 1182 of file Compadre_GMLS.hpp.

◆ setCurvatureWeightingType() [2/2]

void Compadre::GMLS::setCurvatureWeightingType ( const WeightingFunctionType  wt)
inline

Type for weighting kernel for curvature.

Definition at line 1203 of file Compadre_GMLS.hpp.

◆ setDenseSolverType()

void Compadre::GMLS::setDenseSolverType ( const DenseSolverType  dst)
inline

Set dense solver type type.

Definition at line 1143 of file Compadre_GMLS.hpp.

◆ setDimensionOfQuadraturePoints()

void Compadre::GMLS::setDimensionOfQuadraturePoints ( int  dim)
inline

Dimensions of quadrature points to use.

Definition at line 1254 of file Compadre_GMLS.hpp.

◆ setNeighborLists() [1/3]

template<typename view_type >
std::enable_if<view_type::rank==2, void>::type Compadre::GMLS::setNeighborLists ( view_type  neighbor_lists)
inline

Sets neighbor list information.

Should be # targets x maximum number of neighbors for any target + 1. first entry in ever row should be the number of neighbors for the corresponding target.

Definition at line 907 of file Compadre_GMLS.hpp.

◆ setNeighborLists() [2/3]

template<typename view_type >
std::enable_if<view_type::rank==1&&std::is_same<neighbor_lists_type::internal_view_type,view_type>::value==1, void>::type Compadre::GMLS::setNeighborLists ( view_type  neighbor_lists,
view_type  number_of_neighbors_list 
)
inline

Sets neighbor list information from compressed row neighborhood lists data (if same view_type).

Definition at line 881 of file Compadre_GMLS.hpp.

◆ setNeighborLists() [3/3]

template<typename view_type >
std::enable_if<view_type::rank==1&&std::is_same<neighbor_lists_type::internal_view_type,view_type>::value==0, void>::type Compadre::GMLS::setNeighborLists ( view_type  neighbor_lists,
view_type  number_of_neighbors_list 
)
inline

Sets neighbor list information from compressed row neighborhood lists data (if different view_type).

Definition at line 891 of file Compadre_GMLS.hpp.

◆ setOrderOfQuadraturePoints()

void Compadre::GMLS::setOrderOfQuadraturePoints ( int  order)
inline

Number quadrature points to use.

Definition at line 1248 of file Compadre_GMLS.hpp.

◆ setPolynomialOrder()

void Compadre::GMLS::setPolynomialOrder ( const int  poly_order)
inline

Sets basis order to be used when reconstructing any function.

Definition at line 1209 of file Compadre_GMLS.hpp.

◆ setProblemData() [1/2]

template<typename view_type_1 , typename view_type_2 , typename view_type_3 , typename view_type_4 >
void Compadre::GMLS::setProblemData ( view_type_1  cr_neighbor_lists,
view_type_1  number_of_neighbors_list,
view_type_2  source_coordinates,
view_type_3  target_coordinates,
view_type_4  epsilons 
)
inline

Sets basic problem data (neighbor lists data, number of neighbors list, source coordinates, and target coordinates)

Definition at line 846 of file Compadre_GMLS.hpp.

◆ setProblemData() [2/2]

template<typename view_type_1 , typename view_type_2 , typename view_type_3 , typename view_type_4 >
void Compadre::GMLS::setProblemData ( view_type_1  neighbor_lists,
view_type_2  source_coordinates,
view_type_3  target_coordinates,
view_type_4  epsilons 
)
inline

Sets basic problem data (neighbor lists, source coordinates, and target coordinates)

Definition at line 833 of file Compadre_GMLS.hpp.

◆ setQuadratureType()

void Compadre::GMLS::setQuadratureType ( std::string  quadrature_type)
inline

Type of quadrature points.

Definition at line 1260 of file Compadre_GMLS.hpp.

◆ setReferenceOutwardNormalDirection()

template<typename view_type >
void Compadre::GMLS::setReferenceOutwardNormalDirection ( view_type  outward_normal_directions,
bool  use_to_orient_surface = true 
)
inline

(OPTIONAL) Sets outward normal direction.

For manifolds this may be used for orienting surface. It is also accessible for sampling operators that require a normal direction.

Definition at line 1069 of file Compadre_GMLS.hpp.

◆ setSourceExtraData() [1/2]

template<typename view_type >
void Compadre::GMLS::setSourceExtraData ( decltype(_source_extra_data extra_data)
inline

(OPTIONAL) Sets extra data to be used by sampling functionals in certain instances.

(device)

Definition at line 1112 of file Compadre_GMLS.hpp.

◆ setSourceExtraData() [2/2]

template<typename view_type >
void Compadre::GMLS::setSourceExtraData ( view_type  extra_data)
inline

(OPTIONAL) Sets extra data to be used by sampling functionals in certain instances.

Definition at line 1097 of file Compadre_GMLS.hpp.

◆ setSourceSites() [1/2]

template<typename view_type >
void Compadre::GMLS::setSourceSites ( coordinates_type  source_coordinates)
inline

Sets source coordinate information.

Rows of this 2D-array should correspond to neighbor IDs contained in the entries of the neighbor lists 2D array.

Definition at line 946 of file Compadre_GMLS.hpp.

◆ setSourceSites() [2/2]

template<typename view_type >
void Compadre::GMLS::setSourceSites ( view_type  source_coordinates)
inline

Sets source coordinate information.

Rows of this 2D-array should correspond to neighbor IDs contained in the entries of the neighbor lists 2D array.

Definition at line 917 of file Compadre_GMLS.hpp.

◆ setTangentBundle()

template<typename view_type >
void Compadre::GMLS::setTangentBundle ( view_type  tangent_directions)
inline

(OPTIONAL) Sets orthonormal tangent directions for reconstruction on a manifold.

The first rank of this 2D array corresponds to the target indices, i.e., rows of the neighbor lists 2D array. The second rank is the ordinal of the tangent direction (spatial dimensions-1 are tangent, last one is normal), and the third rank is indices into the spatial dimension.

Definition at line 1032 of file Compadre_GMLS.hpp.

◆ setTargetExtraData() [1/2]

template<typename view_type >
void Compadre::GMLS::setTargetExtraData ( decltype(_target_extra_data extra_data)
inline

(OPTIONAL) Sets extra data to be used by target operations in certain instances.

(device)

Definition at line 1136 of file Compadre_GMLS.hpp.

◆ setTargetExtraData() [2/2]

template<typename view_type >
void Compadre::GMLS::setTargetExtraData ( view_type  extra_data)
inline

(OPTIONAL) Sets extra data to be used by target operations in certain instances.

Definition at line 1121 of file Compadre_GMLS.hpp.

◆ setTargetSites() [1/2]

template<typename view_type >
void Compadre::GMLS::setTargetSites ( coordinates_type  target_coordinates)
inline

Sets target coordinate information. Rows of this 2D-array should correspond to rows of the neighbor lists.

Definition at line 991 of file Compadre_GMLS.hpp.

◆ setTargetSites() [2/2]

template<typename view_type >
void Compadre::GMLS::setTargetSites ( view_type  target_coordinates)
inline

Sets target coordinate information. Rows of this 2D-array should correspond to rows of the neighbor lists.

Definition at line 954 of file Compadre_GMLS.hpp.

◆ setWeightingParameter()

void Compadre::GMLS::setWeightingParameter ( int  wp,
int  index = 0 
)
inline

Parameter for weighting kernel for GMLS problem index = 0 sets p paramater for weighting kernel index = 1 sets n paramater for weighting kernel.

Definition at line 1226 of file Compadre_GMLS.hpp.

◆ setWeightingType() [1/2]

void Compadre::GMLS::setWeightingType ( const std::string &  wt)
inline

Type for weighting kernel for GMLS problem.

Definition at line 1155 of file Compadre_GMLS.hpp.

◆ setWeightingType() [2/2]

void Compadre::GMLS::setWeightingType ( const WeightingFunctionType  wt)
inline

Type for weighting kernel for GMLS problem.

Definition at line 1176 of file Compadre_GMLS.hpp.

◆ setWindowSizes() [1/2]

template<typename view_type >
void Compadre::GMLS::setWindowSizes ( decltype(_epsilons epsilons)
inline

Sets window sizes, also called the support of the kernel (device)

Definition at line 1020 of file Compadre_GMLS.hpp.

◆ setWindowSizes() [2/2]

template<typename view_type >
void Compadre::GMLS::setWindowSizes ( view_type  epsilons)
inline

Sets window sizes, also called the support of the kernel.

Definition at line 1008 of file Compadre_GMLS.hpp.

◆ verifyAdditionalPointConnections()

bool Compadre::GMLS::verifyAdditionalPointConnections ( bool  assert_valid = false)
inline

Verify whether _additional_pc is valid.

Definition at line 1300 of file Compadre_GMLS.hpp.

◆ verifyPointConnections()

bool Compadre::GMLS::verifyPointConnections ( bool  assert_valid = false)
inline

Verify whether _pc is valid.

Definition at line 1284 of file Compadre_GMLS.hpp.

◆ Wab()

static KOKKOS_INLINE_FUNCTION double Compadre::GMLS::Wab ( const double  r,
const double  h,
const WeightingFunctionType weighting_type,
const int  p,
const int  n 
)
inlinestatic

Evaluates the weighting kernel.

Parameters
r[in] - Euclidean distance of relative vector. Euclidean distance of (target - neighbor) in some basis.
h[in] - window size. Kernel is guaranteed to take on a value of zero if it exceeds h.
weighting_type[in] - weighting type to be evaluated as the kernel. e,g. power, Gaussian, etc..
p[in] - parameter to be given to the kernel (see Wab definition for context).
n[in] - parameter to be given to the kernel (see Wab definition for context).

Definition at line 541 of file Compadre_GMLS.hpp.

Friends And Related Function Documentation

◆ Evaluator

friend class Evaluator
friend

Definition at line 34 of file Compadre_GMLS.hpp.

Member Data Documentation

◆ _additional_pc

point_connections_type Compadre::GMLS::_additional_pc
private

(OPTIONAL) connections between additional points and neighbors

Definition at line 109 of file Compadre_GMLS.hpp.

◆ _basis_multiplier

int Compadre::GMLS::_basis_multiplier
private

dimension of the reconstructed function e.g.

reconstruction of vector on a 2D manifold in 3D would have _basis_multiplier of 2

Definition at line 189 of file Compadre_GMLS.hpp.

◆ _constraint_type

ConstraintType Compadre::GMLS::_constraint_type
private

constraint type for GMLS problem

Definition at line 149 of file Compadre_GMLS.hpp.

◆ _curvature_poly_order

int Compadre::GMLS::_curvature_poly_order
private

order of basis for curvature reconstruction

Definition at line 121 of file Compadre_GMLS.hpp.

◆ _curvature_support_operations

Kokkos::View<TargetOperation*> Compadre::GMLS::_curvature_support_operations
private

vector containing target functionals to be applied for curvature

Definition at line 161 of file Compadre_GMLS.hpp.

◆ _curvature_weighting_n

int Compadre::GMLS::_curvature_weighting_n
private

second parameter to be used for weighting kernel for curvature

Definition at line 185 of file Compadre_GMLS.hpp.

◆ _curvature_weighting_p

int Compadre::GMLS::_curvature_weighting_p
private

first parameter to be used for weighting kernel for curvature

Definition at line 182 of file Compadre_GMLS.hpp.

◆ _curvature_weighting_type

WeightingFunctionType Compadre::GMLS::_curvature_weighting_type
private

weighting kernel type for curvature problem

Definition at line 173 of file Compadre_GMLS.hpp.

◆ _d_ss

SolutionSet<device_memory_space> Compadre::GMLS::_d_ss
private

Definition at line 115 of file Compadre_GMLS.hpp.

◆ _data_sampling_functional

SamplingFunctional Compadre::GMLS::_data_sampling_functional
private

generally the same as _polynomial_sampling_functional, but can differ if specified at

GMLS class instantiation

Definition at line 158 of file Compadre_GMLS.hpp.

◆ _data_sampling_multiplier

int Compadre::GMLS::_data_sampling_multiplier
private

effective dimension of the data sampling functional e.g.

in 3D, a scalar will be 1, a vector will be 3, and a vector of reused scalars will be 3

Definition at line 198 of file Compadre_GMLS.hpp.

◆ _dense_solver_type

DenseSolverType Compadre::GMLS::_dense_solver_type
private

solver type for GMLS problem - can be QR, SVD or LU

Definition at line 142 of file Compadre_GMLS.hpp.

◆ _dimension_of_quadrature_points

int Compadre::GMLS::_dimension_of_quadrature_points
private

dimension of quadrature rule

Definition at line 229 of file Compadre_GMLS.hpp.

◆ _dimensions

int Compadre::GMLS::_dimensions
private

dimension of the problem, set at class instantiation only

Definition at line 133 of file Compadre_GMLS.hpp.

◆ _entire_batch_computed_at_once

bool Compadre::GMLS::_entire_batch_computed_at_once
private

whether entire calculation was computed at once the alternative is that it was broken up over many smaller groups, in which case this is false, and so the _RHS matrix can not be stored or requested

Definition at line 214 of file Compadre_GMLS.hpp.

◆ _epsilons

Kokkos::View<double*> Compadre::GMLS::_epsilons
private

h supports determined through neighbor search (device)

Definition at line 98 of file Compadre_GMLS.hpp.

◆ _global_dimensions

int Compadre::GMLS::_global_dimensions
private

spatial dimension of the points, set at class instantiation only

Definition at line 127 of file Compadre_GMLS.hpp.

◆ _h_ss

SolutionSet<host_memory_space> Compadre::GMLS::_h_ss
private

Solution Set (contains all alpha values from solution and alpha layout methods)

Definition at line 114 of file Compadre_GMLS.hpp.

◆ _host_operations

Kokkos::View<TargetOperation*>::HostMirror Compadre::GMLS::_host_operations
private

vector containing target functionals to be applied for reconstruction problem (host)

Definition at line 167 of file Compadre_GMLS.hpp.

◆ _host_prestencil_weights

Kokkos::View<const double*****, layout_right>::HostMirror Compadre::GMLS::_host_prestencil_weights
private

generated weights for nontraditional samples required to transform data into expected sampling functional form (host)

Definition at line 106 of file Compadre_GMLS.hpp.

◆ _host_ref_N

Kokkos::View<double*>::HostMirror Compadre::GMLS::_host_ref_N
private

reference outward normal vectors information (host)

Definition at line 77 of file Compadre_GMLS.hpp.

◆ _host_T

Kokkos::View<double*>::HostMirror Compadre::GMLS::_host_T
private

tangent vectors information (host)

Definition at line 74 of file Compadre_GMLS.hpp.

◆ _initial_index_for_batch

int Compadre::GMLS::_initial_index_for_batch
private

initial index for current batch

Definition at line 220 of file Compadre_GMLS.hpp.

◆ _local_dimensions

int Compadre::GMLS::_local_dimensions
private

dimension of the problem, set at class instantiation only. For manifolds, generally _global_dimensions-1

Definition at line 130 of file Compadre_GMLS.hpp.

◆ _manifold_curvature_coefficients

Kokkos::View<double*> Compadre::GMLS::_manifold_curvature_coefficients
private

curvature polynomial coefficients for all problems

Definition at line 83 of file Compadre_GMLS.hpp.

◆ _manifold_curvature_gradient

Kokkos::View<double*> Compadre::GMLS::_manifold_curvature_gradient
private

_dimension-1 gradient values for curvature for all problems

Definition at line 86 of file Compadre_GMLS.hpp.

◆ _manifold_metric_tensor_inverse

Kokkos::View<double*> Compadre::GMLS::_manifold_metric_tensor_inverse
private

metric tensor inverse for all problems

Definition at line 80 of file Compadre_GMLS.hpp.

◆ _NP

int Compadre::GMLS::_NP
private

dimension of basis for polynomial reconstruction

Definition at line 124 of file Compadre_GMLS.hpp.

◆ _operations

Kokkos::View<TargetOperation*> Compadre::GMLS::_operations
private

vector containing target functionals to be applied for reconstruction problem (device)

Definition at line 164 of file Compadre_GMLS.hpp.

◆ _order_of_quadrature_points

int Compadre::GMLS::_order_of_quadrature_points
private

order of exact polynomial integration for quadrature rule

Definition at line 226 of file Compadre_GMLS.hpp.

◆ _orthonormal_tangent_space_provided

bool Compadre::GMLS::_orthonormal_tangent_space_provided
private

whether or not the orthonormal tangent directions were provided by the user.

If they are not, then for the case of calculations on manifolds, a GMLS approximation of the tangent space will be made and stored for use.

Definition at line 203 of file Compadre_GMLS.hpp.

◆ _P

Kokkos::View<double*> Compadre::GMLS::_P
private

P*sqrt(w) matrix for all problems.

Definition at line 56 of file Compadre_GMLS.hpp.

◆ _pc

point_connections_type Compadre::GMLS::_pc
private

connections between points and neighbors

Definition at line 95 of file Compadre_GMLS.hpp.

◆ _pm

ParallelManager Compadre::GMLS::_pm
private

determines scratch level spaces and is used to call kernels

Definition at line 223 of file Compadre_GMLS.hpp.

◆ _poly_order

int Compadre::GMLS::_poly_order
private

order of basis for polynomial reconstruction

Definition at line 118 of file Compadre_GMLS.hpp.

◆ _polynomial_sampling_functional

SamplingFunctional Compadre::GMLS::_polynomial_sampling_functional
private

polynomial sampling functional used to construct P matrix, set at GMLS class instantiation
NOTE: can only be set at object instantiation

Definition at line 153 of file Compadre_GMLS.hpp.

◆ _prestencil_weights

Kokkos::View<double*****, layout_right> Compadre::GMLS::_prestencil_weights
private

generated weights for nontraditional samples required to transform data into expected sampling functional form (device).


Definition at line 102 of file Compadre_GMLS.hpp.

◆ _problem_type

ProblemType Compadre::GMLS::_problem_type
private

problem type for GMLS problem, can also be set to STANDARD for normal or MANIFOLD for manifold problems
NOTE: can only be set at object instantiation

Definition at line 146 of file Compadre_GMLS.hpp.

◆ _qm

Quadrature Compadre::GMLS::_qm
private

manages and calculates quadrature

Definition at line 235 of file Compadre_GMLS.hpp.

◆ _quadrature_type

std::string Compadre::GMLS::_quadrature_type
private

quadrature rule type

Definition at line 232 of file Compadre_GMLS.hpp.

◆ _reconstruction_space

ReconstructionSpace Compadre::GMLS::_reconstruction_space
private

reconstruction space for GMLS problems, set at GMLS class instantiation

Definition at line 136 of file Compadre_GMLS.hpp.

◆ _reconstruction_space_rank

int Compadre::GMLS::_reconstruction_space_rank
private

actual rank of reconstruction basis

Definition at line 139 of file Compadre_GMLS.hpp.

◆ _ref_N

Kokkos::View<double*> Compadre::GMLS::_ref_N
private

Rank 2 tensor for high order approximation of tangent vectors for all problems.

First rank is for the target index, the second is for the spatial dimension (_dimensions)

Definition at line 71 of file Compadre_GMLS.hpp.

◆ _reference_outward_normal_direction_provided

bool Compadre::GMLS::_reference_outward_normal_direction_provided
private

whether or not the reference outward normal directions were provided by the user.

Definition at line 206 of file Compadre_GMLS.hpp.

◆ _RHS

Kokkos::View<double*> Compadre::GMLS::_RHS
private

sqrt(w)*Identity matrix for all problems, later holds polynomial coefficients for all problems

Definition at line 59 of file Compadre_GMLS.hpp.

◆ _sampling_multiplier

int Compadre::GMLS::_sampling_multiplier
private

actual dimension of the sampling functional e.g.

reconstruction of vector on a 2D manifold in 3D would have _basis_multiplier of 2 e.g. in 3D, a scalar will be 1, a vector will be 3, and a vector of reused scalars will be 1

Definition at line 194 of file Compadre_GMLS.hpp.

◆ _source_extra_data

Kokkos::View<double**, layout_right> Compadre::GMLS::_source_extra_data
private

Extra data available to basis functions (optional)

Definition at line 89 of file Compadre_GMLS.hpp.

◆ _store_PTWP_inv_PTW

bool Compadre::GMLS::_store_PTWP_inv_PTW
private

whether polynomial coefficients were requested to be stored (in a state not yet applied to data)

Definition at line 217 of file Compadre_GMLS.hpp.

◆ _T

Kokkos::View<double*> Compadre::GMLS::_T
private

Rank 3 tensor for high order approximation of tangent vectors for all problems.

First rank is for the target index, the second is for the local direction to the manifolds 0..(_dimensions-1) are tangent, _dimensions is the normal, and the third is for the spatial dimension (_dimensions)

Definition at line 67 of file Compadre_GMLS.hpp.

◆ _target_extra_data

Kokkos::View<double**, layout_right> Compadre::GMLS::_target_extra_data
private

Extra data available to target operations (optional)

Definition at line 92 of file Compadre_GMLS.hpp.

◆ _use_reference_outward_normal_direction_provided_to_orient_surface

bool Compadre::GMLS::_use_reference_outward_normal_direction_provided_to_orient_surface
private

whether or not to use reference outward normal directions to orient the surface in a manifold problem.

Definition at line 209 of file Compadre_GMLS.hpp.

◆ _w

Kokkos::View<double*> Compadre::GMLS::_w
private

contains weights for all problems

Definition at line 53 of file Compadre_GMLS.hpp.

◆ _weighting_n

int Compadre::GMLS::_weighting_n
private

second parameter to be used for weighting kernel

Definition at line 179 of file Compadre_GMLS.hpp.

◆ _weighting_p

int Compadre::GMLS::_weighting_p
private

first parameter to be used for weighting kernel

Definition at line 176 of file Compadre_GMLS.hpp.

◆ _weighting_type

WeightingFunctionType Compadre::GMLS::_weighting_type
private

weighting kernel type for GMLS

Definition at line 170 of file Compadre_GMLS.hpp.

◆ _Z

Kokkos::View<double*> Compadre::GMLS::_Z
private

stores evaluations of targets applied to basis

Definition at line 62 of file Compadre_GMLS.hpp.


The documentation for this class was generated from the following files: