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 : FullNLProblem(forms),
119 full_size_(full_size),
120 t_(t),
121 penalty_forms_(penalty_forms),
122 solver_(solver)
123 {
126 }
127
129 {
130 return 1;
131
132 // double total_weight = 0;
133 // for (const auto &f : forms_)
134 // total_weight += f->weight();
135 // if (full_size() == current_size())
136 // {
137 // for (const auto &f : penalty_forms_)
138 // total_weight += f->weight() * f->lagrangian_weight();
139 // }
140
141 // logger().debug("Normalizing forms with scale: {}", total_weight);
142
143 // for (auto &f : forms_)
144 // f->set_scale(total_weight);
145 // for (auto &f : penalty_forms_)
146 // f->set_scale(total_weight);
147
148 // return total_weight;
149 }
150
152 {
153 if (penalty_forms_.empty())
154 {
156 return;
157 }
158 igl::Timer timer;
159
160 if (penalty_forms_.size() == 1 && penalty_forms_.front()->has_projection())
161 {
162 Q2_ = penalty_forms_.front()->constraint_projection_matrix();
163 Q2t_ = Q2_.transpose();
164
165 reduced_size_ = Q2_.cols();
166 if (reduced_size_ == 0)
167 return;
169
170 timer.start();
171 StiffnessMatrix Q2tQ2 = Q2t_ * Q2_;
172 solver_->analyze_pattern(Q2tQ2, Q2tQ2.rows());
173 solver_->factorize(Q2tQ2);
174 timer.stop();
175 logger().debug("Factorization and computation of Q2tQ2 took: {}", timer.getElapsedTime());
176
177 std::vector<std::shared_ptr<Form>> tmp;
178 tmp.insert(tmp.end(), penalty_forms_.begin(), penalty_forms_.end());
179 penalty_problem_ = std::make_shared<FullNLProblem>(tmp);
180
182
183 return;
184 }
185
186 std::vector<Eigen::Triplet<double>> Ae;
187 int index = 0;
188 for (const auto &f : penalty_forms_)
189 {
190 const auto &tmp = f->constraint_matrix();
191 for (int i = 0; i < tmp.outerSize(); i++)
192 {
193 for (typename StiffnessMatrix::InnerIterator it(tmp, i); it; ++it)
194 {
195 Ae.emplace_back(index + it.row(), it.col(), it.value());
196 }
197 }
198 index += tmp.rows();
199 }
200 StiffnessMatrix A(index, full_size_);
201 A.setFromTriplets(Ae.begin(), Ae.end());
202 A.makeCompressed();
203
204 int constraint_size = A.rows();
205 num_penalty_constraints_ = A.rows();
206 Eigen::SparseMatrix<double, Eigen::ColMajor, long> At = A.transpose();
207 At.makeCompressed();
208
209 logger().debug("Constraint size: {} x {}", A.rows(), A.cols());
210
211#ifdef POLYSOLVE_WITH_SPQR
212 timer.start();
213
214 cholmod_common cc;
215 cholmod_l_start(&cc); // start CHOLMOD
216
217 // const int ordering = 0; // all, except 3:given treated as 0:fixed
218 const int ordering = SPQR_ORDERING_DEFAULT; // all, except 3:given treated as 0:fixed
219 const double tol = SPQR_DEFAULT_TOL;
220 SuiteSparse_long econ = At.rows();
221 SuiteSparse_long *E; // permutation of 0:n-1, NULL if identity
222
223 cholmod_sparse Ac, *Qc, *Rc;
224
225 fill_cholmod(At, Ac);
226
227 const auto rank = SuiteSparseQR<double>(ordering, tol, econ, &Ac,
228 // outputs
229 &Qc, &Rc, &E, &cc);
230
231 if (!Rc)
232 log_and_throw_error("Failed to factorize constraints matrix");
233
234 const auto n = Rc->ncol;
235 P_.resize(n);
236 if (E)
237 {
238 for (long j = 0; j < n; j++)
239 P_.indices()(j) = E[j];
240 }
241 else
242 P_.setIdentity();
243
244 if (Qc->stype != 0 || Qc->sorted != 1 || Qc->packed != 1 || Rc->stype != 0 || Rc->sorted != 1 || Rc->packed != 1)
245 log_and_throw_error("Q and R must be unsymmetric sorted and packed");
246
247 const StiffnessMatrix Q = Eigen::Map<Eigen::SparseMatrix<double, Eigen::ColMajor, long>>(
248 Qc->nrow, Qc->ncol, Qc->nzmax,
249 static_cast<long *>(Qc->p), static_cast<long *>(Qc->i), static_cast<double *>(Qc->x));
250
251 const StiffnessMatrix R = Eigen::Map<Eigen::SparseMatrix<double, Eigen::ColMajor, long>>(
252 Rc->nrow, Rc->ncol, Rc->nzmax,
253 static_cast<long *>(Rc->p), static_cast<long *>(Rc->i), static_cast<double *>(Rc->x));
254
255 cholmod_l_free_sparse(&Qc, &cc);
256 cholmod_l_free_sparse(&Rc, &cc);
257 std::free(E);
258 cholmod_l_finish(&cc);
259
260 timer.stop();
261 logger().debug("QR took: {}", timer.getElapsedTime());
262#else
263 timer.start();
264
265 // Eigen::SparseQR<StiffnessMatrix, Eigen::NaturalOrdering<int>> QR(At);
266 Eigen::SparseQR<StiffnessMatrix, Eigen::COLAMDOrdering<int>> QR(At);
267
268 timer.stop();
269 logger().debug("QR took: {}", timer.getElapsedTime());
270
271 if (QR.info() != Eigen::Success)
272 log_and_throw_error("Failed to factorize constraints matrix");
273
274 timer.start();
276 Q = QR.matrixQ();
277 timer.stop();
278 logger().debug("Computation of Q took: {}", timer.getElapsedTime());
279
280 const Eigen::SparseMatrix<double, Eigen::RowMajor> R = QR.matrixR();
281
282 P_ = QR.colsPermutation();
283#endif
284
285 for (; constraint_size >= 0; --constraint_size)
286 {
287 const StiffnessMatrix tmp = R.row(constraint_size);
288 if (tmp.nonZeros() != 0)
289 {
290 constraint_size++;
291 break;
292 }
293 }
294 if (constraint_size != num_penalty_constraints_)
295 logger().warn("Matrix A is not full rank, constraint size: {} instead of {}", constraint_size, num_penalty_constraints_);
296
297 reduced_size_ = full_size_ - constraint_size;
298
299 timer.start();
300
301 Q1_ = Q.leftCols(constraint_size);
302 assert(Q1_.rows() == full_size_);
303 assert(Q1_.cols() == constraint_size);
304
305 Q2_ = Q.rightCols(reduced_size_);
306 Q2t_ = Q2_.transpose();
307
308 assert(Q2_.rows() == full_size_);
309 assert(Q2_.cols() == reduced_size_);
310
311 R1_ = R.topRows(constraint_size);
312 assert(R1_.rows() == constraint_size);
313 assert(R1_.cols() == num_penalty_constraints_);
314
315 assert((Q1_.transpose() * Q2_).norm() < 1e-10);
316
317 timer.stop();
318 logger().debug("Getting Q1 Q2, R1 took: {}", timer.getElapsedTime());
319
320 timer.start();
321
322 // arma::sp_mat q2a = fill_arma(Q2_);
323 // arma::sp_mat q2tq2 = q2a.t() * q2a;
324 // const StiffnessMatrix Q2tQ2 = fill_eigen(q2tq2);
325 StiffnessMatrix Q2tQ2 = Q2t_ * Q2_;
326 timer.stop();
327 logger().debug("Getting Q2'*Q2, took: {}", timer.getElapsedTime());
328
329 timer.start();
330 solver_->analyze_pattern(Q2tQ2, Q2tQ2.rows());
331 solver_->factorize(Q2tQ2);
332 timer.stop();
333 logger().debug("Factorization of Q2'*Q2 took: {}", timer.getElapsedTime());
334
335#ifndef NDEBUG
336 StiffnessMatrix test = R.bottomRows(reduced_size_);
337 assert(test.nonZeros() == 0);
338
339 StiffnessMatrix test1 = R1_.row(R1_.rows() - 1);
340 assert(test1.nonZeros() != 0);
341#endif
342
343 // assert((Q1_ * R1_ - At * P_).norm() < 1e-10);
344
345 std::vector<std::shared_ptr<Form>> tmp;
346 tmp.insert(tmp.end(), penalty_forms_.begin(), penalty_forms_.end());
347 penalty_problem_ = std::make_shared<FullNLProblem>(tmp);
348
350 }
351
353 {
354 if (penalty_forms_.size() == 1 && penalty_forms_.front()->has_projection())
355 {
356 Q1R1iTb_ = penalty_forms_.front()->constraint_projection_vector();
357 return;
358 }
359
360 igl::Timer timer;
361 timer.start();
362 // x = Q1 * R1^(-T) * P^T b + Q2 * y
363 int index = 0;
364 TVector constraint_values(num_penalty_constraints_);
365 for (const auto &f : penalty_forms_)
366 {
367 constraint_values.segment(index, f->constraint_value().rows()) = f->constraint_value();
368 index += f->constraint_value().rows();
369 }
370 constraint_values = P_.transpose() * constraint_values;
371
372 Eigen::VectorXd sol;
373
374 if (R1_.rows() == R1_.cols())
375 {
376 sol = R1_.transpose().triangularView<Eigen::Lower>().solve(constraint_values);
377 }
378 else
379 {
380
381#ifdef POLYSOLVE_WITH_SPQR
382 Eigen::SparseMatrix<double, Eigen::ColMajor, long> R1t = R1_.transpose();
383 cholmod_common cc;
384 cholmod_l_start(&cc); // start CHOLMOD
385 cholmod_sparse R1tc;
386 fill_cholmod(R1t, R1tc);
387
388 cholmod_dense b;
389 b.nrow = constraint_values.size();
390 b.ncol = 1;
391 b.nzmax = constraint_values.size();
392 b.d = constraint_values.size();
393 b.x = constraint_values.data();
394 b.z = 0;
395 b.xtype = CHOLMOD_REAL;
396 b.dtype = 0;
397
398 const int ordering = SPQR_ORDERING_DEFAULT; // all, except 3:given treated as 0:fixed
399 const double tol = SPQR_DEFAULT_TOL;
400
401 cholmod_dense *solc = SuiteSparseQR<double>(ordering, tol, &R1tc, &b, &cc);
402
403 sol = Eigen::Map<Eigen::VectorXd>(static_cast<double *>(solc->x), solc->nrow);
404
405 cholmod_l_free_dense(&solc, &cc);
406 cholmod_l_finish(&cc);
407#else
408 Eigen::SparseQR<StiffnessMatrix, Eigen::COLAMDOrdering<int>> solver;
409 solver.compute(R1_.transpose());
410 if (solver.info() != Eigen::Success)
411 {
412 log_and_throw_error("Failed to factorize R1^T");
413 }
414 sol = solver.solve(constraint_values);
415#endif
416 }
417
418 assert((R1_.transpose() * sol - constraint_values).norm() < 1e-10);
419
420 Q1R1iTb_ = Q1_ * sol;
421
422 timer.stop();
423 logger().debug("Computing Q1R1iTb took: {}", timer.getElapsedTime());
424 }
425
426 void NLProblem::init_lagging(const TVector &x)
427 {
429
431 penalty_problem_->init_lagging(x);
432 }
433
434 void NLProblem::update_lagging(const TVector &x, const int iter_num)
435 {
437
439 penalty_problem_->update_lagging(x, iter_num);
440 }
441
442 void NLProblem::update_quantities(const double t, const TVector &x)
443 {
444 t_ = t;
445 const TVector full = reduced_to_full(x);
446 assert(full.size() == full_size_);
447 for (auto &f : forms_)
448 f->update_quantities(t, full);
449
450 for (auto &f : penalty_forms_)
451 f->update_quantities(t, x);
452
454 }
455
456 void NLProblem::line_search_begin(const TVector &x0, const TVector &x1)
457 {
459
461 penalty_problem_->line_search_begin(x0, x1);
462 }
463
464 double NLProblem::max_step_size(const TVector &x0, const TVector &x1)
465 {
467
469 max_step = std::min(max_step, penalty_problem_->max_step_size(x0, x1));
470
471 return max_step;
472 }
473
474 bool NLProblem::is_step_valid(const TVector &x0, const TVector &x1)
475 {
477
478 if (penalty_problem_ && valid && full_size() == current_size())
479 valid = penalty_problem_->is_step_valid(x0, x1);
480
481 return valid;
482 }
483
484 bool NLProblem::is_step_collision_free(const TVector &x0, const TVector &x1)
485 {
487
488 if (penalty_problem_ && free && full_size() == current_size())
489 free = penalty_problem_->is_step_collision_free(x0, x1);
490
491 return free;
492 }
493
494 double NLProblem::value(const TVector &x)
495 {
496 // TODO: removed fearure const bool only_elastic
498
500 {
501 res += penalty_problem_->value(x);
502 }
503
504 return res;
505 }
506
507 void NLProblem::gradient(const TVector &x, TVector &grad)
508 {
510
511 if (full_size() != current_size())
512 {
513 if (penalty_forms_.size() == 1 && penalty_forms_.front()->can_project())
514 penalty_forms_.front()->project_gradient(grad);
515 else
516 grad = Q2t_ * grad;
517 }
518 else if (penalty_problem_)
519 {
520 TVector tmp;
521 penalty_problem_->gradient(x, tmp);
522 grad += tmp;
523 }
524 }
525
526 void NLProblem::hessian(const TVector &x, THessian &hessian)
527 {
529
530 if (full_size() != current_size())
531 {
533 }
534 else if (penalty_problem_)
535 {
536 THessian tmp;
537 penalty_problem_->hessian(x, tmp);
538 hessian += tmp;
539 }
540 }
541
542 void NLProblem::solution_changed(const TVector &newX)
543 {
545
547 penalty_problem_->solution_changed(newX);
548 }
549
550 void NLProblem::post_step(const polysolve::nonlinear::PostStepData &data)
551 {
552 FullNLProblem::post_step(polysolve::nonlinear::PostStepData(data.iter_num, data.solver_info, reduced_to_full(data.x), reduced_to_full(data.grad)));
553
555 penalty_problem_->post_step(data);
556
557 // TODO: add me back
558 static int nsolves = 0;
559 if (data.iter_num == 0)
560 nsolves++;
561 // if (state && state->args["output"]["advanced"]["save_nl_solve_sequence"])
562 // {
563 // const Eigen::MatrixXd displacements = utils::unflatten(reduced_to_full(data.x), state->mesh->dimension());
564 // io::OBJWriter::write(
565 // state->resolve_output_path(fmt::format("nonlinear_solve{:04d}_iter{:04d}.obj", nsolves, data.iter_num)),
566 // state->collision_mesh.displace_vertices(displacements),
567 // state->collision_mesh.edges(), state->collision_mesh.faces());
568 // }
569 }
570
571 NLProblem::TVector NLProblem::full_to_reduced(const TVector &full) const
572 {
573 // Reduced is already at the full size
574 if (full_size() == current_size() || full.size() == current_size())
575 {
576 return full;
577 }
578
579 TVector reduced(reduced_size());
580 const TVector k = full - Q1R1iTb_;
581 const TVector rhs = Q2t_ * k;
582 solver_->solve(rhs, reduced);
583
584#ifndef NDEBUG
585 StiffnessMatrix Q2tQ2 = Q2t_ * Q2_;
586 // std::cout << "err " << (Q2tQ2 * reduced - rhs).norm() << std::endl;
587 assert((Q2tQ2 * reduced - rhs).norm() < 1e-8);
588#endif
589
590 return reduced;
591 }
592
593 NLProblem::TVector NLProblem::full_to_reduced_grad(const TVector &full) const
594 {
595 TVector grad = full;
596 if (penalty_forms_.size() == 1 && penalty_forms_.front()->can_project())
597 penalty_forms_.front()->project_gradient(grad);
598 else
599 grad = Q2t_ * grad;
600
601 return grad;
602 }
603
604 NLProblem::TVector NLProblem::reduced_to_full(const TVector &reduced) const
605 {
606 // Full is already at the reduced size
607 if (full_size() == current_size() || full_size() == reduced.size())
608 {
609 return reduced;
610 }
611
612 // x = Q1 * R1^(-T) * P^T b + Q2 * y
613
614 const TVector full = Q1R1iTb_ + Q2_ * reduced;
615
616#ifndef NDEBUG
617 for (const auto &f : penalty_forms_)
618 {
619 // std::cout << f->compute_error(full) << std::endl;
620 assert(f->compute_error(full) < 1e-8);
621 }
622#endif
623
624 return full;
625 }
626
628 {
629 if (penalty_forms_.size() == 1 && penalty_forms_.front()->can_project())
630 penalty_forms_.front()->project_hessian(hessian);
631 else
632 {
633 // arma::sp_mat q2a = fill_arma(Q2_);
634 // arma::sp_mat ha = fill_arma(hessian);
635 // arma::sp_mat q2thq2 = q2a.t() * ha * q2a;
636 // hessian = fill_eigen(q2thq2);
637 hessian = Q2t_ * hessian * Q2_;
638 // remove numerical zeros
639 hessian.prune([](const Eigen::Index &row, const Eigen::Index &col, const Scalar &value) {
640 return std::abs(value) > 1e-10;
641 });
642 }
643 }
644} // 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:71
StiffnessMatrix R1_
R1 block of the QR decomposition of the constraints matrix.
Definition NLProblem.hpp:93
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.
Definition NLProblem.hpp:94
double normalize_forms() 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
int reduced_size_
Size of the reduced problem.
Definition NLProblem.hpp:72
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
StiffnessMatrix Q2t_
Q2 transpose.
Definition NLProblem.hpp:92
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
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_
Definition NLProblem.hpp:88
virtual double value(const TVector &x) override
StiffnessMatrix Q1_
Q1 block of the QR decomposition of the constraints matrix.
Definition NLProblem.hpp:90
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_
Definition NLProblem.hpp:96
void solution_changed(const TVector &new_x) override
std::shared_ptr< FullNLProblem > penalty_problem_
Definition NLProblem.hpp:98
TVector Q1R1iTb_
Q1_ * (R1_.transpose().triangularView<Eigen::Upper>().solve(constraint_values_))
Definition NLProblem.hpp:95
void full_hessian_to_reduced_hessian(StiffnessMatrix &hessian) const
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