PolyFEM
|
#include <RBFWithQuadraticLagrange.hpp>
Public Member Functions | |
RBFWithQuadraticLagrange (const assembler::LinearAssembler &assembler, const Eigen::MatrixXd ¢ers, const Eigen::MatrixXd &collocation_points, const Eigen::MatrixXd &local_basis_integral, const quadrature::Quadrature &quadr, Eigen::MatrixXd &rhs, bool with_constraints=true) | |
Initialize RBF functions over a polytope element. | |
void | basis (const int local_index, const Eigen::MatrixXd &uv, Eigen::MatrixXd &val) const |
Evaluates one RBF function over a list of coordinates. | |
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 | bases_values (const Eigen::MatrixXd &samples, Eigen::MatrixXd &val) const |
Batch evaluates the RBF + polynomials on a set of sample points. | |
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. | |
Private Member Functions | |
bool | is_volume () const |
void | compute_kernels_matrix (const Eigen::MatrixXd &samples, Eigen::MatrixXd &A) const |
void | compute_constraints_matrix_2d_old (const int num_bases, const quadrature::Quadrature &quadr, Eigen::MatrixXd &C) const |
void | compute_constraints_matrix_2d (const assembler::LinearAssembler &assembler, const int num_bases, const quadrature::Quadrature &quadr, Eigen::MatrixXd &C) const |
void | compute_constraints_matrix_3d (const int num_bases, const quadrature::Quadrature &quadr, Eigen::MatrixXd &C) const |
void | compute_weights (const assembler::LinearAssembler &assembler, const Eigen::MatrixXd &collocation_points, const Eigen::MatrixXd &local_basis_integral, const quadrature::Quadrature &quadr, Eigen::MatrixXd &rhs, bool with_constraints) |
Private Attributes | |
Eigen::MatrixXd | centers_ |
Eigen::MatrixXd | weights_ |
Definition at line 15 of file RBFWithQuadraticLagrange.hpp.
RBFWithQuadraticLagrange::RBFWithQuadraticLagrange | ( | const assembler::LinearAssembler & | assembler, |
const Eigen::MatrixXd & | centers, | ||
const Eigen::MatrixXd & | collocation_points, | ||
const Eigen::MatrixXd & | local_basis_integral, | ||
const quadrature::Quadrature & | quadr, | ||
Eigen::MatrixXd & | rhs, | ||
bool | with_constraints = true |
||
) |
Initialize RBF functions over a polytope element.
[in] | centers | #C x dim positions of the kernels used to define functions over the polytope. The centers are placed at a small offset distance from the boundary of the element, due to the singularity at the centers |
[in] | collocation_points | #S x dim positions of the collocation points, used to approximate the RBF functions over the boundary of the element |
[in] | local_basis_integral | #B x dim+dim*(dim+1)/2 of the constant right-hand side for the integral constraint for each basis over the polytope |
[in] | quadr | Quadrature points and weights inside the polytope |
[in] | rhs | #S x #B of boundary conditions. Each column defines how the i-th basis of the mesh should evaluate on the collocation points sampled on the boundary of the polytope |
[in] | with_constraints | Impose integral constraints to guarantee linear reproduction for the Poisson equation |
Definition at line 79 of file RBFWithQuadraticLagrange.cpp.
References compute_weights().
void RBFWithQuadraticLagrange::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.
[in] | axis | The axis (0, 1, 2) with respect to which to compute the gradient |
[in] | uv | #uv x dim matrix of points to evaluate |
[out] | val | #uv x n_loc_bases of bases gradient wrt axis over the sample points |
Definition at line 130 of file RBFWithQuadraticLagrange.cpp.
References centers_, is_volume(), val, weights_, and x.
Referenced by grad().
void RBFWithQuadraticLagrange::bases_values | ( | const Eigen::MatrixXd & | samples, |
Eigen::MatrixXd & | val | ||
) | const |
Batch evaluates the RBF + polynomials on a set of sample points.
[in] | uv | #uv x dim matrix of points to evaluate |
[out] | val | #uv x n_loc_bases of bases values over the sample points |
Definition at line 118 of file RBFWithQuadraticLagrange.cpp.
References compute_kernels_matrix(), val, and weights_.
Referenced by basis().
void RBFWithQuadraticLagrange::basis | ( | const int | local_index, |
const Eigen::MatrixXd & | uv, | ||
Eigen::MatrixXd & | val | ||
) | const |
Evaluates one RBF function over a list of coordinates.
[in] | local_index | i-th RBF function to evaluate |
[in] | uv | #uv x dim matrix of coordinates to evaluate (in object domain) |
[out] | val | #uv x 1 matrix of computed values |
Definition at line 95 of file RBFWithQuadraticLagrange.cpp.
References bases_values(), and val.
Referenced by compute_weights().
|
private |
Definition at line 264 of file RBFWithQuadraticLagrange.cpp.
References polyfem::assembler::LinearAssembler::assemble(), polyfem::assembler::ElementAssemblyValues::basis_values, centers_, polyfem::assembler::ElementAssemblyValues::has_parameterization, polyfem::basis::RBFWithQuadratic::index_mapping(), polyfem::assembler::Assembler::is_tensor(), is_volume(), polyfem::quadrature::Quadrature::points, polyfem::basis::RBFWithQuadratic::setup_monomials_strong_2d(), polyfem::basis::RBFWithQuadratic::setup_monomials_vals_2d(), and polyfem::quadrature::Quadrature::weights.
Referenced by compute_weights().
|
private |
Definition at line 194 of file RBFWithQuadraticLagrange.cpp.
References centers_, is_volume(), polyfem::quadrature::Quadrature::points, polyfem::utils::show_matrix_stats(), and polyfem::quadrature::Quadrature::weights.
|
private |
Definition at line 371 of file RBFWithQuadraticLagrange.cpp.
References centers_, is_volume(), polyfem::quadrature::Quadrature::points, polyfem::utils::show_matrix_stats(), and polyfem::quadrature::Quadrature::weights.
Referenced by compute_weights().
|
private |
Definition at line 165 of file RBFWithQuadraticLagrange.cpp.
References centers_, is_volume(), and x.
Referenced by bases_values(), and compute_weights().
|
private |
Definition at line 458 of file RBFWithQuadraticLagrange.cpp.
References basis(), centers_, compute_constraints_matrix_2d(), compute_constraints_matrix_3d(), compute_kernels_matrix(), grad(), is_volume(), polyfem::logger(), polyfem::quadrature::Quadrature::points, val, polyfem::quadrature::Quadrature::weights, weights_, and x.
Referenced by RBFWithQuadraticLagrange().
void RBFWithQuadraticLagrange::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.
[in] | local_index | i-th RBF function to evaluate |
[in] | uv | #uv x dim matrix of coordinates to evaluate (in object domain) |
[out] | val | #uv x dim matrix of computed gradients |
Definition at line 104 of file RBFWithQuadraticLagrange.cpp.
References bases_grads(), centers_, and val.
Referenced by compute_weights().
|
inlineprivate |
Definition at line 68 of file RBFWithQuadraticLagrange.hpp.
References centers_.
Referenced by bases_grads(), compute_constraints_matrix_2d(), compute_constraints_matrix_2d_old(), compute_constraints_matrix_3d(), compute_kernels_matrix(), and compute_weights().
|
private |
Definition at line 87 of file RBFWithQuadraticLagrange.hpp.
Referenced by bases_grads(), compute_constraints_matrix_2d(), compute_constraints_matrix_2d_old(), compute_constraints_matrix_3d(), compute_kernels_matrix(), compute_weights(), grad(), and is_volume().
|
private |
Definition at line 90 of file RBFWithQuadraticLagrange.hpp.
Referenced by bases_grads(), bases_values(), and compute_weights().