PolyFEM
|
Form of the augmented lagrangian for bc constraints. More...
#include <BCLagrangianForm.hpp>
Public Member Functions | |
BCLagrangianForm (const int ndof, const std::vector< int > &boundary_nodes, const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< mesh::LocalBoundary > &local_neumann_boundary, const int n_boundary_samples, const StiffnessMatrix &mass, const assembler::RhsAssembler &rhs_assembler, const size_t obstacle_ndof, const bool is_time_dependent, const double t) | |
Construct a new BCLagrangianForm object with a time dependent Dirichlet boundary. | |
std::string | name () const override |
BCLagrangianForm (const int ndof, const std::vector< int > &boundary_nodes, const StiffnessMatrix &mass, const size_t obstacle_ndof, const Eigen::MatrixXd &target_x) | |
Construct a new BCLagrangianForm object with a fixed Dirichlet boundary. | |
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. | |
void | update_quantities (const double t, const Eigen::VectorXd &x) override |
Update time dependent quantities. | |
void | update_lagrangian (const Eigen::VectorXd &x, const double k_al) override |
double | compute_error (const Eigen::VectorXd &x) const override |
virtual bool | can_project () const override |
virtual void | project_gradient (Eigen::VectorXd &grad) const override |
virtual void | project_hessian (StiffnessMatrix &hessian) const override |
![]() | |
AugmentedLagrangianForm () | |
virtual | ~AugmentedLagrangianForm () |
void | set_initial_weight (const double k_al) |
double | lagrangian_weight () const |
const StiffnessMatrix & | constraint_matrix () const |
const Eigen::MatrixXd & | constraint_value () const |
const StiffnessMatrix & | constraint_projection_matrix () const |
const Eigen::MatrixXd & | constraint_projection_vector () const |
bool | has_projection () const |
void | set_scale (const double scale) override |
sets the scale for the form | |
![]() | |
virtual | ~Form () |
virtual void | init (const Eigen::VectorXd &x) |
Initialize the 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 | 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 | 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) |
Private Member Functions | |
void | init_masked_lumped_mass (const StiffnessMatrix &mass, const size_t obstacle_ndof) |
Initialize the masked lumped mass matrix. | |
void | update_target (const double t) |
Update target x to the Dirichlet boundary values at time t. | |
Private Attributes | |
const std::vector< int > & | boundary_nodes_ |
const std::vector< mesh::LocalBoundary > * | local_boundary_ |
const std::vector< mesh::LocalBoundary > * | local_neumann_boundary_ |
const int | n_boundary_samples_ |
const int | n_dofs_ |
Eigen::VectorXi | constraints_ |
Constraints. | |
Eigen::VectorXi | not_constraints_ |
Not Constraints. | |
const assembler::RhsAssembler * | rhs_assembler_ |
Reference to the RHS assembler. | |
const bool | is_time_dependent_ |
StiffnessMatrix | masked_lumped_mass_sqrt_ |
sqrt mass matrix masked by the AL dofs | |
StiffnessMatrix | masked_lumped_mass_ |
mass matrix masked by the AL dofs | |
Additional Inherited Members | |
![]() | |
double | L_weight () const |
double | A_weight () const |
![]() | |
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. | |
![]() | |
double | k_al_ |
penalty parameter | |
Eigen::VectorXd | lagr_mults_ |
vector of lagrange multipliers | |
StiffnessMatrix | A_ |
Constraints matrix. | |
Eigen::MatrixXd | b_ |
Constraints value. | |
StiffnessMatrix | A_proj_ |
Constraints projection matrix. | |
Eigen::MatrixXd | b_proj_ |
Constraints projection value. | |
![]() | |
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 of the augmented lagrangian for bc constraints.
Definition at line 12 of file BCLagrangianForm.hpp.
polyfem::solver::BCLagrangianForm::BCLagrangianForm | ( | const int | ndof, |
const std::vector< int > & | boundary_nodes, | ||
const std::vector< mesh::LocalBoundary > & | local_boundary, | ||
const std::vector< mesh::LocalBoundary > & | local_neumann_boundary, | ||
const int | n_boundary_samples, | ||
const StiffnessMatrix & | mass, | ||
const assembler::RhsAssembler & | rhs_assembler, | ||
const size_t | obstacle_ndof, | ||
const bool | is_time_dependent, | ||
const double | t | ||
) |
Construct a new BCLagrangianForm object with a time dependent Dirichlet boundary.
ndof | Number of degrees of freedom |
boundary_nodes | DoFs that are part of the Dirichlet boundary |
local_boundary | |
local_neumann_boundary | |
n_boundary_samples | |
mass | Mass matrix |
rhs_assembler | Right hand side assembler |
obstacle_ndof | Obstacle's number of DOF |
is_time_dependent | Whether the problem is time dependent |
t | Current time |
Definition at line 8 of file BCLagrangianForm.cpp.
References init_masked_lumped_mass(), and update_target().
polyfem::solver::BCLagrangianForm::BCLagrangianForm | ( | const int | ndof, |
const std::vector< int > & | boundary_nodes, | ||
const StiffnessMatrix & | mass, | ||
const size_t | obstacle_ndof, | ||
const Eigen::MatrixXd & | target_x | ||
) |
Construct a new BCLagrangianForm object with a fixed Dirichlet boundary.
ndof | Number of degrees of freedom |
boundary_nodes | DoFs that are part of the Dirichlet boundary |
mass | Mass matrix |
obstacle_ndof | Obstacle's number of DOF |
target_x | Target values for the boundary DoFs |
Definition at line 30 of file BCLagrangianForm.cpp.
References polyfem::solver::AugmentedLagrangianForm::b_, polyfem::solver::AugmentedLagrangianForm::b_proj_, constraints_, and init_masked_lumped_mass().
|
overridevirtual |
Reimplemented from polyfem::solver::AugmentedLagrangianForm.
Definition at line 145 of file BCLagrangianForm.cpp.
|
overridevirtual |
Implements polyfem::solver::AugmentedLagrangianForm.
Definition at line 186 of file BCLagrangianForm.cpp.
References polyfem::solver::AugmentedLagrangianForm::A_, polyfem::solver::AugmentedLagrangianForm::b_, and x.
|
overridevirtual |
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 170 of file BCLagrangianForm.cpp.
References polyfem::solver::AugmentedLagrangianForm::A_, polyfem::solver::AugmentedLagrangianForm::A_weight(), polyfem::solver::AugmentedLagrangianForm::b_, polyfem::solver::AugmentedLagrangianForm::L_weight(), polyfem::solver::AugmentedLagrangianForm::lagr_mults_, masked_lumped_mass_, masked_lumped_mass_sqrt_, and x.
|
private |
Initialize the masked lumped mass matrix.
mass | Mass matrix |
obstacle_ndof | Obstacle's number of DOF |
Definition at line 50 of file BCLagrangianForm.cpp.
References polyfem::solver::AugmentedLagrangianForm::A_, polyfem::solver::AugmentedLagrangianForm::A_proj_, boundary_nodes_, constraints_, polyfem::solver::AugmentedLagrangianForm::lagr_mults_, polyfem::logger(), polyfem::utils::lump_matrix(), masked_lumped_mass_, masked_lumped_mass_sqrt_, n_dofs_, not_constraints_, and polyfem::utils::sparse_identity().
Referenced by BCLagrangianForm(), and BCLagrangianForm().
|
inlineoverridevirtual |
Implements polyfem::solver::Form.
Definition at line 37 of file BCLagrangianForm.hpp.
|
overridevirtual |
Reimplemented from polyfem::solver::AugmentedLagrangianForm.
Definition at line 147 of file BCLagrangianForm.cpp.
References not_constraints_.
|
overridevirtual |
Reimplemented from polyfem::solver::AugmentedLagrangianForm.
Definition at line 155 of file BCLagrangianForm.cpp.
References not_constraints_.
|
overridevirtual |
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 175 of file BCLagrangianForm.cpp.
References polyfem::solver::AugmentedLagrangianForm::A_, polyfem::solver::AugmentedLagrangianForm::A_weight(), and masked_lumped_mass_.
|
overridevirtual |
Implements polyfem::solver::AugmentedLagrangianForm.
Definition at line 206 of file BCLagrangianForm.cpp.
References polyfem::solver::AugmentedLagrangianForm::A_, polyfem::solver::AugmentedLagrangianForm::b_, polyfem::solver::AugmentedLagrangianForm::k_al_, polyfem::solver::AugmentedLagrangianForm::lagr_mults_, masked_lumped_mass_sqrt_, and x.
|
overridevirtual |
Update time dependent quantities.
t | New time |
x | Solution at time t |
Reimplemented from polyfem::solver::Form.
Definition at line 180 of file BCLagrangianForm.cpp.
References is_time_dependent_, and update_target().
|
private |
Update target x to the Dirichlet boundary values at time t.
t | Current time |
Definition at line 193 of file BCLagrangianForm.cpp.
References polyfem::solver::AugmentedLagrangianForm::b_, polyfem::solver::AugmentedLagrangianForm::b_proj_, boundary_nodes_, constraints_, local_boundary_, local_neumann_boundary_, n_boundary_samples_, n_dofs_, rhs_assembler_, and polyfem::assembler::RhsAssembler::set_bc().
Referenced by BCLagrangianForm(), and update_quantities().
|
overridevirtual |
Compute the value of the form.
x | Current solution |
Implements polyfem::solver::Form.
Definition at line 161 of file BCLagrangianForm.cpp.
References polyfem::solver::AugmentedLagrangianForm::A_, polyfem::solver::AugmentedLagrangianForm::A_weight(), polyfem::solver::AugmentedLagrangianForm::b_, polyfem::solver::AugmentedLagrangianForm::L_weight(), polyfem::solver::AugmentedLagrangianForm::lagr_mults_, masked_lumped_mass_, masked_lumped_mass_sqrt_, and x.
|
private |
Definition at line 84 of file BCLagrangianForm.hpp.
Referenced by init_masked_lumped_mass(), and update_target().
|
private |
Constraints.
Definition at line 89 of file BCLagrangianForm.hpp.
Referenced by BCLagrangianForm(), init_masked_lumped_mass(), and update_target().
|
private |
Definition at line 93 of file BCLagrangianForm.hpp.
Referenced by update_quantities().
|
private |
Definition at line 85 of file BCLagrangianForm.hpp.
Referenced by update_target().
|
private |
Definition at line 86 of file BCLagrangianForm.hpp.
Referenced by update_target().
|
private |
mass matrix masked by the AL dofs
Definition at line 96 of file BCLagrangianForm.hpp.
Referenced by first_derivative_unweighted(), init_masked_lumped_mass(), second_derivative_unweighted(), and value_unweighted().
|
private |
sqrt mass matrix masked by the AL dofs
Definition at line 95 of file BCLagrangianForm.hpp.
Referenced by first_derivative_unweighted(), init_masked_lumped_mass(), update_lagrangian(), and value_unweighted().
|
private |
Definition at line 87 of file BCLagrangianForm.hpp.
Referenced by update_target().
|
private |
Definition at line 88 of file BCLagrangianForm.hpp.
Referenced by init_masked_lumped_mass(), and update_target().
|
private |
Not Constraints.
Definition at line 90 of file BCLagrangianForm.hpp.
Referenced by init_masked_lumped_mass(), project_gradient(), and project_hessian().
|
private |
Reference to the RHS assembler.
Definition at line 92 of file BCLagrangianForm.hpp.
Referenced by update_target().