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

Janus provides a high-level C++ interface to constrained nonlinear optimization through janus::Opti, wrapping CasADi/IPOPT. It allows you to define optimization variables, write objectives and constraints using standard C++ syntax, and reuse your physicist-written simulation models directly in the optimization loop. This works in symbolic mode, with seamless interop with numeric constants.

Quick Start

#include <janus/janus.hpp>
auto x = opti.variable(0.0);
auto y = opti.variable(0.0);
// Minimize the Rosenbrock function
auto f = janus::pow(1 - x, 2) + 100 * janus::pow(y - janus::pow(x, 2), 2);
opti.minimize(f);
auto sol = opti.solve({.max_iter = 500, .verbose = false});
std::cout << "x=" << sol.value(x) << ", y=" << sol.value(y) << std::endl;
double value(const SymbolicScalar &var) const
Extract scalar value at optimum.
Definition OptiSol.hpp:45
Main optimization environment class.
Definition Opti.hpp:167
void minimize(const SymbolicScalar &objective)
Set objective to minimize.
Definition Opti.hpp:555
OptiSol solve(const OptiOptions &options={})
Solve the optimization problem.
Definition Opti.hpp:603
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
Umbrella header that includes the entire Janus public API.
T pow(const T &base, const T &exponent)
Computes the power function: base^exponent.
Definition Arithmetic.hpp:72

Core API

  • janus::Opti: The optimization problem builder.
  • opti.variable(initial_guess): Create a scalar decision variable.
  • opti.variable(n, initial_guess): Create a vector decision variable of length n.
  • opti.parameter(value): Create a parameter that can be changed between solves.
  • opti.minimize(expr): Set the objective function to minimize.
  • opti.subject_to(constraint): Add an equality (==) or inequality (<=, >=) constraint.
  • opti.subject_to_bounds(var, lb, ub): Set variable bounds (more efficient than general constraints).
  • opti.subject_to_lower(var, lb) / opti.subject_to_upper(var, ub): One-sided bounds.
  • opti.solve(options): Execute IPOPT and return a solution object.
  • sol.value(var): Retrieve the optimal value of a variable from the solution.

Usage Patterns

The Rosenbrock Benchmark

The Rosenbrock function (or "banana function") is a classic test for optimization algorithms. Code reference: examples/optimization/rosenbrock.cpp

auto x = opti.variable(0.0);
auto y = opti.variable(0.0);
auto f = janus::pow(1 - x, 2) + 100 * janus::pow(y - janus::pow(x, 2), 2);
opti.minimize(f);
// Add constraints
opti.subject_to(janus::pow(x, 2) + janus::pow(y, 2) <= 2.0);
auto sol = opti.solve({.max_iter = 500, .verbose = false});
double x_star = sol.value(x);
double y_star = sol.value(y);
std::cout << "Optimal Solution: x=" << x_star << ", y=" << y_star << std::endl;
void subject_to(const SymbolicScalar &constraint)
Add a scalar constraint.
Definition Opti.hpp:396

"Write Once, Use Everywhere" (C++20 Style)

The true power of Janus is reusing your existing physics code. With C++20 Abbreviated Function Templates (auto parameters), this is seamless.

Code reference: examples/optimization/drag_optimization.cpp

The shared physics function works for double, SymbolicScalar, or a mix of both:

auto compute_drag(auto rho, auto v, auto S,
auto Cd0, auto k, auto Cl, auto Cl0) {
auto q = 0.5 * rho * janus::pow(v, 2.0);
auto Cd = Cd0 + k * janus::pow(Cl - Cl0, 2.0);
return q * S * Cd;
}

Using it in optimization requires no casting or explicit templates:

auto V = opti.variable(50.0);
auto Cl = opti.variable(0.5);
const double rho = 1.225;
const double S = 16.0;
auto D = compute_drag(rho, V, S, 0.02, 0.04, Cl, 0.1);
opti.minimize(D);
opti.subject_to_bounds(V, 20.0, 150.0);
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

Why this matters:

  • Zero Boilerplate: No template <typename T> syntax required for users.
  • Type Safety: Still statically checked at compile time.
  • Performance: Compiles down to optimized machine code (Numeric mode) or efficient graph construction (Symbolic mode).

Bounds Handling

Janus provides explicit helpers for variable bounds, which is more efficient than general constraints for the solver.

// Scalar variable bounds
auto x = opti.variable(0.5);
opti.subject_to_bounds(x, 0.0, 1.0); // 0 <= x <= 1
opti.subject_to_lower(x, 0.0); // x >= 0
opti.subject_to_upper(x, 1.0); // x <= 1
// Vector variable bounds (applies to all elements)
auto vec = opti.variable(10, 0.0);
opti.subject_to_bounds(vec, -5.0, 5.0);
void subject_to_upper(const SymbolicScalar &scalar, double upper_bound)
Apply upper bound to a scalar variable.
Definition Opti.hpp:481
void subject_to_lower(const SymbolicScalar &scalar, double lower_bound)
Apply lower bound to a scalar variable.
Definition Opti.hpp:472

Parametric Sweeps

Run the same optimization across a range of parameter values with automatic warm-starting:

auto rho = opti.parameter(1.225);
auto V = opti.variable(50.0);
const double S = 16.0; // Reference area [m^2]
const double Cd0 = 0.02; // Zero-lift drag coefficient
auto drag = 0.5 * rho * V * V * S * Cd0;
opti.minimize(drag);
opti.subject_to(V >= 10.0);
std::vector<double> rho_values = {1.2, 1.0, 0.8, 0.6};
auto result = opti.solve_sweep(rho, rho_values);
for (size_t i = 0; i < result.size(); ++i) {
if (result.converged[i]) {
std::cout << "rho=" << result.param_values[i]
<< " V*=" << result.solutions[i]->value(V) << "\n";
}
}
SweepResult solve_sweep(const SymbolicScalar &param, const std::vector< double > &values, const OptiOptions &options={})
Perform parametric sweep over a parameter.
Definition Opti.hpp:662
SymbolicScalar parameter(double value)
Create a scalar parameter.
Definition Opti.hpp:360

See the full example: examples/optimization/parametric_sweep.cpp

Advanced Usage

Scaling Diagnostics and Nondimensionalization

Poor scaling is a common failure mode for nonlinear programs. janus::Opti exposes three practical tools:

  • Variable scales can still be supplied explicitly through variable(..., scale, ...).
  • If the initial guess is zero, finite bounds are now used to infer a more sensible default scale.
  • Objectives and constraints can be scaled explicitly, and analyze_scaling() reports suspicious magnitudes before solve.
auto x = opti.variable(0.0, std::nullopt, -1e6, 1e6);
opti.subject_to(x == 1e6, 1e6);
opti.minimize(janus::pow(x - 1e6, 2), 1e12);
auto report = opti.analyze_scaling();
if (report.has_issues()) {
for (const auto& issue : report.issues) {
std::cout << issue.label << ": " << issue.message << "\n";
}
}
ScalingReport analyze_scaling(const ScalingAnalysisOptions &opts={}) const
Diagnose variable, objective, and constraint scaling at the current initial point.
Definition Opti.hpp:827

The report summarizes:

  • variable block scales and normalized initial guesses
  • constraint row magnitudes versus applied linear scales
  • objective magnitude versus applied objective scale
  • Jacobian sparsity density for the current NLP

This is intended as a pre-solve diagnostic pass, not a full automatic reformulation engine.

Solution Persistence & Warm Starting

Janus allows you to save optimization results to JSON and use them to warm-start subsequent runs. This is crucial for complex problems where a good initial guess can significantly reduce solve time.

Saving Results:

auto sol = opti.solve();
std::map<std::string, janus::SymbolicScalar> vars;
vars["x"] = x;
vars["y"] = y;
sol.save("solution.json", vars);
void save(const std::string &filename, const std::map< std::string, SymbolicScalar > &named_vars) const
Save solution to JSON file.
Definition OptiSol.hpp:132

Loading & Warm Starting:

try {
auto cache = janus::OptiSol::load("solution.json");
double x_init = cache.count("x") ? cache["x"][0] : 0.0;
auto x = opti.variable(x_init);
} catch (...) {
// Handle cold start (file not found, etc.)
}
static std::map< std::string, std::vector< double > > load(const std::string &filename)
Load solution data from JSON file.
Definition OptiSol.hpp:173

See Also