Skip to main content

ferray_ufunc/ops/
arithmetic.rs

1// ferray-ufunc: Arithmetic functions
2//
3// add, subtract, multiply, divide, true_divide, floor_divide, power,
4// remainder, mod_, fmod, divmod, absolute, fabs, sign, negative, positive,
5// reciprocal, sqrt, cbrt, square, heaviside, gcd, lcm
6//
7// Cumulative: cumsum, cumprod, nancumsum, nancumprod
8// Differences: diff, ediff1d, gradient
9// Products: cross
10// Integration: trapezoid
11//
12// Reduction: add_reduce, add_accumulate, multiply_outer
13
14use ferray_core::Array;
15use ferray_core::dimension::{Dimension, Ix1, IxDyn};
16use ferray_core::dtype::Element;
17use ferray_core::error::{FerrayError, FerrayResult};
18use num_traits::Float;
19
20use crate::helpers::{binary_broadcast_op, binary_float_op, unary_float_op};
21
22// ---------------------------------------------------------------------------
23// Basic arithmetic (binary, same-shape)
24// ---------------------------------------------------------------------------
25
26/// Elementwise addition.
27pub fn add<T, D>(a: &Array<T, D>, b: &Array<T, D>) -> FerrayResult<Array<T, D>>
28where
29    T: Element + std::ops::Add<Output = T> + Copy,
30    D: Dimension,
31{
32    if a.shape() != b.shape() {
33        return Err(FerrayError::shape_mismatch(format!(
34            "add: shapes {:?} and {:?} do not match",
35            a.shape(),
36            b.shape()
37        )));
38    }
39    let data: Vec<T> = a.iter().zip(b.iter()).map(|(&x, &y)| x + y).collect();
40    Array::from_vec(a.dim().clone(), data)
41}
42
43/// Elementwise subtraction.
44pub fn subtract<T, D>(a: &Array<T, D>, b: &Array<T, D>) -> FerrayResult<Array<T, D>>
45where
46    T: Element + std::ops::Sub<Output = T> + Copy,
47    D: Dimension,
48{
49    if a.shape() != b.shape() {
50        return Err(FerrayError::shape_mismatch(format!(
51            "subtract: shapes {:?} and {:?} do not match",
52            a.shape(),
53            b.shape()
54        )));
55    }
56    let data: Vec<T> = a.iter().zip(b.iter()).map(|(&x, &y)| x - y).collect();
57    Array::from_vec(a.dim().clone(), data)
58}
59
60/// Elementwise multiplication.
61pub fn multiply<T, D>(a: &Array<T, D>, b: &Array<T, D>) -> FerrayResult<Array<T, D>>
62where
63    T: Element + std::ops::Mul<Output = T> + Copy,
64    D: Dimension,
65{
66    if a.shape() != b.shape() {
67        return Err(FerrayError::shape_mismatch(format!(
68            "multiply: shapes {:?} and {:?} do not match",
69            a.shape(),
70            b.shape()
71        )));
72    }
73    let data: Vec<T> = a.iter().zip(b.iter()).map(|(&x, &y)| x * y).collect();
74    Array::from_vec(a.dim().clone(), data)
75}
76
77/// Elementwise division.
78pub fn divide<T, D>(a: &Array<T, D>, b: &Array<T, D>) -> FerrayResult<Array<T, D>>
79where
80    T: Element + std::ops::Div<Output = T> + Copy,
81    D: Dimension,
82{
83    if a.shape() != b.shape() {
84        return Err(FerrayError::shape_mismatch(format!(
85            "divide: shapes {:?} and {:?} do not match",
86            a.shape(),
87            b.shape()
88        )));
89    }
90    let data: Vec<T> = a.iter().zip(b.iter()).map(|(&x, &y)| x / y).collect();
91    Array::from_vec(a.dim().clone(), data)
92}
93
94/// Alias for [`divide`] — true division (float).
95pub fn true_divide<T, D>(a: &Array<T, D>, b: &Array<T, D>) -> FerrayResult<Array<T, D>>
96where
97    T: Element + Float,
98    D: Dimension,
99{
100    binary_float_op(a, b, |x, y| x / y)
101}
102
103/// Floor division: floor(a / b).
104pub fn floor_divide<T, D>(a: &Array<T, D>, b: &Array<T, D>) -> FerrayResult<Array<T, D>>
105where
106    T: Element + Float,
107    D: Dimension,
108{
109    binary_float_op(a, b, |x, y| (x / y).floor())
110}
111
112/// Elementwise power: a^b.
113pub fn power<T, D>(a: &Array<T, D>, b: &Array<T, D>) -> FerrayResult<Array<T, D>>
114where
115    T: Element + Float,
116    D: Dimension,
117{
118    binary_float_op(a, b, |x, y| x.powf(y))
119}
120
121/// Elementwise remainder (Python-style modulo).
122pub fn remainder<T, D>(a: &Array<T, D>, b: &Array<T, D>) -> FerrayResult<Array<T, D>>
123where
124    T: Element + Float,
125    D: Dimension,
126{
127    let z = <T as Element>::zero();
128    binary_float_op(a, b, |x, y| {
129        let r = x % y;
130        // Python/NumPy mod: result has same sign as divisor
131        if (r < z && y > z) || (r > z && y < z) {
132            r + y
133        } else {
134            r
135        }
136    })
137}
138
139/// Alias for [`remainder`].
140pub fn mod_<T, D>(a: &Array<T, D>, b: &Array<T, D>) -> FerrayResult<Array<T, D>>
141where
142    T: Element + Float,
143    D: Dimension,
144{
145    remainder(a, b)
146}
147
148/// C-style fmod (remainder has same sign as dividend).
149pub fn fmod<T, D>(a: &Array<T, D>, b: &Array<T, D>) -> FerrayResult<Array<T, D>>
150where
151    T: Element + Float,
152    D: Dimension,
153{
154    binary_float_op(a, b, |x, y| x % y)
155}
156
157/// Return (floor_divide, remainder) as a tuple of arrays.
158pub fn divmod<T, D>(a: &Array<T, D>, b: &Array<T, D>) -> FerrayResult<(Array<T, D>, Array<T, D>)>
159where
160    T: Element + Float,
161    D: Dimension,
162{
163    Ok((floor_divide(a, b)?, remainder(a, b)?))
164}
165
166// ---------------------------------------------------------------------------
167// Unary arithmetic
168// ---------------------------------------------------------------------------
169
170/// Elementwise absolute value.
171///
172/// Uses hardware SIMD for contiguous f64 arrays.
173pub fn absolute<T, D>(input: &Array<T, D>) -> FerrayResult<Array<T, D>>
174where
175    T: Element + Float,
176    D: Dimension,
177{
178    if let Some(r) = crate::helpers::try_simd_f64_unary(input, crate::dispatch::simd_abs_f64) {
179        return r;
180    }
181    unary_float_op(input, T::abs)
182}
183
184/// Alias for [`absolute`] — float abs.
185pub fn fabs<T, D>(input: &Array<T, D>) -> FerrayResult<Array<T, D>>
186where
187    T: Element + Float,
188    D: Dimension,
189{
190    absolute(input)
191}
192
193/// Elementwise sign: -1 for negative, 0 for zero, +1 for positive.
194pub fn sign<T, D>(input: &Array<T, D>) -> FerrayResult<Array<T, D>>
195where
196    T: Element + Float,
197    D: Dimension,
198{
199    unary_float_op(input, |x| {
200        if x.is_nan() {
201            <T as Float>::nan()
202        } else if x > <T as Element>::zero() {
203            <T as Element>::one()
204        } else if x < <T as Element>::zero() {
205            -<T as Element>::one()
206        } else {
207            <T as Element>::zero()
208        }
209    })
210}
211
212/// Elementwise negation.
213///
214/// Uses hardware SIMD for contiguous f64 arrays.
215pub fn negative<T, D>(input: &Array<T, D>) -> FerrayResult<Array<T, D>>
216where
217    T: Element + Float,
218    D: Dimension,
219{
220    if let Some(r) = crate::helpers::try_simd_f64_unary(input, crate::dispatch::simd_neg_f64) {
221        return r;
222    }
223    unary_float_op(input, |x| -x)
224}
225
226/// Elementwise positive (identity for numeric types).
227pub fn positive<T, D>(input: &Array<T, D>) -> FerrayResult<Array<T, D>>
228where
229    T: Element + Float,
230    D: Dimension,
231{
232    unary_float_op(input, |x| x)
233}
234
235/// Elementwise reciprocal: 1/x.
236///
237/// Uses hardware SIMD for contiguous f64 arrays.
238pub fn reciprocal<T, D>(input: &Array<T, D>) -> FerrayResult<Array<T, D>>
239where
240    T: Element + Float,
241    D: Dimension,
242{
243    if let Some(r) = crate::helpers::try_simd_f64_unary(input, crate::dispatch::simd_reciprocal_f64)
244    {
245        return r;
246    }
247    unary_float_op(input, T::recip)
248}
249
250/// Elementwise square root.
251///
252/// Uses hardware SIMD (`vsqrtpd`) for contiguous f64 arrays.
253pub fn sqrt<T, D>(input: &Array<T, D>) -> FerrayResult<Array<T, D>>
254where
255    T: Element + Float,
256    D: Dimension,
257{
258    if let Some(r) = crate::helpers::try_simd_f64_unary(input, crate::dispatch::simd_sqrt_f64) {
259        return r;
260    }
261    unary_float_op(input, T::sqrt)
262}
263
264/// Elementwise cube root.
265pub fn cbrt<T, D>(input: &Array<T, D>) -> FerrayResult<Array<T, D>>
266where
267    T: Element + Float + crate::cr_math::CrMath,
268    D: Dimension,
269{
270    unary_float_op(input, T::cr_cbrt)
271}
272
273/// Elementwise square: x^2.
274///
275/// Uses hardware SIMD for contiguous f64 arrays.
276pub fn square<T, D>(input: &Array<T, D>) -> FerrayResult<Array<T, D>>
277where
278    T: Element + Float,
279    D: Dimension,
280{
281    if let Some(r) = crate::helpers::try_simd_f64_unary(input, crate::dispatch::simd_square_f64) {
282        return r;
283    }
284    unary_float_op(input, |x| x * x)
285}
286
287/// Heaviside step function.
288///
289/// `heaviside(x, h0)` returns 0 for x < 0, h0 for x == 0, and 1 for x > 0.
290pub fn heaviside<T, D>(x: &Array<T, D>, h0: &Array<T, D>) -> FerrayResult<Array<T, D>>
291where
292    T: Element + Float,
293    D: Dimension,
294{
295    binary_float_op(x, h0, |xi, h0i| {
296        if xi < <T as Element>::zero() {
297            <T as Element>::zero()
298        } else if xi == <T as Element>::zero() {
299            h0i
300        } else {
301            <T as Element>::one()
302        }
303    })
304}
305
306/// Integer GCD (works on float representations of integers).
307pub fn gcd<T, D>(a: &Array<T, D>, b: &Array<T, D>) -> FerrayResult<Array<T, D>>
308where
309    T: Element + Float,
310    D: Dimension,
311{
312    binary_float_op(a, b, |mut x, mut y| {
313        x = x.abs();
314        y = y.abs();
315        while y != <T as Element>::zero() {
316            let t = y;
317            y = x % y;
318            x = t;
319        }
320        x
321    })
322}
323
324/// Integer LCM (works on float representations of integers).
325pub fn lcm<T, D>(a: &Array<T, D>, b: &Array<T, D>) -> FerrayResult<Array<T, D>>
326where
327    T: Element + Float,
328    D: Dimension,
329{
330    binary_float_op(a, b, |x, y| {
331        let ax = x.abs();
332        let ay = y.abs();
333        if ax == <T as Element>::zero() || ay == <T as Element>::zero() {
334            return <T as Element>::zero();
335        }
336        // lcm = |a*b| / gcd(a,b)
337        let mut gx = ax;
338        let mut gy = ay;
339        while gy != <T as Element>::zero() {
340            let t = gy;
341            gy = gx % gy;
342            gx = t;
343        }
344        ax / gx * ay
345    })
346}
347
348// ---------------------------------------------------------------------------
349// Broadcasting binary arithmetic
350// ---------------------------------------------------------------------------
351
352/// Elementwise addition with broadcasting.
353pub fn add_broadcast<T, D1, D2>(a: &Array<T, D1>, b: &Array<T, D2>) -> FerrayResult<Array<T, IxDyn>>
354where
355    T: Element + std::ops::Add<Output = T> + Copy,
356    D1: Dimension,
357    D2: Dimension,
358{
359    binary_broadcast_op(a, b, |x, y| x + y)
360}
361
362/// Elementwise subtraction with broadcasting.
363pub fn subtract_broadcast<T, D1, D2>(
364    a: &Array<T, D1>,
365    b: &Array<T, D2>,
366) -> FerrayResult<Array<T, IxDyn>>
367where
368    T: Element + std::ops::Sub<Output = T> + Copy,
369    D1: Dimension,
370    D2: Dimension,
371{
372    binary_broadcast_op(a, b, |x, y| x - y)
373}
374
375/// Elementwise multiplication with broadcasting.
376pub fn multiply_broadcast<T, D1, D2>(
377    a: &Array<T, D1>,
378    b: &Array<T, D2>,
379) -> FerrayResult<Array<T, IxDyn>>
380where
381    T: Element + std::ops::Mul<Output = T> + Copy,
382    D1: Dimension,
383    D2: Dimension,
384{
385    binary_broadcast_op(a, b, |x, y| x * y)
386}
387
388/// Elementwise division with broadcasting.
389pub fn divide_broadcast<T, D1, D2>(
390    a: &Array<T, D1>,
391    b: &Array<T, D2>,
392) -> FerrayResult<Array<T, IxDyn>>
393where
394    T: Element + std::ops::Div<Output = T> + Copy,
395    D1: Dimension,
396    D2: Dimension,
397{
398    binary_broadcast_op(a, b, |x, y| x / y)
399}
400
401// ---------------------------------------------------------------------------
402// Reductions
403// ---------------------------------------------------------------------------
404
405/// Reduce by addition along an axis (column sums, row sums, etc.).
406///
407/// AC-2: `add_reduce` computes correct column sums.
408pub fn add_reduce<T, D>(input: &Array<T, D>, axis: usize) -> FerrayResult<Array<T, IxDyn>>
409where
410    T: Element + std::ops::Add<Output = T> + Copy,
411    D: Dimension,
412{
413    let ndim = input.ndim();
414    if axis >= ndim {
415        return Err(FerrayError::axis_out_of_bounds(axis, ndim));
416    }
417    let shape = input.shape().to_vec();
418    let axis_len = shape[axis];
419
420    // Compute output shape: remove the axis
421    let mut out_shape: Vec<usize> = Vec::with_capacity(ndim - 1);
422    for (i, &s) in shape.iter().enumerate() {
423        if i != axis {
424            out_shape.push(s);
425        }
426    }
427    let out_size: usize = out_shape.iter().product();
428
429    // Compute strides for input traversal
430    let mut stride = 1usize;
431    for d in shape.iter().skip(axis + 1) {
432        stride *= d;
433    }
434    let outer_size: usize = shape[..axis].iter().product();
435    let inner_size = stride;
436
437    let data: Vec<T> = input.iter().copied().collect();
438    let mut result = vec![<T as Element>::zero(); out_size];
439
440    for outer in 0..outer_size {
441        for inner in 0..inner_size {
442            let mut acc = <T as Element>::zero();
443            for k in 0..axis_len {
444                let idx = outer * axis_len * inner_size + k * inner_size + inner;
445                acc = acc + data[idx];
446            }
447            result[outer * inner_size + inner] = acc;
448        }
449    }
450
451    if out_shape.is_empty() {
452        out_shape.push(1);
453    }
454    Array::from_vec(IxDyn::from(&out_shape[..]), result)
455}
456
457/// Running (cumulative) addition along an axis.
458///
459/// AC-2: `add_accumulate` produces running sums.
460pub fn add_accumulate<T, D>(input: &Array<T, D>, axis: usize) -> FerrayResult<Array<T, D>>
461where
462    T: Element + std::ops::Add<Output = T> + Copy,
463    D: Dimension,
464{
465    cumsum(input, Some(axis))
466}
467
468/// Outer product: multiply_outer(a, b)[i, j] = a[i] * b[j].
469///
470/// AC-3: multiply_outer produces correct outer product.
471pub fn multiply_outer<T>(a: &Array<T, Ix1>, b: &Array<T, Ix1>) -> FerrayResult<Array<T, IxDyn>>
472where
473    T: Element + std::ops::Mul<Output = T> + Copy,
474{
475    use ferray_core::dimension::Ix2;
476    let m = a.size();
477    let n = b.size();
478    let a_data: Vec<T> = a.iter().copied().collect();
479    let b_data: Vec<T> = b.iter().copied().collect();
480    let mut data = Vec::with_capacity(m * n);
481    for &ai in &a_data {
482        for &bj in &b_data {
483            data.push(ai * bj);
484        }
485    }
486    let arr = Array::from_vec(Ix2::new([m, n]), data)?;
487    // Convert to IxDyn
488    let dyn_data: Vec<T> = arr.iter().copied().collect();
489    Array::from_vec(IxDyn::from(&[m, n][..]), dyn_data)
490}
491
492// ---------------------------------------------------------------------------
493// Cumulative operations
494// ---------------------------------------------------------------------------
495
496/// Cumulative sum along an axis (or flattened if axis is None).
497///
498/// AC-11: `cumsum([1,2,3,4]) == [1,3,6,10]`.
499pub fn cumsum<T, D>(input: &Array<T, D>, axis: Option<usize>) -> FerrayResult<Array<T, D>>
500where
501    T: Element + std::ops::Add<Output = T> + Copy,
502    D: Dimension,
503{
504    if let Some(ax) = axis {
505        if ax >= input.ndim() {
506            return Err(FerrayError::axis_out_of_bounds(ax, input.ndim()));
507        }
508        // Work along the given axis
509        let shape = input.shape().to_vec();
510        let data: Vec<T> = input.iter().copied().collect();
511        let mut result = data.clone();
512        // Compute strides manually
513        let mut stride = 1usize;
514        for d in shape.iter().skip(ax + 1) {
515            stride *= d;
516        }
517        let axis_len = shape[ax];
518        let outer_size: usize = shape[..ax].iter().product();
519        let inner_size = stride;
520
521        for outer in 0..outer_size {
522            for inner in 0..inner_size {
523                let base = outer * axis_len * inner_size + inner;
524                for k in 1..axis_len {
525                    let prev = base + (k - 1) * inner_size;
526                    let curr = base + k * inner_size;
527                    result[curr] = result[prev] + result[curr];
528                }
529            }
530        }
531        Array::from_vec(input.dim().clone(), result)
532    } else {
533        // Flatten and cumsum
534        let mut data: Vec<T> = input.iter().copied().collect();
535        for i in 1..data.len() {
536            data[i] = data[i - 1] + data[i];
537        }
538        Array::from_vec(input.dim().clone(), data)
539    }
540}
541
542/// Cumulative product along an axis (or flattened if axis is None).
543pub fn cumprod<T, D>(input: &Array<T, D>, axis: Option<usize>) -> FerrayResult<Array<T, D>>
544where
545    T: Element + std::ops::Mul<Output = T> + Copy,
546    D: Dimension,
547{
548    if let Some(ax) = axis {
549        if ax >= input.ndim() {
550            return Err(FerrayError::axis_out_of_bounds(ax, input.ndim()));
551        }
552        let shape = input.shape().to_vec();
553        let data: Vec<T> = input.iter().copied().collect();
554        let mut result = data.clone();
555        let mut stride = 1usize;
556        for d in shape.iter().skip(ax + 1) {
557            stride *= d;
558        }
559        let axis_len = shape[ax];
560        let outer_size: usize = shape[..ax].iter().product();
561        let inner_size = stride;
562
563        for outer in 0..outer_size {
564            for inner in 0..inner_size {
565                let base = outer * axis_len * inner_size + inner;
566                for k in 1..axis_len {
567                    let prev = base + (k - 1) * inner_size;
568                    let curr = base + k * inner_size;
569                    result[curr] = result[prev] * result[curr];
570                }
571            }
572        }
573        Array::from_vec(input.dim().clone(), result)
574    } else {
575        let mut data: Vec<T> = input.iter().copied().collect();
576        for i in 1..data.len() {
577            data[i] = data[i - 1] * data[i];
578        }
579        Array::from_vec(input.dim().clone(), data)
580    }
581}
582
583/// Cumulative sum ignoring NaNs.
584pub fn nancumsum<T, D>(input: &Array<T, D>, axis: Option<usize>) -> FerrayResult<Array<T, D>>
585where
586    T: Element + Float,
587    D: Dimension,
588{
589    // Replace NaN with zero, then cumsum
590    let cleaned: Vec<T> = input
591        .iter()
592        .map(|&x| {
593            if x.is_nan() {
594                <T as Element>::zero()
595            } else {
596                x
597            }
598        })
599        .collect();
600    let arr = Array::from_vec(input.dim().clone(), cleaned)?;
601    cumsum(&arr, axis)
602}
603
604/// Cumulative product ignoring NaNs.
605pub fn nancumprod<T, D>(input: &Array<T, D>, axis: Option<usize>) -> FerrayResult<Array<T, D>>
606where
607    T: Element + Float,
608    D: Dimension,
609{
610    let cleaned: Vec<T> = input
611        .iter()
612        .map(|&x| if x.is_nan() { <T as Element>::one() } else { x })
613        .collect();
614    let arr = Array::from_vec(input.dim().clone(), cleaned)?;
615    cumprod(&arr, axis)
616}
617
618// ---------------------------------------------------------------------------
619// Differences
620// ---------------------------------------------------------------------------
621
622/// Compute the n-th discrete difference along the given axis.
623///
624/// AC-11: `diff([1,3,6,10], 1) == [2,3,4]`.
625pub fn diff<T>(input: &Array<T, Ix1>, n: usize) -> FerrayResult<Array<T, Ix1>>
626where
627    T: Element + std::ops::Sub<Output = T> + Copy,
628{
629    let mut data: Vec<T> = input.iter().copied().collect();
630    for _ in 0..n {
631        if data.len() <= 1 {
632            data.clear();
633            break;
634        }
635        let mut new_data = Vec::with_capacity(data.len() - 1);
636        for i in 1..data.len() {
637            new_data.push(data[i] - data[i - 1]);
638        }
639        data = new_data;
640    }
641    Array::from_vec(Ix1::new([data.len()]), data)
642}
643
644/// Differences between consecutive elements of an array, with optional
645/// prepend/append values.
646pub fn ediff1d<T>(
647    input: &Array<T, Ix1>,
648    to_end: Option<&[T]>,
649    to_begin: Option<&[T]>,
650) -> FerrayResult<Array<T, Ix1>>
651where
652    T: Element + std::ops::Sub<Output = T> + Copy,
653{
654    let data: Vec<T> = input.iter().copied().collect();
655    let mut result = Vec::new();
656
657    if let Some(begin) = to_begin {
658        result.extend_from_slice(begin);
659    }
660
661    for i in 1..data.len() {
662        result.push(data[i] - data[i - 1]);
663    }
664
665    if let Some(end) = to_end {
666        result.extend_from_slice(end);
667    }
668
669    Array::from_vec(Ix1::new([result.len()]), result)
670}
671
672/// Compute the gradient of a 1-D array using central differences.
673///
674/// Edge values use forward/backward differences.
675pub fn gradient<T>(input: &Array<T, Ix1>, spacing: Option<T>) -> FerrayResult<Array<T, Ix1>>
676where
677    T: Element + Float,
678{
679    let data: Vec<T> = input.iter().copied().collect();
680    let n = data.len();
681    if n == 0 {
682        return Array::from_vec(Ix1::new([0]), vec![]);
683    }
684    let h = spacing.unwrap_or_else(|| <T as Element>::one());
685    let two = <T as Element>::one() + <T as Element>::one();
686    let mut result = Vec::with_capacity(n);
687
688    if n == 1 {
689        result.push(<T as Element>::zero());
690    } else {
691        // Forward difference for first element
692        result.push((data[1] - data[0]) / h);
693        // Central differences for interior
694        for i in 1..n - 1 {
695            result.push((data[i + 1] - data[i - 1]) / (two * h));
696        }
697        // Backward difference for last element
698        result.push((data[n - 1] - data[n - 2]) / h);
699    }
700
701    Array::from_vec(Ix1::new([n]), result)
702}
703
704// ---------------------------------------------------------------------------
705// Cross product
706// ---------------------------------------------------------------------------
707
708/// Cross product of two 3-element 1-D arrays.
709pub fn cross<T>(a: &Array<T, Ix1>, b: &Array<T, Ix1>) -> FerrayResult<Array<T, Ix1>>
710where
711    T: Element + std::ops::Mul<Output = T> + std::ops::Sub<Output = T> + Copy,
712{
713    if a.size() != 3 || b.size() != 3 {
714        return Err(FerrayError::invalid_value(
715            "cross product requires 3-element vectors",
716        ));
717    }
718    let ad: Vec<T> = a.iter().copied().collect();
719    let bd: Vec<T> = b.iter().copied().collect();
720    let result = vec![
721        ad[1] * bd[2] - ad[2] * bd[1],
722        ad[2] * bd[0] - ad[0] * bd[2],
723        ad[0] * bd[1] - ad[1] * bd[0],
724    ];
725    Array::from_vec(Ix1::new([3]), result)
726}
727
728// ---------------------------------------------------------------------------
729// Integration
730// ---------------------------------------------------------------------------
731
732/// Integrate using the trapezoidal rule.
733///
734/// If `dx` is provided, it is the spacing between sample points.
735/// If `x` is provided, it gives the sample point coordinates.
736pub fn trapezoid<T>(y: &Array<T, Ix1>, x: Option<&Array<T, Ix1>>, dx: Option<T>) -> FerrayResult<T>
737where
738    T: Element + Float,
739{
740    let ydata: Vec<T> = y.iter().copied().collect();
741    let n = ydata.len();
742    if n < 2 {
743        return Ok(<T as Element>::zero());
744    }
745
746    let two = <T as Element>::one() + <T as Element>::one();
747    let mut total = <T as Element>::zero();
748
749    if let Some(xarr) = x {
750        let xdata: Vec<T> = xarr.iter().copied().collect();
751        if xdata.len() != n {
752            return Err(FerrayError::shape_mismatch(
753                "x and y must have the same length for trapezoid",
754            ));
755        }
756        for i in 1..n {
757            total = total + (ydata[i] + ydata[i - 1]) / two * (xdata[i] - xdata[i - 1]);
758        }
759    } else {
760        let h = dx.unwrap_or_else(|| <T as Element>::one());
761        for i in 1..n {
762            total = total + (ydata[i] + ydata[i - 1]) / two * h;
763        }
764    }
765
766    Ok(total)
767}
768
769// ---------------------------------------------------------------------------
770// f16 variants (f32-promoted)
771// ---------------------------------------------------------------------------
772
773/// Elementwise absolute value for f16 arrays via f32 promotion.
774#[cfg(feature = "f16")]
775pub fn absolute_f16<D>(input: &Array<half::f16, D>) -> FerrayResult<Array<half::f16, D>>
776where
777    D: Dimension,
778{
779    crate::helpers::unary_f16_op(input, f32::abs)
780}
781
782/// Elementwise negation for f16 arrays via f32 promotion.
783#[cfg(feature = "f16")]
784pub fn negative_f16<D>(input: &Array<half::f16, D>) -> FerrayResult<Array<half::f16, D>>
785where
786    D: Dimension,
787{
788    crate::helpers::unary_f16_op(input, |x| -x)
789}
790
791/// Elementwise square root for f16 arrays via f32 promotion.
792#[cfg(feature = "f16")]
793pub fn sqrt_f16<D>(input: &Array<half::f16, D>) -> FerrayResult<Array<half::f16, D>>
794where
795    D: Dimension,
796{
797    crate::helpers::unary_f16_op(input, f32::sqrt)
798}
799
800/// Elementwise cube root for f16 arrays via f32 promotion.
801#[cfg(feature = "f16")]
802pub fn cbrt_f16<D>(input: &Array<half::f16, D>) -> FerrayResult<Array<half::f16, D>>
803where
804    D: Dimension,
805{
806    crate::helpers::unary_f16_op(input, f32::cbrt)
807}
808
809/// Elementwise square for f16 arrays via f32 promotion.
810#[cfg(feature = "f16")]
811pub fn square_f16<D>(input: &Array<half::f16, D>) -> FerrayResult<Array<half::f16, D>>
812where
813    D: Dimension,
814{
815    crate::helpers::unary_f16_op(input, |x| x * x)
816}
817
818/// Elementwise reciprocal for f16 arrays via f32 promotion.
819#[cfg(feature = "f16")]
820pub fn reciprocal_f16<D>(input: &Array<half::f16, D>) -> FerrayResult<Array<half::f16, D>>
821where
822    D: Dimension,
823{
824    crate::helpers::unary_f16_op(input, f32::recip)
825}
826
827/// Elementwise sign for f16 arrays via f32 promotion.
828#[cfg(feature = "f16")]
829pub fn sign_f16<D>(input: &Array<half::f16, D>) -> FerrayResult<Array<half::f16, D>>
830where
831    D: Dimension,
832{
833    crate::helpers::unary_f16_op(input, |x| {
834        if x.is_nan() {
835            f32::NAN
836        } else if x > 0.0 {
837            1.0
838        } else if x < 0.0 {
839            -1.0
840        } else {
841            0.0
842        }
843    })
844}
845
846/// Elementwise addition for f16 arrays via f32 promotion.
847#[cfg(feature = "f16")]
848pub fn add_f16<D>(
849    a: &Array<half::f16, D>,
850    b: &Array<half::f16, D>,
851) -> FerrayResult<Array<half::f16, D>>
852where
853    D: Dimension,
854{
855    crate::helpers::binary_f16_op(a, b, |x, y| x + y)
856}
857
858/// Elementwise subtraction for f16 arrays via f32 promotion.
859#[cfg(feature = "f16")]
860pub fn subtract_f16<D>(
861    a: &Array<half::f16, D>,
862    b: &Array<half::f16, D>,
863) -> FerrayResult<Array<half::f16, D>>
864where
865    D: Dimension,
866{
867    crate::helpers::binary_f16_op(a, b, |x, y| x - y)
868}
869
870/// Elementwise multiplication for f16 arrays via f32 promotion.
871#[cfg(feature = "f16")]
872pub fn multiply_f16<D>(
873    a: &Array<half::f16, D>,
874    b: &Array<half::f16, D>,
875) -> FerrayResult<Array<half::f16, D>>
876where
877    D: Dimension,
878{
879    crate::helpers::binary_f16_op(a, b, |x, y| x * y)
880}
881
882/// Elementwise division for f16 arrays via f32 promotion.
883#[cfg(feature = "f16")]
884pub fn divide_f16<D>(
885    a: &Array<half::f16, D>,
886    b: &Array<half::f16, D>,
887) -> FerrayResult<Array<half::f16, D>>
888where
889    D: Dimension,
890{
891    crate::helpers::binary_f16_op(a, b, |x, y| x / y)
892}
893
894/// Elementwise power for f16 arrays via f32 promotion.
895#[cfg(feature = "f16")]
896pub fn power_f16<D>(
897    a: &Array<half::f16, D>,
898    b: &Array<half::f16, D>,
899) -> FerrayResult<Array<half::f16, D>>
900where
901    D: Dimension,
902{
903    crate::helpers::binary_f16_op(a, b, f32::powf)
904}
905
906/// Floor division for f16 arrays via f32 promotion.
907#[cfg(feature = "f16")]
908pub fn floor_divide_f16<D>(
909    a: &Array<half::f16, D>,
910    b: &Array<half::f16, D>,
911) -> FerrayResult<Array<half::f16, D>>
912where
913    D: Dimension,
914{
915    crate::helpers::binary_f16_op(a, b, |x, y| (x / y).floor())
916}
917
918/// Elementwise remainder for f16 arrays via f32 promotion.
919#[cfg(feature = "f16")]
920pub fn remainder_f16<D>(
921    a: &Array<half::f16, D>,
922    b: &Array<half::f16, D>,
923) -> FerrayResult<Array<half::f16, D>>
924where
925    D: Dimension,
926{
927    crate::helpers::binary_f16_op(a, b, |x, y| {
928        let r = x % y;
929        if (r < 0.0 && y > 0.0) || (r > 0.0 && y < 0.0) {
930            r + y
931        } else {
932            r
933        }
934    })
935}
936
937#[cfg(test)]
938mod tests {
939    use super::*;
940    use ferray_core::dimension::Ix2;
941
942    fn arr1(data: Vec<f64>) -> Array<f64, Ix1> {
943        let n = data.len();
944        Array::from_vec(Ix1::new([n]), data).unwrap()
945    }
946
947    fn arr1_i32(data: Vec<i32>) -> Array<i32, Ix1> {
948        let n = data.len();
949        Array::from_vec(Ix1::new([n]), data).unwrap()
950    }
951
952    #[test]
953    fn test_add() {
954        let a = arr1(vec![1.0, 2.0, 3.0]);
955        let b = arr1(vec![4.0, 5.0, 6.0]);
956        let r = add(&a, &b).unwrap();
957        assert_eq!(r.as_slice().unwrap(), &[5.0, 7.0, 9.0]);
958    }
959
960    #[test]
961    fn test_subtract() {
962        let a = arr1(vec![5.0, 7.0, 9.0]);
963        let b = arr1(vec![1.0, 2.0, 3.0]);
964        let r = subtract(&a, &b).unwrap();
965        assert_eq!(r.as_slice().unwrap(), &[4.0, 5.0, 6.0]);
966    }
967
968    #[test]
969    fn test_multiply() {
970        let a = arr1(vec![2.0, 3.0, 4.0]);
971        let b = arr1(vec![5.0, 6.0, 7.0]);
972        let r = multiply(&a, &b).unwrap();
973        assert_eq!(r.as_slice().unwrap(), &[10.0, 18.0, 28.0]);
974    }
975
976    #[test]
977    fn test_divide() {
978        let a = arr1(vec![10.0, 20.0, 30.0]);
979        let b = arr1(vec![2.0, 4.0, 5.0]);
980        let r = divide(&a, &b).unwrap();
981        assert_eq!(r.as_slice().unwrap(), &[5.0, 5.0, 6.0]);
982    }
983
984    #[test]
985    fn test_floor_divide() {
986        let a = arr1(vec![7.0, -7.0]);
987        let b = arr1(vec![2.0, 2.0]);
988        let r = floor_divide(&a, &b).unwrap();
989        assert_eq!(r.as_slice().unwrap(), &[3.0, -4.0]);
990    }
991
992    #[test]
993    fn test_power() {
994        let a = arr1(vec![2.0, 3.0]);
995        let b = arr1(vec![3.0, 2.0]);
996        let r = power(&a, &b).unwrap();
997        assert_eq!(r.as_slice().unwrap(), &[8.0, 9.0]);
998    }
999
1000    #[test]
1001    fn test_remainder() {
1002        let a = arr1(vec![7.0, -7.0]);
1003        let b = arr1(vec![3.0, 3.0]);
1004        let r = remainder(&a, &b).unwrap();
1005        let s = r.as_slice().unwrap();
1006        assert!((s[0] - 1.0).abs() < 1e-12);
1007        assert!((s[1] - 2.0).abs() < 1e-12);
1008    }
1009
1010    #[test]
1011    fn test_fmod() {
1012        let a = arr1(vec![7.0, -7.0]);
1013        let b = arr1(vec![3.0, 3.0]);
1014        let r = fmod(&a, &b).unwrap();
1015        let s = r.as_slice().unwrap();
1016        assert!((s[0] - 1.0).abs() < 1e-12);
1017        assert!((s[1] - (-1.0)).abs() < 1e-12);
1018    }
1019
1020    #[test]
1021    fn test_absolute() {
1022        let a = arr1(vec![-1.0, 2.0, -3.0]);
1023        let r = absolute(&a).unwrap();
1024        assert_eq!(r.as_slice().unwrap(), &[1.0, 2.0, 3.0]);
1025    }
1026
1027    #[test]
1028    fn test_sign() {
1029        let a = arr1(vec![-5.0, 0.0, 3.0]);
1030        let r = sign(&a).unwrap();
1031        assert_eq!(r.as_slice().unwrap(), &[-1.0, 0.0, 1.0]);
1032    }
1033
1034    #[test]
1035    fn test_negative() {
1036        let a = arr1(vec![1.0, -2.0, 3.0]);
1037        let r = negative(&a).unwrap();
1038        assert_eq!(r.as_slice().unwrap(), &[-1.0, 2.0, -3.0]);
1039    }
1040
1041    #[test]
1042    fn test_sqrt() {
1043        let a = arr1(vec![1.0, 4.0, 9.0, 16.0]);
1044        let r = sqrt(&a).unwrap();
1045        assert_eq!(r.as_slice().unwrap(), &[1.0, 2.0, 3.0, 4.0]);
1046    }
1047
1048    #[test]
1049    fn test_cbrt() {
1050        let a = arr1(vec![8.0, 27.0]);
1051        let r = cbrt(&a).unwrap();
1052        let s = r.as_slice().unwrap();
1053        assert!((s[0] - 2.0).abs() < 1e-12);
1054        assert!((s[1] - 3.0).abs() < 1e-12);
1055    }
1056
1057    #[test]
1058    fn test_square() {
1059        let a = arr1(vec![2.0, 3.0, 4.0]);
1060        let r = square(&a).unwrap();
1061        assert_eq!(r.as_slice().unwrap(), &[4.0, 9.0, 16.0]);
1062    }
1063
1064    #[test]
1065    fn test_reciprocal() {
1066        let a = arr1(vec![2.0, 4.0, 5.0]);
1067        let r = reciprocal(&a).unwrap();
1068        assert_eq!(r.as_slice().unwrap(), &[0.5, 0.25, 0.2]);
1069    }
1070
1071    #[test]
1072    fn test_heaviside() {
1073        let x = arr1(vec![-1.0, 0.0, 1.0]);
1074        let h0 = arr1(vec![0.5, 0.5, 0.5]);
1075        let r = heaviside(&x, &h0).unwrap();
1076        assert_eq!(r.as_slice().unwrap(), &[0.0, 0.5, 1.0]);
1077    }
1078
1079    #[test]
1080    fn test_gcd() {
1081        let a = arr1(vec![12.0, 15.0]);
1082        let b = arr1(vec![8.0, 25.0]);
1083        let r = gcd(&a, &b).unwrap();
1084        assert_eq!(r.as_slice().unwrap(), &[4.0, 5.0]);
1085    }
1086
1087    #[test]
1088    fn test_lcm() {
1089        let a = arr1(vec![4.0, 6.0]);
1090        let b = arr1(vec![6.0, 8.0]);
1091        let r = lcm(&a, &b).unwrap();
1092        assert_eq!(r.as_slice().unwrap(), &[12.0, 24.0]);
1093    }
1094
1095    #[test]
1096    fn test_cumsum_ac11() {
1097        // AC-11: cumsum([1,2,3,4]) == [1,3,6,10]
1098        let a = arr1(vec![1.0, 2.0, 3.0, 4.0]);
1099        let r = cumsum(&a, None).unwrap();
1100        assert_eq!(r.as_slice().unwrap(), &[1.0, 3.0, 6.0, 10.0]);
1101    }
1102
1103    #[test]
1104    fn test_cumsum_i32() {
1105        let a = arr1_i32(vec![1, 2, 3, 4]);
1106        let r = cumsum(&a, None).unwrap();
1107        assert_eq!(r.as_slice().unwrap(), &[1, 3, 6, 10]);
1108    }
1109
1110    #[test]
1111    fn test_cumprod() {
1112        let a = arr1(vec![1.0, 2.0, 3.0, 4.0]);
1113        let r = cumprod(&a, None).unwrap();
1114        assert_eq!(r.as_slice().unwrap(), &[1.0, 2.0, 6.0, 24.0]);
1115    }
1116
1117    #[test]
1118    fn test_diff_ac11() {
1119        // AC-11: diff([1,3,6,10], 1) == [2,3,4]
1120        let a = arr1(vec![1.0, 3.0, 6.0, 10.0]);
1121        let r = diff(&a, 1).unwrap();
1122        assert_eq!(r.as_slice().unwrap(), &[2.0, 3.0, 4.0]);
1123    }
1124
1125    #[test]
1126    fn test_diff_n2() {
1127        let a = arr1(vec![1.0, 3.0, 6.0, 10.0]);
1128        let r = diff(&a, 2).unwrap();
1129        assert_eq!(r.as_slice().unwrap(), &[1.0, 1.0]);
1130    }
1131
1132    #[test]
1133    fn test_ediff1d() {
1134        let a = arr1(vec![1.0, 2.0, 4.0, 7.0]);
1135        let r = ediff1d(&a, None, None).unwrap();
1136        assert_eq!(r.as_slice().unwrap(), &[1.0, 2.0, 3.0]);
1137    }
1138
1139    #[test]
1140    fn test_gradient() {
1141        let a = arr1(vec![1.0, 2.0, 4.0, 7.0, 11.0]);
1142        let r = gradient(&a, None).unwrap();
1143        let s = r.as_slice().unwrap();
1144        // forward: 2-1=1, central: (4-1)/2=1.5, (7-2)/2=2.5, (11-4)/2=3.5, backward: 11-7=4
1145        assert!((s[0] - 1.0).abs() < 1e-12);
1146        assert!((s[1] - 1.5).abs() < 1e-12);
1147        assert!((s[2] - 2.5).abs() < 1e-12);
1148        assert!((s[3] - 3.5).abs() < 1e-12);
1149        assert!((s[4] - 4.0).abs() < 1e-12);
1150    }
1151
1152    #[test]
1153    fn test_cross() {
1154        let a = arr1(vec![1.0, 0.0, 0.0]);
1155        let b = arr1(vec![0.0, 1.0, 0.0]);
1156        let r = cross(&a, &b).unwrap();
1157        assert_eq!(r.as_slice().unwrap(), &[0.0, 0.0, 1.0]);
1158    }
1159
1160    #[test]
1161    fn test_trapezoid() {
1162        // Integrate y=x from 0 to 4: area = 8
1163        let y = arr1(vec![0.0, 1.0, 2.0, 3.0, 4.0]);
1164        let r = trapezoid(&y, None, Some(1.0)).unwrap();
1165        assert!((r - 8.0).abs() < 1e-12);
1166    }
1167
1168    #[test]
1169    fn test_trapezoid_with_x() {
1170        let y = arr1(vec![0.0, 1.0, 4.0]);
1171        let x = arr1(vec![0.0, 1.0, 2.0]);
1172        let r = trapezoid(&y, Some(&x), None).unwrap();
1173        // (0+1)/2*1 + (1+4)/2*1 = 0.5 + 2.5 = 3.0
1174        assert!((r - 3.0).abs() < 1e-12);
1175    }
1176
1177    #[test]
1178    fn test_add_reduce_ac2() {
1179        // AC-2: add_reduce computes correct column sums
1180        let a = Array::<f64, Ix2>::from_vec(Ix2::new([2, 3]), vec![1.0, 2.0, 3.0, 4.0, 5.0, 6.0])
1181            .unwrap();
1182        let r = add_reduce(&a, 0).unwrap();
1183        assert_eq!(r.shape(), &[3]);
1184        let s: Vec<f64> = r.iter().copied().collect();
1185        assert_eq!(s, vec![5.0, 7.0, 9.0]);
1186    }
1187
1188    #[test]
1189    fn test_add_accumulate_ac2() {
1190        let a = arr1(vec![1.0, 2.0, 3.0, 4.0]);
1191        let r = add_accumulate(&a, 0).unwrap();
1192        assert_eq!(r.as_slice().unwrap(), &[1.0, 3.0, 6.0, 10.0]);
1193    }
1194
1195    #[test]
1196    fn test_multiply_outer_ac3() {
1197        // AC-3: multiply_outer produces correct outer product
1198        let a = arr1(vec![1.0, 2.0, 3.0]);
1199        let b = arr1(vec![4.0, 5.0]);
1200        let r = multiply_outer(&a, &b).unwrap();
1201        assert_eq!(r.shape(), &[3, 2]);
1202        let s: Vec<f64> = r.iter().copied().collect();
1203        assert_eq!(s, vec![4.0, 5.0, 8.0, 10.0, 12.0, 15.0]);
1204    }
1205
1206    #[test]
1207    fn test_nancumsum() {
1208        let a = arr1(vec![1.0, f64::NAN, 3.0, 4.0]);
1209        let r = nancumsum(&a, None).unwrap();
1210        let s = r.as_slice().unwrap();
1211        assert_eq!(s[0], 1.0);
1212        assert_eq!(s[1], 1.0); // NaN treated as 0
1213        assert_eq!(s[2], 4.0);
1214        assert_eq!(s[3], 8.0);
1215    }
1216
1217    #[test]
1218    fn test_nancumprod() {
1219        let a = arr1(vec![1.0, f64::NAN, 3.0, 4.0]);
1220        let r = nancumprod(&a, None).unwrap();
1221        let s = r.as_slice().unwrap();
1222        assert_eq!(s[0], 1.0);
1223        assert_eq!(s[1], 1.0); // NaN treated as 1
1224        assert_eq!(s[2], 3.0);
1225        assert_eq!(s[3], 12.0);
1226    }
1227
1228    #[test]
1229    fn test_add_broadcast() {
1230        let a = Array::<f64, Ix2>::from_vec(Ix2::new([2, 1]), vec![1.0, 2.0]).unwrap();
1231        let b = Array::<f64, Ix1>::from_vec(Ix1::new([3]), vec![10.0, 20.0, 30.0]).unwrap();
1232        let r = add_broadcast(&a, &b).unwrap();
1233        assert_eq!(r.shape(), &[2, 3]);
1234    }
1235
1236    #[test]
1237    fn test_divmod() {
1238        let a = arr1(vec![7.0, -7.0]);
1239        let b = arr1(vec![3.0, 3.0]);
1240        let (q, r) = divmod(&a, &b).unwrap();
1241        assert_eq!(q.as_slice().unwrap(), &[2.0, -3.0]);
1242        let rs = r.as_slice().unwrap();
1243        assert!((rs[0] - 1.0).abs() < 1e-12);
1244        assert!((rs[1] - 2.0).abs() < 1e-12);
1245    }
1246
1247    #[test]
1248    fn test_positive() {
1249        let a = arr1(vec![-1.0, 2.0]);
1250        let r = positive(&a).unwrap();
1251        assert_eq!(r.as_slice().unwrap(), &[-1.0, 2.0]);
1252    }
1253
1254    #[test]
1255    fn test_true_divide() {
1256        let a = arr1(vec![10.0, 20.0]);
1257        let b = arr1(vec![3.0, 7.0]);
1258        let r = true_divide(&a, &b).unwrap();
1259        let s = r.as_slice().unwrap();
1260        assert!((s[0] - 10.0 / 3.0).abs() < 1e-12);
1261        assert!((s[1] - 20.0 / 7.0).abs() < 1e-12);
1262    }
1263
1264    #[cfg(feature = "f16")]
1265    mod f16_tests {
1266        use super::*;
1267
1268        fn arr1_f16(data: &[f32]) -> Array<half::f16, Ix1> {
1269            let n = data.len();
1270            let vals: Vec<half::f16> = data.iter().map(|&x| half::f16::from_f32(x)).collect();
1271            Array::from_vec(Ix1::new([n]), vals).unwrap()
1272        }
1273
1274        #[test]
1275        fn test_add_f16() {
1276            let a = arr1_f16(&[1.0, 2.0, 3.0]);
1277            let b = arr1_f16(&[4.0, 5.0, 6.0]);
1278            let r = add_f16(&a, &b).unwrap();
1279            let s = r.as_slice().unwrap();
1280            assert!((s[0].to_f32() - 5.0).abs() < 0.01);
1281            assert!((s[1].to_f32() - 7.0).abs() < 0.01);
1282            assert!((s[2].to_f32() - 9.0).abs() < 0.01);
1283        }
1284
1285        #[test]
1286        fn test_multiply_f16() {
1287            let a = arr1_f16(&[2.0, 3.0]);
1288            let b = arr1_f16(&[4.0, 5.0]);
1289            let r = multiply_f16(&a, &b).unwrap();
1290            let s = r.as_slice().unwrap();
1291            assert!((s[0].to_f32() - 8.0).abs() < 0.01);
1292            assert!((s[1].to_f32() - 15.0).abs() < 0.1);
1293        }
1294
1295        #[test]
1296        fn test_sqrt_f16() {
1297            let a = arr1_f16(&[1.0, 4.0, 9.0, 16.0]);
1298            let r = sqrt_f16(&a).unwrap();
1299            let s = r.as_slice().unwrap();
1300            assert!((s[0].to_f32() - 1.0).abs() < 0.01);
1301            assert!((s[1].to_f32() - 2.0).abs() < 0.01);
1302            assert!((s[2].to_f32() - 3.0).abs() < 0.01);
1303            assert!((s[3].to_f32() - 4.0).abs() < 0.01);
1304        }
1305
1306        #[test]
1307        fn test_absolute_f16() {
1308            let a = arr1_f16(&[-1.0, 2.0, -3.0]);
1309            let r = absolute_f16(&a).unwrap();
1310            let s = r.as_slice().unwrap();
1311            assert!((s[0].to_f32() - 1.0).abs() < 0.01);
1312            assert!((s[1].to_f32() - 2.0).abs() < 0.01);
1313            assert!((s[2].to_f32() - 3.0).abs() < 0.01);
1314        }
1315
1316        #[test]
1317        fn test_power_f16() {
1318            let a = arr1_f16(&[2.0, 3.0]);
1319            let b = arr1_f16(&[3.0, 2.0]);
1320            let r = power_f16(&a, &b).unwrap();
1321            let s = r.as_slice().unwrap();
1322            assert!((s[0].to_f32() - 8.0).abs() < 0.1);
1323            assert!((s[1].to_f32() - 9.0).abs() < 0.1);
1324        }
1325
1326        #[test]
1327        fn test_divide_f16() {
1328            let a = arr1_f16(&[10.0, 20.0]);
1329            let b = arr1_f16(&[2.0, 4.0]);
1330            let r = divide_f16(&a, &b).unwrap();
1331            let s = r.as_slice().unwrap();
1332            assert!((s[0].to_f32() - 5.0).abs() < 0.01);
1333            assert!((s[1].to_f32() - 5.0).abs() < 0.01);
1334        }
1335    }
1336}