PolyFEM
Loading...
Searching...
No Matches
LaggedRegForm.cpp
Go to the documentation of this file.
1#include "LaggedRegForm.hpp"
2
4
6
7namespace polyfem::solver
8{
9 LaggedRegForm::LaggedRegForm(const int n_lagging_iters)
10 : n_lagging_iters_(n_lagging_iters < 0 ? std::numeric_limits<int>::max() : n_lagging_iters)
11 {
12 if (n_lagging_iters_ < 1)
13 disable();
14 }
15
16 double LaggedRegForm::value_unweighted(const Eigen::VectorXd &x) const
17 {
18 return 0.5 * (x - x_lagged_).squaredNorm();
19 }
20
21 void LaggedRegForm::first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
22 {
23 gradv = (x - x_lagged_);
24 }
25
26 void LaggedRegForm::second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const
27 {
28 hessian.resize(x.size(), x.size());
29 hessian.setIdentity();
30 }
31
32 void LaggedRegForm::init_lagging(const Eigen::VectorXd &x)
33 {
34 update_lagging(x, 0);
35 }
36
37 void LaggedRegForm::update_lagging(const Eigen::VectorXd &x, const int iter_num)
38 {
39 x_lagged_ = x;
40
41 const bool enabled_before = enabled();
42 set_enabled(iter_num >= 0 && iter_num < n_lagging_iters_);
43 if (!enabled_before && enabled())
44 logger().debug("Enabling lagged regularization");
45 else if (enabled_before && !enabled())
46 logger().debug("Disabling lagged regularization");
47 }
48} // namespace polyfem::solver
int x
void disable()
Disable the form.
Definition Form.hpp:117
bool enabled() const
Determine if the form is enabled.
Definition Form.hpp:123
void set_enabled(const bool enabled)
Set if the form is enabled.
Definition Form.hpp:119
void init_lagging(const Eigen::VectorXd &x) override
Initialize lagged fields.
void update_lagging(const Eigen::VectorXd &x, const int iter_num) override
Update lagged fields.
void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
LaggedRegForm(const int n_lagging_iters)
Construct a new Lagged Regularization Form object.
int n_lagging_iters_
Number of iterations to lag for.
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.
Eigen::VectorXd x_lagged_
The full variables from the previous lagging solve.
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22