Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
MacroStrainLagrangianForm.cpp
Go to the documentation of this file.
3
4namespace polyfem::solver
5{
7 : macro_strain_constraint_(macro_strain_constraint)
8 {
10 }
11
12 double MacroStrainLagrangianForm::value_unweighted(const Eigen::VectorXd &x) const
13 {
14 const Eigen::VectorXi indices = macro_strain_constraint_.get_fixed_entry().array() + (x.size() - macro_strain_constraint_.dim() * macro_strain_constraint_.dim());
15 const double larg = lagr_mults_.transpose() * (x(indices) - values);
16 const double pen = (x(indices) - values).squaredNorm() / 2.0;
17
18 return L_weight() * larg + A_weight() * pen;
19 }
20
21 void MacroStrainLagrangianForm::first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
22 {
23 const Eigen::VectorXi indices = macro_strain_constraint_.get_fixed_entry().array() + (x.size() - macro_strain_constraint_.dim() * macro_strain_constraint_.dim());
24 gradv.setZero(x.size());
25 gradv(indices) = L_weight() * lagr_mults_ + A_weight() * (x(indices) - values);
26 }
27
29 {
30 const Eigen::VectorXi indices = macro_strain_constraint_.get_fixed_entry().array() + (x.size() - macro_strain_constraint_.dim() * macro_strain_constraint_.dim());
31 hessian.resize(x.size(), x.size());
32 hessian.setZero();
33 for (int i = 0; i < indices.size(); i++)
34 hessian.coeffRef(indices(i), indices(i)) += A_weight();
35 }
36
37 void MacroStrainLagrangianForm::update_quantities(const double t, const Eigen::VectorXd &)
38 {
41 }
42
43 double MacroStrainLagrangianForm::compute_error(const Eigen::VectorXd &x) const
44 {
45 const Eigen::VectorXi indices = macro_strain_constraint_.get_fixed_entry().array() + (x.size() - macro_strain_constraint_.dim() * macro_strain_constraint_.dim());
46 return (x(indices).array() - values.array()).matrix().squaredNorm();
47 }
48
49 void MacroStrainLagrangianForm::update_lagrangian(const Eigen::VectorXd &x, const double k_al)
50 {
51 k_al_ = k_al;
52
53 const Eigen::VectorXi indices = macro_strain_constraint_.get_fixed_entry().array() + (x.size() - macro_strain_constraint_.dim() * macro_strain_constraint_.dim());
54 lagr_mults_ += k_al * (x(indices) - values);
55 }
56} // namespace polyfem::solver
int x
Eigen::MatrixXd eval(const double t) const
const Eigen::VectorXi & get_fixed_entry() const
Eigen::VectorXd lagr_mults_
vector of lagrange multipliers
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.
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the contact barrier potential value.
void update_lagrangian(const Eigen::VectorXd &x, const double k_al) override
const assembler::MacroStrainValue & macro_strain_constraint_
double compute_error(const Eigen::VectorXd &x) const override
void update_quantities(const double t, const Eigen::VectorXd &x) override
Update time-dependent fields.
MacroStrainLagrangianForm(const assembler::MacroStrainValue &macro_strain_constraint)
Construct a new MacroStrainLagrangianForm object.
Eigen::VectorXd flatten(const Eigen::MatrixXd &X)
Flatten rowwises.
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22