PolyFEM
|
Form of the penalty in augmented lagrangian. More...
#include <BCPenaltyForm.hpp>
Public Member Functions | |
BCPenaltyForm (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 BCPenaltyForm object with a time dependent Dirichlet boundary. | |
std::string | name () const override |
BCPenaltyForm (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 BCPenaltyForm 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. | |
StiffnessMatrix & | mask () |
const StiffnessMatrix & | mask () const |
Eigen::VectorXd | target () const |
double | compute_error (const Eigen::VectorXd &x) const |
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 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 void | set_apply_DBC (const Eigen::VectorXd &x, bool apply_DBC) |
Set if the Dirichlet boundary conditions should be enforced. | |
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 int ndof, 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 assembler::RhsAssembler * | rhs_assembler_ |
Reference to the RHS assembler. | |
const bool | is_time_dependent_ |
StiffnessMatrix | masked_lumped_mass_ |
mass matrix masked by the AL dofs | |
StiffnessMatrix | mask_ |
identity matrix masked by the AL dofs | |
Eigen::MatrixXd | target_x_ |
actually a vector with the same size as x with target nodal positions | |
Additional Inherited Members | |
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. | |
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 of the penalty in augmented lagrangian.
Definition at line 12 of file BCPenaltyForm.hpp.
polyfem::solver::BCPenaltyForm::BCPenaltyForm | ( | 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 BCPenaltyForm 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 7 of file BCPenaltyForm.cpp.
References init_masked_lumped_mass(), and update_target().
polyfem::solver::BCPenaltyForm::BCPenaltyForm | ( | 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 BCPenaltyForm 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 28 of file BCPenaltyForm.cpp.
References init_masked_lumped_mass().
double polyfem::solver::BCPenaltyForm::compute_error | ( | const Eigen::VectorXd & | x | ) | const |
Definition at line 104 of file BCPenaltyForm.cpp.
References mask(), target(), 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 88 of file BCPenaltyForm.cpp.
References masked_lumped_mass_, target_x_, and x.
|
private |
Initialize the masked lumped mass matrix.
ndof | Number of degrees of freedom |
mass | Mass matrix |
obstacle | Obstacles |
Definition at line 44 of file BCPenaltyForm.cpp.
References boundary_nodes_, polyfem::utils::lump_matrix(), mask_, masked_lumped_mass_, polyfem::utils::sparse_identity(), and polyfem::solver::Form::value().
Referenced by BCPenaltyForm(), and BCPenaltyForm().
|
inline |
Definition at line 72 of file BCPenaltyForm.hpp.
References mask_.
Referenced by compute_error().
|
inline |
Definition at line 73 of file BCPenaltyForm.hpp.
References mask_.
|
inlineoverridevirtual |
Implements polyfem::solver::Form.
Definition at line 37 of file BCPenaltyForm.hpp.
|
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 93 of file BCPenaltyForm.cpp.
References masked_lumped_mass_.
|
inline |
Definition at line 74 of file BCPenaltyForm.hpp.
References target_x_.
Referenced by compute_error().
|
overridevirtual |
Update time dependent quantities.
t | New time |
x | Solution at time t |
Reimplemented from polyfem::solver::Form.
Definition at line 98 of file BCPenaltyForm.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 109 of file BCPenaltyForm.cpp.
References boundary_nodes_, local_boundary_, local_neumann_boundary_, masked_lumped_mass_, n_boundary_samples_, rhs_assembler_, polyfem::assembler::RhsAssembler::set_bc(), and target_x_.
Referenced by BCPenaltyForm(), and update_quantities().
|
overridevirtual |
Compute the value of the form.
x | Current solution |
Implements polyfem::solver::Form.
Definition at line 81 of file BCPenaltyForm.cpp.
References masked_lumped_mass_, target_x_, and x.
|
private |
Definition at line 79 of file BCPenaltyForm.hpp.
Referenced by init_masked_lumped_mass(), and update_target().
|
private |
Definition at line 85 of file BCPenaltyForm.hpp.
Referenced by update_quantities().
|
private |
Definition at line 80 of file BCPenaltyForm.hpp.
Referenced by update_target().
|
private |
Definition at line 81 of file BCPenaltyForm.hpp.
Referenced by update_target().
|
private |
identity matrix masked by the AL dofs
Definition at line 88 of file BCPenaltyForm.hpp.
Referenced by init_masked_lumped_mass(), mask(), and mask().
|
private |
mass matrix masked by the AL dofs
Definition at line 87 of file BCPenaltyForm.hpp.
Referenced by first_derivative_unweighted(), init_masked_lumped_mass(), second_derivative_unweighted(), update_target(), and value_unweighted().
|
private |
Definition at line 82 of file BCPenaltyForm.hpp.
Referenced by update_target().
|
private |
Reference to the RHS assembler.
Definition at line 84 of file BCPenaltyForm.hpp.
Referenced by update_target().
|
private |
actually a vector with the same size as x with target nodal positions
Definition at line 89 of file BCPenaltyForm.hpp.
Referenced by first_derivative_unweighted(), target(), update_target(), and value_unweighted().