PolyFEM
Loading...
Searching...
No Matches
RBFWithLinear.hpp
Go to the documentation of this file.
1#pragma once
2
4#include <Eigen/Dense>
5
6namespace polyfem
7{
8 namespace basis
9 {
11 {
12 public:
23 RBFWithLinear(const Eigen::MatrixXd &centers, const Eigen::MatrixXd &collocation_points,
24 const Eigen::MatrixXd &local_basis_integral, const quadrature::Quadrature &quadr,
25 Eigen::MatrixXd &rhs, bool with_constraints = true);
26
34 void basis(const int local_index, const Eigen::MatrixXd &uv, Eigen::MatrixXd &val) const;
35
43 void grad(const int local_index, const Eigen::MatrixXd &uv, Eigen::MatrixXd &val) const;
44
51 void bases_values(const Eigen::MatrixXd &samples, Eigen::MatrixXd &val) const;
52
60 void bases_grads(const int axis, const Eigen::MatrixXd &samples, Eigen::MatrixXd &val) const;
61
62 private:
63 bool is_volume() const { return centers_.cols() == 3; }
64
65 // Computes the matrix that evaluates the kernels + polynomial terms on the given sample points
66 void compute_kernels_matrix(const Eigen::MatrixXd &samples, Eigen::MatrixXd &A) const;
67
68 // Computes the relationship w = L v + t between the unknowns (v) and the weights w
69 void compute_constraints_matrix(const int num_bases, const quadrature::Quadrature &quadr,
70 const Eigen::MatrixXd &local_basis_integral, Eigen::MatrixXd &L, Eigen::MatrixXd &t) const;
71
72 // Computes the weights by solving a (possibly constrained) linear least square
73 void compute_weights(const Eigen::MatrixXd &collocation_points,
74 const Eigen::MatrixXd &local_basis_integral, const quadrature::Quadrature &quadr,
75 Eigen::MatrixXd &rhs, bool with_constraints);
76
77 private:
78 // #C x dim matrix of kernel center positions
79 Eigen::MatrixXd centers_;
80
81 // (#C + dim + 1) x #B matrix of weights extending the #B bases that are non-vanishing on the polytope
82 Eigen::MatrixXd weights_;
83 };
84 } // namespace basis
85} // namespace polyfem
double val
Definition Assembler.cpp:86
void compute_weights(const Eigen::MatrixXd &collocation_points, const Eigen::MatrixXd &local_basis_integral, const quadrature::Quadrature &quadr, Eigen::MatrixXd &rhs, bool with_constraints)
void bases_values(const Eigen::MatrixXd &samples, Eigen::MatrixXd &val) const
Batch evaluates the RBF + polynomials on a set of sample points.
void compute_constraints_matrix(const int num_bases, const quadrature::Quadrature &quadr, const Eigen::MatrixXd &local_basis_integral, Eigen::MatrixXd &L, Eigen::MatrixXd &t) const
void compute_kernels_matrix(const Eigen::MatrixXd &samples, Eigen::MatrixXd &A) const
void grad(const int local_index, const Eigen::MatrixXd &uv, Eigen::MatrixXd &val) const
Evaluates the gradient of one RBF function over a list of coordinates.
void basis(const int local_index, const Eigen::MatrixXd &uv, Eigen::MatrixXd &val) const
Evaluates one RBF function over a list of coordinates.
void bases_grads(const int axis, const Eigen::MatrixXd &samples, Eigen::MatrixXd &val) const
Batch evaluates the gradient of the RBF + polynomials on a set of sample points.