PolyFEM
|
#include <InversionBarrierForm.hpp>
Public Member Functions | |
InversionBarrierForm (const Eigen::MatrixXd &rest_positions, const Eigen::MatrixXi &elements, const int dim, const double vhat) | |
std::string | name () const override |
bool | is_step_valid (const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override |
Determine if a step from solution x0 to solution x1 is allowed. | |
Public Member Functions inherited from polyfem::solver::Form | |
virtual | ~Form () |
virtual void | init (const Eigen::VectorXd &x) |
Initialize the form. | |
virtual double | value (const Eigen::VectorXd &x) const |
Compute the value of the form multiplied with the weigth. | |
Eigen::VectorXd | value_per_element (const Eigen::VectorXd &x) const |
Compute the value of the form multiplied with the weigth. | |
virtual void | first_derivative (const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const |
Compute the first derivative of the value wrt x multiplied with the weigth. | |
void | second_derivative (const Eigen::VectorXd &x, StiffnessMatrix &hessian) const |
Compute the second derivative of the value wrt x multiplied with the weigth. | |
virtual double | max_step_size (const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const |
Determine the maximum step size allowable between the current and next solution. | |
virtual void | line_search_begin (const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) |
Initialize variables used during the line search. | |
virtual void | line_search_end () |
Clear variables used during the line search. | |
virtual void | post_step (const polysolve::nonlinear::PostStepData &data) |
Update fields after a step in the optimization. | |
virtual void | solution_changed (const Eigen::VectorXd &new_x) |
Update cached fields upon a change in the solution. | |
virtual void | update_quantities (const double t, const Eigen::VectorXd &x) |
Update time-dependent fields. | |
virtual void | init_lagging (const Eigen::VectorXd &x) |
Initialize lagged fields TODO: more than one step. | |
virtual void | update_lagging (const Eigen::VectorXd &x, const int iter_num) |
Update lagged fields. | |
virtual int | max_lagging_iterations () const |
Get the maximum number of lagging iteration allowable. | |
virtual bool | uses_lagging () const |
Does this form require lagging? | |
void | set_project_to_psd (bool val) |
Set project to psd. | |
bool | is_project_to_psd () const |
Get if the form's second derivative is projected to psd. | |
void | enable () |
Enable the form. | |
void | disable () |
Disable the form. | |
void | set_enabled (const bool enabled) |
Set if the form is enabled. | |
bool | enabled () const |
Determine if the form is enabled. | |
virtual double | weight () const |
Get the form's multiplicative constant weight. | |
void | set_weight (const double weight) |
Set the form's multiplicative constant weight. | |
virtual bool | is_step_collision_free (const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const |
Checks if the step is collision free. | |
void | set_output_dir (const std::string &output_dir) |
Protected Member Functions | |
double | value_unweighted (const Eigen::VectorXd &x) const override |
Compute the value of the form. | |
void | first_derivative_unweighted (const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override |
Compute the first derivative of the value wrt x. | |
void | second_derivative_unweighted (const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override |
Compute the second derivative of the value wrt x. | |
Protected Member Functions inherited from polyfem::solver::Form | |
std::string | resolve_output_path (const std::string &path) const |
virtual Eigen::VectorXd | value_per_element_unweighted (const Eigen::VectorXd &x) const |
Compute the value of the form multiplied per element. | |
Static Protected Member Functions | |
static double | element_volume (const Eigen::MatrixXd &element_vertices) |
static Eigen::VectorXd | element_volume_gradient (const Eigen::MatrixXd &element_vertices) |
static Eigen::MatrixXd | element_volume_hessian (const Eigen::MatrixXd &element_vertices) |
Private Attributes | |
Eigen::MatrixXd | rest_positions_ |
Eigen::MatrixXi | elements_ |
int | dim_ |
double | vhat_ |
Additional Inherited Members | |
Protected Attributes inherited from polyfem::solver::Form | |
bool | project_to_psd_ = false |
If true, the form's second derivative is projected to be positive semidefinite. | |
double | weight_ = 1 |
weight of the form (e.g., AL penalty weight or Δt²) | |
bool | enabled_ = true |
If true, the form is enabled. | |
std::string | output_dir_ |
Definition at line 10 of file InversionBarrierForm.hpp.
polyfem::solver::InversionBarrierForm::InversionBarrierForm | ( | const Eigen::MatrixXd & | rest_positions, |
const Eigen::MatrixXi & | elements, | ||
const int | dim, | ||
const double | vhat | ||
) |
Definition at line 14 of file InversionBarrierForm.cpp.
|
staticprotected |
Definition at line 114 of file InversionBarrierForm.cpp.
References polyfem::utils::tetrahedron_volume(), and polyfem::utils::triangle_area().
Referenced by is_step_valid().
|
staticprotected |
Definition at line 123 of file InversionBarrierForm.cpp.
References polyfem::utils::tetrahedron_volume_gradient(), and polyfem::utils::triangle_area_2D_gradient().
|
staticprotected |
Definition at line 148 of file InversionBarrierForm.cpp.
References polyfem::utils::tetrahedron_volume_hessian(), and polyfem::utils::triangle_area_2D_hessian().
|
overrideprotectedvirtual |
Compute the first derivative of the value wrt x.
[in] | x | Current solution |
[out] | gradv | Output gradient of the value wrt x |
Implements polyfem::solver::Form.
Definition at line 44 of file InversionBarrierForm.cpp.
References dim_, elements_, polyfem::utils::maybe_parallel_for(), rest_positions_, polyfem::utils::unflatten(), V, vhat_, and x.
|
overridevirtual |
Determine if a step from solution x0 to solution x1 is allowed.
x0 | Current solution |
x1 | Proposed next solution |
Reimplemented from polyfem::solver::Form.
Definition at line 173 of file InversionBarrierForm.cpp.
References dim_, element_volume(), elements_, rest_positions_, polyfem::utils::unflatten(), and V.
|
inlineoverridevirtual |
Implements polyfem::solver::Form.
Definition at line 16 of file InversionBarrierForm.hpp.
|
overrideprotectedvirtual |
Compute the second derivative of the value wrt x.
[in] | x | Current solution |
[out] | hessian | Output Hessian of the value wrt x |
Implements polyfem::solver::Form.
Definition at line 72 of file InversionBarrierForm.cpp.
References polyfem::utils::create_thread_storage(), dim_, elements_, polyfem::utils::maybe_parallel_for(), rest_positions_, polyfem::utils::unflatten(), V, vhat_, and x.
|
overrideprotectedvirtual |
Compute the value of the form.
x | Current solution |
Implements polyfem::solver::Form.
Definition at line 20 of file InversionBarrierForm.cpp.
References dim_, elements_, polyfem::utils::maybe_parallel_for(), rest_positions_, polyfem::utils::unflatten(), V, vhat_, and x.
|
private |
Definition at line 49 of file InversionBarrierForm.hpp.
Referenced by first_derivative_unweighted(), is_step_valid(), second_derivative_unweighted(), and value_unweighted().
|
private |
Definition at line 48 of file InversionBarrierForm.hpp.
Referenced by first_derivative_unweighted(), is_step_valid(), second_derivative_unweighted(), and value_unweighted().
|
private |
Definition at line 47 of file InversionBarrierForm.hpp.
Referenced by first_derivative_unweighted(), is_step_valid(), second_derivative_unweighted(), and value_unweighted().
|
private |
Definition at line 50 of file InversionBarrierForm.hpp.
Referenced by first_derivative_unweighted(), second_derivative_unweighted(), and value_unweighted().