rgsl/types/
complex.rs

1//
2// A rust binding for the GSL library by Guillaume Gomez (guillaume1.gomez@gmail.com)
3//
4
5// TODO : port to Rust type : http://doc.rust-lang.org/num/complex/struct.Complex.html
6
7use std::default::Default;
8use std::fmt;
9use std::fmt::{Debug, Formatter};
10
11#[doc(hidden)]
12#[allow(clippy::upper_case_acronyms)]
13pub trait CFFI<T> {
14    fn wrap(s: T) -> Self;
15    fn unwrap(self) -> T;
16}
17
18#[doc(hidden)]
19#[allow(clippy::upper_case_acronyms)]
20pub trait FFFI<T> {
21    fn wrap(self) -> T;
22    fn unwrap(t: T) -> Self;
23}
24
25#[repr(C)]
26#[derive(Clone, Copy, PartialEq)]
27pub struct ComplexF64 {
28    pub dat: [f64; 2],
29}
30
31impl ComplexF64 {
32    /// This function uses the rectangular Cartesian components (x,y) to return the complex number
33    /// z = x + i y.
34    #[doc(alias = "gsl_complex_rect")]
35    pub fn rect(x: f64, y: f64) -> ComplexF64 {
36        unsafe { ::sys::gsl_complex_rect(x, y).wrap() }
37    }
38
39    /// This function returns the complex number z = r \exp(i \theta) = r (\cos(\theta) + i
40    /// \sin(\theta)) from the polar representation (r,theta).
41    #[doc(alias = "gsl_complex_polar")]
42    pub fn polar(r: f64, theta: f64) -> ComplexF64 {
43        unsafe { ::sys::gsl_complex_polar(r, theta).wrap() }
44    }
45
46    /// This function returns the argument of the complex number z, \arg(z), where -\pi < \arg(z)
47    /// <= \pi.
48    #[doc(alias = "gsl_complex_arg")]
49    pub fn arg(&self) -> f64 {
50        unsafe { ::sys::gsl_complex_arg(self.unwrap()) }
51    }
52
53    /// This function returns the magnitude of the complex number z, |z|.
54    #[doc(alias = "gsl_complex_abs")]
55    pub fn abs(&self) -> f64 {
56        unsafe { ::sys::gsl_complex_abs(self.unwrap()) }
57    }
58
59    /// This function returns the squared magnitude of the complex number z, |z|^2.
60    #[doc(alias = "gsl_complex_abs2")]
61    pub fn abs2(&self) -> f64 {
62        unsafe { ::sys::gsl_complex_abs2(self.unwrap()) }
63    }
64
65    /// This function returns the natural logarithm of the magnitude of the complex number z,
66    /// \log|z|.
67    ///
68    /// It allows an accurate evaluation of \log|z| when |z| is close to one.
69    ///
70    /// The direct evaluation of log(gsl_complex_abs(z)) would lead to a loss of precision in this
71    /// case.
72    #[doc(alias = "gsl_complex_logabs")]
73    pub fn logabs(&self) -> f64 {
74        unsafe { ::sys::gsl_complex_logabs(self.unwrap()) }
75    }
76
77    /// This function returns the sum of the complex numbers a and b, z=a+b.
78    #[doc(alias = "gsl_complex_add")]
79    pub fn add(&self, other: &ComplexF64) -> ComplexF64 {
80        unsafe { ::sys::gsl_complex_add(self.unwrap(), other.unwrap()).wrap() }
81    }
82
83    /// This function returns the difference of the complex numbers a and b, z=a-b.
84    #[doc(alias = "gsl_complex_sub")]
85    pub fn sub(&self, other: &ComplexF64) -> ComplexF64 {
86        unsafe { ::sys::gsl_complex_sub(self.unwrap(), other.unwrap()).wrap() }
87    }
88
89    /// This function returns the product of the complex numbers a and b, z=ab.
90    #[doc(alias = "gsl_complex_mul")]
91    pub fn mul(&self, other: &ComplexF64) -> ComplexF64 {
92        unsafe { ::sys::gsl_complex_mul(self.unwrap(), other.unwrap()).wrap() }
93    }
94
95    /// This function returns the quotient of the complex numbers a and b, z=a/b.
96    #[doc(alias = "gsl_complex_div")]
97    pub fn div(&self, other: &ComplexF64) -> ComplexF64 {
98        unsafe { ::sys::gsl_complex_div(self.unwrap(), other.unwrap()).wrap() }
99    }
100
101    /// This function returns the sum of the complex number a and the real number x, z=a+x.
102    #[doc(alias = "gsl_complex_add_real")]
103    pub fn add_real(&self, x: f64) -> ComplexF64 {
104        unsafe { ::sys::gsl_complex_add_real(self.unwrap(), x).wrap() }
105    }
106
107    /// This function returns the difference of the complex number a and the real number x, z=a-x.
108    #[doc(alias = "gsl_complex_sub_real")]
109    pub fn sub_real(&self, x: f64) -> ComplexF64 {
110        unsafe { ::sys::gsl_complex_sub_real(self.unwrap(), x).wrap() }
111    }
112
113    /// This function returns the product of the complex number a and the real number x, z=ax.
114    #[doc(alias = "gsl_complex_mul_real")]
115    pub fn mul_real(&self, x: f64) -> ComplexF64 {
116        unsafe { ::sys::gsl_complex_mul_real(self.unwrap(), x).wrap() }
117    }
118
119    /// This function returns the quotient of the complex number a and the real number x, z=a/x.
120    #[doc(alias = "gsl_complex_div_real")]
121    pub fn div_real(&self, x: f64) -> ComplexF64 {
122        unsafe { ::sys::gsl_complex_div_real(self.unwrap(), x).wrap() }
123    }
124
125    /// This function returns the sum of the complex number a and the imaginary number iy, z=a+iy.
126    #[doc(alias = "gsl_complex_add_imag")]
127    pub fn add_imag(&self, x: f64) -> ComplexF64 {
128        unsafe { ::sys::gsl_complex_add_imag(self.unwrap(), x).wrap() }
129    }
130
131    /// This function returns the difference of the complex number a and the imaginary number iy,
132    /// z=a-iy.
133    #[doc(alias = "gsl_complex_sub_imag")]
134    pub fn sub_imag(&self, x: f64) -> ComplexF64 {
135        unsafe { ::sys::gsl_complex_sub_imag(self.unwrap(), x).wrap() }
136    }
137
138    /// This function returns the product of the complex number a and the imaginary number iy,
139    /// z=a*(iy).
140    #[doc(alias = "gsl_complex_mul_imag")]
141    pub fn mul_imag(&self, x: f64) -> ComplexF64 {
142        unsafe { ::sys::gsl_complex_mul_imag(self.unwrap(), x).wrap() }
143    }
144
145    /// This function returns the quotient of the complex number a and the imaginary number iy,
146    /// z=a/(iy).
147    #[doc(alias = "gsl_complex_div_imag")]
148    pub fn div_imag(&self, x: f64) -> ComplexF64 {
149        unsafe { ::sys::gsl_complex_div_imag(self.unwrap(), x).wrap() }
150    }
151
152    /// This function returns the complex conjugate of the complex number z, z^* = x - i y.
153    #[doc(alias = "gsl_complex_conjugate")]
154    pub fn conjugate(&self) -> ComplexF64 {
155        unsafe { ::sys::gsl_complex_conjugate(self.unwrap()).wrap() }
156    }
157
158    /// This function returns the inverse, or reciprocal, of the complex number z, 1/z = (x - i y)/
159    /// (x^2 + y^2).
160    #[doc(alias = "gsl_complex_inverse")]
161    pub fn inverse(&self) -> ComplexF64 {
162        unsafe { ::sys::gsl_complex_inverse(self.unwrap()).wrap() }
163    }
164
165    /// This function returns the negative of the complex number z, -z = (-x) + i(-y).
166    #[doc(alias = "gsl_complex_negative")]
167    pub fn negative(&self) -> ComplexF64 {
168        unsafe { ::sys::gsl_complex_negative(self.unwrap()).wrap() }
169    }
170
171    /// This function returns the square root of the complex number z, \sqrt z.
172    ///
173    /// The branch cut is the negative real axis. The result always lies in the right half of the
174    /// omplex plane.
175    #[doc(alias = "gsl_complex_sqrt")]
176    pub fn sqrt(&self) -> ComplexF64 {
177        unsafe { ::sys::gsl_complex_sqrt(self.unwrap()).wrap() }
178    }
179
180    /// This function returns the complex square root of the real number x, where x may be negative.
181    #[doc(alias = "gsl_complex_sqrt_real")]
182    pub fn sqrt_real(x: f64) -> ComplexF64 {
183        unsafe { ::sys::gsl_complex_sqrt_real(x).wrap() }
184    }
185
186    /// The function returns the complex number z raised to the complex power a, z^a.
187    /// This is computed as \exp(\log(z)*a) using complex logarithms and complex exponentials.
188    #[doc(alias = "gsl_complex_pow")]
189    pub fn pow(&self, other: &ComplexF64) -> ComplexF64 {
190        unsafe { ::sys::gsl_complex_pow(self.unwrap(), other.unwrap()).wrap() }
191    }
192
193    /// This function returns the complex number z raised to the real power x, z^x.
194    #[doc(alias = "gsl_complex_pow_real")]
195    pub fn pow_real(&self, x: f64) -> ComplexF64 {
196        unsafe { ::sys::gsl_complex_pow_real(self.unwrap(), x).wrap() }
197    }
198
199    /// This function returns the complex exponential of the complex number z, \exp(z).
200    #[doc(alias = "gsl_complex_exp")]
201    pub fn exp(&self) -> ComplexF64 {
202        unsafe { ::sys::gsl_complex_exp(self.unwrap()).wrap() }
203    }
204
205    /// This function returns the complex natural logarithm (base e) of the complex number z,
206    /// \log(z).
207    ///
208    /// The branch cut is the negative real axis.
209    #[doc(alias = "gsl_complex_log")]
210    pub fn log(&self) -> ComplexF64 {
211        unsafe { ::sys::gsl_complex_log(self.unwrap()).wrap() }
212    }
213
214    /// This function returns the complex base-10 logarithm of the complex number z, \log_10 (z).
215    #[doc(alias = "gsl_complex_log10")]
216    pub fn log10(&self) -> ComplexF64 {
217        unsafe { ::sys::gsl_complex_log10(self.unwrap()).wrap() }
218    }
219
220    /// This function returns the complex base-b logarithm of the complex number z, \log_b(z).
221    /// This quantity is computed as the ratio \log(z)/\log(b).
222    #[doc(alias = "gsl_complex_log_b")]
223    pub fn log_b(&self, other: &ComplexF64) -> ComplexF64 {
224        unsafe { ::sys::gsl_complex_log_b(self.unwrap(), other.unwrap()).wrap() }
225    }
226
227    /// This function returns the complex sine of the complex number z, \sin(z) = (\exp(iz) -
228    /// \exp(-iz))/(2i).
229    #[doc(alias = "gsl_complex_sin")]
230    pub fn sin(&self) -> ComplexF64 {
231        unsafe { ::sys::gsl_complex_sin(self.unwrap()).wrap() }
232    }
233
234    /// This function returns the complex cosine of the complex number z, \cos(z) = (\exp(iz) +
235    /// \exp(-iz))/2.
236    #[doc(alias = "gsl_complex_cos")]
237    pub fn cos(&self) -> ComplexF64 {
238        unsafe { ::sys::gsl_complex_cos(self.unwrap()).wrap() }
239    }
240
241    /// This function returns the complex tangent of the complex number z, \tan(z) =
242    /// \sin(z)/\cos(z).
243    #[doc(alias = "gsl_complex_tan")]
244    pub fn tan(&self) -> ComplexF64 {
245        unsafe { ::sys::gsl_complex_tan(self.unwrap()).wrap() }
246    }
247
248    /// This function returns the complex secant of the complex number z, \sec(z) = 1/\cos(z).
249    #[doc(alias = "gsl_complex_sec")]
250    pub fn sec(&self) -> ComplexF64 {
251        unsafe { ::sys::gsl_complex_sec(self.unwrap()).wrap() }
252    }
253
254    /// This function returns the complex cosecant of the complex number z, \csc(z) = 1/\sin(z).
255    #[doc(alias = "gsl_complex_csc")]
256    pub fn csc(&self) -> ComplexF64 {
257        unsafe { ::sys::gsl_complex_csc(self.unwrap()).wrap() }
258    }
259
260    /// This function returns the complex cotangent of the complex number z, \cot(z) = 1/\tan(z).
261    #[doc(alias = "gsl_complex_cot")]
262    pub fn cot(&self) -> ComplexF64 {
263        unsafe { ::sys::gsl_complex_cot(self.unwrap()).wrap() }
264    }
265
266    /// This function returns the complex arcsine of the complex number z, \arcsin(z).
267    /// The branch cuts are on the real axis, less than -1 and greater than 1.
268    #[doc(alias = "gsl_complex_arcsin")]
269    pub fn arcsin(&self) -> ComplexF64 {
270        unsafe { ::sys::gsl_complex_arcsin(self.unwrap()).wrap() }
271    }
272
273    /// This function returns the complex arcsine of the real number z, \arcsin(z).
274    ///
275    /// * For z between -1 and 1, the function returns a real value in the range [-\pi/2,\pi/2].
276    /// * For z less than -1 the result has a real part of -\pi/2 and a positive imaginary part.
277    /// * For z greater than 1 the result has a real part of \pi/2 and a negative imaginary part.
278    #[doc(alias = "gsl_complex_arcsin_real")]
279    pub fn arcsin_real(z: f64) -> ComplexF64 {
280        unsafe { ::sys::gsl_complex_arcsin_real(z).wrap() }
281    }
282
283    /// This function returns the complex arccosine of the complex number z, \arccos(z).
284    /// The branch cuts are on the real axis, less than -1 and greater than 1.
285    #[doc(alias = "gsl_complex_arccos")]
286    pub fn arccos(&self) -> ComplexF64 {
287        unsafe { ::sys::gsl_complex_arccos(self.unwrap()).wrap() }
288    }
289
290    /// This function returns the complex arccosine of the real number z, \arccos(z).
291    ///
292    /// * For z between -1 and 1, the function returns a real value in the range [0,\pi].
293    /// * For z less than -1 the result has a real part of \pi and a negative imaginary part.
294    /// * For z greater than 1 the result is purely imaginary and positive.
295    #[doc(alias = "gsl_complex_arccos_real")]
296    pub fn arccos_real(z: f64) -> ComplexF64 {
297        unsafe { ::sys::gsl_complex_arccos_real(z).wrap() }
298    }
299
300    /// This function returns the complex arctangent of the complex number z, \arctan(z).
301    /// The branch cuts are on the imaginary axis, below -i and above i.
302    #[doc(alias = "gsl_complex_arctan")]
303    pub fn arctan(&self) -> ComplexF64 {
304        unsafe { ::sys::gsl_complex_arctan(self.unwrap()).wrap() }
305    }
306
307    /// This function returns the complex arcsecant of the complex number z, \arcsec(z) =
308    /// \arccos(1/z).
309    #[doc(alias = "gsl_complex_arcsec")]
310    pub fn arcsec(&self) -> ComplexF64 {
311        unsafe { ::sys::gsl_complex_arcsec(self.unwrap()).wrap() }
312    }
313
314    /// This function returns the complex arcsecant of the real number z, \arcsec(z) = \arccos(1/z).
315    #[doc(alias = "gsl_complex_arcsec_real")]
316    pub fn arcsec_real(z: f64) -> ComplexF64 {
317        unsafe { ::sys::gsl_complex_arcsec_real(z).wrap() }
318    }
319
320    /// This function returns the complex arccosecant of the complex number z, \arccsc(z) =
321    /// \arcsin(1/z).
322    #[doc(alias = "gsl_complex_arccsc")]
323    pub fn arccsc(&self) -> ComplexF64 {
324        unsafe { ::sys::gsl_complex_arccsc(self.unwrap()).wrap() }
325    }
326
327    /// This function returns the complex arccosecant of the real number z, \arccsc(z) =
328    /// \arcsin(1/z).
329    #[doc(alias = "gsl_complex_arccsc_real")]
330    pub fn arccsc_real(z: f64) -> ComplexF64 {
331        unsafe { ::sys::gsl_complex_arccsc_real(z).wrap() }
332    }
333
334    /// This function returns the complex arccotangent of the complex number z, \arccot(z) =
335    /// \arctan(1/z).
336    #[doc(alias = "gsl_complex_arccot")]
337    pub fn arccot(&self) -> ComplexF64 {
338        unsafe { ::sys::gsl_complex_arccot(self.unwrap()).wrap() }
339    }
340
341    /// This function returns the complex hyperbolic sine of the complex number z, \sinh(z) =
342    /// (\exp(z) - \exp(-z))/2.
343    #[doc(alias = "gsl_complex_sinh")]
344    pub fn sinh(&self) -> ComplexF64 {
345        unsafe { ::sys::gsl_complex_sinh(self.unwrap()).wrap() }
346    }
347
348    /// This function returns the complex hyperbolic cosine of the complex number z, \cosh(z) =
349    /// (\exp(z) + \exp(-z))/2.
350    #[doc(alias = "gsl_complex_cosh")]
351    pub fn cosh(&self) -> ComplexF64 {
352        unsafe { ::sys::gsl_complex_cosh(self.unwrap()).wrap() }
353    }
354
355    /// This function returns the complex hyperbolic tangent of the complex number z, \tanh(z) =
356    /// \sinh(z)/\cosh(z).
357    #[doc(alias = "gsl_complex_tanh")]
358    pub fn tanh(&self) -> ComplexF64 {
359        unsafe { ::sys::gsl_complex_tanh(self.unwrap()).wrap() }
360    }
361
362    /// This function returns the complex hyperbolic secant of the complex number z, \sech(z) =
363    /// 1/\cosh(z).
364    #[doc(alias = "gsl_complex_sech")]
365    pub fn sech(&self) -> ComplexF64 {
366        unsafe { ::sys::gsl_complex_sech(self.unwrap()).wrap() }
367    }
368
369    /// This function returns the complex hyperbolic cosecant of the complex number z, \csch(z) =
370    /// 1/\sinh(z).
371    #[doc(alias = "gsl_complex_csch")]
372    pub fn csch(&self) -> ComplexF64 {
373        unsafe { ::sys::gsl_complex_csch(self.unwrap()).wrap() }
374    }
375
376    /// This function returns the complex hyperbolic cotangent of the complex number z, \coth(z) =
377    /// 1/\tanh(z).
378    #[doc(alias = "gsl_complex_coth")]
379    pub fn coth(&self) -> ComplexF64 {
380        unsafe { ::sys::gsl_complex_coth(self.unwrap()).wrap() }
381    }
382
383    /// This function returns the complex hyperbolic arcsine of the complex number z, \arcsinh(z).
384    /// The branch cuts are on the imaginary axis, below -i and above i.
385    #[doc(alias = "gsl_complex_arcsinh")]
386    pub fn arcsinh(&self) -> ComplexF64 {
387        unsafe { ::sys::gsl_complex_arcsinh(self.unwrap()).wrap() }
388    }
389
390    /// This function returns the complex hyperbolic arccosine of the complex number z, \arccosh(z).
391    /// The branch cut is on the real axis, less than 1.
392    /// Note that in this case we use the negative square root in formula 4.6.21 of Abramowitz &
393    /// Stegun giving \arccosh(z)=\log(z-\sqrt{z^2-1}).
394    #[doc(alias = "gsl_complex_arccosh")]
395    pub fn arccosh(&self) -> ComplexF64 {
396        unsafe { ::sys::gsl_complex_arccosh(self.unwrap()).wrap() }
397    }
398
399    /// This function returns the complex hyperbolic arccosine of the real number z, \arccosh(z).
400    #[doc(alias = "gsl_complex_arccosh_real")]
401    pub fn arccosh_real(z: f64) -> ComplexF64 {
402        unsafe { ::sys::gsl_complex_arccosh_real(z).wrap() }
403    }
404
405    /// This function returns the complex hyperbolic arctangent of the complex number z,
406    /// \arctanh(z).
407    ///
408    /// The branch cuts are on the real axis, less than -1 and greater than 1.
409    #[doc(alias = "gsl_complex_arctanh")]
410    pub fn arctanh(&self) -> ComplexF64 {
411        unsafe { ::sys::gsl_complex_arctanh(self.unwrap()).wrap() }
412    }
413
414    /// This function returns the complex hyperbolic arctangent of the real number z, \arctanh(z).
415    #[doc(alias = "gsl_complex_arctanh_real")]
416    pub fn arctanh_real(z: f64) -> ComplexF64 {
417        unsafe { ::sys::gsl_complex_arctanh_real(z).wrap() }
418    }
419
420    /// This function returns the complex hyperbolic arcsecant of the complex number z, \arcsech(z)
421    /// = \arccosh(1/z).
422    #[doc(alias = "gsl_complex_arcsech")]
423    pub fn arcsech(&self) -> ComplexF64 {
424        unsafe { ::sys::gsl_complex_arcsech(self.unwrap()).wrap() }
425    }
426
427    /// This function returns the complex hyperbolic arccosecant of the complex number z,
428    /// \arccsch(z) = \arcsin(1/z).
429    #[doc(alias = "gsl_complex_arccsch")]
430    pub fn arccsch(&self) -> ComplexF64 {
431        unsafe { ::sys::gsl_complex_arccsch(self.unwrap()).wrap() }
432    }
433
434    /// This function returns the complex hyperbolic arccotangent of the complex number z,
435    /// \arccoth(z) = \arctanh(1/z).
436    #[doc(alias = "gsl_complex_arccoth")]
437    pub fn arccoth(&self) -> ComplexF64 {
438        unsafe { ::sys::gsl_complex_arccoth(self.unwrap()).wrap() }
439    }
440
441    pub fn real(&self) -> f64 {
442        self.dat[0]
443    }
444
445    pub fn imaginary(&self) -> f64 {
446        self.dat[1]
447    }
448}
449
450impl Debug for ComplexF64 {
451    fn fmt(&self, f: &mut Formatter) -> fmt::Result {
452        write!(f, "[{}, {}]", self.dat[0], self.dat[1])
453    }
454}
455
456impl Default for ComplexF64 {
457    fn default() -> ComplexF64 {
458        ComplexF64 { dat: [0f64, 0f64] }
459    }
460}
461
462impl CFFI<sys::gsl_complex> for ComplexF64 {
463    fn wrap(t: sys::gsl_complex) -> ComplexF64 {
464        unsafe { std::mem::transmute(t) }
465    }
466
467    fn unwrap(self) -> sys::gsl_complex {
468        unsafe { std::mem::transmute(self) }
469    }
470}
471
472impl CFFI<sys::gsl_complex_float> for ComplexF64 {
473    fn wrap(t: sys::gsl_complex_float) -> ComplexF64 {
474        ComplexF64 {
475            dat: [t.dat[0] as f64, t.dat[1] as f64],
476        }
477    }
478
479    fn unwrap(self) -> sys::gsl_complex_float {
480        sys::gsl_complex_float {
481            dat: [self.dat[0] as f32, self.dat[1] as f32],
482        }
483    }
484}
485
486impl FFFI<ComplexF32> for sys::gsl_complex {
487    fn wrap(self) -> ComplexF32 {
488        ComplexF32 {
489            dat: [self.dat[0] as f32, self.dat[1] as f32],
490        }
491    }
492
493    fn unwrap(t: ComplexF32) -> sys::gsl_complex {
494        sys::gsl_complex {
495            dat: [t.dat[0] as f64, t.dat[1] as f64],
496        }
497    }
498}
499
500impl FFFI<ComplexF64> for sys::gsl_complex {
501    fn wrap(self) -> ComplexF64 {
502        unsafe { std::mem::transmute(self) }
503    }
504
505    fn unwrap(t: ComplexF64) -> sys::gsl_complex {
506        unsafe { std::mem::transmute(t) }
507    }
508}
509
510#[repr(C)]
511#[derive(Clone, Copy, PartialEq)]
512pub struct ComplexF32 {
513    pub dat: [f32; 2],
514}
515
516impl ComplexF32 {
517    /// This function uses the rectangular Cartesian components (x,y) to return the complex number
518    /// z = x + i y.
519    #[doc(alias = "gsl_complex_rect")]
520    pub fn rect(x: f32, y: f32) -> ComplexF32 {
521        unsafe { ::sys::gsl_complex_rect(x as f64, y as f64).wrap() }
522    }
523
524    /// This function returns the complex number z = r \exp(i \theta) = r (\cos(\theta) + i
525    /// \sin(\theta)) from the polar representation (r,theta).
526    #[doc(alias = "gsl_complex_polar")]
527    pub fn polar(r: f32, theta: f32) -> ComplexF32 {
528        unsafe { ::sys::gsl_complex_polar(r as f64, theta as f64).wrap() }
529    }
530
531    /// This function returns the argument of the complex number z, \arg(z), where -\pi < \arg(z)
532    /// <= \pi.
533    #[doc(alias = "gsl_complex_arg")]
534    pub fn arg(&self) -> f32 {
535        unsafe { ::sys::gsl_complex_arg(self.unwrap()) as f32 }
536    }
537
538    /// This function returns the magnitude of the complex number z, |z|.
539    #[doc(alias = "gsl_complex_abs")]
540    pub fn abs(&self) -> f32 {
541        unsafe { ::sys::gsl_complex_abs(self.unwrap()) as f32 }
542    }
543
544    /// This function returns the squared magnitude of the complex number z, |z|^2.
545    #[doc(alias = "gsl_complex_abs2")]
546    pub fn abs2(&self) -> f32 {
547        unsafe { ::sys::gsl_complex_abs2(self.unwrap()) as f32 }
548    }
549
550    /// This function returns the natural logarithm of the magnitude of the complex number z,
551    /// \log|z|.
552    ///
553    /// It allows an accurate evaluation of \log|z| when |z| is close to one.
554    /// The direct evaluation of log(gsl_complex_abs(z)) would lead to a loss of precision in
555    /// this case.
556    #[doc(alias = "gsl_complex_logabs")]
557    pub fn logabs(&self) -> f32 {
558        unsafe { ::sys::gsl_complex_logabs(self.unwrap()) as f32 }
559    }
560
561    /// This function returns the sum of the complex numbers a and b, z=a+b.
562    #[doc(alias = "gsl_complex_add")]
563    pub fn add(&self, other: &ComplexF32) -> ComplexF32 {
564        unsafe { ::sys::gsl_complex_add(self.unwrap(), other.unwrap()).wrap() }
565    }
566
567    /// This function returns the difference of the complex numbers a and b, z=a-b.
568    #[doc(alias = "gsl_complex_sub")]
569    pub fn sub(&self, other: &ComplexF32) -> ComplexF32 {
570        unsafe { ::sys::gsl_complex_sub(self.unwrap(), other.unwrap()).wrap() }
571    }
572
573    /// This function returns the product of the complex numbers a and b, z=ab.
574    #[doc(alias = "gsl_complex_mul")]
575    pub fn mul(&self, other: &ComplexF32) -> ComplexF32 {
576        unsafe { ::sys::gsl_complex_mul(self.unwrap(), other.unwrap()).wrap() }
577    }
578
579    /// This function returns the quotient of the complex numbers a and b, z=a/b.
580    #[doc(alias = "gsl_complex_div")]
581    pub fn div(&self, other: &ComplexF32) -> ComplexF32 {
582        unsafe { ::sys::gsl_complex_div(self.unwrap(), other.unwrap()).wrap() }
583    }
584
585    /// This function returns the sum of the complex number a and the real number x, z=a+x.
586    #[doc(alias = "gsl_complex_add_real")]
587    pub fn add_real(&self, x: f32) -> ComplexF32 {
588        unsafe { ::sys::gsl_complex_add_real(self.unwrap(), x as f64).wrap() }
589    }
590
591    /// This function returns the difference of the complex number a and the real number x, z=a-x.
592    #[doc(alias = "gsl_complex_sub_real")]
593    pub fn sub_real(&self, x: f32) -> ComplexF32 {
594        unsafe { ::sys::gsl_complex_sub_real(self.unwrap(), x as f64).wrap() }
595    }
596
597    /// This function returns the product of the complex number a and the real number x, z=ax.
598    #[doc(alias = "gsl_complex_mul_real")]
599    pub fn mul_real(&self, x: f32) -> ComplexF32 {
600        unsafe { ::sys::gsl_complex_mul_real(self.unwrap(), x as f64).wrap() }
601    }
602
603    /// This function returns the quotient of the complex number a and the real number x, z=a/x.
604    #[doc(alias = "gsl_complex_div_real")]
605    pub fn div_real(&self, x: f32) -> ComplexF32 {
606        unsafe { ::sys::gsl_complex_div_real(self.unwrap(), x as f64).wrap() }
607    }
608
609    /// This function returns the sum of the complex number a and the imaginary number iy, z=a+iy.
610    #[doc(alias = "gsl_complex_add_imag")]
611    pub fn add_imag(&self, x: f32) -> ComplexF32 {
612        unsafe { ::sys::gsl_complex_add_imag(self.unwrap(), x as f64).wrap() }
613    }
614
615    /// This function returns the difference of the complex number a and the imaginary number iy, z=a-iy.
616    #[doc(alias = "gsl_complex_sub_imag")]
617    pub fn sub_imag(&self, x: f32) -> ComplexF32 {
618        unsafe { ::sys::gsl_complex_sub_imag(self.unwrap(), x as f64).wrap() }
619    }
620
621    /// This function returns the product of the complex number a and the imaginary number iy, z=a*(iy).
622    #[doc(alias = "gsl_complex_mul_imag")]
623    pub fn mul_imag(&self, x: f32) -> ComplexF32 {
624        unsafe { ::sys::gsl_complex_mul_imag(self.unwrap(), x as f64).wrap() }
625    }
626
627    /// This function returns the quotient of the complex number a and the imaginary number iy, z=a/(iy).
628    #[doc(alias = "gsl_complex_div_imag")]
629    pub fn div_imag(&self, x: f32) -> ComplexF32 {
630        unsafe { ::sys::gsl_complex_div_imag(self.unwrap(), x as f64).wrap() }
631    }
632
633    /// This function returns the complex conjugate of the complex number z, z^* = x - i y.
634    #[doc(alias = "gsl_complex_conjugate")]
635    pub fn conjugate(&self) -> ComplexF32 {
636        unsafe { ::sys::gsl_complex_conjugate(self.unwrap()).wrap() }
637    }
638
639    /// This function returns the inverse, or reciprocal, of the complex number z, 1/z = (x - i y)/
640    /// (x^2 + y^2).
641    #[doc(alias = "gsl_complex_inverse")]
642    pub fn inverse(&self) -> ComplexF32 {
643        unsafe { ::sys::gsl_complex_inverse(self.unwrap()).wrap() }
644    }
645
646    /// This function returns the negative of the complex number z, -z = (-x) + i(-y).
647    #[doc(alias = "gsl_complex_negative")]
648    pub fn negative(&self) -> ComplexF32 {
649        unsafe { ::sys::gsl_complex_negative(self.unwrap()).wrap() }
650    }
651
652    /// This function returns the square root of the complex number z, \sqrt z.
653    ///
654    /// The branch cut is the negative real axis. The result always lies in the right half of the
655    /// complex plane.
656    #[doc(alias = "gsl_complex_sqrt")]
657    pub fn sqrt(&self) -> ComplexF32 {
658        unsafe { ::sys::gsl_complex_sqrt(self.unwrap()).wrap() }
659    }
660
661    /// This function returns the complex square root of the real number x, where x may be negative.
662    #[doc(alias = "gsl_complex_sqrt_real")]
663    pub fn sqrt_real(x: f32) -> ComplexF32 {
664        unsafe { ::sys::gsl_complex_sqrt_real(x as f64).wrap() }
665    }
666
667    /// The function returns the complex number z raised to the complex power a, z^a.
668    ///
669    /// This is computed as \exp(\log(z)*a) using complex logarithms and complex exponentials.
670    #[doc(alias = "gsl_complex_pow")]
671    pub fn pow(&self, other: &ComplexF32) -> ComplexF32 {
672        unsafe { ::sys::gsl_complex_pow(self.unwrap(), other.unwrap()).wrap() }
673    }
674
675    /// This function returns the complex number z raised to the real power x, z^x.
676    #[doc(alias = "gsl_complex_pow_real")]
677    pub fn pow_real(&self, x: f32) -> ComplexF32 {
678        unsafe { ::sys::gsl_complex_pow_real(self.unwrap(), x as f64).wrap() }
679    }
680
681    /// This function returns the complex exponential of the complex number z, \exp(z).
682    #[doc(alias = "gsl_complex_exp")]
683    pub fn exp(&self) -> ComplexF32 {
684        unsafe { ::sys::gsl_complex_exp(self.unwrap()).wrap() }
685    }
686
687    /// This function returns the complex natural logarithm (base e) of the complex number z, \log(z).
688    /// The branch cut is the negative real axis.
689    #[doc(alias = "gsl_complex_log")]
690    pub fn log(&self) -> ComplexF32 {
691        unsafe { ::sys::gsl_complex_log(self.unwrap()).wrap() }
692    }
693
694    /// This function returns the complex base-10 logarithm of the complex number z, \log_10 (z).
695    #[doc(alias = "gsl_complex_log10")]
696    pub fn log10(&self) -> ComplexF32 {
697        unsafe { ::sys::gsl_complex_log10(self.unwrap()).wrap() }
698    }
699
700    /// This function returns the complex base-b logarithm of the complex number z, \log_b(z).
701    /// This quantity is computed as the ratio \log(z)/\log(b).
702    #[doc(alias = "gsl_complex_log_b")]
703    pub fn log_b(&self, other: &ComplexF32) -> ComplexF32 {
704        unsafe { ::sys::gsl_complex_log_b(self.unwrap(), other.unwrap()).wrap() }
705    }
706
707    /// This function returns the complex sine of the complex number z, \sin(z) = (\exp(iz) -
708    /// \exp(-iz))/(2i).
709    #[doc(alias = "gsl_complex_sin")]
710    pub fn sin(&self) -> ComplexF32 {
711        unsafe { ::sys::gsl_complex_sin(self.unwrap()).wrap() }
712    }
713
714    /// This function returns the complex cosine of the complex number z, \cos(z) = (\exp(iz) +
715    /// \exp(-iz))/2.
716    #[doc(alias = "gsl_complex_cos")]
717    pub fn cos(&self) -> ComplexF32 {
718        unsafe { ::sys::gsl_complex_cos(self.unwrap()).wrap() }
719    }
720
721    /// This function returns the complex tangent of the complex number z, \tan(z) =
722    /// \sin(z)/\cos(z).
723    #[doc(alias = "gsl_complex_tan")]
724    pub fn tan(&self) -> ComplexF32 {
725        unsafe { ::sys::gsl_complex_tan(self.unwrap()).wrap() }
726    }
727
728    /// This function returns the complex secant of the complex number z, \sec(z) = 1/\cos(z).
729    #[doc(alias = "gsl_complex_sec")]
730    pub fn sec(&self) -> ComplexF32 {
731        unsafe { ::sys::gsl_complex_sec(self.unwrap()).wrap() }
732    }
733
734    /// This function returns the complex cosecant of the complex number z, \csc(z) = 1/\sin(z).
735    #[doc(alias = "gsl_complex_csc")]
736    pub fn csc(&self) -> ComplexF32 {
737        unsafe { ::sys::gsl_complex_csc(self.unwrap()).wrap() }
738    }
739
740    /// This function returns the complex cotangent of the complex number z, \cot(z) = 1/\tan(z).
741    #[doc(alias = "gsl_complex_cot")]
742    pub fn cot(&self) -> ComplexF32 {
743        unsafe { ::sys::gsl_complex_cot(self.unwrap()).wrap() }
744    }
745
746    /// This function returns the complex arcsine of the complex number z, \arcsin(z).
747    /// The branch cuts are on the real axis, less than -1 and greater than 1.
748    #[doc(alias = "gsl_complex_arcsin")]
749    pub fn arcsin(&self) -> ComplexF32 {
750        unsafe { ::sys::gsl_complex_arcsin(self.unwrap()).wrap() }
751    }
752
753    /// This function returns the complex arcsine of the real number z, \arcsin(z).
754    ///
755    /// * For z between -1 and 1, the function returns a real value in the range [-\pi/2,\pi/2].
756    /// * For z less than -1 the result has a real part of -\pi/2 and a positive imaginary part.
757    /// * For z greater than 1 the result has a real part of \pi/2 and a negative imaginary part.
758    #[doc(alias = "gsl_complex_arcsin_real")]
759    pub fn arcsin_real(z: f32) -> ComplexF32 {
760        unsafe { ::sys::gsl_complex_arcsin_real(z as f64).wrap() }
761    }
762
763    /// This function returns the complex arccosine of the complex number z, \arccos(z).
764    /// The branch cuts are on the real axis, less than -1 and greater than 1.
765    #[doc(alias = "gsl_complex_arccos")]
766    pub fn arccos(&self) -> ComplexF32 {
767        unsafe { ::sys::gsl_complex_arccos(self.unwrap()).wrap() }
768    }
769
770    /// This function returns the complex arccosine of the real number z, \arccos(z).
771    ///
772    /// * For z between -1 and 1, the function returns a real value in the range [0,\pi].
773    /// * For z less than -1 the result has a real part of \pi and a negative imaginary part.
774    /// * For z greater than 1 the result is purely imaginary and positive.
775    #[doc(alias = "gsl_complex_arccos_real")]
776    pub fn arccos_real(z: f32) -> ComplexF32 {
777        unsafe { ::sys::gsl_complex_arccos_real(z as f64).wrap() }
778    }
779
780    /// This function returns the complex arctangent of the complex number z, \arctan(z).
781    /// The branch cuts are on the imaginary axis, below -i and above i.
782    #[doc(alias = "gsl_complex_arctan")]
783    pub fn arctan(&self) -> ComplexF32 {
784        unsafe { ::sys::gsl_complex_arctan(self.unwrap()).wrap() }
785    }
786
787    /// This function returns the complex arcsecant of the complex number z, \arcsec(z) =
788    /// \arccos(1/z).
789    #[doc(alias = "gsl_complex_arcsec")]
790    pub fn arcsec(&self) -> ComplexF32 {
791        unsafe { ::sys::gsl_complex_arcsec(self.unwrap()).wrap() }
792    }
793
794    /// This function returns the complex arcsecant of the real number z, \arcsec(z) = \arccos(1/z).
795    #[doc(alias = "gsl_complex_arcsec_real")]
796    pub fn arcsec_real(z: f32) -> ComplexF32 {
797        unsafe { ::sys::gsl_complex_arcsec_real(z as f64).wrap() }
798    }
799
800    /// This function returns the complex arccosecant of the complex number z, \arccsc(z) =
801    /// \arcsin(1/z).
802    #[doc(alias = "gsl_complex_arccsc")]
803    pub fn arccsc(&self) -> ComplexF32 {
804        unsafe { ::sys::gsl_complex_arccsc(self.unwrap()).wrap() }
805    }
806
807    /// This function returns the complex arccosecant of the real number z, \arccsc(z) =
808    /// \arcsin(1/z).
809    #[doc(alias = "gsl_complex_arccsc_real")]
810    pub fn arccsc_real(z: f32) -> ComplexF32 {
811        unsafe { ::sys::gsl_complex_arccsc_real(z as f64).wrap() }
812    }
813
814    /// This function returns the complex arccotangent of the complex number z, \arccot(z) =
815    /// \arctan(1/z).
816    #[doc(alias = "gsl_complex_arccot")]
817    pub fn arccot(&self) -> ComplexF32 {
818        unsafe { ::sys::gsl_complex_arccot(self.unwrap()).wrap() }
819    }
820
821    /// This function returns the complex hyperbolic sine of the complex number z, \sinh(z) =
822    /// (\exp(z) - \exp(-z))/2.
823    #[doc(alias = "gsl_complex_sinh")]
824    pub fn sinh(&self) -> ComplexF32 {
825        unsafe { ::sys::gsl_complex_sinh(self.unwrap()).wrap() }
826    }
827
828    /// This function returns the complex hyperbolic cosine of the complex number z, \cosh(z) =
829    /// (\exp(z) + \exp(-z))/2.
830    #[doc(alias = "gsl_complex_cosh")]
831    pub fn cosh(&self) -> ComplexF32 {
832        unsafe { ::sys::gsl_complex_cosh(self.unwrap()).wrap() }
833    }
834
835    /// This function returns the complex hyperbolic tangent of the complex number z, \tanh(z) =
836    /// \sinh(z)/\cosh(z).
837    #[doc(alias = "gsl_complex_tanh")]
838    pub fn tanh(&self) -> ComplexF32 {
839        unsafe { ::sys::gsl_complex_tanh(self.unwrap()).wrap() }
840    }
841
842    /// This function returns the complex hyperbolic secant of the complex number z, \sech(z) =
843    /// 1/\cosh(z).
844    #[doc(alias = "gsl_complex_sech")]
845    pub fn sech(&self) -> ComplexF32 {
846        unsafe { ::sys::gsl_complex_sech(self.unwrap()).wrap() }
847    }
848
849    /// This function returns the complex hyperbolic cosecant of the complex number z, \csch(z) =
850    /// 1/\sinh(z).
851    #[doc(alias = "gsl_complex_csch")]
852    pub fn csch(&self) -> ComplexF32 {
853        unsafe { ::sys::gsl_complex_csch(self.unwrap()).wrap() }
854    }
855
856    /// This function returns the complex hyperbolic cotangent of the complex number z, \coth(z) =
857    /// 1/\tanh(z).
858    #[doc(alias = "gsl_complex_coth")]
859    pub fn coth(&self) -> ComplexF32 {
860        unsafe { ::sys::gsl_complex_coth(self.unwrap()).wrap() }
861    }
862
863    /// This function returns the complex hyperbolic arcsine of the complex number z, \arcsinh(z).
864    /// The branch cuts are on the imaginary axis, below -i and above i.
865    #[doc(alias = "gsl_complex_arcsinh")]
866    pub fn arcsinh(&self) -> ComplexF32 {
867        unsafe { ::sys::gsl_complex_arcsinh(self.unwrap()).wrap() }
868    }
869
870    /// This function returns the complex hyperbolic arccosine of the complex number z, \arccosh(z).
871    ///
872    /// The branch cut is on the real axis, less than 1.
873    ///
874    /// Note that in this case we use the negative square root in formula 4.6.21 of Abramowitz &
875    /// Stegun giving \arccosh(z)=\log(z-\sqrt{z^2-1}).
876    #[doc(alias = "gsl_complex_arccosh")]
877    pub fn arccosh(&self) -> ComplexF32 {
878        unsafe { ::sys::gsl_complex_arccosh(self.unwrap()).wrap() }
879    }
880
881    /// This function returns the complex hyperbolic arccosine of the real number z, \arccosh(z).
882    #[doc(alias = "gsl_complex_arccosh_real")]
883    pub fn arccosh_real(z: f32) -> ComplexF32 {
884        unsafe { ::sys::gsl_complex_arccosh_real(z as f64).wrap() }
885    }
886
887    /// This function returns the complex hyperbolic arctangent of the complex number z,
888    /// arctanh(z).
889    ///
890    /// The branch cuts are on the real axis, less than -1 and greater than 1.
891    #[doc(alias = "gsl_complex_arctanh")]
892    pub fn arctanh(&self) -> ComplexF32 {
893        unsafe { ::sys::gsl_complex_arctanh(self.unwrap()).wrap() }
894    }
895
896    /// This function returns the complex hyperbolic arctangent of the real number z, \arctanh(z).
897    #[doc(alias = "gsl_complex_arctanh_real")]
898    pub fn arctanh_real(z: f32) -> ComplexF32 {
899        unsafe { ::sys::gsl_complex_arctanh_real(z as f64).wrap() }
900    }
901
902    /// This function returns the complex hyperbolic arcsecant of the complex number z, \arcsech(z)
903    /// = \arccosh(1/z).
904    #[doc(alias = "gsl_complex_arcsech")]
905    pub fn arcsech(&self) -> ComplexF32 {
906        unsafe { ::sys::gsl_complex_arcsech(self.unwrap()).wrap() }
907    }
908
909    /// This function returns the complex hyperbolic arccosecant of the complex number z,
910    /// \arccsch(z) = \arcsin(1/z).
911    #[doc(alias = "gsl_complex_arccsch")]
912    pub fn arccsch(&self) -> ComplexF32 {
913        unsafe { ::sys::gsl_complex_arccsch(self.unwrap()).wrap() }
914    }
915
916    /// This function returns the complex hyperbolic arccotangent of the complex number z,
917    /// \arccoth(z) = \arctanh(1/z).
918    #[doc(alias = "gsl_complex_arccoth")]
919    pub fn arccoth(&self) -> ComplexF32 {
920        unsafe { ::sys::gsl_complex_arccoth(self.unwrap()).wrap() }
921    }
922
923    pub fn real(&self) -> f32 {
924        self.dat[0]
925    }
926
927    pub fn imaginary(&self) -> f32 {
928        self.dat[1]
929    }
930}
931
932impl Debug for ComplexF32 {
933    fn fmt(&self, f: &mut Formatter) -> fmt::Result {
934        write!(f, "[{}, {}]", self.dat[0], self.dat[1])
935    }
936}
937
938impl Default for ComplexF32 {
939    fn default() -> ComplexF32 {
940        ComplexF32 { dat: [0f32, 0f32] }
941    }
942}
943
944impl CFFI<sys::gsl_complex> for ComplexF32 {
945    fn wrap(s: sys::gsl_complex) -> ComplexF32 {
946        ComplexF32 {
947            dat: [s.dat[0] as f32, s.dat[1] as f32],
948        }
949    }
950
951    fn unwrap(self) -> sys::gsl_complex {
952        sys::gsl_complex {
953            dat: [self.dat[0] as f64, self.dat[1] as f64],
954        }
955    }
956}
957
958impl CFFI<sys::gsl_complex_float> for ComplexF32 {
959    fn wrap(s: sys::gsl_complex_float) -> ComplexF32 {
960        unsafe { std::mem::transmute(s) }
961    }
962
963    fn unwrap(self) -> sys::gsl_complex_float {
964        unsafe { std::mem::transmute(self) }
965    }
966}
967
968// All these tests have been tested against the following C code:
969//
970// ```ignore
971// #include <gsl/gsl_complex.h>
972// #include <gsl/gsl_complex_math.h>
973// #include <gsl/gsl_inline.h>
974// #include <gsl/gsl_complex.h>
975//
976// void print_complex(gsl_complex *c, const char *text) {
977//   printf("%s: %f %f\n", text, c->dat[0], c->dat[1]);
978// }
979//
980// int main (void)
981// {
982//   gsl_complex c = gsl_complex_rect(10., 10.);
983//   gsl_complex c2 = gsl_complex_rect(1., -1.);
984//   print_complex(&c, "rect");
985//   print_complex(&c2, "rect");
986//   gsl_complex c3 = gsl_complex_polar(5., 7.);
987//   print_complex(&c3, "polar");
988//
989//   printf("-> %f\n", gsl_complex_arg(c3));
990//   printf("-> %f\n", gsl_complex_abs(c3));
991//   printf("-> %f\n", gsl_complex_abs2(c3));
992//   printf("-> %f\n", gsl_complex_logabs(c3));
993//
994//   {
995//     gsl_complex c4 = gsl_complex_add(c3, c2);
996//     print_complex(&c4, "\nadd");
997//   }
998//   {
999//     gsl_complex c4 = gsl_complex_sub(c3, c2);
1000//     print_complex(&c4, "sub");
1001//   }
1002//   {
1003//     gsl_complex c4 = gsl_complex_mul(c3, c2);
1004//     print_complex(&c4, "mul");
1005//   }
1006//   {
1007//     gsl_complex c4 = gsl_complex_div(c3, c2);
1008//     print_complex(&c4, "div");
1009//   }
1010//   {
1011//     gsl_complex c4 = gsl_complex_add_real(c3, 5.);
1012//     print_complex(&c4, "add_real");
1013//   }
1014//   {
1015//     gsl_complex c4 = gsl_complex_sub_real(c3, 5.);
1016//     print_complex(&c4, "sub_real");
1017//   }
1018//   {
1019//     gsl_complex c4 = gsl_complex_mul_real(c3, 5.);
1020//     print_complex(&c4, "mul_real");
1021//   }
1022//   {
1023//     gsl_complex c4 = gsl_complex_div_real(c3, 5.);
1024//     print_complex(&c4, "div_real");
1025//   }
1026//
1027//
1028//   {
1029//     gsl_complex c4 = gsl_complex_add_imag(c3, 5.);
1030//     print_complex(&c4, "\nadd_imag");
1031//   }
1032//   {
1033//     gsl_complex c4 = gsl_complex_sub_imag(c3, 5.);
1034//     print_complex(&c4, "sub_imag");
1035//   }
1036//   {
1037//     gsl_complex c4 = gsl_complex_mul_imag(c3, 5.);
1038//     print_complex(&c4, "mul_imag");
1039//   }
1040//   {
1041//     gsl_complex c4 = gsl_complex_div_imag(c3, 5.);
1042//     print_complex(&c4, "div_imag");
1043//   }
1044//
1045//
1046//   {
1047//     gsl_complex c4 = gsl_complex_conjugate(c3);
1048//     print_complex(&c4, "\nconjugate");
1049//   }
1050//   {
1051//     gsl_complex c4 = gsl_complex_inverse(c3);
1052//     print_complex(&c4, "inverse");
1053//   }
1054//   {
1055//     gsl_complex c4 = gsl_complex_negative(c3);
1056//     print_complex(&c4, "negative");
1057//   }
1058//   {
1059//     gsl_complex c4 = gsl_complex_sqrt(c3);
1060//     print_complex(&c4, "sqrt");
1061//   }
1062//   {
1063//     gsl_complex c4 = gsl_complex_sqrt_real(5.);
1064//     print_complex(&c4, "sqrt_real");
1065//   }
1066//
1067//
1068//   {
1069//     gsl_complex c4 = gsl_complex_pow(c3, c2);
1070//     print_complex(&c4, "\npow");
1071//   }
1072//   {
1073//     gsl_complex c4 = gsl_complex_pow_real(c3, 5.);
1074//     print_complex(&c4, "pow_real");
1075//   }
1076//   {
1077//     gsl_complex c4 = gsl_complex_exp(c3);
1078//     print_complex(&c4, "exp");
1079//   }
1080//   {
1081//     gsl_complex c4 = gsl_complex_log(c3);
1082//     print_complex(&c4, "log");
1083//   }
1084//   {
1085//     gsl_complex c4 = gsl_complex_log10(c3);
1086//     print_complex(&c4, "log10");
1087//   }
1088//   {
1089//     gsl_complex c4 = gsl_complex_log_b(c3, c2);
1090//     print_complex(&c4, "log_b");
1091//   }
1092//   {
1093//     gsl_complex c4 = gsl_complex_sin(c3);
1094//     print_complex(&c4, "sin");
1095//   }
1096//   {
1097//     gsl_complex c4 = gsl_complex_cos(c3);
1098//     print_complex(&c4, "cos");
1099//   }
1100//   {
1101//     gsl_complex c4 = gsl_complex_tan(c3);
1102//     print_complex(&c4, "tan");
1103//   }
1104//
1105//
1106//   {
1107//     gsl_complex c4 = gsl_complex_sec(c3);
1108//     print_complex(&c4, "\nsec");
1109//   }
1110//   {
1111//     gsl_complex c4 = gsl_complex_csc(c3);
1112//     print_complex(&c4, "csc");
1113//   }
1114//   {
1115//     gsl_complex c4 = gsl_complex_cot(c3);
1116//     print_complex(&c4, "cot");
1117//   }
1118//   {
1119//     gsl_complex c4 = gsl_complex_arcsin(c3);
1120//     print_complex(&c4, "arcsin");
1121//   }
1122//   {
1123//     gsl_complex c4 = gsl_complex_arcsin_real(5.);
1124//     print_complex(&c4, "arcsin_real");
1125//   }
1126//   {
1127//     gsl_complex c4 = gsl_complex_arccos(c3);
1128//     print_complex(&c4, "arccos");
1129//   }
1130//   {
1131//     gsl_complex c4 = gsl_complex_arccos_real(5.);
1132//     print_complex(&c4, "arccos_real");
1133//   }
1134//   {
1135//     gsl_complex c4 = gsl_complex_arctan(c3);
1136//     print_complex(&c4, "arctan");
1137//   }
1138//   {
1139//     gsl_complex c4 = gsl_complex_arcsec(c3);
1140//     print_complex(&c4, "arcsec");
1141//   }
1142//   {
1143//     gsl_complex c4 = gsl_complex_arcsec_real(5.);
1144//     print_complex(&c4, "arcsec_real");
1145//   }
1146//   {
1147//     gsl_complex c4 = gsl_complex_arccsc(c3);
1148//     print_complex(&c4, "arccsc");
1149//   }
1150//   {
1151//     gsl_complex c4 = gsl_complex_arccsc_real(5.);
1152//     print_complex(&c4, "arccsc_real");
1153//   }
1154//   {
1155//     gsl_complex c4 = gsl_complex_arccot(c3);
1156//     print_complex(&c4, "arccot");
1157//   }
1158//   {
1159//     gsl_complex c4 = gsl_complex_sinh(c3);
1160//     print_complex(&c4, "sinh");
1161//   }
1162//   {
1163//     gsl_complex c4 = gsl_complex_cosh(c3);
1164//     print_complex(&c4, "cosh");
1165//   }
1166//   {
1167//     gsl_complex c4 = gsl_complex_tanh(c3);
1168//     print_complex(&c4, "tanh");
1169//   }
1170//   {
1171//     gsl_complex c4 = gsl_complex_sech(c3);
1172//     print_complex(&c4, "sech");
1173//   }
1174//   {
1175//     gsl_complex c4 = gsl_complex_csch(c3);
1176//     print_complex(&c4, "csch");
1177//   }
1178//   {
1179//     gsl_complex c4 = gsl_complex_coth(c3);
1180//     print_complex(&c4, "coth");
1181//   }
1182//   {
1183//     gsl_complex c4 = gsl_complex_arcsinh(c3);
1184//     print_complex(&c4, "arcsinh");
1185//   }
1186//   {
1187//     gsl_complex c4 = gsl_complex_arccosh(c3);
1188//     print_complex(&c4, "arccosh");
1189//   }
1190//   {
1191//     gsl_complex c4 = gsl_complex_arccosh_real(5.);
1192//     print_complex(&c4, "arccosh_real");
1193//   }
1194//   {
1195//     gsl_complex c4 = gsl_complex_arctanh(c3);
1196//     print_complex(&c4, "arctanh");
1197//   }
1198//   {
1199//     gsl_complex c4 = gsl_complex_arctanh_real(5.);
1200//     print_complex(&c4, "arctanh_real");
1201//   }
1202//   {
1203//     gsl_complex c4 = gsl_complex_arcsech(c3);
1204//     print_complex(&c4, "arcsech");
1205//   }
1206//   {
1207//     gsl_complex c4 = gsl_complex_arccsch(c3);
1208//     print_complex(&c4, "arccsch");
1209//   }
1210//   {
1211//     gsl_complex c4 = gsl_complex_arccoth(c3);
1212//     print_complex(&c4, "arccoth");
1213//   }
1214//   return 0;
1215// }
1216// ```
1217#[test]
1218fn complex_f64() {
1219    let v = ComplexF64::rect(10., 10.);
1220    assert_eq!(v, ComplexF64 { dat: [10., 10.] });
1221    let v2 = ComplexF64::rect(1., -1.);
1222    assert_eq!(v2, ComplexF64 { dat: [1., -1.] });
1223    let v = ComplexF64::polar(5., 7.);
1224    assert_eq!(
1225        format!("{:.4} {:.4}", v.dat[0], v.dat[1]),
1226        "3.7695 3.2849".to_owned()
1227    );
1228
1229    let arg = v.arg();
1230    assert_eq!(format!("{:.4}", arg), "0.7168".to_owned());
1231    let arg = v.abs();
1232    assert_eq!(format!("{:.3}", arg), "5.000".to_owned());
1233    let arg = v.abs2();
1234    assert_eq!(format!("{:.3}", arg), "25.000".to_owned());
1235    let arg = v.logabs();
1236    assert_eq!(format!("{:.4}", arg), "1.6094".to_owned());
1237
1238    let v3 = v.add(&v2);
1239    assert_eq!(
1240        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1241        "4.7695 2.2849".to_owned()
1242    );
1243    let v3 = v.sub(&v2);
1244    assert_eq!(
1245        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1246        "2.7695 4.2849".to_owned()
1247    );
1248    let v3 = v.mul(&v2);
1249    assert_eq!(
1250        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1251        "7.0544 -0.4846".to_owned()
1252    );
1253    let v3 = v.div(&v2);
1254    assert_eq!(
1255        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1256        "0.2423 3.5272".to_owned()
1257    );
1258    let v3 = v.add_real(5.);
1259    assert_eq!(
1260        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1261        "8.7695 3.2849".to_owned()
1262    );
1263    let v3 = v.sub_real(5.);
1264    assert_eq!(
1265        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1266        "-1.2305 3.2849".to_owned()
1267    );
1268    let v3 = v.mul_real(5.);
1269    assert_eq!(
1270        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1271        "18.8476 16.4247".to_owned()
1272    );
1273    let v3 = v.div_real(5.);
1274    assert_eq!(
1275        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1276        "0.7539 0.6570".to_owned()
1277    );
1278
1279    let v3 = v.add_imag(5.);
1280    assert_eq!(
1281        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1282        "3.7695 8.2849".to_owned()
1283    );
1284    let v3 = v.sub_imag(5.);
1285    assert_eq!(
1286        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1287        "3.7695 -1.7151".to_owned()
1288    );
1289    let v3 = v.mul_imag(5.);
1290    assert_eq!(
1291        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1292        "-16.4247 18.8476".to_owned()
1293    );
1294    let v3 = v.div_imag(5.);
1295    assert_eq!(
1296        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1297        "0.6570 -0.7539".to_owned()
1298    );
1299
1300    let v3 = v.conjugate();
1301    assert_eq!(
1302        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1303        "3.7695 -3.2849".to_owned()
1304    );
1305    let v3 = v.inverse();
1306    assert_eq!(
1307        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1308        "0.1508 -0.1314".to_owned()
1309    );
1310    let v3 = v.negative();
1311    assert_eq!(
1312        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1313        "-3.7695 -3.2849".to_owned()
1314    );
1315    let v3 = v.sqrt();
1316    assert_eq!(
1317        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1318        "2.0940 0.7844".to_owned()
1319    );
1320    let v3 = ComplexF64::sqrt_real(5.);
1321    assert_eq!(
1322        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1323        "2.2361 0.0000".to_owned()
1324    );
1325
1326    let v3 = v.pow(&v2);
1327    assert_eq!(
1328        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1329        "6.4240 -7.9737".to_owned()
1330    );
1331    let v3 = v.pow_real(5.);
1332    assert_eq!(
1333        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1334        "-2824.0381 -1338.0708".to_owned()
1335    );
1336    let v3 = v.exp();
1337    assert_eq!(
1338        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1339        "-42.9142 -6.1938".to_owned()
1340    );
1341    let v3 = v.log();
1342    assert_eq!(
1343        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1344        "1.6094 0.7168".to_owned()
1345    );
1346    let v3 = v.log10();
1347    assert_eq!(
1348        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1349        "0.6990 0.3113".to_owned()
1350    );
1351    let v3 = v.log_b(&v2);
1352    assert_eq!(
1353        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1354        "-0.0071 2.0523".to_owned()
1355    );
1356    let v3 = v.sin();
1357    assert_eq!(
1358        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1359        "-7.8557 -10.7913".to_owned()
1360    );
1361    let v3 = v.cos();
1362    assert_eq!(
1363        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1364        "-10.8216 7.8337".to_owned()
1365    );
1366    let v3 = v.tan();
1367    assert_eq!(
1368        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1369        "0.0027 0.9991".to_owned()
1370    );
1371
1372    let v3 = v.sec();
1373    assert_eq!(
1374        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1375        "-0.0606 -0.0439".to_owned()
1376    );
1377    let v3 = v.csc();
1378    assert_eq!(
1379        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1380        "-0.0441 0.0606".to_owned()
1381    );
1382    let v3 = v.cot();
1383    assert_eq!(
1384        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1385        "0.0027 -1.0009".to_owned()
1386    );
1387    let v3 = v.arcsin();
1388    assert_eq!(
1389        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1390        "0.8440 2.3014".to_owned()
1391    );
1392    let v3 = ComplexF64::arcsin_real(5.);
1393    assert_eq!(
1394        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1395        "1.5708 -2.2924".to_owned()
1396    );
1397    let v3 = v.arccos();
1398    assert_eq!(
1399        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1400        "0.7268 -2.3014".to_owned()
1401    );
1402    let v3 = ComplexF64::arccos_real(5.);
1403    assert_eq!(
1404        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1405        "0.0000 2.2924".to_owned()
1406    );
1407    let v3 = v.arctan();
1408    assert_eq!(
1409        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1410        "1.4186 0.1291".to_owned()
1411    );
1412    let v3 = v.arcsec();
1413    assert_eq!(
1414        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1415        "1.4208 0.1325".to_owned()
1416    );
1417    let v3 = ComplexF64::arcsec_real(5.);
1418    assert_eq!(
1419        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1420        "1.3694 0.0000".to_owned()
1421    );
1422    let v3 = v.arccsc();
1423    assert_eq!(
1424        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1425        "0.1500 -0.1325".to_owned()
1426    );
1427    let v3 = ComplexF64::arccsc_real(5.);
1428    assert_eq!(
1429        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1430        "0.2014 0.0000".to_owned()
1431    );
1432    let v3 = v.arccot();
1433    assert_eq!(
1434        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1435        "0.1522 -0.1291".to_owned()
1436    );
1437    let v3 = v.sinh();
1438    assert_eq!(
1439        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1440        "-21.4457 -3.0986".to_owned()
1441    );
1442    let v3 = v.cosh();
1443    assert_eq!(
1444        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1445        "-21.4685 -3.0953".to_owned()
1446    );
1447    let v3 = v.tanh();
1448    assert_eq!(
1449        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1450        "0.9990 0.0003".to_owned()
1451    );
1452    let v3 = v.sech();
1453    assert_eq!(
1454        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1455        "-0.0456 0.0066".to_owned()
1456    );
1457    let v3 = v.csch();
1458    assert_eq!(
1459        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1460        "-0.0457 0.0066".to_owned()
1461    );
1462    let v3 = v.coth();
1463    assert_eq!(
1464        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1465        "1.0010 -0.0003".to_owned()
1466    );
1467    let v3 = v.arcsinh();
1468    assert_eq!(
1469        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1470        "2.3041 0.7070".to_owned()
1471    );
1472    let v3 = v.arccosh();
1473    assert_eq!(
1474        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1475        "2.3014 0.7268".to_owned()
1476    );
1477    let v3 = ComplexF64::arccosh_real(5.);
1478    assert_eq!(
1479        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1480        "2.2924 0.0000".to_owned()
1481    );
1482    let v3 = v.arctanh();
1483    assert_eq!(
1484        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1485        "0.1493 1.4372".to_owned()
1486    );
1487    let v3 = ComplexF64::arctanh_real(5.);
1488    assert_eq!(
1489        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1490        "0.2027 -1.5708".to_owned()
1491    );
1492    let v3 = v.arcsech();
1493    assert_eq!(
1494        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1495        "0.1325 -1.4208".to_owned()
1496    );
1497    let v3 = v.arccsch();
1498    assert_eq!(
1499        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1500        "0.1515 -0.1303".to_owned()
1501    );
1502    let v3 = v.arccoth();
1503    assert_eq!(
1504        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1505        "0.1493 -0.1336".to_owned()
1506    );
1507}
1508
1509#[test]
1510fn complex_f32() {
1511    let v = ComplexF32::rect(10., 10.);
1512    assert_eq!(v, ComplexF32 { dat: [10., 10.] });
1513    let v2 = ComplexF32::rect(1., -1.);
1514    assert_eq!(v2, ComplexF32 { dat: [1., -1.] });
1515    let v = ComplexF32::polar(5., 7.);
1516    assert_eq!(
1517        format!("{:.4} {:.4}", v.dat[0], v.dat[1]),
1518        "3.7695 3.2849".to_owned()
1519    );
1520
1521    let arg = v.arg();
1522    assert_eq!(format!("{:.4}", arg), "0.7168".to_owned());
1523    let arg = v.abs();
1524    assert_eq!(format!("{:.3}", arg), "5.000".to_owned());
1525    let arg = v.abs2();
1526    assert_eq!(format!("{:.3}", arg), "25.000".to_owned());
1527    let arg = v.logabs();
1528    assert_eq!(format!("{:.4}", arg), "1.6094".to_owned());
1529
1530    let v3 = v.add(&v2);
1531    assert_eq!(
1532        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1533        "4.7695 2.2849".to_owned()
1534    );
1535    let v3 = v.sub(&v2);
1536    assert_eq!(
1537        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1538        "2.7695 4.2849".to_owned()
1539    );
1540    let v3 = v.mul(&v2);
1541    assert_eq!(
1542        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1543        "7.0544 -0.4846".to_owned()
1544    );
1545    let v3 = v.div(&v2);
1546    assert_eq!(
1547        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1548        "0.2423 3.5272".to_owned()
1549    );
1550    let v3 = v.add_real(5.);
1551    assert_eq!(
1552        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1553        "8.7695 3.2849".to_owned()
1554    );
1555    let v3 = v.sub_real(5.);
1556    assert_eq!(
1557        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1558        "-1.2305 3.2849".to_owned()
1559    );
1560    let v3 = v.mul_real(5.);
1561    assert_eq!(
1562        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1563        "18.8476 16.4247".to_owned()
1564    );
1565    let v3 = v.div_real(5.);
1566    assert_eq!(
1567        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1568        "0.7539 0.6570".to_owned()
1569    );
1570
1571    let v3 = v.add_imag(5.);
1572    assert_eq!(
1573        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1574        "3.7695 8.2849".to_owned()
1575    );
1576    let v3 = v.sub_imag(5.);
1577    assert_eq!(
1578        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1579        "3.7695 -1.7151".to_owned()
1580    );
1581    let v3 = v.mul_imag(5.);
1582    assert_eq!(
1583        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1584        "-16.4247 18.8476".to_owned()
1585    );
1586    let v3 = v.div_imag(5.);
1587    assert_eq!(
1588        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1589        "0.6570 -0.7539".to_owned()
1590    );
1591
1592    let v3 = v.conjugate();
1593    assert_eq!(
1594        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1595        "3.7695 -3.2849".to_owned()
1596    );
1597    let v3 = v.inverse();
1598    assert_eq!(
1599        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1600        "0.1508 -0.1314".to_owned()
1601    );
1602    let v3 = v.negative();
1603    assert_eq!(
1604        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1605        "-3.7695 -3.2849".to_owned()
1606    );
1607    let v3 = v.sqrt();
1608    assert_eq!(
1609        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1610        "2.0940 0.7844".to_owned()
1611    );
1612    let v3 = ComplexF32::sqrt_real(5.);
1613    assert_eq!(
1614        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1615        "2.2361 0.0000".to_owned()
1616    );
1617
1618    let v3 = v.pow(&v2);
1619    assert_eq!(
1620        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1621        "6.4240 -7.9737".to_owned()
1622    );
1623    let v3 = v.pow_real(5.);
1624    assert_eq!(
1625        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1626        "-2824.0381 -1338.0712".to_owned()
1627    );
1628    let v3 = v.exp();
1629    assert_eq!(
1630        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1631        "-42.9142 -6.1938".to_owned()
1632    );
1633    let v3 = v.log();
1634    assert_eq!(
1635        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1636        "1.6094 0.7168".to_owned()
1637    );
1638    let v3 = v.log10();
1639    assert_eq!(
1640        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1641        "0.6990 0.3113".to_owned()
1642    );
1643    let v3 = v.log_b(&v2);
1644    assert_eq!(
1645        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1646        "-0.0071 2.0523".to_owned()
1647    );
1648    let v3 = v.sin();
1649    assert_eq!(
1650        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1651        "-7.8557 -10.7913".to_owned()
1652    );
1653    let v3 = v.cos();
1654    assert_eq!(
1655        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1656        "-10.8216 7.8337".to_owned()
1657    );
1658    let v3 = v.tan();
1659    assert_eq!(
1660        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1661        "0.0027 0.9991".to_owned()
1662    );
1663
1664    let v3 = v.sec();
1665    assert_eq!(
1666        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1667        "-0.0606 -0.0439".to_owned()
1668    );
1669    let v3 = v.csc();
1670    assert_eq!(
1671        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1672        "-0.0441 0.0606".to_owned()
1673    );
1674    let v3 = v.cot();
1675    assert_eq!(
1676        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1677        "0.0027 -1.0009".to_owned()
1678    );
1679    let v3 = v.arcsin();
1680    assert_eq!(
1681        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1682        "0.8440 2.3014".to_owned()
1683    );
1684    let v3 = ComplexF32::arcsin_real(5.);
1685    assert_eq!(
1686        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1687        "1.5708 -2.2924".to_owned()
1688    );
1689    let v3 = v.arccos();
1690    assert_eq!(
1691        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1692        "0.7268 -2.3014".to_owned()
1693    );
1694    let v3 = ComplexF32::arccos_real(5.);
1695    assert_eq!(
1696        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1697        "0.0000 2.2924".to_owned()
1698    );
1699    let v3 = v.arctan();
1700    assert_eq!(
1701        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1702        "1.4186 0.1291".to_owned()
1703    );
1704    let v3 = v.arcsec();
1705    assert_eq!(
1706        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1707        "1.4208 0.1325".to_owned()
1708    );
1709    let v3 = ComplexF32::arcsec_real(5.);
1710    assert_eq!(
1711        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1712        "1.3694 0.0000".to_owned()
1713    );
1714    let v3 = v.arccsc();
1715    assert_eq!(
1716        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1717        "0.1500 -0.1325".to_owned()
1718    );
1719    let v3 = ComplexF32::arccsc_real(5.);
1720    assert_eq!(
1721        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1722        "0.2014 0.0000".to_owned()
1723    );
1724    let v3 = v.arccot();
1725    assert_eq!(
1726        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1727        "0.1522 -0.1291".to_owned()
1728    );
1729    let v3 = v.sinh();
1730    assert_eq!(
1731        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1732        "-21.4457 -3.0986".to_owned()
1733    );
1734    let v3 = v.cosh();
1735    assert_eq!(
1736        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1737        "-21.4685 -3.0953".to_owned()
1738    );
1739    let v3 = v.tanh();
1740    assert_eq!(
1741        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1742        "0.9990 0.0003".to_owned()
1743    );
1744    let v3 = v.sech();
1745    assert_eq!(
1746        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1747        "-0.0456 0.0066".to_owned()
1748    );
1749    let v3 = v.csch();
1750    assert_eq!(
1751        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1752        "-0.0457 0.0066".to_owned()
1753    );
1754    let v3 = v.coth();
1755    assert_eq!(
1756        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1757        "1.0010 -0.0003".to_owned()
1758    );
1759    let v3 = v.arcsinh();
1760    assert_eq!(
1761        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1762        "2.3041 0.7070".to_owned()
1763    );
1764    let v3 = v.arccosh();
1765    assert_eq!(
1766        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1767        "2.3014 0.7268".to_owned()
1768    );
1769    let v3 = ComplexF32::arccosh_real(5.);
1770    assert_eq!(
1771        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1772        "2.2924 0.0000".to_owned()
1773    );
1774    let v3 = v.arctanh();
1775    assert_eq!(
1776        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1777        "0.1493 1.4372".to_owned()
1778    );
1779    let v3 = ComplexF32::arctanh_real(5.);
1780    assert_eq!(
1781        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1782        "0.2027 -1.5708".to_owned()
1783    );
1784    let v3 = v.arcsech();
1785    assert_eq!(
1786        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1787        "0.1325 -1.4208".to_owned()
1788    );
1789    let v3 = v.arccsch();
1790    assert_eq!(
1791        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1792        "0.1515 -0.1303".to_owned()
1793    );
1794    let v3 = v.arccoth();
1795    assert_eq!(
1796        format!("{:.4} {:.4}", v3.dat[0], v3.dat[1]),
1797        "0.1493 -0.1336".to_owned()
1798    );
1799}