PolyFEM
Loading...
Searching...
No Matches
CompositeForm.cpp
Go to the documentation of this file.
2
3#include <polyfem/State.hpp>
5
6#include <Eigen/Core>
7
8#include <algorithm>
9
10namespace polyfem::solver
11{
12 Eigen::MatrixXd CompositeForm::compute_reduced_adjoint_rhs(const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const
13 {
14 Eigen::VectorXd composite_grad = compose_grad(get_inputs(x));
15
16 Eigen::MatrixXd term;
17 Eigen::VectorXd tmp_grad;
18 for (int i = 0; i < forms_.size(); i++)
19 {
20 if (i == 0)
21 term = composite_grad(i) * forms_[i]->compute_reduced_adjoint_rhs(x, state, diff_cache);
22 else
23 term += composite_grad(i) * forms_[i]->compute_reduced_adjoint_rhs(x, state, diff_cache);
24 }
25
26 return term * weight();
27 }
28
29 void CompositeForm::compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
30 {
31 Eigen::VectorXd composite_grad = compose_grad(get_inputs(x));
32
33 gradv.setZero(x.size());
34 Eigen::VectorXd tmp_grad;
35 for (int i = 0; i < forms_.size(); i++)
36 {
37 forms_[i]->compute_partial_gradient(x, tmp_grad);
38 gradv += composite_grad(i) * tmp_grad;
39 }
40 gradv *= weight();
41 }
42
43 Eigen::VectorXd CompositeForm::get_inputs(const Eigen::VectorXd &x) const
44 {
45 Eigen::VectorXd values;
46 values.setZero(forms_.size());
47
48 for (int i = 0; i < forms_.size(); i++)
49 values(i) = forms_[i]->value(x);
50
51 return values;
52 }
53
54 double CompositeForm::value_unweighted(const Eigen::VectorXd &x) const
55 {
56 return compose(get_inputs(x));
57 }
58
59 void CompositeForm::init(const Eigen::VectorXd &x)
60 {
61 for (const auto &f : forms_)
62 f->init(x);
63 }
64
65 bool CompositeForm::is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
66 {
67 for (const auto &f : forms_)
68 {
69 if (f->enabled() && !f->is_step_valid(x0, x1))
70 return false;
71 }
72 return true;
73 }
74
75 double CompositeForm::max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
76 {
77 double step = 1;
78 for (const auto &f : forms_)
79 if (f->enabled())
80 step = std::min(step, f->max_step_size(x0, x1));
81
82 return step;
83 }
84
85 void CompositeForm::line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1)
86 {
87 for (const auto &f : forms_)
88 f->line_search_begin(x0, x1);
89 }
90
92 {
93 for (const auto &f : forms_)
94 f->line_search_end();
95 }
96
97 void CompositeForm::post_step(const polysolve::nonlinear::PostStepData &data)
98 {
99 for (const auto &f : forms_)
100 f->post_step(data);
101 }
102
103 void CompositeForm::solution_changed(const Eigen::VectorXd &new_x)
104 {
106 for (const auto &f : forms_)
107 f->solution_changed(new_x);
108 }
109 bool CompositeForm::is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
110 {
111 for (const auto &f : forms_)
112 {
113 if (f->enabled() && !f->is_step_collision_free(x0, x1))
114 return false;
115 }
116 return true;
117 }
118} // namespace polyfem::solver
int x
Storage for additional data required by differntial code.
Definition DiffCache.hpp:21
main class that contains the polyfem solver and all its state
Definition State.hpp:113
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 Eigen::MatrixXd compute_reduced_adjoint_rhs(const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const override final
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 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