1 #ifndef _COMPADRE_GMLS_HPP_ 2 #define _COMPADRE_GMLS_HPP_ 4 #include "Compadre_Config.h" 23 struct GMLSSolutionData;
39 Kokkos::View<double**, layout_right>,
53 Kokkos::View<double*>
_w;
56 Kokkos::View<double*>
_P;
59 Kokkos::View<double*>
_RHS;
62 Kokkos::View<double*>
_Z;
67 Kokkos::View<double*>
_T;
74 Kokkos::View<double*>::HostMirror
_host_T;
248 template <
typename view_type>
251 auto additional_evaluation_coordinates =
coordinates_type(
"device additional evaluation coordinates",
252 evaluation_coordinates.extent(0), evaluation_coordinates.extent(1));
254 typedef typename view_type::memory_space input_array_memory_space;
255 if (std::is_same<input_array_memory_space, device_memory_space>::value) {
260 Kokkos::deep_copy(additional_evaluation_coordinates, evaluation_coordinates);
265 auto host_additional_evaluation_coordinates = Kokkos::create_mirror_view(additional_evaluation_coordinates);
266 Kokkos::deep_copy(host_additional_evaluation_coordinates, evaluation_coordinates);
268 Kokkos::deep_copy(additional_evaluation_coordinates, host_additional_evaluation_coordinates);
278 template <
typename view_type>
286 template <
typename view_type>
287 typename std::enable_if<view_type::rank==1&&std::is_same<neighbor_lists_type::internal_view_type,view_type>::value==1,
void>::type
298 template <
typename view_type>
299 typename std::enable_if<view_type::rank==1&&std::is_same<neighbor_lists_type::internal_view_type,view_type>::value==0,
void>::type
303 gmls_view_type d_additional_evaluation_indices(
"compressed row additional evaluation indices lists data", additional_evaluation_indices.extent(0));
304 gmls_view_type d_number_of_neighbors_list(
"number of additional evaluation indices", number_of_neighbors_list.extent(0));
305 Kokkos::deep_copy(d_additional_evaluation_indices, additional_evaluation_indices);
306 Kokkos::deep_copy(d_number_of_neighbors_list, number_of_neighbors_list);
317 template <
typename view_type>
320 auto additional_nla = Convert2DToCompressedRowNeighborLists<decltype(additional_evaluation_indices), Kokkos::View<int*> >(additional_evaluation_indices);
336 std::string solver_type_to_lower = dense_solver_type;
337 transform(solver_type_to_lower.begin(), solver_type_to_lower.end(), solver_type_to_lower.begin(), ::tolower);
338 if (solver_type_to_lower ==
"lu") {
347 std::string problem_type_to_lower = problem_type;
348 transform(problem_type_to_lower.begin(), problem_type_to_lower.end(), problem_type_to_lower.begin(), ::tolower);
349 if (problem_type_to_lower ==
"standard") {
351 }
else if (problem_type_to_lower ==
"manifold") {
360 std::string constraint_type_to_lower = constraint_type;
361 transform(constraint_type_to_lower.begin(), constraint_type_to_lower.end(), constraint_type_to_lower.begin(), ::tolower);
362 if (constraint_type_to_lower ==
"none") {
364 }
else if (constraint_type_to_lower ==
"neumann_grad_scalar") {
384 const int poly_order,
385 const int dimensions,
389 const int manifold_curvature_poly_order) :
411 (
"operations needed for manifold gradient reconstruction", 1);
412 auto curvature_support_operations_mirror =
414 curvature_support_operations_mirror(0) =
463 const int poly_order,
464 const int dimensions = 3,
465 const std::string dense_solver_type = std::string(
"QR"),
466 const std::string problem_type = std::string(
"STANDARD"),
467 const std::string constraint_type = std::string(
"NO_CONSTRAINT"),
468 const int manifold_curvature_poly_order = 2)
474 const int dimensions = 3,
475 const std::string dense_solver_type = std::string(
"QR"),
476 const std::string problem_type = std::string(
"STANDARD"),
477 const std::string constraint_type = std::string(
"NO_CONSTRAINT"),
478 const int manifold_curvature_poly_order = 2)
485 const int poly_order,
486 const int dimensions = 3,
487 const std::string dense_solver_type = std::string(
"QR"),
488 const std::string problem_type = std::string(
"STANDARD"),
489 const std::string constraint_type = std::string(
"NO_CONSTRAINT"),
490 const int manifold_curvature_poly_order = 2)
491 :
GMLS(reconstruction_space, dual_sampling_strategy, dual_sampling_strategy, poly_order, dimensions, dense_solver_type, problem_type, constraint_type, manifold_curvature_poly_order) {}
503 KOKKOS_INLINE_FUNCTION
515 KOKKOS_INLINE_FUNCTION
518 const int np =
getNP(m, dimension, r_space);
522 nn = np * (1.7 + m*0.1);
525 nn = np * (1.4 + m*0.03);
540 KOKKOS_INLINE_FUNCTION
551 double abs_r_over_h_to_n = std::abs(r/h);
552 if (n>1) abs_r_over_h_to_n = std::pow(abs_r_over_h_to_n, n);
553 return std::pow(1.0-abs_r_over_h_to_n, p) * double(1.0-abs_r_over_h_to_n>0.0);
557 double x = std::abs(r/h);
558 return ((1-x)+x*(1-x)*(1-2*x)) * double(x<=1.0);
561 double pi = 3.14159265358979323846;
562 double abs_r_over_h_to_n = std::abs(r/h);
563 return std::cos(0.5*pi*r/h) * double(1.0-abs_r_over_h_to_n>0.0);
568 double h_over_p = h/p;
569 double abs_r_over_h_to_n = std::abs(r/h);
570 return double(1.0-abs_r_over_h_to_n>0.0)/( h_over_p * 2.5066282746310002416124 ) * std::exp(-.5*r*r/(h_over_p*h_over_p));
574 double abs_r_over_h_to_n = std::abs(r/h);
575 return double(1.0-abs_r_over_h_to_n>0.0) / (std::exp(p*r) + std::exp(-p*r) + n);
609 &&
"Entire batch not computed at once, so getFullPolynomialCoefficientsBasis() can not be called.");
611 &&
"generateAlphas() called with keep_coefficients set to false.");
707 double getTangentBundle(
const int target_index,
const int direction,
const int component)
const {
710 scratch_matrix_right_type::HostMirror
712 return T(direction, component);
718 "getRefenceNormalDirection called, but reference outwrad normal directions were never provided.");
719 scratch_vector_type::HostMirror
721 return ref_N(component);
745 &&
"Entire batch not computed at once, so getFullPolynomialCoefficientsBasis() can not be called.");
747 &&
"generateAlphas() called with keep_coefficients set to false.");
770 if (for_target)
return 0;
else return 1;
776 const int target_index_in_weights =
780 const int neighbor_index_in_weights =
785 output_component, input_component);
796 "getSolutionSetHost() called with invalid alpha values.");
803 "getSolutionSetDevice() called with invalid alpha values.");
825 if (
_RHS.extent(0) > 0)
826 _RHS = Kokkos::View<double*>(
"RHS",0);
832 template<
typename view_type_1,
typename view_type_2,
typename view_type_3,
typename view_type_4>
834 view_type_1 neighbor_lists,
835 view_type_2 source_coordinates,
836 view_type_3 target_coordinates,
837 view_type_4 epsilons) {
838 this->setNeighborLists<view_type_1>(neighbor_lists);
839 this->setSourceSites<view_type_2>(source_coordinates);
840 this->setTargetSites<view_type_3>(target_coordinates);
841 this->setWindowSizes<view_type_4>(epsilons);
845 template<
typename view_type_1,
typename view_type_2,
typename view_type_3,
typename view_type_4>
847 view_type_1 cr_neighbor_lists,
848 view_type_1 number_of_neighbors_list,
849 view_type_2 source_coordinates,
850 view_type_3 target_coordinates,
851 view_type_4 epsilons) {
852 this->setNeighborLists<view_type_1>(cr_neighbor_lists, number_of_neighbors_list);
853 this->setSourceSites<view_type_2>(source_coordinates);
854 this->setTargetSites<view_type_3>(target_coordinates);
855 this->setWindowSizes<view_type_4>(epsilons);
859 template<
typename view_type_1,
typename view_type_2>
861 view_type_1 additional_evaluation_indices,
862 view_type_2 additional_evaluation_coordinates) {
863 this->setAuxiliaryEvaluationIndicesLists<view_type_1>(additional_evaluation_indices);
864 this->setAuxiliaryEvaluationCoordinates<view_type_2>(additional_evaluation_coordinates);
868 template<
typename view_type_1,
typename view_type_2>
870 view_type_1 cr_additional_evaluation_indices,
871 view_type_1 number_of_additional_evaluation_indices,
872 view_type_2 additional_evaluation_coordinates) {
873 this->setAuxiliaryEvaluationIndicesLists<view_type_1>(cr_additional_evaluation_indices,
874 number_of_additional_evaluation_indices);
875 this->setAuxiliaryEvaluationCoordinates<view_type_2>(additional_evaluation_coordinates);
879 template <
typename view_type>
880 typename std::enable_if<view_type::rank==1&&std::is_same<neighbor_lists_type::internal_view_type,view_type>::value==1,
void>::type
889 template <
typename view_type>
890 typename std::enable_if<view_type::rank==1&&std::is_same<neighbor_lists_type::internal_view_type,view_type>::value==0,
void>::type
894 gmls_view_type d_neighbor_lists(
"compressed row neighbor lists data", neighbor_lists.extent(0));
895 gmls_view_type d_number_of_neighbors_list(
"number of neighbors list", number_of_neighbors_list.extent(0));
896 Kokkos::deep_copy(d_neighbor_lists, neighbor_lists);
897 Kokkos::deep_copy(d_number_of_neighbors_list, number_of_neighbors_list);
906 template <
typename view_type>
907 typename std::enable_if<view_type::rank==2, void>::type
setNeighborLists(view_type neighbor_lists) {
909 auto nla = Convert2DToCompressedRowNeighborLists<decltype(neighbor_lists), Kokkos::View<int*> >(neighbor_lists);
916 template<
typename view_type>
921 source_coordinates.extent(0), source_coordinates.extent(1));
923 typedef typename view_type::memory_space input_array_memory_space;
924 if (std::is_same<input_array_memory_space, device_memory_space>::value) {
929 Kokkos::deep_copy(sc, source_coordinates);
934 auto host_source_coordinates = Kokkos::create_mirror_view(sc);
935 Kokkos::deep_copy(host_source_coordinates, source_coordinates);
937 Kokkos::deep_copy(sc, host_source_coordinates);
945 template<
typename view_type>
953 template<
typename view_type>
957 target_coordinates.extent(0), target_coordinates.extent(1));
959 typedef typename view_type::memory_space input_array_memory_space;
960 if (std::is_same<input_array_memory_space, device_memory_space>::value) {
965 Kokkos::deep_copy(tc, target_coordinates);
970 auto host_target_coordinates = Kokkos::create_mirror_view(tc);
971 Kokkos::deep_copy(host_target_coordinates, target_coordinates);
973 Kokkos::deep_copy(tc, host_target_coordinates);
977 Kokkos::View<int*>(),
978 Kokkos::View<int*>(
"number of additional evaluation indices",
979 target_coordinates.extent(0))
990 template<
typename view_type>
994 Kokkos::View<int*>(),
995 Kokkos::View<int*>(
"number of additional evaluation indices",
996 target_coordinates.extent(0))
1007 template<
typename view_type>
1019 template<
typename view_type>
1031 template<
typename view_type>
1041 "Memory space does not match between _T and tangent_directions");
1046 Kokkos::parallel_for(
"copy tangent vectors", Kokkos::RangePolicy<device_execution_space>(0,
_pc.
_target_coordinates.extent(0)), KOKKOS_LAMBDA(
const int i) {
1048 for (
int j=0; j<this_dimensions; ++j) {
1049 for (
int k=0; k<this_dimensions; ++k) {
1050 T(j,k) = tangent_directions(i, j, k);
1057 _host_T = Kokkos::create_mirror_view(
_T);
1065 template<
typename view_type>
1072 auto this_ref_N = this->
_ref_N;
1076 Kokkos::parallel_for(
"copy normal vectors", Kokkos::RangePolicy<device_execution_space>(0,
_pc.
_target_coordinates.extent(0)), KOKKOS_LAMBDA(
const int i) {
1077 for (
int j=0; j<this_dimensions; ++j) {
1078 this_ref_N(i*this_dimensions + j) = outward_normal_directions(i, j);
1093 template<
typename view_type>
1100 Kokkos::deep_copy(host_extra_data, extra_data);
1108 template<
typename view_type>
1117 template<
typename view_type>
1124 Kokkos::deep_copy(host_extra_data, extra_data);
1132 template<
typename view_type>
1153 std::string wt_to_lower = wt;
1154 transform(wt_to_lower.begin(), wt_to_lower.end(), wt_to_lower.begin(), ::tolower);
1155 if (wt_to_lower ==
"power") {
1157 }
else if (wt_to_lower ==
"gaussian") {
1159 }
else if (wt_to_lower ==
"cubicspline") {
1161 }
else if (wt_to_lower ==
"cosine") {
1163 }
else if (wt_to_lower ==
"sigmoid") {
1180 std::string wt_to_lower = wt;
1181 transform(wt_to_lower.begin(), wt_to_lower.end(), wt_to_lower.begin(), ::tolower);
1182 if (wt_to_lower ==
"power") {
1184 }
else if (wt_to_lower ==
"gaussian") {
1186 }
else if (wt_to_lower ==
"cubicspline") {
1188 }
else if (wt_to_lower ==
"cosine") {
1190 }
else if (wt_to_lower ==
"sigmoid") {
1284 "Target coordinates and neighbor lists have different size.");
1288 "Source coordinates and target coordinates have different dimensions.");
1292 "Source coordinates not set in GMLS class before calling generatePolynomialCoefficients.");
1300 "Target coordinates and additional evaluation indices have different size.");
1304 "Target coordinates and additional evaluation indices have different size.");
1309 "Target coordinates and additional evaluation coordinates have different dimensions.");
1329 const bool clear_cache =
true);
1344 void generateAlphas(
const int number_of_batches = 1,
const bool keep_coefficients =
false,
1345 const bool clear_cache =
true);
Kokkos::View< double * > _manifold_metric_tensor_inverse
metric tensor inverse for all problems
Divergence-free vector polynomial basis.
decltype(_RHS) getFullPolynomialCoefficientsBasis() const
Get a view (device) of all polynomial coefficients basis.
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 ...
KOKKOS_INLINE_FUNCTION int getMaxNumNeighbors() const
Get the maximum number of neighbors of all targets' neighborhoods (host/device)
Standard GMLS problem type.
const GMLSSolutionData extractSolutionData() const
Get GMLS solution data.
Kokkos::View< const double *****, layout_right >::HostMirror _host_prestencil_weights
generated weights for nontraditional samples required to transform data into expected sampling functi...
Kokkos::View< TargetOperation * > _operations
vector containing target functionals to be applied for reconstruction problem (device) ...
void setCurvatureWeightingParameter(int wp, int index=0)
Parameter for weighting kernel for curvature index = 0 sets p paramater for weighting kernel index = ...
bool _entire_batch_computed_at_once
whether entire calculation was computed at once the alternative is that it was broken up over many sm...
void setNeighborLists(nla_type nla)
Update only target coordinates.
int _weighting_p
first parameter to be used for weighting kernel
Kokkos::View< double * > _manifold_curvature_coefficients
curvature polynomial coefficients for all problems
SolutionSet< host_memory_space > _h_ss
Solution Set (contains all alpha values from solution and alpha layout methods)
decltype(_additional_pc) * getAdditionalPointConnections()
(OPTIONAL) Get a view (device) of all additional evaluation point connection info ...
Kokkos::View< double * > _ref_N
Rank 2 tensor for high order approximation of tangent vectors for all problems.
KOKKOS_INLINE_FUNCTION int getSize(const int degree, const int dimension)
Returns size of basis.
Lightweight Evaluator Helper This class is a lightweight wrapper for extracting and applying all rele...
#define compadre_kernel_assert_release(condition)
compadre_kernel_assert_release is similar to compadre_assert_release, but is a call on the device...
constexpr SamplingFunctional ManifoldVectorPointSample
Point evaluations of the entire vector source function (but on a manifold, so it includes a transform...
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 fun...
bool _reference_outward_normal_direction_provided
whether or not the reference outward normal directions were provided by the user. ...
int _curvature_weighting_p
first parameter to be used for weighting kernel for curvature
int getOrderOfQuadraturePoints() const
Order of quadrature points.
void setCurvaturePolynomialOrder(const int curvature_poly_order)
Sets basis order to be used when reconstruction curvature.
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.
ConstraintType getConstraintType() const
Get constraint type.
Neumann Gradient Scalar Type.
int _NP
dimension of basis for polynomial reconstruction
Kokkos::View< double * > _epsilons
h supports determined through neighbor search (device)
bool _orthonormal_tangent_space_provided
whether or not the orthonormal tangent directions were provided by the user.
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)
void setQuadratureType(std::string quadrature_type)
Type of quadrature points.
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 ...
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...
Kokkos::View< double * > _RHS
sqrt(w)*Identity matrix for all problems, later holds polynomial coefficients for all problems ...
device_mirror_source_view_type _source_coordinates
Kokkos::View< double * > _w
contains weights for all problems
void setWindowSizes(view_type epsilons)
Sets window sizes, also called the support of the kernel.
int _sampling_multiplier
actual dimension of the sampling functional e.g.
int _reconstruction_space_rank
actual rank of reconstruction basis
int _curvature_poly_order
order of basis for curvature reconstruction
Kokkos::View< int * > internal_view_type
decltype(_d_ss) * getSolutionSetDevice(bool alpha_validity_check=true)
Get solution set on device.
Each target applies a different data transform, but the same to each neighbor.
SamplingFunctional _polynomial_sampling_functional
polynomial sampling functional used to construct P matrix, set at GMLS class instantiation NOTE: ca...
std::enable_if< view_type::rank==2, void >::type setNeighborLists(view_type neighbor_lists)
Sets neighbor list information.
int _data_sampling_multiplier
effective dimension of the data sampling functional e.g.
KOKKOS_INLINE_FUNCTION int getActualReconstructionSpaceRank(const int &index)
Number of actual components in the ReconstructionSpace.
int getNumberOfQuadraturePoints() const
Number of quadrature points.
Each target applies a different transform for each neighbor.
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...
Quadrature _qm
manages and calculates quadrature
Kokkos::View< double * >::HostMirror _host_T
tangent vectors information (host)
void clearTargets()
Empties the vector of target functionals to apply to the reconstruction.
void setTargetExtraData(view_type extra_data)
(OPTIONAL) Sets extra data to be used by target operations in certain instances.
void setReferenceOutwardNormalDirection(view_type outward_normal_directions, bool use_to_orient_surface=true)
(OPTIONAL) Sets outward normal direction.
void setDimensionOfQuadraturePoints(int dim)
Dimensions of quadrature points to use.
Scalar polynomial basis centered at the target site and scaled by sum of basis powers e...
void addTargets(TargetOperation lro)
Adds a target to the vector of target functional to be applied to the reconstruction.
Solve GMLS problem on a manifold (will use QR or SVD to solve the resultant GMLS problem dependent on...
SamplingFunctional getPolynomialSamplingFunctional() const
Get the polynomial sampling functional specified at instantiation.
int getLocalDimensions() const
Local dimension of the GMLS problem (less than global dimension if on a manifold), set only at class instantiation.
WeightingFunctionType _curvature_weighting_type
weighting kernel type for curvature problem
WeightingFunctionType _weighting_type
weighting kernel type for GMLS
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.
int _weighting_n
second parameter to be used for weighting kernel
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.
void clearTargets()
Empties the vector of target functionals to apply to the reconstruction.
bool _contains_valid_alphas
whether internal alpha values are valid (set externally on a solve)
static ConstraintType parseConstraintType(const std::string &constraint_type)
Parses a string to determine constraint type.
bool verifyPointConnections(bool assert_valid=false)
Verify whether _pc is valid.
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 ...
DenseSolverType _dense_solver_type
solver type for GMLS problem - can be QR, SVD or LU
static ProblemType parseProblemType(const std::string &problem_type)
Parses a string to determine problem type.
QR+Pivoting factorization performed on P*sqrt(w) matrix.
void addTargets(TargetOperation lro)
Adds a target to the vector of target functional to be applied to the reconstruction.
Kokkos::View< double * > _manifold_curvature_gradient
_dimension-1 gradient values for curvature for all problems
void setSourceExtraData(view_type extra_data)
(OPTIONAL) Sets extra data to be used by sampling functionals in certain instances.
ProblemType
Problem type, that optionally can handle manifolds.
Kokkos::View< int *, host_execution_space > host_managed_local_index_type
device_mirror_target_view_type _target_coordinates
Kokkos::View< double *****, layout_right > _prestencil_weights
generated weights for nontraditional samples required to transform data into expected sampling functi...
void addTargets(std::vector< TargetOperation > lro)
Adds a vector of target functionals to the vector of target functionals already to be applied to the ...
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)...
KOKKOS_INLINE_FUNCTION int getSize(const int degree, const int dimension)
Returns size of basis.
void setCurvatureWeightingType(const std::string &wt)
Type for weighting kernel for curvature.
int _initial_index_for_batch
initial index for current batch
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...
SamplingFunctional getDataSamplingFunctional() const
Get the data sampling functional specified at instantiation (often the same as the polynomial samplin...
Kokkos::View< double **, layout_right > _source_extra_data
Extra data available to basis functions (optional)
Scalar basis reused as many times as there are components in the vector resulting in a much cheaper p...
bool containsValidAlphas() const
Check if GMLS solution set contains valid alpha values (has generateAlphas been called) ...
void setAdditionalEvaluationSitesData(view_type_1 additional_evaluation_indices, view_type_2 additional_evaluation_coordinates)
(OPTIONAL) Sets additional evaluation sites for each target site
std::string getQuadratureType() const
Type of quadrature points.
Bernstein polynomial basis.
ReconstructionSpace
Space in which to reconstruct polynomial.
neighbor_lists_type * getAdditionalEvaluationIndices() const
(OPTIONAL) Get additional evaluation sites neighbor list-like accessor
void copyAlphas(SolutionSet< other_memory_space > &other)
Copies alphas between two instances of SolutionSet Copying of alphas is intentionally omitted in copy...
NeighborLists< Kokkos::View< int * > > neighbor_lists_type
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 use...
KOKKOS_INLINE_FUNCTION void getPDims(DenseSolverType dense_solver_type, ConstraintType constraint_type, ReconstructionSpace reconstruction_space, const int dimension, const int M, const int N, int &out_row, int &out_col)
host_managed_local_index_type getPolynomialCoefficientsMemorySize() const
Returns 2D array size in memory on which coefficients are stored.
Kokkos::View< double **, layout_right, Kokkos::MemoryTraits< Kokkos::Unmanaged > > scratch_matrix_right_type
const GMLSBasisData extractBasisData() const
Get GMLS basis data.
KOKKOS_INLINE_FUNCTION int getNumberOfQuadraturePoints() const
void setWeightingType(const WeightingFunctionType wt)
Type for weighting kernel for GMLS problem.
void setSourceSites(coordinates_type source_coordinates)
Sets source coordinate information.
void setDenseSolverType(const DenseSolverType dst)
Set dense solver type type.
WeightingFunctionType
Available weighting kernel function types.
int getManifoldWeightingParameter(const int index=0) const
Get parameter for weighting kernel for curvature.
point_connections_type _additional_pc
(OPTIONAL) connections between additional points and neighbors
Kokkos::View< double **, layout_right > _target_extra_data
Extra data available to target operations (optional)
int _dimension_of_quadrature_points
dimension of quadrature rule
int getPolynomialOrder() const
Get basis order used for reconstruction.
double getReferenceNormalDirection(const int target_index, const int component) const
Get component of tangent or normal directions for manifold problems.
Kokkos::View< double * >::HostMirror _host_ref_N
reference outward normal vectors information (host)
SamplingFunctional _data_sampling_functional
generally the same as _polynomial_sampling_functional, but can differ if specified at ...
Point evaluation of the gradient of a scalar.
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.
Kokkos::View< double **, layout_right > coordinates_type
Combines NeighborLists with the PointClouds from which it was derived Assumed that memory_space is th...
ParallelManager _pm
determines scratch level spaces and is used to call kernels
host_managed_local_index_type getPolynomialCoefficientsDomainRangeSize() const
Returns (size of the basis used in instance's polynomial reconstruction) x (data input dimension) ...
Kokkos::View< double * > _Z
stores evaluations of targets applied to basis
void setWindowSizes(decltype(_epsilons) epsilons)
Sets window sizes, also called the support of the kernel (device)
int getDimensions() const
Dimension of the GMLS problem, set only at class instantiation.
int getPolynomialCoefficientsSize() const
Returns size of the basis used in instance's polynomial reconstruction.
Kokkos::View< TargetOperation * > _curvature_support_operations
vector containing target functionals to be applied for curvature
Generalized Moving Least Squares (GMLS)
void setTargetExtraData(decltype(_target_extra_data) extra_data)
(OPTIONAL) Sets extra data to be used by target operations in certain instances.
ProblemType _problem_type
problem type for GMLS problem, can also be set to STANDARD for normal or MANIFOLD for manifold proble...
int _basis_multiplier
dimension of the reconstructed function e.g.
decltype(_T) * getTangentDirections()
Get a view (device) of all tangent direction bundles.
void setTargetSites(coordinates_type target_coordinates)
Sets target coordinate information. Rows of this 2D-array should correspond to rows of the neighbor l...
WeightingFunctionType getWeightingType() const
Type for weighting kernel for GMLS problem.
int _global_dimensions
spatial dimension of the points, set at class instantiation only
SolutionSet< device_memory_space > _d_ss
DenseSolverType
Dense solver type.
Kokkos::View< double * > _T
Rank 3 tensor for high order approximation of tangent vectors for all problems.
double getTangentBundle(const int target_index, const int direction, const int component) const
Get component of tangent or normal directions for manifold problems.
decltype(_epsilons) * getWindowSizes()
Get a view (device) of all window sizes.
void resetCoefficientData()
void setWeightingParameter(int wp, int index=0)
Parameter for weighting kernel for GMLS problem index = 0 sets p paramater for weighting kernel index...
int getCurvaturePolynomialOrder() const
Get basis order used for curvature reconstruction.
int getDimensionOfQuadraturePoints() const
Dimensions of quadrature points.
ReconstructionSpace getReconstructionSpace() const
Get the reconstruction space specified at instantiation.
void setTargetCoordinates(view_type_1 target_coordinates)
Update only target coordinates.
void setSourceSites(view_type source_coordinates)
Sets source coordinate information.
KOKKOS_INLINE_FUNCTION int getSize(const int degree, const int dimension)
Returns size of basis.
bool verifyAdditionalPointConnections(bool assert_valid=false)
Verify whether _additional_pc is valid.
LU factorization performed on P^T*W*P matrix.
PointConnections< Kokkos::View< double **, layout_right >, Kokkos::View< double **, layout_right >, NeighborLists< Kokkos::View< int * > > > point_connections_type
bool _store_PTWP_inv_PTW
whether polynomial coefficients were requested to be stored (in a state not yet applied to data) ...
int transform_type
Describes the SamplingFunction relationship to targets, neighbors.
void setPolynomialOrder(const int poly_order)
Sets basis order to be used when reconstructing any function.
void setTargetSites(view_type target_coordinates)
Sets target coordinate information. Rows of this 2D-array should correspond to rows of the neighbor l...
ReconstructionSpace _reconstruction_space
reconstruction space for GMLS problems, set at GMLS class instantiation
void setAuxiliaryEvaluationCoordinates(view_type evaluation_coordinates)
(OPTIONAL) Sets additional points for evaluation of target operation on polynomial reconstruction...
Kokkos::View< double * > _P
P*sqrt(w) matrix for all problems.
int _curvature_weighting_n
second parameter to be used for weighting kernel for curvature
void setWeightingType(const std::string &wt)
Type for weighting kernel for GMLS problem.
KOKKOS_INLINE_FUNCTION void getRHSDims(DenseSolverType dense_solver_type, ConstraintType constraint_type, ReconstructionSpace reconstruction_space, const int dimension, const int M, const int N, int &RHS_row, int &RHS_col)
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) ...
constexpr SamplingFunctional PointSample
Available sampling functionals.
int _dimensions
dimension of the problem, set at class instantiation only
decltype(_pc) * getPointConnections()
Get a view (device) of all point connection info.
std::string _quadrature_type
quadrature rule type
DenseSolverType getDenseSolverType() const
Get dense solver type.
int getWeightingParameter(const int index=0) const
Get parameter for weighting kernel for GMLS problem.
Kokkos::View< TargetOperation * >::HostMirror _host_operations
vector containing target functionals to be applied for reconstruction problem (host) ...
decltype(_ref_N) * getReferenceNormalDirections()
Get a view (device) of all reference outward normal directions.
void setConstraintType(const ConstraintType ct)
Set constraint type.
#define compadre_assert_debug(condition)
compadre_assert_debug is used for assertions that are checked in loops, as these significantly impact...
void setAuxiliaryEvaluationCoordinates(coordinates_type evaluation_coordinates)
(OPTIONAL) Sets additional points for evaluation of target operation on polynomial reconstruction...
void setCurvatureWeightingType(const WeightingFunctionType wt)
Type for weighting kernel for curvature.
int getGlobalDimensions() const
Dimension of the GMLS problem's point data (spatial description of points in ambient space)...
decltype(_h_ss) * getSolutionSetHost(bool alpha_validity_check=true)
Get solution set on host.
point_connections_type _pc
connections between points and neighbors
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
ConstraintType _constraint_type
constraint type for GMLS problem
ProblemType getProblemType() const
Get problem type.
int _local_dimensions
dimension of the problem, set at class instantiation only. For manifolds, generally _global_dimension...
NeighborLists assists in accessing entries of compressed row neighborhood lists.
int _max_evaluation_sites_per_target
maximum number of evaluation sites for each target (includes target site)
int _order_of_quadrature_points
order of exact polynomial integration for quadrature rule
std::enable_if< view_type::rank==2, void >::type setAuxiliaryEvaluationIndicesLists(view_type additional_evaluation_indices)
(OPTIONAL) Sets the additional target evaluation indices list information.
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 polynomi...
#define compadre_assert_release(condition)
compadre_assert_release is used for assertions that should always be checked, but generally are not e...
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 provi...
WeightingFunctionType getManifoldWeightingType() const
Type for weighting kernel for curvature.
void setSourceCoordinates(view_type_2 source_coordinates)
Update only source coordinates.
constexpr SamplingFunctional VectorPointSample
Point evaluations of the entire vector source function.
void setOrderOfQuadraturePoints(int order)
Number quadrature points to use.
static DenseSolverType parseSolverType(const std::string &dense_solver_type)
Parses a string to determine solver type.
int _poly_order
order of basis for polynomial reconstruction
TargetOperation
Available target functionals.
void setSourceExtraData(decltype(_source_extra_data) extra_data)
(OPTIONAL) Sets extra data to be used by sampling functionals in certain instances.
void setTangentBundle(view_type tangent_directions)
(OPTIONAL) Sets orthonormal tangent directions for reconstruction on a manifold.
neighbor_lists_type * getNeighborLists() const
Get neighbor list accessor.
ConstraintType
Constraint 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)...