SciPy.jl
A Julia interface for SciPy using PyCall.jl.
You can use many useful scientific functions of SciPy from Julia codes.
You can know which kind of functions are available in the SciPy Reference Guide.
Requirements
Julia 1.4.x or higher.
Install
The package is tested against, and being developed for, Julia 1.4 and above on Linux, macOS, and Windows.
You can install this package with:
using Pkg;Pkg.add("SciPy")
and then just import it with using SciPy
.
If you want use latest development code, clone this repo.
Then, command this in Pkg mode:
Pkg> dev /path/to/module
How to use
You can access each SciPy modules as belows:
scipy.cluster
SciPy.cluster
— Constantscipy.cluster module
Examples
You can do k-means clustering using this module:
julia> features = [[ 1.9 2.3];
[ 1.5 2.5];
[ 0.8 0.6];
[ 0.4 1.8];
[ 0.1 0.1];
[ 0.2 1.8];
[ 2.0 0.5];
[ 0.3 1.5];
[ 1.0 1.0]]
9×2 Array{Float64,2}:
1.9 2.3
1.5 2.5
0.8 0.6
0.4 1.8
0.1 0.1
0.2 1.8
2.0 0.5
0.3 1.5
1.0 1.0
julia> whitened = SciPy.cluster.vq.whiten(features)
9×2 Array{Float64,2}:
2.7396 2.91001
2.16284 3.16306
1.15351 0.759134
0.576757 2.2774
0.144189 0.126522
0.288379 2.2774
2.88379 0.632612
0.432568 1.89784
1.44189 1.26522
julia> SciPy.cluster.vq.kmeans(whitened, [whitened[1,:] whitened[3,:]] )
([1.1174670798453024 1.8345740800894272; 2.8837860125040065 0.6326117517549749], 1.073399
3090584457)
scipy.constants
SciPy.constants
— Constantscipy.constants module
Examples
You can access each constants:
julia> SciPy.constants.golden
1.618033988749895
julia> SciPy.constants.physical_constants["electron mass"]
(9.10938356e-31, "kg", 1.1e-38)
julia> SciPy.constants.convert_temperature([-40, 40.0], "Celsius", "Kelvin")
2-element Array{Float64,1}:
233.14999999999998
313.15
scipy.fft
SciPy.fft
— Constantscipy.fft module
Examples
You can use FFT (Fast Fourier Transform):
julia> SciPy.fft.fft(exp.(π/8 * collect(1:8)))
8-element Array{Complex{Float64},1}:
68.17385416403044 - 0.0im
1.408601300061675 + 31.248171041435185im
-10.268363617931092 + 15.207165888808841im
-12.695027025520982 + 6.493878653648949im
-13.216494113330363 - 0.0im
-12.695027025520982 - 6.493878653648949im
-10.268363617931092 - 15.207165888808841im
1.408601300061675 - 31.248171041435185im
scipy.interpolate
SciPy.interpolate
— Constantscipy.interpolate module
Examples
You can interpolate 1D data:
julia> x = collect(0:10);
julia> y = exp.(-x/3.0);
julia> f = SciPy.interpolate.interp1d(x, y);
julia> f(0.5)
0-dimensional Array{Float64,0}:
0.8582656552868946
scipy.io
SciPy.io
— Constantscipy.io module
Examples
You can save a MATLAB-style .mat file:
julia> mdic = Dict([("a", 100), ("label", "experiment")])
Dict{String,Any} with 2 entries:
"label" => "experiment"
"a" => 100
julia> SciPy.io.savemat("sample_data.mat", mdic)
scipy.linalg
SciPy.linalg
— Constantscipy.linalg module
Examples
You can compute the matrix exponential using Pade approximation:
julia> SciPy.linalg.expm(zeros((2,2)))
2×2 Array{Float64,2}:
1.0 0.0
0.0 1.0
scipy.ndimage
SciPy.ndimage
— Constantscipy.ndimage module
Examples
You can compute a multidimensional convolution:
julia> k = [[1 1 1];[1 1 0];[1 0 0]]
3×3 Array{Int64,2}:
1 1 1
1 1 0
1 0 0
julia> a = [[1 2 0 0];
[5 3 0 4];
[0 0 0 7];
[9 3 0 0]]
4×4 Array{Int64,2}:
1 2 0 0
5 3 0 4
0 0 0 7
9 3 0 0
julia> SciPy.ndimage.convolve(a, k, mode="constant", cval=0.0)
4×4 Array{Int64,2}:
11 10 7 4
10 3 11 11
15 12 14 7
12 3 7 0
scipy.odr
SciPy.odr
— Constantscipy.odr module
Examples
You can calculate orthogonal distance regression with an exponential model:
julia> x = collect(0.0:5.0);
julia> y = -10.0 .+ exp.(0.5*x);
julia> data = SciPy.odr.Data(x, y)
PyObject <scipy.odr.odrpack.Data object at 0x7fe5fda4ccc0>
julia> data = SciPy.odr.Data(x, y);
julia> odr_obj = SciPy.odr.ODR(data, SciPy.odr.exponential);
julia> output = odr_obj.run();
julia> println(output.beta)
[-10.0, 0.5]
scipy.optimize
SciPy.optimize
— Constantscipy.optimize module
Examples
You can optimize the Rosen function:
julia> x0 = [1.3, 0.7, 0.8, 1.9, 1.2];
julia> res = SciPy.optimize.minimize(SciPy.optimize.rosen, x0, method="Nelder-Mead", tol=1e-6)
Dict{Any,Any} with 8 entries:
"fun" => 1.94206e-13
"nit" => 295
"nfev" => 494
"status" => 0
"success" => true
"message" => "Optimization terminated successfully."
"x" => [1.0, 1.0, 1.0, 1.0, 1.0]
"final_simplex" => ([1.0 1.0 … 1.0 1.0; 1.0 1.0 … 1.0 1.0; … ; 1.0 1.0 … 1.0 1.0; 1.0 1.0 … 1.0 1.0], [1…
scipy.signal
SciPy.signal
— Constantscipy.signal module
Examples
You can compute the Kaiser parameter beta, given the attenuation a:
julia> SciPy.signal.kaiser_beta(65)
6.20426
scipy.sparse
SciPy.sparse
— Constantscipy.sparse module
Examples
You can do sparse matrix calculation:
julia> A = SciPy.sparse.csc_matrix([[1.0 0.0];
[1.0 2.0]]);
julia> Ainv = SciPy.sparse.linalg.inv(A);
julia> A.dot(Ainv).todense()
2×2 Array{Float64,2}:
1.0 0.0
0.0 1.0
scipy.spatial
SciPy.spatial
— Constantscipy.spatial module
Examples
You can calculate several rotation representations:
julia> R = SciPy.spatial.transform.Rotation;
julia> r = R.from_quat([0, 0, sin(π/4.0), cos(π/4)]);
julia> r.as_matrix()
3×3 Array{Float64,2}:
2.22045e-16 -1.0 0.0
1.0 2.22045e-16 0.0
0.0 0.0 1.0
julia> r.as_rotvec()
3-element Array{Float64,1}:
0.0
0.0
1.5707963267948963
julia> r.as_euler("zyx", degrees=true)
3-element Array{Float64,1}:
90.0
0.0
0.0
scipy.stats
SciPy.stats
— Constantscipy.stats module
Examples
You can calculate Pearson correlation coefficient and p-value:
julia> a = [0, 0, 0, 1, 1, 1, 1];
julia> b = collect(0:6);
julia> SciPy.stats.pearsonr(a, b)
(0.8660254037844386, 0.011724811003954649)
Accessing docstring
You can access docstring of a SciPy function:
help?> SciPy.io.savemat
Save a dictionary of names and arrays into a MATLAB-style .mat file.
...
Index
SciPy.SciPy
SciPy.cluster
SciPy.constants
SciPy.fft
SciPy.integrate
SciPy.interpolate
SciPy.io
SciPy.linalg
SciPy.ndimage
SciPy.odr
SciPy.optimize
SciPy.signal
SciPy.sparse
SciPy.spatial
SciPy.stats
SciPy.__init__
SciPy.print_configulations
API Reference
SciPy.SciPy
— ModuleA Julia interface module for SciPy
SciPy.integrate
— Constantscipy.integrate module
Examples
You can compute a definite integral:
julia> f(x) = x^2
f (generic function with 1 method)
julia> SciPy.integrate.quad(f, 0, 4)
(21.333333333333336, 2.368475785867001e-13)
SciPy.print_configulations
— MethodPrint configulations:
- Julia version
- Python version
- Python path
- scipy version
SciPy.__init__
— MethodModule initialization function