Compadre
1.5.9
|
Generalized Moving Least Squares (GMLS) More...
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 40 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_right > | coordinates_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_type * | getNeighborLists () const |
Get neighbor list accessor. More... | |
decltype(_pc)* | getPointConnections () |
Get a view (device) of all point connection info. More... | |
neighbor_lists_type * | getAdditionalEvaluationIndices () 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_weights) | getPrestencilWeights () 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(_RHS) | getFullPolynomialCoefficientsBasis () 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 |
typedef Kokkos::View<double**, layout_right> Compadre::GMLS::coordinates_type |
Definition at line 53 of file Compadre_GMLS.hpp.
typedef NeighborLists<Kokkos::View<int*> > Compadre::GMLS::neighbor_lists_type |
Definition at line 51 of file Compadre_GMLS.hpp.
typedef PointConnections<Kokkos::View<double**, layout_right>, Kokkos::View<double**, layout_right>, NeighborLists<Kokkos::View<int*> > > Compadre::GMLS::point_connections_type |
Definition at line 49 of file Compadre_GMLS.hpp.
|
inline |
Maximal constructor, not intended for users.
Definition at line 389 of file Compadre_GMLS.hpp.
|
inline |
Maximal constructor, but with string arguments.
Definition at line 468 of file Compadre_GMLS.hpp.
|
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 481 of file Compadre_GMLS.hpp.
|
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 491 of file Compadre_GMLS.hpp.
|
inline |
Adds a target to the vector of target functional to be applied to the reconstruction.
Definition at line 1286 of file Compadre_GMLS.hpp.
|
inline |
Adds a vector of target functionals to the vector of target functionals already to be applied to the reconstruction.
Definition at line 1292 of file Compadre_GMLS.hpp.
|
inline |
Empties the vector of target functionals to apply to the reconstruction.
Definition at line 1298 of file Compadre_GMLS.hpp.
|
inline |
Check if GMLS solution set contains valid alpha values (has generateAlphas been called)
Definition at line 824 of file Compadre_GMLS.hpp.
const GMLSBasisData Compadre::GMLS::extractBasisData | ( | ) | const |
Get GMLS basis data.
Definition at line 529 of file Compadre_GMLS.cpp.
const GMLSSolutionData Compadre::GMLS::extractSolutionData | ( | ) | const |
Get GMLS solution data.
Definition at line 481 of file Compadre_GMLS.cpp.
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.
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 475 of file Compadre_GMLS.cpp.
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.
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 14 of file Compadre_GMLS.cpp.
|
inline |
(OPTIONAL) Get additional evaluation sites neighbor list-like accessor
Definition at line 706 of file Compadre_GMLS.hpp.
|
inline |
(OPTIONAL) Get a view (device) of all additional evaluation point connection info
Definition at line 711 of file Compadre_GMLS.hpp.
|
inline |
Get constraint type.
Definition at line 653 of file Compadre_GMLS.hpp.
|
inline |
Get basis order used for curvature reconstruction.
Definition at line 659 of file Compadre_GMLS.hpp.
|
inline |
Get the data sampling functional specified at instantiation (often the same as the polynomial sampling functional)
Definition at line 775 of file Compadre_GMLS.hpp.
|
inline |
Get dense solver type.
Definition at line 647 of file Compadre_GMLS.hpp.
|
inline |
Dimensions of quadrature points.
Definition at line 692 of file Compadre_GMLS.hpp.
|
inline |
Dimension of the GMLS problem, set only at class instantiation.
Definition at line 638 of file Compadre_GMLS.hpp.
Get a view (device) of all polynomial coefficients basis.
Definition at line 759 of file Compadre_GMLS.hpp.
|
inline |
Dimension of the GMLS problem's point data (spatial description of points in ambient space), set only at class instantiation.
Definition at line 641 of file Compadre_GMLS.hpp.
|
inline |
Local dimension of the GMLS problem (less than global dimension if on a manifold), set only at class instantiation.
Definition at line 644 of file Compadre_GMLS.hpp.
|
inline |
Get parameter for weighting kernel for curvature.
Definition at line 677 of file Compadre_GMLS.hpp.
|
inline |
Type for weighting kernel for curvature.
Definition at line 665 of file Compadre_GMLS.hpp.
|
inline |
Get neighbor list accessor.
Definition at line 698 of file Compadre_GMLS.hpp.
|
inlinestatic |
Returns number of neighbors needed for unisolvency for a given basis order and dimension.
Definition at line 524 of file Compadre_GMLS.hpp.
|
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.
Definition at line 512 of file Compadre_GMLS.hpp.
|
inline |
Number of quadrature points.
Definition at line 686 of file Compadre_GMLS.hpp.
|
inline |
Order of quadrature points.
Definition at line 689 of file Compadre_GMLS.hpp.
Get a view (device) of all point connection info.
Definition at line 703 of file Compadre_GMLS.hpp.
|
inline |
Returns (size of the basis used in instance's polynomial reconstruction) x (data input dimension)
Definition at line 608 of file Compadre_GMLS.hpp.
|
inline |
Returns 2D array size in memory on which coefficients are stored.
Definition at line 622 of file Compadre_GMLS.hpp.
|
inline |
Returns size of the basis used in instance's polynomial reconstruction.
Definition at line 616 of file Compadre_GMLS.hpp.
|
inline |
Get basis order used for reconstruction.
Definition at line 656 of file Compadre_GMLS.hpp.
|
inline |
Get the polynomial sampling functional specified at instantiation.
Definition at line 772 of file Compadre_GMLS.hpp.
|
inline |
Returns a stencil to transform data from its existing state into the input expected for some sampling functionals.
Definition at line 782 of file Compadre_GMLS.hpp.
|
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 754 of file Compadre_GMLS.hpp.
|
inline |
Get problem type.
Definition at line 650 of file Compadre_GMLS.hpp.
|
inline |
Type of quadrature points.
Definition at line 695 of file Compadre_GMLS.hpp.
|
inline |
Get the reconstruction space specified at instantiation.
Definition at line 778 of file Compadre_GMLS.hpp.
|
inline |
Get component of tangent or normal directions for manifold problems.
Definition at line 732 of file Compadre_GMLS.hpp.
Get a view (device) of all reference outward normal directions.
Definition at line 720 of file Compadre_GMLS.hpp.
Get solution set on device.
Definition at line 817 of file Compadre_GMLS.hpp.
Get solution set on host.
Definition at line 805 of file Compadre_GMLS.hpp.
|
inline |
Get component of tangent or normal directions for manifold problems.
Definition at line 723 of file Compadre_GMLS.hpp.
Get a view (device) of all tangent direction bundles.
Definition at line 717 of file Compadre_GMLS.hpp.
|
inline |
Get parameter for weighting kernel for GMLS problem.
Definition at line 668 of file Compadre_GMLS.hpp.
|
inline |
Type for weighting kernel for GMLS problem.
Definition at line 662 of file Compadre_GMLS.hpp.
Get a view (device) of all window sizes.
Definition at line 714 of file Compadre_GMLS.hpp.
|
inlinestaticprivate |
Parses a string to determine constraint type.
Definition at line 367 of file Compadre_GMLS.hpp.
|
inlinestaticprivate |
Parses a string to determine problem type.
Definition at line 354 of file Compadre_GMLS.hpp.
|
inlinestaticprivate |
Parses a string to determine solver type.
Definition at line 343 of file Compadre_GMLS.hpp.
|
inline |
Definition at line 840 of file Compadre_GMLS.hpp.
|
inline |
(OPTIONAL) Sets additional evaluation sites for each target site
Definition at line 876 of file Compadre_GMLS.hpp.
|
inline |
(OPTIONAL) Sets additional evaluation sites for each target site
Definition at line 885 of file Compadre_GMLS.hpp.
|
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 257 of file Compadre_GMLS.hpp.
|
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 287 of file Compadre_GMLS.hpp.
|
inlineprivate |
(OPTIONAL) Sets the additional target evaluation indices list information from compressed row format (if same view_type)
Definition at line 296 of file Compadre_GMLS.hpp.
|
inlineprivate |
(OPTIONAL) Sets the additional target evaluation indices list information from compressed row format (if different view_type)
Definition at line 308 of file Compadre_GMLS.hpp.
|
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 326 of file Compadre_GMLS.hpp.
|
inline |
Set constraint type.
Definition at line 1165 of file Compadre_GMLS.hpp.
|
inline |
Sets basis order to be used when reconstruction curvature.
Definition at line 1237 of file Compadre_GMLS.hpp.
|
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 1258 of file Compadre_GMLS.hpp.
|
inline |
Type for weighting kernel for curvature.
Definition at line 1200 of file Compadre_GMLS.hpp.
|
inline |
Type for weighting kernel for curvature.
Definition at line 1223 of file Compadre_GMLS.hpp.
|
inline |
Set dense solver type type.
Definition at line 1159 of file Compadre_GMLS.hpp.
|
inline |
Dimensions of quadrature points to use.
Definition at line 1274 of file Compadre_GMLS.hpp.
|
inline |
Sets neighbor list information from compressed row neighborhood lists data (if same view_type).
Definition at line 897 of file Compadre_GMLS.hpp.
|
inline |
Sets neighbor list information from compressed row neighborhood lists data (if different view_type).
Definition at line 907 of file Compadre_GMLS.hpp.
|
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 923 of file Compadre_GMLS.hpp.
|
inline |
Number quadrature points to use.
Definition at line 1268 of file Compadre_GMLS.hpp.
|
inline |
Sets basis order to be used when reconstructing any function.
Definition at line 1229 of file Compadre_GMLS.hpp.
|
inline |
Sets basic problem data (neighbor lists, source coordinates, and target coordinates)
Definition at line 849 of file Compadre_GMLS.hpp.
|
inline |
Sets basic problem data (neighbor lists data, number of neighbors list, source coordinates, and target coordinates)
Definition at line 862 of file Compadre_GMLS.hpp.
|
inline |
Type of quadrature points.
Definition at line 1280 of file Compadre_GMLS.hpp.
|
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 1085 of file Compadre_GMLS.hpp.
|
inline |
(OPTIONAL) Sets extra data to be used by sampling functionals in certain instances.
Definition at line 1113 of file Compadre_GMLS.hpp.
|
inline |
(OPTIONAL) Sets extra data to be used by sampling functionals in certain instances.
(device)
Definition at line 1128 of file Compadre_GMLS.hpp.
|
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 933 of file Compadre_GMLS.hpp.
|
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 962 of file Compadre_GMLS.hpp.
|
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 1048 of file Compadre_GMLS.hpp.
|
inline |
(OPTIONAL) Sets extra data to be used by target operations in certain instances.
Definition at line 1137 of file Compadre_GMLS.hpp.
|
inline |
(OPTIONAL) Sets extra data to be used by target operations in certain instances.
(device)
Definition at line 1152 of file Compadre_GMLS.hpp.
|
inline |
Sets target coordinate information. Rows of this 2D-array should correspond to rows of the neighbor lists.
Definition at line 970 of file Compadre_GMLS.hpp.
|
inline |
Sets target coordinate information. Rows of this 2D-array should correspond to rows of the neighbor lists.
Definition at line 1007 of file Compadre_GMLS.hpp.
|
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 1246 of file Compadre_GMLS.hpp.
|
inline |
Type for weighting kernel for GMLS problem.
Definition at line 1171 of file Compadre_GMLS.hpp.
|
inline |
Type for weighting kernel for GMLS problem.
Definition at line 1194 of file Compadre_GMLS.hpp.
|
inline |
Sets window sizes, also called the support of the kernel.
Definition at line 1024 of file Compadre_GMLS.hpp.
|
inline |
Sets window sizes, also called the support of the kernel (device)
Definition at line 1036 of file Compadre_GMLS.hpp.
|
inline |
Verify whether _additional_pc is valid.
Definition at line 1320 of file Compadre_GMLS.hpp.
|
inline |
Verify whether _pc is valid.
Definition at line 1304 of file Compadre_GMLS.hpp.
|
inlinestatic |
Evaluates the weighting kernel.
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 549 of file Compadre_GMLS.hpp.
|
friend |
Definition at line 42 of file Compadre_GMLS.hpp.
|
private |
(OPTIONAL) connections between additional points and neighbors
Definition at line 117 of file Compadre_GMLS.hpp.
|
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 197 of file Compadre_GMLS.hpp.
|
private |
constraint type for GMLS problem
Definition at line 157 of file Compadre_GMLS.hpp.
|
private |
order of basis for curvature reconstruction
Definition at line 129 of file Compadre_GMLS.hpp.
|
private |
vector containing target functionals to be applied for curvature
Definition at line 169 of file Compadre_GMLS.hpp.
|
private |
second parameter to be used for weighting kernel for curvature
Definition at line 193 of file Compadre_GMLS.hpp.
|
private |
first parameter to be used for weighting kernel for curvature
Definition at line 190 of file Compadre_GMLS.hpp.
|
private |
weighting kernel type for curvature problem
Definition at line 181 of file Compadre_GMLS.hpp.
|
private |
Definition at line 123 of file Compadre_GMLS.hpp.
|
private |
generally the same as _polynomial_sampling_functional, but can differ if specified at
GMLS class instantiation
Definition at line 166 of file Compadre_GMLS.hpp.
|
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 206 of file Compadre_GMLS.hpp.
|
private |
solver type for GMLS problem - can be QR, SVD or LU
Definition at line 150 of file Compadre_GMLS.hpp.
|
private |
dimension of quadrature rule
Definition at line 237 of file Compadre_GMLS.hpp.
|
private |
dimension of the problem, set at class instantiation only
Definition at line 141 of file Compadre_GMLS.hpp.
|
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 222 of file Compadre_GMLS.hpp.
|
private |
h supports determined through neighbor search (device)
Definition at line 106 of file Compadre_GMLS.hpp.
|
private |
spatial dimension of the points, set at class instantiation only
Definition at line 135 of file Compadre_GMLS.hpp.
|
private |
Solution Set (contains all alpha values from solution and alpha layout methods)
Definition at line 122 of file Compadre_GMLS.hpp.
|
private |
vector containing target functionals to be applied for reconstruction problem (host)
Definition at line 175 of file Compadre_GMLS.hpp.
|
private |
generated weights for nontraditional samples required to transform data into expected sampling functional form (host)
Definition at line 114 of file Compadre_GMLS.hpp.
|
private |
reference outward normal vectors information (host)
Definition at line 85 of file Compadre_GMLS.hpp.
|
private |
tangent vectors information (host)
Definition at line 82 of file Compadre_GMLS.hpp.
|
private |
initial index for current batch
Definition at line 228 of file Compadre_GMLS.hpp.
|
private |
dimension of the problem, set at class instantiation only. For manifolds, generally _global_dimensions-1
Definition at line 138 of file Compadre_GMLS.hpp.
|
private |
curvature polynomial coefficients for all problems
Definition at line 91 of file Compadre_GMLS.hpp.
|
private |
_dimension-1 gradient values for curvature for all problems
Definition at line 94 of file Compadre_GMLS.hpp.
|
private |
metric tensor inverse for all problems
Definition at line 88 of file Compadre_GMLS.hpp.
|
private |
dimension of basis for polynomial reconstruction
Definition at line 132 of file Compadre_GMLS.hpp.
|
private |
vector containing target functionals to be applied for reconstruction problem (device)
Definition at line 172 of file Compadre_GMLS.hpp.
|
private |
order of exact polynomial integration for quadrature rule
Definition at line 234 of file Compadre_GMLS.hpp.
|
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 211 of file Compadre_GMLS.hpp.
|
private |
P*sqrt(w) matrix for all problems.
Definition at line 64 of file Compadre_GMLS.hpp.
|
private |
connections between points and neighbors
Definition at line 103 of file Compadre_GMLS.hpp.
|
private |
determines scratch level spaces and is used to call kernels
Definition at line 231 of file Compadre_GMLS.hpp.
|
private |
order of basis for polynomial reconstruction
Definition at line 126 of file Compadre_GMLS.hpp.
|
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 161 of file Compadre_GMLS.hpp.
|
private |
generated weights for nontraditional samples required to transform data into expected sampling functional form (device).
Definition at line 110 of file Compadre_GMLS.hpp.
|
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 154 of file Compadre_GMLS.hpp.
|
private |
manages and calculates quadrature
Definition at line 243 of file Compadre_GMLS.hpp.
|
private |
quadrature rule type
Definition at line 240 of file Compadre_GMLS.hpp.
|
private |
reconstruction space for GMLS problems, set at GMLS class instantiation
Definition at line 144 of file Compadre_GMLS.hpp.
|
private |
actual rank of reconstruction basis
Definition at line 147 of file Compadre_GMLS.hpp.
|
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 79 of file Compadre_GMLS.hpp.
|
private |
whether or not the reference outward normal directions were provided by the user.
Definition at line 214 of file Compadre_GMLS.hpp.
|
private |
sqrt(w)*Identity matrix for all problems, later holds polynomial coefficients for all problems
Definition at line 67 of file Compadre_GMLS.hpp.
|
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 202 of file Compadre_GMLS.hpp.
|
private |
Extra data available to basis functions (optional)
Definition at line 97 of file Compadre_GMLS.hpp.
|
private |
whether polynomial coefficients were requested to be stored (in a state not yet applied to data)
Definition at line 225 of file Compadre_GMLS.hpp.
|
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 75 of file Compadre_GMLS.hpp.
|
private |
Extra data available to target operations (optional)
Definition at line 100 of file Compadre_GMLS.hpp.
|
private |
whether or not to use reference outward normal directions to orient the surface in a manifold problem.
Definition at line 217 of file Compadre_GMLS.hpp.
|
private |
contains weights for all problems
Definition at line 61 of file Compadre_GMLS.hpp.
|
private |
second parameter to be used for weighting kernel
Definition at line 187 of file Compadre_GMLS.hpp.
|
private |
first parameter to be used for weighting kernel
Definition at line 184 of file Compadre_GMLS.hpp.
|
private |
weighting kernel type for GMLS
Definition at line 178 of file Compadre_GMLS.hpp.
|
private |
stores evaluations of targets applied to basis
Definition at line 70 of file Compadre_GMLS.hpp.