API#
SCP algorithms and utility functions#
- class scocp.SCvxStar(problem, tol_opt=1e-06, tol_feas=1e-06, rho0=0.0, rho1=0.25, rho2=0.7, alpha1=2.0, alpha2=3.0, beta=2.0, gamma=0.9, r_bounds=[1e-08, 10.0], steps_minimum_trust_region=10, weight_max=1e+16)#
SCvx* algorithm for optimal control problems
Hyperparameters are defined according to Oguri, 2023 (doi: 10.1109/CDC49753.2023.10383462).
- Parameters:
problem (ContinuousControlSCOCP or ImpulsiveControlSCOCP) – SCOCP instance
tol_opt (float) – optimality tolerance
tol_feas (float) – feasibility tolerance
rho0 (float) – step acceptance criterion parameter
rho1 (float) – trust-region radius contraction threshold
rho2 (float) – trust-region radius expansion threshold
alpha1 (float) – trust-region radius contraction factor s.t. r_k+1 = max(r_k/alpha1, r_bounds[0])
alpha2 (float) – trust-region radius expansion factor s.t. r_k+1 = min(r_k*alpha2, r_bounds[1])
beta (float) – weight update factor
gamma (float) – update factor for Lagrange multiplier update criterion delta
r_bounds (list) – trust region bounds
weight_max (float) – maximum weight
- evaluate_penalty(gdyn, g, h)#
Evaluate penalty function according to Augmented Lagrangian formulation
- Parameters:
gdyn (np.array) – (N-1)-by-nx array of nonlinear dynamics constraints violations
g (np.array) – ng-by-1 array of nonlinear equality constraints violations
h (np.array) – nh-by-1 array of nonlinear inequality constraints violations
- plot_DeltaJ(axis, summary_dict: dict, s=5)#
Plot iterations of DeltaJ
- plot_J0(axis, summary_dict: dict, s=5)#
Plot iterations of J0
- plot_chi(axis, summary_dict: dict, s=5)#
Plot iterations of chi
- plot_trust_region_radius_u(axis, summary_dict: dict, s=5)#
Plot iterations of trust region radius u
- plot_trust_region_radius_x(axis, summary_dict: dict, s=5)#
Plot iterations of trust region radius x
- plot_w(axis, summary_dict: dict, s=5)#
Plot iterations of weight
- solve(xbar, ubar, vbar=None, ybar=None, maxiter: int = 10, verbose: bool = True, feasability_norm=inf, debug=False, exit_condition_rho=False, save_all_iterations=False)#
Solve optimal control problem via SCvx* algorithm
- Parameters:
xbar (np.array) – N-by-nx array of reference states
ubar (np.array) – N-by-nu array of reference controls
vbar (np.array) – N-by-1 array of reference constraints
maxiter (int) – maximum number of iterations
verbose (bool) – whether to print verbose output
feasability_norm (str) – norm to use for feasibility evaluation
debug (bool) – whether to print debug output
exit_condition_rho (bool) – whether to exit if rho < rho0
save_all_iterations (bool) – whether to save all iterations
- Returns:
solution object
- Return type:
SCPSolution
- scocp.get_augmented_lagrangian_penalty(weight, xi_dyn, lmb_dyn, xi=None, lmb_eq=None, zeta=None, lmb_ineq=None)#
Evaluate augmented Lagrangian penalty function
- Parameters:
weight (float) – weight of the penalty function
xi_dyn (cp.Variable) – slack variable for dynamics
lmb_dyn (cp.Parameter) – multiplier for dynamics
xi (cp.Variable, optional) – slack variable for equality constraints
lmb_eq (cp.Parameter, optional) – multiplier for equality constraints
zeta (cp.Variable, optional) – slack variable for inequality constraints
lmb_ineq (cp.Parameter, optional) – multiplier for inequality constraints
- Returns:
augmented Lagrangian penalty function
- Return type:
(cp.Expression)
- class scocp.MovingTarget(eval_state: Callable, eval_state_derivative: Callable)#
Define moving target for rendezvous problem
We assume the target state equality constraint is of the form
- g(r_N,v_N,t_N) = [ r_N - r_ref(t_N)
v_N - v_ref(t_N) ]
where r_ref(t_N) and v_ref(t_N) are the position and velocity of the target at time t_N.
- Parameters:
eval_state (function) – function to evaluate target state
eval_state_derivative (function) – function to evaluate derivative of target state
- target_state(t: float) ndarray #
Get target state at time t
- target_state_derivative(t: float) ndarray #
Get target state derivative w.r.t. time
Continuous control problems#
Base class for continuous control problems#
- class scocp.ContinuousControlSCOCP(integrator, times, ng: int = 0, nh: int = 0, ny: int = 0, augment_Gamma: bool = False, weight: float = 100.0, trust_region_radius_x: float = 0.1, trust_region_radius_u: float = None, solver='CLARABEL', verbose_solver: bool = False, impulsive_B=None)#
Sequentially convexified optimal control problem (SCOCP) for continuous dynamics
- Parameters:
integrator (obj) – integrator object
times (np.array) – time grid
ng (int) – number of nonlinear equality constraints, excluding dynamics constraints
nh (int) – number of nonlinear inequality constraints
ny (int) – number of other variables
augment_Gamma (bool) – whether to augment the control with the constraint vector when integrating the dynamics
weight (float) – weight of the objective function
trust_region_radius_x (float) – trust region radius for state
trust_region_radius_u (float) – trust region radius for control
solver (str) – solver to use
verbose_solver (bool) – whether to print verbose output
- build_linear_model(xbar, ubar, vbar)#
Construct linear model for dynamics multiple-shooting constraints within SCP algorithm This function computes the Phi_A, Phi_B, and Phi_c matrices and stores them in the class attributes.
- Parameters:
xbar (np.array) – (N, self.integrator.nx) array of state history
ubar (np.array) – (N-1, self.integrator.nu) array of control history
vbar (np.array) – (N-1, self.integrator.nv) array of constraint history
- evaluate_nonlinear_constraints(xs, us, gs, ys=None)#
Evaluate nonlinear constraints
- Returns:
tuple of 1D arrays of nonlinear equality and inequality constraints
- Return type:
(tuple)
- evaluate_nonlinear_dynamics(xs, us, vs, stm=False, steps=None)#
Evaluate nonlinear dynamics along given state and control history
- Parameters:
integrator (obj) – integrator object
times (np.array) – time grid
vs (np.array) – state history
ubar (np.array) – control history
stm (bool) – whether to propagate STMs, defaults to False
- evaluate_objective(xs, us, vs, ys=None)#
Evaluate the objective function
- reset()#
Reset problem parameters and storages
- solve_convex_problem(xbar, ubar, vbar, ybar=None)#
Solve the convex subproblem
- Parameters:
xbar (np.array) – (N, self.integrator.nx) array of reference state history
ubar (np.array) – (N-1, self.integrator.nu) array of reference control history
vbar (np.array) – (N-1, self.integrator.nv) array of reference constraint history
ybar (np.array) – (self.ny,) other reference variables
- Returns:
np.array values of xs, us, gs, xi_dyn, xi_eq, zeta_ineq
- Return type:
(tuple)
Example continuous control problems#
- class scocp.FixedTimeContinuousRdv(x0, xf, umax, *args, **kwargs)#
Fixed-time continuous rendezvous problem class
- evaluate_objective(xs, us, vs, ys=None)#
Evaluate the objective function
- solve_convex_problem(xbar, ubar, vbar, ybar=None)#
Solve the convex subproblem
- Parameters:
xbar (np.array) – (N, self.integrator.nx) array of reference state history
ubar (np.array) – (N-1, self.integrator.nu) array of reference control history
vbar (np.array) – (N-1, self.integrator.nv) array of reference constraint history
- Returns:
np.array values of xs, us, vs, xi_dyn, xi_eq, zeta_ineq
- Return type:
(tuple)
- class scocp.FixedTimeContinuousRdvLogMass(x0, xf, Tmax, N, *args, **kwargs)#
Fixed-time continuous rendezvous problem class with log-mass dynamics
- evaluate_nonlinear_constraints(xs, us, vs, ys=None)#
Evaluate nonlinear constraints
- Returns:
tuple of 1D arrays of nonlinear equality and inequality constraints
- Return type:
(tuple)
- evaluate_objective(xs, us, vs, ys=None)#
Evaluate the objective function
- solve_convex_problem(xbar, ubar, vbar, ybar=None)#
Solve the convex subproblem
- Parameters:
xbar (np.array) – (N, self.integrator.nx) array of reference state history
ubar (np.array) – (N-1, self.integrator.nu) array of reference control history
vbar (np.array) – (N-1, self.integrator.nv) array of reference constraint history
- Returns:
np.array values of xs, us, gs, xi_dyn, xi_eq, zeta_ineq
- Return type:
(tuple)
- class scocp.FreeTimeContinuousRdv(x0, xf, umax, tf_bounds, s_bounds, *args, **kwargs)#
Free-time continuous rendezvous problem class A good initial guess for the dilation factor s is the guessed TOF, since dt/dtau = s.
- Parameters:
x0 (np.array) – initial state
xf (np.array) – final state
umax (float) – maximum control magnitude
tf_bounds (list) – bounds on the final time, given as [tf_min, tf_max]
s_bounds (list) – bounds on time dilation factor, given as [s_min, s_max]
- evaluate_objective(xs, us, vs, ys=None)#
Evaluate the objective function
- solve_convex_problem(xbar, ubar, vbar, ybar=None)#
Solve the convex subproblem
- Parameters:
xbar (np.array) – (N, self.integrator.nx) array of reference state history
ubar (np.array) – (N-1, self.integrator.nu) array of reference control history
vbar (np.array) – (N-1, self.integrator.nv) array of reference constraint history
- Returns:
np.array values of xs, us, gs, xi_dyn, xi_eq, zeta_ineq
- Return type:
(tuple)
- class scocp.FreeTimeContinuousRdvLogMass(x0, xf, Tmax, tf_bounds, s_bounds, N, *args, **kwargs)#
Free-time continuous rendezvous problem class with log-mass dynamics
Note the ordering expected for the state and the control vectors:
state = [x,y,z,vx,vy,vz,log(mass),t] u = [ax, ay, az, s, Gamma] where s is the dilation factor, Gamma is the control magnitude (at convergence)
- evaluate_nonlinear_constraints(xs, us, gs, ys=None)#
Evaluate nonlinear constraints
- Returns:
tuple of 1D arrays of nonlinear equality and inequality constraints
- Return type:
(tuple)
- evaluate_objective(xs, us, vs, ys=None)#
Evaluate the objective function
- solve_convex_problem(xbar, ubar, vbar, ybar=None)#
Solve the convex subproblem
- Parameters:
xbar (np.array) – (N, self.integrator.nx) array of reference state history
ubar (np.array) – (N-1, self.integrator.nu) array of reference control history
vbar (np.array) – (N-1, self.integrator.nv) array of reference constraint history
- Returns:
np.array values of xs, us, gs, xi_dyn, xi_eq, zeta_ineq
- Return type:
(tuple)
- class scocp.FreeTimeContinuousMovingTargetRdvLogMass(x0, target: MovingTarget, Tmax, tf_bounds, s_bounds, N, *args, **kwargs)#
Free-time continuous moving target rendezvous problem class with log-mass dynamics
Note the ordering expected for the state and the control vectors:
state = [x,y,z,vx,vy,vz,log(mass),t] u = [ax, ay, az, s, Gamma] where s is the dilation factor, Gamma is the control magnitude (at convergence)
The target state is given by a MovingTarget object, which result in 6 non-convex equality constraints.
- evaluate_nonlinear_constraints(xs, us, vs, ys=None)#
Evaluate nonlinear constraints
- Returns:
tuple of 1D arrays of nonlinear equality and inequality constraints
- Return type:
(tuple)
- evaluate_objective(xs, us, vs, ys=None)#
Evaluate the objective function
- solve_convex_problem(xbar, ubar, vbar, ybar=None)#
Solve the convex subproblem
- Parameters:
xbar (np.array) – (N, self.integrator.nx) array of reference state history
ubar (np.array) – (N-1, self.integrator.nu) array of reference control history
vbar (np.array) – (N-1, self.integrator.nv) array of reference constraint history
- Returns:
np.array values of xs, us, gs, xi_dyn, xi_eq, zeta_ineq
- Return type:
(tuple)
- class scocp.FreeTimeContinuousMovingTargetRdvMass(x0, target: MovingTarget, Tmax, tf_bounds, s_bounds, N, *args, **kwargs)#
Free-time continuous moving target rendezvous problem class with mass dynamics
Note the ordering expected for the state and the control vectors:
state = [x,y,z,vx,vy,vz,mass,t] u = [ux, uy, uz, s, Gamma] where s is the dilation factor, Gamma is the control magnitude (at convergence)
The target state is given by a MovingTarget object, which result in 6 non-convex equality constraints.
- evaluate_nonlinear_constraints(xs, us, vs, ys=None)#
Evaluate nonlinear constraints
- Returns:
tuple of 1D arrays of nonlinear equality and inequality constraints
- Return type:
(tuple)
- evaluate_objective(xs, us, vs, ys=None)#
Evaluate the objective function
- solve_convex_problem(xbar, ubar, vbar, ybar=None)#
Solve the convex subproblem
- Parameters:
xbar (np.array) – (N, self.integrator.nx) array of reference state history
ubar (np.array) – (N-1, self.integrator.nu) array of reference control history
vbar (np.array) – (N-1, self.integrator.nv) array of reference constraint history
- Returns:
np.array values of xs, us, gs, xi_dyn, xi_eq, zeta_ineq
- Return type:
(tuple)
Impulsive control problems#
Base class for impulsive control problems#
- class scocp.ImpulsiveControlSCOCP(integrator, times, ng: int = 0, nh: int = 0, ny: int = 0, augment_Gamma: bool = False, B=None, weight: float = 100.0, trust_region_radius_x: float = 0.1, trust_region_radius_u: float = None, solver='CLARABEL', verbose_solver: bool = False)#
Sequentially convexified optimal control problem (SCOCP) for impulsive dynamics
- Parameters:
integrator (obj) – integrator object
times (np.array) – time grid
ng (int) – number of nonlinear equality constraints, excluding dynamics constraints
nh (int) – number of nonlinear inequality constraints
ny (int) – number of other variables
augment_Gamma (bool) – whether to augment the control with the constraint vector when integrating the dynamics
B (np.array) – control matrix
weight (float) – weight of the objective function
trust_region_radius_x (float) – trust region radius for state
trust_region_radius_u (float) – trust region radius for control
solver (str) – solver to use
verbose_solver (bool) – whether to print verbose output
- evaluate_nonlinear_constraints(xs, us, gs, ys=None)#
Evaluate nonlinear constraints
- Returns:
tuple of 1D arrays of nonlinear equality and inequality constraints
- Return type:
(tuple)
- evaluate_nonlinear_dynamics(xbar, ubar, stm=False, steps=None)#
Evaluate nonlinear dynamics along given state and control history
- Parameters:
integrator (obj) – integrator object
times (np.array) – time grid
xbar (np.array) – state history
ubar (np.array) – control history
stm (bool) – whether to propagate STMs, defaults to False
- evaluate_objective(xs, us, gs, ys=None)#
Evaluate the objective function
- reset()#
Reset problem parameters and storages
- solve_convex_problem(xbar, ubar, vbar, ybar=None)#
Solve the convex subproblem
Example impulsive control problems#
- class scocp.FixedTimeImpulsiveRdv(x0, xf, *args, **kwargs)#
Fixed-time impulsive rendezvous problem class
- evaluate_objective(xs, us, vs, ys=None)#
Evaluate the objective function
- solve_convex_problem(xbar, ubar, vbar, ybar=None)#
Solve the convex subproblem
- Parameters:
xbar (np.array) – (N, self.integrator.nx) array of reference state history
ubar (np.array) – (N-1, self.integrator.nu) array of reference control history
vbar (np.array) – (N-1, self.integrator.nv) array of reference constraint history
- Returns:
np.array values of xs, us, vs, ys, xi_dyn, xi_eq, zeta_ineq
- Return type:
(tuple)
Integrators#
- class scocp.ScipyIntegrator(nx, nu, rhs: Callable, rhs_stm: Callable, impulsive=True, nv=0, method='RK45', reltol=1e-12, abstol=1e-12, args=None)#
Integrator class using scipy.integrate.solve_ivp
- Parameters:
nx (int) – dimension of states
nu (int) – dimension of controls
rhs (function) – right-hand side function
rhs_stm (function) – right-hand side function for state transition matrix
impulsive (bool) – whether the dynamics are impulsive
nv (int) – dimensions corresponding to control norms to be augmented
method (str) – integration method
reltol (float) – relative tolerance
abstol (float) – absolute tolerance
args (tuple) – additional arguments for the right-hand side function
- solve(tspan, x0, u=None, stm=False, t_eval=None, args=None, get_ODESolution=False, dense_output=False)#
Solve initial value problem
- Parameters:
tspan (tuple) – time span
x0 (np.array) – initial state
stm (bool) – whether to solve for state transition matrix
t_eval (np.array) – evaluation times
args (tuple) – additional arguments for the right-hand side function
get_ODESolution (bool) – whether to return an ODESolution object
dense_output (bool) – whether to return a dense output
- Returns:
if get_ODESolution is False, return a tuple of times and state with shape N-by-nx if get_ODESolution is True, return an ODESolution object
- Return type:
(tuple or ODESolution)