PolyFEM
Loading...
Searching...
No Matches
NLProblem.cpp
Go to the documentation of this file.
1#include "NLProblem.hpp"
2
6
8
9#include <polysolve/linear/Solver.hpp>
10#ifdef POLYSOLVE_WITH_SPQR
11#include <Eigen/SPQRSupport>
12#include <SuiteSparseQR.hpp>
13#endif
14
15#ifdef POLYFEM_WITH_ARMADILLO
16#include <armadillo>
17#endif
18
19#include <igl/cat.h>
20#include <igl/Timer.h>
21
22#include <set>
23
24/*
25m \frac{\partial^2 u}{\partial t^2} = \psi = \text{div}(\sigma[u])\newline
26u^{t+1} = u(t+\Delta t)\approx u(t) + \Delta t \dot u + \frac{\Delta t^2} 2 \ddot u \newline
27= u(t) + \Delta t \dot u + \frac{\Delta t^2}{2} \psi\newline
28M u^{t+1}_h \approx M u^t_h + \Delta t M v^t_h + \frac{\Delta t^2} {2} A u^{t+1}_h \newline
29%
30M (u^{t+1}_h - (u^t_h + \Delta t v^t_h)) - \frac{\Delta t^2} {2} A u^{t+1}_h
31*/
32// mü = ψ = div(σ[u])
33// uᵗ⁺¹ = u(t + Δt) ≈ u(t) + Δtu̇ + ½Δt²ü = u(t) + Δtu̇ + ½Δt²ψ
34// Muₕᵗ⁺¹ ≈ Muₕᵗ + ΔtMvₕᵗ ½Δt²Auₕᵗ⁺¹
35// Root-finding form:
36// M(uₕᵗ⁺¹ - (uₕᵗ + Δtvₕᵗ)) - ½Δt²Auₕᵗ⁺¹ = 0
37
38namespace polyfem::solver
39{
40
41 namespace
42 {
43#ifdef POLYSOLVE_WITH_SPQR
44 void fill_cholmod(Eigen::SparseMatrix<double, Eigen::ColMajor, long> &mat, cholmod_sparse &cmat)
45 {
46 long *p = mat.outerIndexPtr();
47
48 cmat.nzmax = mat.nonZeros();
49 cmat.nrow = mat.rows();
50 cmat.ncol = mat.cols();
51 cmat.p = p;
52 cmat.i = mat.innerIndexPtr();
53 cmat.x = mat.valuePtr();
54 cmat.z = 0;
55 cmat.sorted = 1;
56 cmat.packed = 1;
57 cmat.nz = 0;
58 cmat.dtype = 0;
59 cmat.stype = -1;
60 cmat.xtype = CHOLMOD_REAL;
61 cmat.dtype = CHOLMOD_DOUBLE;
62 cmat.stype = 0;
63 cmat.itype = CHOLMOD_LONG;
64 }
65#endif
66
67#ifdef POLYFEM_WITH_ARMADILLO
68 arma::sp_mat fill_arma(const StiffnessMatrix &mat)
69 {
70 std::vector<unsigned long long> rowind_vect(mat.innerIndexPtr(), mat.innerIndexPtr() + mat.nonZeros());
71 std::vector<unsigned long long> colptr_vect(mat.outerIndexPtr(), mat.outerIndexPtr() + mat.outerSize() + 1);
72 std::vector<double> values_vect(mat.valuePtr(), mat.valuePtr() + mat.nonZeros());
73
74 arma::dvec values(values_vect.data(), values_vect.size(), false);
75 arma::uvec rowind(rowind_vect.data(), rowind_vect.size(), false);
76 arma::uvec colptr(colptr_vect.data(), colptr_vect.size(), false);
77
78 arma::sp_mat amat(rowind, colptr, values, mat.rows(), mat.cols(), false);
79
80 return amat;
81 }
82
83 StiffnessMatrix fill_eigen(const arma::sp_mat &mat)
84 {
85 // convert to eigen sparse
86 std::vector<long> outerIndexPtr(mat.row_indices, mat.row_indices + mat.n_nonzero);
87 std::vector<long> innerIndexPtr(mat.col_ptrs, mat.col_ptrs + mat.n_cols + 1);
88
89 const StiffnessMatrix out = Eigen::Map<const Eigen::SparseMatrix<double, Eigen::ColMajor, long>>(
90 mat.n_rows, mat.n_cols, mat.n_nonzero, innerIndexPtr.data(), outerIndexPtr.data(), mat.values);
91 return out;
92 }
93#endif
94
95 } // namespace
97 const int full_size,
98 const std::vector<std::shared_ptr<Form>> &forms,
99 const std::vector<std::shared_ptr<AugmentedLagrangianForm>> &penalty_forms,
100 const std::shared_ptr<polysolve::linear::Solver> &solver)
101 : FullNLProblem(forms),
102 full_size_(full_size),
103 t_(0),
104 penalty_forms_(penalty_forms),
105 solver_(solver)
106 {
109 }
110
112 const int full_size,
113 const std::shared_ptr<utils::PeriodicBoundary> &periodic_bc,
114 const double t,
115 const std::vector<std::shared_ptr<Form>> &forms,
116 const std::vector<std::shared_ptr<AugmentedLagrangianForm>> &penalty_forms,
117 const std::shared_ptr<polysolve::linear::Solver> &solver,
118 const double char_length,
119 const double char_force,
120 StiffnessMatrix lumped_mass,
121 const int dimension)
122 : FullNLProblem(forms),
123 full_size_(full_size),
124 t_(t),
125 penalty_forms_(penalty_forms),
126 solver_(solver),
127 L(char_length),
128 F0(char_force),
129 lumped_mass_(lumped_mass.diagonal().asDiagonal()),
130 dim(dimension)
131 {
134
135 double total_lumped_mass = 0;
136 int num_nonzero_mass_entries = 0;
137 for (int i = 0; i < lumped_mass_.diagonal().size(); i++)
138 {
139 if (lumped_mass_.diagonal()[i] > 0)
140 {
141 total_lumped_mass += lumped_mass_.diagonal()[i];
142 num_nonzero_mass_entries++;
143 }
144 }
145 const double avg_lumped_mass = total_lumped_mass / num_nonzero_mass_entries;
146 for (int i = 0; i < lumped_mass_.diagonal().size(); i++)
147 {
148 if (lumped_mass_.diagonal()[i] == 0)
149 {
150 lumped_mass_.diagonal()[i] = avg_lumped_mass;
151 }
152 }
153 }
154
156 {
157 return 1;
158
159 // double total_weight = 0;
160 // for (const auto &f : forms_)
161 // total_weight += f->weight();
162 // if (full_size() == current_size())
163 // {
164 // for (const auto &f : penalty_forms_)
165 // total_weight += f->weight() * f->lagrangian_weight();
166 // }
167
168 // logger().debug("Normalizing forms with scale: {}", total_weight);
169
170 // for (auto &f : forms_)
171 // f->set_scale(total_weight);
172 // for (auto &f : penalty_forms_)
173 // f->set_scale(total_weight);
174
175 // return total_weight;
176 }
177
178 double NLProblem::grad_norm_rescaling(const polysolve::nonlinear::NormType norm_type) const
179 {
180 switch (norm_type)
181 {
182 case polysolve::nonlinear::NormType::EUCLIDEAN:
183 return 1;
184 case polysolve::nonlinear::NormType::L2:
185 return F0 * (dim == 2 ? L : std::pow(L, 1.5));
186 case polysolve::nonlinear::NormType::Linf:
187 return F0;
188 }
189 log_and_throw_error("Unrecognized norm type!");
190 }
191
192 double NLProblem::step_norm_rescaling(const polysolve::nonlinear::NormType norm_type) const
193 {
194 switch (norm_type)
195 {
196 case polysolve::nonlinear::NormType::EUCLIDEAN:
197 return 1;
198 case polysolve::nonlinear::NormType::L2:
199 return dim == 2 ? L * L : std::pow(L, 2.5);
200 case polysolve::nonlinear::NormType::Linf:
201 return L;
202 }
203 log_and_throw_error("Unrecognized norm type!");
204 }
205
206 double NLProblem::energy_norm_rescaling(const polysolve::nonlinear::NormType norm_type) const
207 {
208 const double density_scale = dim == 2 ? L * L : L * L * L;
209 return F0 * density_scale * L;
210 }
211
212 double NLProblem::grad_norm(const TVector &grad, const polysolve::nonlinear::NormType norm_type) const
213 {
214 switch (norm_type)
215 {
216 case polysolve::nonlinear::NormType::EUCLIDEAN:
217 return grad.norm();
218 case polysolve::nonlinear::NormType::L2:
219 return sqrt(grad.transpose() * current_lumped_mass().inverse() * grad);
220 case polysolve::nonlinear::NormType::Linf:
221 return (current_lumped_mass().inverse() * grad).cwiseAbs().maxCoeff();
222 }
223 log_and_throw_error("Unrecognized norm type!");
224 }
225
226 double NLProblem::step_norm(const TVector &x, const polysolve::nonlinear::NormType norm_type) const
227 {
228 switch (norm_type)
229 {
230 case polysolve::nonlinear::NormType::EUCLIDEAN:
231 return x.norm();
232 case polysolve::nonlinear::NormType::L2:
233 return sqrt(x.transpose() * current_lumped_mass() * x);
234 case polysolve::nonlinear::NormType::Linf:
235 return x.cwiseAbs().maxCoeff();
236 }
237 log_and_throw_error("Unrecognized norm type!");
238 }
239
241 {
242 if (penalty_forms_.empty())
243 {
245 return;
246 }
247 igl::Timer timer;
248
249 if (penalty_forms_.size() == 1 && penalty_forms_.front()->has_projection())
250 {
251 Q2_ = penalty_forms_.front()->constraint_projection_matrix();
252 Q2t_ = Q2_.transpose();
253
254 reduced_size_ = Q2_.cols();
255 if (reduced_size_ == 0)
256 return;
258
259 timer.start();
260 StiffnessMatrix Q2tQ2 = Q2t_ * Q2_;
261 solver_->analyze_pattern(Q2tQ2, Q2tQ2.rows());
262 solver_->factorize(Q2tQ2);
263 timer.stop();
264 logger().debug("Factorization and computation of Q2tQ2 took: {}", timer.getElapsedTime());
265
266 std::vector<std::shared_ptr<Form>> tmp;
267 tmp.insert(tmp.end(), penalty_forms_.begin(), penalty_forms_.end());
268 penalty_problem_ = std::make_shared<FullNLProblem>(tmp);
269
271
272 return;
273 }
274
275 std::vector<Eigen::Triplet<double>> Ae;
276 int index = 0;
277 for (const auto &f : penalty_forms_)
278 {
279 const auto &tmp = f->constraint_matrix();
280 for (int i = 0; i < tmp.outerSize(); i++)
281 {
282 for (typename StiffnessMatrix::InnerIterator it(tmp, i); it; ++it)
283 {
284 Ae.emplace_back(index + it.row(), it.col(), it.value());
285 }
286 }
287 index += tmp.rows();
288 }
289 StiffnessMatrix A(index, full_size_);
290 A.setFromTriplets(Ae.begin(), Ae.end());
291 A.makeCompressed();
292
293 int constraint_size = A.rows();
294 num_penalty_constraints_ = A.rows();
295 Eigen::SparseMatrix<double, Eigen::ColMajor, long> At = A.transpose();
296 At.makeCompressed();
297
298 logger().debug("Constraint size: {} x {}", A.rows(), A.cols());
299
300#ifdef POLYSOLVE_WITH_SPQR
301 timer.start();
302
303 cholmod_common cc;
304 cholmod_l_start(&cc); // start CHOLMOD
305
306 // const int ordering = 0; // all, except 3:given treated as 0:fixed
307 const int ordering = SPQR_ORDERING_DEFAULT; // all, except 3:given treated as 0:fixed
308 const double tol = SPQR_DEFAULT_TOL;
309 SuiteSparse_long econ = At.rows();
310 SuiteSparse_long *E; // permutation of 0:n-1, NULL if identity
311
312 cholmod_sparse Ac, *Qc, *Rc;
313
314 fill_cholmod(At, Ac);
315
316 const auto rank = SuiteSparseQR<double>(ordering, tol, econ, &Ac,
317 // outputs
318 &Qc, &Rc, &E, &cc);
319
320 if (!Rc)
321 log_and_throw_error("Failed to factorize constraints matrix");
322
323 const auto n = Rc->ncol;
324 P_.resize(n);
325 if (E)
326 {
327 for (long j = 0; j < n; j++)
328 P_.indices()(j) = E[j];
329 }
330 else
331 P_.setIdentity();
332
333 if (Qc->stype != 0 || Qc->sorted != 1 || Qc->packed != 1 || Rc->stype != 0 || Rc->sorted != 1 || Rc->packed != 1)
334 log_and_throw_error("Q and R must be unsymmetric sorted and packed");
335
336 const StiffnessMatrix Q = Eigen::Map<Eigen::SparseMatrix<double, Eigen::ColMajor, long>>(
337 Qc->nrow, Qc->ncol, Qc->nzmax,
338 static_cast<long *>(Qc->p), static_cast<long *>(Qc->i), static_cast<double *>(Qc->x));
339
340 const StiffnessMatrix R = Eigen::Map<Eigen::SparseMatrix<double, Eigen::ColMajor, long>>(
341 Rc->nrow, Rc->ncol, Rc->nzmax,
342 static_cast<long *>(Rc->p), static_cast<long *>(Rc->i), static_cast<double *>(Rc->x));
343
344 cholmod_l_free_sparse(&Qc, &cc);
345 cholmod_l_free_sparse(&Rc, &cc);
346 std::free(E);
347 cholmod_l_finish(&cc);
348
349 timer.stop();
350 logger().debug("QR took: {}", timer.getElapsedTime());
351#else
352 timer.start();
353
354 // Eigen::SparseQR<StiffnessMatrix, Eigen::NaturalOrdering<int>> QR(At);
355 Eigen::SparseQR<StiffnessMatrix, Eigen::COLAMDOrdering<int>> QR(At);
356
357 timer.stop();
358 logger().debug("QR took: {}", timer.getElapsedTime());
359
360 if (QR.info() != Eigen::Success)
361 log_and_throw_error("Failed to factorize constraints matrix");
362
363 timer.start();
365 Q = QR.matrixQ();
366 timer.stop();
367 logger().debug("Computation of Q took: {}", timer.getElapsedTime());
368
369 const Eigen::SparseMatrix<double, Eigen::RowMajor> R = QR.matrixR();
370
371 P_ = QR.colsPermutation();
372#endif
373
374 for (; constraint_size >= 0; --constraint_size)
375 {
376 const StiffnessMatrix tmp = R.row(constraint_size);
377 if (tmp.nonZeros() != 0)
378 {
379 constraint_size++;
380 break;
381 }
382 }
383 if (constraint_size != num_penalty_constraints_)
384 logger().warn("Matrix A is not full rank, constraint size: {} instead of {}", constraint_size, num_penalty_constraints_);
385
386 reduced_size_ = full_size_ - constraint_size;
387
388 timer.start();
389
390 Q1_ = Q.leftCols(constraint_size);
391 assert(Q1_.rows() == full_size_);
392 assert(Q1_.cols() == constraint_size);
393
394 Q2_ = Q.rightCols(reduced_size_);
395 Q2t_ = Q2_.transpose();
396
397 assert(Q2_.rows() == full_size_);
398 assert(Q2_.cols() == reduced_size_);
399
400 R1_ = R.topRows(constraint_size);
401 assert(R1_.rows() == constraint_size);
402 assert(R1_.cols() == num_penalty_constraints_);
403
404 assert((Q1_.transpose() * Q2_).norm() < 1e-10);
405
406 timer.stop();
407 logger().debug("Getting Q1 Q2, R1 took: {}", timer.getElapsedTime());
408
409 timer.start();
410
411 // arma::sp_mat q2a = fill_arma(Q2_);
412 // arma::sp_mat q2tq2 = q2a.t() * q2a;
413 // const StiffnessMatrix Q2tQ2 = fill_eigen(q2tq2);
414 StiffnessMatrix Q2tQ2 = Q2t_ * Q2_;
415 timer.stop();
416 logger().debug("Getting Q2'*Q2, took: {}", timer.getElapsedTime());
417
418 timer.start();
419 solver_->analyze_pattern(Q2tQ2, Q2tQ2.rows());
420 solver_->factorize(Q2tQ2);
421 timer.stop();
422 logger().debug("Factorization of Q2'*Q2 took: {}", timer.getElapsedTime());
423
424#ifndef NDEBUG
425 StiffnessMatrix test = R.bottomRows(reduced_size_);
426 assert(test.nonZeros() == 0);
427
428 StiffnessMatrix test1 = R1_.row(R1_.rows() - 1);
429 assert(test1.nonZeros() != 0);
430#endif
431
432 // assert((Q1_ * R1_ - At * P_).norm() < 1e-10);
433
434 std::vector<std::shared_ptr<Form>> tmp;
435 tmp.insert(tmp.end(), penalty_forms_.begin(), penalty_forms_.end());
436 penalty_problem_ = std::make_shared<FullNLProblem>(tmp);
437
439 }
440
442 {
443 if (penalty_forms_.size() == 1 && penalty_forms_.front()->has_projection())
444 {
445 Q1R1iTb_ = penalty_forms_.front()->constraint_projection_vector();
446 return;
447 }
448
449 igl::Timer timer;
450 timer.start();
451 // x = Q1 * R1^(-T) * P^T b + Q2 * y
452 int index = 0;
453 TVector constraint_values(num_penalty_constraints_);
454 for (const auto &f : penalty_forms_)
455 {
456 constraint_values.segment(index, f->constraint_value().rows()) = f->constraint_value();
457 index += f->constraint_value().rows();
458 }
459 constraint_values = P_.transpose() * constraint_values;
460
461 Eigen::VectorXd sol;
462
463 if (R1_.rows() == R1_.cols())
464 {
465 sol = R1_.transpose().triangularView<Eigen::Lower>().solve(constraint_values);
466 }
467 else
468 {
469
470#ifdef POLYSOLVE_WITH_SPQR
471 Eigen::SparseMatrix<double, Eigen::ColMajor, long> R1t = R1_.transpose();
472 cholmod_common cc;
473 cholmod_l_start(&cc); // start CHOLMOD
474 cholmod_sparse R1tc;
475 fill_cholmod(R1t, R1tc);
476
477 cholmod_dense b;
478 b.nrow = constraint_values.size();
479 b.ncol = 1;
480 b.nzmax = constraint_values.size();
481 b.d = constraint_values.size();
482 b.x = constraint_values.data();
483 b.z = 0;
484 b.xtype = CHOLMOD_REAL;
485 b.dtype = 0;
486
487 const int ordering = SPQR_ORDERING_DEFAULT; // all, except 3:given treated as 0:fixed
488 const double tol = SPQR_DEFAULT_TOL;
489
490 cholmod_dense *solc = SuiteSparseQR<double>(ordering, tol, &R1tc, &b, &cc);
491
492 sol = Eigen::Map<Eigen::VectorXd>(static_cast<double *>(solc->x), solc->nrow);
493
494 cholmod_l_free_dense(&solc, &cc);
495 cholmod_l_finish(&cc);
496#else
497 Eigen::SparseQR<StiffnessMatrix, Eigen::COLAMDOrdering<int>> solver;
498 solver.compute(R1_.transpose());
499 if (solver.info() != Eigen::Success)
500 {
501 log_and_throw_error("Failed to factorize R1^T");
502 }
503 sol = solver.solve(constraint_values);
504#endif
505 }
506
507 assert((R1_.transpose() * sol - constraint_values).norm() < 1e-10);
508
509 Q1R1iTb_ = Q1_ * sol;
510
511 timer.stop();
512 logger().debug("Computing Q1R1iTb took: {}", timer.getElapsedTime());
513 }
514
515 void NLProblem::init_lagging(const TVector &x)
516 {
518
520 penalty_problem_->init_lagging(x);
521 }
522
523 void NLProblem::update_lagging(const TVector &x, const int iter_num)
524 {
526
528 penalty_problem_->update_lagging(x, iter_num);
529 }
530
531 void NLProblem::update_quantities(const double t, const TVector &x)
532 {
533 t_ = t;
534 const TVector full = reduced_to_full(x);
535 assert(full.size() == full_size_);
536 for (auto &f : forms_)
537 f->update_quantities(t, full);
538
539 for (auto &f : penalty_forms_)
540 f->update_quantities(t, x);
541
543 }
544
545 void NLProblem::line_search_begin(const TVector &x0, const TVector &x1)
546 {
548
550 penalty_problem_->line_search_begin(x0, x1);
551 }
552
553 double NLProblem::max_step_size(const TVector &x0, const TVector &x1)
554 {
556
558 max_step = std::min(max_step, penalty_problem_->max_step_size(x0, x1));
559
560 return max_step;
561 }
562
563 bool NLProblem::is_step_valid(const TVector &x0, const TVector &x1)
564 {
566
567 if (penalty_problem_ && valid && full_size() == current_size())
568 valid = penalty_problem_->is_step_valid(x0, x1);
569
570 return valid;
571 }
572
573 bool NLProblem::is_step_collision_free(const TVector &x0, const TVector &x1)
574 {
576
577 if (penalty_problem_ && free && full_size() == current_size())
578 free = penalty_problem_->is_step_collision_free(x0, x1);
579
580 return free;
581 }
582
583 double NLProblem::value(const TVector &x)
584 {
585 // TODO: removed fearure const bool only_elastic
587
589 {
590 res += penalty_problem_->value(x);
591 }
592
593 return res;
594 }
595
596 void NLProblem::gradient(const TVector &x, TVector &grad)
597 {
599
600 if (full_size() != current_size())
601 {
602 if (penalty_forms_.size() == 1 && penalty_forms_.front()->can_project())
603 penalty_forms_.front()->project_gradient(grad);
604 else
605 grad = Q2t_ * grad;
606 }
607 else if (penalty_problem_)
608 {
609 TVector tmp;
610 penalty_problem_->gradient(x, tmp);
611 grad += tmp;
612 }
613 }
614
615 void NLProblem::hessian(const TVector &x, THessian &hessian)
616 {
618
619 if (full_size() != current_size())
620 {
622 }
623 else if (penalty_problem_)
624 {
625 THessian tmp;
626 penalty_problem_->hessian(x, tmp);
627 hessian += tmp;
628 }
629 }
630
631 void NLProblem::solution_changed(const TVector &newX)
632 {
634
636 penalty_problem_->solution_changed(newX);
637 }
638
639 void NLProblem::post_step(const polysolve::nonlinear::PostStepData &data)
640 {
641 FullNLProblem::post_step(polysolve::nonlinear::PostStepData(data.iter_num, data.solver_info, reduced_to_full(data.x), reduced_to_full(data.grad)));
642
644 penalty_problem_->post_step(data);
645
646 // TODO: add me back
647 static int nsolves = 0;
648 if (data.iter_num == 0)
649 nsolves++;
650 // if (state && state->args["output"]["advanced"]["save_nl_solve_sequence"])
651 // {
652 // const Eigen::MatrixXd displacements = utils::unflatten(reduced_to_full(data.x), state->mesh->dimension());
653 // io::OBJWriter::write(
654 // state->resolve_output_path(fmt::format("nonlinear_solve{:04d}_iter{:04d}.obj", nsolves, data.iter_num)),
655 // state->collision_mesh.displace_vertices(displacements),
656 // state->collision_mesh.edges(), state->collision_mesh.faces());
657 // }
658 }
659
660 NLProblem::TVector NLProblem::full_to_reduced(const TVector &full) const
661 {
662 // Reduced is already at the full size
663 if (full_size() == current_size() || full.size() == current_size())
664 {
665 return full;
666 }
667
668 TVector reduced(reduced_size());
669 const TVector k = full - Q1R1iTb_;
670 const TVector rhs = Q2t_ * k;
671 solver_->solve(rhs, reduced);
672
673#ifndef NDEBUG
674 StiffnessMatrix Q2tQ2 = Q2t_ * Q2_;
675 // std::cout << "err " << (Q2tQ2 * reduced - rhs).norm() << std::endl;
676 assert((Q2tQ2 * reduced - rhs).norm() < 1e-8);
677#endif
678
679 return reduced;
680 }
681
682 NLProblem::TVector NLProblem::full_to_reduced_grad(const TVector &full) const
683 {
684 TVector grad = full;
685 if (penalty_forms_.size() == 1 && penalty_forms_.front()->can_project())
686 penalty_forms_.front()->project_gradient(grad);
687 else
688 grad = Q2t_ * grad;
689
690 return grad;
691 }
692
693 NLProblem::TVector NLProblem::full_to_reduced_diag(const TVector &full_diag) const
694 {
695 if (full_size() == current_size() || full_diag.size() == current_size())
696 {
697 return full_diag;
698 }
699
700 TVector diag = full_diag;
701 if (penalty_forms_.size() == 1 && penalty_forms_.front()->can_project())
702 penalty_forms_.front()->project_diag(diag);
703 else
704 {
705 Eigen::SparseMatrix<double> reduced_mat = Q2t_ * diag.asDiagonal() * Q2_;
706 diag = reduced_mat.diagonal();
707 }
708
709 return diag;
710 }
711
712 NLProblem::TVector NLProblem::reduced_to_full(const TVector &reduced) const
713 {
714 // Full is already at the reduced size
715 if (full_size() == current_size() || full_size() == reduced.size())
716 {
717 return reduced;
718 }
719
720 // x = Q1 * R1^(-T) * P^T b + Q2 * y
721
722 const TVector full = Q1R1iTb_ + Q2_ * reduced;
723
724#ifndef NDEBUG
725 for (const auto &f : penalty_forms_)
726 {
727 // std::cout << f->compute_error(full) << std::endl;
728 assert(f->compute_error(full) < 1e-8);
729 }
730#endif
731
732 return full;
733 }
734
736 {
737 if (penalty_forms_.size() == 1 && penalty_forms_.front()->can_project())
738 penalty_forms_.front()->project_hessian(hessian);
739 else
740 {
741 // arma::sp_mat q2a = fill_arma(Q2_);
742 // arma::sp_mat ha = fill_arma(hessian);
743 // arma::sp_mat q2thq2 = q2a.t() * ha * q2a;
744 // hessian = fill_eigen(q2thq2);
745 hessian = Q2t_ * hessian * Q2_;
746 // remove numerical zeros
747 hessian.prune([](const Eigen::Index &row, const Eigen::Index &col, const Scalar &value) {
748 return std::abs(value) > 1e-10;
749 });
750 }
751 }
752} // namespace polyfem::solver
int x
virtual double max_step_size(const TVector &x0, const TVector &x1) override
virtual void init_lagging(const TVector &x)
virtual void hessian(const TVector &x, THessian &hessian) override
virtual bool is_step_collision_free(const TVector &x0, const TVector &x1)
std::vector< std::shared_ptr< Form > > forms_
virtual void post_step(const polysolve::nonlinear::PostStepData &data) override
virtual void update_lagging(const TVector &x, const int iter_num)
virtual bool is_step_valid(const TVector &x0, const TVector &x1) override
virtual double value(const TVector &x) override
virtual void solution_changed(const TVector &new_x) override
virtual void gradient(const TVector &x, TVector &gradv) override
virtual void line_search_begin(const TVector &x0, const TVector &x1) override
const int full_size_
Size of the full problem.
Definition NLProblem.hpp:84
virtual double grad_norm_rescaling(const polysolve::nonlinear::NormType norm_type) const override
StiffnessMatrix R1_
R1 block of the QR decomposition of the constraints matrix.
void line_search_begin(const TVector &x0, const TVector &x1) override
Eigen::PermutationMatrix< Eigen::Dynamic, Eigen::Dynamic > P_
Permutation matrix of the QR decomposition of the constraints matrix.
virtual double energy_norm_rescaling(const polysolve::nonlinear::NormType norm_type) const override
double normalize_forms() 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
int reduced_size_
Size of the reduced problem.
Definition NLProblem.hpp:85
virtual double step_norm_rescaling(const polysolve::nonlinear::NormType norm_type) const override
virtual void post_step(const polysolve::nonlinear::PostStepData &data) override
NLProblem(const int full_size, const std::vector< std::shared_ptr< Form > > &forms, const std::vector< std::shared_ptr< AugmentedLagrangianForm > > &penalty_forms, const std::shared_ptr< polysolve::linear::Solver > &solver)
Definition NLProblem.cpp:96
virtual TVector full_to_reduced_grad(const TVector &full) const
virtual TVector full_to_reduced_diag(const TVector &full_diag) const
virtual double step_norm(const TVector &x, const polysolve::nonlinear::NormType norm_type) const override
StiffnessMatrix Q2t_
Q2 transpose.
TVector full_to_reduced(const TVector &full) const
virtual void update_quantities(const double t, const TVector &x)
virtual void gradient(const TVector &x, TVector &gradv) override
void init_lagging(const TVector &x) override
Eigen::DiagonalMatrix< double, Eigen::Dynamic > lumped_mass_
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
std::vector< std::shared_ptr< AugmentedLagrangianForm > > penalty_forms_
Eigen::DiagonalMatrix< double, Eigen::Dynamic > current_lumped_mass() const
Definition NLProblem.hpp:98
virtual double value(const TVector &x) override
StiffnessMatrix Q1_
Q1 block of the QR decomposition of the constraints matrix.
virtual double max_step_size(const TVector &x0, const TVector &x1) override
virtual void hessian(const TVector &x, THessian &hessian) override
std::shared_ptr< polysolve::linear::Solver > solver_
void solution_changed(const TVector &new_x) override
virtual double grad_norm(const TVector &grad, const polysolve::nonlinear::NormType norm_type) const override
std::shared_ptr< FullNLProblem > penalty_problem_
TVector Q1R1iTb_
Q1_ * (R1_.transpose().triangularView<Eigen::Upper>().solve(constraint_values_))
void full_hessian_to_reduced_hessian(StiffnessMatrix &hessian) const
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > inverse(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > &mat)
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:44
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:73
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24