The JiTCODE module ================== **Note and remember** *that some relevant information can be found in the* `common JiTC*DE documentation`_. *This includes installation instructions, compiler issues and optimisation, general performance considerations, how to implement network dynamics and conditional statements, and a small FAQ.* Overview -------- JiTCODE (just-in-time compilation for ordinary differential equations) is an extension of SciPy’s `ODE`_ (`scipy.integrate.ode`) or `Solve IVP`_ (`scipy.integrate.solve_ivp`). Where the latter take a Python function as an argument, JiTCODE takes an iterable (or generator function or dictionary) of symbolic expressions, which it translates to C code, compiles on the fly, and uses as the function to feed into SciPy’s ODE or Solve IVP. Symbolic expressions are mostly handled by `SymEngine`_, `SymPy`_’s compiled-backend-to-be (see `SymPy vs. SymEngine`_ for details). This design has the following advantages: * **Speed boost through compilation** There are two main computational burdens when numerically solving ODEs. The first are simple vector operations for the integration procedure, which can be efficiently handled by NumPy or similar. The second is calculating the derivative (and Jacobian, if desired), which can be efficiently handled with just-in-time compilation. Note that for smaller ODEs and small integration times, the Python overhead may still take a considerable toll in comparison to pure compiled code. * **Speed boost through symbolic optimisation** If your derivative is automatically generated by some routine, you can simplify it symbolically to boost the speed. In fact, blatant optimisations such as :math:`y·(x-x)=0` are done on the fly by SymEngine. This is for example interesting if you want to simulate dynamics on a sparse network, as non-existing links are not taken into account when calculating the derivative when integrating. * **Automatically calculated Jacobian and Lyapunov exponents** As the derivative is provided symbolically, SymEngines’s automatic differentiation routines can be employed to calculate the Jacobian desired by some integrators as well as the ODE for the tangent vector required for calculating the Lyapunov exponents (see `lyapunov`). * **Symbolic interface** You can enter your differential equations almost like you would on paper. Also, if you are working with SymPy or SymEngine anyway – e.g., to calculate fixed points –, you do not need to bother much with translating your equations. If compilation fails to work for whatever reason, pure Python functions can be employed as a fallback. For most applications, the only difference to SciPy’s ODE in terms of handling is that you have to supply the derivative in the correct format – JiTCODE automatically takes care of the compilation. However, you can also tweak this step by step, if you desire (see `tweak` and the `notes on compilation`_). As with SciPy’s ODE, this documentation assumes that the differential equation you want to solve is: .. math:: \dot{y} = f(t,y) .. _example: A quick example --------------- *This example is particularly intended for readers familiar with SciPy’s* `ODE`_. *If you aren’t, skip to the next section.* .. automodule:: double_fhn_example .. literalinclude:: ../examples/double_fhn_example.py :linenos: :dedent: 1 :start-after: example-start An explicit example ------------------- *This example is particularly intended for readers who are using a numerical integrator for the first time.* .. automodule:: lotka_volterra .. _tweak: Details of the building process ------------------------------- To generate the functions needed by SciPy’s ODE or Solve IVP, JiTCODE executes a series of distinct processing steps, each of which can be invoked through a command and tweaked as needed. These commands execute previous steps as needed, i.e., if they have not been performed yet. If you are happy with the default options, however, you do not need to bother with this and can just use the commands at the very end of the chain, namely `set_integrator`, `set_initial_value`, or `save_compiled`. The following diagram details which command calls which other command when needed: .. digraph:: foo bgcolor = "transparent" edge [arrowsize=0.5]; node [shape=box fontname="Monospace" fontsize=10 margin=0.05 fillcolor="#C8D5E3" penwidth=0 style=filled]; {rank = same; "save_compiled"; "set_integrator"} "set_initial_value" -> "set_integrator"; "save_compiled" -> "compile_C"; "set_integrator" -> "generate_functions" -> "compile_C" -> "generate_f_C" -> "generate_helpers_C"; "compile_C" -> "generate_jac_C" -> "generate_jac_sym"; "generate_jac_C" -> "generate_helpers_C"; "generate_functions" -> "generate_lambdas" -> "generate_f_lambda"; "generate_lambdas" -> "generate_jac_lambda"-> "generate_jac_sym"; .. _example_2: A more complicated example -------------------------- .. automodule:: SW_of_Roesslers .. _lyapunov: Calculating Lyapunov exponents with `jitcode_lyap` -------------------------------------------------- `jitcode_lyap` is a simple extension of `jitcode` that almost automatically handles calculating Lyapunov exponents by evolving tangent vectors [BGGS80]_. It works just like `jitcode`, except that it generates and integrates additional differential equations for the tangent vectors. After every call of `integrate`, the tangent vectors are orthonormalised, and the “local” Lyapunov exponents for this integration step are returned alongside with the system’s state. These can then be further processed to obtain the Lyapunov exponents. The tangent vectors are initialised with random vectors, and you have to take care of the preiterations that the tangent vectors require to align themselves. You also have to take care not to `integrate` for too long to avoid a numerical overflow in the tangent vectors. Estimates for the Lyapunov vectors are returned as well. .. automodule:: double_fhn_lyapunov .. literalinclude:: ../examples/double_fhn_lyapunov.py :emphasize-lines: 20-21,26,28-36 :start-after: example-start :dedent: 1 :linenos: .. _transversal: Calculating transversal Lyapunov exponents with `jitcode_transversal_lyap` -------------------------------------------------------------------------- `jitcode_transversal_lyap` is a variant of `jitcode_lyap` that calculates Lyapunov exponents transversal to a user-defined synchronisation manifold. It automatically conflates the differential equations of a group of synchronised components into one equation on the synchronisation manifold. Moreover, it transforms the equations for the tangent vectors such that the tangent vector is automatically orthogonal to the synchronisation manifold. If you are interested in the mathematical details, please read `the accompanying paper`_. .. automodule:: double_fhn_transversal_lyap Other things you may want to do ------------------------------- * `JiTCDDE `_ is like JiTCODE for delay differential equations. * `JiTCSDE `_ is like JiTCODE for stochastics differential equations. * If you want to have input that cannot be expressed in a simple function, you can use `jitcdde_input `_ or use a callback (see the next point). * If you want to call a Python function within the derivative, use the `callback_functions` argument. * The `zvode` integrator and complex arithmetics in general do not work, as they are not that easily implementable in C. Command reference ----------------- .. automodule:: _jitcode :members: :exclude-members: jitcode, jitcode_lyap, jitcode_restricted_lyap, jitcode_transversal_lyap .. autoclass:: integrator_tools.UnsuccessfulIntegration The main class ^^^^^^^^^^^^^^ .. autoclass:: jitcode :members: :inherited-members: Lyapunov exponents ^^^^^^^^^^^^^^^^^^ .. autoclass:: jitcode_lyap :members: .. autoclass:: jitcode_transversal_lyap :members: .. autoclass:: jitcode_restricted_lyap :members: .. _reference: References ---------- .. _notes on compilation: https://jitcde-common.readthedocs.io .. _common JiTC*DE documentation: https://jitcde-common.readthedocs.io .. _SymPy vs. SymEngine: https://jitcde-common.readthedocs.io/#sympy-vs-symengine .. [BGGS80] G. Benettin, L. Galgani, A. Giorgilli, and J.-M. Strelcyn: Lyapunov Characteristic Exponents for smooth dynamical systems and for Hamiltonian systems; A method for computing all of them. Meccanica 15, pp. 9–30 (1980), `10.1007/BF02128236 `_. .. _ODE: http://docs.scipy.org/doc/scipy/reference/generated/scipy.integrate.ode.html .. _Solve IVP: http://docs.scipy.org/doc/scipy/reference/generated/scipy.integrate.solve_ivp.html .. _SymPy Issue 8997: https://github.com/sympy/sympy/issues/8997 .. _SymEngine: https://github.com/symengine/symengine .. _SymPy: http://www.sympy.org/ .. _the accompanying paper: http://arxiv.org/abs/1711.09886