PolyFEM
|
Form representing the contact potential and forces. More...
#include <NormalAdhesionForm.hpp>
Public Member Functions | |
NormalAdhesionForm (const ipc::CollisionMesh &collision_mesh, const double dhat_p, const double dhat_a, const double Y, 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 NormalAdhesion Form object. | |
virtual | ~NormalAdhesionForm ()=default |
std::string | name () const override |
void | init (const Eigen::VectorXd &x) override |
Initialize the form. | |
virtual void | force_shape_derivative (const ipc::NormalCollisions &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. | |
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. | |
Eigen::MatrixXd | compute_displaced_surface (const Eigen::VectorXd &x) const |
Compute the displaced positions of the surface nodes. | |
bool | enable_shape_derivatives () const |
double | dhat_a () const |
double | dhat_p () const |
double | Y () const |
const ipc::NormalCollisions & | collision_set () const |
const ipc::NormalAdhesionPotential & | normal_adhesion_potential () const |
![]() | |
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 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 | 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) |
virtual void | set_scale (const double scale) |
sets the scale for the form | |
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 normal adhesion 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. | |
![]() | |
std::string | resolve_output_path (const std::string &path) const |
Protected Attributes | |
const ipc::CollisionMesh & | collision_mesh_ |
Collision mesh. | |
const double | dhat_p_ |
Maximum adhesion strength distance. | |
const double | dhat_a_ |
Adhesion activation distance. | |
const double | Y_ |
Adhesion strength. | |
const double | dmin_ = 0 |
Minimum distance between elements. | |
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 ipc::TightInclusionCCD | tight_inclusion_ccd_ |
Continuous collision detection specification object. | |
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::NormalCollisions | collision_set_ |
Cached constraint set for the current solution. | |
ipc::Candidates | candidates_ |
Cached candidate set for the current solution. | |
const ipc::NormalAdhesionPotential | normal_adhesion_potential_ |
![]() | |
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 18 of file NormalAdhesionForm.hpp.
polyfem::solver::NormalAdhesionForm::NormalAdhesionForm | ( | const ipc::CollisionMesh & | collision_mesh, |
const double | dhat_p, | ||
const double | dhat_a, | ||
const double | Y, | ||
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 NormalAdhesion Form object.
collision_mesh | Reference to the collision mesh |
dhat_p | Distance of largest adhesion force |
dhat_a | Adhesion activation distance |
Y | Adhesion strength |
is_time_dependent | Is the simulation time dependent? |
enable_shape_derivatives | Enable shape derivatives |
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 NormalAdhesionForm.cpp.
References collision_set_, dhat_a(), dhat_p(), enable_shape_derivatives(), and prev_distance_.
|
virtualdefault |
|
inline |
Definition at line 105 of file NormalAdhesionForm.hpp.
References collision_set_.
Referenced by force_shape_derivative().
Eigen::MatrixXd polyfem::solver::NormalAdhesionForm::compute_displaced_surface | ( | const Eigen::VectorXd & | x | ) | const |
Compute the displaced positions of the surface nodes.
Definition at line 67 of file NormalAdhesionForm.cpp.
References collision_mesh_, polyfem::utils::unflatten(), and x.
Referenced by polyfem::solver::TangentialAdhesionForm::compute_displaced_surface(), first_derivative_unweighted(), force_shape_derivative(), init(), line_search_begin(), post_step(), second_derivative_unweighted(), solution_changed(), update_quantities(), value_per_element_unweighted(), and value_unweighted().
|
inline |
Definition at line 102 of file NormalAdhesionForm.hpp.
References dhat_a_.
Referenced by NormalAdhesionForm(), and polyfem::solver::TangentialAdhesionForm::update_lagging().
|
inline |
Definition at line 103 of file NormalAdhesionForm.hpp.
References dhat_p_.
Referenced by NormalAdhesionForm().
|
inline |
Definition at line 97 of file NormalAdhesionForm.hpp.
References enable_shape_derivatives_.
Referenced by NormalAdhesionForm().
|
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 143 of file NormalAdhesionForm.cpp.
References collision_mesh_, collision_set_, compute_displaced_surface(), normal_adhesion_potential_, and x.
|
virtual |
Definition at line 52 of file NormalAdhesionForm.cpp.
References collision_mesh_, collision_set(), compute_displaced_surface(), and normal_adhesion_potential_.
|
overridevirtual |
Initialize the form.
x | Current solution |
Reimplemented from polyfem::solver::Form.
Definition at line 47 of file NormalAdhesionForm.cpp.
References compute_displaced_surface(), update_collision_set(), and x.
|
overridevirtual |
Initialize variables used during the line search.
x0 | Current solution |
x1 | Next solution |
Reimplemented from polyfem::solver::Form.
Definition at line 170 of file NormalAdhesionForm.cpp.
References broad_phase_method_, candidates_, collision_mesh_, compute_displaced_surface(), dhat_a_, and use_cached_candidates_.
|
overridevirtual |
Clear variables used during the line search.
Reimplemented from polyfem::solver::Form.
Definition at line 182 of file NormalAdhesionForm.cpp.
References candidates_, and use_cached_candidates_.
|
inlineoverridevirtual |
Implements polyfem::solver::Form.
Definition at line 42 of file NormalAdhesionForm.hpp.
|
inline |
Definition at line 106 of file NormalAdhesionForm.hpp.
References normal_adhesion_potential_.
Referenced by polyfem::solver::TangentialAdhesionForm::force_shape_derivative(), and polyfem::solver::TangentialAdhesionForm::update_lagging().
|
overridevirtual |
Update fields after a step in the optimization.
iter_num | Optimization iteration number |
x | Current solution |
Reimplemented from polyfem::solver::Form.
Definition at line 188 of file NormalAdhesionForm.cpp.
References collision_mesh_, collision_set_, compute_displaced_surface(), and prev_distance_.
|
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.
Definition at line 149 of file NormalAdhesionForm.cpp.
References collision_mesh_, collision_set_, compute_displaced_surface(), normal_adhesion_potential_, POLYFEM_SCOPED_TIMER, polyfem::solver::Form::project_to_psd_, and x.
|
overridevirtual |
Update cached fields upon a change in the solution.
new_x | New solution |
Reimplemented from polyfem::solver::Form.
Definition at line 165 of file NormalAdhesionForm.cpp.
References compute_displaced_surface(), and update_collision_set().
|
protected |
Update the cached candidate set for the current solution.
displaced_surface | Vertex positions displaced by the current solution |
Definition at line 72 of file NormalAdhesionForm.cpp.
References broad_phase_method_, candidates_, collision_mesh_, collision_set_, dhat_a_, dmin_, and use_cached_candidates_.
Referenced by init(), solution_changed(), and update_quantities().
|
overridevirtual |
Update time-dependent fields.
t | Current time |
x | Current solution at time t |
Reimplemented from polyfem::solver::Form.
Definition at line 62 of file NormalAdhesionForm.cpp.
References compute_displaced_surface(), update_collision_set(), and x.
|
overrideprotectedvirtual |
Compute the value of the form multiplied per element.
x | Current solution |
Reimplemented from polyfem::solver::Form.
Definition at line 93 of file NormalAdhesionForm.cpp.
References collision_mesh_, collision_set_, compute_displaced_surface(), polyfem::F, polyfem::utils::maybe_parallel_for(), V, and x.
|
overrideprotectedvirtual |
Compute the normal adhesion potential value.
x | Current solution |
Implements polyfem::solver::Form.
Definition at line 88 of file NormalAdhesionForm.cpp.
References collision_mesh_, collision_set_, compute_displaced_surface(), normal_adhesion_potential_, and x.
|
inline |
Definition at line 104 of file NormalAdhesionForm.hpp.
References Y_.
|
protected |
Broad phase method to use for distance and CCD evaluations.
Definition at line 135 of file NormalAdhesionForm.hpp.
Referenced by line_search_begin(), and update_collision_set().
|
protected |
Cached candidate set for the current solution.
Definition at line 147 of file NormalAdhesionForm.hpp.
Referenced by line_search_begin(), line_search_end(), and update_collision_set().
|
protected |
Collision mesh.
Definition at line 114 of file NormalAdhesionForm.hpp.
Referenced by compute_displaced_surface(), first_derivative_unweighted(), force_shape_derivative(), line_search_begin(), post_step(), second_derivative_unweighted(), update_collision_set(), value_per_element_unweighted(), and value_unweighted().
|
protected |
Cached constraint set for the current solution.
Definition at line 145 of file NormalAdhesionForm.hpp.
Referenced by collision_set(), first_derivative_unweighted(), NormalAdhesionForm(), post_step(), second_derivative_unweighted(), update_collision_set(), value_per_element_unweighted(), and value_unweighted().
|
protected |
Adhesion activation distance.
Definition at line 120 of file NormalAdhesionForm.hpp.
Referenced by dhat_a(), line_search_begin(), and update_collision_set().
|
protected |
Maximum adhesion strength distance.
Definition at line 117 of file NormalAdhesionForm.hpp.
Referenced by dhat_p().
|
protected |
Minimum distance between elements.
Definition at line 126 of file NormalAdhesionForm.hpp.
Referenced by update_collision_set().
|
protected |
Enable shape derivatives computation.
Definition at line 132 of file NormalAdhesionForm.hpp.
Referenced by enable_shape_derivatives().
|
protected |
Is the simulation time dependent?
Definition at line 129 of file NormalAdhesionForm.hpp.
|
protected |
Definition at line 149 of file NormalAdhesionForm.hpp.
Referenced by first_derivative_unweighted(), force_shape_derivative(), normal_adhesion_potential(), second_derivative_unweighted(), and value_unweighted().
|
protected |
Previous minimum distance between all elements.
Definition at line 140 of file NormalAdhesionForm.hpp.
Referenced by NormalAdhesionForm(), and post_step().
bool polyfem::solver::NormalAdhesionForm::save_ccd_debug_meshes = false |
If true, output debug files.
Definition at line 100 of file NormalAdhesionForm.hpp.
|
protected |
Continuous collision detection specification object.
Definition at line 137 of file NormalAdhesionForm.hpp.
|
protected |
If true, use the cached candidate set for the current solution.
Definition at line 143 of file NormalAdhesionForm.hpp.
Referenced by line_search_begin(), line_search_end(), and update_collision_set().
|
protected |