# NumbaCS
**Documentation:** [https://numbacs.readthedocs.io/](https://numbacs.readthedocs.io/en/latest/)
**Source Code (MPL-2.0):** [https://github.com/alb3rtjarvis/numbacs](https://github.com/alb3rtjarvis/numbacs)
NumbaCS (Numba Coherent Structures) is a Python package for performing coherent structure analysis in a user-friendly and efficient manner. It implements methods referred to as "coherent structure methods", an umbrella term for any method that can be used to infer or extract Lagrangian and Eulerian coherent structures. These are tools for performing analysis of time-dependent dynamical systems, mainly with a focus on material transport in fluid flows. While this package can be used with any (incompressible) flow, methods are provided which make it simple for the user to take advantage of these tools for large-scale geophysical flows.
## Features
NumbaCS currently implements the following features:
### Diagnostics:
* **FTLE** -- finite time Lyapunov exponent
* **iLE** -- instantaneous Lyapunov exponent
* **LAVD** -- Lagrangian averaged vorticity deviation
* **IVD** -- instantaneous vorticity deviation
### Feature extraction:
* **LCS** -- hyperbolic and elliptic Lagrangian coherent structures
* **OECS** -- hyperbolic and elliptic objective Eulerian coherent structures
* **Ridge extraction** -- FTLE and iLE ridge extraction
## Installation
Conda:
```
conda install -c conda-forge numbacs
```
Pip:
```
pip install numbacs
```
**_NOTE:_** It is strongly recommended to use conda for installation
due to some reported issues when using pip for one of the dependencies
(see issues [#20](https://github.com/Nicholaswogan/numbalsoda/issues/20)
and [#29](https://github.com/Nicholaswogan/numbalsoda/issues/29) of the
numbalsoda package -- these seem to primarily affect Windows installations).
## Basic usage
Below we show basic demonstrations for computing FTLE. See the [Examples Gallery](https://numbacs.readthedocs.io/en/latest/auto_examples/index.html) for more examples showcasing functionality.
### Predefined flow
```python
import numpy as np
from math import copysign
from numbacs.flows import get_predefined_flow
from numbacs.integration import flowmap_grid_2D
from numbacs.diagnostics import ftle_grid_2D
import matplotlib.pyplot as plt
# set initial time, integration span, and integration direction
t0 = 0.
T = -10.
int_direction = copysign(1, T)
# get ode to be used by 'flowmap_grid_2D'
funcptr, params, domain = get_predefined_flow('double_gyre', int_direction=int_direction)
# set up domain
nx,ny = 401,201
x = np.linspace(domain[0][0], domain[0][1], nx)
y = np.linspace(domain[1][0], domain[1][1], ny)
dx = abs(x[1] - x[0])
dy = abs(y[1] - y[0])
# computes final position of particle trajectories over grid
flowmap = flowmap_grid_2D(funcptr, t0, T, x, y, params)
# compute FTLE over grid
ftle = ftle_grid_2D(flowmap, T, dx, dy)
# plot FTLE
fig, ax = plt.subplots(dpi=200)
ax.contourf(x, y, ftle.T, levels=50)
ax.set_aspect('equal')
plt.show()
```

### Numerical velocity data
```python
import numpy as np
from math import copysign
from numbacs.flows import get_interp_arrays_2D, get_flow_2D
from numbacs.integration import flowmap_grid_2D
from numbacs.diagnostics import ftle_grid_2D
import matplotlib.pyplot as plt
# load in atmospheric data
dates = np.load('../data/merra_june2020/dates.npy')
dt = (dates[1] - dates[0]).astype('timedelta64[h]').astype(int)
t = np.arange(0, len(dates)*dt, dt, np.float64)
lon = np.load('../data/merra_june2020/lon.npy')
lat = np.load('../data/merra_june2020/lat.npy')
# NumbaCS uses 'ij' indexing, most geophysical data uses 'xy'
# indexing for the spatial coordintes. We need to switch axes and
# scale by 3.6 since velocity data is in m/s and we want km/hr.
u = np.moveaxis(np.load('../data/merra_june2020/u_500_800hPa.npy'),1, 2)*3.6
v = np.moveaxis(np.load('../data/merra_june2020/v_500_800hPa.npy'),1, 2)*3.6
nt,nx,ny = u.shape
# set domain on which ftle will be computed
dx = 0.2
dy = 0.2
lonf = np.arange(-100, 35+dx, dx)
latf = np.arange(-5, 45+dy, dy)
# set integration span and integration direction
day = 16
t0_date = np.datetime64("2020-06-{:02d}".format(day))
t0 = t[np.nonzero(dates == t0_date)[0][0]]
T = -72.0
params = np.array([copysign(1, T)])
# get interpolant arrays of velocity field
grid_vel, C_eval_u, C_eval_v = get_interp_arrays_2D(t, lon, lat, u, v)
# retrieve flow and set spherical = 1 since flow is on spherical domain
# and lon is from [-180, 180)
funcptr = get_flow_2D(grid_vel, C_eval_u, C_eval_v, spherical=1)
# compute final position of particle trajectories over grid
flowmap = flowmap_grid_2D(funcptr, t0, T, lonf, latf, params)
# compute FTLE over grid
ftle = ftle_grid_2D(flowmap, T, dx, dy)
# plot coastlines and FTLE
coastlines = np.load('../data/merra_june2020/coastlines.npy')
fig,ax = plt.subplots(dpi=200)
ax.scatter(
coastlines[:,0], coastlines[:,1], 1, 'k', marker='.', edgecolors=None, linewidths=0
)
ax.contourf(lonf, latf, ftle.T, levels=80, zorder=0)
ax.set_xlim([lonf[0], lonf[-1]])
ax.set_ylim([latf[0], latf[-1]])
ax.set_aspect('equal')
plt.show()
```

## Key dependencies
NumbaCS is built on top of three main packages: [Numba](https://numba.pydata.org), [numbalsoda](https://github.com/Nicholaswogan/numbalsoda), and [interpolation](https://www.econforge.org/interpolation.py/). Numba is a JIT-compiler for Python array and numerical functions that generates optimized machine code just-in-time (using the LLVM compiler library) to significantly speed up numerical operations in Python. Numbalsoda is a Python wrapper to ODE solvers in both C++ (LSODA) and FORTRAN (DOP853) that is compatible with Numba (standard Python ODE solvers cannot be used within Numba functions) and bypasses the Python interpreter, speeding up the most computationally expensive piece of Lagrangian coherent structure methods, particle integration. The interpolation package is a Numba-compatible package which optimizes interpolation in Python. It is used in NumbaCS to generate JIT-compiled interpolant functions of numerical velocity fields which can be fed into solvers from the numbalsoda package, resulting in efficient particle integration for flows defined by numerical velocity data. All of this interfacing, which is done behind the scenes through modules the user can call in a straightforward manner, is how NumbaCS achieves impressive speeds while maintaining a relatively simple user experience. We are grateful to the creators, maintainers, and contributors of each of these packages, as well as the other dependencies which NumbaCS relies on ([NumPy](https://numpy.org/), [SciPy](https://scipy.org/), and [ContourPy](https://contourpy.readthedocs.io/en/v1.3.0/)).
## Roadmap
Future releases aim to extend certain methods to higher dimensions, implement new features that should be straightforward within this framework (shape coherent sets, lobe dynamics, etc.), and further streamline and optimize the process for large-scale geophysical applications.
## Contributing
Pull requests are welcome. For major changes, please open an issue first
to discuss what you would like to change.
Please make sure to update tests as appropriate. See the [Contributing guide](https://numbacs.readthedocs.io/en/latest/contributing.html)
for more details.
## Similar software
This section lists similar packages, their functionality, what ODE solvers are available (and what language they are implemented in), and the available interpolation routines. For performance comparisons of *some* packages on core functionality, see the [Benchmarks](https://github.com/alb3rtjarvis/coherent_benchmarks) repository.
---
[`Lagrangian`](https://lagrangian.readthedocs.io/en/latest/index.html) --
Python package for computing FSLE, FTLE, and eigenvectors of Cauchy-Green tensor with a
focus on geophysical flows. Only works with NetCDF files for velocity data.
Largely written in C++ with pybind11 used for binding, resulting in fast
runtimes.
- *Features*: FTLE, FSLE, Cauchy Green eigenvectors
- *Integration*: RK4 (C++)
- *Interpolation*: Linear
[`Dynlab`](https://github.com/hokiepete/dynlab) --
Object oriented Python package
which computes Lagrangian and Eulerian diagnostics along with ridge extraction.
Provides a large collection of predefined flows and is very user friendly.
- *Features*: FTLE, iLE, Trajectory repulsion rate, FTLE ridge extraction, iLE ridge extraction
- *Integration*: LSODA (Python)
- *Interpolation*: N/A
[`TBarrier`](https://github.com/haller-group/TBarrier) --
Collection of Jupyter
notebooks accompanying the book *Transport Barriers and Coherent Structures
in Flow Data -- Advective, Diffusive, Stochastic, and Active methods by George
Haller*. Python code which implements a wide variety of Lagrangian
and Eulerian diagnostics and extraction methods for a variety of different
transport settings, in both 2 and 3 dimensions (NumbaCS currently only
implements purely advective methods in 2D).
- *Features* (both 2D and 3D): FTLE, iLE, variational hyperbolic LCS and OECS, variational elliptic LCS and OECS, variational parabolic LCS, active hyperbolic LCS and OECS, active elliptic LCS and OECS, DBS, diffusive and stochastic elliptic LCS and OECS
- *Integration*: RK4 (Python)
- *Interpolation*: Linear in time, cubic in space
[`Newman`](https://github.com/RossDynamics/Newmanv3.1) --
Fast C++ code for computing FTLE which works with geophysical flows and storm
tracking. Velocity data must be raw binary or ASCII format. No longer maintained.
- *Features* (both 2D and 3D): FTLE, Cauchy Green eigenvectors
- *Integration*: RK45, RK4, Euler (C++)
- *Interpolation*: Linear
[`Aquila-LCS`](https://github.com/ChristianLagares/Aquila-LCS) --
Python code designed to compute FTLE for high-speed turbulent boundary layers in 3D.
Utilizes Numba to implement GPU and CPU versions of the code for fast runtimes.
- *Features* (both 2D and 3D): FTLE, FSLE
- *Integration*: Euler (Python/Numba)
- *Interpolation*: Linear
[`CoherentStructures.jl`](https://coherentstructures.github.io/CoherentStructures.jl/stable/) --
Julia toolbox for computing LCS/FTCS in aperiodic flows. Implements elliptic
LCS methods, FEM-based methods (FEM approximation of dynamic Laplacian for FTCS
extraction), and Graph Laplacian-based methods (spectral clustering and
diffusion maps for coherent sets).
- *Features*: FTLE, finite time coherent sets (via both FEM approximation of dynamic Laplacian and graph Laplacian-based methods), variational elliptic LCS and OECS, diffusive and stochastic elliptic LCS and OECS
- *Integration*: DifferentialEquations.jl, a very advanced and efficient suite of DE solvers (Julia)
- *Interpolation*: Linear, cubic, B-Spline
[`LCS Tool`](https://github.com/haller-group/LCStool) --
MATLAB code used to compute elliptic LCS, hyperbolic LCS, and FTLE.
- *Features*: FTLE, variational hyperbolic LCS, variational elliptic LCS
- *Integration*: ode45 - based off of RK5(4) due to Dormand and Prince (MATLAB)
- *Interpolation*: Linear, cubic
[`LCS MATLAB Kit`](https://dabirilab.com/software) --
MATLAB GUI for computing FTLE from a time series of 2D velocity data. Has FTLE
implementation for intertial particles as well (iFTLE).
- *Features*: FTLE, iFTLE
- *Integration*: Version 1.0 -- RK4 (MATLAB), Version 2.3 -- Euler (MATLAB)
- *Interpolation* Version 1.0 -- Cubic, Version 2.3 -- Linear
[`NumbaCS`](https://numbacs.readthedocs.io/en/latest/) --
Numba accelerated Python package which efficiently computes a variety of
coherent structure methods.
- *Features*: FTLE, iLE, FTLE ridge extraction, variational hyperbolic LCS and OECS, LAVD-based elliptic LCS, IVD-based elliptic OECS, flow map composition
- *Integration*: DOP853 (FORTRAN), LSODA (C++)
- *Interpolation*: Linear, cubic
## Acknowledgments
This work was partially supported by the National Science Foundation (NSF) under
grant number 1821145 and the National Aeronautics and Space Administration
(NASA) under grant number 80NSSC20K1532 issued through the Interdisciplinary
Research in Earth Science (IDS) and Biological Diversity & Ecological
Conservation programs.
Raw data
{
"_id": null,
"home_page": null,
"name": "numbacs",
"maintainer": null,
"docs_url": null,
"requires_python": "<3.12,>=3.10",
"maintainer_email": null,
"keywords": "applied-mathematics, dynamical-systems, fluid-dynamics, ftle, lagrangian-coherent-structures, lcs, numba",
"author": null,
"author_email": "alb3rtjarvis <ajarvis@vt.edu>",
"download_url": "https://files.pythonhosted.org/packages/8c/95/51bd62cc5c3f691993bc2a1d5d323ada0f348942ff2d929fd05c6f9f0b06/numbacs-0.1.2.tar.gz",
"platform": null,
"description": "# NumbaCS\n\n**Documentation:** [https://numbacs.readthedocs.io/](https://numbacs.readthedocs.io/en/latest/)\n\n**Source Code (MPL-2.0):** [https://github.com/alb3rtjarvis/numbacs](https://github.com/alb3rtjarvis/numbacs)\n\nNumbaCS (Numba Coherent Structures) is a Python package for performing coherent structure analysis in a user-friendly and efficient manner. It implements methods referred to as \"coherent structure methods\", an umbrella term for any method that can be used to infer or extract Lagrangian and Eulerian coherent structures. These are tools for performing analysis of time-dependent dynamical systems, mainly with a focus on material transport in fluid flows. While this package can be used with any (incompressible) flow, methods are provided which make it simple for the user to take advantage of these tools for large-scale geophysical flows.\n\n## Features\n\nNumbaCS currently implements the following features:\n\n### Diagnostics:\n\n* **FTLE** -- finite time Lyapunov exponent\n* **iLE** -- instantaneous Lyapunov exponent\n* **LAVD** -- Lagrangian averaged vorticity deviation\n* **IVD** -- instantaneous vorticity deviation\n\n### Feature extraction:\n\n* **LCS** -- hyperbolic and elliptic Lagrangian coherent structures\n* **OECS** -- hyperbolic and elliptic objective Eulerian coherent structures\n* **Ridge extraction** -- FTLE and iLE ridge extraction\n\n## Installation\n\nConda:\n```\nconda install -c conda-forge numbacs\n```\nPip:\n```\npip install numbacs\n```\n\n**_NOTE:_** It is strongly recommended to use conda for installation\ndue to some reported issues when using pip for one of the dependencies\n(see issues [#20](https://github.com/Nicholaswogan/numbalsoda/issues/20)\nand [#29](https://github.com/Nicholaswogan/numbalsoda/issues/29) of the\nnumbalsoda package -- these seem to primarily affect Windows installations).\n\n## Basic usage\n\nBelow we show basic demonstrations for computing FTLE. See the [Examples Gallery](https://numbacs.readthedocs.io/en/latest/auto_examples/index.html) for more examples showcasing functionality.\n\n### Predefined flow\n\n```python\nimport numpy as np\nfrom math import copysign\nfrom numbacs.flows import get_predefined_flow\nfrom numbacs.integration import flowmap_grid_2D\nfrom numbacs.diagnostics import ftle_grid_2D\nimport matplotlib.pyplot as plt\n\n# set initial time, integration span, and integration direction\nt0 = 0.\nT = -10.\nint_direction = copysign(1, T)\n\n# get ode to be used by 'flowmap_grid_2D'\nfuncptr, params, domain = get_predefined_flow('double_gyre', int_direction=int_direction)\n\n# set up domain\nnx,ny = 401,201\nx = np.linspace(domain[0][0], domain[0][1], nx)\ny = np.linspace(domain[1][0], domain[1][1], ny)\ndx = abs(x[1] - x[0])\ndy = abs(y[1] - y[0])\n\n# computes final position of particle trajectories over grid\nflowmap = flowmap_grid_2D(funcptr, t0, T, x, y, params)\n\n# compute FTLE over grid\nftle = ftle_grid_2D(flowmap, T, dx, dy)\n\n# plot FTLE\nfig, ax = plt.subplots(dpi=200)\nax.contourf(x, y, ftle.T, levels=50)\nax.set_aspect('equal')\nplt.show()\n```\n\n\n### Numerical velocity data\n\n```python\nimport numpy as np\nfrom math import copysign\nfrom numbacs.flows import get_interp_arrays_2D, get_flow_2D\nfrom numbacs.integration import flowmap_grid_2D\nfrom numbacs.diagnostics import ftle_grid_2D\nimport matplotlib.pyplot as plt\n\n# load in atmospheric data\ndates = np.load('../data/merra_june2020/dates.npy')\ndt = (dates[1] - dates[0]).astype('timedelta64[h]').astype(int)\nt = np.arange(0, len(dates)*dt, dt, np.float64)\nlon = np.load('../data/merra_june2020/lon.npy')\nlat = np.load('../data/merra_june2020/lat.npy')\n\n# NumbaCS uses 'ij' indexing, most geophysical data uses 'xy'\n# indexing for the spatial coordintes. We need to switch axes and\n# scale by 3.6 since velocity data is in m/s and we want km/hr.\nu = np.moveaxis(np.load('../data/merra_june2020/u_500_800hPa.npy'),1, 2)*3.6\nv = np.moveaxis(np.load('../data/merra_june2020/v_500_800hPa.npy'),1, 2)*3.6\nnt,nx,ny = u.shape\n\n# set domain on which ftle will be computed\ndx = 0.2\ndy = 0.2\nlonf = np.arange(-100, 35+dx, dx)\nlatf = np.arange(-5, 45+dy, dy)\n\n# set integration span and integration direction\nday = 16\nt0_date = np.datetime64(\"2020-06-{:02d}\".format(day))\nt0 = t[np.nonzero(dates == t0_date)[0][0]]\nT = -72.0\nparams = np.array([copysign(1, T)])\n\n# get interpolant arrays of velocity field\ngrid_vel, C_eval_u, C_eval_v = get_interp_arrays_2D(t, lon, lat, u, v)\n\n# retrieve flow and set spherical = 1 since flow is on spherical domain\n# and lon is from [-180, 180)\nfuncptr = get_flow_2D(grid_vel, C_eval_u, C_eval_v, spherical=1)\n\n# compute final position of particle trajectories over grid\nflowmap = flowmap_grid_2D(funcptr, t0, T, lonf, latf, params)\n\n# compute FTLE over grid\nftle = ftle_grid_2D(flowmap, T, dx, dy)\n\n# plot coastlines and FTLE\ncoastlines = np.load('../data/merra_june2020/coastlines.npy')\nfig,ax = plt.subplots(dpi=200)\nax.scatter(\n coastlines[:,0], coastlines[:,1], 1, 'k', marker='.', edgecolors=None, linewidths=0\n)\nax.contourf(lonf, latf, ftle.T, levels=80, zorder=0)\nax.set_xlim([lonf[0], lonf[-1]])\nax.set_ylim([latf[0], latf[-1]])\nax.set_aspect('equal')\nplt.show()\n```\n\n\n\n## Key dependencies\n\nNumbaCS is built on top of three main packages: [Numba](https://numba.pydata.org), [numbalsoda](https://github.com/Nicholaswogan/numbalsoda), and [interpolation](https://www.econforge.org/interpolation.py/). Numba is a JIT-compiler for Python array and numerical functions that generates optimized machine code just-in-time (using the LLVM compiler library) to significantly speed up numerical operations in Python. Numbalsoda is a Python wrapper to ODE solvers in both C++ (LSODA) and FORTRAN (DOP853) that is compatible with Numba (standard Python ODE solvers cannot be used within Numba functions) and bypasses the Python interpreter, speeding up the most computationally expensive piece of Lagrangian coherent structure methods, particle integration. The interpolation package is a Numba-compatible package which optimizes interpolation in Python. It is used in NumbaCS to generate JIT-compiled interpolant functions of numerical velocity fields which can be fed into solvers from the numbalsoda package, resulting in efficient particle integration for flows defined by numerical velocity data. All of this interfacing, which is done behind the scenes through modules the user can call in a straightforward manner, is how NumbaCS achieves impressive speeds while maintaining a relatively simple user experience. We are grateful to the creators, maintainers, and contributors of each of these packages, as well as the other dependencies which NumbaCS relies on ([NumPy](https://numpy.org/), [SciPy](https://scipy.org/), and [ContourPy](https://contourpy.readthedocs.io/en/v1.3.0/)).\n\n## Roadmap\n\nFuture releases aim to extend certain methods to higher dimensions, implement new features that should be straightforward within this framework (shape coherent sets, lobe dynamics, etc.), and further streamline and optimize the process for large-scale geophysical applications.\n\n## Contributing\n\nPull requests are welcome. For major changes, please open an issue first\nto discuss what you would like to change.\n\nPlease make sure to update tests as appropriate. See the [Contributing guide](https://numbacs.readthedocs.io/en/latest/contributing.html)\nfor more details.\n\n## Similar software\n\nThis section lists similar packages, their functionality, what ODE solvers are available (and what language they are implemented in), and the available interpolation routines. For performance comparisons of *some* packages on core functionality, see the [Benchmarks](https://github.com/alb3rtjarvis/coherent_benchmarks) repository.\n\n---\n\n[`Lagrangian`](https://lagrangian.readthedocs.io/en/latest/index.html) --\nPython package for computing FSLE, FTLE, and eigenvectors of Cauchy-Green tensor with a\nfocus on geophysical flows. Only works with NetCDF files for velocity data.\nLargely written in C++ with pybind11 used for binding, resulting in fast\nruntimes.\n- *Features*: FTLE, FSLE, Cauchy Green eigenvectors\n- *Integration*: RK4 (C++)\n- *Interpolation*: Linear\n\n[`Dynlab`](https://github.com/hokiepete/dynlab) --\nObject oriented Python package\nwhich computes Lagrangian and Eulerian diagnostics along with ridge extraction.\nProvides a large collection of predefined flows and is very user friendly.\n- *Features*: FTLE, iLE, Trajectory repulsion rate, FTLE ridge extraction, iLE ridge extraction\n- *Integration*: LSODA (Python)\n- *Interpolation*: N/A\n\n[`TBarrier`](https://github.com/haller-group/TBarrier) --\nCollection of Jupyter\nnotebooks accompanying the book *Transport Barriers and Coherent Structures\nin Flow Data -- Advective, Diffusive, Stochastic, and Active methods by George\nHaller*. Python code which implements a wide variety of Lagrangian\nand Eulerian diagnostics and extraction methods for a variety of different\ntransport settings, in both 2 and 3 dimensions (NumbaCS currently only\nimplements purely advective methods in 2D).\n- *Features* (both 2D and 3D): FTLE, iLE, variational hyperbolic LCS and OECS, variational elliptic LCS and OECS, variational parabolic LCS, active hyperbolic LCS and OECS, active elliptic LCS and OECS, DBS, diffusive and stochastic elliptic LCS and OECS\n- *Integration*: RK4 (Python)\n- *Interpolation*: Linear in time, cubic in space\n\n[`Newman`](https://github.com/RossDynamics/Newmanv3.1) --\nFast C++ code for computing FTLE which works with geophysical flows and storm\ntracking. Velocity data must be raw binary or ASCII format. No longer maintained.\n- *Features* (both 2D and 3D): FTLE, Cauchy Green eigenvectors\n- *Integration*: RK45, RK4, Euler (C++)\n- *Interpolation*: Linear\n\n[`Aquila-LCS`](https://github.com/ChristianLagares/Aquila-LCS) --\nPython code designed to compute FTLE for high-speed turbulent boundary layers in 3D.\nUtilizes Numba to implement GPU and CPU versions of the code for fast runtimes.\n- *Features* (both 2D and 3D): FTLE, FSLE\n- *Integration*: Euler (Python/Numba)\n- *Interpolation*: Linear\n\n[`CoherentStructures.jl`](https://coherentstructures.github.io/CoherentStructures.jl/stable/) --\nJulia toolbox for computing LCS/FTCS in aperiodic flows. Implements elliptic\nLCS methods, FEM-based methods (FEM approximation of dynamic Laplacian for FTCS\nextraction), and Graph Laplacian-based methods (spectral clustering and\ndiffusion maps for coherent sets).\n- *Features*: FTLE, finite time coherent sets (via both FEM approximation of dynamic Laplacian and graph Laplacian-based methods), variational elliptic LCS and OECS, diffusive and stochastic elliptic LCS and OECS\n- *Integration*: DifferentialEquations.jl, a very advanced and efficient suite of DE solvers (Julia)\n- *Interpolation*: Linear, cubic, B-Spline\n\n[`LCS Tool`](https://github.com/haller-group/LCStool) --\nMATLAB code used to compute elliptic LCS, hyperbolic LCS, and FTLE.\n- *Features*: FTLE, variational hyperbolic LCS, variational elliptic LCS\n- *Integration*: ode45 - based off of RK5(4) due to Dormand and Prince (MATLAB)\n- *Interpolation*: Linear, cubic\n\n[`LCS MATLAB Kit`](https://dabirilab.com/software) --\nMATLAB GUI for computing FTLE from a time series of 2D velocity data. Has FTLE\nimplementation for intertial particles as well (iFTLE).\n- *Features*: FTLE, iFTLE\n- *Integration*: Version 1.0 -- RK4 (MATLAB), Version 2.3 -- Euler (MATLAB)\n- *Interpolation* Version 1.0 -- Cubic, Version 2.3 -- Linear\n\n[`NumbaCS`](https://numbacs.readthedocs.io/en/latest/) --\nNumba accelerated Python package which efficiently computes a variety of\ncoherent structure methods.\n- *Features*: FTLE, iLE, FTLE ridge extraction, variational hyperbolic LCS and OECS, LAVD-based elliptic LCS, IVD-based elliptic OECS, flow map composition\n- *Integration*: DOP853 (FORTRAN), LSODA (C++)\n- *Interpolation*: Linear, cubic\n\n## Acknowledgments\n\nThis work was partially supported by the National Science Foundation (NSF) under\ngrant number 1821145 and the National Aeronautics and Space Administration\n(NASA) under grant number 80NSSC20K1532 issued through the Interdisciplinary\nResearch in Earth Science (IDS) and Biological Diversity & Ecological\nConservation programs.\n",
"bugtrack_url": null,
"license": null,
"summary": "Numba-accelerated coherent structure package.",
"version": "0.1.2",
"project_urls": {
"Documentation": "https://numbacs.readthedocs.io/",
"Issues": "https://github.com/alb3rtjarvis/numbacs/issues",
"Source": "https://github.com/alb3rtjarvis/numbacs"
},
"split_keywords": [
"applied-mathematics",
" dynamical-systems",
" fluid-dynamics",
" ftle",
" lagrangian-coherent-structures",
" lcs",
" numba"
],
"urls": [
{
"comment_text": null,
"digests": {
"blake2b_256": "b96074faee854a5a09aadfaa3898d9fae11219b5751d672f7609ad21ac7971b9",
"md5": "126c5ffd6491bdb2a03d2ea8df479c46",
"sha256": "10efbcaa80d7d0e4bdeb08d2b926169d9bba5b48b498a79ef2a9fbfd96756151"
},
"downloads": -1,
"filename": "numbacs-0.1.2-py3-none-any.whl",
"has_sig": false,
"md5_digest": "126c5ffd6491bdb2a03d2ea8df479c46",
"packagetype": "bdist_wheel",
"python_version": "py3",
"requires_python": "<3.12,>=3.10",
"size": 45658,
"upload_time": "2025-08-16T18:22:27",
"upload_time_iso_8601": "2025-08-16T18:22:27.386142Z",
"url": "https://files.pythonhosted.org/packages/b9/60/74faee854a5a09aadfaa3898d9fae11219b5751d672f7609ad21ac7971b9/numbacs-0.1.2-py3-none-any.whl",
"yanked": false,
"yanked_reason": null
},
{
"comment_text": null,
"digests": {
"blake2b_256": "8c9551bd62cc5c3f691993bc2a1d5d323ada0f348942ff2d929fd05c6f9f0b06",
"md5": "c01b373f3339dd7946bee19214ddb214",
"sha256": "6f03cd777d26ce145e01994494c370d509d790431c26b7bffcdfd4a7ac88bac6"
},
"downloads": -1,
"filename": "numbacs-0.1.2.tar.gz",
"has_sig": false,
"md5_digest": "c01b373f3339dd7946bee19214ddb214",
"packagetype": "sdist",
"python_version": "source",
"requires_python": "<3.12,>=3.10",
"size": 122411,
"upload_time": "2025-08-16T18:22:28",
"upload_time_iso_8601": "2025-08-16T18:22:28.446968Z",
"url": "https://files.pythonhosted.org/packages/8c/95/51bd62cc5c3f691993bc2a1d5d323ada0f348942ff2d929fd05c6f9f0b06/numbacs-0.1.2.tar.gz",
"yanked": false,
"yanked_reason": null
}
],
"upload_time": "2025-08-16 18:22:28",
"github": true,
"gitlab": false,
"bitbucket": false,
"codeberg": false,
"github_user": "alb3rtjarvis",
"github_project": "numbacs",
"travis_ci": false,
"coveralls": true,
"github_actions": true,
"lcname": "numbacs"
}