PolyFEM
Loading...
Searching...
No Matches
CompositeForm.cpp
Go to the documentation of this file.
1#include "CompositeForm.hpp"
2#include <polyfem/State.hpp>
3
4namespace polyfem::solver
5{
6 Eigen::MatrixXd CompositeForm::compute_reduced_adjoint_rhs(const Eigen::VectorXd &x, const State &state) const
7 {
8 Eigen::VectorXd composite_grad = compose_grad(get_inputs(x));
9
10 Eigen::MatrixXd term;
11 Eigen::VectorXd tmp_grad;
12 for (int i = 0; i < forms_.size(); i++)
13 {
14 if (i == 0)
15 term = composite_grad(i) * forms_[i]->compute_reduced_adjoint_rhs(x, state);
16 else
17 term += composite_grad(i) * forms_[i]->compute_reduced_adjoint_rhs(x, state);
18 }
19
20 return term * weight();
21 }
22
23 void CompositeForm::compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
24 {
25 Eigen::VectorXd composite_grad = compose_grad(get_inputs(x));
26
27 gradv.setZero(x.size());
28 Eigen::VectorXd tmp_grad;
29 for (int i = 0; i < forms_.size(); i++)
30 {
31 forms_[i]->compute_partial_gradient(x, tmp_grad);
32 gradv += composite_grad(i) * tmp_grad;
33 }
34 gradv *= weight();
35 }
36
37 Eigen::VectorXd CompositeForm::get_inputs(const Eigen::VectorXd &x) const
38 {
39 Eigen::VectorXd values;
40 values.setZero(forms_.size());
41
42 for (int i = 0; i < forms_.size(); i++)
43 values(i) = forms_[i]->value(x);
44
45 return values;
46 }
47
48 double CompositeForm::value_unweighted(const Eigen::VectorXd &x) const
49 {
50 return compose(get_inputs(x));
51 }
52
53 void CompositeForm::init(const Eigen::VectorXd &x)
54 {
55 for (const auto &f : forms_)
56 f->init(x);
57 }
58
59 bool CompositeForm::is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
60 {
61 for (const auto &f : forms_)
62 {
63 if (f->enabled() && !f->is_step_valid(x0, x1))
64 return false;
65 }
66 return true;
67 }
68
69 double CompositeForm::max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
70 {
71 double step = 1;
72 for (const auto &f : forms_)
73 if (f->enabled())
74 step = std::min(step, f->max_step_size(x0, x1));
75
76 return step;
77 }
78
79 void CompositeForm::line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1)
80 {
81 for (const auto &f : forms_)
82 f->line_search_begin(x0, x1);
83 }
84
86 {
87 for (const auto &f : forms_)
88 f->line_search_end();
89 }
90
91 void CompositeForm::post_step(const polysolve::nonlinear::PostStepData &data)
92 {
93 for (const auto &f : forms_)
94 f->post_step(data);
95 }
96
97 void CompositeForm::solution_changed(const Eigen::VectorXd &new_x)
98 {
100 for (const auto &f : forms_)
101 f->solution_changed(new_x);
102 }
103 bool CompositeForm::is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
104 {
105 for (const auto &f : forms_)
106 {
107 if (f->enabled() && !f->is_step_collision_free(x0, x1))
108 return false;
109 }
110 return true;
111 }
112} // namespace polyfem::solver
int x
main class that contains the polyfem solver and all its state
Definition State.hpp:79
virtual void solution_changed(const Eigen::VectorXd &new_x) override
Update cached fields upon a change in the solution.
virtual void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override final
std::vector< std::shared_ptr< AdjointForm > > forms_
virtual bool is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const final override
Determine if a step from solution x0 to solution x1 is allowed.
virtual double compose(const Eigen::VectorXd &inputs) const =0
virtual void init(const Eigen::VectorXd &x) final override
Initialize the form.
virtual double value_unweighted(const Eigen::VectorXd &x) const final override
Compute the value of the form.
virtual double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const final override
Determine the maximum step size allowable between the current and next solution.
Eigen::VectorXd get_inputs(const Eigen::VectorXd &x) const
virtual Eigen::MatrixXd compute_reduced_adjoint_rhs(const Eigen::VectorXd &x, const State &state) const override final
virtual void post_step(const polysolve::nonlinear::PostStepData &data) final override
Update fields after a step in the optimization.
virtual Eigen::VectorXd compose_grad(const Eigen::VectorXd &inputs) const =0
virtual bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const final override
Checks if the step is collision free.
virtual void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) final override
Initialize variables used during the line search.
virtual void solution_changed(const Eigen::VectorXd &new_x) final override
Update cached fields upon a change in the solution.
virtual void line_search_end() final override
Clear variables used during the line search.
virtual double weight() const
Get the form's multiplicative constant weight.
Definition Form.hpp:128