PolyFEM
Loading...
Searching...
No Matches
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 return lagr_mults_.transpose() * (x(indices) - values);
16 }
17
18 void MacroStrainLagrangianForm::first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
19 {
20 const Eigen::VectorXi indices = macro_strain_constraint_.get_fixed_entry().array() + (x.size() - macro_strain_constraint_.dim() * macro_strain_constraint_.dim());
21 gradv.setZero(x.size());
22 gradv(indices) = lagr_mults_;
23 }
24
26 {
27 hessian.resize(x.size(), x.size());
28 hessian.setZero();
29 }
30
31 void MacroStrainLagrangianForm::update_quantities(const double t, const Eigen::VectorXd &)
32 {
35 }
36
37 void MacroStrainLagrangianForm::update_lagrangian(const Eigen::VectorXd &x, const double k_al)
38 {
39 const Eigen::VectorXi indices = macro_strain_constraint_.get_fixed_entry().array() + (x.size() - macro_strain_constraint_.dim() * macro_strain_constraint_.dim());
40 lagr_mults_ += k_al * (x(indices) - values);
41 }
42}
int x
Eigen::MatrixXd eval(const double t) const
const Eigen::VectorXi & get_fixed_entry() const
void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Compute the first derivative of the value wrt x.
void update_lagrangian(const Eigen::VectorXd &x, const double k_al)
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.
Eigen::VectorXd lagr_mults_
vector of lagrange multipliers
const assembler::MacroStrainValue & macro_strain_constraint_
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