Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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#include <armadillo>
16
17#include <igl/cat.h>
18#include <igl/Timer.h>
19
20#include <set>
21
22/*
23m \frac{\partial^2 u}{\partial t^2} = \psi = \text{div}(\sigma[u])\newline
24u^{t+1} = u(t+\Delta t)\approx u(t) + \Delta t \dot u + \frac{\Delta t^2} 2 \ddot u \newline
25= u(t) + \Delta t \dot u + \frac{\Delta t^2}{2} \psi\newline
26M 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
27%
28M (u^{t+1}_h - (u^t_h + \Delta t v^t_h)) - \frac{\Delta t^2} {2} A u^{t+1}_h
29*/
30// mü = ψ = div(σ[u])
31// uᵗ⁺¹ = u(t + Δt) ≈ u(t) + Δtu̇ + ½Δt²ü = u(t) + Δtu̇ + ½Δt²ψ
32// Muₕᵗ⁺¹ ≈ Muₕᵗ + ΔtMvₕᵗ ½Δt²Auₕᵗ⁺¹
33// Root-finding form:
34// M(uₕᵗ⁺¹ - (uₕᵗ + Δtvₕᵗ)) - ½Δt²Auₕᵗ⁺¹ = 0
35
36namespace polyfem::solver
37{
38
39 namespace
40 {
41#ifdef POLYSOLVE_WITH_SPQR
42 void fill_cholmod(Eigen::SparseMatrix<double, Eigen::ColMajor, long> &mat, cholmod_sparse &cmat)
43 {
44 long *p = mat.outerIndexPtr();
45
46 cmat.nzmax = mat.nonZeros();
47 cmat.nrow = mat.rows();
48 cmat.ncol = mat.cols();
49 cmat.p = p;
50 cmat.i = mat.innerIndexPtr();
51 cmat.x = mat.valuePtr();
52 cmat.z = 0;
53 cmat.sorted = 1;
54 cmat.packed = 1;
55 cmat.nz = 0;
56 cmat.dtype = 0;
57 cmat.stype = -1;
58 cmat.xtype = CHOLMOD_REAL;
59 cmat.dtype = CHOLMOD_DOUBLE;
60 cmat.stype = 0;
61 cmat.itype = CHOLMOD_LONG;
62 }
63#endif
64 arma::sp_mat fill_arma(const StiffnessMatrix &mat)
65 {
66 std::vector<unsigned long long> rowind_vect(mat.innerIndexPtr(), mat.innerIndexPtr() + mat.nonZeros());
67 std::vector<unsigned long long> colptr_vect(mat.outerIndexPtr(), mat.outerIndexPtr() + mat.outerSize() + 1);
68 std::vector<double> values_vect(mat.valuePtr(), mat.valuePtr() + mat.nonZeros());
69
70 arma::dvec values(values_vect.data(), values_vect.size(), false);
71 arma::uvec rowind(rowind_vect.data(), rowind_vect.size(), false);
72 arma::uvec colptr(colptr_vect.data(), colptr_vect.size(), false);
73
74 arma::sp_mat amat(rowind, colptr, values, mat.rows(), mat.cols(), false);
75
76 return amat;
77 }
78
79 StiffnessMatrix fill_eigen(const arma::sp_mat &mat)
80 {
81 // convert to eigen sparse
82 std::vector<long> outerIndexPtr(mat.row_indices, mat.row_indices + mat.n_nonzero);
83 std::vector<long> innerIndexPtr(mat.col_ptrs, mat.col_ptrs + mat.n_cols + 1);
84
85 const StiffnessMatrix out = Eigen::Map<const Eigen::SparseMatrix<double, Eigen::ColMajor, long>>(
86 mat.n_rows, mat.n_cols, mat.n_nonzero, innerIndexPtr.data(), outerIndexPtr.data(), mat.values);
87 return out;
88 }
89
90 } // namespace
92 const int full_size,
93 const std::vector<std::shared_ptr<Form>> &forms,
94 const std::vector<std::shared_ptr<AugmentedLagrangianForm>> &penalty_forms,
95 const std::shared_ptr<polysolve::linear::Solver> &solver)
96 : FullNLProblem(forms),
97 full_size_(full_size),
98 t_(0),
99 penalty_forms_(penalty_forms),
100 solver_(solver)
101 {
104 }
105
107 const int full_size,
108 const std::shared_ptr<utils::PeriodicBoundary> &periodic_bc,
109 const double t,
110 const std::vector<std::shared_ptr<Form>> &forms,
111 const std::vector<std::shared_ptr<AugmentedLagrangianForm>> &penalty_forms,
112 const std::shared_ptr<polysolve::linear::Solver> &solver)
113 : FullNLProblem(forms),
114 full_size_(full_size),
115 t_(t),
116 penalty_forms_(penalty_forms),
117 solver_(solver)
118 {
121 }
122
124 {
125 return 1;
126
127 // double total_weight = 0;
128 // for (const auto &f : forms_)
129 // total_weight += f->weight();
130 // if (full_size() == current_size())
131 // {
132 // for (const auto &f : penalty_forms_)
133 // total_weight += f->weight() * f->lagrangian_weight();
134 // }
135
136 // logger().debug("Normalizing forms with scale: {}", total_weight);
137
138 // for (auto &f : forms_)
139 // f->set_scale(total_weight);
140 // for (auto &f : penalty_forms_)
141 // f->set_scale(total_weight);
142
143 // return total_weight;
144 }
145
147 {
148 if (penalty_forms_.empty())
149 {
151 return;
152 }
153 igl::Timer timer;
154
155 if (penalty_forms_.size() == 1 && penalty_forms_.front()->has_projection())
156 {
157 Q2_ = penalty_forms_.front()->constraint_projection_matrix();
158 Q2t_ = Q2_.transpose();
159
160 reduced_size_ = Q2_.cols();
162
163 timer.start();
164 StiffnessMatrix Q2tQ2 = Q2t_ * Q2_;
165 solver_->analyze_pattern(Q2tQ2, Q2tQ2.rows());
166 solver_->factorize(Q2tQ2);
167 timer.stop();
168 logger().debug("Factorization and computation of Q2tQ2 took: {}", timer.getElapsedTime());
169
170 std::vector<std::shared_ptr<Form>> tmp;
171 tmp.insert(tmp.end(), penalty_forms_.begin(), penalty_forms_.end());
172 penalty_problem_ = std::make_shared<FullNLProblem>(tmp);
173
175
176 return;
177 }
178
179 std::vector<Eigen::Triplet<double>> Ae;
180 int index = 0;
181 for (const auto &f : penalty_forms_)
182 {
183 const auto &tmp = f->constraint_matrix();
184 for (int i = 0; i < tmp.outerSize(); i++)
185 {
186 for (typename StiffnessMatrix::InnerIterator it(tmp, i); it; ++it)
187 {
188 Ae.emplace_back(index + it.row(), it.col(), it.value());
189 }
190 }
191 index += tmp.rows();
192 }
193 StiffnessMatrix A(index, full_size_);
194 A.setFromTriplets(Ae.begin(), Ae.end());
195 A.makeCompressed();
196
197 int constraint_size = A.rows();
198 num_penalty_constraints_ = A.rows();
199 Eigen::SparseMatrix<double, Eigen::ColMajor, long> At = A.transpose();
200 At.makeCompressed();
201
202 logger().debug("Constraint size: {} x {}", A.rows(), A.cols());
203
204#ifdef POLYSOLVE_WITH_SPQR
205 timer.start();
206
207 cholmod_common cc;
208 cholmod_l_start(&cc); // start CHOLMOD
209
210 // const int ordering = 0; // all, except 3:given treated as 0:fixed
211 const int ordering = SPQR_ORDERING_DEFAULT; // all, except 3:given treated as 0:fixed
212 const double tol = SPQR_DEFAULT_TOL;
213 SuiteSparse_long econ = At.rows();
214 SuiteSparse_long *E; // permutation of 0:n-1, NULL if identity
215
216 cholmod_sparse Ac, *Qc, *Rc;
217
218 fill_cholmod(At, Ac);
219
220 const auto rank = SuiteSparseQR<double>(ordering, tol, econ, &Ac,
221 // outputs
222 &Qc, &Rc, &E, &cc);
223
224 if (!Rc)
225 log_and_throw_error("Failed to factorize constraints matrix");
226
227 const auto n = Rc->ncol;
228 P_.resize(n);
229 if (E)
230 {
231 for (long j = 0; j < n; j++)
232 P_.indices()(j) = E[j];
233 }
234 else
235 P_.setIdentity();
236
237 if (Qc->stype != 0 || Qc->sorted != 1 || Qc->packed != 1 || Rc->stype != 0 || Rc->sorted != 1 || Rc->packed != 1)
238 log_and_throw_error("Q and R must be unsymmetric sorted and packed");
239
240 const StiffnessMatrix Q = Eigen::Map<Eigen::SparseMatrix<double, Eigen::ColMajor, long>>(
241 Qc->nrow, Qc->ncol, Qc->nzmax,
242 static_cast<long *>(Qc->p), static_cast<long *>(Qc->i), static_cast<double *>(Qc->x));
243
244 const StiffnessMatrix R = Eigen::Map<Eigen::SparseMatrix<double, Eigen::ColMajor, long>>(
245 Rc->nrow, Rc->ncol, Rc->nzmax,
246 static_cast<long *>(Rc->p), static_cast<long *>(Rc->i), static_cast<double *>(Rc->x));
247
248 cholmod_l_free_sparse(&Qc, &cc);
249 cholmod_l_free_sparse(&Rc, &cc);
250 std::free(E);
251 cholmod_l_finish(&cc);
252
253 timer.stop();
254 logger().debug("QR took: {}", timer.getElapsedTime());
255#else
256 timer.start();
257
258 // Eigen::SparseQR<StiffnessMatrix, Eigen::NaturalOrdering<int>> QR(At);
259 Eigen::SparseQR<StiffnessMatrix, Eigen::COLAMDOrdering<int>> QR(At);
260
261 timer.stop();
262 logger().debug("QR took: {}", timer.getElapsedTime());
263
264 if (QR.info() != Eigen::Success)
265 log_and_throw_error("Failed to factorize constraints matrix");
266
267 timer.start();
269 Q = QR.matrixQ();
270 timer.stop();
271 logger().debug("Computation of Q took: {}", timer.getElapsedTime());
272
273 const Eigen::SparseMatrix<double, Eigen::RowMajor> R = QR.matrixR();
274
275 P_ = QR.colsPermutation();
276#endif
277
278 for (; constraint_size >= 0; --constraint_size)
279 {
280 const StiffnessMatrix tmp = R.row(constraint_size);
281 if (tmp.nonZeros() != 0)
282 {
283 constraint_size++;
284 break;
285 }
286 }
287 if (constraint_size != num_penalty_constraints_)
288 logger().warn("Matrix A is not full rank, constraint size: {} instead of {}", constraint_size, num_penalty_constraints_);
289
290 reduced_size_ = full_size_ - constraint_size;
291
292 timer.start();
293
294 Q1_ = Q.leftCols(constraint_size);
295 assert(Q1_.rows() == full_size_);
296 assert(Q1_.cols() == constraint_size);
297
298 Q2_ = Q.rightCols(reduced_size_);
299 Q2t_ = Q2_.transpose();
300
301 assert(Q2_.rows() == full_size_);
302 assert(Q2_.cols() == reduced_size_);
303
304 R1_ = R.topRows(constraint_size);
305 assert(R1_.rows() == constraint_size);
306 assert(R1_.cols() == num_penalty_constraints_);
307
308 assert((Q1_.transpose() * Q2_).norm() < 1e-10);
309
310 timer.stop();
311 logger().debug("Getting Q1 Q2, R1 took: {}", timer.getElapsedTime());
312
313 timer.start();
314
315 // arma::sp_mat q2a = fill_arma(Q2_);
316 // arma::sp_mat q2tq2 = q2a.t() * q2a;
317 // const StiffnessMatrix Q2tQ2 = fill_eigen(q2tq2);
318 StiffnessMatrix Q2tQ2 = Q2t_ * Q2_;
319 timer.stop();
320 logger().debug("Getting Q2'*Q2, took: {}", timer.getElapsedTime());
321
322 timer.start();
323 solver_->analyze_pattern(Q2tQ2, Q2tQ2.rows());
324 solver_->factorize(Q2tQ2);
325 timer.stop();
326 logger().debug("Factorization of Q2'*Q2 took: {}", timer.getElapsedTime());
327
328#ifndef NDEBUG
329 StiffnessMatrix test = R.bottomRows(reduced_size_);
330 assert(test.nonZeros() == 0);
331
332 StiffnessMatrix test1 = R1_.row(R1_.rows() - 1);
333 assert(test1.nonZeros() != 0);
334#endif
335
336 // assert((Q1_ * R1_ - At * P_).norm() < 1e-10);
337
338 std::vector<std::shared_ptr<Form>> tmp;
339 tmp.insert(tmp.end(), penalty_forms_.begin(), penalty_forms_.end());
340 penalty_problem_ = std::make_shared<FullNLProblem>(tmp);
341
343 }
344
346 {
347 if (penalty_forms_.size() == 1 && penalty_forms_.front()->has_projection())
348 {
349 Q1R1iTb_ = penalty_forms_.front()->constraint_projection_vector();
350 return;
351 }
352
353 igl::Timer timer;
354 timer.start();
355 // x = Q1 * R1^(-T) * P^T b + Q2 * y
356 int index = 0;
357 TVector constraint_values(num_penalty_constraints_);
358 for (const auto &f : penalty_forms_)
359 {
360 constraint_values.segment(index, f->constraint_value().rows()) = f->constraint_value();
361 index += f->constraint_value().rows();
362 }
363 constraint_values = P_.transpose() * constraint_values;
364
365 Eigen::VectorXd sol;
366
367 if (R1_.rows() == R1_.cols())
368 {
369 sol = R1_.transpose().triangularView<Eigen::Lower>().solve(constraint_values);
370 }
371 else
372 {
373
374#ifdef POLYSOLVE_WITH_SPQR
375 Eigen::SparseMatrix<double, Eigen::ColMajor, long> R1t = R1_.transpose();
376 cholmod_common cc;
377 cholmod_l_start(&cc); // start CHOLMOD
378 cholmod_sparse R1tc;
379 fill_cholmod(R1t, R1tc);
380
381 cholmod_dense b;
382 b.nrow = constraint_values.size();
383 b.ncol = 1;
384 b.nzmax = constraint_values.size();
385 b.d = constraint_values.size();
386 b.x = constraint_values.data();
387 b.z = 0;
388 b.xtype = CHOLMOD_REAL;
389 b.dtype = 0;
390
391 const int ordering = SPQR_ORDERING_DEFAULT; // all, except 3:given treated as 0:fixed
392 const double tol = SPQR_DEFAULT_TOL;
393
394 cholmod_dense *solc = SuiteSparseQR<double>(ordering, tol, &R1tc, &b, &cc);
395
396 sol = Eigen::Map<Eigen::VectorXd>(static_cast<double *>(solc->x), solc->nrow);
397
398 cholmod_l_free_dense(&solc, &cc);
399 cholmod_l_finish(&cc);
400#else
401 Eigen::SparseQR<StiffnessMatrix, Eigen::COLAMDOrdering<int>> solver;
402 solver.compute(R1_.transpose());
403 if (solver.info() != Eigen::Success)
404 {
405 log_and_throw_error("Failed to factorize R1^T");
406 }
407 sol = solver.solve(constraint_values);
408#endif
409 }
410
411 assert((R1_.transpose() * sol - constraint_values).norm() < 1e-10);
412
413 Q1R1iTb_ = Q1_ * sol;
414
415 timer.stop();
416 logger().debug("Computing Q1R1iTb took: {}", timer.getElapsedTime());
417 }
418
419 void NLProblem::init_lagging(const TVector &x)
420 {
422
424 penalty_problem_->init_lagging(x);
425 }
426
427 void NLProblem::update_lagging(const TVector &x, const int iter_num)
428 {
430
432 penalty_problem_->update_lagging(x, iter_num);
433 }
434
435 void NLProblem::update_quantities(const double t, const TVector &x)
436 {
437 t_ = t;
438 const TVector full = reduced_to_full(x);
439 assert(full.size() == full_size_);
440 for (auto &f : forms_)
441 f->update_quantities(t, full);
442
443 for (auto &f : penalty_forms_)
444 f->update_quantities(t, x);
445
447 }
448
449 void NLProblem::line_search_begin(const TVector &x0, const TVector &x1)
450 {
452
454 penalty_problem_->line_search_begin(x0, x1);
455 }
456
457 double NLProblem::max_step_size(const TVector &x0, const TVector &x1)
458 {
460
462 max_step = std::min(max_step, penalty_problem_->max_step_size(x0, x1));
463
464 return max_step;
465 }
466
467 bool NLProblem::is_step_valid(const TVector &x0, const TVector &x1)
468 {
470
471 if (penalty_problem_ && valid && full_size() == current_size())
472 valid = penalty_problem_->is_step_valid(x0, x1);
473
474 return valid;
475 }
476
477 bool NLProblem::is_step_collision_free(const TVector &x0, const TVector &x1)
478 {
480
481 if (penalty_problem_ && free && full_size() == current_size())
482 free = penalty_problem_->is_step_collision_free(x0, x1);
483
484 return free;
485 }
486
487 double NLProblem::value(const TVector &x)
488 {
489 // TODO: removed fearure const bool only_elastic
491
493 {
494 res += penalty_problem_->value(x);
495 }
496
497 return res;
498 }
499
500 void NLProblem::gradient(const TVector &x, TVector &grad)
501 {
503
504 if (full_size() != current_size())
505 {
506 if (penalty_forms_.size() == 1 && penalty_forms_.front()->can_project())
507 penalty_forms_.front()->project_gradient(grad);
508 else
509 grad = Q2t_ * grad;
510 }
511 else if (penalty_problem_)
512 {
513 TVector tmp;
514 penalty_problem_->gradient(x, tmp);
515 grad += tmp;
516 }
517 }
518
519 void NLProblem::hessian(const TVector &x, THessian &hessian)
520 {
522
523 if (full_size() != current_size())
524 {
526 }
527 else if (penalty_problem_)
528 {
529 THessian tmp;
530 penalty_problem_->hessian(x, tmp);
531 hessian += tmp;
532 }
533 }
534
535 void NLProblem::solution_changed(const TVector &newX)
536 {
538
540 penalty_problem_->solution_changed(newX);
541 }
542
543 void NLProblem::post_step(const polysolve::nonlinear::PostStepData &data)
544 {
545 FullNLProblem::post_step(polysolve::nonlinear::PostStepData(data.iter_num, data.solver_info, reduced_to_full(data.x), reduced_to_full(data.grad)));
546
548 penalty_problem_->post_step(data);
549
550 // TODO: add me back
551 static int nsolves = 0;
552 if (data.iter_num == 0)
553 nsolves++;
554 // if (state && state->args["output"]["advanced"]["save_nl_solve_sequence"])
555 // {
556 // const Eigen::MatrixXd displacements = utils::unflatten(reduced_to_full(data.x), state->mesh->dimension());
557 // io::OBJWriter::write(
558 // state->resolve_output_path(fmt::format("nonlinear_solve{:04d}_iter{:04d}.obj", nsolves, data.iter_num)),
559 // state->collision_mesh.displace_vertices(displacements),
560 // state->collision_mesh.edges(), state->collision_mesh.faces());
561 // }
562 }
563
564 NLProblem::TVector NLProblem::full_to_reduced(const TVector &full) const
565 {
566 // Reduced is already at the full size
567 if (full_size() == current_size() || full.size() == current_size())
568 {
569 return full;
570 }
571
572 TVector reduced(reduced_size());
573 const TVector k = full - Q1R1iTb_;
574 const TVector rhs = Q2t_ * k;
575 solver_->solve(rhs, reduced);
576
577#ifndef NDEBUG
578 StiffnessMatrix Q2tQ2 = Q2t_ * Q2_;
579 // std::cout << "err " << (Q2tQ2 * reduced - rhs).norm() << std::endl;
580 assert((Q2tQ2 * reduced - rhs).norm() < 1e-8);
581#endif
582
583 return reduced;
584 }
585
586 NLProblem::TVector NLProblem::full_to_reduced_grad(const TVector &full) const
587 {
588 TVector grad = full;
589 if (penalty_forms_.size() == 1 && penalty_forms_.front()->can_project())
590 penalty_forms_.front()->project_gradient(grad);
591 else
592 grad = Q2t_ * grad;
593
594 return grad;
595 }
596
597 NLProblem::TVector NLProblem::reduced_to_full(const TVector &reduced) const
598 {
599 // Full is already at the reduced size
600 if (full_size() == current_size() || full_size() == reduced.size())
601 {
602 return reduced;
603 }
604
605 // x = Q1 * R1^(-T) * P^T b + Q2 * y
606
607 const TVector full = Q1R1iTb_ + Q2_ * reduced;
608
609#ifndef NDEBUG
610 for (const auto &f : penalty_forms_)
611 {
612 // std::cout << f->compute_error(full) << std::endl;
613 assert(f->compute_error(full) < 1e-8);
614 }
615#endif
616
617 return full;
618 }
619
621 {
622 if (penalty_forms_.size() == 1 && penalty_forms_.front()->can_project())
623 penalty_forms_.front()->project_hessian(hessian);
624 else
625 {
626 // arma::sp_mat q2a = fill_arma(Q2_);
627 // arma::sp_mat ha = fill_arma(hessian);
628 // arma::sp_mat q2thq2 = q2a.t() * ha * q2a;
629 // hessian = fill_eigen(q2thq2);
630 hessian = Q2t_ * hessian * Q2_;
631 // remove numerical zeros
632 hessian.prune([](const Eigen::Index &row, const Eigen::Index &col, const Scalar &value) {
633 return std::abs(value) > 1e-10;
634 });
635 }
636 }
637} // 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:91
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:22