nutopy.ivp.dexp¶
- dexp(df, tf, tfd, t0, t0d, x0, x0d, pars=None, parsd=None, *, options=None, time_steps=None, args=(), kwargs={})[source]¶
The
dexp
function computes a numerical approximation of the derivative of the solution of the initial value problem\[\frac{d}{d t} x = \dot{x} = f(t, x, pars), \quad x(t_0) = x_0,\]where \(pars \in \mathbf{R}^k\) is a vector of
parameters
, \(t\) is thetime
and \(x\) is thestate
variable. The time \(t_0\) is refered as theinitial time
, \(x_0\) as theinitial condition
and \(f\) as thedynamics
. The solution of the initial value problem at thefinal time
\(t_f\) is denoted \(x(t_f, t_0, x_0, pars)\).Hence,
dexp
approximates\[\begin{split}x'(t_f, t_0, x_0, pars) \cdot (\delta t_f, \delta t_0, \delta x_0, \delta pars) =& \phantom{+} \frac{\partial x}{\partial t_f}(t_f, t_0, x_0, pars) \cdot \delta t_f \\ & + \frac{\partial x}{\partial t_0}(t_f, t_0, x_0, pars) \cdot \delta t_0 \\ & + \frac{\partial x}{\partial x_0}(t_f, t_0, x_0, pars) \cdot \delta x_0 \\ & + \frac{\partial x}{\partial pars}(t_f, t_0, x_0, pars) \cdot \delta pars.\end{split}\]Let us introduce \(\delta x(t) = x'(t, t_0, x_0, pars) \cdot (0, \delta t_0, \delta x_0, \delta pars)\). Then, \(\delta x(t_f)\) is computed solving the following variational equations:
\[\begin{split}\dot{x} &= f(t, x, pars), \\ \dot{\delta x} &= \frac{\partial f}{\partial x}(t, x, pars) \cdot \delta x + \frac{\partial f}{\partial pars}(t, x, pars) \cdot \delta pars, \\ x(t_0) &= x_0, \\ \delta x(t_0) &= \delta x_0 - f(t_0, x_0, pars)\, \delta t_0.\end{split}\]The missing part is given by:
\[\frac{\partial x}{\partial t_f}(t_f, t_0, x_0, pars) \cdot \delta t_f = f(t_f, x(t_f, t_0, x_0, pars), pars)\, \delta t_f.\]- Parameters
df (callable) –
A vector function
y, yd = df(t, x, xd)
if(pars, parsd)
isNone
, wherey = f(t, x)
andyd = df/dx(t, x) . xd
or
y, yd = df(t, x, xd, pars, parsd)
if(pars, parsd)
is notNone
, wherey = f(t, x, pars)
andyd = df/dx(t, x, pars) . xd + df/dpars(t, x, pars) . parsd
.tf (
float
) – Final timetfd (
float
) – Increment \(\delta t_f\)t0 (
float
) – Initial timet0d (
float
) – Increment \(\delta t_0\)x0 (float vector) – Initial condition
x0d (float vector) – Increment \(\delta x_0\)
pars (float vector, optional, positional) – Parameters
parsd (float vector, optional, positional) – Increment \(\delta pars\)
options (
nutopy.ivp.Options
, optional, key only) – A dictionnary of solver options. The integrator scheme is chosen settingSolverMethod
.time_steps (float vector, optional, key only) –
If
time_steps
is notNone
, then it is a vector of time steps, that issol.tout = time_steps
. Iftime_steps
isNone
, thensol.tout = [t0, tf]
if a fixed-step integrator scheme is used whilesol.tout
is given by the integrator if a variable step-size scheme is used.Remark. Note that in the case of a variable step-size integrator scheme,
time_steps
is used only for output, that is the steps are computed by the integrator but dense output is used to provide by interpolation the solution x at times intime_steps
.args (tuple, key only) – optional arguments for
df
kwargs (dictionnary, key only) – keywords optional arguments for
df
- Returns
sol – A dictionary/struct of outputs:
sol.xf
orsol.get('xf')
: final pointx(tf, t0, x0, pars)
;sol.xfd
: derivativex'(tf, t0, x0, pars) . (tfd, t0d, x0d, parsd)
;sol.tout
: the integration time-steps:tout=time_steps
if provided;sol.xout
: the integration points:xout[i] = x(tout[i], t0, x0, pars)
;sol.xdout
: derivative at integration time-steps:xdout[i] = x'(tout[i], t0, x0, pars) . (tfd, t0d, x0d, parsd)
;sol.success
: whether or not the integrator exited successfully;sol.status
: termination status of the integrator;sol.message
: description of the cause of the termination;sol.nfev
: number of evaluations of the dynamics;sol.nsteps
: number of integration time-steps.
Notes
Use
numpy.ndarray
to define float vectors.nutopy.ivp.dexp
is the derivative ofnutopy.ivp.exp
Examples
Example: dx/dx2, x = (x1, x2)
>>> import numpy as np >>> from nutopy import ivp
>>> def df(t, x, xd): ... y = np.zeros(x.size) ... y[0] = -x[0] + x[1] ... y[1] = x[1] ... yd = np.zeros(x.size) ... yd[0] = -xd[0] + xd[1] ... yd[1] = xd[1] ... return (y, yd)
>>> tf = 1.0 >>> t0 = 0.0 >>> x0 = np.array([-1.5, 0.5]) >>> tfd = 0.0 >>> t0d = 0.0 >>> x0d = np.array([0.0, 1.0]) >>> sol = ivp.dexp(df, tf, tfd, t0, t0d, x0, x0d) >>> print('dx/dx2(tf, t0, x0) = ', sol.xfd) dx/dx2(tf, t0, x0) = [1.17520119 2.71828183]
- Raises
ArgumentTypeError – Variable tf must be a float or cast to float.
ArgumentTypeError – Variable tfd must be a float or cast to float.
ArgumentTypeError – Variable t0 must be a float or cast to float.
ArgumentTypeError – Variable t0d must be a float or cast to float.
ArgumentTypeError – Variable options must be a Options object.
ArgumentDimensionError – Variable x0 must be a one dimensional array.
ArgumentDimensionError – Variable x0d must be a one dimensional array.
ArgumentDimensionError – Variable x0d must have the same length than x0.
ArgumentDimensionError – Variable pars must be a one dimensional array.
ArgumentDimensionError – Variable parsd must be a one dimensional array.
ArgumentDimensionError – Variable parsd must have the same length than pars.
InputArgumentError – Either give pars and parsd or any.
InputArgumentError – You must have time_steps[0]=t0 and time_steps[-1]=tf.
InputArgumentError – The argument time_steps cannot be of size 1.