Janus 2.0.0
High-performance C++20 dual-mode numerical framework
Loading...
Searching...
No Matches
janus::Opti Class Reference

Main optimization environment class. More...

#include <Opti.hpp>

Public Member Functions

 Opti ()
 Construct a new optimization environment.
 Opti (const std::vector< std::string > &categories_to_freeze)
 Construct with frozen categories.
 ~Opti ()=default
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.
SymbolicScalar variable (double init_guess, const VariableOptions &opts, 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 with options.
SymbolicVector variable (int n_vars, 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 vector of decision variables with scalar init guess.
SymbolicVector variable (const NumericVector &init_guess, std::optional< double > scale=std::nullopt, std::optional< double > lower_bound=std::nullopt, std::optional< double > upper_bound=std::nullopt)
 Create a vector of decision variables with per-element init guess.
SymbolicScalar parameter (double value)
 Create a scalar parameter.
SymbolicVector parameter (const NumericVector &value)
 Create a vector parameter.
void subject_to (const SymbolicScalar &constraint)
 Add a scalar constraint.
void subject_to (const SymbolicScalar &constraint, double linear_scale)
 Add a scaled scalar constraint.
void subject_to (const std::vector< SymbolicScalar > &constraints)
 Add multiple constraints.
void subject_to (const std::vector< SymbolicScalar > &constraints, double linear_scale)
 Add multiple constraints with one shared linear scale.
void subject_to (std::initializer_list< SymbolicScalar > constraints)
 Add constraints from initializer list.
void subject_to (std::initializer_list< SymbolicScalar > constraints, double linear_scale)
 Add an initializer-list of constraints with one shared linear scale.
void subject_to_lower (const SymbolicScalar &scalar, double lower_bound)
 Apply lower bound to a scalar variable.
void subject_to_upper (const SymbolicScalar &scalar, double upper_bound)
 Apply upper bound to a scalar variable.
void subject_to_bounds (const SymbolicScalar &scalar, double lower_bound, double upper_bound)
 Apply both lower and upper bounds to a scalar variable.
void subject_to_lower (const SymbolicVector &vec, double lower_bound)
 Apply lower bound to all elements of a vector.
void subject_to_upper (const SymbolicVector &vec, double upper_bound)
 Apply upper bound to all elements of a vector.
void subject_to_bounds (const SymbolicVector &vec, double lower_bound, double upper_bound)
 Apply both lower and upper bounds to all elements.
void minimize (const SymbolicScalar &objective)
 Set objective to minimize.
void minimize (const SymbolicScalar &objective, double objective_scale)
 Set objective to minimize with explicit objective scaling.
void maximize (const SymbolicScalar &objective)
 Set objective to maximize.
void maximize (const SymbolicScalar &objective, double objective_scale)
 Set objective to maximize with explicit objective scaling.
OptiSol solve (const OptiOptions &options={})
 Solve the optimization problem.
SweepResult solve_sweep (const SymbolicScalar &param, const std::vector< double > &values, const OptiOptions &options={})
 Perform parametric sweep over a parameter.
SymbolicVector derivative_of (const SymbolicVector &var, const NumericVector &with_respect_to, double derivative_init_guess, const std::string &method="trapezoidal")
 Create a derivative variable constrained by integration.
void constrain_derivative (const SymbolicVector &derivative, const SymbolicVector &var, const NumericVector &with_respect_to, const std::string &method="trapezoidal")
 Constrain an existing variable to be a derivative.
void constrain_derivative (const SymbolicVector &derivative, const SymbolicVector &var, const NumericVector &with_respect_to, IntegrationMethod method)
 Constrain an existing variable to be a derivative (enum version).
ScalingReport analyze_scaling (const ScalingAnalysisOptions &opts={}) const
 Diagnose variable, objective, and constraint scaling at the current initial point.
const casadi::Opti & casadi_opti () const
 Access the underlying CasADi Opti object.
void set_initial (const casadi::MX &x, const casadi::DM &v)
 Set initial guess for a symbolic expression.
std::vector< SymbolicScalarget_category (const std::string &category) const
 Get all variables in a category.
std::vector< std::string > get_category_names () const
 Get all registered category names.
bool is_category_frozen (const std::string &category) const
 Check if a category is marked for freezing.

Detailed Description

Main optimization environment class.

Wraps CasADi's Opti interface to provide Janus-native types and a clean C++ API for nonlinear programming with IPOPT backend.

Example:

auto x = opti.variable(0.0); // scalar, init_guess=0
auto y = opti.variable(0.0);
opti.minimize((1 - x) * (1 - x) + 100 * (y - x * x) * (y - x * x));
auto sol = opti.solve();
double x_opt = sol.value(x); // ~1.0
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

Rosenbrock Benchmark:

auto x = opti.variable(10, 0.0); // 10 variables
opti.subject_to(x >= 0);
// minimize sum(100*(x[i+1] - x[i]^2)^2 + (1 - x[i])^2)
for (int i = 0; i < 9; ++i) {
obj = obj + 100 * janus::pow(x(i+1) - x(i)*x(i), 2)
+ janus::pow(1 - x(i), 2);
}
opti.minimize(obj);
auto sol = opti.solve(); // All x[i] ~= 1.0
void subject_to(const SymbolicScalar &constraint)
Add a scalar constraint.
Definition Opti.hpp:396
T pow(const T &base, const T &exponent)
Computes the power function: base^exponent.
Definition Arithmetic.hpp:72
casadi::MX SymbolicScalar
CasADi MX symbolic scalar.
Definition JanusTypes.hpp:70
See also
OptiSol for extracting results
OptiOptions for solver configuration
SweepResult for parametric sweep results
ScalingReport for scaling diagnostics

Constructor & Destructor Documentation

◆ Opti() [1/2]

janus::Opti::Opti ( )
inline

Construct a new optimization environment.

◆ Opti() [2/2]

janus::Opti::Opti ( const std::vector< std::string > & categories_to_freeze)
inlineexplicit

Construct with frozen categories.

Variables created with a category in categories_to_freeze will be automatically frozen at their initial guess values.

Parameters
categories_to_freezeList of category names to freeze

◆ ~Opti()

janus::Opti::~Opti ( )
default

Member Function Documentation

◆ analyze_scaling()

ScalingReport janus::Opti::analyze_scaling ( const ScalingAnalysisOptions & opts = {}) const
inline

Diagnose variable, objective, and constraint scaling at the current initial point.

Parameters
optsthresholds controlling warning/critical classification
Returns
scaling report with per-variable and per-constraint metadata
See also
ScalingReport for the returned structure

◆ casadi_opti()

const casadi::Opti & janus::Opti::casadi_opti ( ) const
inline

Access the underlying CasADi Opti object.

For advanced usage when CasADi-specific features are needed.

Returns
reference to the CasADi Opti object

◆ constrain_derivative() [1/2]

void janus::Opti::constrain_derivative ( const SymbolicVector & derivative,
const SymbolicVector & var,
const NumericVector & with_respect_to,
const std::string & method = "trapezoidal" )
inline

Constrain an existing variable to be a derivative.

Adds constraints: d(variable)/d(with_respect_to) == derivative

Parameters
derivativeThe derivative variable
varThe variable being differentiated
with_respect_toIndependent variable (e.g., time)
methodIntegration method (string or IntegrationMethod enum)

◆ constrain_derivative() [2/2]

void janus::Opti::constrain_derivative ( const SymbolicVector & derivative,
const SymbolicVector & var,
const NumericVector & with_respect_to,
IntegrationMethod method )
inline

Constrain an existing variable to be a derivative (enum version).

◆ derivative_of()

SymbolicVector janus::Opti::derivative_of ( const SymbolicVector & var,
const NumericVector & with_respect_to,
double derivative_init_guess,
const std::string & method = "trapezoidal" )
inline

Create a derivative variable constrained by integration.

Returns a new variable that is constrained to be the derivative of variable with respect to with_respect_to. This is the core mechanism for trajectory optimization via direct collocation.

Example:

NumericVector time = janus::linspace(0, 1, 100);
auto position = opti.variable(100, 0.0);
auto velocity = opti.derivative_of(position, time, 0.0);
// velocity is now constrained: d(position)/dt = velocity
JanusVector< NumericScalar > NumericVector
Eigen::VectorXd equivalent.
Definition JanusTypes.hpp:67
JanusVector< T > linspace(const T &start, const T &end, int n)
Generates linearly spaced vector.
Definition Spacing.hpp:26
Parameters
varThe quantity to differentiate
with_respect_toIndependent variable (e.g., time array)
derivative_init_guessInitial guess for derivative values
methodIntegration method: "trapezoidal", "forward_euler", "backward_euler"
Returns
Symbolic vector representing the derivative

◆ get_category()

std::vector< SymbolicScalar > janus::Opti::get_category ( const std::string & category) const
inline

Get all variables in a category.

Parameters
categoryCategory name
Returns
Vector of symbolic scalars in that category (empty if not found)

◆ get_category_names()

std::vector< std::string > janus::Opti::get_category_names ( ) const
inline

Get all registered category names.

Returns
Vector of category names

◆ is_category_frozen()

bool janus::Opti::is_category_frozen ( const std::string & category) const
inline

Check if a category is marked for freezing.

Parameters
categoryCategory name
Returns
True if the category was specified in constructor to be frozen

◆ maximize() [1/2]

void janus::Opti::maximize ( const SymbolicScalar & objective)
inline

Set objective to maximize.

Parameters
objectiveSymbolic expression to maximize

◆ maximize() [2/2]

void janus::Opti::maximize ( const SymbolicScalar & objective,
double objective_scale )
inline

Set objective to maximize with explicit objective scaling.

The solver sees -objective / objective_scale.

Parameters
objectivesymbolic expression to maximize
objective_scalepositive scale factor

◆ minimize() [1/2]

void janus::Opti::minimize ( const SymbolicScalar & objective)
inline

Set objective to minimize.

Parameters
objectiveSymbolic expression to minimize

◆ minimize() [2/2]

void janus::Opti::minimize ( const SymbolicScalar & objective,
double objective_scale )
inline

Set objective to minimize with explicit objective scaling.

The solver sees objective / objective_scale.

Parameters
objectivesymbolic expression to minimize
objective_scalepositive scale factor

◆ parameter() [1/2]

SymbolicVector janus::Opti::parameter ( const NumericVector & value)
inline

Create a vector parameter.

Parameters
valueVector of parameter values
Returns
Symbolic vector representing the parameters

◆ parameter() [2/2]

SymbolicScalar janus::Opti::parameter ( double value)
inline

Create a scalar parameter.

Parameters are fixed values that can be changed between solves.

Parameters
valueParameter value
Returns
Symbolic scalar representing the parameter

◆ set_initial()

void janus::Opti::set_initial ( const casadi::MX & x,
const casadi::DM & v )
inline

Set initial guess for a symbolic expression.

Forwards to CasADi's set_initial for element-level initial guesses (e.g. matrix slices from transcription methods).

Parameters
xSymbolic expression
vInitial value

◆ solve()

OptiSol janus::Opti::solve ( const OptiOptions & options = {})
inline

Solve the optimization problem.

Uses the solver specified in options (default: IPOPT).

Parameters
optionsSolver configuration (solver, max_iter, verbose, etc.)
Returns
OptiSol containing optimized values
Exceptions
RuntimeErrorif selected solver unavailable or fails

◆ solve_sweep()

SweepResult janus::Opti::solve_sweep ( const SymbolicScalar & param,
const std::vector< double > & values,
const OptiOptions & options = {} )
inline

Perform parametric sweep over a parameter.

Solves the optimization problem for each value in the sweep range, automatically warm-starting subsequent solves from the previous solution.

Example:

auto rho = opti.parameter(1.225); // Air density
// ... define problem using rho ...
auto result = opti.solve_sweep(rho, {1.0, 1.1, 1.2, 1.3});
for (size_t i = 0; i < result.size(); ++i) {
if (result.converged[i]) {
std::cout << "rho=" << result.param_values[i]
<< " x*=" << result.solutions[i]->value(x) << "\n";
}
}
Parameters
paramParameter to sweep (from opti.parameter())
valuesVector of parameter values to try
optionsSolver options (applied to all solves)
Returns
SweepResult containing all solutions

◆ subject_to() [1/6]

void janus::Opti::subject_to ( const std::vector< SymbolicScalar > & constraints)
inline

Add multiple constraints.

Parameters
constraintsVector of constraint expressions

◆ subject_to() [2/6]

void janus::Opti::subject_to ( const std::vector< SymbolicScalar > & constraints,
double linear_scale )
inline

Add multiple constraints with one shared linear scale.

Parameters
constraintsvector of constraint expressions
linear_scalepositive scale factor applied to all

◆ subject_to() [3/6]

void janus::Opti::subject_to ( const SymbolicScalar & constraint)
inline

Add a scalar constraint.

Supports both equality and inequality constraints:

opti.subject_to(x >= 0); // inequality
opti.subject_to(x == 1); // equality
opti.subject_to(x * x + y * y <= 1); // nonlinear
Parameters
constraintSymbolic inequality/equality expression

◆ subject_to() [4/6]

void janus::Opti::subject_to ( const SymbolicScalar & constraint,
double linear_scale )
inline

Add a scaled scalar constraint.

The provided scale is forwarded to CasADi so the solver sees constraint / linear_scale.

Parameters
constraintsymbolic inequality/equality expression
linear_scalepositive scale factor

◆ subject_to() [5/6]

void janus::Opti::subject_to ( std::initializer_list< SymbolicScalar > constraints)
inline

Add constraints from initializer list.

Example:

opti.subject_to({
x >= 0,
y >= 0,
x + y <= 10
});

◆ subject_to() [6/6]

void janus::Opti::subject_to ( std::initializer_list< SymbolicScalar > constraints,
double linear_scale )
inline

Add an initializer-list of constraints with one shared linear scale.

Parameters
constraintsinitializer list of constraint expressions
linear_scalepositive scale factor applied to all

◆ subject_to_bounds() [1/2]

void janus::Opti::subject_to_bounds ( const SymbolicScalar & scalar,
double lower_bound,
double upper_bound )
inline

Apply both lower and upper bounds to a scalar variable.

Parameters
scalarSymbolic scalar
lower_boundLower bound
upper_boundUpper bound

◆ subject_to_bounds() [2/2]

void janus::Opti::subject_to_bounds ( const SymbolicVector & vec,
double lower_bound,
double upper_bound )
inline

Apply both lower and upper bounds to all elements.

Example:

auto x = opti.variable(10, 0.5);
opti.subject_to_bounds(x, 0.0, 1.0); // 0 <= x[i] <= 1
Parameters
vecSymbolic vector
lower_boundLower bound for all elements
upper_boundUpper bound for all elements

◆ subject_to_lower() [1/2]

void janus::Opti::subject_to_lower ( const SymbolicScalar & scalar,
double lower_bound )
inline

Apply lower bound to a scalar variable.

Parameters
scalarSymbolic scalar
lower_boundLower bound

◆ subject_to_lower() [2/2]

void janus::Opti::subject_to_lower ( const SymbolicVector & vec,
double lower_bound )
inline

Apply lower bound to all elements of a vector.

Example:

auto x = opti.variable(10, 0.0);
opti.subject_to_lower(x, 0.0); // All x[i] >= 0
Parameters
vecSymbolic vector
lower_boundLower bound for all elements

◆ subject_to_upper() [1/2]

void janus::Opti::subject_to_upper ( const SymbolicScalar & scalar,
double upper_bound )
inline

Apply upper bound to a scalar variable.

Parameters
scalarSymbolic scalar
upper_boundUpper bound

◆ subject_to_upper() [2/2]

void janus::Opti::subject_to_upper ( const SymbolicVector & vec,
double upper_bound )
inline

Apply upper bound to all elements of a vector.

Parameters
vecSymbolic vector
upper_boundUpper bound for all elements

◆ variable() [1/4]

SymbolicVector janus::Opti::variable ( const NumericVector & init_guess,
std::optional< double > scale = std::nullopt,
std::optional< double > lower_bound = std::nullopt,
std::optional< double > upper_bound = std::nullopt )
inline

Create a vector of decision variables with per-element init guess.

Parameters
init_guessVector of initial guesses (size determines n_vars)
scaleOptional scale for numerical conditioning
lower_boundOptional lower bound for all elements
upper_boundOptional upper bound for all elements
Returns
Symbolic vector of optimization variables

◆ variable() [2/4]

SymbolicScalar janus::Opti::variable ( double init_guess,
const VariableOptions & opts,
std::optional< double > scale = std::nullopt,
std::optional< double > lower_bound = std::nullopt,
std::optional< double > upper_bound = std::nullopt )
inline

Create a scalar decision variable with options.

Parameters
init_guessInitial guess value (where optimizer starts)
optsVariable options (category, freeze)
scaleOptional scale for numerical conditioning (default: |init_guess| or 1)
lower_boundOptional lower bound constraint
upper_boundOptional upper bound constraint
Returns
Symbolic scalar representing the variable (or frozen parameter)

◆ variable() [3/4]

SymbolicScalar janus::Opti::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 )
inline

Create a scalar decision variable.

Parameters
init_guessInitial guess value (where optimizer starts)
scaleOptional scale for numerical conditioning (default: |init_guess| or 1)
lower_boundOptional lower bound constraint
upper_boundOptional upper bound constraint
Returns
Symbolic scalar representing the optimization variable

◆ variable() [4/4]

SymbolicVector janus::Opti::variable ( int n_vars,
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 )
inline

Create a vector of decision variables with scalar init guess.

Parameters
n_varsNumber of variables
init_guessInitial guess applied to all elements
scaleOptional scale for numerical conditioning
lower_boundOptional lower bound for all elements
upper_boundOptional upper bound for all elements
Returns
Symbolic vector (column vector) of optimization variables

The documentation for this class was generated from the following file: