PolyFEM
Loading...
Searching...
No Matches
GenericLagrangianForm.cpp
Go to the documentation of this file.
2
4
5namespace polyfem::solver
6{
7 GenericLagrangianForm::GenericLagrangianForm(const std::vector<int> &constraint_nodes,
8 const StiffnessMatrix &A,
9 const Eigen::VectorXd &b)
10 : AugmentedLagrangianForm(constraint_nodes), A(A), b(b)
11 {
12 }
13
14 double GenericLagrangianForm::value_unweighted(const Eigen::VectorXd &x) const
15 {
16 const Eigen::VectorXd res = A * x - b;
17 const double L_penalty = lagr_mults_.transpose() * res;
18 const double A_penalty = res.squaredNorm() / 2;
19 return L_penalty + k_al_ * A_penalty;
20 }
21
22 void GenericLagrangianForm::first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
23 {
24 gradv = A * lagr_mults_ + k_al_ * (A * x - b);
25 }
26
27 void GenericLagrangianForm::second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const
28 {
29 hessian = A;
30 }
31
32 void GenericLagrangianForm::update_lagrangian(const Eigen::VectorXd &x, const double k_al)
33 {
34 lagr_mults_ += k_al * (A * x - b);
35 }
36
37 double GenericLagrangianForm::compute_error(const Eigen::VectorXd &x) const
38 {
39 const Eigen::VectorXd res = A * x - b;
40 return res.norm();
41 }
42} // namespace polyfem::solver
int x
Eigen::VectorXd lagr_mults_
vector of lagrange multipliers
void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form.
double compute_error(const Eigen::VectorXd &x) const override
void update_lagrangian(const Eigen::VectorXd &x, const double k_al) override
GenericLagrangianForm(const std::vector< int > &constraint_nodes, const StiffnessMatrix &A, const Eigen::VectorXd &b)
Construct a new GenericLagrangianForm object for the constraints Ax = b.
void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Compute the first derivative of the value wrt x.
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22