Skip to content

RK4 Integrator

Classical 4th-order Runge-Kutta integrator for fixed-step integration.

RK4Integrator

RK4Integrator(dimension: int, dynamics_fn: callable, jacobian: Union[DAnalyticJacobian, DNumericalJacobian] = None, sensitivity: Any = None, control_fn: Any = None, config: IntegratorConfig = None)

4th-order Runge-Kutta fixed-step integrator.

Classical RK4 method with fixed timesteps. Provides good accuracy for most problems with 4th-order error convergence.

Parameters:

Name Type Description Default
dimension int

State vector dimension

required
dynamics_fn callable

Dynamics function with signature (t: float, state: ndarray) -> ndarray

required
jacobian DAnalyticJacobian or DNumericalJacobian

Jacobian provider for variational matrix propagation

None
config IntegratorConfig

Integration configuration

None
Example
import brahe as bh
import numpy as np

# Define dynamics: simple harmonic oscillator
def dynamics(t, state):
    return np.array([state[1], -state[0]])

# Create integrator for 2D system
integrator = bh.RK4Integrator(2, dynamics)

# Integrate one step
state = np.array([1.0, 0.0])
new_state = integrator.step(0.0, state, 0.01)

# With variational matrix (state transition matrix)
def jacobian_fn(t, state):
    return np.array([[0.0, 1.0], [-1.0, 0.0]])

jacobian = bh.DAnalyticJacobian(jacobian_fn)
integrator_varmat = bh.RK4Integrator(2, dynamics, jacobian)

phi = np.eye(2)
new_state, new_phi = integrator_varmat.step_with_varmat(0.0, state, phi, 0.01)

Initialize instance.

dimension property

dimension: Any

Get state vector dimension.

step method descriptor

step(t: float, state: ndarray, dt: float = None) -> ndarray

Perform one integration step.

Parameters:

Name Type Description Default
t float

Current time

required
state ndarray

State vector at time t

required
dt float

Timestep. If None, uses the step size from configuration.

None

Returns:

Name Type Description
ndarray ndarray

State vector at time t + dt

Example
1
2
3
4
5
6
7
# Using explicit dt
new_state = integrator.step(0.0, state, 0.01)

# Using config-based dt
config = bh.IntegratorConfig.fixed_step(0.01)
integrator = bh.RK4Integrator(2, dynamics, config=config)
new_state = integrator.step(0.0, state)

step_with_sensmat method descriptor

step_with_sensmat(t: float, state: ndarray, sens: ndarray, params: ndarray, dt: float = None) -> Tuple

Advance state and sensitivity matrix by one timestep.

Propagates the sensitivity matrix S that maps parameter uncertainties to state uncertainties. The sensitivity evolves according to dS/dt = A*S + B.

Parameters:

Name Type Description Default
t float

Current time

required
state ndarray

State vector at time t

required
sens ndarray

Sensitivity matrix at time t (state_dim x param_dim)

required
params ndarray

Parameter vector

required
dt float

Integration timestep

None

Returns:

Name Type Description
tuple Tuple

(new_state, new_sensitivity)

Example
import brahe as bh
import numpy as np

# Setup integrator with Jacobian and sensitivity providers
integrator = bh.RK4Integrator(6, dynamics, jacobian, sensitivity)
state = np.array([...])
sens = np.zeros((6, 2))  # 6 state dims, 2 params
params = np.array([1.0, 2.0])

new_state, new_sens = integrator.step_with_sensmat(0.0, state, sens, params, 1.0)

step_with_varmat method descriptor

step_with_varmat(t: float, state: ndarray, phi: ndarray, dt: float = None) -> Tuple

Perform one integration step with variational matrix propagation.

Parameters:

Name Type Description Default
t float

Current time

required
state ndarray

State vector at time t

required
phi ndarray

State transition matrix at time t (dimension × dimension)

required
dt float

Timestep. If None, uses the step size from configuration.

None

Returns:

Name Type Description
tuple Tuple

(new_state, new_phi) - State vector and STM at time t + dt

Example
1
2
3
4
5
6
phi = np.eye(2)
# Using explicit dt
new_state, new_phi = integrator.step_with_varmat(0.0, state, phi, 0.01)

# Using config-based dt
new_state, new_phi = integrator.step_with_varmat(0.0, state, phi)

step_with_varmat_sensmat method descriptor

step_with_varmat_sensmat(t: float, state: ndarray, phi: ndarray, sens: ndarray, params: ndarray, dt: float = None) -> Tuple

Advance state, variational matrix (STM), and sensitivity matrix by one timestep.

Propagates both matrices simultaneously for complete uncertainty quantification: - STM (Phi): Maps initial state uncertainties to current state uncertainties - Sensitivity (S): Maps parameter uncertainties to state uncertainties

Parameters:

Name Type Description Default
t float

Current time

required
state ndarray

State vector at time t

required
phi ndarray

State transition matrix at time t (state_dim x state_dim)

required
sens ndarray

Sensitivity matrix at time t (state_dim x param_dim)

required
params ndarray

Parameter vector

required
dt float

Integration timestep

None

Returns:

Name Type Description
tuple Tuple

(new_state, new_phi, new_sensitivity)

Example
import brahe as bh
import numpy as np

# Setup integrator with Jacobian and sensitivity providers
integrator = bh.RK4Integrator(6, dynamics, jacobian, sensitivity)
state = np.array([...])
phi = np.eye(6)
sens = np.zeros((6, 2))  # 6 state dims, 2 params
params = np.array([1.0, 2.0])

new_state, new_phi, new_sens = integrator.step_with_varmat_sensmat(
    0.0, state, phi, sens, params, 1.0
)

See Also