Legend:
Library
Module
Module type
Parameter
Class
Class type
N-dimensional array module: including creation, manipulation, and various vectorised mathematical operations.
About the comparison of two complex numbers x and y, Owl uses the following conventions: 1) x and y are equal iff both real and imaginary parts are equal; 2) x is less than y if the magnitude of x is less than the magnitude of y; in case both x and y have the same magnitudes, x is less than y if the phase of x is less than the phase of y; 3) less or equal, greater, greater or equal relation can be further defined atop of the aforementioned conventions.
The generic module supports operations for the following Bigarry element types: Int8_signed, Int8_unsigned, Int16_signed, Int16_unsigned, Int32, Int64, Float32, Float64, Complex32, Complex64.
Type definition
type('a, 'b) t = ('a, 'b, Stdlib.Bigarray.c_layout)Stdlib.Bigarray.Genarray.t
N-dimensional array type, i.e. Bigarray Genarray type.
type('a, 'b) kind = ('a, 'b)Stdlib.Bigarray.kind
Type of the ndarray, e.g., Bigarray.Float32, Bigarray.Complex64, and etc.
empty Bigarray.Float64 [|3;4;5|] creates a three diemensional array of Bigarray.Float64 type. Each dimension has the following size: 3, 4, and 5. The elements in the array are not initialised, they can be any value. empty is faster than zeros to create a ndarray.
The module only supports the following four types of ndarray: Bigarray.Float32, Bigarray.Float64, Bigarray.Complex32, and Bigarray.Complex64.
val create : ('a, 'b)kind->int array->'a->('a, 'b)t
create Bigarray.Float64 [|3;4;5|] 2. creates a three-diemensional array of Bigarray.Float64 type. Each dimension has the following size: 3, 4, and 5. The elements in the array are initialised to 2.
val init : ('a, 'b)kind->int array->(int ->'a)->('a, 'b)t
init Bigarray.Float64 d f creates a ndarray x of shape d, then using f to initialise the elements in x. The input of f is 1-dimensional index of the ndarray. You need to explicitly convert it if you need N-dimensional index. The function ind can help you.
val init_nd : ('a, 'b)kind->int array->(int array->'a)->('a, 'b)t
init_nd is almost the same as init but f receives n-dimensional index as input. It is more convenient since you don't have to convert the index by yourself, but this also means init_nd is slower than init.
zeros Bigarray.Complex32 [|3;4;5|] creates a three-diemensional array of Bigarray.Complex32 type. Each dimension has the following size: 3, 4, and 5. The elements in the array are initialised to "zero". Depending on the kind, zero can be 0. or Complex.zero.
ones Bigarray.Complex32 [|3;4;5|] creates a three-diemensional array of Bigarray.Complex32 type. Each dimension has the following size: 3, 4, and 5. The elements in the array are initialised to "one". Depending on the kind, one can be 1. or Complex.one.
val uniform : ('a, 'b)kind->?a:'a->?b:'a->int array->('a, 'b)t
uniform Bigarray.Float64 [|3;4;5|] creates a three-diemensional array of type Bigarray.Float64. Each dimension has the following size: 3, 4, and 5. The elements in the array follow a uniform distribution 0,1.
val gaussian : ('a, 'b)kind->?mu:'a->?sigma:'a->int array->('a, 'b)t
gaussian Float64 [|3;4;5|] ...
val poisson : ('a, 'b)kind->mu:float ->int array->('a, 'b)t
poisson Float64 [|3;4;5|] ...
val sequential : ('a, 'b)kind->?a:'a->?step:'a->int array->('a, 'b)t
sequential Bigarray.Float64 [|3;4;5|] 2. creates a three-diemensional array of type Bigarray.Float64. Each dimension has the following size: 3, 4, and 5. The elements in the array are assigned sequential values.
?a specifies the starting value and the default value is zero; whilst ?step specifies the step size with default value one.
val linspace : ('a, 'b)kind->'a->'a->int ->('a, 'b)t
linspace k 0. 9. 10 ...
val logspace : ('a, 'b)kind->?base:float ->'a->'a->int ->('a, 'b)t
logspace k 0. 9. 10 ...
val bernoulli : ('a, 'b)kind->?p:float ->int array->('a, 'b)t
bernoulli k ~p:0.3 [|2;3;4|]
val complex :
('a, 'b)kind->('c, 'd)kind->('a, 'b)t->('a, 'b)t->('c, 'd)t
complex re im constructs a complex ndarray/matrix from re and im. re and im contain the real and imaginary part of x respectively.
Note that both re and im can be complex but must have same type. The real part of re will be the real part of x and the imaginary part of im will be the imaginary part of x.
val polar :
('a, 'b)kind->('c, 'd)kind->('a, 'b)t->('a, 'b)t->('c, 'd)t
complex rho theta constructs a complex ndarray/matrix from polar coordinates rho and theta. rho contains the magnitudes and theta contains phase angles. Note that the behaviour is undefined if rho has negative elelments or theta has infinity elelments.
val unit_basis : ('a, 'b)kind->int ->int ->('a, 'b)t
unit_basis k n i returns a unit basis vector with ith element set to 1.
same_data x y checks whether x and y share the same underlying data in the memory. Namely, both variables point to the same memory address. This is done by checking the Data pointer in the Bigarray structure.
This function is very useful for avoiding unnecessary copying between two ndarrays especially if one has been reshaped or sliced.
kind x returns the type of ndarray x. It is one of the four possible values: Bigarray.Float32, Bigarray.Float64, Bigarray.Complex32, and Bigarray.Complex64.
val get_index : ('a, 'b)t->int array array->'a array
get_index i x returns an array of element values specified by the indices i. The length of array i equals the number of dimensions of x. The arrays in i must have the same length, and each represents the indices in that dimension.
E.g., [| [|1;2|]; [|3;4|] |] returns the value of elements at position (1,3) and (2,4) respectively.
val set_index : ('a, 'b)t->int array array->'a array-> unit
set_index i x a sets the value of elements in x according to the indices specified by i. The length of array i equals the number of dimensions of x. The arrays in i must have the same length, and each represents the indices in that dimension.
If the length of a equals to the length of i, then each element will be assigned by the value in the corresponding position in x. If the length of a equals to one, then all the elements will be assigned the same value.
val get_fancy : Owl_types.index list->('a, 'b)t->('a, 'b)t
get_fancy s x returns a copy of the slice in x. The slice is defined by a which is an int option array. E.g., for a ndarray x of dimension [|2; 2; 3|], slice [0] x takes the following slices of index \(0,*,*\), i.e., [|0;0;0|], [|0;0;1|], [|0;0;2|] ... Also note that if the length of s is less than the number of dimensions of x, slice function will append slice definition to higher diemensions by assuming all the elements in missing dimensions will be taken.
Basically, slice function offers very much the same semantic as that in numpy, i.e., start:stop:step grammar, so if you how to index and slice ndarray in numpy, you should not find it difficult to use this function. Please just refer to numpy documentation or my tutorial.
There are two differences between slice_left and slice: slice_left does not make a copy but simply moving the pointer; slice_left can only make a slice from left-most axis whereas slice is much more flexible and can work on arbitrary axis which need not start from left-most side.
val set_fancy : Owl_types.index list->('a, 'b)t->('a, 'b)t-> unit
set_fancy axis x y set the slice defined by axis in x according to the values in y. y must have the same shape as the one defined by axis.
About the slice definition of axis, please refer to get_fancy function.
val get_fancy_ext : Owl_types.index array->('a, 'b)t->('a, 'b)t
This function is used for extended indexing operator since ocaml 4.10.0. The indexing and slicing syntax become much ligher.
val set_fancy_ext : Owl_types.index array->('a, 'b)t->('a, 'b)t-> unit
This function is used for extended indexing operator since ocaml 4.10.0. The indexing and slicing syntax become much ligher.
val get_slice : int list list->('a, 'b)t->('a, 'b)t
get_slice axis x aims to provide a simpler version of get_fancy. This function assumes that every list element in the passed in int list list represents a range, i.e., R constructor.
E.g., [[];[0;3];[0]] is equivalent to [R []; R [0;3]; R [0]].
val set_slice : int list list->('a, 'b)t->('a, 'b)t-> unit
set_slice axis x y aims to provide a simpler version of set_fancy. This function assumes that every list element in the passed in int list list represents a range, i.e., R constructor.
E.g., [[];[0;3];[0]] is equivalent to [R []; R [0;3]; R [0]].
val get_slice_ext : int list array->('a, 'b)t->('a, 'b)t
get_slice_ext axis x is used for extended indexing operator since ocaml 4.10.0. The indexing and slicing syntax become much ligher.
E.g., x.%{0;1;2}.
val set_slice_ext : int list array->('a, 'b)t->('a, 'b)t-> unit
Similar to get_slice_ext axis x, this function is used for extended indexing operator since ocaml 4.10.0. The indexing and slicing syntax become much ligher.
Some as Bigarray.sub_left, please refer to Bigarray documentation.
val sub_ndarray : int array->('a, 'b)t->('a, 'b)t array
sub_ndarray parts x is similar to Bigarray.sub_left. It splits the passed in ndarray x along the axis 0 according to parts. The elelments in parts do not need to be equal but they must sum up to the dimension along axis zero.
The returned sub-ndarrays share the same memory as x. Because there is no copies made, this function is much faster than using `split` function to divide the lowest dimensionality of x.
val resize : ?head:bool ->('a, 'b)t->int array->('a, 'b)t
resize ~head x d resizes the ndarray x. If there are less number of elelments in the new shape than the old one, the new ndarray shares part of the memory with the old x. head indicates the alignment between the new and old data, either from head or from tail. Note the data is flattened before the operation.
If there are more elements in the new shape d. Then new memory space will be allocated and the content of x will be copied to the new memory. The rest of the allocated space will be filled with zeros. The default value of head is true.
reshape x d transforms x into a new shape definted by d. Note the reshape function will not make a copy of x, the returned ndarray shares the same memory with the original x.
One shape dimension (only one) can be set to -1. In this case, the value is inferred from the length of the array and remaining dimensions.
rotate x d rotates x clockwise d degrees. d must be multiple times of 90, otherwise the function will fail. If x is an n-dimensional array, then the function rotates the plane formed by the first and second dimensions.
val transpose : ?axis:int array->('a, 'b)t->('a, 'b)t
transpose ~axis x makes a copy of x, then transpose it according to ~axis. ~axis must be a valid permutation of x dimension indices. E.g., for a three-dimensional ndarray, it can be [2;1;0], [0;2;1], [1;2;0], and etc.
tile x a tiles the data in x according the repetition specified by a. This function provides the exact behaviour as numpy.tile, please refer to the numpy's online documentation for details.
repeat x a repeats the elements of x according the repetition specified by a. The i-th element of a specifies the number of times that the individual entries of the i-th dimension of x should be repeated.
val concat_vertical : ('a, 'b)t->('a, 'b)t->('a, 'b)t
concat_vertical x y concatenates two ndarray x and y vertically. This is just a convenient function for concatenating two ndarrays along their lowest dimension, i.e. 0.
The associated operator is @||, please refer to :doc:`owl_operator`.
val concat_horizontal : ('a, 'b)t->('a, 'b)t->('a, 'b)t
concat_horizontal x y concatenates two ndarrays x and y horizontally. This is just a convenient function for concatenating two ndarrays along their highest dimension.
The associated operator is @=, please refer to :doc:`owl_operator`.
concat_vh is used to assemble small parts of matrices into a bigger one. E.g. In [| [|a; b; c|]; [|d; e; f|]; [|g; h; i|] |], wherein `a, b, c ... i` are matrices of different shapes. They will be concatenated into a big matrix as follows.
.. math:: \beginmatrix a & b & c \\ d & e & f \\ g & h & i \endmatrix
This is achieved by first concatenating along axis:1 for each element in the array, then concatenating along axis:0. The number of elements in each array needs not to be equal as long as the aggregated dimensions match. E.g., please check the following example.
.. code-block:: ocaml
let a00 = Mat.sequential 2 3 in let a01 = Mat.sequential 2 2 in let a02 = Mat.sequential 2 1 in let a10 = Mat.sequential 3 3 in let a11 = Mat.sequential 3 3 in Mat.concat_vh | [|a00; a01; a02|]; [|a10; a11|] |;;
val concatenate : ?axis:int ->('a, 'b)t array->('a, 'b)t
concatenate ~axis:2 x concatenates an array of ndarrays along the third dimension. For the ndarrays in x, they must have the same shape except the dimension specified by axis. The default value of axis is 0, i.e., the lowest dimension of a matrix/ndarray.
val stack : ?axis:int ->('a, 'b)t array->('a, 'b)t
stack ~axis x stacks an array of ndarrays along the axis dimension. For example, if x contains K ndarrays of shape |2;3|, then stack ~axis:1 x will return an ndarray of dimensions |2;K;3|. The ndarrays in x, they must all have the same shape. The default value of axis is 0.
val split : ?axis:int ->int array->('a, 'b)t->('a, 'b)t array
split ~axis parts x splits an ndarray x into parts along the specified axis. This function is the inverse operation of concatenate. The elements in x must sum up to the dimension in the specified axis.
split_vh parts x splits a passed in ndarray x along the first two dimensions, i.e. axis 0 and axis 1. This is the inverse operation of concat_vh function, and the function is very useful in dividing a big matrix into smaller (especially heterogeneous) parts.
For example, given a matrix x of shape [|8;10|], it is possible to split in the following ways.
val squeeze : ?axis:int array->('a, 'b)t->('a, 'b)t
squeeze ~axis x removes single-dimensional entries from the shape of x.
val expand : ?hi:bool ->('a, 'b)t->int ->('a, 'b)t
expand x d reshapes x by increasing its rank from num_dims x to d. The opposite operation is squeeze x. The hi parameter is used to specify whether the expandsion is along high dimension (by setting true), or along the low dimension (by setting false). The default value is false.
val pad : ?v:'a->int list list->('a, 'b)t->('a, 'b)t
pad ~v p x pads a ndarray x with a constant value v. The padding index p is a list of lists of 2 integers. These two integers denote padding width at both edges of one dimension of x.
top x n returns the indices of n greatest values of x. The indices are arranged according to the corresponding element values, from the greatest one to the smallest one.
bottom x n returns the indices of n smallest values of x. The indices are arranged according to the corresponding element values, from the smallest one to the greatest one.
sort x performs quicksort of the elelments in x. A new copy is returned as result, the original x remains intact. If you want to perform in-place sorting, please use `sort_` instead.
val argsort : ('a, 'b)t->(int64, Stdlib.Bigarray.int64_elt)t
argsort x returns the indices with which the elements in x are sorted in increasing order. Note that the returned index ndarray has the same shape as that of x, and the indices are 1D indices.
val draw : ?axis:int ->('a, 'b)t->int ->('a, 'b)t * int array
draw ~axis x n draws n samples from x along the specified axis, with replacement. axis is set to zero by default. The return is a tuple of both samples and the indices of the selected samples.
val mmap :
Unix.file_descr ->?pos:int64 ->('a, 'b)kind->bool ->int array->('a, 'b)t
map f x is similar to mapi f x except the index is not passed.
val foldi :
?axis:int ->(int ->'a->'a->'a)->'a->('a, 'b)t->('a, 'b)t
foldi ~axis f a x folds (or reduces) the elements in x from left along the specified axis using passed in function f. a is the initial element and in f i acc bacc is the accumulater and b is one of the elements in x along the same axis. Note that i is 1d index of b.
val fold : ?axis:int ->('a->'a->'a)->'a->('a, 'b)t->('a, 'b)t
Similar to foldi, except that the index of an element is not passed to f.
val scani : ?axis:int ->(int ->'a->'a->'a)->('a, 'b)t->('a, 'b)t
scan ~axis f x scans the x along the specified axis using passed in function f. f acc a b returns an updated acc which will be passed in the next call to f i acc a. This function can be used to implement accumulative operations such as sum and prod functions. Note that the i is 1d index of a in x.
val scan : ?axis:int ->('a->'a->'a)->('a, 'b)t->('a, 'b)t
Similar to scani, except that the index of an element is not passed to f.
val filteri : (int ->'a-> bool)->('a, 'b)t->int array
filteri f x uses f to filter out certain elements in x. An element will be included if f returns true. The returned result is an array of 1-dimensional indices of the selected elements. To obtain the n-dimensional indices, you need to convert it manually with Owl's helper function.
Similar to filteri, but the indices are not passed to f.
val iter2i : (int ->'a->'b-> unit)->('a, 'c)t->('b, 'd)t-> unit
Similar to iteri but applies to two N-dimensional arrays x and y. Both x and y must have the same shape.
val iter2 : ('a->'b-> unit)->('a, 'c)t->('b, 'd)t-> unit
Similar to iter2i, except that the index not passed to f.
val map2i : (int ->'a->'a->'a)->('a, 'b)t->('a, 'b)t->('a, 'b)t
map2i f x y applies f to two elements of the same position in both x and y. Note that 1d index is passed to function f.
val map2 : ('a->'a->'a)->('a, 'b)t->('a, 'b)t->('a, 'b)t
map2 f x y is similar to map2i f x y except the index is not passed.
val iteri_nd : (int array->'a-> unit)->('a, 'b)t-> unit
Similar to iteri but n-d indices are passed to the user function.
val mapi_nd : (int array->'a->'a)->('a, 'b)t->('a, 'b)t
Similar to mapi but n-d indices are passed to the user function.
val foldi_nd :
?axis:int ->(int array->'a->'a->'a)->'a->('a, 'b)t->('a, 'b)t
Similar to foldi but n-d indices are passed to the user function.
val scani_nd :
?axis:int ->(int array->'a->'a->'a)->('a, 'b)t->('a, 'b)t
Similar to scani but n-d indices are passed to the user function.
val filteri_nd : (int array->'a-> bool)->('a, 'b)t->int array array
Similar to filteri but n-d indices are returned.
val iter2i_nd :
(int array->'a->'c-> unit)->('a, 'b)t->('c, 'd)t->
unit
Similar to iter2i but n-d indices are passed to the user function.
val map2i_nd :
(int array->'a->'a->'a)->('a, 'b)t->('a, 'b)t->('a, 'b)t
Similar to map2i but n-d indices are passed to the user function.
val iteri_slice :
?axis:int ->(int ->('a, 'b)t-> unit)->('a, 'b)t->
unit
iteri_slice ~axis f x iterates the slices along the specified axis in x and applies the function f. The 1-d index of of the slice is passed in. By default, the axis is 0. Setting axis to the highest dimension is not allowed because in that case you can just use `iteri` to iterate all the elements in x which is more efficient.
Note that the slice is obtained by slicing left (due to Owl's C-layout ndarray) a sub-array out of x. E.g., if x has shape [|3;4;5|], setting axis=0 will iterate three 4 x 5 matrices. The slice shares the same memory with x so no copy is made.
val iter_slice : ?axis:int ->(('a, 'b)t-> unit)->('a, 'b)t-> unit
Similar to iteri_slice but slice index is not passed in.
val mapi_slice :
?axis:int ->(int ->('a, 'b)t->'c)->('a, 'b)t->'c array
mapi_slice ~axis f x maps the slices along the specified axis in x and applies the function f. By default, axis is 0. The index of of the slice is passed in.
Please refer to iteri_slice for more details.
val map_slice : ?axis:int ->(('a, 'b)t->'c)->('a, 'b)t->'c array
Similar to mapi_slice but slice index is not passed in.
elt_equal x y performs element-wise = comparison of x and y. Assume that a is from x and b is the corresponding element of a from y of the same position. The function returns another binary (0 and 1) ndarray/matrix wherein 1 indicates a = b.
The function supports broadcast operation.
val elt_not_equal : ('a, 'b)t->('a, 'b)t->('a, 'b)t
elt_not_equal x y performs element-wise != comparison of x and y. Assume that a is from x and b is the corresponding element of a from y of the same position. The function returns another binary (0 and 1) ndarray/matrix wherein 1 indicates a <> b.
elt_less x y performs element-wise < comparison of x and y. Assume that a is from x and b is the corresponding element of a from y of the same position. The function returns another binary (0 and 1) ndarray/matrix wherein 1 indicates a < b.
elt_greater x y performs element-wise > comparison of x and y. Assume that a is from x and b is the corresponding element of a from y of the same position. The function returns another binary (0 and 1) ndarray/matrix wherein 1 indicates a > b.
The function supports broadcast operation.
val elt_less_equal : ('a, 'b)t->('a, 'b)t->('a, 'b)t
elt_less_equal x y performs element-wise <= comparison of x and y. Assume that a is from x and b is the corresponding element of a from y of the same position. The function returns another binary (0 and 1) ndarray/matrix wherein 1 indicates a <= b.
The function supports broadcast operation.
val elt_greater_equal : ('a, 'b)t->('a, 'b)t->('a, 'b)t
elt_greater_equal x y performs element-wise >= comparison of x and y. Assume that a is from x and b is the corresponding element of a from y of the same position. The function returns another binary (0 and 1) ndarray/matrix wherein 1 indicates a >= b.
elt_equal_scalar x a performs element-wise = comparison of x and a. Assume that b is one element from x The function returns another binary (0 and 1) ndarray/matrix wherein 1 of the corresponding position indicates a = b, otherwise 0.
val elt_not_equal_scalar : ('a, 'b)t->'a->('a, 'b)t
elt_not_equal_scalar x a performs element-wise != comparison of x and a. Assume that b is one element from x The function returns another binary (0 and 1) ndarray/matrix wherein 1 of the corresponding position indicates a <> b, otherwise 0.
elt_less_scalar x a performs element-wise < comparison of x and a. Assume that b is one element from x The function returns another binary (0 and 1) ndarray/matrix wherein 1 of the corresponding position indicates a < b, otherwise 0.
elt_greater_scalar x a performs element-wise > comparison of x and a. Assume that b is one element from x The function returns another binary (0 and 1) ndarray/matrix wherein 1 of the corresponding position indicates a > b, otherwise 0.
val elt_less_equal_scalar : ('a, 'b)t->'a->('a, 'b)t
elt_less_equal_scalar x a performs element-wise <= comparison of x and a. Assume that b is one element from x The function returns another binary (0 and 1) ndarray/matrix wherein 1 of the corresponding position indicates a <= b, otherwise 0.
val elt_greater_equal_scalar : ('a, 'b)t->'a->('a, 'b)t
elt_greater_equal_scalar x a performs element-wise >= comparison of x and a. Assume that b is one element from x The function returns another binary (0 and 1) ndarray/matrix wherein 1 of the corresponding position indicates a >= b, otherwise 0.
val approx_equal : ?eps:float ->('a, 'b)t->('a, 'b)t-> bool
approx_equal ~eps x y returns true if x and y are approximately equal, i.e., for any two elements a from x and b from y, we have abs (a - b) < eps. For complex numbers, the eps applies to both real and imaginary part.
Note: the threshold check is exclusive for passed in eps, i.e., the threshold interval is (a-eps, a+eps).
val approx_equal_scalar : ?eps:float ->('a, 'b)t->'a-> bool
approx_equal_scalar ~eps x a returns true all the elements in x are approximately equal to a, i.e., abs (x - a) < eps. For complex numbers, the eps applies to both real and imaginary part.
Note: the threshold check is exclusive for the passed in eps.
val approx_elt_equal : ?eps:float ->('a, 'b)t->('a, 'b)t->('a, 'b)t
approx_elt_equal ~eps x y compares the element-wise equality of x and y, then returns another binary (i.e., 0 and 1) ndarray/matrix wherein 1 indicates that two corresponding elements a from x and b from y are considered as approximately equal, namely abs (a - b) < eps.
val approx_elt_equal_scalar : ?eps:float ->('a, 'b)t->'a->('a, 'b)t
approx_elt_equal_scalar ~eps x a compares all the elements of x to a scalar value a, then returns another binary (i.e., 0 and 1) ndarray/matrix wherein 1 indicates that the element b from x is considered as approximately equal to a, namely abs (a - b) < eps.
Input/Output functions
val of_array : ('a, 'b)kind->'a array->int array->('a, 'b)t
of_array k x d takes an array x and converts it into an ndarray of type k and shape d.
to_array x converts an ndarray x to OCaml's array type. Note that the ndarray x is flattened before conversion.
val print :
?max_row:int ->?max_col:int ->?header:bool ->?fmt:('a-> string)->('a, 'b)t->
unit
print x prints all the elements in x as well as their indices. max_row and max_col specify the maximum number of rows and columns to display. header specifies whether or not to print out the headers. fmt is the function to format every element into string.
val pp_dsnda : Stdlib.Format.formatter ->('a, 'b)t-> unit
pp_dsnda x prints x in OCaml toplevel. If the ndarray is too long, pp_dsnda only prints out parts of the ndarray.
load_npy file load a npy file into a matrix of type k. If the matrix is in the file is not of type k, fails with [file]: incorrect format. This function is implemented using npy-ocaml https://github.com/LaurentMazare/npy-ocaml.
Unary math operators
val re_c2s :
(Stdlib.Complex.t, Stdlib.Bigarray.complex32_elt)t->(float, Stdlib.Bigarray.float32_elt)t
re_c2s x returns all the real components of x in a new ndarray of same shape.
val re_z2d :
(Stdlib.Complex.t, Stdlib.Bigarray.complex64_elt)t->(float, Stdlib.Bigarray.float64_elt)t
re_d2z x returns all the real components of x in a new ndarray of same shape.
val im_c2s :
(Stdlib.Complex.t, Stdlib.Bigarray.complex32_elt)t->(float, Stdlib.Bigarray.float32_elt)t
im_c2s x returns all the imaginary components of x in a new ndarray of same shape.
val im_z2d :
(Stdlib.Complex.t, Stdlib.Bigarray.complex64_elt)t->(float, Stdlib.Bigarray.float64_elt)t
im_d2z x returns all the imaginary components of x in a new ndarray of same shape.
val sum : ?axis:int ->?keep_dims:bool ->('a, 'b)t->('a, 'b)t
sum ~axis x sums the elements in x along specified axis.
sem' x calculates the standard error of mean of all the elements in x.
val min : ?axis:int ->?keep_dims:bool ->('a, 'b)t->('a, 'b)t
min x returns the minimum of all elements in x along specified axis. If no axis is specified, x will be flattened and the minimum of all the elements will be returned. For two complex numbers, the one with the smaller magnitude will be selected. If two magnitudes are the same, the one with the smaller phase will be selected.
min' x is similar to min but returns the minimum of all elements in x in scalar value.
val max : ?axis:int ->?keep_dims:bool ->('a, 'b)t->('a, 'b)t
max x returns the maximum of all elements in x along specified axis. If no axis is specified, x will be flattened and the maximum of all the elements will be returned. For two complex numbers, the one with the greater magnitude will be selected. If two magnitudes are the same, the one with the greater phase will be selected.
max_i x returns the maximum of all elements in x as well as its index.
val minmax_i : ('a, 'b)t->('a * int array) * ('a * int array)
minmax_i x returns ((min_v,min_i), (max_v,max_i)) where (min_v,min_i) is the minimum value in x along with its index while (max_v,max_i) is the maximum value along its index.
reci_tol ~tol x computes the reciprocal of every element in x. Different from reci, reci_tol sets the elements whose abs value smaller than tol to zeros. If tol is not specified, the default Owl_utils.eps Float32 will be used. For complex numbers, refer to Owl's doc to see how to compare.
fix x rounds each element of x to the nearest integer toward zero. For positive elements, the behavior is the same as floor. For negative ones, the behavior is the same as ceil.
modf x performs modf over all the elements in x, the fractal part is saved in the first element of the returned tuple whereas the integer part is saved in the second element.
relu x computes the rectified linear unit function max(x, 0) of the elements in x and returns the result in a new ndarray.
val elu : ?alpha:float ->(float, 'a)t->(float, 'a)t
elu alpha x computes the exponential linear unit function x >= 0. ? x : (alpha * (exp(x) - 1)) of the elements in x and returns the result in a new ndarray.
val leaky_relu : ?alpha:float ->(float, 'a)t->(float, 'a)t
leaky_relu alpha x computes the leaky rectified linear unit function x >= 0. ? x : (alpha * x) of the elements in x and returns the result in a new ndarray.
softsign x computes the softsign function x / (1 + abs(x)) of the elements in x and returns the result in a new ndarray.
val softmax : ?axis:int ->(float, 'a)t->(float, 'a)t
softmax x computes the softmax functions (exp x) / (sum (exp x)) of all the elements along the specified axis in x and returns the result in a new ndarray.
By default, axis = -1, i.e. along the highest dimension.
l2norm_sqr x calculates the square of l2-norm (or l2norm, Euclidean norm) of all elements in x. The function uses conjugate transpose in the product, hence it always returns a float number.
val vecnorm :
?axis:int ->?p:float ->?keep_dims:bool ->('a, 'b)t->('a, 'b)t
vecnorm ~axis ~p x calculates the generalised vector p-norm along the specified axis. The generalised p-norm is defined as below.
cumsum ~axis x : performs cumulative sum of the elements along the given axis ~axis. If ~axis is None, then the cumsum is performed along the lowest dimension. The returned result however always remains the same shape.
cummax ~axis x : performs cumulative max along axis dimension.
val diff : ?axis:int ->?n:int ->('a, 'b)t->('a, 'b)t
diff ~axis ~n x calculates the n-th difference of x along the specified axis.
Parameters: * axis: axis to calculate the difference. The default value is the highest dimension. * n: how many times to calculate the difference. The default value is 1.
Return: * The difference ndarray y. Note that the shape of y 1 less than that of x along specified axis.
val angle : (Stdlib.Complex.t, 'a)t->(Stdlib.Complex.t, 'a)t
angle x calculates the phase angle of all complex numbers in x.
val proj : (Stdlib.Complex.t, 'a)t->(Stdlib.Complex.t, 'a)t
proj x computes the projection on Riemann sphere of all elelments in x.
add x y adds all the elements in x and y elementwise, and returns the result in a new ndarray.
General broadcast operation is automatically applied to add/sub/mul/div, etc. The function compares the dimension element-wise from the highest to the lowest with the following broadcast rules (same as numpy): 1. equal; 2. either is 1.
ssqr x a computes the sum of squared differences of all the elements in x from constant a. This function only computes the square of each element rather than the conjugate transpose as l2norm_sqr does.
ssqr_diff x y computes the sum of squared differences of every elements in x and its corresponding element in y.
val cross_entropy' : (float, 'a)t->(float, 'a)t-> float
cross_entropy x y calculates the cross entropy between x and y using base e.
val clip_by_value : ?amin:'a->?amax:'a->('a, 'b)t->('a, 'b)t
clip_by_value ~amin ~amax x clips the elements in x based on amin and amax. The elements smaller than amin will be set to amin, and the elements greater than amax will be set to amax.
clip_by_l2norm t x clips the x according to the threshold set by t.
val fma : ('a, 'b)t->('a, 'b)t->('a, 'b)t->('a, 'b)t
fma x y z calculates the `fused multiply add`, i.e. (x * y) + z.
Tensor Calculus
val contract1 : (int * int) array->('a, 'b)t->('a, 'b)t
contract1 index_pairs x performs indices contraction (a.k.a tensor contraction) on x. index_pairs is an array of contracted indices.
Caveat: Not well tested yet, use with care! Also, consider to use TTGT in future for better performance.
val contract2 : (int * int) array->('a, 'b)t->('a, 'b)t->('a, 'b)t
contract2 index_pairs x y performs indices contraction (a.k.a tensor contraction) on two ndarrays x and y. index_pairs is an array of contracted indices, the first element is the index of x, the second is that of y.
Caveat: Not well tested yet, use with care! Also, consider to use TTGT in future for better performance.
cast kind x casts x of type ('c, 'd) t to type ('a, 'b) t specify by the passed in kind parameter. This function is a generalisation of the other casting functions such as cast_s2d, cast_c2z, and etc.
val cast_s2d :
(float, Stdlib.Bigarray.float32_elt)t->(float, Stdlib.Bigarray.float64_elt)t
cast_s2d x casts x from float32 to float64.
val cast_d2s :
(float, Stdlib.Bigarray.float64_elt)t->(float, Stdlib.Bigarray.float32_elt)t
cast_d2s x casts x from float64 to float32.
val cast_c2z :
(Stdlib.Complex.t, Stdlib.Bigarray.complex32_elt)t->(Stdlib.Complex.t, Stdlib.Bigarray.complex64_elt)t
cast_c2z x casts x from complex32 to complex64.
val cast_z2c :
(Stdlib.Complex.t, Stdlib.Bigarray.complex64_elt)t->(Stdlib.Complex.t, Stdlib.Bigarray.complex32_elt)t
cast_z2c x casts x from complex64 to complex32.
val cast_s2c :
(float, Stdlib.Bigarray.float32_elt)t->(Stdlib.Complex.t, Stdlib.Bigarray.complex32_elt)t
cast_s2c x casts x from float32 to complex32.
val cast_d2z :
(float, Stdlib.Bigarray.float64_elt)t->(Stdlib.Complex.t, Stdlib.Bigarray.complex64_elt)t
cast_d2z x casts x from float64 to complex64.
val cast_s2z :
(float, Stdlib.Bigarray.float32_elt)t->(Stdlib.Complex.t, Stdlib.Bigarray.complex64_elt)t
cast_s2z x casts x from float32 to complex64.
val cast_d2c :
(float, Stdlib.Bigarray.float64_elt)t->(Stdlib.Complex.t, Stdlib.Bigarray.complex32_elt)t
cast_d2c x casts x from float64 to complex32.
Neural network related
val conv1d :
?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t
TODO
val conv2d :
?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t
TODO
val conv3d :
?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t
TODO
val dilated_conv1d :
?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->int array->('a, 'b)t
TODO
val dilated_conv2d :
?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->int array->('a, 'b)t
TODO
val dilated_conv3d :
?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->int array->('a, 'b)t
TODO
val transpose_conv1d :
?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t
TODO
val transpose_conv2d :
?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t
TODO
val transpose_conv3d :
?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t
TODO
val max_pool1d :
?padding:Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t
TODO
val max_pool2d :
?padding:Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t
TODO
val max_pool3d :
?padding:Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t
TODO
val avg_pool1d :
?padding:Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t
TODO
val avg_pool2d :
?padding:Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t
TODO
val avg_pool3d :
?padding:Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t
val upsampling2d : ('a, 'b)t->int array->('a, 'b)t
TODO
val conv1d_backward_input :
('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->('a, 'b)t
TODO
val conv1d_backward_kernel :
('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->('a, 'b)t
TODO
val conv2d_backward_input :
('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->('a, 'b)t
TODO
val conv2d_backward_kernel :
('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->('a, 'b)t
TODO
val conv3d_backward_input :
('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->('a, 'b)t
TODO
val conv3d_backward_kernel :
('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->('a, 'b)t
TODO
val dilated_conv1d_backward_input :
('a, 'b)t->('a, 'b)t->int array->int array->('a, 'b)t->('a, 'b)t
TODO
val dilated_conv1d_backward_kernel :
('a, 'b)t->('a, 'b)t->int array->int array->('a, 'b)t->('a, 'b)t
TODO
val dilated_conv2d_backward_input :
('a, 'b)t->('a, 'b)t->int array->int array->('a, 'b)t->('a, 'b)t
TODO
val dilated_conv2d_backward_kernel :
('a, 'b)t->('a, 'b)t->int array->int array->('a, 'b)t->('a, 'b)t
TODO
val dilated_conv3d_backward_input :
('a, 'b)t->('a, 'b)t->int array->int array->('a, 'b)t->('a, 'b)t
TODO
val dilated_conv3d_backward_kernel :
('a, 'b)t->('a, 'b)t->int array->int array->('a, 'b)t->('a, 'b)t
TODO
val transpose_conv1d_backward_input :
('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->('a, 'b)t
TODO
val transpose_conv1d_backward_kernel :
('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->('a, 'b)t
TODO
val transpose_conv2d_backward_input :
('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->('a, 'b)t
TODO
val transpose_conv2d_backward_kernel :
('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->('a, 'b)t
TODO
val transpose_conv3d_backward_input :
('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->('a, 'b)t
TODO
val transpose_conv3d_backward_kernel :
('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->('a, 'b)t
TODO
val max_pool1d_backward :
Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t->('a, 'b)t
TODO
val max_pool2d_backward :
Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t->('a, 'b)t
TODO
val max_pool3d_backward :
Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t->('a, 'b)t
TODO
val avg_pool1d_backward :
Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t->('a, 'b)t
TODO
val avg_pool2d_backward :
Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t->('a, 'b)t
TODO
val avg_pool3d_backward :
Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t->('a, 'b)t
TODO
val upsampling2d_backward : ('a, 'b)t->int array->('a, 'b)t->('a, 'b)t
TODO
Helper functions
The following functions are helper functions for some other functions in both Ndarray and Ndview modules. In general, you are not supposed to use these functions directly.
one_hot idx depth creates one-hot vectors according to the indices ndarray and the specified depth. If idx is rank N, then the return is rank N+1. More specifically, if idx is of shape [|a;b;c|], the return is of shape [|a;b;c;depth|].
sum_slices ~axis:2 x for x of [|2;3;4;5|], it returns an ndarray of shape [|4;5|]. Currently, the operation is done using gemm, it is fast but consumes more memory.
val slide :
?axis:int ->?ofs:int ->?step:int ->window:int ->('a, 'b)t->('a, 'b)t
slide ~axis ~window x generates a new ndarray by sliding a window along specified axis in x. E.g., if x has shape [|a;b;c|] and axis = 1, then [|a; number of windows; window; c|] is the shape of the returned ndarray.
Parameters: * axis is the axis for sliding, the default is -1, i.e. highest dimension. * ofs is the starting position of the sliding window. The default is 0. * step is the step size, the default is 1. * window is the size of the sliding window.
sigmoid_ x is similar to sigmoid but output is written to x
val softmax_ : ?out:('a, 'b)t->?axis:int ->('a, 'b)t-> unit
softmax_ x is similar to softmax but output is written to x
val cumsum_ : ?out:('a, 'b)t->?axis:int ->('a, 'b)t-> unit
cumsum_ x is similar to cumsum but output is written to x
val cumprod_ : ?out:('a, 'b)t->?axis:int ->('a, 'b)t-> unit
cumprod_ x is similar to cumprod but output is written to x
val cummin_ : ?out:('a, 'b)t->?axis:int ->('a, 'b)t-> unit
cummin_ x is similar to cummin but output is written to x
val cummax_ : ?out:('a, 'b)t->?axis:int ->('a, 'b)t-> unit
cummax_ x is similar to cummax but output is written to x
val dropout_ : ?out:('a, 'b)t->?rate:float ->('a, 'b)t-> unit
dropout_ x is similar to dropout but output is written to x
val elt_equal_ : ?out:('a, 'b)t->('a, 'b)t->('a, 'b)t-> unit
elt_equal_ x y is similar to elt_equal function but the output is written to out. You need to make sure out is big enough to hold the output result.
val elt_not_equal_ : ?out:('a, 'b)t->('a, 'b)t->('a, 'b)t-> unit
elt_not_equal_ x y is similar to elt_not_equal function but the output is written to out. You need to make sure out is big enough to hold the output result.
val elt_less_ : ?out:('a, 'b)t->('a, 'b)t->('a, 'b)t-> unit
elt_less_ x y is similar to elt_less function but the output is written to out. You need to make sure out is big enough to hold the output result.
val elt_greater_ : ?out:('a, 'b)t->('a, 'b)t->('a, 'b)t-> unit
elt_greater_ x y is similar to elt_greater function but the output is written to out. You need to make sure out is big enough to hold the output result.
val elt_less_equal_ : ?out:('a, 'b)t->('a, 'b)t->('a, 'b)t-> unit
elt_less_equal_ x y is similar to elt_less_equal function but the output is written to out. You need to make sure out is big enough to hold the output result.
val elt_greater_equal_ : ?out:('a, 'b)t->('a, 'b)t->('a, 'b)t-> unit
elt_greater_equal_ x y is similar to elt_greater_equal function but the output is written to out. You need to make sure out is big enough to hold the output result.
val elt_equal_scalar_ : ?out:('a, 'b)t->('a, 'b)t->'a-> unit
elt_equal_scalar_ x a is similar to elt_equal_scalar function but the output is written to x.
val elt_not_equal_scalar_ : ?out:('a, 'b)t->('a, 'b)t->'a-> unit
elt_not_equal_scalar_ x a is similar to elt_not_equal_scalar function but the output is written to x.
val elt_less_scalar_ : ?out:('a, 'b)t->('a, 'b)t->'a-> unit
elt_less_scalar_ x a is similar to elt_less_scalar function but the output is written to x.
val elt_greater_scalar_ : ?out:('a, 'b)t->('a, 'b)t->'a-> unit
elt_greater_scalar_ x a is similar to elt_greater_scalar function but the output is written to x.
val elt_less_equal_scalar_ : ?out:('a, 'b)t->('a, 'b)t->'a-> unit
elt_less_equal_scalar_ x a is similar to elt_less_equal_scalar function but the output is written to x.
val elt_greater_equal_scalar_ : ?out:('a, 'b)t->('a, 'b)t->'a-> unit
elt_greater_equal_scalar_ x a is similar to elt_greater_equal_scalar function but the output is written to x.
val conv1d_ :
out:('a, 'b)t->?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->
unit
TODO
val conv2d_ :
out:('a, 'b)t->?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->
unit
TODO
val conv3d_ :
out:('a, 'b)t->?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->
unit
TODO
val dilated_conv1d_ :
out:('a, 'b)t->?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->int array->
unit
TODO
val dilated_conv2d_ :
out:('a, 'b)t->?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->int array->
unit
TODO
val dilated_conv3d_ :
out:('a, 'b)t->?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->int array->
unit
TODO
val transpose_conv1d_ :
out:('a, 'b)t->?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->
unit
TODO
val transpose_conv2d_ :
out:('a, 'b)t->?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->
unit
TODO
val transpose_conv3d_ :
out:('a, 'b)t->?padding:Owl_types.padding ->('a, 'b)t->('a, 'b)t->int array->
unit
TODO
val max_pool1d_ :
out:('a, 'b)t->?padding:Owl_types.padding ->('a, 'b)t->int array->int array->
unit
TODO
val max_pool2d_ :
out:('a, 'b)t->?padding:Owl_types.padding ->('a, 'b)t->int array->int array->
unit
TODO
val max_pool3d_ :
out:('a, 'b)t->?padding:Owl_types.padding ->('a, 'b)t->int array->int array->
unit
TODO
val avg_pool1d_ :
out:('a, 'b)t->?padding:Owl_types.padding ->('a, 'b)t->int array->int array->
unit
TODO
val avg_pool2d_ :
out:('a, 'b)t->?padding:Owl_types.padding ->('a, 'b)t->int array->int array->
unit
TODO
val avg_pool3d_ :
out:('a, 'b)t->?padding:Owl_types.padding ->('a, 'b)t->int array->int array->
unit
TODO
val upsampling2d_ : out:('a, 'b)t->('a, 'b)t->int array-> unit
TODO
val conv1d_backward_input_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->
unit
TODO
val conv1d_backward_kernel_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->
unit
TODO
val conv2d_backward_input_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->
unit
TODO
val conv2d_backward_kernel_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->
unit
TODO
val conv3d_backward_input_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->
unit
TODO
val conv3d_backward_kernel_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->
unit
TODO
val dilated_conv1d_backward_input_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->int array->('a, 'b)t->
unit
TODO
val dilated_conv1d_backward_kernel_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->int array->('a, 'b)t->
unit
TODO
val dilated_conv2d_backward_input_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->int array->('a, 'b)t->
unit
TODO
val dilated_conv2d_backward_kernel_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->int array->('a, 'b)t->
unit
TODO
val dilated_conv3d_backward_input_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->int array->('a, 'b)t->
unit
TODO
val dilated_conv3d_backward_kernel_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->int array->('a, 'b)t->
unit
TODO
val transpose_conv1d_backward_input_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->
unit
TODO
val transpose_conv1d_backward_kernel_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->
unit
TODO
val transpose_conv2d_backward_input_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->
unit
TODO
val transpose_conv2d_backward_kernel_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->
unit
TODO
val transpose_conv3d_backward_input_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->
unit
TODO
val transpose_conv3d_backward_kernel_ :
out:('a, 'b)t->('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->
unit
TODO
val max_pool1d_backward_ :
out:('a, 'b)t->Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t->
unit
TODO
val max_pool2d_backward_ :
out:('a, 'b)t->Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t->
unit
TODO
val max_pool3d_backward_ :
out:('a, 'b)t->Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t->
unit
TODO
val avg_pool1d_backward_ :
out:('a, 'b)t->Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t->
unit
TODO
val avg_pool2d_backward_ :
out:('a, 'b)t->Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t->
unit
TODO
val avg_pool3d_backward_ :
out:('a, 'b)t->Owl_types.padding ->('a, 'b)t->int array->int array->('a, 'b)t->
unit
TODO
val upsampling2d_backward_ :
out:('a, 'b)t->('a, 'b)t->int array->('a, 'b)t->
unit
TODO
val fused_adagrad_ : ?out:('a, 'b)t->rate:'a->eps:'a->('a, 'b)t-> unit