File Coverage

complex.pd
Criterion Covered Total %
statement 104 110 94.5
branch 40 58 68.9
condition 6 11 54.5
subroutine 33 35 94.2
pod 11 19 57.8
total 194 233 83.2


line stmt bran cond sub pod time code
1             use strict;
2             use warnings;
3             use PDL::Types qw(ppdefs ppdefs_complex types);
4             my $R = [ppdefs()];
5             my $F = [map $_->ppsym, grep $_->real && !$_->integer, types()];
6             my $C = [ppdefs_complex()];
7              
8             pp_core_importList('()');
9             pp_beginwrap; # required for overload to work
10              
11             # pp_def functions go into the PDL::Complex namespace
12             # to avoid clashing with PDL::FFTW funcs of the same name that go
13             # into the PDL namespace
14             # it should be of no effect to the user of the module but you
15             # never know....
16             pp_bless('PDL::Complex');
17              
18             pp_addpm {At => 'Top'}, <<'EOD';
19 1     1   519 use strict;
  1         1  
  1         38  
20 1     1   5 use warnings;
  1         1  
  1         45  
21 1     1   5 use Carp;
  1         2  
  1         200  
22             our $VERSION = '2.011';
23              
24             =encoding utf8
25              
26             =head1 NAME
27              
28             PDL::Complex - handle complex numbers (DEPRECATED - use native complex)
29              
30             =head1 SYNOPSIS
31              
32             use PDL;
33             use PDL::Complex;
34              
35             =head1 DESCRIPTION
36              
37             This module is deprecated in favour of using "native complex" data types, e.g.:
38              
39             use PDL;
40             my $complex_pdl = cdouble('[1+3i]');
41             print $complex_pdl * pdl('i'); # [-3+i]
42              
43             This module features a growing number of functions manipulating complex
44             numbers. These are usually represented as a pair C<[ real imag ]> or
45             C<[ magnitude phase ]>. If not explicitly mentioned, the functions can work
46             inplace (not yet implemented!!!) and require rectangular form.
47              
48             While there is a procedural interface available (C<< $x/$y*$c <=> Cmul
49             (Cdiv ($x, $y), $c) >>), you can also opt to cast your pdl's into the
50             C datatype, which works just like your normal ndarrays, but
51             with all the normal perl operators overloaded.
52              
53             The latter means that C will be evaluated using the
54             normal rules of complex numbers, while other pdl functions (like C)
55             just treat the ndarray as a real-valued ndarray with a lowest dimension of
56             size 2, so C will return the maximum of all real and imaginary parts,
57             not the "highest" (for some definition)
58              
59             =head2 Native complex support
60              
61             2.027 added changes in complex number handling, with support for C99
62             complex floating-point types, and most functions and modules in the core
63             distribution support these as well.
64              
65             PDL can now handle complex numbers natively as scalars. This has
66             the advantage that real and complex valued ndarrays have the same
67             dimensions. Consider this when writing code in the future.
68              
69             See L, L, L, L,
70             L for more.
71              
72             =head1 TIPS, TRICKS & CAVEATS
73              
74             =over 4
75              
76             =item *
77              
78             C is a function (not, as of 2.047, a constant) exported by this module,
79             which represents C<-1**0.5>, i.e. the imaginary unit. it can be used to
80             quickly and conveniently write complex constants like this: C<4+3*i>.
81              
82             B This will override the PDL::Core function of the same name, which
83             returns a native complex value.
84              
85             =item *
86              
87             Use C to convert from real to complex, as in C<$r
88             = Cpow $cplx, r2C 2>. The overloaded operators automatically do that for
89             you, all the other functions, do not. So C will return all
90             the fifths roots of 1+1*i (due to broadcasting).
91              
92             =item *
93              
94             use C to cast from normal ndarrays into the
95             complex datatype. Use C to cast back. This
96             requires a copy, though.
97              
98             =back
99              
100             =head1 EXAMPLE WALK-THROUGH
101              
102             The complex constant five is equal to C:
103              
104             pdl> p $x = r2C 5
105             5 +0i
106              
107             Now calculate the three cubic roots of five:
108              
109             pdl> p $r = Croots $x, 3
110             [1.70998 +0i -0.854988 +1.48088i -0.854988 -1.48088i]
111              
112             Check that these really are the roots:
113              
114             pdl> p $r ** 3
115             [5 +0i 5 -1.22465e-15i 5 -7.65714e-15i]
116              
117             Duh! Could be better. Now try by multiplying C<$r> three times with itself:
118              
119             pdl> p $r*$r*$r
120             [5 +0i 5 -4.72647e-15i 5 -7.53694e-15i]
121              
122             Well... maybe C (which is used by the C<**> operator) isn't as
123             bad as I thought. Now multiply by C and negate, then take the complex
124             conjugate, which is just a very expensive way of swapping real and
125             imaginary parts.
126              
127             pdl> p Cconj(-($r*i))
128             [0 +1.70998i 1.48088 -0.854988i -1.48088 -0.854988i]
129              
130             Now plot the magnitude of (part of) the complex sine. First generate the
131             coefficients:
132              
133             pdl> $sin = i * zeroes(50)->xlinvals(2,4) + zeroes(50)->xlinvals(0,7)
134              
135             Now plot the imaginary part, the real part and the magnitude of the sine
136             into the same diagram:
137              
138             pdl> use PDL::Graphics::Gnuplot
139             pdl> gplot( with => 'lines',
140             PDL::cat(im ( sin $sin ),
141             re ( sin $sin ),
142             abs( sin $sin ) ))
143              
144             An ASCII version of this plot looks like this:
145              
146             30 ++-----+------+------+------+------+------+------+------+------+-----++
147             + + + + + + + + + + +
148             | $$|
149             | $ |
150             25 ++ $$ ++
151             | *** |
152             | ** *** |
153             | $$* *|
154             20 ++ $** ++
155             | $$$* #|
156             | $$$ * # |
157             | $$ * # |
158             15 ++ $$$ * # ++
159             | $$$ ** # |
160             | $$$$ * # |
161             | $$$$ * # |
162             10 ++ $$$$$ * # ++
163             | $$$$$ * # |
164             | $$$$$$$ * # |
165             5 ++ $$$############ * # ++
166             |*****$$$### ### * # |
167             * #***** # * # |
168             | ### *** ### ** # |
169             0 ## *** # * # ++
170             | * # * # |
171             | *** # ** # |
172             | * # * # |
173             -5 ++ ** # * # ++
174             | *** ## ** # |
175             | * #* # |
176             | **** ***## # |
177             -10 ++ **** # # ++
178             | # # |
179             | ## ## |
180             + + + + + + + ### + ### + + +
181             -15 ++-----+------+------+------+------+------+-----###-----+------+-----++
182             0 5 10 15 20 25 30 35 40 45 50
183              
184              
185             =head1 OPERATORS
186              
187             The following operators are overloaded:
188              
189             =over 4
190              
191             =item +, += (addition)
192              
193             =item -, -= (subtraction)
194              
195             =item *, *= (multiplication; L)
196              
197             =item /, /= (division; L)
198              
199             =item **, **= (exponentiation; L)
200              
201             =item atan2 (4-quadrant arc tangent)
202              
203             =item sin (L)
204              
205             =item cos (L)
206              
207             =item exp (L)
208              
209             =item abs (L)
210              
211             =item log (L)
212              
213             =item sqrt (L)
214              
215             =item ++, -- (increment, decrement; they affect the real part of the complex number only)
216              
217             =item "" (stringification)
218              
219             =back
220              
221             Comparing complex numbers other than for equality is a fatal error.
222              
223             =cut
224              
225 1     1   9 my $i;
226             BEGIN { $i = bless PDL->pdl(0,1) }
227 1     1   174 {
  1         5  
  1         1490  
228 49 100   49 0 12737 no warnings 'redefine';
229             sub i { $i->copy + (@_ ? $_[0] : 0) };
230             }
231              
232             # sensible aliases from PDL::LinearAlgebra
233             *r2p = \&Cr2p;
234             *p2r = \&Cp2r;
235             *conj = \&Cconj;
236             *abs = \&Cabs;
237             *abs2 = \&Cabs2;
238             *arg = \&Carg;
239             *tan = \&Ctan;
240             *proj = \&Cproj;
241             *asin = \&Casin;
242             *acos = \&Cacos;
243             *atan = \&Catan;
244             *sinh = \&Csinh;
245             *cosh = \&Ccosh;
246             *tanh = \&Ctanh;
247             *asinh = \&Casinh;
248             *acosh = \&Cacosh;
249             *atanh = \&Catanh;
250             *tricpy = \&Ctricpy;
251             *mstack = \&Cmstack;
252             *augment = \&Caugment;
253             EOD
254              
255             for (qw(Ctan Catan re im i cplx real)) {
256             pp_add_exported '', $_;
257             }
258              
259             pp_addhdr <<'EOH';
260              
261             #include
262              
263             #ifndef M_PI
264             # define M_PI 3.1415926535897932384626433832795029
265             #endif
266             #ifndef M_2PI
267             # define M_2PI (2. * M_PI)
268             #endif
269              
270             #if __GLIBC__ > 1 && (defined __USE_MISC || defined __USE_XOPEN || defined __USE_ISOC9X)
271             # define CABS(r,i) hypot (r, i)
272             #else
273             static double
274             CABS (double r, double i)
275             {
276             double t;
277              
278             if (r < 0) r = - r;
279             if (i < 0) i = - i;
280              
281             if (i > r)
282             {
283             t = r; r = i; i = t;
284             }
285              
286             if (r + i == r)
287             return r;
288              
289             t = i / r;
290             return r * sqrt (1 + t*t);
291             }
292             #endif
293              
294             #if __GLIBC__ >= 2 && __GLIBC_MINOR__ >= 1 && defined __USE_GNU
295             # define SINCOS(x,s,c) sincos ((x), &(s), &(c))
296             #else
297             # define SINCOS(x,s,c) \
298             (s) = sin (x); \
299             (c) = cos (x);
300             #endif
301              
302              
303             #define CSQRT(type,ar,ai,cr,ci) \
304             type mag = CABS ((ar), (ai)); \
305             type t; \
306             \
307             if (mag == 0) \
308             (cr) = (ci) = 0; \
309             else if ((ar) > 0) \
310             { \
311             t = sqrt (0.5 * (mag + (ar))); \
312             (cr) = t; \
313             (ci) = 0.5 * (ai) / t; \
314             } \
315             else \
316             { \
317             t = sqrt (0.5 * (mag - (ar))); \
318             \
319             if ((ai) < 0) \
320             t = -t; \
321             \
322             (cr) = 0.5 * (ai) / t; \
323             (ci) = t; \
324             }
325              
326              
327             #define CLOG(ar,ai,cr,ci) \
328             (cr) = log (CABS ((ar), (ai))); \
329             (ci) = atan2 ((ai), (ar));
330              
331             EOH
332              
333             pp_addpm <<'EOP';
334              
335             =head2 from_native
336              
337             =for ref
338              
339             Class method to convert a native-complex ndarray to a PDL::Complex object.
340              
341             =for usage
342              
343             PDL::Complex->from_native($native_complex_ndarray)
344              
345             =cut
346              
347 8     8 1 18919 sub from_native {
348 8 100       37 my ($class, $ndarray) = @_;
349 5 50       15 return $ndarray if UNIVERSAL::isa($ndarray,'PDL::Complex'); # NOOP if P:C
350 5 50       84 croak "not an ndarray" if !UNIVERSAL::isa($ndarray,'PDL');
351 5         198 croak "not a native complex ndarray" if $ndarray->type->real;
352             bless PDL::append($ndarray->re->dummy(0),$ndarray->im->dummy(0)), $class;
353             }
354              
355             =head2 as_native
356              
357             =for ref
358              
359             Object method to convert a PDL::Complex object to a native-complex ndarray.
360              
361             =for usage
362              
363             $pdl_complex_obj->as_native
364              
365             =cut
366              
367 11     11 1 1893 sub as_native {
368             PDL::Ops::czip(map $_[0]->slice("($_)"), 0..1);
369             }
370              
371             =head2 cplx
372              
373             =for ref
374              
375             Cast a real-valued ndarray to the complex datatype.
376              
377             The first dimension of the ndarray must be of size 2. After this the
378             usual (complex) arithmetic operators are applied to this pdl, rather
379             than the normal elementwise pdl operators. Dataflow to the complex
380             parent works. Use C on the result if you don't want this.
381              
382             =for usage
383              
384             cplx($real_valued_pdl)
385              
386             =head2 complex
387              
388             =for ref
389              
390             Cast a real-valued ndarray to the complex datatype I dataflow
391             and I.
392              
393             Achieved by merely reblessing an ndarray. The first dimension of the
394             ndarray must be of size 2.
395              
396             =for usage
397              
398             complex($real_valued_pdl)
399              
400             =head2 real
401              
402             =for ref
403              
404             Cast a complex valued pdl back to the "normal" pdl datatype.
405              
406             Afterwards the normal elementwise pdl operators are used in
407             operations. Dataflow to the real parent works. Use C on the
408             result if you don't want this.
409              
410             =for usage
411              
412             real($cplx_valued_pdl)
413              
414             =cut
415              
416 10 50   10 1 7184 sub cplx($) {
417 10 50 33     84 return $_[0] if UNIVERSAL::isa($_[0],'PDL::Complex'); # NOOP if just ndarray
418 10         53 croak "first dimsize must be 2" unless $_[0]->dims > 0 && $_[0]->dim(0) == 2;
419             bless $_[0]->slice('');
420             }
421              
422 12 100   12 1 7977 sub complex($) {
423 10 50 33     87 return $_[0] if UNIVERSAL::isa($_[0],'PDL::Complex'); # NOOP if just ndarray
424 10         43 croak "first dimsize must be 2" unless $_[0]->dims > 0 && $_[0]->dim(0) == 2;
425             bless $_[0];
426             }
427              
428             *PDL::cplx = \&cplx;
429             *PDL::complex = \&complex;
430              
431 63 50   63 1 6181 sub real($) {
432 63         145 return $_[0] unless UNIVERSAL::isa($_[0],'PDL::Complex'); # NOOP unless complex
433             bless $_[0]->slice(''), 'PDL';
434             }
435              
436             =head2 t
437              
438             =for usage
439              
440             $pdl = $pdl->t(SCALAR(conj))
441             conj : Conjugate Transpose = 1 | Transpose = 0, default = 0;
442              
443             =for ref
444              
445             Convenient function for transposing real or complex 2D array(s).
446             For complex data, if conj is true returns conjugate transposed array(s).
447             Supports broadcasting. Not exported.
448              
449             Originally by Grégory Vanuxem.
450              
451             =cut
452              
453 4     4 1 344 sub t {
454 4         9 my ($m, $conj) = @_;
455 4 100       18 my $ndims = $m->dims;
    100          
456             my $r = $ndims > 2 ? $m->xchg(1,2) :
457             $ndims > 1 ? $m->dummy(1) :
458 4 100       47 $m->dummy(1)->dummy(1);
459             $conj ? $r->conj : $r;
460             }
461              
462             EOP
463              
464             pp_def 'r2C',
465             Pars => 'r(); [o]c(m=2)',
466             Doc => 'convert real to complex, assuming an imaginary part of zero',
467             PMCode => << 'EOPM',
468              
469             undef &PDL::r2C;
470             *PDL::r2C = \&PDL::Complex::r2C;
471             sub PDL::Complex::r2C {
472             return $_[0] if UNIVERSAL::isa($_[0],'PDL::Complex');
473             my $r = __PACKAGE__->initialize;
474             &PDL::Complex::_r2C_int($_[0], $r);
475             $r }
476              
477             EOPM
478             Code => q!
479             $c(m=>0) = $r();
480             $c(m=>1) = 0;
481             !
482             ;
483              
484             pp_def 'i2C',
485             Pars => 'r(); [o]c(m=2)',
486             Doc => 'convert imaginary to complex, assuming a real part of zero',
487             PMCode => 'undef &PDL::i2C; *PDL::i2C = \&PDL::Complex::i2C; sub PDL::Complex::i2C { my $r = __PACKAGE__->initialize; &PDL::Complex::_i2C_int($_[0], $r); $r }',
488             Code => q!
489             $c(m=>0) = 0;
490             $c(m=>1) = $r();
491             !
492             ;
493              
494             pp_def 'Cr2p',
495             Pars => 'r(m=2); [o]p(m=2)',
496             Inplace => 1,
497             GenericTypes => $F,
498             Doc => 'convert complex numbers in rectangular form to polar (mod,arg) form. Works inplace',
499             Code => q!
500             $GENERIC() x = $r(m=>0);
501             $GENERIC() y = $r(m=>1);
502             $p(m=>0) = CABS (x, y);
503             $p(m=>1) = atan2 (y, x);
504             !
505             ;
506              
507             pp_def 'Cp2r',
508             Pars => 'r(m=2); [o]p(m=2)',
509             Inplace => 1,
510             GenericTypes => $F,
511             Doc => 'convert complex numbers in polar (mod,arg) form to rectangular form. Works inplace',
512             Code => q!
513             $GENERIC() m = $r(m=>0);
514             $GENERIC() a = $r(m=>1);
515             double s, c;
516              
517             SINCOS (a, s, c);
518             $p(m=>0) = c * m;
519             $p(m=>1) = s * m;
520             !
521             ;
522              
523             pp_def 'Cadd', # this is here for a) completeness and b) not having to mess with PDL::Ops
524             Pars => 'a(m=2); b(m=2); [o]c(m=2)',
525             Doc => undef,
526             Code => q^
527             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
528             $GENERIC() br = $b(m=>0), bi = $b(m=>1);
529             $c(m=>0) = ar + br;
530             $c(m=>1) = ai + bi;
531             ^
532             ;
533              
534             pp_def 'Csub', # this is here for a) completeness and b) not having to mess with PDL::Ops
535             Pars => 'a(m=2); b(m=2); [o]c(m=2)',
536             Doc => undef,
537             Code => q^
538             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
539             $GENERIC() br = $b(m=>0), bi = $b(m=>1);
540             $c(m=>0) = ar - br;
541             $c(m=>1) = ai - bi;
542             ^
543             ;
544              
545             pp_def 'Cmul',
546             Pars => 'a(m=2); b(m=2); [o]c(m=2)',
547             Doc => 'complex multiplication',
548             Code => q^
549             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
550             $GENERIC() br = $b(m=>0), bi = $b(m=>1);
551             $c(m=>0) = ar*br - ai*bi;
552             $c(m=>1) = ar*bi + ai*br;
553             ^
554             ;
555              
556             pp_def 'Cprodover',
557             Pars => 'a(m=2,n); [o]c(m=2)',
558             Doc => 'Project via product to N-1 dimension',
559             Code => q^
560             $GENERIC() br, bi, cr, ci,tmp;
561             cr = $a(m=>0,n=>0);
562             ci = $a(m=>1,n=>0);
563             loop(n=1) %{
564             br = $a(m=>0);
565             bi = $a(m=>1);
566             tmp = cr*bi + ci*br;
567             cr = cr*br - ci*bi;
568             ci = tmp;
569             %}
570             $c(m=>0) = cr;
571             $c(m=>1) = ci;
572             ^
573             ;
574              
575             pp_def 'Cscale',
576             Pars => 'a(m=2); b(); [o]c(m=2)',
577             Doc => 'mixed complex/real multiplication',
578             Code => q^
579             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
580             $c(m=>0) = ar * $b();
581             $c(m=>1) = ai * $b();
582             ^
583             ;
584              
585             pp_def 'Cdiv',
586             Pars => 'a(m=2); b(m=2); [o]c(m=2)',
587             GenericTypes => $F,
588             Doc => 'complex division',
589             Code => q^
590             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
591             $GENERIC() br = $b(m=>0), bi = $b(m=>1);
592              
593             if (fabsl (br) > fabsl (bi))
594             {
595             $GENERIC() tt = bi / br;
596             $GENERIC() dn = br + tt * bi;
597             $c(m=>0) = (ar + tt * ai) / dn;
598             $c(m=>1) = (ai - tt * ar) / dn;
599             }
600             else
601             {
602             $GENERIC() tt = br / bi;
603             $GENERIC() dn = br * tt + bi;
604             $c(m=>0) = (ar * tt + ai) / dn;
605             $c(m=>1) = (ai * tt - ar) / dn;
606             }
607             ^
608             ;
609              
610             pp_def 'Ceq',
611             Pars => 'a(m=2); b(m=2); [o]c()',
612             GenericTypes => $F,
613             Doc => "=for ref\n\nComplex equality operator.",
614             Code => q^
615             $c() = (($a(m=>0) == $b(m=>0)) && ($a(m=>1) == $b(m=>1)));
616             ^,
617             PMCode => <<'EOF',
618             sub PDL::Complex::Ceq {
619             my @args = !$_[2] ? @_[1,0] : @_[0,1];
620             $args[1] = r2C($args[1]) if ref $args[1] ne __PACKAGE__;
621             PDL::Complex::_Ceq_int($args[0], $args[1], my $r = PDL->null);
622             $r;
623             }
624             EOF
625             ;
626              
627             pp_def 'Cconj',
628             Pars => 'a(m=2); [o]c(m=2)',
629             Inplace => 1,
630             Doc => 'complex conjugation. Works inplace',
631             Code => q^
632             $c(m=>0) = $a(m=>0);
633             $c(m=>1) = -$a(m=>1);
634             ^
635             ;
636              
637             pp_def 'Cabs',
638             Pars => 'a(m=2); [o]c()',
639             GenericTypes => $F,
640             Doc => 'complex C (also known as I)',
641             PMCode => q^sub PDL::Complex::Cabs($) {
642             my $pdl= shift;
643             my $abs = PDL->null;
644             &PDL::Complex::_Cabs_int($pdl, $abs);
645             $abs;
646             }^,
647             Code => q^
648             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
649             $c() = CABS (ar, ai);
650             ^
651             ;
652              
653             pp_def 'Cabs2',
654             Pars => 'a(m=2); [o]c()',
655             Doc => 'complex squared C (also known I)',
656             PMCode => q^sub PDL::Complex::Cabs2($) {
657             my $pdl= shift;
658             my $abs2 = PDL->null;
659             &PDL::Complex::_Cabs2_int($pdl, $abs2);
660             $abs2;
661             }^,
662             Code => q^
663             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
664             $c() = ar*ar + ai*ai;
665             ^
666             ;
667              
668             pp_def 'Carg',
669             Pars => 'a(m=2); [o]c()',
670             GenericTypes => $F,
671             Doc => 'complex argument function ("angle")',
672             PMCode => q^sub PDL::Complex::Carg($) {
673             my $pdl= shift;
674             my $arg = PDL->null;
675             &PDL::Complex::_Carg_int($pdl, $arg);
676             $arg;
677             }^,
678             Code => q^
679             $c() = atan2 ($a(m=>1), $a(m=>0));
680             ^
681             ;
682              
683             pp_def 'Csin',
684             Pars => 'a(m=2); [o]c(m=2)',
685             Inplace => 1,
686             GenericTypes => $F,
687             Doc => ' sin (a) = 1/(2*i) * (exp (a*i) - exp (-a*i)). Works inplace',
688             Code => q^
689             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
690             double s, c;
691              
692             SINCOS (ar, s, c);
693             $c(m=>0) = s * cosh (ai);
694             $c(m=>1) = c * sinh (ai);
695             ^
696             ;
697              
698             pp_def 'Ccos',
699             Pars => 'a(m=2); [o]c(m=2)',
700             Inplace => 1,
701             GenericTypes => $F,
702             Doc => ' cos (a) = 1/2 * (exp (a*i) + exp (-a*i)). Works inplace',
703             Code => q^
704             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
705             double s, c;
706              
707             SINCOS (ar, s, c);
708             $c(m=>0) = c * cosh (ai);
709             $c(m=>1) = - s * sinh (ai);
710             ^
711             ;
712              
713             pp_addpm <<'EOD';
714              
715             =head2 Ctan
716              
717             =for ref
718              
719             Complex tangent
720              
721             tan (a) = -i * (exp (a*i) - exp (-a*i)) / (exp (a*i) + exp (-a*i))
722              
723             Does not work inplace.
724              
725             =cut
726 5     5 1 125  
727             sub Ctan($) { Csin($_[0]) / Ccos($_[0]) }
728              
729              
730             EOD
731              
732             pp_def 'Cexp',
733             Pars => 'a(m=2); [o]c(m=2)',
734             Inplace => 1,
735             GenericTypes => $F,
736             Doc => ' exp (a) = exp (real (a)) * (cos (imag (a)) + i * sin (imag (a))). Works inplace',
737             Code => q^
738             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
739             $GENERIC() ex = exp (ar);
740             double s, c;
741              
742             SINCOS (ai, s, c);
743             $c(m=>0) = ex * c;
744             $c(m=>1) = ex * s;
745             ^
746             ;
747              
748             pp_def 'Clog',
749             Pars => 'a(m=2); [o]c(m=2)',
750             Inplace => 1,
751             GenericTypes => $F,
752             Doc => ' log (a) = log (cabs (a)) + i * carg (a). Works inplace',
753             Code => q^
754             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
755              
756             CLOG (ar, ai, $c(m=>0), $c(m=>1));
757             ^
758             ;
759              
760             pp_def 'Cpow',
761             Pars => 'a(m=2); b(m=2); [o]c(m=2)',
762             Inplace => ['a'],
763             GenericTypes => $F,
764             Doc => 'complex C (C<**>-operator)',
765             Code => q^
766             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
767             $GENERIC() br = $b(m=>0), bi = $b(m=>1);
768              
769             /* real ndarray (scalar or 1-ndarray) */
770             if($PDL(b)->dims[0]==0)
771             bi = 0;
772             double logr, logi, x, y;
773             double s, c;
774              
775             if(ar == 0 && ai == 0){
776             if(br == 0 && bi == 0) {
777             $c(m=>0) = 1;
778             $c(m=>1) = 0;
779             }
780             else {
781             $c(m=>0) = 0;
782             $c(m=>1) = 0;
783             }
784             }
785             else {
786             CLOG (ar, ai, logr, logi);
787             x = exp (logr*br - logi*bi);
788             y = logr*bi + logi*br;
789              
790             SINCOS (y, s, c);
791              
792             $c(m=>0) = x * c;
793             if(ai == 0 && bi == 0) $c(m=>1) = 0;
794             else $c(m=>1) = x * s;
795             }
796             ^
797             ;
798              
799             pp_def 'Csqrt',
800             Pars => 'a(m=2); [o]c(m=2)',
801             Inplace => 1,
802             GenericTypes => $F,
803             Doc => 'Works inplace',
804             Code => q^
805             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
806              
807             CSQRT ($GENERIC(), ar, ai, $c(m=>0), $c(m=>1));
808             ^
809             ;
810              
811             pp_def 'Casin',
812             Pars => 'a(m=2); [o]c(m=2)',
813             Inplace => 1,
814             GenericTypes => $F,
815             Doc => 'Works inplace',
816             Code => q^
817             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
818              
819             $GENERIC() t1 = sqrt ((ar+1)*(ar+1) + ai*ai);
820             $GENERIC() t2 = sqrt ((ar-1)*(ar-1) + ai*ai);
821             $GENERIC() alpha = (t1+t2)*0.5;
822             $GENERIC() beta = (t1-t2)*0.5;
823              
824             if (alpha < 1) alpha = 1;
825             if (beta > 1) beta = 1;
826             else if (beta < -1) beta = -1;
827              
828             $c(m=>0) = atan2 (beta, sqrt (1-beta*beta));
829             $c(m=>1) = - log (alpha + sqrt (alpha*alpha-1));
830             if (ai > 0 || (ai == 0 && ar < -1))
831             $c(m=>1) = - $c(m=>1);
832             ^
833             ;
834              
835             pp_def 'Cacos',
836             Pars => 'a(m=2); [o]c(m=2)',
837             Inplace => 1,
838             GenericTypes => $F,
839             Doc => 'Works inplace',
840             Code => q^
841             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
842              
843             $GENERIC() t1 = sqrt ((ar+1)*(ar+1) + ai*ai);
844             $GENERIC() t2 = sqrt ((ar-1)*(ar-1) + ai*ai);
845             $GENERIC() alpha = (t1+t2)*0.5;
846             $GENERIC() beta = (t1-t2)*0.5;
847              
848             if (alpha < 1) alpha = 1;
849             if (beta > 1) beta = 1;
850             else if (beta < -1) beta = -1;
851              
852             $c(m=>0) = atan2 (sqrt (1-beta*beta), beta);
853             $c(m=>1) = log (alpha + sqrt (alpha*alpha-1));
854             if (ai > 0 || (ai == 0 && ar < -1))
855             $c(m=>1) = - $c(m=>1);
856             ^
857             ;
858              
859             pp_addpm <<'EOD';
860              
861             =head2 Catan
862              
863             =for ref
864              
865             Return the complex C.
866              
867             Does not work inplace.
868              
869             =cut
870              
871 2     2 1 4 sub Catan($) {
872 2         6 my $z = shift;
873             Cmul Clog(Cdiv (PDL::Complex::i()+$z, PDL::Complex::i()-$z)), PDL->pdl(0, 0.5);
874             }
875              
876             EOD
877              
878             pp_def 'Csinh',
879             Pars => 'a(m=2); [o]c(m=2)',
880             Inplace => 1,
881             GenericTypes => $F,
882             Doc => ' sinh (a) = (exp (a) - exp (-a)) / 2. Works inplace',
883             Code => q^
884             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
885             double s, c;
886              
887             SINCOS (ai, s, c);
888             $c(m=>0) = sinh (ar) * c;
889             $c(m=>1) = cosh (ar) * s;
890             ^
891             ;
892              
893             pp_def 'Ccosh',
894             Pars => 'a(m=2); [o]c(m=2)',
895             Inplace => 1,
896             GenericTypes => $F,
897             Doc => ' cosh (a) = (exp (a) + exp (-a)) / 2. Works inplace',
898             Code => q^
899             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
900             double s, c;
901              
902             SINCOS (ai, s, c);
903             $c(m=>0) = cosh (ar) * c;
904             $c(m=>1) = sinh (ar) * s;
905             ^
906             ;
907              
908             pp_def 'Ctanh',
909             Pars => 'a(m=2); [o]c(m=2)',
910             Inplace => 1,
911             GenericTypes => $F,
912             Doc => 'Works inplace',
913             Code => q^
914             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
915             double den;
916             double s, c;
917              
918             SINCOS (2*ai, s, c);
919             den = cosh (2*ar) + c;
920              
921             $c(m=>0) = sinh (2*ar) / den;
922             $c(m=>1) = s / den;
923             ^
924             ;
925              
926             pp_def 'Casinh',
927             Pars => 'a(m=2); [o]c(m=2)',
928             Inplace => 1,
929             GenericTypes => $F,
930             Doc => 'Works inplace',
931             Code => q^
932             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
933             $GENERIC() yr = (ar-ai) * (ar+ai) + 1;
934             $GENERIC() yi = 2*ar*ai;
935              
936             CSQRT ($GENERIC(), yr, yi, yr, yi)
937              
938             yr += ar;
939             yi += ai;
940              
941             CLOG (yr, yi, $c(m=>0), $c(m=>1));
942             ^
943             ;
944              
945             pp_def 'Cacosh',
946             Pars => 'a(m=2); [o]c(m=2)',
947             Inplace => 1,
948             GenericTypes => $F,
949             Doc => 'Works inplace',
950             Code => q^
951             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
952              
953             $GENERIC() yr = (ar-ai) * (ar+ai) - 1;
954             $GENERIC() yi = 2*ar*ai;
955              
956             CSQRT ($GENERIC(), yr, yi, yr, yi)
957              
958             yr += ar;
959             yi += ai;
960              
961             CLOG (yr, yi, $c(m=>0), $c(m=>1));
962             ^
963             ;
964              
965             pp_def 'Catanh',
966             Pars => 'a(m=2); [o]c(m=2)',
967             Inplace => 1,
968             GenericTypes => $F,
969             Doc => 'Works inplace',
970             Code => q^
971             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
972              
973             double i2 = ai*ai;
974             double num = i2 + (1+ar) * (1+ar);
975             double den = i2 + (1-ar) * (1-ar);
976              
977             $c(m=>0) = 0.25 * (log(num) - log(den));
978             $c(m=>1) = 0.5 * atan2 (2*ai, 1 - ar*ar - i2);
979             ^
980             ;
981              
982             pp_def 'Cproj',
983             Pars => 'a(m=2); [o]c(m=2)',
984             Inplace => 1,
985             GenericTypes => $F,
986             Doc => 'compute the projection of a complex number to the riemann sphere. Works inplace',
987             Code => q^
988             $GENERIC() ar = $a(m=>0), ai = $a(m=>1);
989              
990             double den = ar*ar + ai*ai + 1;
991              
992             $c(m=>0) = 2*ar / den;
993             $c(m=>1) = 2*ai / den;
994             ^
995             ;
996              
997             pp_def 'Croots',
998             Pars => 'a(m=2); [o]c(m=2,n)',
999             OtherPars => 'int n => n',
1000             GenericTypes => $F,
1001             Doc => 'Compute the C roots of C. C must be a positive integer. The result will always be a complex type!',
1002             PMCode => q^sub PDL::Complex::Croots($$) {
1003             my ($pdl, $n) = @_;
1004             my $r = PDL->null;
1005             &PDL::Complex::_Croots_int($pdl, $r, $n);
1006             bless $r;
1007             }^,
1008             Code => q^
1009             double s, c;
1010             double ar = $a(m=>0), ai = $a(m=>1),
1011             n1 = 1 / (double)$COMP(n),
1012             rr = pow (CABS (ar, ai), n1), /* do not optimize the sqrt out of this expr! */
1013             at = atan2 (ai, ar) * n1,
1014             ti = M_2PI * n1;
1015              
1016             loop(n) %{
1017             SINCOS (at, s, c);
1018              
1019             $c(m=>0) = rr * c;
1020             $c(m=>1) = rr * s;
1021              
1022             at += ti;
1023             %}
1024             ^
1025             ;
1026              
1027             pp_addpm <<'EOD';
1028              
1029             =head2 re, im
1030              
1031             Return the real or imaginary part of the complex number(s) given.
1032              
1033             These are slicing operators, so data flow works. The real and
1034             imaginary parts are returned as ndarrays (ref eq PDL).
1035              
1036             =cut
1037 21     21 1 207  
1038 21     21 1 574 sub re($) { $_[0]->slice("(0)") }
1039             sub im($) { $_[0]->slice("(1)") }
1040              
1041 1     1   7 {
  1         2  
  1         350  
1042             no warnings 'redefine';
1043             # if the argument does anything other than pass through 0-th dim, re-bless
1044 137 50   137 0 645 sub slice :lvalue {
1045 137 100 100     821 my $first = ref $_[1] ? $_[1][0] : (split ',', $_[1])[0];
1046 137         516 my $class = ($first//'') =~ /^[:x]?$/i ? ref($_[0]) : 'PDL';
1047 137         1432 my $ret = bless $_[0]->SUPER::slice(@_[1..$#_]), $class;
1048             $ret;
1049             }
1050             }
1051              
1052             EOD
1053              
1054             pp_def 'rCpolynomial',
1055             Pars => 'coeffs(n); x(c=2,m); [o]out(c=2,m)',
1056             Doc => 'evaluate the polynomial with (real) coefficients C at the (complex) position(s) C. C is the constant term.',
1057             GenericTypes => $F,
1058             PMCode=> q!
1059             sub rCpolynomial {
1060             my $coeffs = shift;
1061             my $x = shift;
1062             my $out = $x->copy;
1063             _rCpolynomial_int($coeffs,$x,$out);
1064             return PDL::complex($out);
1065             }
1066             !,
1067             Code => q!
1068             loop(m) %{
1069             double xr = 1;
1070             double xi = 0;
1071             double or = 0;
1072             double oi = 0;
1073             double Xr;
1074             loop(n) %{
1075             or += $coeffs() * xr;
1076             oi += $coeffs() * xi;
1077             Xr = xr;
1078             xr = Xr * $x(c=>0) - xi * $x(c=>1);
1079             xi = xi * $x(c=>0) + Xr * $x(c=>1);
1080             %}
1081             $out(c=>0) = or;
1082             $out(c=>1) = oi;
1083             %}
1084             !
1085             ;
1086              
1087             pp_def('Ctricpy',
1088             Pars => 'A(c=2,m,n);[o] C(c=2,m,n)',
1089             OtherPars => 'int uplo',
1090             OtherParsDefaults => {uplo => 0},
1091             ArgOrder => [qw(A uplo C)],
1092             Code => '
1093             if ($COMP(uplo)) { /* lower */
1094             broadcastloop %{ loop(n,m=:n+1,c) %{ $C() = $A(); %} %}
1095             } else {
1096             broadcastloop %{ loop(n,m=n,c) %{ $C() = $A(); %} %}
1097             }
1098             ',
1099             Doc => <<'EOT',
1100             =for usage
1101              
1102             tricpy(PDL(A), int(uplo), PDL(C))
1103              
1104             =for example
1105              
1106             $c = $a->tricpy($uplo); # explicit uplo
1107             $c = $a->tricpy; # default upper
1108              
1109             or
1110              
1111             tricpy($a, $uplo, $c); # modify c
1112              
1113             =for ref
1114              
1115             Copy triangular part to another matrix. If uplo == 0 copy upper triangular
1116             part.
1117              
1118             Originally by Grégory Vanuxem.
1119             EOT
1120             );
1121              
1122             pp_def('Cmstack',
1123             DefaultFlow => 1,
1124             TwoWay => 1,
1125             Pars => 'x(c=2,n,m);y(c,n,p);[o]out(c,n,q=CALC($SIZE(m)+$SIZE(p)));',
1126             Code => '
1127             loop (m,n,c) %{ $out(q=>m) = $x(); %}
1128             loop (q=$SIZE(m),n,c) %{ $out() = $y(p=>q-$SIZE(m)); %}
1129             ',
1130             BackCode => '
1131             loop (m,n,c) %{ $x() = $out(q=>m); %}
1132             loop (q=$SIZE(m),n,c) %{ $y(p=>q-$SIZE(m)) = $out(); %}
1133             ',
1134             Doc => <
1135             =for ref
1136              
1137             Combine two 2D ndarrays into a single ndarray, along the second
1138             ("vertical") dim.
1139             This routine does backward and forward dataflow automatically.
1140              
1141             Originally by Grégory Vanuxem.
1142             EOT
1143             );
1144              
1145             pp_def('Caugment',
1146             DefaultFlow => 1,
1147             TwoWay => 1,
1148             Pars => 'x(c=2,n);y(c,p);[o]out(c,q=CALC($SIZE(n)+$SIZE(p)));',
1149             Code => '
1150             loop (q=:$SIZE(n),c) %{ $out() = $x(n=>q); %}
1151             loop (q=$SIZE(n),c) %{ $out() = $y(p=>q-$SIZE(n)); %}
1152             ',
1153             BackCode => '
1154             loop (q=:$SIZE(n),c) %{ $x(n=>q) = $out(); %}
1155             loop (q=$SIZE(n),c) %{ $y(p=>q-$SIZE(n)) = $out(); %}
1156             ',
1157             Doc => <
1158             =for ref
1159              
1160             Combine two ndarrays into a single ndarray along the 0-th ("horizontal") dim.
1161             This routine does backward and forward dataflow automatically.
1162              
1163             Originally by Grégory Vanuxem.
1164             EOT
1165             );
1166              
1167             pp_add_isa 'PDL';
1168              
1169             pp_addpm {At => 'Bot'}, <<'EOD';
1170              
1171 1     1 0 3 # undocumented compatibility functions (thanks to Luis Mochan!)
1172 0     0 1 0 sub Catan2 { Clog( $_[1] + i()*$_[0])/i }
1173             sub atan2 { Clog( $_[1] + i()*$_[0])/i }
1174              
1175             =begin comment
1176              
1177             In _gen_biop, the '+' or '-' between the operator (e.g., '*') and the
1178             function that it overloads (e.g., 'Cmul') flags whether the operation
1179             is ('+') or is not ('-') commutative. See the discussion of argument
1180             swapping in the section "Calling Conventions and Magic Autogeneration"
1181             in "perldoc overload".
1182              
1183             =end comment
1184              
1185             =cut
1186 1     1   321  
1187             my %NO_MUTATE; BEGIN { @NO_MUTATE{qw(atan2 .= ==)} = (); }
1188 8     8   15 sub _gen_biop {
1189 8         11 local $_ = shift;
1190 8 50       56 my $sub;
1191 8 100   8   1053 die if !(my ($op, $commutes, $func) = /(\S+)([-+])(\w+)/);
  8 50       677  
  8 50       39  
  8 50       19  
  35 50       934  
  35 100       164  
  35         140  
  79         2878  
  79         362  
  79         304  
  10         322  
  10         39  
  10         43  
  5         1360  
  5         31  
  5         30  
  3         101  
  3         16  
  3         13  
  1         6  
  1         3  
  1         3  
  13         662  
  13         57  
  13         58  
1192             $sub = eval 'sub {
1193             my ($x, $y) = '.($commutes eq '+' ? '' : '$_[2] ? @_[1,0] : ').'@_[0,1];
1194             $_ = r2C $_ for grep ref $_ ne __PACKAGE__, $x, $y;
1195             '.$func.'($x, $y);
1196 8 50       28 }'; #need to swap?
1197 8 100       45 die if $@;
1198             ($op, $sub, exists $NO_MUTATE{$op} ? () : ("$op=", $sub));
1199             }
1200              
1201 6     6   19 sub _gen_unop {
1202 1     1   8 my ($op, $func) = split '@', $_[0];
  1         2  
  1         429  
1203 6 50       49 no strict 'refs';
1204 6     0   339 *$op = \&$func if $op =~ /\w+/; # create an alias
  0         0  
  2         71  
  0         0  
  0         0  
  0         0  
  0         0  
1205             ($op, eval 'sub { '.$func.' $_[0] }');
1206             }
1207              
1208             sub initialize {
1209             # Bless a null PDL into the supplied 1st arg package
1210 661   66 661 0 382300 # If 1st arg is a ref, get the package from it
1211             bless PDL->null, ref($_[0]) || $_[0];
1212             }
1213              
1214             # so broadcasting doesn't also assign the real value into the imaginary
1215 3 50   3 0 13 sub Cassgn {
1216 3 50       12 my @args = !$_[2] ? @_[1,0] : @_[0,1];
1217 3         122 $args[1] = r2C($args[1]) if ref $args[1] ne __PACKAGE__;
1218 3         15 PDL::Ops::assgn(@args);
1219             $args[1];
1220             }
1221              
1222             use overload
1223             (map _gen_biop($_), qw(++Cadd --Csub *+Cmul /-Cdiv **-Cpow atan2-Catan2 ==+Ceq .=-Cassgn)),
1224 4     4   1429 (map _gen_unop($_), qw(sin@Csin cos@Ccos exp@Cexp abs@Cabs log@Clog sqrt@Csqrt)),
1225 2     2   36 (map +($_ => sub { confess "Can't compare complex numbers" }), qw(< > <= >=)),
1226 10 100   10   1673 "!=" => sub { !($_[0] == $_[1]) },
1227 1     1   34 '""' => sub { $_[0]->isnull ? "PDL::Complex->null" : $_[0]->as_native->string },
  1         3  
  1         5  
1228             ;
1229              
1230 5     5 0 684 sub sum {
1231 5 50       27 my($x) = @_;
1232 5         37 return $x if $x->dims==1;
1233 5         213 my $tmp = $x->mv(0,-1)->clump(-2)->mv(1,0)->sumover;
1234             return $tmp;
1235             }
1236              
1237 14     14 0 256 sub sumover{
1238 14         58 my $m = shift;
1239             PDL::Ufunc::sumover($m->transpose);
1240             }
1241              
1242             *PDL::Complex::Csumover=\&sumover; # define through alias
1243              
1244             *PDL::Complex::prodover=\&Cprodover; # define through alias
1245              
1246 4     4 0 18 sub prod {
1247 4 50       22 my($x) = @_;
1248 4         52 return $x if $x->dims==1;
1249 4         155 my $tmp = $x->mv(0,-1)->clump(-2)->mv(1,0)->prodover;
1250             return $tmp;
1251             }
1252              
1253             =head1 AUTHOR
1254              
1255             Copyright (C) 2000 Marc Lehmann .
1256             All rights reserved. There is no warranty. You are allowed
1257             to redistribute this software / documentation as described
1258             in the file COPYING in the PDL distribution.
1259              
1260             =head1 SEE ALSO
1261              
1262             perl(1), L.
1263              
1264             =cut
1265             EOD
1266              
1267             pp_done;
1268