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