.. automodule:: nbodypy Nbody Class ***************************************** This module provides the basic functions to build objects for the simulation of a N-body problem. Build a Nbody instance ====================== To build a `Nbody` instance, we simply call the constructor. .. code-block:: python :caption: Simple constructor :linenos: import nbodypy system = nbodypy.Nbody() Here the documentation of the constructor. .. automethod:: nbodypy.Nbody.__init__ Let us show a more complete example. .. code-block:: python :caption: Complete example :linenos: import nbodypy import numpy zinit = numpy.array([0.165649854848924, 3.0549363634996e-151 , 0.209141794901608, 8.86909514289322e-29 , 2.03225904939571, -9.31938237331704e-29, -0.0828249274244622 , 0.160070428973807, -0.104570897450804, -2.04906075007578 , -1.01612952469786 , 1.15731474026243, -0.0828249274244618 , -0.160070428973807 , -0.104570897450804, 2.04906075007578, -1.01612952469785, -1.15731474026243]) Arotating = numpy.array([[0.0, 5.13873206001209, 0.0], [-5.13873206001209, 0.0, 0.0], [0.0,0.0, 0.0]]) System = nbodypy.Nbody(N=3, dim=6, init = zinit, rotating=Arotating) With a Json File ---------------- For the same example, one can build an instance of Nbody class with a Json file. .. code-block:: json :caption: `example.json` :linenos: { "N" : "3", "dim" : "6", "init" : "[0.165649854848924, 3.0549363634996e-151 , 0.209141794901608,8.86909514289322e-29 , 2.03225904939571, -9.31938237331704e-29,-0.0828249274244622 , 0.160070428973807, -0.104570897450804,-2.04906075007578 , -1.01612952469786 , 1.15731474026243,-0.0828249274244618 , -0.160070428973807 , -0.104570897450804,2.04906075007578, -1.01612952469785, -1.15731474026243]", "rotating" : "[[0.0, 5.13873206001209, 0.0],[-5.13873206001209, 0.0, 0.0],[0.0,0.0, 0.0]]" } .. code-block:: python :caption: Complete example :linenos: import nbodypy import numpy System = nbodypy.Nbody(jsonFile="example.json") Access to the Parameters ========================= Get the Number of Bodies ------------------------ .. automethod:: nbodypy.Nbody.get_N .. code-block:: python :linenos: import nbodypy system = nbodypy.Nbody() Nbodies = system.get_N() Get the Dimension ----------------- .. automethod:: nbodypy.Nbody.get_dim The dimension parametrer corresponds to the dimension of the phase space for *one* body. .. code-block:: python :linenos: import nbodypy system = nbodypy.Nbody() Nbodies = system.get_dim() Get the Current State --------------------- .. automethod:: nbodypy.Nbody.get_z To get the *current* state of the system, a `numpy.array` of dimension the number of bodies times the dimension. .. code-block:: python :linenos: import nbodypy system = nbodypy.Nbody() Nbodies = system.get_z() One can simply use `system.z`. Get the Position and the Velocity --------------------------------- .. automethod:: nbodypy.Nbody.get_r .. automethod:: nbodypy.Nbody.get_v We provide two methods to get the position and the velocity of one chosen body. Of course, this depends on the dimension of the system. .. code-block:: python :linenos: import nbodypy import numpy zinit = numpy.array([0, 0.316046051487446, 0.864077674448677, 0, 0.100941149773621, -0.158023025743723, -0.432038837224337, 2.02906597904133, -0.100941149773621 , -0.158023025743723, -0.43203883722434, -2.02906597904133]) system = nbodypy.Nbody(N=3,init=zinit) R1 = system.get_r(0) # get the position of the first body V3 = system.get_v(2) # get the velocity of the third body Integrate a solution ==================== The dynamics ------------ This part of the `Nbody` class depends on the function: .. autofunction:: nbodypy.dynamics.system Normally, you do not have to use this function, but if you want to, you have to import it as follows: .. code-block:: python import nbodypy.dynamics The `integrate` method ------------------------- The method to integrate a solution is: .. automethod:: nbodypy.Nbody.integrate For instance, the following exemple: .. code-block:: python :caption: Example of integration :linenos: import nbodypy import numpy zinit = numpy.array([0, 0.316046051487446, 0.864077674448677, 0, 0.100941149773621, -0.158023025743723, -0.432038837224337, 2.02906597904133, -0.100941149773621 , -0.158023025743723, -0.43203883722434, -2.02906597904133]) system = nbodypy.Nbody(N=3,init=zinit) t = numpy.linspace(0, 1.0/3, 100) solution = system.integrate(t) Save the solution in a file ########################### One can save the integration of a solution in a file. To do that, one just has to use the optional parameter `fileName`. .. code-block:: python :caption: Example of integration :linenos: import nbodypy import numpy zinit = numpy.array([0, 0.316046051487446, 0.864077674448677, 0, 0.100941149773621, -0.158023025743723, -0.432038837224337, 2.02906597904133, -0.100941149773621 , -0.158023025743723, -0.43203883722434, -2.02906597904133]) system = nbodypy.Nbody(N=3,init=zinit) t = numpy.linspace(0, 1.0/3, 100) solution = system.integrate(t,fileName="exportSol.dat") This can be use without getting the solution in a variable: .. code-block:: python :caption: Example of integration :linenos: system.integrate(t,fileName="exportSol.dat") The file contains in its first column all the value of the time at each step, and the others contain the phase space states of each body. For example, in the planar case for :math:`N` bodies: .. math:: t_i, x_{i,1}^1, x_{i,1}^2, v_{i,1}^1, v_{i,1}^2,\dots, x_{i,N}^1, x_{i,N}^2, v_{i,N}^1, v_{i,N}^2 Move the Initial Condition -------------------------- You can also "move" the initial position `system.z`. With the `move` membre you can integrate the system during a time `t`, and modify the current state `system.z`. .. automethod:: nbodypy.Nbody.move .. code-block:: python :caption: Example of moving :linenos: import nbodypy import numpy zinit = numpy.array([0, 0.316046051487446, 0.864077674448677, 0, 0.100941149773621, -0.158023025743723, -0.432038837224337, 2.02906597904133, -0.100941149773621 , -0.158023025743723, -0.43203883722434, -2.02906597904133]) system = nbodypy.Nbody(N=3,init=zinit) system.move(0.3) The Linearized Dynamics ======================= This part of the `Nbody` class depends on the function: .. autofunction:: nbodypy.dynamics.Dsystem Normally, you do not have to use this function, but if you want to, you have to import it as follows: .. code-block:: python import nbodypy.dynamics The `Dintegrate` Method ------------------------- The method to integrate a solution and the linearized system around this solution is: .. automethod:: nbodypy.Nbody.Dintegrate For instance, the following exemple: .. code-block:: python :caption: Example of Dintegration :linenos: import nbodypy import numpy zinit = numpy.array([0, 0.316046051487446, 0.864077674448677, 0, 0.100941149773621, -0.158023025743723, -0.432038837224337, 2.02906597904133, -0.100941149773621 , -0.158023025743723, -0.43203883722434, -2.02906597904133]) system = nbodypy.Nbody(N=3,init=zinit) t = numpy.linspace(0, 1.0/3, 100) e1 = numpy.zeros(3*4, dtype=float) e1[0] = 1.0 system.zdz = numpy.hstack((system.z, e1)) solution = system.Dintegrate(t) Compute the Monodromy Matrix ---------------------------- The monodromy matrix gives us information about the influence of perturbations around the periodic orbit. Let us consider a periodic solution :math:`x^{*}(\cdot)` of period :math:`T`, the flow from the initial point :math:`x_{0}^{*}=x^{*}(0)` is denoted by :math:`\varphi(t,x_{0}^{*})`. The matrix .. math:: M=\dfrac{\partial\varphi(T,x_{0}^{*})}{\partial x_{0}} determines whether initial perturbations :math:`x_{0}^{*}` from the periodic orbit :math:`x^{*}(\cdot)` decay or grow. This matrix is called the *monodromy* matrix. Integrating the linearized system ################################# One way to compute it is to solve a certain dynamical system. Let us consider the dynamical system :math:`\dot x=f(x)` and a periodic solution :math:`x^{*}`. The monodromy matrix :math:`M` is the matrix :math:`\Phi(T)` where :math:`\Phi` is a solution of the linearized system .. math:: \left\{ \begin{array}{l} \dfrac{\partial\Phi}{\partial t}=\dfrac{\partial f(x^*)}{\partial x}\Phi,\\ \Phi(0)=Id. \end{array} \right. The function to compute the monodromy matrix following this way is: .. automethod:: nbodypy.Nbody.monodromy This can be used as follow: .. code-block:: python :caption: Example of computing the monodromy matrix :linenos: import nbodypy import numpy zinit = numpy.array([0, 0.316046051487446, 0.864077674448677, 0, 0.100941149773621, -0.158023025743723, -0.432038837224337, 2.02906597904133, -0.100941149773621 , -0.158023025743723, -0.43203883722434, -2.02906597904133]) system = nbodypy.Nbody(N=3,init=zinit) M = system.monodromy(1.0) By Finite Differences ##################### There is another function to compute the monodromy matrix using finite differences. It is faster but less precise. .. automethod:: nbodypy.Nbody.monodromyDF .. code-block:: python :caption: Example of computing the monodromy matrix :linenos: import nbodypy import numpy zinit = numpy.array([0, 0.316046051487446, 0.864077674448677, 0, 0.100941149773621, -0.158023025743723, -0.432038837224337, 2.02906597904133, -0.100941149773621 , -0.158023025743723, -0.43203883722434, -2.02906597904133]) system = nbodypy.Nbody(N=3,init=zinit) M = system.monodromyDF(1.0) .. automodule:: nbodypy :members: