PolyFEM
Loading...
Searching...
No Matches
polyfem::basis::RBFWithLinear Class Reference

#include <RBFWithLinear.hpp>

Public Member Functions

 RBFWithLinear (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.
 
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 (const int num_bases, const quadrature::Quadrature &quadr, const Eigen::MatrixXd &local_basis_integral, Eigen::MatrixXd &L, Eigen::MatrixXd &t) const
 
void compute_weights (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_
 

Detailed Description

Definition at line 10 of file RBFWithLinear.hpp.

Constructor & Destructor Documentation

◆ RBFWithLinear()

RBFWithLinear::RBFWithLinear ( 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.

Parameters
[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 containing the expected value of the integral of each basis over the polytope
[in]quadrQuadrature 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_constraintsImpose integral constraints to guarantee linear reproduction for the Poisson equation

Definition at line 62 of file RBFWithLinear.cpp.

References compute_weights().

Here is the call graph for this function:

Member Function Documentation

◆ bases_grads()

void RBFWithLinear::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.

Parameters
[in]axisThe 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 111 of file RBFWithLinear.cpp.

References centers_, is_volume(), val, weights_, and x.

Referenced by grad().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ bases_values()

void RBFWithLinear::bases_values ( const Eigen::MatrixXd &  samples,
Eigen::MatrixXd &  val 
) const

Batch evaluates the RBF + polynomials on a set of sample points.

Parameters
[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 99 of file RBFWithLinear.cpp.

References compute_kernels_matrix(), val, and weights_.

Referenced by basis().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ basis()

void RBFWithLinear::basis ( const int  local_index,
const Eigen::MatrixXd &  uv,
Eigen::MatrixXd &  val 
) const

Evaluates one RBF function over a list of coordinates.

Parameters
[in]local_indexi-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 76 of file RBFWithLinear.cpp.

References bases_values(), and val.

Here is the call graph for this function:

◆ compute_constraints_matrix()

void RBFWithLinear::compute_constraints_matrix ( const int  num_bases,
const quadrature::Quadrature quadr,
const Eigen::MatrixXd &  local_basis_integral,
Eigen::MatrixXd &  L,
Eigen::MatrixXd &  t 
) const
private

Definition at line 151 of file RBFWithLinear.cpp.

References centers_, is_volume(), polyfem::quadrature::Quadrature::points, polyfem::quadrature::Quadrature::weights, and x.

Referenced by compute_weights().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_kernels_matrix()

void RBFWithLinear::compute_kernels_matrix ( const Eigen::MatrixXd &  samples,
Eigen::MatrixXd &  A 
) const
private

Definition at line 134 of file RBFWithLinear.cpp.

References centers_, is_volume(), and x.

Referenced by bases_values(), and compute_weights().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_weights()

void RBFWithLinear::compute_weights ( const Eigen::MatrixXd &  collocation_points,
const Eigen::MatrixXd &  local_basis_integral,
const quadrature::Quadrature quadr,
Eigen::MatrixXd &  rhs,
bool  with_constraints 
)
private

Definition at line 191 of file RBFWithLinear.cpp.

References centers_, compute_constraints_matrix(), compute_kernels_matrix(), polyfem::logger(), polyfem::quadrature::Quadrature::weights, and weights_.

Referenced by RBFWithLinear().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ grad()

void RBFWithLinear::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.

Parameters
[in]local_indexi-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 85 of file RBFWithLinear.cpp.

References bases_grads(), centers_, and val.

Here is the call graph for this function:

◆ is_volume()

bool polyfem::basis::RBFWithLinear::is_volume ( ) const
inlineprivate

Definition at line 63 of file RBFWithLinear.hpp.

References centers_.

Referenced by bases_grads(), compute_constraints_matrix(), and compute_kernels_matrix().

Here is the caller graph for this function:

Member Data Documentation

◆ centers_

Eigen::MatrixXd polyfem::basis::RBFWithLinear::centers_
private

◆ weights_

Eigen::MatrixXd polyfem::basis::RBFWithLinear::weights_
private

Definition at line 82 of file RBFWithLinear.hpp.

Referenced by bases_grads(), bases_values(), and compute_weights().


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