arr_rs/math/operations/trigonometric.rs
1use crate::{
2 core::prelude::*,
3 errors::prelude::*,
4 extensions::prelude::*,
5 numeric::prelude::*,
6 validators::prelude::*,
7};
8
9/// `ArrayTrait` - Array Trigonometric functions
10pub trait ArrayTrigonometric<N: NumericOps> where Self: Sized + Clone {
11
12 /// Compute the sine of array elements
13 ///
14 /// # Examples
15 ///
16 /// ```
17 /// use arr_rs::prelude::*;
18 ///
19 /// let arr = Array::flat(vec![-1., 0., 1.]).unwrap();
20 /// assert_eq!(Array::flat(vec![-0.8414709848078965, 0., 0.8414709848078965]), arr.sin());
21 /// ```
22 ///
23 /// # Errors
24 ///
25 /// may returns `ArrayError`
26 fn sin(&self) -> Result<Array<N>, ArrayError>;
27
28 /// Compute the cosine of array elements
29 ///
30 /// # Examples
31 ///
32 /// ```
33 /// use arr_rs::prelude::*;
34 ///
35 /// let arr = Array::flat(vec![-1., 0., 1.]).unwrap();
36 /// assert_eq!(Array::flat(vec![0.5403023058681398, 1., 0.5403023058681398]), arr.cos());
37 /// ```
38 ///
39 /// # Errors
40 ///
41 /// may returns `ArrayError`
42 fn cos(&self) -> Result<Array<N>, ArrayError>;
43
44 /// Compute the tangent of array elements
45 ///
46 /// # Examples
47 ///
48 /// ```
49 /// use arr_rs::prelude::*;
50 ///
51 /// let arr = Array::flat(vec![-1., 0., 1.]).unwrap();
52 /// assert_eq!(Array::flat(vec![-1.5574077246549023, 0., 1.5574077246549023]), arr.tan());
53 /// ```
54 ///
55 /// # Errors
56 ///
57 /// may returns `ArrayError`
58 fn tan(&self) -> Result<Array<N>, ArrayError>;
59
60 /// Compute the inverse sine of array elements
61 ///
62 /// # Examples
63 ///
64 /// ```
65 /// use arr_rs::prelude::*;
66 ///
67 /// let arr = Array::flat(vec![-1., 0., 1.]).unwrap();
68 /// assert_eq!(Array::flat(vec![-std::f64::consts::FRAC_PI_2, 0., std::f64::consts::FRAC_PI_2]), arr.asin());
69 /// ```
70 ///
71 /// # Errors
72 ///
73 /// may returns `ArrayError`
74 fn asin(&self) -> Result<Array<N>, ArrayError>;
75
76 /// Compute the inverse cosine of array elements
77 ///
78 /// # Examples
79 ///
80 /// ```
81 /// use arr_rs::prelude::*;
82 ///
83 /// let arr = Array::flat(vec![-1., 0., 1.]).unwrap();
84 /// assert_eq!(Array::flat(vec![std::f64::consts::PI, std::f64::consts::FRAC_PI_2, 0.]), arr.acos());
85 /// ```
86 ///
87 /// # Errors
88 ///
89 /// may returns `ArrayError`
90 fn acos(&self) -> Result<Array<N>, ArrayError>;
91
92 /// Compute the inverse tangent of array elements
93 ///
94 /// # Examples
95 ///
96 /// ```
97 /// use arr_rs::prelude::*;
98 ///
99 /// let arr = Array::flat(vec![-1., 0., 1.]).unwrap();
100 /// assert_eq!(Array::flat(vec![-std::f64::consts::FRAC_PI_4, 0., std::f64::consts::FRAC_PI_4]), arr.atan());
101 /// ```
102 ///
103 /// # Errors
104 ///
105 /// may returns `ArrayError`
106 fn atan(&self) -> Result<Array<N>, ArrayError>;
107
108 /// Compute the inverse tangent of x1/x2 choosing the quadrant correctly
109 ///
110 /// # Arguments
111 ///
112 /// * `other` - array to perform the operation with
113 ///
114 /// # Examples
115 ///
116 /// ```
117 /// use arr_rs::prelude::*;
118 ///
119 /// let arr = Array::flat(vec![-1., 0., 1.]).unwrap();
120 /// let other = Array::flat(vec![2., 4., 6.]).unwrap();
121 /// assert_eq!(Array::flat(vec![-0.4636476090008061,0., 0.16514867741462683]), arr.atan2(&other));
122 /// ```
123 ///
124 /// # Errors
125 ///
126 /// may returns `ArrayError`
127 fn atan2(&self, other: &Array<N>) -> Result<Array<N>, ArrayError>;
128
129 /// Given the “legs” of a right triangle, return its hypotenuse
130 ///
131 /// # Arguments
132 ///
133 /// * `other` - array to perform the operation with
134 ///
135 /// # Examples
136 ///
137 /// ```
138 /// use arr_rs::prelude::*;
139 ///
140 /// let arr = Array::full(vec![3, 3], 3).unwrap();
141 /// let other = Array::full(vec![3, 3], 4).unwrap();
142 /// assert_eq!(Array::full(vec![3, 3], 5), arr.hypot(&other));
143 /// ```
144 ///
145 /// # Errors
146 ///
147 /// may returns `ArrayError`
148 fn hypot(&self, other: &Array<N>) -> Result<Array<N>, ArrayError>;
149
150 /// Convert angles from radians to degrees
151 ///
152 /// # Examples
153 ///
154 /// ```
155 /// use arr_rs::prelude::*;
156 ///
157 /// let arr = Array::flat(vec![std::f64::consts::FRAC_PI_6, std::f64::consts::FRAC_PI_4, std::f64::consts::FRAC_PI_3, std::f64::consts::FRAC_PI_2]).unwrap();
158 /// let expected = Array::flat(vec![30., 45., 60., 90.]);
159 /// assert_eq!(format!("{expected:.8?}"), format!("{:.8?}", arr.degrees()));
160 /// ```
161 ///
162 /// # Errors
163 ///
164 /// may returns `ArrayError`
165 fn degrees(&self) -> Result<Array<N>, ArrayError>;
166
167 /// Convert angles from radians to degrees. alias on `degrees`
168 ///
169 /// # Examples
170 ///
171 /// ```
172 /// use arr_rs::prelude::*;
173 ///
174 /// let arr = Array::flat(vec![std::f64::consts::FRAC_PI_6, std::f64::consts::FRAC_PI_4, std::f64::consts::FRAC_PI_3, std::f64::consts::FRAC_PI_2]).unwrap();
175 /// let expected = Array::flat(vec![30., 45., 60., 90.]);
176 /// assert_eq!(format!("{expected:.8?}"), format!("{:.8?}", arr.rad2deg()));
177 /// ```
178 ///
179 /// # Errors
180 ///
181 /// may returns `ArrayError`
182 fn rad2deg(&self) -> Result<Array<N>, ArrayError>;
183
184 /// Convert angles from degrees to radians
185 ///
186 /// # Examples
187 ///
188 /// ```
189 /// use arr_rs::prelude::*;
190 ///
191 /// let arr = Array::flat(vec![30., 45., 60., 90.]).unwrap();
192 /// let expected = Array::flat(vec![std::f64::consts::FRAC_PI_6, std::f64::consts::FRAC_PI_4, std::f64::consts::FRAC_PI_3, std::f64::consts::FRAC_PI_2]);
193 /// assert_eq!(format!("{expected:.8?}"), format!("{:.8?}", arr.radians()));
194 /// ```
195 ///
196 /// # Errors
197 ///
198 /// may returns `ArrayError`
199 fn radians(&self) -> Result<Array<N>, ArrayError>;
200
201 /// Convert angles from degrees to radians. alias on `radians`
202 ///
203 /// # Examples
204 ///
205 /// ```
206 /// use arr_rs::prelude::*;
207 ///
208 /// let arr = Array::flat(vec![30., 45., 60., 90.]).unwrap();
209 /// let expected = Array::flat(vec![std::f64::consts::FRAC_PI_6, std::f64::consts::FRAC_PI_4, std::f64::consts::FRAC_PI_3, std::f64::consts::FRAC_PI_2]);
210 /// assert_eq!(format!("{expected:.8?}"), format!("{:.8?}", arr.deg2rad()));
211 /// ```
212 ///
213 /// # Errors
214 ///
215 /// may returns `ArrayError`
216 fn deg2rad(&self) -> Result<Array<N>, ArrayError>;
217
218 /// Unwrap by taking the complement of large deltas with respect to the period
219 ///
220 /// # Arguments
221 ///
222 /// * `discont` - Maximum discontinuity between values, default is period/2. Values below period/2 are treated as if they were period/2.
223 /// * `axis` - Axis along which unwrap will operate, default is the last axis.
224 /// * `period` - Size of the range over which the input wraps, default is 2 pi.
225 ///
226 /// # Examples
227 ///
228 /// ```
229 /// use arr_rs::prelude::*;
230 ///
231 /// let arr = Array::flat(vec![0.5, 1.5, 3.2, 6.8, 5.9]).unwrap();
232 /// let expected = Array::flat(vec![0.5, 1.5, 3.2, 0.51681469, -0.38318531]);
233 /// assert_eq!(format!("{expected:.8?}"), format!("{:.8?}", arr.unwrap_phase(None, None, None)));
234 /// ```
235 ///
236 /// # Errors
237 ///
238 /// may returns `ArrayError`
239 fn unwrap_phase(&self, discont: Option<Array<f64>>, axis: Option<isize>, period: Option<Array<f64>>) -> Result<Array<N>, ArrayError>;
240}
241
242impl <N: NumericOps> ArrayTrigonometric<N> for Array<N> {
243
244 fn sin(&self) -> Result<Self, ArrayError> {
245 self.map(|i| N::from(i.to_f64().sin()))
246 }
247
248 fn cos(&self) -> Result<Self, ArrayError> {
249 self.map(|i| N::from(i.to_f64().cos()))
250 }
251
252 fn tan(&self) -> Result<Self, ArrayError> {
253 self.map(|i| N::from(i.to_f64().tan()))
254 }
255
256 fn asin(&self) -> Result<Self, ArrayError> {
257 self.map(|i| N::from(i.to_f64().asin()))
258 }
259
260 fn acos(&self) -> Result<Self, ArrayError> {
261 self.map(|i| N::from(i.to_f64().acos()))
262 }
263
264 fn atan(&self) -> Result<Self, ArrayError> {
265 self.map(|i| N::from(i.to_f64().atan()))
266 }
267
268 fn atan2(&self, other: &Self) -> Result<Self, ArrayError> {
269 let broadcasted = self.broadcast(other)?;
270 let elements = broadcasted.clone().into_iter()
271 .map(|tuple| N::from(tuple.0.to_f64().atan2(tuple.1.to_f64())))
272 .collect();
273 Self::new(elements, broadcasted.get_shape()?)
274 }
275
276 fn hypot(&self, other: &Self) -> Result<Self, ArrayError> {
277 let broadcasted = self.broadcast(other)?;
278 let elements = broadcasted.clone().into_iter()
279 .map(|tuple| N::from(tuple.0.to_f64().hypot(tuple.1.to_f64())))
280 .collect();
281 Self::new(elements, broadcasted.get_shape()?)
282 }
283
284 fn degrees(&self) -> Result<Self, ArrayError> {
285 self.map(|i| N::from(i.to_f64().to_degrees()))
286 }
287
288 fn rad2deg(&self) -> Result<Self, ArrayError> {
289 self.degrees()
290 }
291
292 fn radians(&self) -> Result<Self, ArrayError> {
293 self.map(|i| N::from(i.to_f64().to_radians()))
294 }
295
296 fn deg2rad(&self) -> Result<Self, ArrayError> {
297 self.radians()
298 }
299
300 fn unwrap_phase(&self, discont: Option<Array<f64>>, axis: Option<isize>, period: Option<Array<f64>>) -> Result<Self, ArrayError> {
301
302 fn parse_parameter<N: Numeric>(array: &Array<N>, parameter: &Array<N>) -> Result<Array<N>, ArrayError> {
303 let self_len = array.len()? - 1;
304 parameter.len()?.is_one_of(vec![&1, &self_len])?;
305 let result =
306 if parameter.len()? == self_len { parameter.clone() }
307 else { Array::flat(vec![parameter[0]; self_len])? };
308 Ok(result)
309 }
310
311 let period = period.unwrap_or(Array::single(std::f64::consts::PI * 2.)?);
312 let discont = discont.unwrap_or((period.clone() / 2.)?);
313 let (mut discont, mut period) = (discont.to_array_num()?, period.to_array_num()?);
314
315 if self.ndim()? == 1 {
316 period = parse_parameter(self, &period)?;
317 discont = parse_parameter(self, &discont)?;
318
319 let mut unwrapped_phase = self.get_elements()?;
320 unwrapped_phase.clone()
321 .windows(2)
322 .enumerate()
323 .filter(|(idx, window)| N::from((window[1] - window[0]).to_f64().abs()) > discont[*idx])
324 .for_each(|(i, _)| {
325 let diff = (unwrapped_phase[i + 1] - unwrapped_phase[i]).to_f64();
326 unwrapped_phase[i + 1..].iter_mut().for_each(|val| *val -= N::from(diff.signum() * ((diff.abs() / period[i].to_f64()).ceil() * period[i].to_f64())));
327 });
328 Self::new(unwrapped_phase, self.get_shape()?)
329 } else {
330 let axis = axis.unwrap_or(-1);
331 let axis = self.normalize_axis(axis);
332
333 let b_shape = self.get_shape()?.update_at(axis, 1);
334 let period = period.broadcast_to(b_shape.clone())?;
335 let discont = discont.broadcast_to(b_shape)?;
336
337 let parts = self.get_shape()?.remove_at(axis).into_iter().product();
338 let partial = self
339 .moveaxis(vec![axis.to_isize()], vec![self.ndim()?.to_isize()])?
340 .ravel().split(parts, None)?;
341
342 let parameter_parts = period.get_shape()?.remove_at(axis).into_iter().product();
343 let period_partial = period
344 .moveaxis(vec![axis.to_isize()], vec![self.ndim()?.to_isize()])?
345 .ravel().split(parameter_parts, None)?;
346 let discont_partial = discont
347 .moveaxis(vec![axis.to_isize()], vec![self.ndim()?.to_isize()])?
348 .ravel().split(parameter_parts, None)?;
349
350 let mut results = vec![];
351 for i in 0..partial.len() {
352 let discont_p = Some(discont_partial[i].to_array_f64()?);
353 let period_p = Some(period_partial[i].to_array_f64()?);
354 results.push(partial[i].unwrap_phase(discont_p, None, period_p)?);
355 };
356 let mut tmp_shape = self.get_shape()?;
357 let tmp_to_push = tmp_shape.remove(axis);
358 tmp_shape.push(tmp_to_push);
359
360 let result = results.into_iter()
361 .flatten()
362 .collect::<Self>()
363 .reshape(&tmp_shape);
364 if axis == 0 { result.rollaxis((self.ndim()? - 1).to_isize(), None) }
365 else { result.moveaxis(vec![axis.to_isize()], vec![self.ndim()?.to_isize()]) }
366 .reshape(&self.get_shape()?)
367 }
368 }
369}
370
371impl <N: NumericOps> ArrayTrigonometric<N> for Result<Array<N>, ArrayError> {
372
373 fn sin(&self) -> Self {
374 self.clone()?.sin()
375 }
376
377 fn cos(&self) -> Self {
378 self.clone()?.cos()
379 }
380
381 fn tan(&self) -> Self {
382 self.clone()?.tan()
383 }
384
385 fn asin(&self) -> Self {
386 self.clone()?.asin()
387 }
388
389 fn acos(&self) -> Self {
390 self.clone()?.acos()
391 }
392
393 fn atan(&self) -> Self {
394 self.clone()?.atan()
395 }
396
397 fn atan2(&self, other: &Array<N>) -> Self {
398 self.clone()?.atan2(other)
399 }
400
401 fn hypot(&self, other: &Array<N>) -> Self {
402 self.clone()?.hypot(other)
403 }
404
405 fn degrees(&self) -> Self {
406 self.clone()?.degrees()
407 }
408
409 fn rad2deg(&self) -> Self {
410 self.clone()?.rad2deg()
411 }
412
413 fn radians(&self) -> Self {
414 self.clone()?.radians()
415 }
416
417 fn deg2rad(&self) -> Self {
418 self.clone()?.deg2rad()
419 }
420
421 fn unwrap_phase(&self, discont: Option<Array<f64>>, axis: Option<isize>, period: Option<Array<f64>>) -> Self {
422 self.clone()?.unwrap_phase(discont, axis, period)
423 }
424}