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 const std::shared_ptr<polysolve::linear::Solver> &solver,
18 const double char_length,
19 const double char_force,
20 StiffnessMatrix lumped_mass,
21 const int dimension)
22 : NLProblem(full_size, state.periodic_bc, t, forms, penalty_forms, solver, char_length, char_force, lumped_mass, dimension),
23 state_(state),
24 only_symmetric(solve_symmetric_macro_strain),
25 macro_strain_constraint_(macro_strain_constraint)
26 {
28 }
29
31 {
32 const int dim = state_.mesh->dimension();
34 {
35 macro_mid_to_full_.setZero(dim * dim, (dim * (dim + 1)) / 2);
36 macro_full_to_mid_.setZero((dim * (dim + 1)) / 2, dim * dim);
37 for (int i = 0, idx = 0; i < dim; i++)
38 {
39 for (int j = i; j < dim; j++)
40 {
41 macro_full_to_mid_(idx, i * dim + j) = 1;
42
43 macro_mid_to_full_(j * dim + i, idx) = 1;
44 macro_mid_to_full_(i * dim + j, idx) = 1;
45
46 idx++;
47 }
48 }
49 }
50 else
51 {
52 macro_mid_to_full_.setIdentity(dim * dim, dim * dim);
53 macro_full_to_mid_.setIdentity(dim * dim, dim * dim);
54 }
56 }
57
58 Eigen::VectorXd NLHomoProblem::extended_to_reduced(const Eigen::VectorXd &extended) const
59 {
60 const int dim = state_.mesh->dimension();
61 const int dof2 = macro_reduced_size();
62 const int dof1 = reduced_size();
63
64 Eigen::VectorXd reduced(dof1 + dof2);
65 reduced.head(dof1) = NLProblem::full_to_reduced(extended.head(extended.size() - dim * dim));
66 reduced.tail(dof2) = macro_full_to_reduced(extended.tail(dim * dim));
67
68 return reduced;
69 }
70
71 Eigen::VectorXd NLHomoProblem::reduced_to_extended(const Eigen::VectorXd &reduced, bool homogeneous) const
72 {
73 const int dim = state_.mesh->dimension();
74 const int dof2 = macro_reduced_size();
75 const int dof1 = reduced_size();
76 assert(reduced.size() == dof1 + dof2);
77
78 Eigen::VectorXd full = NLProblem::reduced_to_full(reduced.head(dof1));
79 Eigen::VectorXd disp_grad = macro_reduced_to_full(reduced.tail(dof2), homogeneous);
80 Eigen::VectorXd extended(full.size() + disp_grad.size());
81 extended << full, disp_grad;
82
83 return extended;
84 }
85
86 Eigen::VectorXd NLHomoProblem::extended_to_reduced_grad(const Eigen::VectorXd &extended) const
87 {
88 const int dim = state_.mesh->dimension();
89 const int dof2 = macro_reduced_size();
90 const int dof1 = reduced_size();
91
92 Eigen::VectorXd grad(dof1 + dof2);
93 grad.head(dof1) = NLProblem::full_to_reduced_grad(extended.head(extended.size() - dim * dim));
94 grad.tail(dof2) = macro_full_to_reduced_grad(extended.tail(dim * dim));
95
96 return grad;
97 }
98
99 double NLHomoProblem::value(const TVector &x)
100 {
102
104 {
105 val += penalty_problem_->value(x);
106 }
107
108 for (auto &form : homo_forms)
109 if (form->enabled())
110 val += form->value(reduced_to_extended(x));
111
112 return val;
113 }
114 void NLHomoProblem::gradient(const TVector &x, TVector &gradv)
115 {
117
118 if (full_size() != current_size())
119 {
120 gradv = full_to_reduced_grad(gradv);
121 }
122 else if (penalty_problem_)
123 {
124 TVector tmp;
125 penalty_problem_->gradient(x, tmp);
126 gradv += tmp;
127 }
128
129 for (auto &form : homo_forms)
130 if (form->enabled())
131 {
132 Eigen::VectorXd grad_extended;
133 form->first_derivative(reduced_to_extended(x), grad_extended);
134 gradv += extended_to_reduced_grad(grad_extended);
135 }
136 }
137 void NLHomoProblem::extended_hessian_to_reduced_hessian(const THessian &extended, THessian &reduced) const
138 {
139 const int dim = state_.mesh->dimension();
140 const int dof2 = macro_reduced_size();
141 const int dof1 = reduced_size();
142
143 Eigen::MatrixXd A12, A22;
144 {
145 Eigen::MatrixXd tmp = Eigen::MatrixXd(extended.rightCols(dim * dim));
146 A12 = macro_full_to_reduced_grad(tmp.topRows(tmp.rows() - dim * dim).transpose()).transpose();
147 A22 = macro_full_to_reduced_grad(macro_full_to_reduced_grad(tmp.bottomRows(dim * dim)).transpose());
148 }
149
150 std::vector<Eigen::Triplet<double>> entries;
151 entries.reserve(extended.nonZeros() + A12.size() * 2 + A22.size());
152
153 for (int k = 0; k < full_size_; ++k)
154 for (StiffnessMatrix::InnerIterator it(extended, k); it; ++it)
155 if (it.row() < full_size_ && it.col() < full_size_)
156 entries.emplace_back(it.row(), it.col(), it.value());
157
158 for (int i = 0; i < A12.rows(); i++)
159 {
160 for (int j = 0; j < A12.cols(); j++)
161 {
162 entries.emplace_back(i, full_size_ + j, A12(i, j));
163 entries.emplace_back(full_size_ + j, i, A12(i, j));
164 }
165 }
166
167 for (int i = 0; i < A22.rows(); i++)
168 for (int j = 0; j < A22.cols(); j++)
169 entries.emplace_back(i + full_size_, j + full_size_, A22(i, j));
170
171 // Eigen::VectorXd tmp = macro_full_to_reduced_grad(Eigen::VectorXd::Ones(dim*dim));
172 // for (int i = 0; i < tmp.size(); i++)
173 // entries.emplace_back(i + full_size_, i + full_size_, 1 - tmp(i));
174
175 reduced.resize(0, 0);
176 reduced.resize(full_size_ + dof2, full_size_ + dof2);
177 reduced.setFromTriplets(entries.begin(), entries.end());
178
179 {
180 assert(dof1 == Q2_.cols());
181 StiffnessMatrix Q2_extended(Q2_.rows() + dof2, Q2_.cols() + dof2);
182
183 {
184 entries.clear();
185 for (int k = 0; k < Q2_.cols(); ++k)
186 for (StiffnessMatrix::InnerIterator it(Q2_, k); it; ++it)
187 {
188 entries.emplace_back(it.row(), it.col(), it.value());
189 }
190
191 for (int k = 0; k < dof2; k++)
192 {
193 entries.emplace_back(Q2_.rows() + k, dof1 + k, 1);
194 }
195
196 Q2_extended.setFromTriplets(entries.begin(), entries.end());
197 }
198
199 reduced = Q2_extended.transpose() * reduced * Q2_extended;
200 // remove numerical zeros
201 reduced.prune([](const Eigen::Index &row, const Eigen::Index &col, const Scalar &value) {
202 return std::abs(value) > 1e-10;
203 });
204 }
205 }
206 Eigen::MatrixXd NLHomoProblem::reduced_to_disp_grad(const TVector &reduced, bool homogeneous) const
207 {
208 const int dim = state_.mesh->dimension();
209 const int dof2 = macro_reduced_size();
210 const int dof1 = reduced_size();
211
212 return utils::unflatten(macro_reduced_to_full(reduced.tail(dof2), homogeneous), dim);
213 }
214 void NLHomoProblem::hessian(const TVector &x, THessian &hessian)
215 {
217
219
220 for (auto &form : homo_forms)
221 if (form->enabled())
222 {
223 THessian hess_extended;
224 form->second_derivative(reduced_to_extended(x), hess_extended);
225
226 THessian hess;
227 extended_hessian_to_reduced_hessian(hess_extended, hess);
228 hessian += hess;
229 }
230 }
231
232 void NLHomoProblem::set_fixed_entry(const Eigen::VectorXi &fixed_entry)
233 {
234 const int dim = state_.mesh->dimension();
235
236 Eigen::VectorXd fixed_mask;
237 fixed_mask.setZero(dim * dim);
238 fixed_mask(fixed_entry.array()).setOnes();
239 fixed_mask = (macro_full_to_mid_ * fixed_mask).eval();
240
241 fixed_mask_.setZero(fixed_mask.size());
242 for (int i = 0; i < fixed_mask.size(); i++)
243 if (abs(fixed_mask(i)) > 1e-8)
244 fixed_mask_(i) = true;
245
246 const int new_reduced_size = fixed_mask_.size() - fixed_mask_.sum();
247 macro_mid_to_reduced_.setZero(new_reduced_size, fixed_mask_.size());
248 for (int i = 0, j = 0; i < fixed_mask_.size(); i++)
249 if (!fixed_mask_(i))
250 macro_mid_to_reduced_(j++, i) = 1;
251 }
252
254 {
255 const int dim = state_.mesh->dimension();
256 const int dof2 = macro_reduced_size();
257 const int dof1 = reduced_size();
258 const int full_size = hessian.rows();
259
260 Eigen::MatrixXd tmp = constraint_grad();
261 Eigen::MatrixXd A12 = hessian * tmp.transpose();
262 Eigen::MatrixXd A22 = tmp * A12;
263
264 std::vector<Eigen::Triplet<double>> entries;
265 entries.reserve(hessian.nonZeros() + A12.size() * 2 + A22.size());
266
267 for (int k = 0; k < hessian.outerSize(); ++k)
268 for (StiffnessMatrix::InnerIterator it(hessian, k); it; ++it)
269 entries.emplace_back(it.row(), it.col(), it.value());
270
271 for (int i = 0; i < A12.rows(); i++)
272 for (int j = 0; j < A12.cols(); j++)
273 {
274 entries.emplace_back(i, full_size + j, A12(i, j));
275 entries.emplace_back(full_size + j, i, A12(i, j));
276 }
277
278 for (int i = 0; i < A22.rows(); i++)
279 for (int j = 0; j < A22.cols(); j++)
280 entries.emplace_back(i + full_size, j + full_size, A22(i, j));
281
282 hessian.resize(0, 0);
283 hessian.resize(full_size + dof2, full_size + dof2);
284 hessian.setFromTriplets(entries.begin(), entries.end());
285 // NLProblem::full_hessian_to_reduced_hessian(hessian);
286
287 {
288 assert(dof1 == Q2_.cols());
289 StiffnessMatrix Q2_extended(Q2_.rows() + dof2, Q2_.cols() + dof2);
290
291 {
292 entries.clear();
293 for (int k = 0; k < Q2_.cols(); ++k)
294 for (StiffnessMatrix::InnerIterator it(Q2_, k); it; ++it)
295 {
296 entries.emplace_back(it.row(), it.col(), it.value());
297 }
298
299 for (int k = 0; k < dof2; k++)
300 {
301 entries.emplace_back(Q2_.rows() + k, dof1 + k, 1);
302 }
303
304 Q2_extended.setFromTriplets(entries.begin(), entries.end());
305 }
306
307 hessian = Q2_extended.transpose() * hessian * Q2_extended;
308 // remove numerical zeros
309 hessian.prune([](const Eigen::Index &row, const Eigen::Index &col, const Scalar &value) {
310 return std::abs(value) > 1e-10;
311 });
312 }
313 }
314
315 NLHomoProblem::TVector NLHomoProblem::full_to_reduced(const TVector &full) const
316 {
317 log_and_throw_error("Invalid function!");
318 return TVector();
319 }
320
321 NLHomoProblem::TVector NLHomoProblem::full_to_reduced(const TVector &full, const Eigen::MatrixXd &disp_grad) const
322 {
323 const int dim = state_.mesh->dimension();
324 const int dof2 = macro_reduced_size();
325 const int dof1 = reduced_size();
326
327 TVector reduced;
328 reduced.setZero(dof1 + dof2);
329
331 reduced.tail(dof2) = macro_full_to_reduced(utils::flatten(disp_grad));
332
333 return reduced;
334 }
335 NLHomoProblem::TVector NLHomoProblem::full_to_reduced_grad(const TVector &full) const
336 {
337 const int dim = state_.mesh->dimension();
338 const int dof2 = macro_reduced_size();
339 const int dof1 = reduced_size();
340
341 TVector reduced;
342 reduced.setZero(dof1 + dof2);
343 reduced.head(dof1) = NLProblem::full_to_reduced_grad(full);
344 reduced.tail(dof2) = constraint_grad() * full;
345
346 return reduced;
347 }
348 NLHomoProblem::TVector NLHomoProblem::full_to_reduced_diag(const TVector &full_diag) const
349 {
350 const int dof2 = macro_reduced_size();
351 const int dof1 = reduced_size();
352
353 TVector reduced_diag;
354 reduced_diag.setZero(dof1 + dof2);
355 reduced_diag.head(dof1) = NLProblem::full_to_reduced_diag(full_diag);
356 reduced_diag.tail(dof2) = constraint_grad().cwiseAbs2() * full_diag;
357
358 return reduced_diag;
359 }
360 NLHomoProblem::TVector NLHomoProblem::reduced_to_full(const TVector &reduced) const
361 {
362 const int dim = state_.mesh->dimension();
363 const int dof2 = macro_reduced_size();
364 const int dof1 = reduced_size();
365
366 Eigen::MatrixXd disp_grad = utils::unflatten(macro_reduced_to_full(reduced.tail(dof2)), dim);
368 }
369
371 {
372 return macro_mid_to_reduced_.rows();
373 }
374 NLHomoProblem::TVector NLHomoProblem::macro_full_to_reduced(const TVector &full) const
375 {
377 }
378 Eigen::MatrixXd NLHomoProblem::macro_full_to_reduced_grad(const Eigen::MatrixXd &full) const
379 {
380 return macro_mid_to_reduced_ * macro_mid_to_full_.transpose() * full;
381 }
382 NLHomoProblem::TVector NLHomoProblem::macro_reduced_to_full(const TVector &reduced, bool homogeneous) const
383 {
384 TVector mid = macro_mid_to_reduced_.transpose() * reduced;
385 const TVector fixed_values = homogeneous ? TVector::Zero(macro_full_to_mid_.rows()) : TVector(macro_full_to_mid_ * utils::flatten(macro_strain_constraint_.eval(t_)));
386 for (int i = 0; i < fixed_mask_.size(); i++)
387 if (fixed_mask_(i))
388 mid(i) = fixed_values(i);
389
390 return macro_mid_to_full_ * mid;
391 }
392
393 void NLHomoProblem::init(const TVector &x0)
394 {
395 for (auto &form : homo_forms)
396 form->init(reduced_to_extended(x0));
398 }
399
400 bool NLHomoProblem::is_step_valid(const TVector &x0, const TVector &x1)
401 {
402 bool flag = NLProblem::is_step_valid(x0.head(reduced_size()), x1.head(reduced_size()));
403 for (auto &form : homo_forms)
404 if (form->enabled())
405 flag &= form->is_step_valid(reduced_to_extended(x0), reduced_to_extended(x1));
406
407 return flag;
408 }
409 bool NLHomoProblem::is_step_collision_free(const TVector &x0, const TVector &x1)
410 {
411 bool flag = NLProblem::is_step_collision_free(x0.head(reduced_size()), x1.head(reduced_size()));
412 for (auto &form : homo_forms)
413 if (form->enabled())
414 flag &= form->is_step_collision_free(reduced_to_extended(x0), reduced_to_extended(x1));
415
416 return flag;
417 }
418 double NLHomoProblem::max_step_size(const TVector &x0, const TVector &x1)
419 {
420 double size = NLProblem::max_step_size(x0.head(reduced_size()), x1.head(reduced_size()));
421 for (auto &form : homo_forms)
422 if (form->enabled())
423 size = std::min(size, form->max_step_size(reduced_to_extended(x0), reduced_to_extended(x1)));
424 return size;
425 }
426
427 void NLHomoProblem::line_search_begin(const TVector &x0, const TVector &x1)
428 {
430 for (auto &form : homo_forms)
431 form->line_search_begin(reduced_to_extended(x0), reduced_to_extended(x1));
432 }
433 void NLHomoProblem::post_step(const polysolve::nonlinear::PostStepData &data)
434 {
435 NLProblem::post_step(polysolve::nonlinear::PostStepData(data.iter_num, data.solver_info, reduced_to_full(data.x.head(reduced_size())), reduced_to_full(data.grad.head(reduced_size()))));
436 for (auto &form : homo_forms)
437 form->post_step(polysolve::nonlinear::PostStepData(
438 data.iter_num, data.solver_info, reduced_to_extended(data.x), reduced_to_extended(data.grad)));
439 }
440
441 void NLHomoProblem::solution_changed(const TVector &new_x)
442 {
444 for (auto &form : homo_forms)
445 form->solution_changed(reduced_to_extended(new_x));
446 }
447
448 void NLHomoProblem::init_lagging(const TVector &x)
449 {
451 for (auto &form : homo_forms)
452 form->init_lagging(reduced_to_extended(x));
453 }
454 void NLHomoProblem::update_lagging(const TVector &x, const int iter_num)
455 {
456 NLProblem::update_lagging(x.head(reduced_size()), iter_num);
457 for (auto &form : homo_forms)
458 form->update_lagging(reduced_to_extended(x), iter_num);
459 }
460
461 void NLHomoProblem::update_quantities(const double t, const TVector &x)
462 {
464 for (auto &form : homo_forms)
465 form->update_quantities(t, reduced_to_extended(x));
466 }
467
468 Eigen::MatrixXd NLHomoProblem::constraint_grad() const
469 {
470 const int dim = state_.mesh->dimension();
471 Eigen::MatrixXd jac; // (dim*dim) x (dim*n_bases)
472
474
475 jac.setZero(dim * dim, full_size_);
476 for (int i = 0; i < X.rows(); i++)
477 for (int j = 0; j < dim; j++)
478 for (int k = 0; k < dim; k++)
479 jac(j * dim + k, i * dim + j) = X(i, k);
480
481 return macro_full_to_reduced_grad(jac);
482 }
483} // 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:114
int n_bases
number of bases
Definition State.hpp:214
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:594
std::shared_ptr< polyfem::mesh::MeshNodes > mesh_nodes
Mapping from input nodes to FE nodes.
Definition State.hpp:229
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 hessian(const TVector &x, THessian &hessian) override
virtual double value(const TVector &x) override
virtual void init(const TVector &x0) override
virtual void gradient(const TVector &x, TVector &gradv) 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
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, const std::shared_ptr< polysolve::linear::Solver > &solver, const double char_length, const double char_force, StiffnessMatrix lumped_mass, const int dimension)
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
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 full_to_reduced_diag(const TVector &full_diag) const override
TVector reduced_to_extended(const TVector &reduced, bool homogeneous=false) const
TVector reduced_to_full(const TVector &reduced) const
TVector extended_to_reduced(const TVector &extended) const
std::vector< std::shared_ptr< Form > > homo_forms
Eigen::MatrixXd constraint_grad() const
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 full_hessian_to_reduced_hessian(THessian &hessian) const
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:84
void line_search_begin(const TVector &x0, const TVector &x1) override
StiffnessMatrix Q2_
Q2 block of the QR decomposition of the constraints matrix.
virtual bool is_step_valid(const TVector &x0, const TVector &x1) override
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_diag(const TVector &full_diag) const
TVector full_to_reduced(const TVector &full) const
virtual void update_quantities(const double t, const TVector &x)
void init_lagging(const TVector &x) override
virtual bool is_step_collision_free(const TVector &x0, const TVector &x1) override
TVector reduced_to_full(const TVector &reduced) const
void update_lagging(const TVector &x, const int iter_num) override
virtual double max_step_size(const TVector &x0, const TVector &x1) override
void solution_changed(const TVector &new_x) override
std::shared_ptr< FullNLProblem > penalty_problem_
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:73
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24