PolyFEM
Loading...
Searching...
No Matches
TransientForm.cpp
Go to the documentation of this file.
2
3#include <polyfem/State.hpp>
6
7#include <Eigen/Core>
8
9#include <cassert>
10#include <vector>
11
12namespace polyfem::solver
13{
15 {
16 std::vector<double> weights;
17 weights.assign(time_steps_ + 1, dt_);
18 if (transient_integral_type_ == "uniform")
19 {
20 weights[0] = 0;
21 }
22 else if (transient_integral_type_ == "trapezoidal")
23 {
24 weights[0] = dt_ / 2.;
25 weights[weights.size() - 1] = dt_ / 2.;
26 }
27 else if (transient_integral_type_ == "simpson")
28 {
29 weights[0] = dt_ / 3.;
30 weights[weights.size() - 1] = dt_ / 3.;
31 for (int i = 1; i < weights.size() - 1; i++)
32 {
33 if (i % 2)
34 weights[i] = dt_ * 4. / 3.;
35 else
36 weights[i] = dt_ * 2. / 4.;
37 }
38 }
39 else if (transient_integral_type_ == "final")
40 {
41 weights.assign(time_steps_ + 1, 0);
42 weights[time_steps_] = 1;
43 }
44 else if (transient_integral_type_ == "steps")
45 {
46 weights.assign(time_steps_ + 1, 0);
47 for (const int step : steps_)
48 {
49 assert(step > 0 && step < weights.size());
50 weights[step] += 1. / steps_.size();
51 }
52 }
53 else
54 assert(false);
55
56 return weights;
57 }
58
59 double TransientForm::value_unweighted(const Eigen::VectorXd &x) const
60 {
61 std::vector<double> weights = get_transient_quadrature_weights();
62
63 double value = 0;
64 for (int i = 0; i < time_steps_ + 1; i++)
65 {
66 if (weights[i] == 0)
67 continue;
68 const double tmp = obj_->value_unweighted_step(i, x);
69 value += (weights[i] * obj_->weight()) * tmp;
70 }
71
72 return value;
73 }
74 Eigen::MatrixXd TransientForm::compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const
75 {
76 Eigen::MatrixXd terms;
77 terms.setZero(state.ndof(), time_steps_ + 1);
78 std::vector<double> weights = get_transient_quadrature_weights();
79
80 for (int i = 0; i < time_steps_ + 1; i++)
81 {
82 if (weights[i] == 0)
83 continue;
84 terms.col(i) = weights[i] * obj_->compute_adjoint_rhs_step(i, x, state, diff_cache);
85 if (obj_->depends_on_step_prev() && i > 0)
86 terms.col(i - 1) = weights[i] * obj_->compute_adjoint_rhs_step_prev(i, x, state, diff_cache);
87 }
88
89 return terms * weight();
90 }
91 void TransientForm::compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
92 {
93 gradv.setZero(x.size());
94 std::vector<double> weights = get_transient_quadrature_weights();
95
96 Eigen::VectorXd tmp;
97 for (int i = 0; i < time_steps_ + 1; i++)
98 {
99 if (weights[i] == 0)
100 continue;
101 obj_->compute_partial_gradient_step(i, x, tmp);
102 gradv += weights[i] * tmp;
103 }
104
105 gradv *= weight();
106 }
107
108 void TransientForm::init(const Eigen::VectorXd &x)
109 {
110 obj_->init(x);
111 }
112
113 bool TransientForm::is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
114 {
115 return obj_->is_step_valid(x0, x1);
116 }
117
118 double TransientForm::max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
119 {
120 return obj_->max_step_size(x0, x1);
121 }
122
123 void TransientForm::line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1)
124 {
125 obj_->line_search_begin(x0, x1);
126 }
127
129 {
130 obj_->line_search_end();
131 }
132
133 void TransientForm::post_step(const polysolve::nonlinear::PostStepData &data)
134 {
135 obj_->post_step(data);
136 }
137
138 void TransientForm::solution_changed(const Eigen::VectorXd &new_x)
139 {
141 // obj_->solution_changed(new_x);
142 std::vector<double> weights = get_transient_quadrature_weights();
143 for (int i = 0; i <= time_steps_; i++)
144 {
145 if (weights[i] == 0)
146 continue;
147 obj_->solution_changed_step(i, new_x);
148 }
149 }
150 bool TransientForm::is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
151 {
152 return obj_->is_step_collision_free(x0, x1);
153 }
154
155 double ProxyTransientForm::value_unweighted(const Eigen::VectorXd &x) const
156 {
157 Eigen::VectorXd vals(steps_.size());
158 int j = 0;
159 for (int i : steps_)
160 {
161 vals(j++) = obj_->value_unweighted_step(i, x);
162 }
163
164 return eval(vals);
165 }
166 Eigen::MatrixXd ProxyTransientForm::compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const
167 {
168 Eigen::VectorXd vals(steps_.size());
169 Eigen::MatrixXd terms;
170 terms.setZero(state.ndof(), steps_.size());
171
172 int j = 0;
173 for (int i : steps_)
174 {
175 vals(j) = obj_->value_unweighted_step(i, x);
176 terms.col(j++) = obj_->compute_adjoint_rhs_step(i, x, state, diff_cache);
177 }
178
179 const Eigen::VectorXd g = eval_grad(vals);
180 Eigen::MatrixXd out;
181 out.setZero(state.ndof(), time_steps_ + 1);
182 j = 0;
183 for (int i : steps_)
184 {
185 out.col(i) = terms.col(j) * (g(j) * weight());
186 j++;
187 }
188
189 return out;
190 }
191 void ProxyTransientForm::compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
192 {
193 Eigen::VectorXd vals(steps_.size());
194 Eigen::MatrixXd terms;
195 terms.setZero(x.size(), steps_.size());
196
197 int j = 0;
198 Eigen::VectorXd tmp;
199 for (int i : steps_)
200 {
201 vals(j) = obj_->value_unweighted_step(i, x);
202 obj_->compute_partial_gradient_step(i, x, tmp);
203 terms.col(j++) = tmp;
204 }
205
206 gradv = terms * eval_grad(vals) * weight();
207 }
208
209 double ProxyTransientForm::eval(const Eigen::VectorXd &y) const
210 {
211 return 1. / y.array().inverse().sum();
212 }
213
214 Eigen::VectorXd ProxyTransientForm::eval_grad(const Eigen::VectorXd &y) const
215 {
216 return y.array().pow(-2.0) * pow(eval(y), 2);
217 }
218} // namespace polyfem::solver
ElementAssemblyValues vals
Definition Assembler.cpp:22
int y
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
int ndof() const
Definition State.hpp:309
virtual void solution_changed(const Eigen::VectorXd &new_x) override
Update cached fields upon a change in the solution.
double value(const Eigen::VectorXd &x) const override
Compute the value of the form multiplied with the weigth.
virtual double weight() const
Get the form's multiplicative constant weight.
Definition Form.hpp:128
Eigen::MatrixXd compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const override
Eigen::VectorXd eval_grad(const Eigen::VectorXd &y) const
double eval(const Eigen::VectorXd &y) const
void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form.
virtual double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form.
void solution_changed(const Eigen::VectorXd &new_x) override
Update cached fields upon a change in the solution.
double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Determine the maximum step size allowable between the current and next solution.
virtual Eigen::MatrixXd compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const override
bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Checks if the step is collision free.
virtual void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
bool is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Determine if a step from solution x0 to solution x1 is allowed.
std::shared_ptr< StaticForm > obj_
void init(const Eigen::VectorXd &x) override
Initialize the form.
void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override
Initialize variables used during the line search.
std::vector< double > get_transient_quadrature_weights() const
void line_search_end() override
Clear variables used during the line search.
void post_step(const polysolve::nonlinear::PostStepData &data) override
Update fields after a step in the optimization.