Janus 2.0.0
High-performance C++20 dual-mode numerical framework
Loading...
Searching...
No Matches
Root Finding

Janus provides nonlinear root-finding utilities in RootFinding.hpp for solving F(x) = 0 systems. The API serves two distinct roles: numeric nonlinear solves with a globalization stack for difficult residual systems, and symbolic implicit solves that stay differentiable inside CasADi graphs. Both modes work with dense column-vector inputs and outputs of matching dimensions.

Quick Start

#include <janus/janus.hpp>
// Define a residual symbolically
auto x = janus::sym("x");
janus::Function F("f", {x}, {x * x - 2.0});
// Solve F(x) = 0 starting from x0 = 1
Eigen::VectorXd x0(1);
x0 << 1.0;
auto result = janus::rootfinder<double>(F, x0);
std::cout << "Root: " << result.x.transpose() << std::endl; // ~1.4142
Wrapper around casadi::Function providing Eigen-native IO.
Definition Function.hpp:46
Umbrella header that includes the entire Janus public API.
SymbolicScalar sym(const std::string &name)
Create a named symbolic scalar variable.
Definition JanusTypes.hpp:90
RootResult< Scalar > rootfinder(const janus::Function &F, const Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > &x0, const RootFinderOptions &opts={})
Solve F(x) = 0 for x given an initial guess.
Definition RootFinding.hpp:849

Core API

  • janus::rootfinder<double>(F, x0, opts): Numeric root-finding with automatic globalization fallback.
  • janus::NewtonSolver(F, opts): Persistent solver for repeated solves of the same residual system.
  • janus::create_implicit_function(G, x_guess, opts): Create a differentiable implicit function x(p) from G(x, p) = 0.
  • janus::rootfinder<janus::SymbolicScalar>(): CasADi-backed symbolic solve node.
  • janus::RootFinderOptions: Configuration for tolerance, max iterations, strategy, and pseudo-transient parameters.
  • janus::RootResult<double>: Result struct with x, converged, method, iterations, residual_norm, step_norm, and message.

Usage Patterns

Numeric Solves

rootfinder<double>() and NewtonSolver::solve() use Janus's own globalization stack. The default is RootSolveStrategy::Auto, which tries:

  1. Trust-region Newton (Levenberg-Marquardt)
  2. Line-search Newton
  3. Quasi-Newton Broyden updates
  4. Pseudo-transient continuation

This is designed to be more robust than a bare Newton step when the initial Jacobian is poor, singular, or badly scaled.

Strategy Selection

@ Auto
Trust-region LM, then line-search Newton, Broyden, pseudo-transient.
Definition RootFinding.hpp:28
Options for root finding algorithms.
Definition RootFinding.hpp:54
RootSolveStrategy strategy
Definition RootFinding.hpp:62

Available numeric strategies:

Strategy Use Case
Auto Default. Best general-purpose option for trim and nonlinear residual solves.
TrustRegionNewton Good first choice when Jacobians are reliable and you want fast local convergence.
LineSearchNewton Useful when exact Newton steps are too aggressive but the Jacobian is still trustworthy.
QuasiNewtonBroyden Useful when recomputing exact Jacobians is expensive and secant updates are acceptable.
PseudoTransientContinuation Last-resort path for highly nonlinear or singular-start problems.

Result Diagnostics

Numeric solves return a RootResult<double> with method and convergence diagnostics:

auto result = janus::rootfinder<double>(F, x0);
if (result.converged) {
std::cout << result.x.transpose() << "\n";
std::cout << result.iterations << "\n";
std::cout << result.residual_norm << "\n";
}

Useful fields:

  • x: final iterate
  • converged: whether the residual tolerance was met
  • method: which stage actually converged
  • iterations: total iterations used across the stack
  • residual_norm: final infinity norm of the residual
  • step_norm: last accepted step infinity norm
  • message: short status string

Automatic Fallback

This example starts from a singular Jacobian, so Auto falls through to pseudo-transient continuation:

auto x = janus::sym("x");
janus::Function F("singular_start", {x}, {x * x - 1.0});
Eigen::VectorXd x0(1);
x0 << 0.0;
opts.max_iter = 60;
auto result = janus::rootfinder<double>(F, x0, opts);
int max_iter
Definition RootFinding.hpp:57
double pseudo_transient_dt0
Definition RootFinding.hpp:70

Persistent Solver Reuse

When you will solve the same residual system multiple times, use NewtonSolver so the residual and Jacobian kernels are compiled once:

janus::NewtonSolver solver(F, opts);
auto res1 = solver.solve(x0_a);
auto res2 = solver.solve(x0_b);
Persistent nonlinear root solver.
Definition RootFinding.hpp:677
@ LineSearchNewton
Exact-Jacobian Newton with backtracking line search.
Definition RootFinding.hpp:30

Differentiable Implicit Solves

create_implicit_function() is the right tool when the solve itself must remain inside a symbolic graph:

// G(x, p) = 0 -> x(p)
auto x = janus::sym("x");
auto p = janus::sym("p");
janus::Function G("G", {x, p}, {x * x - p});
Eigen::VectorXd guess(1);
guess << 1.0;
auto x_of_p = janus::create_implicit_function(G, guess);
auto dxdp = janus::jacobian(x_of_p(p)[0](0), p);
janus::Function create_implicit_function(const janus::Function &G, const Eigen::VectorXd &x_guess, const RootFinderOptions &opts={}, const ImplicitFunctionOptions &implicit_opts={})
Create a differentiable implicit solve wrapper for G(...) = 0.
Definition RootFinding.hpp:870
auto jacobian(const Expr &expression, const Vars &...variables)
Computes Jacobian of an expression with respect to variables.
Definition AutoDiff.hpp:109

Janus keeps this path on CasADi's differentiable newton rootfinder so exact implicit sensitivities are preserved.

Non-Default Implicit Slots

If the unknown state is not the first input or the residual is not the first output, use ImplicitFunctionOptions:

implicit_opts.implicit_input_index = 1;
implicit_opts.implicit_output_index = 1;
auto implicit = janus::create_implicit_function(G, guess, {}, implicit_opts);
Options for building differentiable implicit solve wrappers.
Definition RootFinding.hpp:78
int implicit_output_index
Definition RootFinding.hpp:80
int implicit_input_index
Definition RootFinding.hpp:79

Symbolic rootfinder

rootfinder<janus::SymbolicScalar>() also remains CasADi-rootfinder-backed. Use it when you want a symbolic solve node directly, but for optimization workflows the higher-level create_implicit_function() wrapper is usually cleaner.

See Also