man PDL::Primitive () - primitive operations for pdl

NAME

PDL::Primitive - primitive operations for pdl

DESCRIPTION

This module provides some primitive and useful functions defined using PDL::PP and able to use the new indexing tricks.

See PDL::Indexing for how to use indices creatively. For explanation of the signature format, see PDL::PP.

SYNOPSIS

 use PDL::Primitive;

FUNCTIONS

inner

  Signature: (a(n); b(n); [o]c())

Inner product over one dimension

 c = sum_i a_i * b_i

If CWa() * b() contains only bad data, CWc() is set bad. Otherwise CWc() will have its bad flag cleared, as it will not contain any bad values.

outer

  Signature: (a(n); b(m); [o]c(n,m))

outer product over one dimension

Naturally, it is possible to achieve the effects of outer product simply by threading over the "CW*" operator but this function is provided for convenience.

outer does handle bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

matmult

 Signature: (a(i,x),b(z,i),[o]c(x,z))

Matrix multiplication

We peruse the inner product to define matrix multiplication via a threaded inner product

innerwt

  Signature: (a(n); b(n); c(n); [o]d())

Weighted (i.e. triple) inner product

 d = sum_i a(i) b(i) c(i)

innerwt does handle bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

inner2

  Signature: (a(n); b(n,m); c(m); [o]d())

Inner product of two vectors and a matrix

 d = sum_ij a(i) b(i,j) c(j)

Note that you should probably not thread over CWa and CWc since that would be very wasteful. Instead, you should use a temporary for CWb*c.

inner2 does handle bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

inner2d

  Signature: (a(n,m); b(n,m); [o]c())

Inner product over 2 dimensions.

Equivalent to

 $c = inner($a->clump(2), $b->clump(2))

inner2d does handle bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

inner2t

  Signature: (a(j,n); b(n,m); c(m,k); [t]tmp(n,k); [o]d(j,k)))

Efficient Triple matrix product CWa*b*c

Efficiency comes from by using the temporary CWtmp. This operation only scales as CWN**3 whereas threading using inner2 would scale as CWN**4.

The reason for having this routine is that you do not need to have the same thread-dimensions for CWtmp as for the other arguments, which in case of large numbers of matrices makes this much more memory-efficient.

It is hoped that things like this could be taken care of as a kind of closures at some point.

inner2t does handle bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

crossp

  Signature: (a(tri=3); b(tri); [o] c(tri))

Cross product of two 3D vectors

After

 $c = crossp $a, $b

the inner product CW$c*$a and CW$c*$b will be zero, i.e. CW$c is orthogonal to CW$a and CW$b

crossp does not process bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

norm

  Signature: (vec(n); [o] norm(n))

Normalises a vector to unit Euclidean length

norm does handle bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

indadd

  Signature: (a(); int ind(); [o] sum(m))

Threaded Index Add: Add CWa to the CWind element of CWsum, i.e:

 sum(ind) += a

Simple Example:

  $a = 2;
  $ind = 3;
  $sum = zeroes(1);
  indadd($a,$ind, $sum);
  print $sum
  #Result: ( 2 added to element 3 of $sum)
  # [0 0 0 2 0 0 0 0 0 0]

Threaded Example:

  $a = pdl( 1,2,3);
  $ind = pdl( 1,4,6);
  $sum = zeroes(1);
  indadd($a,$ind, $sum);
  print $sum."\n";
  #Result: ( 1, 2, and 3 added to elements 1,4,6 $sum)
  # [0 1 0 0 2 0 3 0 0 0]

The routine barfs if any of the indices are bad.

conv1d

  Signature: (a(m); kern(p); [o]b(m); int reflect)

1d convolution along first dimension

  $con = conv1d sequence(1), pdl(-1,0,1), {Boundary => 'reflect'};

By default, periodic boundary conditions are assumed (i.e. wrap around). Alternatively, you can request reflective boundary conditions using the CWBoundary option:

  {Boundary => 'reflect'} # case in 'reflect' doesn't matter

The convolution is performed along the first dimension. To apply it across another dimension use the slicing routines, e.g.

  $b = $a->mv(2,0)->conv1d($kernel)->mv(0,2); # along third dim

This function is useful for threaded filtering of 1D signals.

Compare also conv2d, convolve, fftconvolve, fftwconv, rfftwconv

conv1d does not process bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

in

  Signature: (a(); b(n); [o] c())

test if a is in the set of values b

   $goodmsk = $labels->in($goodlabels);
   print pdl(4,3,1)->in(pdl(2,3,3));
  [0 1 0]

CWin is akin to the is an element of of set theory. In priciple, PDL threading could be used to achieve its functionality by using a construct like

   $msk = ($labels->dummy(0) == $goodlabels)->orover;

However, CWin doesn't create a (potentially large) intermediate and is generally faster.

in does not process bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

uniq

return all unique elements of a piddle

The unique elements are returned in ascending order.

  print pdl(2,2,2,4,0,-1,6,6)->uniq;
 [-1 0 2 4 6]

Note: The returned pdl is 1D; any structure of the input piddle is lost.

See uniqind if you need the indices of the unique elements rather than the values.

uniqind

return the indices of all unique elements of a piddle The order is in the order of the values to be consistent with uniq

  print pdl(2,2,2,4,0,-1,6,6)->uniqind;
         [5, 4, 1, 3, 6]

Note: The returned pdl is 1D; any structure of the input piddle is lost.

See uniq if you want the unique values instead of the indices.

uniqvec

return all unique vectors out of a collection

The unique vectors are returned in lexicographically sorted ascending order. The 0th dimension of the input PDL is treated as a dimensional index within each vector, and the 1st and any higher dimensions are taken to run across vectors. The return value is always 2D; any structure of the input PDL (beyond using the 0th dimension for vector index) is lost.

See also uniq for a uniqe list of scalars; and qsortvec for sorting a list of vectors lexicographcally.

hclip

  Signature: (a(); b(); [o] c())

clip CW$a by CW$b (CW$b is upper bound)

hclip does handle bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

lclip

  Signature: (a(); b(); [o] c())

clip CW$a by CW$b (CW$b is lower bound)

lclip does handle bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

clip

Clip a piddle by (optional) upper or lower bounds.

 $b = $a->clip(0,3);
 $c = $a->clip(undef, $x);

clip handles bad values since it is just a wrapper around hclip and lclip.

wtstat

  Signature: (a(n); wt(n); avg(); [o]b(); int deg)

Weighted statistical moment of given degree

This calculates a weighted statistic over the vector CWa. The formula is

 b() = (sum_i wt_i * (a_i ** degree - avg)) / (sum_i wt_i)

Bad values are ignored in any calculation; CW$b will only have its bad flag set if the output contains any bad data.

statsover

  Signature: (a(n); w(n); float+ [o]avg(); float+ [o]prms(); int+ [o]min(); int+ [o]max(); float+ [o]adev(); float+ [o]rms())

Calculate useful statistics over a dimension of a piddle

  ($mean, $rms, $median, $min, $max, $adev) = statover($piddle, $weights);

This utility function calculates various useful quantities of a piddle. These are:

* the mean:
  MEAN = sum (x)/ N
with CWN being the number of elements in x
* RMS deviation from the mean:
  RMS = sqrt(sum( (x-mean(x))^2 )/N)
(also known as the root-mean-square deviation, or the square root of the variance)
* the median
The median is the 50th percentile data value. Median is found by medover, so WEIGHTING IS IGNORED FOR THE MEDIAN CALCULATION.
* the minimum
* the maximum
* the absolute deviation:
  ADEV = sqrt(sum( abs(x-mean(x)) )/N)
(This is also called the standard deviation)
* the population RMS deviation from the mean:
  PRMS = sqrt( sum( (x-mean(x))^2 )/(N-1)
The population deviation is the best-estimate of the deviation of the population from which a sample is drawn.

This operator is a projection operator so the calculation will take place over the final dimension. Thus if the input is N-dimensional each returned value will be N-1 dimensional, to calculate the statistics for the entire piddle either use CWclump(-1) directly on the piddle or call CWstats.

Bad values are simply ignored in the calculation, effectively reducing the sample size. If all data are bad then the output data are marked bad.

stats

Calculates useful statistics on a piddle

 ($mean,$prms,$median,$min,$max,$adev,$rms) = stats($piddle,[$weights]);

This utility calculates all the most useful quantities in one call. It works the same way as statsover, except that the quantities are calculated considering the entire input PDL as a single sample, rather than as a collection of rows.

Bad values are handled; if all input values are bad, then all of the output values are flagged bad.

histogram

  Signature: (in(n); int+[o] hist(m); double step; double min; int msize => m)

Calculates a histogram for given stepsize and minimum.

 $h = histogram($data, $step, $min, $numbins);
 $hist = zeroes $numbins;  # Put histogram in existing piddle.
 histogram($data, $hist, $step, $min, $numbins);

The histogram will contain CW$numbins bins starting from CW$min, each CW$step wide. The value in each bin is the number of values in CW$data that lie within the bin limits.

Data below the lower limit is put in the first bin, and data above the upper limit is put in the last bin.

The output is reset in a different threadloop so that you can take a histogram of CW$a(10,12) into CW$b(1) and get the result you want.

Use hist instead for a high-level interface.

 perldl> p histogram(pdl(1,1,2),1,0,3)
 [0 2 1]

histogram does handle bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

whistogram

  Signature: (in(n); float+ wt(n);float+[o] hist(m); double step; double min; int msize => m)

Calculates a histogram from weighted data for given stepsize and minimum.

 $h = whistogram($data, $weights, $step, $min, $numbins);
 $hist = zeroes $numbins;  # Put histogram in existing piddle.
 whistogram($data, $weights, $hist, $step, $min, $numbins);

The histogram will contain CW$numbins bins starting from CW$min, each CW$step wide. The value in each bin is the sum of the values in CW$weights that correspond to values in CW$data that lie within the bin limits.

Data below the lower limit is put in the first bin, and data above the upper limit is put in the last bin.

The output is reset in a different threadloop so that you can take a histogram of CW$a(10,12) into CW$b(1) and get the result you want.

 perldl> p whistogram(pdl(1,1,2), pdl(0.1,0.1,0.5), 1, 0, 4)
 [0 0.2 0.5 0]

whistogram does handle bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

histogram2d

  Signature: (ina(n); inb(n); int+[o] hist(ma,mb); double stepa; double mina; int masize => ma;
                     double stepb; double minb; int mbsize => mb;)

Calculates a 2d histogram.

 $h = histogram2d($datax, $datay,
       $stepx, $minx, $nbinx, $stepy, $miny, $nbiny);
 $hist = zeroes $nbinx, $nbiny;  # Put histogram in existing piddle.
 histogram2d($datax, $datay, $hist,
       $stepx, $minx, $nbinx, $stepy, $miny, $nbiny);

The histogram will contain CW$nbinx x CW$nbiny bins, with the lower limits of the first one at CW($minx, $miny), and with bin size CW($stepx, $stepy). The value in each bin is the number of values in CW$datax and CW$datay that lie within the bin limits.

Data below the lower limit is put in the first bin, and data above the upper limit is put in the last bin.

 perldl> p histogram2d(pdl(1,1,1,2,2),pdl(2,1,1,1,1),1,0,3,1,0,3)
 [
  [0 0 0]
  [0 2 2]
  [0 1 0]
 ]

histogram2d does handle bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

whistogram2d

  Signature: (ina(n); inb(n); float+ wt(n);float+[o] hist(ma,mb); double stepa; double mina; int masize => ma;
                     double stepb; double minb; int mbsize => mb;)

Calculates a 2d histogram from weighted data.

 $h = whistogram2d($datax, $datay, $weights,
       $stepx, $minx, $nbinx, $stepy, $miny, $nbiny);
 $hist = zeroes $nbinx, $nbiny;  # Put histogram in existing piddle.
 whistogram2d($datax, $datay, $weights, $hist,
       $stepx, $minx, $nbinx, $stepy, $miny, $nbiny);

The histogram will contain CW$nbinx x CW$nbiny bins, with the lower limits of the first one at CW($minx, $miny), and with bin size CW($stepx, $stepy). The value in each bin is the sum of the values in CW$weights that correspond to values in CW$datax and CW$datay that lie within the bin limits.

Data below the lower limit is put in the first bin, and data above the upper limit is put in the last bin.

 perldl> p whistogram2d(pdl(1,1,1,2,2),pdl(2,1,1,1,1),pdl(0.1,0.2,0.3,0.4,0.5),1,0,3,1,0,3)
 [
  [  0   0   0]
  [  0 0.5 0.9]
  [  0 0.1   0]
 ]

whistogram2d does handle bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

fibonacci

  Signature: ([o]x(n))

Constructor - a vector with Fibonacci's sequence

fibonacci does not process bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

append

  Signature: (a(n); b(m); [o] c(mn))

append two or more piddles by concatenating along their first dimensions

 $a = ones(2,4,7);
 $b = sequence 5;
 $c = $a->append($b);  # size of $c is now (7,4,7) (a jumbo-piddle ;)

CWappend appends two piddles along their first dims. Rest of the dimensions must be compatible in the threading sense. Resulting size of first dim is the sum of the sizes of the first dims of the two argument piddles - ie CWn + m.

append does not process bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

glue

  $c = $a->glue(<dim>,$b,...)

Glue two or more PDLs together along an arbitrary dimension (N-D append).

Sticks CW$a, CW$b, and all following arguments together along the specified dimension. All other dimensions must be compatible in the threading sense.

Glue is permissive, in the sense that every PDL is treated as having an infinite number of trivial dimensions of order 1 so CW$a-glue(3,$b)> works, even if CW$a and CW$b are only one dimensional.

If one of the PDLs has no elements, it is ignored. Likewise, if one of them is actually the undefined value, it is treated as if it had no elements.

If the first parameter is a defined perl scalar rather than a pdl, then it is taken as a dimension along which to glue everything else, so you can say CW$cube = PDL::glue(3,@image_list); if you like.

CWglue is implemented in pdl, using a combination of xchg and append. It should probably be updated (one day) to a pure PP function.

axisvalues

  Signature: ([o,nc]a(n))

Internal routine

CWaxisvalues is the internal primitive that implements axisvals and alters its argument.

axisvalues does not process bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

random

Constructor which returns piddle of random numbers

 $a = random([type], $nx, $ny, $nz,...);
 $a = random $b;

etc (see zeroes).

This is the uniform distribution between 0 and 1 (assumedly excluding 1 itself). The arguments are the same as CWzeroes (q.v.) - i.e. one can specify dimensions, types or give a template.

You can use the perl function srand to seed the random generator. For further details consult Perl's srand documentation.

randsym

Constructor which returns piddle of random numbers

 $a = randsym([type], $nx, $ny, $nz,...);
 $a = randsym $b;

etc (see zeroes).

This is the uniform distribution between 0 and 1 (excluding both 0 and 1, cf random). The arguments are the same as CWzeroes (q.v.) - i.e. one can specify dimensions, types or give a template.

You can use the perl function srand to seed the random generator. For further details consult Perl's srand documentation.

grandom

Constructor which returns piddle of Gaussian random numbers

 $a = grandom([type], $nx, $ny, $nz,...);
 $a = grandom $b;

etc (see zeroes).

This is generated using the math library routine CWndtri.

Mean = 0, Stddev = 1

You can use the perl function srand to seed the random generator. For further details consult Perl's srand documentation.

vsearch

  Signature: (i(); x(n); int [o]ip())

routine for searching 1D values i.e. step-function interpolation.

 $inds = vsearch($vals, $xs);

Returns for each value of CW$vals the index of the least larger member of CW$xs (which need to be in increasing order). If the value is larger than any member of CW$xs, the index to the last element of CW$xs is returned.

This function is useful e.g. when you have a list of probabilities for events and want to generate indices to events:

 $a = pdl(.01,.86,.93,1); # Barnsley IFS probabilities cumulatively
 $b = random 20;
 $c = vsearch($b, $a); # Now, $c will have the appropriate distr.

It is possible to use the cumusumover function to obtain cumulative probabilities from absolute probabilities.

needs major (?) work to handles bad values

interpolate

  Signature: (xi(); x(n); y(n); [o] yi(); int [o] err())

routine for 1D linear interpolation

 ( $yi, $err ) = interpolate($xi, $x, $y)

Given a set of points CW($x,$y), use linear interpolation to find the values CW$yi at a set of points CW$xi.

CWinterpolate uses a binary search to find the suspects, er..., interpolation indices and therefore abscissas (ie CW$x) have to be strictly ordered (increasing or decreasing). For interpolation at lots of closely spaced abscissas an approach that uses the last index found as a start for the next search can be faster (compare Numerical Recipes CWhunt routine). Feel free to implement that on top of the binary search if you like. For out of bounds values it just does a linear extrapolation and sets the corresponding element of CW$err to 1, which is otherwise 0.

See also interpol, which uses the same routine, differing only in the handling of extrapolation - an error message is printed rather than returning an error piddle.

needs major (?) work to handles bad values

interpol

 Signature: (xi(); x(n); y(n); [o] yi())

routine for 1D linear interpolation

 $yi = interpol($xi, $x, $y)

CWinterpol uses the same search method as interpolate, hence CW$x must be strictly ordered (either increasing or decreasing). The difference occurs in the handling of out-of-bounds values; here an error message is printed.

interpND

Interpolate values from an N-D piddle, with switchable method

  $source = 10*xvals(10,10) + yvals(10,10);
  $index = pdl([[2.2,3.5],[4.1,5.0]],[[6.0,7.4],[8,9]]);
  print $source->interpND( $index );

InterpND acts like indexND, collapsing CW$index by lookup into CW$source; but it does interpolation rather than direct sampling. The interpolation method and boundary condition are switchable via an options hash.

By default, linear or sample interpolation is used, with constant value outside the boundaries of the source pdl. No dataflow occurs, because in general the output is computed rather than indexed.

All the interpolation methods treat the pixels as value-centered, so the CWsample method will return CW$a->(0) for coordinate values on the set [-0.5,0.5), and all methods will return CW$a->(1) for a coordinate value of exactly 1.

Recognized options:

method
Values can be:
* 0, s, sample, Sample (default for integer source types)
The nearest value is taken. Pixels are regarded as centered on their respective integer coordinates (no offset from the linear case).
* 1, l, linear, Linear (default for floating point source types)
The values are N-linearly interpolated from an N-dimensional cube of size 2.
* 3, c, cube, cubic, Cubic
The values are interpolated using a local cubic fit to the data. The fit is constrained to match the original data and its derivative at the data points. The second derivative of the fit is not continuous at the data points. Multidimensional datasets are interpolated by the successive-collapse method. (Note that the constraint on the first derivative causes a small amount of ringing around sudden features such as step functions).
* f, fft, fourier, Fourier
The source is Fourier transformed, and the interpolated values are explicitly calculated from the coefficients. The boundary condition option is ignored periodic boundaries are imposed. If you pass in the option fft, and it is a list (ARRAY) ref, then it is a stash for the magnitude and phase of the source FFT. If the list has two elements then they are taken as already computed; otherwise they are calculated and put in the stash.
b, bound, boundary, Boundary
This option is passed unmodified into indexND, which is used as the indexing engine for the interpolation. Some current allowed values are 'extend', 'periodic', 'truncate', and 'mirror' (default is 'truncate').
bad
contains the fill value used for 'truncate' boundary. (default 0)
fft
An array ref whose associated list is used to stash the FFT of the source data, for the FFT method.

one2nd

Converts a one dimensional index piddle to a set of ND coordinates

 @coords=one2nd($a, $indices)

returns an array of piddles containing the ND indexes corresponding to the one dimensional list indices. The indices are assumed to correspond to array CW$a clumped using CWclump(-1). This routine is used in whichND, but is useful on its own occasionally.

 perldl> $a=pdl [[[1,2],[-1,1]], [[0,-3],[3,2]]]; $c=$a->clump(-1)
 perldl> $maxind=maximum_ind($c); p $maxind;
 6
 perldl> print one2nd($a, maximum_ind($c))
 0 1 1
 perldl> p $a->at(0,1,1)
 3

which

  Signature: (mask(n); int [o] inds(m))

Returns piddle of indices of non-zero values.

 $i = which($mask);

returns a pdl with indices for all those elements that are nonzero in the mask. Note that the returned indices will be 1D. If you want to index into the original mask or a similar piddle remember to flatten it before calling index:

  $data = random 5, 5;
  $idx = which $data > 0.5; # $idx is now 1D
  $bigsum = $data->flat->index($idx)->sum;  # flatten before indexing

Compare also where for similar functionality.

If you want to return both the indices of non-zero values and the complement, use the function which_both.

 perldl> $x = sequence(1); p $x
 [0 1 2 3 4 5 6 7 8 9]
 perldl> $indx = which($x>6); p $indx
 [7 8 9]

which does handle bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

which_both

  Signature: (mask(n); int [o] inds(m); int [o]notinds(q))

Returns piddle of indices of non-zero values and their complement

 ($i, $c_i) = which_both($mask);

This works just as which, but the complement of CW$i will be in CW$c_i.

 perldl> $x = sequence(1); p $x
 [0 1 2 3 4 5 6 7 8 9]
 perldl> ($small, $big) = which_both ($x >= 5); p "$small\n $big"
 [5 6 7 8 9]
 [0 1 2 3 4]

which_both does handle bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.

where

Return the values from a piddle for whose indices a condition piddle is nonzero.

 $i = $x->where($x+5 > 0); # $i contains those elements of $x
                           # where mask ($x+5 > 0) is 1
 $i .= -5;  # Set those elements (of $x) to -5. Together, these
            # commands clamp $x to a maximum of -5.

It is also possible to use the same mask for several piddles with the same call:

 ($i,$j,$k) = where($x,$y,$z, $x+5>0);

Note: CW$i is always 1-D, even if CW$x is >1-D.

WARNING: The first argument (the values) and the second argument (the mask) currently have to have the exact same dimensions (or horrible things happen). You *cannot* thread over a smaller mask, for example.

whichND

Returns the coordinates for non-zero values.

For historical reasons the return value is different in list and scalar context. In scalar context, you get back a PDL containing coordinates suitable for use in indexND or range; in list context, the coordinates are broken out into separate PDLs.

 $coords = whichND($mask);

returns a PDL containing the coordinates of the elements that are non-zero in CW$mask, suitable for use in indexND. The 0th dimension contains the full coordinate listing of each point; the 1st dimension lists all the points. For example, if CW$mask has rank 4 and 100 matching elements, then CW$coords has dimension 4x100.

 @coords=whichND($mask);

returns a perl list of piddles containing the coordinates of the elements that are non-zero in CW$mask. Each element corresponds to a particular index dimension. For example, if CW$mask has rank 4 and 100 matching elements, then CW@coords has 4 elements, each of which is a pdl of size 100.

 perldl> $a=sequence(10,10,3,4)
 perldl> ($x, $y, $z, $w)=whichND($a == 203); p $x, $y, $z, $w
 [3] [0] [2] [0]
 perldl> print $a->at(list(cat($x,$y,$z,$w)))
 203

AUTHOR

Copyright (C) Tuomas J. Lukka 1997 (lukka@husc.harvard.edu). Contributions by Christian Soeller (c.soeller@auckland.ac.nz), Karl Glazebrook (kgb@aaoepp.aao.gov.au), and Craig DeForest (deforest@boulder.swri.edu). All rights reserved. There is no warranty. You are allowed to redistribute this software / documentation under certain conditions. For details, see the file COPYING in the PDL distribution. If this file is separated from the PDL distribution, the copyright notice should be included in the file.