Integrators¶
RK45Integrator¶
openscvx.integrators.solve_ivp_rk45(f: Callable[[jnp.ndarray, jnp.ndarray, Any], jnp.ndarray], tau_final: float, y_0: jnp.ndarray, args, tau_0: float = 0.0, num_substeps: int = 50, is_not_compiled: bool = False)
¶
Solve an initial-value ODE problem using fixed-step RK45 integration.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
f
|
Callable[[ndarray, ndarray, Any], ndarray]
|
ODE right-hand side; signature f(t, y, *args) -> dy/dt. |
required |
tau_final
|
float
|
Final integration time. |
required |
y_0
|
ndarray
|
Initial state at tau_0. |
required |
args
|
tuple
|
Extra arguments to pass to |
required |
tau_0
|
float
|
Initial time. Defaults to 0.0. |
0.0
|
num_substeps
|
int
|
Number of output time points. Defaults to 50. |
50
|
is_not_compiled
|
bool
|
If True, use Python loop instead of
JAX |
False
|
Returns:
Type | Description |
---|---|
jnp.ndarray: Array of shape (num_substeps, state_dim) with solution at each time. |
openscvx.integrators.rk45_step(f: Callable[[jnp.ndarray, jnp.ndarray, Any], jnp.ndarray], t: jnp.ndarray, y: jnp.ndarray, h: float, *args) -> jnp.ndarray
¶
Perform a single RK45 (Runge-Kutta-Fehlberg) integration step.
This implements the classic Dorman-Prince coefficients for an explicit 4(5) method, returning the fourth-order estimate.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
f
|
Callable[[ndarray, ndarray, Any], ndarray]
|
ODE right-hand side; signature f(t, y, *args) -> dy/dt. |
required |
t
|
ndarray
|
Current time. |
required |
y
|
ndarray
|
Current state vector. |
required |
h
|
float
|
Step size. |
required |
*args
|
Additional arguments passed to |
()
|
Returns:
Type | Description |
---|---|
ndarray
|
jnp.ndarray: Next state estimate at t + h. |
Diffrax Integrators¶
openscvx.integrators.solve_ivp_diffrax(f: Callable[[jnp.ndarray, jnp.ndarray, Any], jnp.ndarray], tau_final: float, y_0: jnp.ndarray, args, tau_0: float = 0.0, num_substeps: int = 50, solver_name: str = 'Dopri8', rtol: float = 0.001, atol: float = 1e-06, extra_kwargs=None)
¶
Solve an initial-value ODE problem using a Diffrax adaptive solver.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
f
|
Callable[[ndarray, ndarray, Any], ndarray]
|
ODE right-hand side; f(t, y, *args). |
required |
tau_final
|
float
|
Final integration time. |
required |
y_0
|
ndarray
|
Initial state at tau_0. |
required |
args
|
tuple
|
Extra arguments to pass to |
required |
tau_0
|
float
|
Initial time. Defaults to 0.0. |
0.0
|
num_substeps
|
int
|
Number of save points between tau_0 and tau_final. Defaults to 50. |
50
|
solver_name
|
str
|
Key into SOLVER_MAP for the Diffrax solver class. Defaults to "Dopri8". |
'Dopri8'
|
rtol
|
float
|
Relative tolerance for adaptive stepping. Defaults to 1e-3. |
0.001
|
atol
|
float
|
Absolute tolerance for adaptive stepping. Defaults to 1e-6. |
1e-06
|
extra_kwargs
|
dict
|
Additional keyword arguments forwarded to |
None
|
Returns:
Type | Description |
---|---|
jnp.ndarray: Solution states at the requested save points, shape (num_substeps, state_dim). |
Raises:
Type | Description |
---|---|
ValueError
|
If |
openscvx.integrators.solve_ivp_diffrax_prop(f: Callable[[jnp.ndarray, jnp.ndarray, Any], jnp.ndarray], tau_final: float, y_0: jnp.ndarray, args, tau_0: float = 0.0, num_substeps: int = 50, solver_name: str = 'Dopri8', rtol: float = 0.001, atol: float = 1e-06, extra_kwargs=None, save_time: jnp.ndarray = None, mask: jnp.ndarray = None)
¶
Solve an initial-value ODE problem using a Diffrax adaptive solver. This function is specifically designed for use in the context of trajectory optimization and handles the nonlinear single-shot propagation of state variables in undilated time.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
f
|
Callable[[ndarray, ndarray, Any], ndarray]
|
ODE right-hand side; signature f(t, y, *args) -> dy/dt. |
required |
tau_final
|
float
|
Final integration time. |
required |
y_0
|
ndarray
|
Initial state at tau_0. |
required |
args
|
tuple
|
Extra arguments to pass to |
required |
tau_0
|
float
|
Initial time. Defaults to 0.0. |
0.0
|
num_substeps
|
int
|
Number of save points between tau_0 and tau_final. Defaults to 50. |
50
|
solver_name
|
str
|
Key into SOLVER_MAP for the Diffrax solver class. Defaults to "Dopri8". |
'Dopri8'
|
rtol
|
float
|
Relative tolerance for adaptive stepping. Defaults to 1e-3. |
0.001
|
atol
|
float
|
Absolute tolerance for adaptive stepping. Defaults to 1e-6. |
1e-06
|
extra_kwargs
|
dict
|
Additional keyword arguments forwarded to |
None
|
save_time
|
ndarray
|
Time points at which to evaluate the solution. Must be provided for export compatibility. |
None
|
mask
|
ndarray
|
Boolean mask for the save_time points. |
None
|
Returns:
Type | Description |
---|---|
jnp.ndarray: Solution states at the requested save points, shape (num_substeps, state_dim). |
Raises:
ValueError: If solver_name
is not in SOLVER_MAP or if save_time is not provided.