NAME
PDL::LinearAlgebra::Trans - Linear Algebra based transcendental functions for PDL
SYNOPSIS
$a = random (100,100);
$sqrt = msqrt( $a );
|
DESCRIPTION
This module provides some transcendental functions for matrices. Moreover it provides sec, asec, sech, asech, cot, acot, acoth, coth, csc, acsc, csch, acsch. Beware, importing this module will overwrite the hidden PDL routine sec. If you need to call it specify its origin module : PDL::Basic::sec(args)
Computes exp(t*A), the matrix exponential of a general matrix, using the irreducible rational Pade approximation to the exponential function exp(x) = r(x) = (+/-)( I + 2*(q(x)/p(x)) ), combined with scaling-and-squaring and optionally normalization of the trace. The algorithm is described in Roger B. Sidje (rbs@maths.uq.edu.au) "EXPOKIT: Software Package for Computing Matrix Exponentials". ACM - Transactions On Mathematical Software, 24(1):130-156, 1998
A: On input argument matrix. On output exp (t *A ).
Use Fortran storage type.
deg: the degre of the diagonal Pade to be used.
a value of 6 is generally satisfactory.
scale: time -scale (can be < 0).
trace: on input, boolean value indicating whether or not perform
a trace normalization. On output value used.
ns: on output number of scaling-squaring used.
info: exit flag.
0 - no problem
> 0 - Singularity in LU factorization when solving
Pade approximation
|
$a = random(5,5);
$trace = pdl(1);
$a ->t->geexp(6,1, $trace , ( $ns = null), ( $info = null));
|
Complex version of geexp. The value used for trace normalization is not returned. The algorithm is described in Roger B. Sidje (rbs@maths.uq.edu.au) "EXPOKIT: Software Package for Computing Matrix Exponentials". ACM - Transactions On Mathematical Software, 24(1):130-156, 1998
', GenericTypes => [D], Code => ' /* ----------------------------------------------------------------------| int dgpadm_(integer *ideg, integer *m, double *t, PDL_CDouble *h__, integer *ldh, PDL_CDouble *wsp, integer *lwsp, integer *ipiv, integer *iexph, integer *ns, integer *iflag)
-----Purpose----------------------------------------------------------|
Computes exp (t *H ), the matrix exponential of a general matrix in
full, using the irreducible rational Pade approximation to the
exponential function exp (x) = r(x) = (+/-)( I + 2*( q(x) /p(x)) ),
combined with scaling-and-squaring.
-----Arguments--------------------------------------------------------|
ideg|deg : (input) the degre of the diagonal Pade to be used.
a value of 6 is generally satisfactory.
m : (input) order of H.
H|A(ldh,m): (input) argument matrix.
t|scale : (input) time -scale (can be < 0).
wsp| exp (lwsp) : (workspace/output) lwsp .ge. 4 *m *m +ideg+1.
ipiv(m) : (workspace)
>>>> iexph : (output) number such that wsp(iexph) points to exp (tH)
i.e., exp (tH) is located at wsp(iexph ... iexph+m *m -1)
^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
NOTE: if the routine was called with wsp(iptr),
then exp (tH) will start at wsp(iptr+iexph-1).
ns : (output) number of scaling-squaring used.
iflag|info : (output) exit flag.
0 - no problem
> 0 - Singularity in LU factorization when solving
Pade approximation
----------------------------------------------------------------------|
Roger B. Sidje (rbs @maths .uq.edu.au)
EXPOKIT: Software Package for Computing Matrix Exponentials.
ACM - Transactions On Mathematical Software, 24(1):130-156, 1998
Gregory Vanuxem (vanuxemg @yahoo .fr)
Minor modifications (Fortran to C with help of f2c, memory allocation,
PDL adaptation, null entry matrix, error return ,
exp (tH) in entry matrix), trace normalization.
----------------------------------------------------------------------|
*/
static PDL_CDouble c_b1 = 0;
static PDL_CDouble c_b2 = 1;
static integer c__2 = 2;
static integer c__1 = 1;
static double c_b19 = 2.;
static double c_b21 = -1.;
integer i__1, i__2, i__3, i__, j, k;
integer ip, mm, iq, ih2, iodd, iget, iput, ifree, iused;
double hnorm,d__1, d__2;
PDL_CDouble scale, tracen, cp, cq, z__1;
PDL_CDouble scale2;
PDL_CDouble *wsp , *h__ ;
extern int FORTRAN(zgemm)(char *, char *, integer *, integer *,
integer *, PDL_CDouble *, PDL_CDouble *, integer *,
PDL_CDouble *, integer *, PDL_CDouble *, PDL_CDouble *,
integer *, ftnlen, ftnlen);
extern int FORTRAN(zgesv)(integer *, integer *, PDL_CDouble *,
integer *, integer *, PDL_CDouble *, integer *, integer *);
extern int FORTRAN(zaxpy)(integer *, PDL_CDouble *,
PDL_CDouble *, integer *, PDL_CDouble *, integer *), FORTRAN(zdscal)(
integer *, double *, PDL_CDouble *, integer *);
extern double FORTRAN(zlange)(char *norm , integer *m , integer *n , PDL_CDouble *a ,
integer *lda , PDL_CDouble *work );
h__ = (PDL_CDouble *) $P (A);
wsp = malloc (( $deg () + 1 + 4 * $SIZE (n) * $SIZE (n)) * sizeof (PDL_CDouble));
mm = $SIZE (n) * $SIZE (n);
$info () = 0;
ih2 = $deg () + 1;
ip = ih2 + mm;
iq = ip + mm;
ifree = iq + mm;
tracen = 0;
if ( $trace ())
{
loop(n0) %{
tracen += h__[n0 + n0 * $SIZE (n)];
%}
tracen /= $SIZE (n);
if (creal(tracen) < 0){
tracen = cimag(tracen);
}
loop(n0) %{
h__[n0 + n0 * $SIZE (n)] -= tracen;
%}
}
// --- scaling: seek ns such that ||t *H /2^ns|| < 1/2;
// and set scale = t/2^ns ...
// Compute Infinity norm
hnorm = FORTRAN(zlange)( "I" , &(integer){ $SIZE (n)}, &(integer){ $SIZE (n)}, h__,
&(integer){ $SIZE (n)}, wsp);
hnorm = abs ( $scale () * hnorm);
if (hnorm == 0.)
{
$ns () = 0;
loop(n0,n1) %{
h__[n1 + n0 * $SIZE (n) ] = (n0 == n1) ? 1 : 0;
%}
}
else
{
i__2 = (integer) (( log (hnorm) / log (2.)) + 2);
$ns () = max(0,i__2);
scale = $scale () / pow(2, (double ) $ns ());
scale2 = scale * scale;
// --- compute Pade coefficients ...
i__ = $deg () + 1;
j = ( $deg () << 1) + 1;
wsp[0] = 1;
for (k = 1; k <= $deg (); ++k) {
d__1 = (double) (i__ - k);
d__2 = (double) (k * (j - k));
wsp[k] = (d__1 * wsp[k - 1]) / d__2;
}
// --- H2 = scale2 *H *H ...
FORTRAN(zgemm)( "n" , "n" , &(integer){ $SIZE (n)}, &(integer){ $SIZE (n)}, &(integer){ $SIZE (n)}, &scale2 , &h__ [0], &(integer){ $SIZE (n)}, &h__ [0],
&(integer){ $SIZE (n)}, &c_b1 , &wsp [ih2], &(integer){ $SIZE (n)}, (ftnlen)1, (ftnlen)1);
// --- initialise p (numerator) and q (denominator) ...
cp = wsp[ $deg () - 1];
cq = wsp[ $deg ()];
loop(n0) %{
loop(n1) %{
i__3 = ip + n0 * $SIZE (n) + n1;
wsp[i__3] = 0.;
i__3 = iq + n0 * $SIZE (n) + n1;
wsp[i__3] = 0.;
%}
i__2 = ip + n0 * ( $SIZE (n) + 1);
wsp[i__2] = cp;
i__2 = iq + n0 * ( $SIZE (n) + 1);
wsp[i__2] = cq;
%}
// --- Apply Horner rule ...
iodd = 1;
k = $deg () - 2;
do
{
iused = iodd * iq + (1 - iodd) * ip;
FORTRAN(zgemm)( "n" , "n" , &(integer){ $SIZE (n)}, &(integer){ $SIZE (n)}, &(integer){ $SIZE (n)}, &c_b2 , &wsp [iused], &(integer){ $SIZE (n)}, &wsp [ih2], &(integer){ $SIZE (n)}, &c_b1 , &
wsp[ifree], &(integer){ $SIZE (n)}, (ftnlen)1, (ftnlen)1);
loop(n0) %{
i__1 = ifree + n0 * ( $SIZE (n) + 1);
i__2 = ifree + n0 * ( $SIZE (n) + 1);
i__3 = k;
wsp[i__1] = wsp[i__2] + wsp[i__3];
%}
ip = (1 - iodd) * ifree + iodd * ip;
iq = iodd * ifree + (1 - iodd) * iq;
ifree = iused;
iodd = 1 - iodd;
--k;
}
while (k >= 0);
// --- Obtain (+/-)(I + 2*(p\ q)) ...
if (iodd != 0)
{
FORTRAN(zgemm)( "n" , "n" , &(integer){ $SIZE (n)}, &(integer){ $SIZE (n)}, &(integer){ $SIZE (n)}, &scale , &wsp [iq], &(integer){ $SIZE (n)}, &h__ [0], &(integer){ $SIZE (n)}, &
c_b1, &wsp [ifree], &(integer){ $SIZE (n)}, (ftnlen)1, (ftnlen)1);
iq = ifree;
} else
{
FORTRAN(zgemm)( "n" , "n" , &(integer){ $SIZE (n)}, &(integer){ $SIZE (n)}, &(integer){ $SIZE (n)}, &scale , &wsp [ip], &(integer){ $SIZE (n)}, &h__ [0], &(integer){ $SIZE (n)}, &
c_b1, &wsp [ifree], &(integer){ $SIZE (n)}, (ftnlen)1, (ftnlen)1);
ip = ifree;
}
z__1 = -1.;
FORTRAN(zaxpy)( &mm , &z__1 , &wsp [ip], &c__1 , &wsp [iq], &c__1 );
FORTRAN(zgesv)(&(integer){ $SIZE (n)}, &(integer){ $SIZE (n)}, &wsp [iq], &(integer){ $SIZE (n)}, $P (ipiv), &wsp [ip], &(integer){ $SIZE (n)}, $P (info));
if ( $info () == 0)
{
FORTRAN(zdscal)( &mm , &c_b19 , &wsp [ip], &c__1 );
loop(n0) %{
i__2 = ip + n0 * ( $SIZE (n) + 1);
i__3 = ip + n0 * ( $SIZE (n) + 1);
wsp[i__2] = wsp[i__3] + 1.;
%}
iput = ip;
if ( $ns () == 0 && iodd != 0)
FORTRAN(zdscal)( &mm , &c_b21 , &wsp [ip], &c__1 );
else
{
// -- squaring : exp (t *H ) = ( exp (t *H ))^(2^ns) ...
iodd = 1;
for (k = 0; k < $ns (); k++)
{
iget = iodd * ip + (1 - iodd) * iq;
iput = (1 - iodd) * ip + iodd * iq;
FORTRAN(zgemm)( "n" , "n" , &(integer){ $SIZE (n)}, &(integer){ $SIZE (n)}, &(integer){ $SIZE (n)}, &c_b2 , &wsp [iget], &(integer){ $SIZE (n)}, &wsp [iget], &(integer){ $SIZE (n)}, &c_b1 ,
&wsp [iput], &(integer){ $SIZE (n)}, (ftnlen)1, (ftnlen)1);
iodd = 1 - iodd;
}
}
}
i__2 = iput + mm;
i__ = 0;
if ( $trace ())
{
tracen = cexp(tracen);
for (i__1 = iput; i__1 < i__2 ; i__1++)
{
h__[i__++] = tracen * wsp[i__1];
}
}
else
{
for (i__1 = iput; i__1 < i__2 ; i__1++)
{
h__[i__++] = wsp[i__1];
}
}
}
free(wsp);
',
);
|
pp_defc('trsqrt', Pars => '[io]A(2,n,n);int uplo();[o] B(2,n,n);int [o]info()', Doc => '
Root square of complex triangular matrix. Uses a recurrence of Bj�rck and Hammarling. (See Nicholas J. Higham. A new sqrtm for MATLAB. Numerical Analysis Report No. 336, Manchester Centre for Computational Mathematics, Manchester, England, January 1999. It\'s available at http://www.ma.man.ac.uk/~higham/pap-mf.html) If uplo is true, A is lower triangular.
', GenericTypes => [D], Code => ' PDL_CDouble *cb, *ca; PDL_CDouble s, snum, sdenum; double *b; double tt, dn; integer i, j, k, l, ind, ind1, ind2;
ca = (PDL_CDouble *) $P (A);
b = (void *) $P (B);
cb = (PDL_CDouble *) $P (B);
|
#define subscr(a_1,a_2)\ ($uplo()) ? (a_2) * $SIZE(n) + (a_1) : (a_1) * $SIZE(n) + (a_2)
$info () = 0;
ind = $SIZE (n) * $SIZE (n) * 2;
for (i = 0; i < ind;i++)
b[i] = 0;
loop(n0) %{
ind = subscr(n0,n0);
cb[ind] = csqrt(ca[ind]);
%}
for (k = 0; k < $SIZE (n)-1; k++)
{
for (i = 0; i < $SIZE (n)-(k+1); i++)
{
j = i + k + 1;
ind = subscr(i,j);
s = 0;
for (l = i+1; l < j; l++)
{
ind1 = subscr(i,l);
ind2 = subscr(l,j);
s += cb[ind1] * cb[ind2];
}
ind1 = subscr(i,i);
ind2 = subscr(j,j);
sdenum = cb[ind1] + cb[ind2];
snum = ca[ind] - s;
if (fabs (creal(sdenum)) > fabs (cimag(sdenum)))
{
if (creal(sdenum) == 0){
$info () = -1;
goto ENDCTRSQRT;
}
tt = cimag(sdenum) / creal(sdenum);
dn = creal(sdenum) + tt * cimag(sdenum);
cb[ind] = (snum + tt * (-I * snum)) / dn;
}
else
{
if (cimag(sdenum) == 0){
$info () = -1;
goto ENDCTRSQRT;
}
tt = creal(sdenum) / cimag(sdenum);
dn = creal(sdenum) * tt + cimag(sdenum);
cb[ind] = (snum * tt + (-I * snum)) / dn;
}
}
}
ENDCTRSQRT:
;
',
);
|
pp_addhdr(' void dfunc_wrapper(PDL_CDouble *p, integer n, SV* dfunc); ');
pp_defc('trfun', Pars => '[io]A(2,n,n);int uplo();[o] B(2,n,n);int [o]info(); complex [t]diag(n);', OtherPars => "SV* func" , Doc => ' =for ref
Apply an arbitrary function to a complex triangular matrix. Uses a recurrence of Parlett. If uplo is true, A is lower triangular. ', GenericTypes => [D], Code => ' PDL_CDouble *cb, *ca; PDL_CDouble s, snum, sdenum; double *b; double tt, dn; integer i, j, k, l, p, ind, ind1, ind2, tmp, c__1; extern double FORTRAN(zdotu)(PDL_CDouble *ret, integer *n, PDL_CDouble *dx, integer *incx, PDL_CDouble *dy, integer *incy);
ca = (PDL_CDouble *) $P (A);
b = (double *) $P (B);
cb = (PDL_CDouble *) $P (B);
c__1 = 1;
|
#define subscr(a_1,a_2) \ ( $uplo() ) ? (a_2)+ $SIZE(n) * (a_1) : (a_1)+ $SIZE(n) * (a_2)
$info () = 0;
ind = $SIZE (n) * $SIZE (n) * 2;
for (i = 0; i < ind;i++)
b[i] = 0;
loop (n) %{
$diag () = ca[subscr(n,n)];
%}
dfunc_wrapper( $P (diag), $SIZE (n), $COMP (func));
loop (n) %{
cb[subscr(n,n)] = $diag ();
%}
for (p = 1; p < $SIZE (n); p++)
{
tmp = $SIZE (n) - p;
for (i = 0; i < tmp; i++)
{
j = i + p;
// $s = $T (,( $j ),( $i ))->Cmul( $F (,( $j ),( $j ))->Csub( $F (,( $i ),( $i ))));
ind1 = subscr(i,i);
ind2 = subscr(j,j);
ind = subscr(j,i);
s = snum = ca[ind] * (cb[ind2] - cb[ind1]);
if (i < (j-1))
{
// $s = $s + $T (, $i +1: $j -1,( $i ))->cdot(1, $F (,( $j ), $i +1: $j -1),1)->Csub( $F (, $i +1: $j -1,( $i ))->cdot(1, $T (,( $j ), $i +1: $j -1),1));
// $T (, $i +1: $j -1,( $i ))->cdot(1, $F (,( $j ), $i +1: $j -1),1)
l = 0;
for (k = i+1; k < j; k++)
$diag ( n =>l++) = cb[subscr(j,k)];
ind = subscr(i+1,i);
// TODO : Humm
FORTRAN(zdotu)( &s , &l , &ca [ind], &c__1 , $P (diag), &c__1 );
snum += s;
// $F (, $i +1: $j -1,( $i ))->cdot(1, $T (,( $j ), $i +1: $j -1),1)
l = 0;
for (k = i+1; k < j; k++)
$diag ( n =>l++) = ca[subscr(j,k)];
ind = subscr(i+1,i);
// TODO : Humm
FORTRAN(zdotu)( &s , &l , &cb [ind], &c__1 , $P (diag), &c__1 );
snum -= s;
ind = subscr(j,i);
}
// $sdenum = $T (,( $j ),( $j ))->Csub( $T (,( $i ),( $i )));
sdenum = ca[ind2] - ca[ind1];
// $s = $s / $sdenum ;
if (fabs (creal(sdenum)) > fabs (cimag(sdenum)))
{
if (creal(sdenum) == 0)
{
$info () = -1;
goto ENDCTRFUN;
}
tt = cimag(sdenum) / creal(sdenum);
dn = creal(sdenum) + tt * cimag(sdenum);
cb[ind] = (snum + tt * (-I * snum)) / dn;
}
else
{
if (cimag(sdenum) == 0)
{
$info () = -1;
goto ENDCTRFUN;
}
tt = creal(sdenum) / cimag(sdenum);
dn = creal(sdenum) * tt + cimag(sdenum);
cb[ind] = (snum * tt + (-I * snum)) / dn;
}
}
}
ENDCTRFUN:
;
',
);
|
pp_addpm(<<'EOD'); my $pi; BEGIN { $pi = pdl(3.1415926535897932384626433832795029) } sub pi () { $pi->copy };
*sec = \&PDL::sec; sub PDL::sec{1/cos($_[0])}
*csc = \&PDL::csc; sub PDL::csc($) {1/sin($_[0])}
*cot = \&PDL::cot; sub PDL::cot($) {1/(sin($_[0])/cos($_[0]))}
*sech = \&PDL::sech; sub PDL::sech($){1/pdl($_[0])->cosh}
*csch = \&PDL::csch; sub PDL::csch($) {1/pdl($_[0])->sinh}
*coth = \&PDL::coth; sub PDL::coth($) {1/pdl($_[0])->tanh}
*asec = \&PDL::asec; sub PDL::asec($) {my $tmp = 1/pdl($_[0]) ; $tmp->acos}
*acsc = \&PDL::acsc; sub PDL::acsc($) {my $tmp = 1/pdl($_[0]) ; $tmp->asin}
*acot = \&PDL::acot; sub PDL::acot($) {my $tmp = 1/pdl($_[0]) ; $tmp->atan}
*asech = \&PDL::asech; sub PDL::asech($) {my $tmp = 1/pdl($_[0]) ; $tmp->acosh}
*acsch = \&PDL::acsch; sub PDL::acsch($) {my $tmp = 1/pdl($_[0]) ; $tmp->asinh}
*acoth = \&PDL::acoth; sub PDL::acoth($) {my $tmp = 1/pdl($_[0]) ; $tmp->atanh}
my $_tol = 9.99999999999999e-15;
sub toreal { $_tol = $_[1] if defined $_[1]; return $_[0] if $_[0]->isempty or $_[0]->type->real; my ($min, $max, $tmp); ($min, $max) = $_[0]->im->minmax; return $_[0]->re->sever unless (abs($min) > $_tol || abs($max) > $_tol); $_[0]; }
mlog
Return matrix logarithm of a square matrix.
my $a = random(10,10);
my $log = mlog( $a );
|
msqrt
Return matrix square root (principal) of a square matrix.
my $a = random(10,10);
my $sqrt = msqrt( $a );
|
mexp
Return matrix exponential of a square matrix.
my $a = random(10,10);
my $exp = mexp( $a );
|
mpow
Return matrix power of a square matrix.
PDL = mpow(PDL(A), SCALAR(exponent))
|
my $a = random(10,10);
my $powered = mpow( $a ,2.5);
|
mcos
Return matrix cosine of a square matrix.
my $a = random(10,10);
my $cos = mcos( $a );
|
macos
Return matrix inverse cosine of a square matrix.
my $a = random(10,10);
my $acos = macos( $a );
|
msin
Return matrix sine of a square matrix.
my $a = random(10,10);
my $sin = msin( $a );
|
masin
Return matrix inverse sine of a square matrix.
my $a = random(10,10);
my $asin = masin( $a );
|
mtan
Return matrix tangent of a square matrix.
my $a = random(10,10);
my $tan = mtan( $a );
|
matan
Return matrix inverse tangent of a square matrix.
my $a = random(10,10);
my $atan = matan( $a );
|
mcot
Return matrix cotangent of a square matrix.
my $a = random(10,10);
my $cot = mcot( $a );
|
macot
Return matrix inverse cotangent of a square matrix.
my $a = random(10,10);
my $acot = macot( $a );
|
msec
Return matrix secant of a square matrix.
my $a = random(10,10);
my $sec = msec( $a );
|
masec
Return matrix inverse secant of a square matrix.
my $a = random(10,10);
my $asec = masec( $a );
|
mcsc
Return matrix cosecant of a square matrix.
my $a = random(10,10);
my $csc = mcsc( $a );
|
macsc
Return matrix inverse cosecant of a square matrix.
my $a = random(10,10);
my $acsc = macsc( $a );
|
mcosh
Return matrix hyperbolic cosine of a square matrix.
my $a = random(10,10);
my $cos = mcosh( $a );
|
macosh
Return matrix hyperbolic inverse cosine of a square matrix.
my $a = random(10,10);
my $acos = macosh( $a );
|
msinh
Return matrix hyperbolic sine of a square matrix.
my $a = random(10,10);
my $sinh = msinh( $a );
|
masinh
Return matrix hyperbolic inverse sine of a square matrix.
my $a = random(10,10);
my $asinh = masinh( $a );
|
mtanh
Return matrix hyperbolic tangent of a square matrix.
my $a = random(10,10);
my $tanh = mtanh( $a );
|
matanh
Return matrix hyperbolic inverse tangent of a square matrix.
my $a = random(10,10);
my $atanh = matanh( $a );
|
mcoth
Return matrix hyperbolic cotangent of a square matrix.
my $a = random(10,10);
my $coth = mcoth( $a );
|
macoth
Return matrix hyperbolic inverse cotangent of a square matrix.
my $a = random(10,10);
my $acoth = macoth( $a );
|
msech
Return matrix hyperbolic secant of a square matrix.
my $a = random(10,10);
my $sech = msech( $a );
|
masech
Return matrix hyperbolic inverse secant of a square matrix.
my $a = random(10,10);
my $asech = masech( $a );
|
mcsch
Return matrix hyperbolic cosecant of a square matrix.
my $a = random(10,10);
my $csch = mcsch( $a );
|
macsch
Return matrix hyperbolic inverse cosecant of a square matrix.
my $a = random(10,10);
my $acsch = macsch( $a );
|
mfun
Return matrix function of second argument of a square matrix. Function will be applied on a complex ndarray.
my $a = random(10,10);
my $fun = mfun( $a , 'cos' );
sub sinbycos2{
$_ [0]->set_inplace(0);
$_ [0] .= $_ [0]->Csin/ $_ [0]->Ccos**2;
}
$fun = mfun( $a , \ &sinbycos2 ,1);
$fun = mfun( $a , \ &sinbycos2 );
scalar msolve( $a ->mcos->mpow(2), $a ->msin);
|
TODO
Improve error return and check singularity. Improve (msqrt,mlog) / r2C
AUTHOR
Copyright (C) Gr�gory Vanuxem 2005-2018.
This library is free software; you can redistribute it and/or modify it under the terms of the Perl Artistic License as in the file Artistic_2 in this distribution.