|
| | BiharmonicEnergy3D (const MatrixXd &pts, const double &M, const VectorXd &L_w, const double weight=1) |
| | The biharmonic energy of a vertex position on a polyline.
|
| |
|
TVector | initial_position () const |
| |
|
double | value (const TVector &x) override |
| |
|
void | gradient (const TVector &x, TVector &gradv) override |
| |
|
void | hessian (const TVector &x, THessian &hessian) override |
| |
|
void | hessian (const TVector &x, MatrixXd &hessian) override |
| |
|
void | solution_changed (const TVector &new_x) override |
| |
|
|
static void | global_mass_and_stiffness (const MatrixXd &pts, const MatrixXi &tris, Eigen::SparseMatrix< double > &M, Eigen::SparseMatrix< double > &L_w) |
| | Construct mass and stiffness for all vertices using IGL.
|
| |
| static void | extract_local_mass_and_stiffness (const size_t vid, const Eigen::SparseMatrix< double > &M_glob, const Eigen::SparseMatrix< double > &L_w_glob, double &M_loc, VectorXd &L_w_loc) |
| | Extract the mass and stiffness for local optimization from the global mass and stiffness matrices.
|
| |
| static void | adjacency_from_stiffness (const size_t vid, const Eigen::SparseMatrix< double > &L_w_glob, std::vector< size_t > &adj) |
| | Get the adjacency for vid from the global stiffness matrix.
|
| |
|
static void | uniform_mass_and_stiffness (const MatrixXd &pts, double &M, VectorXd &L_w) |
| |
◆ BiharmonicEnergy3D()
| wmtk::optimization::BiharmonicEnergy3D::BiharmonicEnergy3D |
( |
const MatrixXd & |
pts, |
|
|
const double & |
M, |
|
|
const VectorXd & |
L_w, |
|
|
const double |
weight = 1 |
|
) |
| |
The biharmonic energy of a vertex position on a polyline.
The energy is defined as: \int_M \Vert \Delta_M p \Vert^2 In descrete form: \sum_{i=1}^n A_i \Vert L p_i \Vert^2
- Parameters
-
| pts | The optimized position (pts.row(0)) and its neighbors. |
◆ adjacency_from_stiffness()
| void wmtk::optimization::BiharmonicEnergy3D::adjacency_from_stiffness |
( |
const size_t |
vid, |
|
|
const Eigen::SparseMatrix< double > & |
L_w_glob, |
|
|
std::vector< size_t > & |
adj |
|
) |
| |
|
static |
Get the adjacency for vid from the global stiffness matrix.
The adjacency does not contain vid itself.
◆ extract_local_mass_and_stiffness()
| void wmtk::optimization::BiharmonicEnergy3D::extract_local_mass_and_stiffness |
( |
const size_t |
vid, |
|
|
const Eigen::SparseMatrix< double > & |
M_glob, |
|
|
const Eigen::SparseMatrix< double > & |
L_w_glob, |
|
|
double & |
M_loc, |
|
|
VectorXd & |
L_w_loc |
|
) |
| |
|
static |
Extract the mass and stiffness for local optimization from the global mass and stiffness matrices.
The local stiffness order is changed to fit to the convention, i.e., L_w_loc[0] = L_w_glob(vid,vid). All other entries are sorted according to their ID.
- Parameters
-
| vid | The vertex ID for the local optimization. |
| M_glob | The global mass matrix. |
| L_w_glob | The global stiffness matrix. |
| M_loc | The local mass. |
| L_w_loc | The local stiffness without 0 entries and the vid entry moved to the front. |
The stiffness matrix is symmetric. It doesn't matter if we extract the row or column for vid.
◆ m_L_w_row
| VectorXd wmtk::optimization::BiharmonicEnergy3D::m_L_w_row |
|
private |
The first row of stiffness matrix L_w which belongs to the optimized vertex. By convention, this needs to be sorted such that the first entry belongs to the optimized vertex.
The full stiffness matrix L_w has shape NxN but the rows of fixed vertices are empty. Therefore, we only need to store the one non-empty row of the optimized vertex. By convention of m_pts, this must be the first one, i.e., L_w.row(0).
◆ m_LTML_row
| VectorXd wmtk::optimization::BiharmonicEnergy3D::m_LTML_row |
|
private |
We only need the first row of LTML, because all other vertices are fixed.
L = M_inv * L_w LTML = L^T * M * L = L_w^T * M_inv * L_w We only want to optimize one vertex, so we only need one row of LTML: LTML.row(0) = M_inv * L_w(0,0) * L_w.row(0) <– this is what is stored in m_LTML_row
◆ m_pts
| MatrixXd wmtk::optimization::BiharmonicEnergy3D::m_pts |
|
private |
The optimized point and its neighbors. For the optimization, we use the convention that m_pts.row(0) is the optimized vertex. The neighbors are fixed.
The documentation for this class was generated from the following files:
- /home/runner/work/wildmeshing-toolkit/wildmeshing-toolkit/src/wmtk/optimization/DirichletEnergy.hpp
- /home/runner/work/wildmeshing-toolkit/wildmeshing-toolkit/src/wmtk/optimization/DirichletEnergy.cpp