Janus 2.0.0
High-performance C++20 dual-mode numerical framework
Loading...
Searching...
No Matches
Linalg.hpp
Go to the documentation of this file.
1#pragma once
7
12#include "janus/math/Logic.hpp"
13#include "janus/math/Trig.hpp"
14#include <Eigen/Dense>
15#include <Eigen/IterativeLinearSolvers>
16#include <Eigen/SparseCholesky>
17#include <Eigen/SparseLU>
18#include <Eigen/SparseQR>
19#include <algorithm>
20#include <array>
21#include <casadi/casadi.hpp>
22#include <cmath>
23#include <complex>
24#include <functional>
25#include <limits>
26#include <numeric>
27#include <string>
28#include <unsupported/Eigen/IterativeSolvers>
29#include <vector>
30
31namespace janus {
32
33// --- Conversion Helpers ---
34// (Moved to JanusTypes.hpp: to_mx, to_eigen, as_vector)
35
44
55
65
73
81
92 double tolerance = 1e-10;
93 int max_iterations = 500;
94 int gmres_restart = 30;
97 casadi::Dict symbolic_options;
98
101 LinearSolvePolicy policy;
103 policy.dense_solver = solver;
104 return policy;
105 }
106
107 static LinearSolvePolicy
114
115 static LinearSolvePolicy
124
126 tolerance = value;
127 return *this;
128 }
129
131 max_iterations = value;
132 return *this;
133 }
134
136 gmres_restart = value;
137 return *this;
138 }
139
142 preconditioner_hook = std::move(hook);
143 return *this;
144 }
145
146 LinearSolvePolicy &set_symbolic_solver(const std::string &solver,
147 const casadi::Dict &opts = casadi::Dict()) {
148 symbolic_linear_solver = solver;
149 symbolic_options = opts;
150 return *this;
151 }
152};
153
154namespace detail {
155
157 public:
158 using Scalar = double;
159 using RealScalar = double;
160 using StorageIndex = int;
161 enum { ColsAtCompileTime = Eigen::Dynamic, MaxColsAtCompileTime = Eigen::Dynamic };
162
164
165 void set_apply(std::function<NumericVector(const NumericVector &)> apply) {
166 apply_ = std::move(apply);
167 }
168
169 Eigen::Index rows() const noexcept { return size_; }
170 Eigen::Index cols() const noexcept { return size_; }
171
172 template <typename MatType> FunctionalPreconditioner &analyzePattern(const MatType &) {
173 return *this;
174 }
175
176 template <typename MatType> FunctionalPreconditioner &factorize(const MatType &mat) {
177 size_ = mat.cols();
178 initialized_ = true;
179 return *this;
180 }
181
182 template <typename MatType> FunctionalPreconditioner &compute(const MatType &mat) {
183 return factorize(mat);
184 }
185
186 template <typename Rhs, typename Dest> void _solve_impl(const Rhs &b, Dest &x) const {
187 NumericVector rhs = b;
188 x = apply_ ? apply_(rhs) : rhs;
189 }
190
191 template <typename Rhs>
192 inline const Eigen::Solve<FunctionalPreconditioner, Rhs>
193 solve(const Eigen::MatrixBase<Rhs> &b) const {
194 eigen_assert(initialized_ && "FunctionalPreconditioner is not initialized.");
195 eigen_assert(size_ == b.rows() && "FunctionalPreconditioner::solve(): invalid rhs size");
196 return Eigen::Solve<FunctionalPreconditioner, Rhs>(*this, b.derived());
197 }
198
199 Eigen::ComputationInfo info() const { return Eigen::Success; }
200
201 private:
202 std::function<NumericVector(const NumericVector &)> apply_;
203 Eigen::Index size_ = 0;
204 bool initialized_ = false;
205};
206
207inline void validate_linear_solve_dims(Eigen::Index a_rows, Eigen::Index a_cols,
208 Eigen::Index b_rows, const std::string &context) {
209 if (a_rows == 0 || a_cols == 0) {
210 throw InvalidArgument(context + ": coefficient matrix must be non-empty");
211 }
212 if (b_rows != a_rows) {
213 throw InvalidArgument(context + ": rhs row count must match coefficient matrix rows");
214 }
215}
216
217inline void validate_square_required(Eigen::Index rows, Eigen::Index cols,
218 const std::string &context, const std::string &solver_name) {
219 if (rows != cols) {
220 throw InvalidArgument(context + ": " + solver_name + " requires a square matrix");
221 }
222}
223
224inline void validate_iterative_policy(const LinearSolvePolicy &policy, const std::string &context) {
225 if (policy.tolerance <= 0.0) {
226 throw InvalidArgument(context + ": tolerance must be positive");
227 }
228 if (policy.max_iterations <= 0) {
229 throw InvalidArgument(context + ": max_iterations must be positive");
230 }
231 if (policy.gmres_restart <= 0) {
232 throw InvalidArgument(context + ": gmres_restart must be positive");
233 }
234}
235
236inline std::function<NumericVector(const NumericVector &)>
238 if (policy.preconditioner_hook) {
239 return policy.preconditioner_hook;
240 }
241
242 switch (policy.iterative_preconditioner) {
244 return [](const NumericVector &rhs) { return rhs; };
246 NumericVector inv_diag(A.rows());
247 const double eps = std::numeric_limits<double>::epsilon();
248 for (int i = 0; i < A.rows(); ++i) {
249 const double diag = A.coeff(i, i);
250 inv_diag(i) = std::abs(diag) > eps ? 1.0 / diag : 1.0;
251 }
252 return [inv_diag](const NumericVector &rhs) { return inv_diag.array() * rhs.array(); };
253 }
254 }
255 throw InvalidArgument("solve: unsupported iterative preconditioner");
256}
257
258template <typename MatrixLike> SparseMatrix dense_to_sparse(const MatrixLike &A) {
259 SparseMatrix sparse = A.sparseView();
260 sparse.makeCompressed();
261 return sparse;
262}
263
264template <typename DerivedB>
265auto solve_sparse_direct_numeric(const SparseMatrix &A, const Eigen::MatrixBase<DerivedB> &b,
266 const LinearSolvePolicy &policy) {
267 using Result = Eigen::Matrix<double, Eigen::Dynamic, DerivedB::ColsAtCompileTime>;
268
269 SparseMatrix matrix = A;
270 matrix.makeCompressed();
271
272 switch (policy.sparse_direct_solver) {
274 Eigen::SparseLU<SparseMatrix> solver;
275 solver.compute(matrix);
276 if (solver.info() != Eigen::Success) {
277 throw InvalidArgument("solve: SparseLU factorization failed");
278 }
279 return solver.solve(b).eval();
280 }
282 Eigen::SparseQR<SparseMatrix, Eigen::COLAMDOrdering<int>> solver;
283 solver.compute(matrix);
284 if (solver.info() != Eigen::Success) {
285 throw InvalidArgument("solve: SparseQR factorization failed");
286 }
287 return solver.solve(b).eval();
288 }
290 validate_square_required(matrix.rows(), matrix.cols(), "solve", "SimplicialLLT");
291 Eigen::SimplicialLLT<SparseMatrix> solver;
292 solver.compute(matrix);
293 if (solver.info() != Eigen::Success) {
294 throw InvalidArgument("solve: SimplicialLLT factorization failed");
295 }
296 return solver.solve(b).eval();
297 }
299 validate_square_required(matrix.rows(), matrix.cols(), "solve", "SimplicialLDLT");
300 Eigen::SimplicialLDLT<SparseMatrix> solver;
301 solver.compute(matrix);
302 if (solver.info() != Eigen::Success) {
303 throw InvalidArgument("solve: SimplicialLDLT factorization failed");
304 }
305 return solver.solve(b).eval();
306 }
307 }
308
309 return Result();
310}
311
312template <typename Solver, typename DerivedB>
313auto solve_iterative_with_solver(Solver &solver, const Eigen::MatrixBase<DerivedB> &b) {
314 using Result = Eigen::Matrix<double, Eigen::Dynamic, DerivedB::ColsAtCompileTime>;
315
316 Result x(b.rows(), b.cols());
317 for (int col = 0; col < b.cols(); ++col) {
318 NumericVector rhs = b.col(col);
319 NumericVector sol = solver.solve(rhs);
320 if (solver.info() != Eigen::Success) {
321 throw InvalidArgument("solve: iterative solver failed to converge");
322 }
323 x.col(col) = sol;
324 }
325 return x;
326}
327
328template <typename DerivedB>
329auto solve_iterative_numeric(const SparseMatrix &A, const Eigen::MatrixBase<DerivedB> &b,
330 const LinearSolvePolicy &policy) {
331 validate_square_required(A.rows(), A.cols(), "solve", "iterative Krylov backends");
332 validate_iterative_policy(policy, "solve");
333
334 FunctionalPreconditioner preconditioner;
335 preconditioner.set_apply(make_preconditioner(A, policy));
336
337 switch (policy.iterative_solver) {
339 Eigen::BiCGSTAB<SparseMatrix, FunctionalPreconditioner> solver;
340 solver.setTolerance(policy.tolerance);
341 solver.setMaxIterations(policy.max_iterations);
342 solver.preconditioner() = preconditioner;
343 solver.compute(A);
344 return solve_iterative_with_solver(solver, b);
345 }
347 Eigen::GMRES<SparseMatrix, FunctionalPreconditioner> solver;
348 solver.setTolerance(policy.tolerance);
349 solver.setMaxIterations(policy.max_iterations);
350 solver.set_restart(policy.gmres_restart);
351 solver.preconditioner() = preconditioner;
352 solver.compute(A);
353 return solve_iterative_with_solver(solver, b);
354 }
355 }
356
357 using Result = Eigen::Matrix<double, Eigen::Dynamic, DerivedB::ColsAtCompileTime>;
358 return Result();
359}
360
361template <typename DerivedA, typename DerivedB>
362auto solve_dense_numeric(const Eigen::MatrixBase<DerivedA> &A, const Eigen::MatrixBase<DerivedB> &b,
363 const LinearSolvePolicy &policy) {
364 switch (policy.dense_solver) {
366 return A.colPivHouseholderQr().solve(b).eval();
368 validate_square_required(A.rows(), A.cols(), "solve", "PartialPivLU");
369 return A.partialPivLu().solve(b).eval();
371 validate_square_required(A.rows(), A.cols(), "solve", "FullPivLU");
372 return A.fullPivLu().solve(b).eval();
374 validate_square_required(A.rows(), A.cols(), "solve", "LLT");
375 Eigen::LLT<NumericMatrix> solver(A.eval());
376 if (solver.info() != Eigen::Success) {
377 throw InvalidArgument("solve: LLT factorization failed");
378 }
379 return solver.solve(b).eval();
380 }
382 validate_square_required(A.rows(), A.cols(), "solve", "LDLT");
383 Eigen::LDLT<NumericMatrix> solver(A.eval());
384 if (solver.info() != Eigen::Success) {
385 throw InvalidArgument("solve: LDLT factorization failed");
386 }
387 return solver.solve(b).eval();
388 }
389 }
390
391 return A.colPivHouseholderQr().solve(b).eval();
392}
393
394} // namespace detail
395
396// --- solve(A, b) ---
407template <typename DerivedA, typename DerivedB>
408auto solve(const Eigen::MatrixBase<DerivedA> &A, const Eigen::MatrixBase<DerivedB> &b) {
409 return solve(A, b, LinearSolvePolicy());
410}
411
415template <typename DerivedA, typename DerivedB>
416auto solve(const Eigen::MatrixBase<DerivedA> &A, const Eigen::MatrixBase<DerivedB> &b,
417 const LinearSolvePolicy &policy) {
418 using Scalar = typename DerivedA::Scalar;
419 detail::validate_linear_solve_dims(A.rows(), A.cols(), b.rows(), "solve");
420
421 if constexpr (std::is_floating_point_v<Scalar>) {
422 switch (policy.backend) {
424 return detail::solve_dense_numeric(A, b, policy);
427 policy);
429 return detail::solve_iterative_numeric(detail::dense_to_sparse(A.eval()), b, policy);
430 }
431 } else {
432 SymbolicScalar A_mx = to_mx(A);
433 SymbolicScalar b_mx = to_mx(b);
434 if (policy.symbolic_linear_solver.empty()) {
435 if (!policy.symbolic_options.empty()) {
436 throw InvalidArgument(
437 "solve: symbolic_options require a non-empty symbolic_linear_solver");
438 }
439 SymbolicScalar x_mx = SymbolicScalar::solve(A_mx, b_mx);
440 return to_eigen(x_mx);
441 }
442 SymbolicScalar x_mx = SymbolicScalar::solve(A_mx, b_mx, policy.symbolic_linear_solver,
443 policy.symbolic_options);
444 return to_eigen(x_mx);
445 }
446}
447
453template <typename DerivedB>
454auto solve(const SparseMatrix &A, const Eigen::MatrixBase<DerivedB> &b) {
456}
457
461template <typename DerivedB>
462auto solve(const SparseMatrix &A, const Eigen::MatrixBase<DerivedB> &b,
463 const LinearSolvePolicy &policy) {
464 detail::validate_linear_solve_dims(A.rows(), A.cols(), b.rows(), "solve");
465
466 switch (policy.backend) {
468 return detail::solve_dense_numeric(NumericMatrix(A), b, policy);
470 return detail::solve_sparse_direct_numeric(A, b, policy);
472 return detail::solve_iterative_numeric(A, b, policy);
473 }
474
475 return detail::solve_sparse_direct_numeric(A, b, policy);
476}
477
478// --- outer(x, y) ---
485template <typename DerivedX, typename DerivedY>
486auto outer(const Eigen::MatrixBase<DerivedX> &x, const Eigen::MatrixBase<DerivedY> &y) {
487 // Eigen's outer product works efficiently for both numeric and symbolic scalars
488 // because MX * MX (scalar mult) creates a standard multiplication node.
489 return x * y.transpose();
490}
491
492// --- Dot Product ---
499template <typename DerivedA, typename DerivedB>
500auto dot(const Eigen::MatrixBase<DerivedA> &a, const Eigen::MatrixBase<DerivedB> &b) {
501 return a.dot(b);
502}
503
504// --- Cross Product ---
512template <typename DerivedA, typename DerivedB>
513auto cross(const Eigen::MatrixBase<DerivedA> &a, const Eigen::MatrixBase<DerivedB> &b) {
514 if (a.size() != 3 || b.size() != 3) {
515 throw InvalidArgument("cross: both vectors must have exactly 3 elements");
516 }
517 // Manual implementation to support Dynamic vectors (Eigen::cross requires fixed size 3)
518 using Scalar = typename DerivedA::Scalar;
519 Eigen::Matrix<Scalar, Eigen::Dynamic, 1> res(3);
520 res(0) = a(1) * b(2) - a(2) * b(1);
521 res(1) = a(2) * b(0) - a(0) * b(2);
522 res(2) = a(0) * b(1) - a(1) * b(0);
523 return res;
524}
525
526// --- Inverse ---
532template <typename Derived> auto inv(const Eigen::MatrixBase<Derived> &A) {
533 using Scalar = typename Derived::Scalar;
534 if constexpr (std::is_floating_point_v<Scalar>) {
535 return A.inverse().eval();
536 } else {
537 SymbolicScalar A_mx = to_mx(A);
538 SymbolicScalar inv_mx = inv(A_mx);
539 return to_eigen(inv_mx);
540 }
541}
542
543// --- Determinant ---
549template <typename Derived> auto det(const Eigen::MatrixBase<Derived> &A) {
550 using Scalar = typename Derived::Scalar;
551 if constexpr (std::is_floating_point_v<Scalar>) {
552 return A.determinant();
553 } else {
554 SymbolicScalar A_mx = to_mx(A);
555 return det(A_mx);
556 }
557}
558
559// --- Inner Product ---
566template <typename DerivedA, typename DerivedB>
567auto inner(const Eigen::MatrixBase<DerivedA> &a, const Eigen::MatrixBase<DerivedB> &b) {
568 return janus::dot(a, b);
569}
570
571// --- Pseudo-Inverse ---
577template <typename Derived> auto pinv(const Eigen::MatrixBase<Derived> &A) {
578 using Scalar = typename Derived::Scalar;
579 if constexpr (std::is_floating_point_v<Scalar>) {
580 return A.completeOrthogonalDecomposition().pseudoInverse();
581 } else {
582 SymbolicScalar A_mx = to_mx(A);
583 SymbolicScalar pinv_mx = SymbolicScalar::pinv(A_mx);
584 return to_eigen(pinv_mx);
585 }
586}
587
588// --- Norm Extensions ---
589
599
606template <typename Derived>
607auto norm(const Eigen::MatrixBase<Derived> &x, NormType type = NormType::L2) {
608 using Scalar = typename Derived::Scalar;
609
610 if constexpr (std::is_floating_point_v<Scalar>) {
611 switch (type) {
612 case NormType::L1:
613 return x.template lpNorm<1>();
614 case NormType::L2:
615 return x.norm();
616 case NormType::Inf:
617 return x.template lpNorm<Eigen::Infinity>();
619 return x.norm(); // Frobenius equals L2 for vectors
620 default:
621 return x.norm();
622 }
623 } else {
624 SymbolicScalar x_mx = to_mx(x);
625 switch (type) {
626 case NormType::L1:
627 return SymbolicScalar::norm_1(x_mx);
628 case NormType::L2:
629 return SymbolicScalar::norm_2(x_mx);
630 case NormType::Inf:
631 return SymbolicScalar::norm_inf(x_mx);
633 return SymbolicScalar::norm_fro(x_mx);
634 default:
635 return SymbolicScalar::norm_2(x_mx);
636 }
637 }
638}
639
640// Backwards compatibility overload for just L2 norm (default handled above really, but if called
641// without args, it works)
642
651
652namespace detail {
653
654template <typename Scalar> JanusVector<Scalar> normalize_vector(const JanusVector<Scalar> &v) {
655 const auto v_norm = janus::norm(v);
656 if constexpr (std::is_floating_point_v<Scalar>) {
657 if (v_norm <= std::numeric_limits<Scalar>::epsilon()) {
658 throw InvalidArgument("eigendecomposition: eigenvector construction failed");
659 }
660 }
661 return v / v_norm;
662}
663
664template <typename Scalar>
666 const auto n01 = janus::dot(cands[0], cands[0]);
667 const auto n02 = janus::dot(cands[1], cands[1]);
668 const auto n12 = janus::dot(cands[2], cands[2]);
669
670 if constexpr (std::is_floating_point_v<Scalar>) {
671 const auto *best = &cands[0];
672 Scalar best_norm = n01;
673 if (n02 > best_norm) {
674 best = &cands[1];
675 best_norm = n02;
676 }
677 if (n12 > best_norm) {
678 best = &cands[2];
679 }
680 return normalize_vector(*best);
681 } else {
682 auto best = janus::detail::select(n02 > n01, cands[1], cands[0]);
683 auto best_norm = janus::where(n02 > n01, n02, n01);
684 best = janus::detail::select(n12 > best_norm, cands[2], best);
685 return normalize_vector(best);
686 }
687}
688
689template <typename Scalar>
691 JanusVector<Scalar> first(2);
692 first << A(0, 1), lambda - A(0, 0);
693
694 JanusVector<Scalar> second(2);
695 second << lambda - A(1, 1), A(1, 0);
696
697 if constexpr (std::is_floating_point_v<Scalar>) {
698 return normalize_vector(janus::dot(first, first) >= janus::dot(second, second) ? first
699 : second);
700 } else {
702 janus::dot(first, first) >= janus::dot(second, second), first, second));
703 }
704}
705
706template <typename Scalar>
708 JanusMatrix<Scalar> shifted = A;
709 shifted(0, 0) = shifted(0, 0) - lambda;
710 shifted(1, 1) = shifted(1, 1) - lambda;
711 shifted(2, 2) = shifted(2, 2) - lambda;
712
713 const JanusVector<Scalar> r0 = shifted.row(0).transpose();
714 const JanusVector<Scalar> r1 = shifted.row(1).transpose();
715 const JanusVector<Scalar> r2 = shifted.row(2).transpose();
716
718 {janus::cross(r0, r1), janus::cross(r0, r2), janus::cross(r1, r2)});
719}
720
721template <typename Scalar> Scalar determinant_3x3(const JanusMatrix<Scalar> &A) {
722 return A(0, 0) * (A(1, 1) * A(2, 2) - A(1, 2) * A(2, 1)) -
723 A(0, 1) * (A(1, 0) * A(2, 2) - A(1, 2) * A(2, 0)) +
724 A(0, 2) * (A(1, 0) * A(2, 1) - A(1, 1) * A(2, 0));
725}
726
727template <typename Scalar>
728void sort_eigenpairs(JanusVector<Scalar> &eigenvalues, JanusMatrix<Scalar> &eigenvectors) {
729 std::vector<Eigen::Index> order(static_cast<size_t>(eigenvalues.size()));
730 std::iota(order.begin(), order.end(), Eigen::Index{0});
731 std::sort(order.begin(), order.end(), [&](Eigen::Index lhs, Eigen::Index rhs) {
732 return eigenvalues(lhs) < eigenvalues(rhs);
733 });
734
735 JanusVector<Scalar> sorted_values(eigenvalues.size());
736 JanusMatrix<Scalar> sorted_vectors(eigenvectors.rows(), eigenvectors.cols());
737 for (Eigen::Index i = 0; i < eigenvalues.size(); ++i) {
738 sorted_values(i) = eigenvalues(order[static_cast<size_t>(i)]);
739 sorted_vectors.col(i) = eigenvectors.col(order[static_cast<size_t>(i)]);
740 }
741
742 eigenvalues = std::move(sorted_values);
743 eigenvectors = std::move(sorted_vectors);
744}
745
746template <typename Scalar>
748 if (A.rows() != A.cols()) {
749 throw InvalidArgument("eig_symmetric: input must be square");
750 }
751
752 if (A.rows() == 1) {
754 result.eigenvalues.resize(1);
755 result.eigenvalues(0) = A(0, 0);
757 return result;
758 }
759
760 if (A.rows() == 2) {
761 const Scalar trace = A(0, 0) + A(1, 1);
762 const Scalar disc = janus::sqrt((A(0, 0) - A(1, 1)) * (A(0, 0) - A(1, 1)) +
763 Scalar(4.0) * A(0, 1) * A(0, 1));
764
766 result.eigenvalues.resize(2);
767 result.eigenvalues(0) = Scalar(0.5) * (trace - disc);
768 result.eigenvalues(1) = Scalar(0.5) * (trace + disc);
769
770 result.eigenvectors.resize(2, 2);
771 result.eigenvectors.col(0) = symmetric_eigenvector_2x2(A, result.eigenvalues(0));
772 result.eigenvectors.col(1) = symmetric_eigenvector_2x2(A, result.eigenvalues(1));
773 return result;
774 }
775
776 if (A.rows() == 3) {
777 const Scalar q = (A(0, 0) + A(1, 1) + A(2, 2)) / Scalar(3.0);
778 const Scalar p1 = A(0, 1) * A(0, 1) + A(0, 2) * A(0, 2) + A(1, 2) * A(1, 2);
779 const Scalar a00 = A(0, 0) - q;
780 const Scalar a11 = A(1, 1) - q;
781 const Scalar a22 = A(2, 2) - q;
782 const Scalar p2 = a00 * a00 + a11 * a11 + a22 * a22 + Scalar(2.0) * p1;
783
784 if constexpr (std::is_floating_point_v<Scalar>) {
785 if (p2 <= std::numeric_limits<Scalar>::epsilon()) {
789 return result;
790 }
791 }
792
793 const auto has_spread = p2 > Scalar(0.0);
794 const Scalar p = janus::where(has_spread, janus::sqrt(p2 / Scalar(6.0)), Scalar(1.0));
795
796 JanusMatrix<Scalar> centered = A;
797 centered(0, 0) = centered(0, 0) - q;
798 centered(1, 1) = centered(1, 1) - q;
799 centered(2, 2) = centered(2, 2) - q;
800 const JanusMatrix<Scalar> B = centered / p;
801
802 const Scalar r = janus::clamp(determinant_3x3(B) / Scalar(2.0), -1.0, 1.0);
803 const Scalar phi = janus::where(has_spread, janus::acos(r) / Scalar(3.0), Scalar(0.0));
804
805 constexpr double kTwoPiOverThree = 2.0943951023931954923;
806 Scalar largest = q + Scalar(2.0) * p * janus::cos(phi);
807 Scalar smallest = q + Scalar(2.0) * p * janus::cos(phi + Scalar(kTwoPiOverThree));
808 Scalar middle = Scalar(3.0) * q - largest - smallest;
809
810 largest = janus::where(has_spread, largest, q);
811 middle = janus::where(has_spread, middle, q);
812 smallest = janus::where(has_spread, smallest, q);
813
815 result.eigenvalues.resize(3);
816 result.eigenvalues << smallest, middle, largest;
817
818 JanusMatrix<Scalar> vectors(3, 3);
819 vectors.col(0) = symmetric_eigenvector_3x3(A, smallest);
820 vectors.col(1) = symmetric_eigenvector_3x3(A, middle);
821 vectors.col(2) = symmetric_eigenvector_3x3(A, largest);
822
824 result.eigenvectors = janus::detail::select(has_spread, vectors, identity);
825 return result;
826 }
827
828 throw InvalidArgument(
829 "eig_symmetric: symbolic support is limited to 1x1, 2x2, and 3x3 matrices");
830}
831
832} // namespace detail
833
834// --- Eigendecomposition ---
842template <typename Derived> auto eig(const Eigen::MatrixBase<Derived> &A) {
843 using Scalar = typename Derived::Scalar;
844
845 if (A.rows() != A.cols()) {
846 throw InvalidArgument("eig: input must be square");
847 }
848
849 if constexpr (std::is_floating_point_v<Scalar>) {
850 using Matrix = JanusMatrix<Scalar>;
851 using Complex = std::complex<Scalar>;
852 using ComplexMatrix = Eigen::Matrix<Complex, Eigen::Dynamic, Eigen::Dynamic>;
853 using ComplexVector = Eigen::Matrix<Complex, Eigen::Dynamic, 1>;
854
855 if (A.rows() == 1) {
857 result.eigenvalues.resize(1);
858 result.eigenvalues(0) = A(0, 0);
859 result.eigenvectors = Matrix::Identity(1, 1);
860 return result;
861 }
862
863 Eigen::EigenSolver<Matrix> solver(A.eval());
864 if (solver.info() != Eigen::Success) {
865 throw InvalidArgument("eig: EigenSolver failed");
866 }
867
868 const ComplexVector values_complex = solver.eigenvalues();
869 const ComplexMatrix vectors_complex = solver.eigenvectors();
870 constexpr Scalar kImagTol = Scalar(1e-10);
871
872 for (Eigen::Index i = 0; i < values_complex.size(); ++i) {
873 if (std::abs(values_complex(i).imag()) > kImagTol) {
874 throw InvalidArgument("eig: complex eigenvalues are not supported");
875 }
876 }
877
878 for (Eigen::Index i = 0; i < vectors_complex.rows(); ++i) {
879 for (Eigen::Index j = 0; j < vectors_complex.cols(); ++j) {
880 if (std::abs(vectors_complex(i, j).imag()) > kImagTol) {
881 throw InvalidArgument("eig: complex eigenvectors are not supported");
882 }
883 }
884 }
885
887 result.eigenvalues = values_complex.real();
888 result.eigenvectors = vectors_complex.real();
890 return result;
891 } else {
892 throw InvalidArgument(
893 "eig: symbolic general eigendecomposition is not supported for CasADi MX; use "
894 "eig_symmetric() for 1x1, 2x2, or 3x3 symmetric matrices");
895 }
896}
897
905template <typename Derived> auto eig_symmetric(const Eigen::MatrixBase<Derived> &A) {
906 using Scalar = typename Derived::Scalar;
907
908 if (A.rows() != A.cols()) {
909 throw InvalidArgument("eig_symmetric: input must be square");
910 }
911
912 if constexpr (std::is_floating_point_v<Scalar>) {
913 if (!A.isApprox(A.transpose(), 1e-12)) {
914 throw InvalidArgument("eig_symmetric: numeric input must be symmetric");
915 }
916
917 using Matrix = JanusMatrix<Scalar>;
918 Eigen::SelfAdjointEigenSolver<Matrix> solver(A.eval());
919 if (solver.info() != Eigen::Success) {
920 throw InvalidArgument("eig_symmetric: SelfAdjointEigenSolver failed");
921 }
922
924 result.eigenvalues = solver.eigenvalues();
925 result.eigenvectors = solver.eigenvectors();
927 return result;
928 } else {
929 return detail::eig_symmetric_symbolic(A.eval());
930 }
931}
932
933// --- Explicit 3x3 Symmetric Inverse (AeroSandbox Helper) ---
950template <typename T>
951std::tuple<T, T, T, T, T, T> inv_symmetric_3x3_explicit(const T &m11, const T &m22, const T &m33,
952 const T &m12, const T &m23, const T &m13) {
953
954 T det = m11 * (m33 * m22 - m23 * m23) - m12 * (m33 * m12 - m23 * m13) +
955 m13 * (m23 * m12 - m22 * m13);
956
957 T inv_det = 1.0 / det;
958
959 T a11 = (m33 * m22 - m23 * m23) * inv_det;
960 T a12 = (m13 * m23 - m33 * m12) * inv_det;
961 T a13 = (m12 * m23 - m13 * m22) * inv_det;
962 T a22 = (m33 * m11 - m13 * m13) * inv_det;
963 T a23 = (m12 * m13 - m11 * m23) * inv_det;
964 T a33 = (m11 * m22 - m12 * m12) * inv_det;
965
966 return {a11, a12, a13, a22, a23, a33};
967}
968
969// --- Sparse Matrix Utilities ---
970// NOTE: Sparse matrices are for NUMERIC data only. CasADi MX handles sparsity internally.
971// For symbolic sparsity analysis, use janus::SparsityPattern.
972
979template <typename Scalar>
980constexpr bool is_numeric_scalar_v = std::is_floating_point_v<Scalar> || std::is_integral_v<Scalar>;
981
999inline SparseMatrix sparse_from_triplets(int rows, int cols,
1000 const std::vector<SparseTriplet> &triplets) {
1001 SparseMatrix m(rows, cols);
1002 m.setFromTriplets(triplets.begin(), triplets.end());
1003 return m;
1004}
1005
1015inline SparseMatrix to_sparse(const NumericMatrix &dense, double tol = 0.0) {
1016 SparseMatrix sparse = dense.sparseView(1.0, tol);
1017 sparse.makeCompressed();
1018 return sparse;
1019}
1020
1027inline NumericMatrix to_dense(const SparseMatrix &sparse) { return NumericMatrix(sparse); }
1028
1036 SparseMatrix I(n, n);
1037 I.setIdentity();
1038 return I;
1039}
1040
1041} // namespace janus
Scalar and element-wise arithmetic functions (abs, sqrt, pow, exp, log, etc.).
C++20 concepts constraining valid Janus scalar types.
Custom exception hierarchy for Janus framework.
Core type aliases for numeric and symbolic Eigen/CasADi interop.
Conditional selection, comparison, and logical operations.
Trigonometric and inverse trigonometric functions.
Input validation failed (e.g., mismatched sizes, invalid parameters).
Definition JanusError.hpp:31
Eigen::ComputationInfo info() const
Definition Linalg.hpp:199
void _solve_impl(const Rhs &b, Dest &x) const
Definition Linalg.hpp:186
void set_apply(std::function< NumericVector(const NumericVector &)> apply)
Definition Linalg.hpp:165
FunctionalPreconditioner & factorize(const MatType &mat)
Definition Linalg.hpp:176
double Scalar
Definition Linalg.hpp:158
const Eigen::Solve< FunctionalPreconditioner, Rhs > solve(const Eigen::MatrixBase< Rhs > &b) const
Definition Linalg.hpp:193
Eigen::Index rows() const noexcept
Definition Linalg.hpp:169
int StorageIndex
Definition Linalg.hpp:160
FunctionalPreconditioner & compute(const MatType &mat)
Definition Linalg.hpp:182
Eigen::Index cols() const noexcept
Definition Linalg.hpp:170
@ MaxColsAtCompileTime
Definition Linalg.hpp:161
@ ColsAtCompileTime
Definition Linalg.hpp:161
FunctionalPreconditioner & analyzePattern(const MatType &)
Definition Linalg.hpp:172
double RealScalar
Definition Linalg.hpp:159
Smooth approximation of ReLU function: softplus(x) = (1/beta) * log(1 + exp(beta * x)).
Definition Diagnostics.hpp:131
void validate_iterative_policy(const LinearSolvePolicy &policy, const std::string &context)
Definition Linalg.hpp:224
JanusVector< Scalar > symmetric_eigenvector_3x3(const JanusMatrix< Scalar > &A, const Scalar &lambda)
Definition Linalg.hpp:707
auto solve_sparse_direct_numeric(const SparseMatrix &A, const Eigen::MatrixBase< DerivedB > &b, const LinearSolvePolicy &policy)
Definition Linalg.hpp:265
void validate_linear_solve_dims(Eigen::Index a_rows, Eigen::Index a_cols, Eigen::Index b_rows, const std::string &context)
Definition Linalg.hpp:207
void validate_square_required(Eigen::Index rows, Eigen::Index cols, const std::string &context, const std::string &solver_name)
Definition Linalg.hpp:217
std::function< NumericVector(const NumericVector &)> make_preconditioner(const SparseMatrix &A, const LinearSolvePolicy &policy)
Definition Linalg.hpp:237
Scalar determinant_3x3(const JanusMatrix< Scalar > &A)
Definition Linalg.hpp:721
JanusVector< Scalar > best_eigenvector_candidate(const std::array< JanusVector< Scalar >, 3 > &cands)
Definition Linalg.hpp:665
auto solve_iterative_with_solver(Solver &solver, const Eigen::MatrixBase< DerivedB > &b)
Definition Linalg.hpp:313
EigenDecomposition< Scalar > eig_symmetric_symbolic(const JanusMatrix< Scalar > &A)
Definition Linalg.hpp:747
auto select(const Cond &cond, const Eigen::MatrixBase< DerivedTrue > &if_true, const Eigen::MatrixBase< DerivedFalse > &if_false)
Definition Logic.hpp:55
SparseMatrix dense_to_sparse(const MatrixLike &A)
Definition Linalg.hpp:258
JanusVector< Scalar > normalize_vector(const JanusVector< Scalar > &v)
Definition Linalg.hpp:654
auto solve_dense_numeric(const Eigen::MatrixBase< DerivedA > &A, const Eigen::MatrixBase< DerivedB > &b, const LinearSolvePolicy &policy)
Definition Linalg.hpp:362
auto solve_iterative_numeric(const SparseMatrix &A, const Eigen::MatrixBase< DerivedB > &b, const LinearSolvePolicy &policy)
Definition Linalg.hpp:329
JanusVector< Scalar > symmetric_eigenvector_2x2(const JanusMatrix< Scalar > &A, const Scalar &lambda)
Definition Linalg.hpp:690
void sort_eigenpairs(JanusVector< Scalar > &eigenvalues, JanusMatrix< Scalar > &eigenvectors)
Definition Linalg.hpp:728
Definition Diagnostics.hpp:19
Eigen::SparseMatrix< double > SparseMatrix
Sparse matrix types for efficient storage of large, sparse numeric data.
Definition JanusTypes.hpp:80
LinearSolveBackend
Backend selection for linear system solves.
Definition Linalg.hpp:39
@ IterativeKrylov
Iterative Krylov subspace method.
Definition Linalg.hpp:42
@ Dense
Dense matrix factorization.
Definition Linalg.hpp:40
@ SparseDirect
Sparse direct factorization.
Definition Linalg.hpp:41
auto where(const Cond &cond, const T1 &if_true, const T2 &if_false)
Select values based on condition (ternary operator) Returns: cond ? if_true : if_false Supports mixed...
Definition Logic.hpp:43
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > JanusMatrix
Dynamic-size matrix for both numeric and symbolic backends.
Definition JanusTypes.hpp:43
Solver
Available NLP solvers.
Definition OptiOptions.hpp:21
NormType
Norm type selection.
Definition Linalg.hpp:593
@ Inf
Infinity (max absolute) norm.
Definition Linalg.hpp:596
@ L2
L2 (Euclidean) norm.
Definition Linalg.hpp:595
@ L1
L1 (Manhattan) norm.
Definition Linalg.hpp:594
@ Frobenius
Frobenius norm.
Definition Linalg.hpp:597
JanusMatrix< NumericScalar > NumericMatrix
Eigen::MatrixXd equivalent.
Definition JanusTypes.hpp:66
auto dot(const Eigen::MatrixBase< DerivedA > &a, const Eigen::MatrixBase< DerivedB > &b)
Computes dot product of two vectors.
Definition Linalg.hpp:500
T sqrt(const T &x)
Computes the square root of a scalar.
Definition Arithmetic.hpp:46
std::tuple< T, T, T, T, T, T > inv_symmetric_3x3_explicit(const T &m11, const T &m22, const T &m33, const T &m12, const T &m23, const T &m13)
Explicit inverse of a symmetric 3x3 matrix.
Definition Linalg.hpp:951
auto clamp(const T &val, const TLow &low, const THigh &high)
Clamps value between low and high.
Definition Logic.hpp:236
auto det(const Eigen::MatrixBase< Derived > &A)
Computes matrix determinant.
Definition Linalg.hpp:549
T acos(const T &x)
Computes arc cosine of x.
Definition Trig.hpp:121
auto eig_symmetric(const Eigen::MatrixBase< Derived > &A)
Computes the eigendecomposition of a symmetric matrix.
Definition Linalg.hpp:905
JanusVector< NumericScalar > NumericVector
Eigen::VectorXd equivalent.
Definition JanusTypes.hpp:67
SparseMatrix sparse_from_triplets(int rows, int cols, const std::vector< SparseTriplet > &triplets)
Create sparse matrix from triplets.
Definition Linalg.hpp:999
auto inv(const Eigen::MatrixBase< Derived > &A)
Computes matrix inverse.
Definition Linalg.hpp:532
casadi::MX to_mx(const Eigen::MatrixBase< Derived > &e)
Convert Eigen matrix of MX (or numeric) to CasADi MX.
Definition JanusTypes.hpp:189
SparseMatrix sparse_identity(int n)
Create identity sparse matrix.
Definition Linalg.hpp:1035
T cos(const T &x)
Computes cosine of x.
Definition Trig.hpp:46
NumericMatrix to_dense(const SparseMatrix &sparse)
Convert sparse matrix to dense.
Definition Linalg.hpp:1027
IterativeKrylovSolver
Iterative Krylov solver algorithm.
Definition Linalg.hpp:69
@ BiCGSTAB
Biconjugate gradient stabilized.
Definition Linalg.hpp:70
@ GMRES
Generalized minimal residual.
Definition Linalg.hpp:71
auto norm(const Eigen::MatrixBase< Derived > &x, NormType type=NormType::L2)
Computes vector/matrix norm.
Definition Linalg.hpp:607
@ None
Definition AutoDiff.hpp:30
Eigen::Matrix< casadi::MX, Eigen::Dynamic, Eigen::Dynamic > to_eigen(const casadi::MX &m)
Convert CasADi MX to Eigen matrix of MX.
Definition JanusTypes.hpp:213
auto pinv(const Eigen::MatrixBase< Derived > &A)
Computes Moore-Penrose pseudo-inverse.
Definition Linalg.hpp:577
const char * solver_name(Solver solver)
Get the CasADi solver name string.
Definition OptiOptions.hpp:38
auto cross(const Eigen::MatrixBase< DerivedA > &a, const Eigen::MatrixBase< DerivedB > &b)
Computes 3D cross product.
Definition Linalg.hpp:513
DenseLinearSolver
Dense linear solver algorithm.
Definition Linalg.hpp:48
@ FullPivLU
Full-pivot LU (square only).
Definition Linalg.hpp:51
@ LLT
Cholesky (SPD only).
Definition Linalg.hpp:52
@ LDLT
LDLT (symmetric only).
Definition Linalg.hpp:53
@ PartialPivLU
Partial-pivot LU (square only).
Definition Linalg.hpp:50
@ ColPivHouseholderQR
Column-pivoted Householder QR (default, general).
Definition Linalg.hpp:49
IterativePreconditioner
Preconditioner for iterative solvers.
Definition Linalg.hpp:77
@ None
No preconditioning.
Definition Linalg.hpp:78
@ Diagonal
Jacobi (diagonal) preconditioner.
Definition Linalg.hpp:79
SparseMatrix to_sparse(const NumericMatrix &dense, double tol=0.0)
Convert dense matrix to sparse.
Definition Linalg.hpp:1015
constexpr bool is_numeric_scalar_v
Compile-time check for numeric scalar types.
Definition Linalg.hpp:980
auto inner(const Eigen::MatrixBase< DerivedA > &a, const Eigen::MatrixBase< DerivedB > &b)
Computes inner product of two vectors (dot product).
Definition Linalg.hpp:567
SparseDirectLinearSolver
Sparse direct solver algorithm.
Definition Linalg.hpp:59
@ SimplicialLDLT
Simplicial LDLT (symmetric only).
Definition Linalg.hpp:63
@ SimplicialLLT
Simplicial Cholesky (SPD only).
Definition Linalg.hpp:62
@ SparseLU
Sparse LU factorization.
Definition Linalg.hpp:60
@ SparseQR
Sparse QR factorization.
Definition Linalg.hpp:61
auto solve(const Eigen::MatrixBase< DerivedA > &A, const Eigen::MatrixBase< DerivedB > &b)
Solves linear system Ax = b using the default backend policy.
Definition Linalg.hpp:408
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > JanusVector
Dynamic-size column vector for both numeric and symbolic backends.
Definition JanusTypes.hpp:49
casadi::MX SymbolicScalar
CasADi MX symbolic scalar.
Definition JanusTypes.hpp:70
auto eig(const Eigen::MatrixBase< Derived > &A)
Computes the eigendecomposition of a square matrix with a real spectrum.
Definition Linalg.hpp:842
auto outer(const Eigen::MatrixBase< DerivedX > &x, const Eigen::MatrixBase< DerivedY > &y)
Computes outer product x * y^T.
Definition Linalg.hpp:486
Result of eigendecomposition: eigenvalues and eigenvectors.
Definition Linalg.hpp:647
JanusMatrix< Scalar > eigenvectors
Eigenvectors as columns.
Definition Linalg.hpp:649
JanusVector< Scalar > eigenvalues
Eigenvalues in ascending order.
Definition Linalg.hpp:648
Configuration for linear system solve backend and algorithm.
Definition Linalg.hpp:86
IterativeKrylovSolver iterative_solver
Definition Linalg.hpp:90
LinearSolveBackend backend
Definition Linalg.hpp:87
int gmres_restart
Definition Linalg.hpp:94
LinearSolvePolicy & set_max_iterations(int value)
Definition Linalg.hpp:130
std::function< NumericVector(const NumericVector &)> preconditioner_hook
Definition Linalg.hpp:95
static LinearSolvePolicy sparse_direct(SparseDirectLinearSolver solver=SparseDirectLinearSolver::SparseLU)
Definition Linalg.hpp:108
LinearSolvePolicy & set_symbolic_solver(const std::string &solver, const casadi::Dict &opts=casadi::Dict())
Definition Linalg.hpp:146
double tolerance
Definition Linalg.hpp:92
static LinearSolvePolicy iterative(IterativeKrylovSolver solver=IterativeKrylovSolver::BiCGSTAB, IterativePreconditioner preconditioner=IterativePreconditioner::Diagonal)
Definition Linalg.hpp:116
int max_iterations
Definition Linalg.hpp:93
SparseDirectLinearSolver sparse_direct_solver
Definition Linalg.hpp:89
LinearSolvePolicy & set_preconditioner_hook(std::function< NumericVector(const NumericVector &)> hook)
Definition Linalg.hpp:141
IterativePreconditioner iterative_preconditioner
Definition Linalg.hpp:91
std::string symbolic_linear_solver
Definition Linalg.hpp:96
DenseLinearSolver dense_solver
Definition Linalg.hpp:88
LinearSolvePolicy & set_tolerance(double value)
Definition Linalg.hpp:125
static LinearSolvePolicy dense(DenseLinearSolver solver=DenseLinearSolver::ColPivHouseholderQR)
Definition Linalg.hpp:100
LinearSolvePolicy & set_gmres_restart(int value)
Definition Linalg.hpp:135
casadi::Dict symbolic_options
Definition Linalg.hpp:97