PolyFEM
|
Form representing the contact potential and forces. More...
#include <ContactForm.hpp>
Public Member Functions | |
ContactForm (const ipc::CollisionMesh &collision_mesh, const double dhat, const double avg_mass, const bool use_convergent_formulation, const bool use_adaptive_barrier_stiffness, const bool is_time_dependent, const bool enable_shape_derivatives, const ipc::BroadPhaseMethod broad_phase_method, const double ccd_tolerance, const int ccd_max_iterations) | |
Construct a new Contact Form object. | |
virtual | ~ContactForm ()=default |
std::string | name () const override |
void | init (const Eigen::VectorXd &x) override |
Initialize the form. | |
virtual void | force_shape_derivative (const ipc::Collisions &collision_set, const Eigen::MatrixXd &solution, const Eigen::VectorXd &adjoint_sol, Eigen::VectorXd &term) |
void | update_quantities (const double t, const Eigen::VectorXd &x) override |
Update time-dependent fields. | |
double | max_step_size (const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override |
Determine the maximum step size allowable between the current and next solution. | |
void | line_search_begin (const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override |
Initialize variables used during the line search. | |
void | line_search_end () override |
Clear variables used during the line search. | |
void | solution_changed (const Eigen::VectorXd &new_x) override |
Update cached fields upon a change in the solution. | |
void | post_step (const polysolve::nonlinear::PostStepData &data) override |
Update fields after a step in the optimization. | |
bool | is_step_collision_free (const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override |
Checks if the step is collision free. | |
virtual void | update_barrier_stiffness (const Eigen::VectorXd &x, const Eigen::MatrixXd &grad_energy) |
Update the barrier stiffness based on the current elasticity energy. | |
Eigen::MatrixXd | compute_displaced_surface (const Eigen::VectorXd &x) const |
Compute the displaced positions of the surface nodes. | |
double | barrier_stiffness () const |
Get the current barrier stiffness. | |
void | set_barrier_stiffness (const double barrier_stiffness) |
Get the current barrier stiffness. | |
bool | use_adaptive_barrier_stiffness () const |
Get use_adaptive_barrier_stiffness. | |
bool | use_convergent_formulation () const |
Get use_convergent_formulation. | |
bool | enable_shape_derivatives () const |
double | weight () const override |
Get the form's multiplicative constant weight. | |
double | dhat () const |
const ipc::Collisions & | collision_set () const |
const ipc::BarrierPotential & | barrier_potential () const |
Public Member Functions inherited from polyfem::solver::Form | |
virtual | ~Form () |
virtual void | finish () |
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 bool | is_step_valid (const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const |
Determine if a step from solution x0 to solution x1 is allowed. | |
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. | |
void | set_weight (const double weight) |
Set the form's multiplicative constant weight. | |
void | set_output_dir (const std::string &output_dir) |
Public Attributes | |
bool | save_ccd_debug_meshes = false |
If true, output debug files. | |
Protected Member Functions | |
virtual double | value_unweighted (const Eigen::VectorXd &x) const override |
Compute the contact barrier potential value. | |
Eigen::VectorXd | value_per_element_unweighted (const Eigen::VectorXd &x) const override |
Compute the value of the form multiplied per element. | |
virtual void | first_derivative_unweighted (const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override |
Compute the first derivative of the value wrt x. | |
virtual void | second_derivative_unweighted (const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override |
Compute the second derivative of the value wrt x. | |
void | update_collision_set (const Eigen::MatrixXd &displaced_surface) |
Update the cached candidate set for the current solution. | |
Protected Member Functions inherited from polyfem::solver::Form | |
std::string | resolve_output_path (const std::string &path) const |
Protected Attributes | |
const ipc::CollisionMesh & | collision_mesh_ |
Collision mesh. | |
const double | dhat_ |
Barrier activation distance. | |
const double | dmin_ = 0 |
Minimum distance between elements. | |
const bool | use_adaptive_barrier_stiffness_ |
If true, use an adaptive barrier stiffness. | |
double | barrier_stiffness_ |
Barrier stiffness. | |
double | max_barrier_stiffness_ |
Maximum barrier stiffness to use when using adaptive barrier stiffness. | |
const double | avg_mass_ |
Average mass of the mesh (used for adaptive barrier stiffness) | |
const bool | is_time_dependent_ |
Is the simulation time dependent? | |
const bool | enable_shape_derivatives_ |
Enable shape derivatives computation. | |
const ipc::BroadPhaseMethod | broad_phase_method_ |
Broad phase method to use for distance and CCD evaluations. | |
const double | ccd_tolerance_ |
Continuous collision detection tolerance. | |
const int | ccd_max_iterations_ |
Continuous collision detection maximum iterations. | |
double | prev_distance_ |
Previous minimum distance between all elements. | |
bool | use_cached_candidates_ = false |
If true, use the cached candidate set for the current solution. | |
ipc::Collisions | collision_set_ |
Cached constraint set for the current solution. | |
ipc::Candidates | candidates_ |
Cached candidate set for the current solution. | |
const ipc::BarrierPotential | barrier_potential_ |
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_ |
Form representing the contact potential and forces.
Definition at line 35 of file ContactForm.hpp.
polyfem::solver::ContactForm::ContactForm | ( | const ipc::CollisionMesh & | collision_mesh, |
const double | dhat, | ||
const double | avg_mass, | ||
const bool | use_convergent_formulation, | ||
const bool | use_adaptive_barrier_stiffness, | ||
const bool | is_time_dependent, | ||
const bool | enable_shape_derivatives, | ||
const ipc::BroadPhaseMethod | broad_phase_method, | ||
const double | ccd_tolerance, | ||
const int | ccd_max_iterations | ||
) |
Construct a new Contact Form object.
collision_mesh | Reference to the collision mesh |
dhat | Barrier activation distance |
avg_mass | Average mass of the mesh |
use_adaptive_barrier_stiffness | If true, use an adaptive barrier stiffness |
is_time_dependent | Is the simulation time dependent? |
broad_phase_method | Broad phase method to use for distance and CCD evaluations |
ccd_tolerance | Continuous collision detection tolerance |
ccd_max_iterations | Continuous collision detection maximum iterations |
Definition at line 20 of file ContactForm.cpp.
References collision_set_, dhat_, enable_shape_derivatives(), prev_distance_, and use_convergent_formulation().
|
virtualdefault |
|
inline |
Definition at line 146 of file ContactForm.hpp.
References barrier_potential_.
Referenced by polyfem::solver::FrictionForm::force_shape_derivative(), and polyfem::solver::FrictionForm::update_lagging().
|
inline |
Get the current barrier stiffness.
Definition at line 129 of file ContactForm.hpp.
References barrier_stiffness_.
Referenced by polyfem::solver::FrictionForm::force_shape_derivative(), force_shape_derivative(), post_step(), set_barrier_stiffness(), update_barrier_stiffness(), and polyfem::solver::FrictionForm::update_lagging().
|
inline |
Definition at line 145 of file ContactForm.hpp.
References collision_set_.
Referenced by force_shape_derivative().
Eigen::MatrixXd polyfem::solver::ContactForm::compute_displaced_surface | ( | const Eigen::VectorXd & | x | ) | const |
Compute the displaced positions of the surface nodes.
Definition at line 69 of file ContactForm.cpp.
References collision_mesh_, polyfem::utils::unflatten(), and x.
Referenced by polyfem::solver::FrictionForm::compute_displaced_surface(), first_derivative_unweighted(), polyfem::solver::PeriodicContactForm::force_periodic_shape_derivative(), force_shape_derivative(), init(), is_step_collision_free(), line_search_begin(), max_step_size(), post_step(), second_derivative_unweighted(), solution_changed(), update_barrier_stiffness(), update_quantities(), value_per_element_unweighted(), and value_unweighted().
|
inline |
Definition at line 144 of file ContactForm.hpp.
References dhat_.
Referenced by polyfem::solver::FrictionForm::update_lagging().
|
inline |
Definition at line 137 of file ContactForm.hpp.
References enable_shape_derivatives_.
Referenced by ContactForm(), and polyfem::solver::FrictionForm::update_lagging().
|
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.
Reimplemented in polyfem::solver::PeriodicContactForm.
Definition at line 200 of file ContactForm.cpp.
References barrier_potential_, collision_mesh_, collision_set_, compute_displaced_surface(), and x.
Referenced by polyfem::solver::PeriodicContactForm::first_derivative_unweighted().
|
virtual |
Definition at line 54 of file ContactForm.cpp.
References barrier_potential_, barrier_stiffness(), collision_mesh_, collision_set(), and compute_displaced_surface().
|
overridevirtual |
Initialize the form.
x | Current solution |
Reimplemented from polyfem::solver::Form.
Reimplemented in polyfem::solver::PeriodicContactForm.
Definition at line 49 of file ContactForm.cpp.
References compute_displaced_surface(), update_collision_set(), and x.
Referenced by polyfem::solver::PeriodicContactForm::init().
|
overridevirtual |
Checks if the step is collision free.
Reimplemented from polyfem::solver::Form.
Reimplemented in polyfem::solver::PeriodicContactForm.
Definition at line 319 of file ContactForm.cpp.
References broad_phase_method_, candidates_, ccd_max_iterations_, ccd_tolerance_, collision_mesh_, compute_displaced_surface(), dmin_, polyfem::utils::is_valid(), and use_cached_candidates_.
Referenced by polyfem::solver::PeriodicContactForm::is_step_collision_free().
|
overridevirtual |
Initialize variables used during the line search.
x0 | Current solution |
x1 | Next solution |
Reimplemented from polyfem::solver::Form.
Reimplemented in polyfem::solver::PeriodicContactForm.
Definition at line 265 of file ContactForm.cpp.
References broad_phase_method_, candidates_, collision_mesh_, compute_displaced_surface(), dhat_, and use_cached_candidates_.
Referenced by polyfem::solver::PeriodicContactForm::line_search_begin().
|
overridevirtual |
Clear variables used during the line search.
Reimplemented from polyfem::solver::Form.
Definition at line 277 of file ContactForm.cpp.
References candidates_, and use_cached_candidates_.
|
overridevirtual |
Determine the maximum step size allowable between the current and next solution.
x0 | Current solution (step size = 0) |
x1 | Next solution (step size = 1) |
Reimplemented from polyfem::solver::Form.
Reimplemented in polyfem::solver::PeriodicContactForm.
Definition at line 218 of file ContactForm.cpp.
References broad_phase_method_, candidates_, ccd_max_iterations_, ccd_tolerance_, collision_mesh_, compute_displaced_surface(), dmin_, polyfem::F, polyfem::log_and_throw_error(), polyfem::logger(), polyfem::solver::Form::resolve_output_path(), save_ccd_debug_meshes, use_cached_candidates_, polyfem::mesh::V0, and polyfem::mesh::V1.
Referenced by polyfem::solver::PeriodicContactForm::max_step_size().
|
inlineoverridevirtual |
Implements polyfem::solver::Form.
Definition at line 59 of file ContactForm.hpp.
|
overridevirtual |
Update fields after a step in the optimization.
iter_num | Optimization iteration number |
x | Current solution |
Reimplemented from polyfem::solver::Form.
Reimplemented in polyfem::solver::PeriodicContactForm.
Definition at line 283 of file ContactForm.cpp.
References barrier_stiffness(), barrier_stiffness_, collision_mesh_, collision_set_, compute_displaced_surface(), is_time_dependent_, polyfem::logger(), max_barrier_stiffness_, prev_distance_, and use_adaptive_barrier_stiffness_.
Referenced by polyfem::solver::PeriodicContactForm::post_step().
|
overrideprotectedvirtual |
Compute the second derivative of the value wrt x.
x | Current solution |
hessian | Output Hessian of the value wrt x |
Implements polyfem::solver::Form.
Reimplemented in polyfem::solver::PeriodicContactForm.
Definition at line 206 of file ContactForm.cpp.
References barrier_potential_, collision_mesh_, collision_set_, compute_displaced_surface(), POLYFEM_SCOPED_TIMER, polyfem::solver::Form::project_to_psd_, and x.
Referenced by polyfem::solver::PeriodicContactForm::force_periodic_shape_derivative(), and polyfem::solver::PeriodicContactForm::second_derivative_unweighted().
|
inline |
Get the current barrier stiffness.
Definition at line 131 of file ContactForm.hpp.
References barrier_stiffness(), and barrier_stiffness_.
|
overridevirtual |
Update cached fields upon a change in the solution.
new_x | New solution |
Reimplemented from polyfem::solver::Form.
Reimplemented in polyfem::solver::PeriodicContactForm.
Definition at line 213 of file ContactForm.cpp.
References compute_displaced_surface(), and update_collision_set().
Referenced by polyfem::solver::PeriodicContactForm::solution_changed().
|
virtual |
Update the barrier stiffness based on the current elasticity energy.
x | Current solution |
Reimplemented in polyfem::solver::PeriodicContactForm.
Definition at line 74 of file ContactForm.cpp.
References avg_mass_, barrier_potential_, barrier_stiffness(), barrier_stiffness_, broad_phase_method_, collision_mesh_, collision_set_, compute_displaced_surface(), dhat_, dmin_, polyfem::logger(), max_barrier_stiffness_, update_collision_set(), use_adaptive_barrier_stiffness(), use_convergent_formulation(), polyfem::solver::Form::weight_, and x.
Referenced by polyfem::solver::PeriodicContactForm::update_barrier_stiffness().
|
protected |
Update the cached candidate set for the current solution.
displaced_surface | Vertex positions displaced by the current solution |
Definition at line 129 of file ContactForm.cpp.
References broad_phase_method_, candidates_, collision_mesh_, collision_set_, dhat_, dmin_, and use_cached_candidates_.
Referenced by init(), solution_changed(), update_barrier_stiffness(), and update_quantities().
|
overridevirtual |
Update time-dependent fields.
t | Current time |
x | Current solution at time t |
Reimplemented from polyfem::solver::Form.
Reimplemented in polyfem::solver::PeriodicContactForm.
Definition at line 64 of file ContactForm.cpp.
References compute_displaced_surface(), update_collision_set(), and x.
Referenced by polyfem::solver::PeriodicContactForm::update_quantities().
|
inline |
Get use_adaptive_barrier_stiffness.
Definition at line 133 of file ContactForm.hpp.
References use_adaptive_barrier_stiffness_.
Referenced by update_barrier_stiffness().
|
inline |
Get use_convergent_formulation.
Definition at line 135 of file ContactForm.hpp.
References collision_set_.
Referenced by ContactForm(), update_barrier_stiffness(), and polyfem::solver::FrictionForm::update_lagging().
|
overrideprotectedvirtual |
Compute the value of the form multiplied per element.
x | Current solution |
Reimplemented from polyfem::solver::Form.
Definition at line 150 of file ContactForm.cpp.
References collision_mesh_, collision_set_, compute_displaced_surface(), polyfem::F, polyfem::utils::maybe_parallel_for(), V, and x.
|
overrideprotectedvirtual |
Compute the contact barrier potential value.
x | Current solution |
Implements polyfem::solver::Form.
Reimplemented in polyfem::solver::PeriodicContactForm.
Definition at line 145 of file ContactForm.cpp.
References barrier_potential_, collision_mesh_, collision_set_, compute_displaced_surface(), and x.
Referenced by polyfem::solver::PeriodicContactForm::value_unweighted().
|
inlineoverridevirtual |
Get the form's multiplicative constant weight.
Reimplemented from polyfem::solver::Form.
Definition at line 139 of file ContactForm.hpp.
References barrier_stiffness_, and polyfem::solver::Form::weight_.
Referenced by polyfem::solver::PeriodicContactForm::force_periodic_shape_derivative().
|
protected |
Average mass of the mesh (used for adaptive barrier stiffness)
Definition at line 170 of file ContactForm.hpp.
Referenced by update_barrier_stiffness().
|
protected |
Definition at line 195 of file ContactForm.hpp.
Referenced by barrier_potential(), first_derivative_unweighted(), polyfem::solver::PeriodicContactForm::force_periodic_shape_derivative(), force_shape_derivative(), second_derivative_unweighted(), update_barrier_stiffness(), and value_unweighted().
|
protected |
Barrier stiffness.
Definition at line 165 of file ContactForm.hpp.
Referenced by barrier_stiffness(), post_step(), set_barrier_stiffness(), update_barrier_stiffness(), and weight().
|
protected |
Broad phase method to use for distance and CCD evaluations.
Definition at line 179 of file ContactForm.hpp.
Referenced by is_step_collision_free(), line_search_begin(), max_step_size(), update_barrier_stiffness(), and update_collision_set().
|
protected |
Cached candidate set for the current solution.
Definition at line 193 of file ContactForm.hpp.
Referenced by is_step_collision_free(), line_search_begin(), line_search_end(), max_step_size(), and update_collision_set().
|
protected |
Continuous collision detection maximum iterations.
Definition at line 183 of file ContactForm.hpp.
Referenced by is_step_collision_free(), and max_step_size().
|
protected |
Continuous collision detection tolerance.
Definition at line 181 of file ContactForm.hpp.
Referenced by is_step_collision_free(), and max_step_size().
|
protected |
Collision mesh.
Definition at line 154 of file ContactForm.hpp.
Referenced by compute_displaced_surface(), first_derivative_unweighted(), polyfem::solver::PeriodicContactForm::force_periodic_shape_derivative(), force_shape_derivative(), is_step_collision_free(), line_search_begin(), max_step_size(), polyfem::solver::PeriodicContactForm::PeriodicContactForm(), post_step(), second_derivative_unweighted(), polyfem::solver::PeriodicContactForm::single_to_tiled(), polyfem::solver::PeriodicContactForm::tiled_to_single_grad(), update_barrier_stiffness(), update_collision_set(), polyfem::solver::PeriodicContactForm::update_projection(), value_per_element_unweighted(), and value_unweighted().
|
protected |
Cached constraint set for the current solution.
Definition at line 191 of file ContactForm.hpp.
Referenced by collision_set(), ContactForm(), first_derivative_unweighted(), post_step(), second_derivative_unweighted(), update_barrier_stiffness(), update_collision_set(), use_convergent_formulation(), value_per_element_unweighted(), and value_unweighted().
|
protected |
Barrier activation distance.
Definition at line 157 of file ContactForm.hpp.
Referenced by ContactForm(), dhat(), line_search_begin(), update_barrier_stiffness(), and update_collision_set().
|
protected |
Minimum distance between elements.
Definition at line 160 of file ContactForm.hpp.
Referenced by is_step_collision_free(), max_step_size(), update_barrier_stiffness(), and update_collision_set().
|
protected |
Enable shape derivatives computation.
Definition at line 176 of file ContactForm.hpp.
Referenced by enable_shape_derivatives().
|
protected |
Is the simulation time dependent?
Definition at line 173 of file ContactForm.hpp.
Referenced by post_step().
|
protected |
Maximum barrier stiffness to use when using adaptive barrier stiffness.
Definition at line 167 of file ContactForm.hpp.
Referenced by post_step(), and update_barrier_stiffness().
|
protected |
Previous minimum distance between all elements.
Definition at line 186 of file ContactForm.hpp.
Referenced by ContactForm(), and post_step().
bool polyfem::solver::ContactForm::save_ccd_debug_meshes = false |
If true, output debug files.
Definition at line 142 of file ContactForm.hpp.
Referenced by max_step_size().
|
protected |
If true, use an adaptive barrier stiffness.
Definition at line 163 of file ContactForm.hpp.
Referenced by post_step(), and use_adaptive_barrier_stiffness().
|
protected |
If true, use the cached candidate set for the current solution.
Definition at line 189 of file ContactForm.hpp.
Referenced by is_step_collision_free(), line_search_begin(), line_search_end(), max_step_size(), and update_collision_set().