Janus 2.0.0
High-performance C++20 dual-mode numerical framework
Loading...
Searching...
No Matches
Janus v2.0.0 Usage Guide

Purpose: Map of the Janus library for AI agents and developers. Each section gives a brief overview and links to the detailed user guide. For the full API, see the Doxygen docs.


Table of Contents

  1. Best Practices
  2. Type System
  3. Module Reference
    • Core Layer
    • Math Layer
    • Optimization Layer
  4. User Guide Index

Best Practices

1. Template-First Design (MANDATORY)

All physics/math functions must be templated on Scalar:

// Correct
template <typename Scalar>
Scalar my_function(const Scalar& x) { ... }
// Wrong -- breaks symbolic mode
double my_function(double x) { ... }

2. Math Dispatch – Use janus:: Namespace

Always use Janus math functions instead of std:::

// Never: std::sin(x), std::pow(x, 2), etc.
T sqrt(const T &x)
Computes the square root of a scalar.
Definition Arithmetic.hpp:46
T sin(const T &x)
Computes sine of x.
Definition Trig.hpp:21
T pow(const T &base, const T &exponent)
Computes the power function: base^exponent.
Definition Arithmetic.hpp:72
T exp(const T &x)
Computes the exponential function e^x.
Definition Arithmetic.hpp:131

3. Branching – Use janus::where(), Never if/else

Scalar result = janus::where(x > 0, x, -x);
// Multi-way branching
Scalar cd = janus::select(
{mach < 0.3, mach < 0.8, mach < 1.2},
{Scalar(0.02), Scalar(0.025), Scalar(0.05)},
Scalar(0.03)); // default
Scalar select(const std::vector< CondType > &conditions, const std::vector< Scalar > &values, const Scalar &default_value)
Multi-way conditional selection (cleaner alternative to nested where).
Definition Logic.hpp:547
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

4. Loops – Structural Bounds Only

// Correct -- structural bound (known at trace time)
for (int i = 0; i < N; ++i) { ... }
// Wrong -- dynamic bound breaks symbolic mode
while (error > tolerance) { ... }

5. Type Aliases – Use Janus Native Types

janus::Vec3<Scalar> // 3D vector
janus::Mat3<Scalar> // 3x3 matrix
janus::VecX<Scalar> // Dynamic vector
janus::MatX<Scalar> // Dynamic matrix
Core type aliases for numeric and symbolic Eigen/CasADi interop.
Eigen::Matrix< Scalar, 3, 3 > Mat3
Definition JanusTypes.hpp:61
Eigen::Matrix< Scalar, 3, 1 > Vec3
Definition JanusTypes.hpp:57

6. Include Convention

#include <janus/janus.hpp> // Everything (recommended for applications)
#include <janus/using.hpp> // Convenience header -- brings common symbols into scope
Umbrella header that includes the entire Janus public API.
Convenience header bringing common Janus symbols into scope.

Type System

Backend Types

Type Numeric Mode Symbolic Mode
Scalar double casadi::MX
Matrix Eigen::MatrixXd Eigen::Matrix<casadi::MX>
Vector Eigen::VectorXd Eigen::Matrix<casadi::MX, Dynamic, 1>

Janus Type Aliases (janus/core/JanusTypes.hpp)

// Symbolic types (for graph building)
janus::SymbolicScalar // casadi::MX
janus::SymbolicMatrix // Eigen::Matrix<casadi::MX, Dynamic, Dynamic>
janus::SymbolicVector // Eigen::Matrix<casadi::MX, Dynamic, 1>
// Numeric types (for evaluation)
janus::NumericMatrix // Eigen::MatrixXd
janus::NumericVector // Eigen::VectorXd
// Fixed-size templated types
janus::VecX<T>, janus::MatX<T>, janus::RowVecX<T>
// Sparse types (numeric only)
janus::SparseMatrix // Eigen::SparseMatrix<double>
janus::SparseTriplet // Eigen::Triplet<double>
Eigen::SparseMatrix< double > SparseMatrix
Sparse matrix types for efficient storage of large, sparse numeric data.
Definition JanusTypes.hpp:80
JanusVector< SymbolicScalar > SymbolicVector
Eigen vector of MX elements.
Definition JanusTypes.hpp:72
Eigen::Matrix< Scalar, 2, 1 > Vec2
Fixed-size vectors and matrices for performance-critical code.
Definition JanusTypes.hpp:56
JanusMatrix< NumericScalar > NumericMatrix
Eigen::MatrixXd equivalent.
Definition JanusTypes.hpp:66
Eigen::Matrix< Scalar, 4, 4 > Mat4
Definition JanusTypes.hpp:62
Eigen::Triplet< double > SparseTriplet
(row, col, value) triplet
Definition JanusTypes.hpp:81
JanusVector< NumericScalar > NumericVector
Eigen::VectorXd equivalent.
Definition JanusTypes.hpp:67
Eigen::Matrix< Scalar, 2, 2 > Mat2
Definition JanusTypes.hpp:60
Eigen::Matrix< Scalar, 4, 1 > Vec4
Definition JanusTypes.hpp:58
JanusMatrix< SymbolicScalar > SymbolicMatrix
Eigen matrix of MX elements.
Definition JanusTypes.hpp:71
casadi::MX SymbolicScalar
CasADi MX symbolic scalar.
Definition JanusTypes.hpp:70

Symbolic Variable Creation

auto x = janus::sym("x"); // Scalar
auto M = janus::sym("M", rows, cols); // Matrix (MX)
auto v = janus::sym_vector("v", size); // SymbolicVector (Eigen)
auto [vec, mx] = janus::sym_vec_pair("state", 3); // Both representations
std::pair< SymbolicVector, SymbolicScalar > sym_vec_pair(const std::string &name, int size)
Create symbolic vector and return both SymbolicVector and underlying MX.
Definition JanusTypes.hpp:155
SymbolicScalar sym(const std::string &name)
Create a named symbolic scalar variable.
Definition JanusTypes.hpp:90
SymbolicVector sym_vector(const std::string &name, int size)
Create a named symbolic vector (returns SymbolicVector).
Definition JanusTypes.hpp:115

Conversion Utilities

janus::to_mx(eigen_matrix) // Eigen -> CasADi MX
janus::to_eigen(casadi_mx) // CasADi MX -> Eigen
janus::as_mx(symbolic_vector) // SymbolicVector -> single MX
janus::as_vector(casadi_mx) // MX -> SymbolicVector
SymbolicVector as_vector(const casadi::MX &m)
Convert CasADi MX vector to SymbolicVector (Eigen container of MX).
Definition JanusTypes.hpp:228
casadi::MX to_mx(const Eigen::MatrixBase< Derived > &e)
Convert Eigen matrix of MX (or numeric) to CasADi MX.
Definition JanusTypes.hpp:189
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
SymbolicScalar as_mx(const SymbolicVector &v)
Get the underlying MX representation of a SymbolicVector.
Definition JanusTypes.hpp:173

Module Reference

Core Layer

File Description
JanusTypes.hpp Type system, aliases, janus::sym(), janus::to_mx(), janus::to_eigen()
JanusConcepts.hpp C++20 concepts: ScalarType, NumericScalar, SymbolicScalar
JanusError.hpp Exception hierarchy: JanusError, InvalidArgument, RuntimeError, IntegrationError, InterpolationError
JanusIO.hpp janus::eval(), janus::print(), janus::to_dot(), janus::graphviz()
Function.hpp janus::Function – compiled symbolic function wrapper
Sparsity.hpp Sparsity patterns, graph coloring, janus::sparse_jacobian(), janus::sparse_hessian()
Diagnostics.hpp Structural observability/identifiability: janus::analyze_structural_observability(), janus::analyze_structural_identifiability()
StructuralTransforms.hpp janus::alias_eliminate(), janus::block_triangularize(), janus::structural_analyze()

See docs/user_guides/sparsity.md, docs/user_guides/structural_diagnostics.md, docs/user_guides/structural_transforms.md.


Math Layer

Arithmetic & Trigonometry

Standard math dispatch (janus::sin, janus::pow, janus::exp, etc.) with scalar and matrix overloads. See docs/user_guides/math_functions.md for the full function table.

Logic & Branching

janus::where(), janus::select(), janus::min(), janus::max(), janus::clamp(), element-wise comparisons, smooth blending (janus::sigmoid_blend, janus::blend). See docs/user_guides/math_functions.md.

Calculus & Autodiff

Gradient, Jacobian, Hessian, Hessian-vector products, Lagrangian second-order adjoints, sensitivity regime selection. See docs/user_guides/symbolic_computing.md.

auto J = janus::jacobian(f, x);
auto H = janus::hessian(f, x);
auto Hv = janus::hessian_vector_product(f, x, direction);
auto jacobian(const Expr &expression, const Vars &...variables)
Computes Jacobian of an expression with respect to variables.
Definition AutoDiff.hpp:109
SymbolicMatrix hessian_vector_product(const SymbolicArg &expr, const SymbolicArg &vars, const SymbolicArg &direction)
Hessian-vector product for a scalar expression without forming the dense Hessian.
Definition AutoDiff.hpp:477
SymbolicMatrix hessian(const SymbolicArg &expr, const SymbolicArg &vars)
Hessian matrix (second-order derivatives).
Definition AutoDiff.hpp:437

Linear Algebra

janus::solve(A, b) with optional LinearSolvePolicy for backend selection:

auto x = janus::solve(A, b, policy);
@ SparseDirect
Sparse direct factorization.
Definition Linalg.hpp:41
@ SparseLU
Sparse LU factorization.
Definition Linalg.hpp:60
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
Configuration for linear system solve backend and algorithm.
Definition Linalg.hpp:86
LinearSolveBackend backend
Definition Linalg.hpp:87
SparseDirectLinearSolver sparse_direct_solver
Definition Linalg.hpp:89

Available backends: Dense (ColPivHouseholderQR, PartialPivLU, FullPivLU, LLT, LDLT), SparseDirect (SparseLU, SparseQR, SimplicialLLT, SimplicialLDLT), IterativeKrylov (BiCGSTAB, GMRES with preconditioner hooks). Also includes janus::dot, janus::cross, janus::norm, janus::inv, janus::det, janus::eye, janus::zeros, janus::ones, janus::block_diag, and more.

See docs/user_guides/math_functions.md.

Interpolation

1D and N-dimensional interpolation with "linear", "cubic", and "monotonic" methods. See docs/user_guides/interpolation.md.

auto y = janus::interp1(x, xp, fp);
auto y = janus::interp_nd(point, table);

Polynomial Chaos Expansions (PCE)

Askey-scheme orthogonal polynomial bases (Hermite, Legendre, Jacobi, Laguerre), total-order and tensor-product truncation, coefficient fitting via projection or regression, symbolic mean/variance extraction.

janus::PolynomialChaosBasis basis(dimensions, order, options);
auto coeffs = janus::pce_projection_coefficients(basis, grid, values);
auto mu = janus::pce_mean(coeffs);
auto var = janus::pce_variance(basis, coeffs);
Multidimensional polynomial chaos basis with fixed truncation/order.
Definition PolynomialChaos.hpp:456
JanusVector< Scalar > pce_projection_coefficients(const PolynomialChaosBasis &basis, const NumericMatrix &samples, const NumericVector &weights, const JanusVector< Scalar > &sample_values)
Compute PCE projection coefficients from weighted samples (vector).
Definition PolynomialChaos.hpp:554
Scalar pce_mean(const JanusVector< Scalar > &coefficients)
Extract PCE mean (zeroth coefficient).
Definition PolynomialChaos.hpp:665
Scalar pce_variance(const PolynomialChaosBasis &basis, const JanusVector< Scalar > &coefficients)
Compute PCE variance from coefficients.
Definition PolynomialChaos.hpp:680

See docs/user_guides/polynomial_chaos.md.

Stochastic Quadrature

Probability-measure quadrature rules, tensor-product grids, and Smolyak sparse grids for high-dimensional integration and PCE projection.

auto rule = janus::stochastic_quadrature_rule(dim, order, family);
auto sparse = janus::smolyak_sparse_grid(dimensions, level, options);
StochasticQuadratureGrid tensor_product_quadrature(const std::vector< UnivariateQuadratureRule > &rules)
Build the tensor-product grid from a list of one-dimensional rules.
Definition Quadrature.hpp:493
UnivariateQuadratureRule stochastic_quadrature_rule(const PolynomialChaosDimension &dimension, int order, StochasticQuadratureRule rule=StochasticQuadratureRule::Gauss)
Build a one-dimensional stochastic quadrature rule with a fixed order.
Definition Quadrature.hpp:426
StochasticQuadratureGrid smolyak_sparse_grid(const std::vector< PolynomialChaosDimension > &dimensions, int level, SmolyakQuadratureOptions options={})
Build a Smolyak sparse grid on probability measures.
Definition Quadrature.hpp:554

See docs/user_guides/stochastic_quadrature.md.

Root Finding

Nonlinear solve for F(x) = 0 with a numeric globalization stack (trust-region Newton, line-search Newton, Broyden, pseudo-transient continuation) and differentiable implicit function wrappers for embedding solves inside symbolic graphs.

auto result = janus::rootfinder(function, x0, opts);
// Differentiable implicit solve for use inside optimization
auto implicit_fn = janus::create_implicit_function(function, x_guess, opts, implicit_opts);
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
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

See docs/user_guides/root_finding.md.

ODE Integration

IVP solvers with multiple steppers (RK4, CVODES, BDF1, RosenbrockEuler), definite integration via Gauss-Kronrod quadrature. See docs/user_guides/integration.md.

auto result = janus::solve_ivp(dynamics, x0, t_span);
double I = janus::quad(f, a, b);
QuadResult< T > quad(Func &&func, T a, T b, double abstol=1e-8, double reltol=1e-6)
Compute definite integral using adaptive quadrature (numeric) or CVODES (symbolic).
Definition Integrate.hpp:354
OdeResult< double > solve_ivp(Func &&fun, std::pair< double, double > t_span, const NumericVector &y0, int n_eval=100, double abstol=1e-8, double reltol=1e-6)
Solve initial value problem: dy/dt = f(t, y), y(t0) = y0.
Definition Integrate.hpp:481

Second-Order Integrators

Dedicated solvers for systems of the form q'' = a(t, q):

auto result = janus::solve_second_order_ivp(accel, q0, v0, t_span);
// Single steps:
janus::stormer_verlet_step(accel, q, v, t, dt); // Symplectic
janus::rkn4_step(accel, q, v, t, dt); // 4th-order RKN
SecondOrderStepResult< Scalar > rkn4_step(AccelFunc &&acceleration, const JanusVector< Scalar > &q, const JanusVector< Scalar > &v, Scalar t, Scalar dt)
Classical 4th-order Runge-Kutta-Nystrom step for q'' = a(t, q).
Definition IntegratorStep.hpp:193
SecondOrderStepResult< Scalar > stormer_verlet_step(AccelFunc &&acceleration, const JanusVector< Scalar > &q, const JanusVector< Scalar > &v, Scalar t, Scalar dt)
Stormer-Verlet / velocity-Verlet step for q'' = a(t, q).
Definition IntegratorStep.hpp:167
SecondOrderOdeResult< double > solve_second_order_ivp(AccelFunc &&acceleration, std::pair< double, double > t_span, const NumericVector &q0, const NumericVector &v0, int n_eval=100, const SecondOrderIvpOptions &opts={})
Solve a second-order IVP q'' = a(t, q), q(t0) = q0, v(t0) = v0.
Definition Integrate.hpp:558

Mass-Matrix Integrators

Native support for stiff systems M(t,y) y' = f(t,y):

auto result = janus::solve_ivp_mass_matrix(rhs, M, x0, t_span); // Numeric
auto result = janus::solve_ivp_mass_matrix_expr(rhs, M, t, y, x0, t_span); // Symbolic (IDAS)
OdeResult< double > solve_ivp_mass_matrix(RhsFunc &&rhs, MassFunc &&mass_matrix, std::pair< double, double > t_span, const NumericVector &y0, int n_eval=100, const MassMatrixIvpOptions &opts={})
Solve M(t, y) y' = f(t, y) with a native numeric stiff integrator.
Definition Integrate.hpp:638
OdeResult< double > solve_ivp_mass_matrix_expr(const SymbolicScalar &rhs_expr, const SymbolicScalar &mass_expr, const SymbolicScalar &t_var, const SymbolicScalar &y_var, std::pair< double, double > t_span, const NumericVector &y0, int n_eval=100, const MassMatrixIvpOptions &opts={})
Solve a symbolic mass-matrix IVP using CasADi IDAS.
Definition Integrate.hpp:718

Sparsity Pipelines

NaN-propagation sparsity detection, graph coloring, and compiled sparse derivative kernels that avoid materializing dense Jacobian/Hessian matrices.

auto J = janus::sparse_jacobian(result, x);
auto H = janus::sparse_hessian(objective, vars);
auto nz = J.values(x_val); // Evaluate only nonzero entries
SparseHessianEvaluator sparse_hessian(const SymbolicArg &expression, const SymbolicArg &variables, const std::string &name="")
Compile a sparse Hessian evaluator from symbolic expressions.
Definition Sparsity.hpp:1044
SparseJacobianEvaluator sparse_jacobian(const SymbolicArg &expression, const SymbolicArg &variables, const std::string &name="")
Compile a sparse Jacobian evaluator from symbolic expressions.
Definition Sparsity.hpp:1005

See docs/user_guides/sparsity.md.

Structural Diagnostics

Preflight checks for structural observability and identifiability before committing to an optimization solve.

auto obs = janus::analyze_structural_observability(measurement_fn, 0);
auto id = janus::analyze_structural_identifiability(measurement_fn, 1);
auto all = janus::analyze_structural_diagnostics(system_fn, options);
StructuralDiagnosticsReport analyze_structural_diagnostics(const Function &fn, const StructuralDiagnosticsOptions &opts)
Run structural observability and identifiability checks together.
Definition Diagnostics.hpp:575
StructuralSensitivityReport analyze_structural_identifiability(const Function &fn, int parameter_input_idx, const StructuralSensitivityOptions &opts={})
Analyze which parameters are structurally identifiable from selected outputs.
Definition Diagnostics.hpp:561
StructuralSensitivityReport analyze_structural_observability(const Function &fn, int state_input_idx=0, const StructuralSensitivityOptions &opts={})
Analyze which states are structurally observable from selected outputs.
Definition Diagnostics.hpp:547

See docs/user_guides/structural_diagnostics.md.

Structural Transforms

Alias elimination, BLT (block lower-triangular) decomposition, and tearing recommendations for large-scale equation systems.

auto alias = janus::alias_eliminate(residual_fn);
auto blt = janus::block_triangularize(residual_fn);
auto analysis = janus::structural_analyze(residual_fn);
BLTDecomposition block_triangularize(const Function &fn, const StructuralTransformOptions &opts={})
Compute a block-triangular decomposition of a selected residual block.
Definition StructuralTransforms.hpp:550
AliasEliminationResult alias_eliminate(const Function &fn, const StructuralTransformOptions &opts={})
Eliminate trivial affine alias rows from a selected residual block.
Definition StructuralTransforms.hpp:402
StructuralAnalysis structural_analyze(const Function &fn, const StructuralTransformOptions &opts={})
Run alias elimination followed by BLT decomposition.
Definition StructuralTransforms.hpp:596

See docs/user_guides/structural_transforms.md.

Other Math Modules


Optimization Layer

Opti Interface (Opti.hpp)

// Variables
auto x = opti.variable(1.0); // Scalar
auto v = opti.variable(3, 0.0); // Vector
auto x = opti.variable(1.0, {.category = "Wing", .freeze = true});
// Parameters (fixed between solves)
auto p = opti.parameter(5.0);
// Objective
opti.minimize(cost_function);
opti.minimize(cost_function, 1e6); // Explicit objective scaling
opti.maximize(profit_function);
// Constraints
opti.subject_to(x >= 0);
opti.subject_to(g == 0, 1e3); // Explicit constraint scaling
opti.subject_to(x * x + y * y <= 1);
opti.subject_to_bounds(x, lower, upper); // Box constraints
Main optimization environment class.
Definition Opti.hpp:167
void minimize(const SymbolicScalar &objective)
Set objective to minimize.
Definition Opti.hpp:555
void subject_to_bounds(const SymbolicScalar &scalar, double lower_bound, double upper_bound)
Apply both lower and upper bounds to a scalar variable.
Definition Opti.hpp:491
void maximize(const SymbolicScalar &objective)
Set objective to maximize.
Definition Opti.hpp:575
SymbolicScalar variable(double init_guess=0.0, std::optional< double > scale=std::nullopt, std::optional< double > lower_bound=std::nullopt, std::optional< double > upper_bound=std::nullopt)
Create a scalar decision variable.
Definition Opti.hpp:200
void subject_to(const SymbolicScalar &constraint)
Add a scalar constraint.
Definition Opti.hpp:396
SymbolicScalar parameter(double value)
Create a scalar parameter.
Definition Opti.hpp:360

See docs/user_guides/optimization.md.

Solving and Options (OptiOptions.hpp)

Options are passed to solve() via OptiOptions:

auto sol = opti.solve(); // Defaults
auto sol = opti.solve({.max_iter = 500, .verbose = false}); // Designated initializers
auto sol = opti.solve(janus::OptiOptions{}.set_tol(1e-10)); // Builder pattern
// Solver selection
auto sol = opti.solve({.solver = janus::Solver::Ipopt}); // Default
auto sol = opti.solve({.solver = janus::Solver::Snopt}); // Requires SNOPT license
OptiSol solve(const OptiOptions &options={})
Solve the optimization problem.
Definition Opti.hpp:603
@ Ipopt
Interior Point OPTimizer (default, always available).
Definition OptiOptions.hpp:22
@ Snopt
Sparse Nonlinear OPTimizer (requires license).
Definition OptiOptions.hpp:23
bool solver_available(Solver solver)
Check if a solver is available in the current CasADi build.
Definition OptiOptions.hpp:58
Options for solving optimization problems.
Definition OptiOptions.hpp:126
OptiOptions & set_tol(double v)
Set convergence tolerance.
Definition OptiOptions.hpp:157

Solution Extraction (OptiSol.hpp)

auto sol = opti.solve();
double x_opt = sol.value(x); // Scalar
janus::NumericVector v_opt = sol.value(v); // Vector
janus::NumericMatrix M_opt = sol.value(M); // Matrix
auto stats = sol.stats(); // Solver statistics
// Save / load
sol.save("result.json", {{"x", x}, {"y", y}});
auto data = janus::OptiSol::load("result.json");
static std::map< std::string, std::vector< double > > load(const std::string &filename)
Load solution data from JSON file.
Definition OptiSol.hpp:173
double value(const SymbolicScalar &var) const
Extract scalar value at optimum.
Definition OptiSol.hpp:45

Scaling Diagnostics

Preflight analysis of variable, constraint, and objective scaling before solving:

auto report = opti.analyze_scaling();
// Report contains ScalingIssue entries with severity, suggested scales, etc.
ScalingReport analyze_scaling(const ScalingAnalysisOptions &opts={}) const
Diagnose variable, objective, and constraint scaling at the current initial point.
Definition Opti.hpp:827

Parametric Sweep (OptiSweep.hpp)

janus::OptiSweep sweep(opti);
auto results = sweep.run(parameter, values);

Trajectory Optimization

Four transcription methods are available, all sharing a common TranscriptionBase interface. See docs/user_guides/transcription_methods.md for comparison.

auto [X, U, tau] = colloc.setup(n_states, n_controls, t0, tf);
colloc.set_dynamics([](const auto& x, const auto& u, const auto& t) {
return dynamics(x, u, t);
});
colloc.add_dynamics_constraints();
colloc.set_initial_state(x0);
colloc.set_final_state(xf);
opti.minimize(objective);
auto sol = opti.solve();
Direct collocation transcription.
Definition Collocation.hpp:38

User Guide Index

Guide File Topics
Numeric Computing docs/user_guides/numeric_computing.md Numeric mode, evaluation
Symbolic Computing docs/user_guides/symbolic_computing.md Symbolic mode, graph building
Math Functions docs/user_guides/math_functions.md Full janus:: math dispatch table
Interpolation docs/user_guides/interpolation.md 1D/ND interpolation
Integration docs/user_guides/integration.md ODE solvers, quadrature
Root Finding docs/user_guides/root_finding.md Nonlinear solves, implicit functions
Polynomial Chaos docs/user_guides/polynomial_chaos.md PCE bases, fitting, moments
Stochastic Quadrature docs/user_guides/stochastic_quadrature.md Quadrature rules, sparse grids
Sparsity docs/user_guides/sparsity.md Sparsity, coloring, sparse derivatives
Structural Diagnostics docs/user_guides/structural_diagnostics.md Observability, identifiability
Structural Transforms docs/user_guides/structural_transforms.md Alias elimination, BLT, tearing
Graph Visualization docs/user_guides/graph_visualization.md DOT export, Graphviz
Optimization docs/user_guides/optimization.md Opti interface, constraints
Collocation docs/user_guides/collocation.md Direct collocation
Multiple Shooting docs/user_guides/multiple_shooting.md Multiple shooting
Pseudospectral docs/user_guides/pseudospectral.md LGL/CGL pseudospectral
Birkhoff Pseudospectral docs/user_guides/birkhoff_pseudospectral.md LGL/CGL Birkhoff
Transcription Methods docs/user_guides/transcription_methods.md Comparison of methods

Summary for Agents

DO NOT Reimplement

The following functionality already exists in Janus – check the relevant user guide before building anything new:

  • All basic math (sin, cos, pow, exp, log, sqrt, etc.)
  • Linear algebra (dot, cross, norm, inv, det, solve with policies)
  • Quaternion algebra and rotations
  • Interpolation (1D, ND, multiple methods)
  • Root finding (globalization stack, implicit function wrappers)
  • ODE integration (solve_ivp, quad, second-order, mass-matrix)
  • Polynomial chaos and stochastic quadrature
  • Sparsity analysis and sparse derivative kernels
  • Structural diagnostics and transforms
  • Discrete integration (trapz, Simpson, etc.)
  • Branching logic (where, select, clamp, min, max)
  • Function compilation
  • Optimization (Opti, IPOPT, SNOPT)
  • Trajectory optimization (collocation, multiple shooting, pseudospectral, Birkhoff)
  • Scaling diagnostics

When Building on Janus

  1. Import via #include <janus/janus.hpp> (includes everything)
  2. Use Janus types (janus::Vec3<Scalar>, janus::SymbolicScalar)
  3. Use Janus math (janus::sin, not std::sin)
  4. Use Janus branching (janus::where, not if/else)
  5. Template everything on Scalar
  6. Test both modes (numeric AND symbolic)