PolyFEM
Loading...
Searching...
No Matches
NLHomoProblem.cpp
Go to the documentation of this file.
1#include "NLHomoProblem.hpp"
2#include <polyfem/State.hpp>
7
8namespace polyfem::solver
9{
10 NLHomoProblem::NLHomoProblem(const int full_size,
11 const assembler::MacroStrainValue &macro_strain_constraint,
12 const State &state,
13 const double t,
14 const std::vector<std::shared_ptr<Form>> &forms,
15 const std::vector<std::shared_ptr<AugmentedLagrangianForm>> &penalty_forms,
16 const bool solve_symmetric_macro_strain)
17 : NLProblem(full_size, state.periodic_bc, t, forms, penalty_forms),
18 state_(state),
19 only_symmetric(solve_symmetric_macro_strain),
20 macro_strain_constraint_(macro_strain_constraint)
21 {
23 }
24
26 {
27 const int dim = state_.mesh->dimension();
29 {
30 macro_mid_to_full_.setZero(dim * dim, (dim * (dim + 1)) / 2);
31 macro_full_to_mid_.setZero((dim * (dim + 1)) / 2, dim * dim);
32 for (int i = 0, idx = 0; i < dim; i++)
33 {
34 for (int j = i; j < dim; j++)
35 {
36 macro_full_to_mid_(idx, i * dim + j) = 1;
37
38 macro_mid_to_full_(j * dim + i, idx) = 1;
39 macro_mid_to_full_(i * dim + j, idx) = 1;
40
41 idx++;
42 }
43 }
44 }
45 else
46 {
47 macro_mid_to_full_.setIdentity(dim * dim, dim * dim);
48 macro_full_to_mid_.setIdentity(dim * dim, dim * dim);
49 }
51 }
52
53 Eigen::VectorXd NLHomoProblem::extended_to_reduced(const Eigen::VectorXd &extended) const
54 {
55 const int dim = state_.mesh->dimension();
56 const int dof2 = macro_reduced_size();
57 const int dof1 = reduced_size();
58
59 Eigen::VectorXd reduced(dof1 + dof2);
60 reduced.head(dof1) = NLProblem::full_to_reduced(extended.head(extended.size() - dim * dim));
61 reduced.tail(dof2) = macro_full_to_reduced(extended.tail(dim * dim));
62
63 return reduced;
64 }
65
66 Eigen::VectorXd NLHomoProblem::reduced_to_extended(const Eigen::VectorXd &reduced, bool homogeneous) const
67 {
68 const int dim = state_.mesh->dimension();
69 const int dof2 = macro_reduced_size();
70 const int dof1 = reduced_size();
71 assert(reduced.size() == dof1 + dof2);
72
73 Eigen::VectorXd fluctuation = NLProblem::reduced_to_full(reduced.head(dof1));
74 Eigen::VectorXd disp_grad = macro_reduced_to_full(reduced.tail(dof2), homogeneous);
75 Eigen::VectorXd extended(fluctuation.size() + disp_grad.size());
76 extended << fluctuation, disp_grad;
77
78 return extended;
79 }
80
81 Eigen::VectorXd NLHomoProblem::extended_to_reduced_grad(const Eigen::VectorXd &extended) const
82 {
83 const int dim = state_.mesh->dimension();
84 const int dof2 = macro_reduced_size();
85 const int dof1 = reduced_size();
86
87 Eigen::VectorXd grad(dof1 + dof2);
88 grad.head(dof1) = NLProblem::full_to_reduced_grad(extended.head(extended.size() - dim * dim));
89 grad.tail(dof2) = macro_full_to_reduced_grad(extended.tail(dim * dim));
90
91 return grad;
92 }
93
94 NLHomoProblem::TVector NLHomoProblem::reduced_to_full_shape_derivative(const Eigen::MatrixXd &disp_grad, const TVector &adjoint_full) const
95 {
96 const int dim = state_.mesh->dimension();
97
98 Eigen::VectorXd term;
99 term.setZero(state_.n_bases * dim);
100 for (int i = 0; i < state_.n_bases; i++)
101 term.segment(i * dim, dim) += disp_grad.transpose() * adjoint_full.segment(i * dim, dim);
102
104 }
105
106 double NLHomoProblem::value(const TVector &x)
107 {
108 double val = NLProblem::value(x);
109 for (auto &form : homo_forms)
110 if (form->enabled())
111 val += form->value(reduced_to_extended(x));
112
113 return val;
114 }
115 void NLHomoProblem::gradient(const TVector &x, TVector &gradv)
116 {
117 NLProblem::gradient(x, gradv);
118
119 for (auto &form : homo_forms)
120 if (form->enabled())
121 {
122 Eigen::VectorXd grad_extended;
123 form->first_derivative(reduced_to_extended(x), grad_extended);
124 gradv += extended_to_reduced_grad(grad_extended);
125 }
126 }
127 void NLHomoProblem::extended_hessian_to_reduced_hessian(const THessian &extended, THessian &reduced) const
128 {
129 const int dim = state_.mesh->dimension();
130 const int dof2 = macro_reduced_size();
131 const int dof1 = reduced_size();
132
133 Eigen::MatrixXd A12, A22;
134 {
135 Eigen::MatrixXd tmp = Eigen::MatrixXd(extended.rightCols(dim * dim));
136 A12 = macro_full_to_reduced_grad(tmp.topRows(tmp.rows() - dim * dim).transpose()).transpose();
137 A22 = macro_full_to_reduced_grad(macro_full_to_reduced_grad(tmp.bottomRows(dim * dim)).transpose());
138 }
139
140 std::vector<Eigen::Triplet<double>> entries;
141 entries.reserve(extended.nonZeros() + A12.size() * 2 + A22.size());
142
143 for (int k = 0; k < full_size_; ++k)
144 for (StiffnessMatrix::InnerIterator it(extended, k); it; ++it)
145 if (it.row() < full_size_ && it.col() < full_size_)
146 entries.emplace_back(it.row(), it.col(), it.value());
147
148 for (int i = 0; i < A12.rows(); i++)
149 {
150 for (int j = 0; j < A12.cols(); j++)
151 {
152 entries.emplace_back(i, full_size_ + j, A12(i, j));
153 entries.emplace_back(full_size_ + j, i, A12(i, j));
154 }
155 }
156
157 for (int i = 0; i < A22.rows(); i++)
158 for (int j = 0; j < A22.cols(); j++)
159 entries.emplace_back(i + full_size_, j + full_size_, A22(i, j));
160
161 // Eigen::VectorXd tmp = macro_full_to_reduced_grad(Eigen::VectorXd::Ones(dim*dim));
162 // for (int i = 0; i < tmp.size(); i++)
163 // entries.emplace_back(i + full_size_, i + full_size_, 1 - tmp(i));
164
165 THessian mid(full_size_ + dof2, full_size_ + dof2);
166 mid.setFromTriplets(entries.begin(), entries.end());
167
169 }
170 Eigen::MatrixXd NLHomoProblem::reduced_to_disp_grad(const TVector &reduced, bool homogeneous) const
171 {
172 const int dim = state_.mesh->dimension();
173 const int dof2 = macro_reduced_size();
174 const int dof1 = reduced_size();
175
176 return utils::unflatten(macro_reduced_to_full(reduced.tail(dof2), homogeneous), dim);
177 }
178 void NLHomoProblem::hessian(const TVector &x, THessian &hessian)
179 {
181
182 for (auto &form : homo_forms)
183 if (form->enabled())
184 {
185 THessian hess_extended;
186 form->second_derivative(reduced_to_extended(x), hess_extended);
187
188 THessian hess;
189 extended_hessian_to_reduced_hessian(hess_extended, hess);
190 hessian += hess;
191 }
192 }
193
194 void NLHomoProblem::set_fixed_entry(const Eigen::VectorXi &fixed_entry)
195 {
196 const int dim = state_.mesh->dimension();
197
198 Eigen::VectorXd fixed_mask;
199 fixed_mask.setZero(dim * dim);
200 fixed_mask(fixed_entry.array()).setOnes();
201 fixed_mask = (macro_full_to_mid_ * fixed_mask).eval();
202
203 fixed_mask_.setZero(fixed_mask.size());
204 for (int i = 0; i < fixed_mask.size(); i++)
205 if (abs(fixed_mask(i)) > 1e-8)
206 fixed_mask_(i) = true;
207
208 const int new_reduced_size = fixed_mask_.size() - fixed_mask_.sum();
209 macro_mid_to_reduced_.setZero(new_reduced_size, fixed_mask_.size());
210 for (int i = 0, j = 0; i < fixed_mask_.size(); i++)
211 if (!fixed_mask_(i))
212 macro_mid_to_reduced_(j++, i) = 1;
213 }
214
215 void NLHomoProblem::full_hessian_to_reduced_hessian(const THessian &full, THessian &reduced) const
216 {
217 const int dim = state_.mesh->dimension();
218 const int dof2 = macro_reduced_size();
219 const int dof1 = reduced_size();
220
221 Eigen::MatrixXd tmp = constraint_grad();
222 Eigen::MatrixXd A12 = full * tmp.transpose();
223 Eigen::MatrixXd A22 = tmp * A12;
224
225 std::vector<Eigen::Triplet<double>> entries;
226 entries.reserve(full.nonZeros() + A12.size() * 2 + A22.size());
227
228 for (int k = 0; k < full.outerSize(); ++k)
229 for (StiffnessMatrix::InnerIterator it(full, k); it; ++it)
230 entries.emplace_back(it.row(), it.col(), it.value());
231
232 for (int i = 0; i < A12.rows(); i++)
233 for (int j = 0; j < A12.cols(); j++)
234 {
235 entries.emplace_back(i, full.cols() + j, A12(i, j));
236 entries.emplace_back(full.rows() + j, i, A12(i, j));
237 }
238
239 for (int i = 0; i < A22.rows(); i++)
240 for (int j = 0; j < A22.cols(); j++)
241 entries.emplace_back(i + full.rows(), j + full.cols(), A22(i, j));
242
243 THessian mid(full.rows() + dof2, full.cols() + dof2);
244 mid.setFromTriplets(entries.begin(), entries.end());
245
247 }
248
249 NLHomoProblem::TVector NLHomoProblem::full_to_reduced(const TVector &full) const
250 {
251 log_and_throw_error("Invalid function!");
252 return TVector();
253 }
254
255 NLHomoProblem::TVector NLHomoProblem::full_to_reduced(const TVector &full, const Eigen::MatrixXd &disp_grad) const
256 {
257 const int dim = state_.mesh->dimension();
258 const int dof2 = macro_reduced_size();
259 const int dof1 = reduced_size();
260
261 TVector reduced;
262 reduced.setZero(dof1 + dof2);
263
265 reduced.tail(dof2) = macro_full_to_reduced(utils::flatten(disp_grad));
266
267 return reduced;
268 }
269 NLHomoProblem::TVector NLHomoProblem::full_to_reduced_grad(const TVector &full) const
270 {
271 const int dim = state_.mesh->dimension();
272 const int dof2 = macro_reduced_size();
273 const int dof1 = reduced_size();
274
275 TVector reduced;
276 reduced.setZero(dof1 + dof2);
277
278 reduced.head(dof1) = NLProblem::full_to_reduced_grad(full);
279 reduced.tail(dof2) = constraint_grad() * full;
280
281 return reduced;
282 }
283 NLHomoProblem::TVector NLHomoProblem::reduced_to_full(const TVector &reduced) const
284 {
285 const int dim = state_.mesh->dimension();
286 const int dof2 = macro_reduced_size();
287 const int dof1 = reduced_size();
288
289 Eigen::MatrixXd disp_grad = utils::unflatten(macro_reduced_to_full(reduced.tail(dof2)), dim);
291 }
292
294 {
295 return macro_mid_to_reduced_.rows();
296 }
297 NLHomoProblem::TVector NLHomoProblem::macro_full_to_reduced(const TVector &full) const
298 {
300 }
301 Eigen::MatrixXd NLHomoProblem::macro_full_to_reduced_grad(const Eigen::MatrixXd &full) const
302 {
303 return macro_mid_to_reduced_ * macro_mid_to_full_.transpose() * full;
304 }
305 NLHomoProblem::TVector NLHomoProblem::macro_reduced_to_full(const TVector &reduced, bool homogeneous) const
306 {
307 TVector mid = macro_mid_to_reduced_.transpose() * reduced;
308 const TVector fixed_values = homogeneous ? TVector::Zero(macro_full_to_mid_.rows()) : TVector(macro_full_to_mid_ * utils::flatten(macro_strain_constraint_.eval(t_)));
309 for (int i = 0; i < fixed_mask_.size(); i++)
310 if (fixed_mask_(i))
311 mid(i) = fixed_values(i);
312
313 return macro_mid_to_full_ * mid;
314 }
315
316 void NLHomoProblem::init(const TVector &x0)
317 {
318 for (auto &form : homo_forms)
319 form->init(reduced_to_extended(x0));
321 }
322
323 bool NLHomoProblem::is_step_valid(const TVector &x0, const TVector &x1)
324 {
325 bool flag = NLProblem::is_step_valid(x0, x1);
326 for (auto &form : homo_forms)
327 if (form->enabled())
328 flag &= form->is_step_valid(reduced_to_extended(x0), reduced_to_extended(x1));
329
330 return flag;
331 }
332 bool NLHomoProblem::is_step_collision_free(const TVector &x0, const TVector &x1)
333 {
334 bool flag = NLProblem::is_step_collision_free(x0, x1);
335 for (auto &form : homo_forms)
336 if (form->enabled())
337 flag &= form->is_step_collision_free(reduced_to_extended(x0), reduced_to_extended(x1));
338
339 return flag;
340 }
341 double NLHomoProblem::max_step_size(const TVector &x0, const TVector &x1)
342 {
343 double size = NLProblem::max_step_size(x0, x1);
344 for (auto &form : homo_forms)
345 if (form->enabled())
346 size = std::min(size, form->max_step_size(reduced_to_extended(x0), reduced_to_extended(x1)));
347 return size;
348 }
349
350 void NLHomoProblem::line_search_begin(const TVector &x0, const TVector &x1)
351 {
353 for (auto &form : homo_forms)
354 form->line_search_begin(reduced_to_extended(x0), reduced_to_extended(x1));
355 }
356 void NLHomoProblem::post_step(const polysolve::nonlinear::PostStepData &data)
357 {
359 for (auto &form : homo_forms)
360 form->post_step(polysolve::nonlinear::PostStepData(
361 data.iter_num, data.solver_info, reduced_to_extended(data.x), reduced_to_extended(data.grad)));
362 }
363
364 void NLHomoProblem::solution_changed(const TVector &new_x)
365 {
367 for (auto &form : homo_forms)
368 form->solution_changed(reduced_to_extended(new_x));
369 }
370
371 void NLHomoProblem::init_lagging(const TVector &x)
372 {
374 for (auto &form : homo_forms)
375 form->init_lagging(reduced_to_extended(x));
376 }
377 void NLHomoProblem::update_lagging(const TVector &x, const int iter_num)
378 {
379 NLProblem::update_lagging(x, iter_num);
380 for (auto &form : homo_forms)
381 form->update_lagging(reduced_to_extended(x), iter_num);
382 }
383
384 void NLHomoProblem::update_quantities(const double t, const TVector &x)
385 {
387 for (auto &form : homo_forms)
388 form->update_quantities(t, reduced_to_extended(x));
389 }
390
391 Eigen::MatrixXd NLHomoProblem::constraint_grad() const
392 {
393 const int dim = state_.mesh->dimension();
394 Eigen::MatrixXd jac; // (dim*dim) x (dim*n_bases)
395
397
398 jac.setZero(dim * dim, full_size_);
399 for (int i = 0; i < X.rows(); i++)
400 for (int j = 0; j < dim; j++)
401 for (int k = 0; k < dim; k++)
402 jac(j * dim + k, i * dim + j) = X(i, k);
403
404 return macro_full_to_reduced_grad(jac);
405 }
406
407 Eigen::MatrixXd NLHomoProblem::constraint_values(const TVector &) const
408 {
409 Eigen::MatrixXd result = Eigen::MatrixXd::Zero(full_size(), 1);
410 return result;
411 }
412} // namespace polyfem::solver
double val
Definition Assembler.cpp:86
std::vector< Eigen::Triplet< double > > entries
int x
main class that contains the polyfem solver and all its state
Definition State.hpp:79
int n_bases
number of bases
Definition State.hpp:178
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:466
std::shared_ptr< polyfem::mesh::MeshNodes > mesh_nodes
Mapping from input nodes to FE nodes.
Definition State.hpp:193
StiffnessMatrix basis_nodes_to_gbasis_nodes
Definition State.hpp:699
Eigen::MatrixXd eval(const double t) const
static Eigen::MatrixXd generate_linear_field(const int n_bases, const std::shared_ptr< mesh::MeshNodes > mesh_nodes, const Eigen::MatrixXd &grad)
static Eigen::MatrixXd get_bases_position(const int n_bases, const std::shared_ptr< mesh::MeshNodes > mesh_nodes)
virtual void init(const TVector &x0) override
Eigen::MatrixXd reduced_to_disp_grad(const TVector &reduced, bool homogeneous=false) const
void gradient(const TVector &x, TVector &gradv) override
void init(const TVector &x0) override
void post_step(const polysolve::nonlinear::PostStepData &data) override
void set_fixed_entry(const Eigen::VectorXi &fixed_entry)
TVector macro_full_to_reduced(const TVector &full) const
bool is_step_valid(const TVector &x0, const TVector &x1) override
Eigen::MatrixXd constraint_values(const TVector &) const override
const assembler::MacroStrainValue & macro_strain_constraint_
void line_search_begin(const TVector &x0, const TVector &x1) override
void update_quantities(const double t, const TVector &x) override
TVector extended_to_reduced_grad(const TVector &extended) const
void hessian(const TVector &x, THessian &hessian) override
void init_lagging(const TVector &x) override
TVector reduced_to_extended(const TVector &reduced, bool homogeneous=false) const
TVector extended_to_reduced(const TVector &extended) const
NLHomoProblem(const int full_size, const assembler::MacroStrainValue &macro_strain_constraint, const State &state, const double t, const std::vector< std::shared_ptr< Form > > &forms, const std::vector< std::shared_ptr< AugmentedLagrangianForm > > &penalty_forms, const bool solve_symmetric_macro_strain)
std::vector< std::shared_ptr< Form > > homo_forms
Eigen::MatrixXd constraint_grad() const
void full_hessian_to_reduced_hessian(const THessian &full, THessian &reduced) const override
TVector reduced_to_full_shape_derivative(const Eigen::MatrixXd &disp_grad, const TVector &adjoint_full) const
TVector reduced_to_full(const TVector &reduced) const override
Eigen::MatrixXd macro_full_to_reduced_grad(const Eigen::MatrixXd &full) const
double value(const TVector &x) override
TVector full_to_reduced(const TVector &full, const Eigen::MatrixXd &disp_grad) const
TVector full_to_reduced_grad(const TVector &full) const override
void solution_changed(const TVector &new_x) override
double max_step_size(const TVector &x0, const TVector &x1) override
void extended_hessian_to_reduced_hessian(const THessian &extended, THessian &reduced) const
bool is_step_collision_free(const TVector &x0, const TVector &x1) override
void update_lagging(const TVector &x, const int iter_num) override
TVector macro_reduced_to_full(const TVector &reduced, bool homogeneous=false) const
const int full_size_
Size of the full problem.
Definition NLProblem.hpp:65
void line_search_begin(const TVector &x0, const TVector &x1) override
Definition NLProblem.cpp:88
virtual bool is_step_valid(const TVector &x0, const TVector &x1) override
Definition NLProblem.cpp:98
virtual void post_step(const polysolve::nonlinear::PostStepData &data) override
virtual TVector full_to_reduced_grad(const TVector &full) const
virtual TVector full_to_reduced(const TVector &full) const
virtual void update_quantities(const double t, const TVector &x)
Definition NLProblem.cpp:80
virtual void gradient(const TVector &x, TVector &gradv) override
void init_lagging(const TVector &x) override
Definition NLProblem.cpp:70
virtual bool is_step_collision_free(const TVector &x0, const TVector &x1) override
virtual TVector reduced_to_full(const TVector &reduced) const
void update_lagging(const TVector &x, const int iter_num) override
Definition NLProblem.cpp:75
virtual void full_hessian_to_reduced_hessian(const THessian &full, THessian &reduced) const
virtual double value(const TVector &x) override
virtual double max_step_size(const TVector &x0, const TVector &x1) override
Definition NLProblem.cpp:93
virtual void hessian(const TVector &x, THessian &hessian) override
void solution_changed(const TVector &new_x) override
Eigen::MatrixXd unflatten(const Eigen::VectorXd &x, int dim)
Unflatten rowwises, so every dim elements in x become a row.
Eigen::VectorXd flatten(const Eigen::MatrixXd &X)
Flatten rowwises.
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71