use Config; my $ptrsize = $Config{'ptrsize'}; my $int_ptr_type = ($ptrsize == 4) ? $Config{i32type} : $Config{i64type}; # pointer type in C code my $pp_ptr_type = ($ptrsize == 4) ? 'int' : 'longlong'; # pointer type in pp_def Pars sections # User can set this global variable to 1 if he wants # to use the normal plplot order of arguments, not the PP-required # order for functions with OtherPars. $PDL::Graphics::PLplot::standard_order = 0; # Read in options determined by Makefile.PL and written to the OPTIONS! file my $plversion = do 'OPTIONS!'; pp_addpm({At => Top}, <<'EOD'); our $VERSION; BEGIN { $VERSION = '0.61' }; =head1 NAME PDL::Graphics::PLplot - Object-oriented interface from perl/PDL to the PLPLOT plotting library =head1 SYNOPSIS use PDL; use PDL::Graphics::PLplot; my $pl = PDL::Graphics::PLplot->new (DEV => "png", FILE => "test.png"); my $x = sequence(10); my $y = $x**2; $pl->xyplot($x, $y); $pl->close; For more information on PLplot, see http://www.plplot.org/ Also see the test file, F in this distribution for some working examples. =head1 LONG NAMES If you are annoyed by the long constructor call, consider installing the L CPAN package. Using C, the above example becomes use PDL; use aliased 'PDL::Graphics::PLplot'; my $pl = PLplot->new (DEV => "png", FILE => "test.png"); my $x = sequence(10); # etc, as above =head1 DESCRIPTION This is the PDL interface to the PLplot graphics library. It provides a familiar 'perlish' Object Oriented interface as well as access to the low-level PLplot commands from the C-API. =head1 OPTIONS The following options are supported. Most options can be used with any function. A few are only supported on the call to 'new'. =head2 Options used upon creation of a PLplot object (with 'new'): =head3 BACKGROUND Set the color for index 0, the plot background =head3 DEV Set the output device type. To see a list of allowed types, try: PDL::Graphics::PLplot->new(); =for example PDL::Graphics::PLplot->new(DEV => 'png', FILE => 'test.png'); =head3 FILE Set the output file or display. For file output devices, sets the output file name. For graphical displays (like C<'xwin'>) sets the name of the display, eg (C<'hostname.foobar.com:0'>) =for example PDL::Graphics::PLplot->new(DEV => 'png', FILE => 'test.png'); PDL::Graphics::PLplot->new(DEV => 'xwin', FILE => ':0'); =head3 OPTS Set plotting options. See the PLplot documentation for the complete listing of available options. The value of C<'OPTS'> must be a hash reference, whose keys are the names of the options. For instance, to obtain PostScript fonts with the ps output device, use: =for example PDL::Graphics::PLplot->new(DEV => 'ps', OPTS => {drvopt => 'text=1'}); =head3 MEM This option is used in conjunction with C<< DEV => 'mem' >>. This option takes as input a PDL image and allows one to 'decorate' it using PLplot. The 'decorated' PDL image can then be written to an image file using, for example, L. This option may not be available if plplot does not include the 'mem' driver. =for example # read in Earth image and draw an equator. my $pl = PDL::Graphics::PLplot->new (MEM => $earth, DEV => 'mem'); my $x = pdl(-180, 180); my $y = zeroes(2); $pl->xyplot($x, $y, BOX => [-180,180,-90,90], VIEWPORT => [0.0, 1.0, 0.0, 1.0], XBOX => '', YBOX => '', PLOTTYPE => 'LINE'); $pl->close; =head3 FRAMECOLOR Set color index 1, the frame color =head3 JUST A flag used to specify equal scale on the axes. If this is not specified, the default is to scale the axes to fit best on the page. =for example PDL::Graphics::PLplot->new(DEV => 'png', FILE => 'test.png', JUST => 1); =head3 ORIENTATION The orientation of the plot: 0 -- 0 degrees (landscape mode) 1 -- 90 degrees (portrait mode) 2 -- 180 degrees (seascape mode) 3 -- 270 degrees (upside-down mode) Intermediate values (0.2) are acceptable if you are feeling daring. =for example # portrait orientation PDL::Graphics::PLplot->new(DEV => 'png', FILE => 'test.png', ORIENTATION => 1); =head3 PAGESIZE Set the size in pixels of the output page. =for example # PNG 500 by 600 pixels PDL::Graphics::PLplot->new(DEV => 'png', FILE => 'test.png', PAGESIZE => [500,600]); =head3 SUBPAGES Set the number of sub pages in the plot, [$nx, $ny] =for example # PNG 300 by 600 pixels # Two subpages stacked on top of one another. PDL::Graphics::PLplot->new(DEV => 'png', FILE => 'test.png', PAGESIZE => [300,600], SUBPAGES => [1,2]); =head2 Options used after initialization (after 'new') =head3 BOX Set the plotting box in world coordinates. Used to explicitly set the size of the plotting area. =for example my $pl = PDL::Graphics::PLplot->new(DEV => 'png', FILE => 'test.png'); $pl->xyplot ($x, $y, BOX => [0,100,0,200]); =head3 CHARSIZE Set the size of text in multiples of the default size. C<< CHARSIZE => 1.5 >> gives characters 1.5 times the normal size. =head3 COLOR Set the current color for plotting and character drawing. Colors are specified not as color indices but as RGB triples. Some pre-defined triples are included: BLACK GREEN WHEAT BLUE RED AQUAMARINE GREY BLUEVIOLET YELLOW PINK BROWN CYAN TURQUOISE MAGENTA SALMON WHITE ROYALBLUE DEEPSKYBLUE VIOLET STEELBLUE1 DEEPPINK MAGENTA DARKORCHID1 PALEVIOLETRED2 TURQUOISE1 LIGHTSEAGREEN SKYBLUE FORESTGREEN CHARTREUSE3 GOLD2 SIENNA1 CORAL HOTPINK LIGHTCORAL LIGHTPINK1 LIGHTGOLDENROD =for example # These two are equivalent: $pl->xyplot ($x, $y, COLOR => 'YELLOW'); $pl->xyplot ($x, $y, COLOR => [0,255,0]); =head3 LINEWIDTH Set the line width for plotting. Values range from 1 to a device dependent maximum. =head3 LINESTYLE Set the line style for plotting. Pre-defined line styles use values 1 to 8, one being a solid line, 2-8 being various dashed patterns. =head3 MAJTICKSIZE Set the length of major ticks as a fraction of the default setting. One (default) means leave these ticks the normal size. =head3 MINTICKSIZE Set the length of minor ticks (and error bar terminals) as a fraction of the default setting. One (default) means leave these ticks the normal size. =head3 NXSUB The number of minor tick marks between each major tick mark on the X axis. Specify zero (default) to let PLplot compute this automatically. =head3 NYSUB The number of minor tick marks between each major tick mark on the Y axis. Specify zero (default) to let PLplot compute this automatically. =head3 PALETTE Load pre-defined color map 1 color ranges. Currently, values include: RAINBOW -- from Red to Violet through the spectrum REVERSERAINBOW -- Violet through Red GREYSCALE -- from black to white via grey. REVERSEGREYSCALE -- from white to black via grey. GREENRED -- from green to red REDGREEN -- from red to green =for example # Plot x/y points with the z axis in color $pl->xyplot ($x, $y, PALETTE => 'RAINBOW', PLOTTYPE => 'POINTS', COLORMAP => $z); =head3 PLOTTYPE Specify which type of XY plot is desired: LINE -- A line POINTS -- A bunch of symbols LINEPOINTS -- both =head3 SUBPAGE Set which subpage to plot on. Subpages are numbered 1 to N. A zero can be specified meaning 'advance to the next subpage' (just a call to L). =for example my $pl = PDL::Graphics::PLplot->new(DEV => 'png', FILE => 'test.png', SUBPAGES => [1,2]); $pl->xyplot ($x, $y, SUBPAGE => 1); $pl->xyplot ($a, $b, SUBPAGE => 2); =head3 SYMBOL Specify which symbol to use when plotting C<< PLOTTYPE => 'POINTS' >>. A large variety of symbols are available, see: http://plplot.sourceforge.net/examples-data/demo07/x07.*.png, where * is 01 - 17. You are most likely to find good plotting symbols in the 800s: http://plplot.sourceforge.net/examples-data/demo07/x07.06.png =head3 SYMBOLSIZE Specify the size of symbols plotted in multiples of the default size (1). Value are real numbers from 0 to large. =head3 TEXTPOSITION Specify the placement of text. Either relative to border, specified as: [$side, $disp, $pos, $just] Where side = 't', 'b', 'l', or 'r' for top, bottom, left and right disp is the number of character heights out from the edge pos is the position along the edge of the viewport, from 0 to 1. just tells where the reference point of the string is: 0 = left, 1 = right, 0.5 = center. or inside the plot window, specified as: [$x, $y, $dx, $dy, $just] Where x = x coordinate of reference point of string. y = y coordinate of reference point of string. dx Together with dy, this specifies the inclination of the string. The baseline of the string is parallel to a line joining (x, y) to (x+dx, y+dy). dy Together with dx, this specifies the inclination of the string. just Specifies the position of the string relative to its reference point. If just=0, the reference point is at the left and if just=1, it is at the right of the string. Other values of just give intermediate justifications. =for example # Plot text on top of plot $pl->text ("Top label", TEXTPOSITION => ['t', 4.0, 0.5, 0.5]); # Plot text in plotting area $pl->text ("Line label", TEXTPOSITION => [50, 60, 5, 5, 0.5]); =head3 TITLE Add a title on top of a plot. =for example # Plot text on top of plot $pl->xyplot ($x, $y, TITLE => 'X vs. Y'); =head3 UNFILLED_BARS For 'bargraph', if set to true then plot the bars as outlines in the current color and not as filled boxes =for example # Plot text on top of plot $pl->bargraph($labels, $values, UNFILLED_BARS => 1); =head3 VIEWPORT Set the location of the plotting window on the page. Takes a four element array ref specifying: xmin -- The coordinate of the left-hand edge of the viewport. (0 to 1) xmax -- The coordinate of the right-hand edge of the viewport. (0 to 1) ymin -- The coordinate of the bottom edge of the viewport. (0 to 1) ymax -- The coordinate of the top edge of the viewport. (0 to 1) You will need to use this to make color keys or insets. =for example # Make a small plotting window in the lower left of the page $pl->xyplot ($x, $y, VIEWPORT => [0.1, 0.5, 0.1, 0.5]); # Also useful in creating color keys: $pl->xyplot ($x, $y, PALETTE => 'RAINBOW', PLOTTYPE => 'POINTS', COLORMAP => $z); $pl->colorkey ($z, 'v', VIEWPORT => [0.93, 0.96, 0.15, 0.85]); # Plot an inset; first the primary data and then the inset. In this # case, the inset contains a selection of the orignal data $pl->xyplot ($x, $y); $pl->xyplot (where($x, $y, $x < 1.2), VIEWPORT => [0.7, 0.9, 0.6, 0.8]); =head3 XBOX Specify how to label the X axis of the plot as a string of option letters: a: Draws axis, X-axis is horizontal line (y=0), and Y-axis is vertical line (x=0). b: Draws bottom (X) or left (Y) edge of frame. c: Draws top (X) or right (Y) edge of frame. f: Always use fixed point numeric labels. g: Draws a grid at the major tick interval. h: Draws a grid at the minor tick interval. i: Inverts tick marks, so they are drawn outwards, rather than inwards. l: Labels axis logarithmically. This only affects the labels, not the data, and so it is necessary to compute the logarithms of data points before passing them to any of the drawing routines. m: Writes numeric labels at major tick intervals in the unconventional location (above box for X, right of box for Y). n: Writes numeric labels at major tick intervals in the conventional location (below box for X, left of box for Y). s: Enables subticks between major ticks, only valid if t is also specified. t: Draws major ticks. The default is C<'BCNST'> which draws lines around the plot, draws major and minor ticks and labels major ticks. =for example # plot two lines in a box with independent X axes labeled # differently on top and bottom $pl->xyplot($x1, $y, XBOX => 'bnst', # bottom line, bottom numbers, ticks, subticks YBOX => 'bnst'); # left line, left numbers, ticks, subticks $pl->xyplot($x2, $y, XBOX => 'cmst', # top line, top numbers, ticks, subticks YBOX => 'cst', # right line, ticks, subticks BOX => [$x2->minmax, $y->minmax]); =head3 XERRORBAR Used only with L. Draws horizontal error bars at all points (C<$x>, C<$y>) in the plot. Specify a PDL containing the same number of points as C<$x> and C<$y> which specifies the width of the error bar, which will be centered at (C<$x>, C<$y>). =head3 XLAB Specify a label for the X axis. =head3 XTICK Interval (in graph units/world coordinates) between major x axis tick marks. Specify zero (default) to allow PLplot to compute this automatically. =head3 YBOX Specify how to label the Y axis of the plot as a string of option letters. See L. =head3 YERRORBAR Used only for xyplot. Draws vertical error bars at all points (C<$x>, C<$y>) in the plot. Specify a PDL containing the same number of points as C<$x> and C<$y> which specifies the width of the error bar, which will be centered at (C<$x>, C<$y>). =head3 YLAB Specify a label for the Y axis. =head3 YTICK Interval (in graph units/world coordinates) between major y axis tick marks. Specify zero (default) to allow PLplot to compute this automatically. =head3 ZRANGE For L (when C is specified), for L and for L. Normally, the range of the Z variable (color) is taken as C<< $z->minmax >>. If a different range is desired, specify it in C, like so: $pl->shadeplot ($z, $nlevels, PALETTE => 'GREENRED', ZRANGE => [0,100]); or $pl->xyplot ($x, $y, PALETTE => 'RAINBOW', PLOTTYPE => 'POINTS', COLORMAP => $z, ZRANGE => [-90,-20]); $pl->colorkey ($z, 'v', VIEWPORT => [0.93, 0.96, 0.13, 0.85], ZRANGE => [-90,-20]); =head1 METHODS These are the high-level, object oriented methods for PLplot. =head2 new =for ref Create an object representing a plot. =for usage Arguments: none. Supported options: BACKGROUND DEV FILE FRAMECOLOR JUST PAGESIZE SUBPAGES =for example my $pl = PDL::Graphics::PLplot->new(DEV => 'png', FILE => 'test.png'); =head2 setparm =for ref Set options for a plot object. =for usage Arguments: none. Supported options: All options except: BACKGROUND DEV FILE FRAMECOLOR JUST PAGESIZE SUBPAGES (These must be set in call to 'new'.) =for example $pl->setparm (TEXTSIZE => 2); =head2 xyplot =for ref Plot XY lines and/or points. Also supports color scales for points. This function works with bad values. If a bad value is specified for a points plot, it is omitted. If a bad value is specified for a line plot, the bad value makes a gap in the line. This is useful for drawing maps; for example C<$x> and C<$y> can be the continent boundary latitude and longitude. =for usage Arguments: $x, $y Supported options: All options except: BACKGROUND DEV FILE FRAMECOLOR JUST PAGESIZE SUBPAGES (These must be set in call to 'new'.) =for example $pl->xyplot($x, $y, PLOTTYPE => 'POINTS', COLOR => 'BLUEVIOLET', SYMBOL => 1, SYMBOLSIZE => 4); $pl->xyplot($x, $y, PLOTTYPE => 'LINEPOINTS', COLOR => [50,230,30]); $pl->xyplot($x, $y, PALETTE => 'RAINBOW', PLOTTYPE => 'POINTS', COLORMAP => $z); =head2 stripplots =for ref Plot a set of strip plots with a common X axis, but with different Y axes. Looks like a stack of long, thin XY plots, all line up on the same X axis. =for usage Arguments: $xs -- 1D PDL with common X axis values, length = N $ys -- reference to a list of 1D PDLs with Y-axis values, length = N or 2D PDL with N x M elements -- OR -- $xs -- reference to a list of 1D PDLs with X-axis values $ys -- reference to a list of 1D PDLs with Y-axis values %opts -- Options hash Supported options: All options except: BACKGROUND DEV FILE FRAMECOLOR JUST PAGESIZE SUBPAGES (These must be set in call to 'new'.) =for example my $x = sequence(20); my $y1 = $x**2; my $y2 = sqrt($x); my $y3 = $x**3; my $y4 = sin(($x/20) * 2 * $pi); $ys = cat($y1, $y2, $y3, $y4); $pl->stripplots($x, $ys, PLOTTYPE => 'LINE', TITLE => 'functions', YLAB => ['x**2', 'sqrt(x)', 'x**3', 'sin(x/20*2pi)'], COLOR => ['GREEN', 'DEEPSKYBLUE', 'DARKORCHID1', 'DEEPPINK'], XLAB => 'X label'); # Equivalent to above: $pl->stripplots($x, [$y1, $y2, $y3, $y4], PLOTTYPE => 'LINE', TITLE => 'functions', YLAB => ['x**2', 'sqrt(x)', 'x**3', 'sin(x/20*2pi)'], COLOR => ['GREEN', 'DEEPSKYBLUE', 'DARKORCHID1', 'DEEPPINK'], XLAB => 'X label'); # Here's something a bit different. Notice that different xs have # different lengths. $x1 = sequence(20); $y1 = $x1**2; $x2 = sequence(18); $y2 = sqrt($x2); $x3 = sequence(24); $y3 = $x3**3; my $x4 = sequence(27); $a = ($x4/20) * 2 * $pi; my $y4 = sin($a); $xs = [$x1, $x2, $x3, $x4]; $ys = [$y1, $y2, $y3, $y4]; $pl->stripplots($xs, $ys, PLOTTYPE => 'LINE', TITLE => 'functions', YLAB => ['x**2', 'sqrt(x)', 'x**3', 'sin(x/20*2pi)'], COLOR => ['GREEN', 'DEEPSKYBLUE', 'DARKORCHID1', 'DEEPPINK'], XLAB => 'X label'); In addition, COLOR may be specified as a reference to a list of colors. If this is done, the colors are applied separately to each plot. Also, the options Y_BASE and Y_GUTTER can be specified. Y_BASE gives the Y offset of the bottom of the lowest plot (0-1, specified like a VIEWPORT, defaults to 0.1) and Y_GUTTER gives the gap between the graphs (0-1, default = 0.02). =head2 colorkey =for ref Plot a color key showing which color represents which value =for usage Arguments: $range : A PDL which tells the range of the color values $orientation : 'v' for vertical color key, 'h' for horizontal Supported options: All options except: BACKGROUND DEV FILE FRAMECOLOR JUST PAGESIZE SUBPAGES (These must be set in call to 'new'.) =for example # Plot X vs. Y with Z shown by the color. Then plot # vertical key to the right of the original plot. $pl->xyplot ($x, $y, PALETTE => 'RAINBOW', PLOTTYPE => 'POINTS', COLORMAP => $z); $pl->colorkey ($z, 'v', VIEWPORT => [0.93, 0.96, 0.15, 0.85]); =head2 shadeplot =for ref Create a shaded contour plot of 2D PDL 'z' with 'nsteps' contour levels. Linear scaling is used to map the coordinates of Z(X, Y) to world coordinates via the L option. =for usage Arguments: $z : A 2D PDL which contains surface values at each XY coordinate. $nsteps : The number of contour levels requested for the plot. Supported options: All options except: BACKGROUND DEV FILE FRAMECOLOR JUST PAGESIZE SUBPAGES (These must be set in call to 'new'.) =for example # vertical key to the right of the original plot. # The BOX must be specified to give real coordinate values to the $z array. $pl->shadeplot ($z, $nsteps, BOX => [-1, 1, -1, 1], PALETTE => 'RAINBOW', ZRANGE => [0,100]); $pl->colorkey ($z, 'v', VIEWPORT => [0.93, 0.96, 0.15, 0.85], ZRANGE => [0,100]); =head2 histogram =for ref Create a histogram of a 1-D variable. =for usage Arguments: $x : A 1D PDL $nbins : The number of bins to use in the histogram. Supported options: All options except: BACKGROUND DEV FILE FRAMECOLOR JUST PAGESIZE SUBPAGES (These must be set in call to 'new'.) =for example $pl->histogram ($x, $nbins, BOX => [$min, $max, 0, 100]); =head2 bargraph =for ref Simple utility to plot a bar chart with labels on the X axis. The usual options can be specified, plus one other: MAXBARLABELS specifies the maximum number of labels to allow on the X axis. The default is 20. If this value is exceeded, then every other label is plotted. If twice MAXBARLABELS is exceeded, then only every third label is printed, and so on. if UNFILLED_BARS is set to true, then plot the bars as outlines and not as filled rectangles. =for usage Arguments: $labels -- A reference to a perl list of strings. $values -- A PDL of values to be plotted. Supported options: All options except: BACKGROUND DEV FILE FRAMECOLOR JUST PAGESIZE SUBPAGES (These must be set in call to 'new'.) =for example $labels = ['one', 'two', 'three']; $values = pdl(1, 2, 3); # Note if TEXTPOSITION is specified, it must be in 4 argument mode (border mode): # [$side, $disp, $pos, $just] # # Where side = 't', 'b', 'l', or 'r' for top, bottom, left and right # 'tv', 'bv', 'lv' or 'rv' for top, bottom, left or right perpendicular to the axis. # # disp is the number of character heights out from the edge # pos is the position along the edge of the viewport, from 0 to 1. # just tells where the reference point of the string is: 0 = left, 1 = right, 0.5 = center. # # The '$pos' entry will be ignored (computed by the bargraph routine) $pl->bargraph($labels, $values, MAXBARLABELS => 30, TEXTPOSITION => ['bv', 0.5, 1.0, 1.0]); =head2 text =for ref Write text on a plot. Text can either be written with respect to the borders or at an arbitrary location and angle (see the L entry). =for usage Arguments: $t : The text. Supported options: All options except: BACKGROUND DEV FILE FRAMECOLOR JUST PAGESIZE SUBPAGES (These must be set in call to 'new'.) =for example $pl->text("Count", COLOR => 'PINK', TEXTPOSITION => ['t', 3, 0.5, 0.5]); # top, 3 units out, string ref. pt in # center of string, middle of axis =head2 close =for ref Close a PLplot object, writing out the file and cleaning up. =for usage Arguments: None Returns: Nothing This closing of the PLplot object can be done explicitly though the 'close' method. Alternatively, a DESTROY block does an automatic close whenever the PLplot object passes out of scope. =for example $pl->close; =cut # pull in low level interface use vars qw(%_constants %_actions); # Colors (from rgb.txt) are stored as RGB triples # with each value from 0-255 sub cc2t { [map {hex} split ' ', shift] } %_constants = ( BLACK => [ 0, 0, 0], RED => [240, 50, 50], YELLOW => [255,255, 0], GREEN => [ 0,255, 0], AQUAMARINE => [127,255,212], PINK => [255,192,203], WHEAT => [245,222,179], GREY => [190,190,190], BROWN => [165, 42, 42], BLUE => [ 0, 0,255], BLUEVIOLET => [138, 43,226], CYAN => [ 0,255,255], TURQUOISE => [ 64,224,208], MAGENTA => [255, 0,255], SALMON => [250,128,114], WHITE => [255,255,255], ROYALBLUE => cc2t('2B 60 DE'), DEEPSKYBLUE => cc2t('3B B9 FF'), VIOLET => cc2t('8D 38 C9'), STEELBLUE1 => cc2t('5C B3 FF'), DEEPPINK => cc2t('F5 28 87'), MAGENTA => cc2t('FF 00 FF'), DARKORCHID1 => cc2t('B0 41 FF'), PALEVIOLETRED2 => cc2t('E5 6E 94'), TURQUOISE1 => cc2t('52 F3 FF'), LIGHTSEAGREEN => cc2t('3E A9 9F'), SKYBLUE => cc2t('66 98 FF'), FORESTGREEN => cc2t('4E 92 58'), CHARTREUSE3 => cc2t('6C C4 17'), GOLD2 => cc2t('EA C1 17'), SIENNA1 => cc2t('F8 74 31'), CORAL => cc2t('F7 65 41'), HOTPINK => cc2t('F6 60 AB'), LIGHTCORAL => cc2t('E7 74 71'), LIGHTPINK1 => cc2t('F9 A7 B0'), LIGHTGOLDENROD => cc2t('EC D8 72'), ); # a hash of subroutines to invoke when certain keywords are specified # These are called with arg(0) = $self (the plot object) # and arg(1) = value specified for keyword %_actions = ( # Set color for index 0, the plot background BACKGROUND => sub { my $self = shift; my $color = _color(shift); $self->{COLORS}[0] = $color; plscolbg (@$color); }, # set plotting box in world coordinates BOX => sub { my $self = shift; my $box = shift; die "Box must be a ref to a four element array" unless (ref($box) =~ /ARRAY/ and @$box == 4); $self->{BOX} = $box; }, CHARSIZE => sub { my $self = shift; $self->{CHARSIZE} = $_[0]; plschr (0, $_[0]) }, # 0 - N COLOR => # maintain color map, set to specified rgb triple sub { my $self = shift; my $color = _color(shift); # init. $self->{COLORS} = [] unless exists($self->{COLORS}); my @idx = @{$self->{COLORS}}; # map of color index (0-15) to RGB triples my $found = 0; for (my $i=2;$i<@idx;$i++) { # map entries 0 and 1 are reserved for BACKGROUND and FRAMECOLOR if (_coloreq ($color, $idx[$i])) { $self->{CURRENT_COLOR_IDX} = $i; $found = 1; plscol0 ($self->{CURRENT_COLOR_IDX}, @$color); } } return if ($found); die "Too many colors used! (max 15)" if (@{$self->{COLORS}} > 14); # add this color as index 2 or greater (entries 0 and 1 reserved) my $idx = (@{$self->{COLORS}} > 1) ? @{$self->{COLORS}} : 2; $self->{COLORS}[$idx] = $color; $self->{CURRENT_COLOR_IDX} = $idx; plscol0 ($self->{CURRENT_COLOR_IDX}, @$color); }, # set output device type DEV => sub { my $self = shift; my $dev = shift; $self->{DEV} = $dev; plsdev ($dev) }, # this must be specified with call to new! # set PDL to plot into (alternative to specifying DEV) MEM => sub { my $self = shift; my $pdl = shift; my $x = $pdl->getdim(1); my $y = $pdl->getdim(2); plsmem ($x, $y, $pdl); }, # set output file FILE => sub { plsfnam ($_[1]) }, # this must be specified with call to new! # set color for index 1, the plot frame and text FRAMECOLOR => # set color index 1, the frame color sub { my $self = shift; my $color = _color(shift); $self->{COLORS}[1] = $color; plscol0 (1, @$color); }, # Set flag for equal scale axes JUST => sub { my $self = shift; my $just = shift; die "JUST must be 0 or 1 (defaults to 0)" unless ($just == 0 or $just == 1); $self->{JUST} = $just; }, LINEWIDTH => sub { my $self = shift; my $wid = shift; die "LINEWIDTH must range from 0 to LARGE8" unless ($wid >= 0); $self->{LINEWIDTH} = $wid; }, LINESTYLE => sub { my $self = shift; my $sty = shift; die "LINESTYLE must range from 1 to 8" unless ($sty >= 1 and $sty <= 8); $self->{LINESTYLE} = $sty; }, MAJTICKSIZE => sub { my $self = shift; my $val = shift; die "MAJTICKSIZE must be greater than or equal to zero" unless ($val >= 0); plsmaj (0, $val); }, MINTICKSIZE => sub { my $self = shift; my $val = shift; die "MINTICKSIZE must be greater than or equal to zero" unless ($val >= 0); plsmin (0, $val); }, NXSUB => sub { my $self = shift; my $val = shift; die "NXSUB must be an integer greater than or equal to zero" unless ($val >= 0 and int($val) == $val); $self->{NXSUB} = $val; }, NYSUB => sub { my $self = shift; my $val = shift; die "NYSUB must be an integer greater than or equal to zero" unless ($val >= 0 and int($val) == $val); $self->{NYSUB} = $val; }, # set driver options, example for ps driver, {text => 1} is accepted OPTS => sub { my $self = shift; my $opts = shift; foreach my $opt (keys %$opts) { plsetopt ($opt, $$opts{$opt}); } }, # set driver options, example for ps driver, {text => 1} is accepted ORIENTATION => sub { my $self = shift; my $orient = shift; die "Orientation must be between 0 and 4" unless ($orient >= 0 and $orient <= 4); $self->{ORIENTATION} = $orient; }, PAGESIZE => # set plot size in mm. Only useful in call to 'new' sub { my $self = shift; my $dims = shift; die "plot size must be a 2 element array ref: X size in pixels, Y size in pixels" if ((ref($dims) !~ /ARRAY/) || @$dims != 2); $self->{PAGESIZE} = $dims; }, PALETTE => # load some pre-done color map 1 setups sub { my $self = shift; my $pal = shift; my %legal = (REVERSERAINBOW => 1, REVERSEGREYSCALE => 1, REDGREEN => 1, RAINBOW => 1, GREYSCALE => 1, GREENRED => 1); if ($legal{$pal}) { $self->{PALETTE} = $pal; if ($pal eq 'RAINBOW') { plscmap1l (0, PDL->new(0,1), PDL->new(0,300), PDL->new(0.5, 0.5), PDL->new(1,1), PDL->new(0,0)); } elsif ($pal eq 'REVERSERAINBOW') { plscmap1l (0, PDL->new(0,1), PDL->new(270,-30), PDL->new(0.5, 0.5), PDL->new(1,1), PDL->new(0,0)); } elsif ($pal eq 'GREYSCALE') { plscmap1l (0, PDL->new(0,1), PDL->new(0,0), PDL->new(0,1), PDL->new(0,0), PDL->new(0,0)); } elsif ($pal eq 'REVERSEGREYSCALE') { plscmap1l (0, PDL->new(0,1), PDL->new(0,0), PDL->new(1,0), PDL->new(0,0), PDL->new(0,0)); } elsif ($pal eq 'GREENRED') { plscmap1l (0, PDL->new(0,1), PDL->new(120,0), PDL->new(0.5, 0.5), PDL->new(1,1), PDL->new(1,1)); } elsif ($pal eq 'REDGREEN') { plscmap1l (0, PDL->new(0,1), PDL->new(0,120), PDL->new(0.5, 0.5), PDL->new(1,1), PDL->new(1,1)); } } else { die "Illegal palette name. Legal names are: " . join (" ", keys %legal); } }, PLOTTYPE => # specify plot type (LINE, POINTS, LINEPOINTS) sub { my $self = shift; my $val = shift; my %legal = (LINE => 1, POINTS => 1, LINEPOINTS => 1); if ($legal{$val}) { $self->{PLOTTYPE} = $val; } else { die "Illegal plot type. Legal options are: " . join (" ", keys %legal); } }, SUBPAGE => # specify which subpage to plot on 1-N or 0 (meaning 'next') sub { my $self = shift; my $val = shift; my $err = "SUBPAGE = \$npage where \$npage = 1-N or 0 (for 'next subpage')"; if ($val >= 0) { $self->{SUBPAGE} = $val; } else { die $err; } }, SUBPAGES => # specify number of sub pages [nx, ny] sub { my $self = shift; my $val = shift; my $err = "SUBPAGES = [\$nx, \$ny] where \$nx and \$ny are between 1 and 127"; if (ref($val) =~ /ARRAY/ and @$val == 2) { my ($nx, $ny) = @$val; if ($nx > 0 and $nx < 128 and $ny > 0 and $ny < 128) { $self->{SUBPAGES} = [$nx, $ny]; } else { die $err; } } else { die $err; } }, SYMBOL => # specify type of symbol to plot sub { my $self = shift; my $val = shift; if ($val >= 0 && $val < 3000) { $self->{SYMBOL} = $val; } else { die "Illegal symbol number. Legal symbols are between 0 and 3000"; } }, SYMBOLSIZE => sub { my ($self, $size) = @_; die "symbol size must be a real number from 0 to (large)" unless ($size >= 0); $self->{SYMBOLSIZE} = $size; }, TEXTPOSITION => # specify placement of text. Either relative to border, specified as: # [$side, $disp, $pos, $just] # or # inside plot window, specified as: # [$x, $y, $dx, $dy, $just] (see POD doc for details) sub { my $self = shift; my $val = shift; die "TEXTPOSITION value must be an array ref with either: [$side, $disp, $pos, $just] or [$x, $y, $dx, $dy, $just]" unless ((ref($val) =~ /ARRAY/) and ((@$val == 4) || (@$val == 5))); if (@$val == 4) { $self->{TEXTMODE} = 'border'; } else { $self->{TEXTMODE} = 'plot'; } $self->{TEXTPOSITION} = $val; }, # draw a title for the graph TITLE => sub { my $self = shift; my $text = shift; $self->{TITLE} = $text; }, # Specify outline bars for bargraph UNFILLED_BARS => sub { my $self = shift; my $val = shift; $self->{UNFILLED_BARS} = $val; }, # set the location of the plotting window on the page VIEWPORT => sub { my $self = shift; my $vp = shift; die "Viewport must be a ref to a four element array" unless (ref($vp) =~ /ARRAY/ and @$vp == 4); $self->{VIEWPORT} = $vp; }, XBOX => # set X axis label options. See pod for definitions. sub { my $self = shift; my $opts = lc shift; my @opts = split '', $opts; map { 'abcdfghilmnst' =~ /$_/i || die "Illegal option $_. Only abcdfghilmnst permitted" } @opts; $self->{XBOX} = $opts; }, # draw an X axis label for the graph XLAB => sub { my $self = shift; my $text = shift; $self->{XLAB} = $text; }, XTICK => sub { my $self = shift; my $val = shift; die "XTICK must be greater than or equal to zero" unless ($val >= 0); $self->{XTICK} = $val; }, YBOX => # set Y axis label options. See pod for definitions. sub { my $self = shift; my $opts = shift; my @opts = split '', $opts; map { 'abcfghilmnstv' =~ /$_/i || die "Illegal option $_. Only abcfghilmnstv permitted" } @opts; $self->{YBOX} = $opts; }, # draw an Y axis label for the graph YLAB => sub { my $self = shift; my $text = shift; $self->{YLAB} = $text; }, YTICK => sub { my $self = shift; my $val = shift; die "YTICK must be greater than or equal to zero" unless ($val >= 0); $self->{YTICK} = $val; }, ZRANGE => sub { my $self = shift; my $val = shift; die "ZRANGE must be a perl array ref with min and max Z values" unless (ref($val) =~ /ARRAY/ && @$val == 2); $self->{ZRANGE} = $val; }, ); # ## Internal utility routines # # handle color as string in _constants hash or [r,g,b] triple # Input: either color name or [r,g,b] array ref # Output: [r,g,b] array ref or exception sub _color { my $c = shift; if (ref($c) =~ /ARRAY/) { return $c; } elsif ($c = $_constants{$c}) { return $c; } else { die "Color $c not defined"; } } # return 1 if input [r,g,b] triples are equal. sub _coloreq { my ($a, $b) = @_; for (my $i=0;$i<3;$i++) { return 0 if ($$a[$i] != $$b[$i]); } return 1; } # Initialize plotting window given the world coordinate box and # a 'justify' flag (for equal axis scales). sub _setwindow { my $self = shift; # choose correct subwindow pladv ($self->{SUBPAGE}) if (exists ($self->{SUBPAGE})); delete ($self->{SUBPAGE}); # get rid of SUBPAGE so future plots will stay on same # page unless user asks for specific page my $box = $self->{BOX} || [0,1,0,1]; # default window sub MAX { ($_[0] > $_[1]) ? $_[0] : $_[1]; } # get subpage offsets from page left/bottom of image my ($spxmin, $spxmax, $spymin, $spymax) = (PDL->new(0),PDL->new(0),PDL->new(0),PDL->new(0)); plgspa($spxmin, $spxmax, $spymin, $spymax); $spxmin = $spxmin->at(0); $spxmax = $spxmax->at(0); $spymin = $spymin->at(0); $spymax = $spymax->at(0); my $xsize = $spxmax - $spxmin; my $ysize = $spymax - $spymin; my @vp = @{$self->{VIEWPORT}}; # view port xmin, xmax, ymin, ymax in fraction of image size # if JUSTify is zero, set to the user specified (or default) VIEWPORT if ($self->{JUST} == 0) { plvpor(@vp); # compute viewport to allow the same scales for both axes } else { my $p_def = PDL->new(0); my $p_ht = PDL->new(0); plgchr ($p_def, $p_ht); $p_def = $p_def->at(0); my $lb = 8.0 * $p_def; my $rb = 5.0 * $p_def; my $tb = 5.0 * $p_def; my $bb = 5.0 * $p_def; my $dx = $$box[1] - $$box[0]; my $dy = $$box[3] - $$box[2]; my $xscale = $dx / ($xsize - $lb - $rb); my $yscale = $dy / ($ysize - $tb - $bb); my $scale = MAX($xscale, $yscale); my $vpxmin = MAX($lb, 0.5 * ($xsize - $dx / $scale)); my $vpxmax = $vpxmin + ($dx / $scale); my $vpymin = MAX($bb, 0.5 * ($ysize - $dy / $scale)); my $vpymax = $vpymin + ($dy / $scale); plsvpa($vpxmin, $vpxmax, $vpymin, $vpymax); $self->{VIEWPORT} = [$vpxmin/$xsize, $vpxmax/$xsize, $vpymin/$ysize, $vpymax/$ysize]; } # set up world coords in window plwind (@$box); } # Add title and axis labels. sub _drawlabels { my $self = shift; plcol0 (1); # set to frame color plmtex (2.5, 0.5, 0.5, 't', $self->{TITLE}) if ($self->{TITLE}); plmtex (3.0, 0.5, 0.5, 'b', $self->{XLAB}) if ($self->{XLAB}); plmtex (3.5, 0.5, 0.5, 'l', $self->{YLAB}) if ($self->{YLAB}); plcol0 ($self->{CURRENT_COLOR_IDX}); # set back } # ## user-visible routines # # Pool of PLplot stream numbers. One of these stream numbers is taken when 'new' is called # and when the corresponding 'close' is called, it is returned to the pool. The pool is # just a queue: 'new' shifts stream numbers from the top of the queue, 'close' pushes them # back on the bottom of the queue. my @plplot_stream_pool = (0..99); # This routine starts out a plot. Generally one specifies # DEV and FILE (device and output file name) as options. sub new { my $type = shift; my $self = {}; # set up object $self->{PLOTTYPE} = 'LINE'; # $self->{CURRENT_COLOR_IDX} = 1; $self->{COLORS} = []; bless $self, $type; # set stream number first $self->{STREAMNUMBER} = shift @plplot_stream_pool; die "No more PLplot streams left, too many open PLplot objects!" if (!defined($self->{STREAMNUMBER})); plsstrm($self->{STREAMNUMBER}); # set background and frame color first $self->setparm(BACKGROUND => 'WHITE', FRAMECOLOR => 'BLACK'); # set defaults, allow input options to override my %opts = ( COLOR => 'BLACK', XBOX => 'BCNST', YBOX => 'BCNST', JUST => 0, SUBPAGES => [1,1], VIEWPORT => [0.1, 0.87, 0.13, 0.82], SUBPAGE => 0, PAGESIZE => [600, 500], LINESTYLE => 1, LINEWIDTH => 0, SYMBOL => 751, # a small square NXSUB => 0, NYSUB => 0, ORIENTATION=> 0, XTICK => 0, YTICK => 0, CHARSIZE => 1, @_); # apply options $self->setparm(%opts); # Do initial setup plspage (0, 0, @{$self->{PAGESIZE}}, 0, 0) if (defined($self->{PAGESIZE})); plssub (@{$self->{SUBPAGES}}); plfontld (1); # extented symbol pages plscmap0n (16); # set up color map 0 to 16 colors. Is this needed? plscmap1n (128); # set map 1 to 128 colors (should work for devices with 256 colors) plinit (); # set page orientation plsdiori ($self->{ORIENTATION}); # set up plotting box $self->_setwindow; return $self; } # set parameters. Called from user directly or from other routines. sub setparm { my $self = shift; my %opts = @_; # Set PLplot to right output stream plsstrm($self->{STREAMNUMBER}); # apply all options OPTION: foreach my $o (keys %opts) { unless (exists($_actions{$o})) { warn "Illegal option $o, ignoring"; next OPTION; } &{$_actions{$o}}($self, $opts{$o}); } } # handle 2D plots sub xyplot { my $self = shift; my $x = shift; my $y = shift; my %opts = @_; # Set PLplot to right output stream plsstrm($self->{STREAMNUMBER}); # only process COLORMAP entries once my $z = $opts{COLORMAP}; delete ($opts{COLORMAP}); # handle ERRORBAR options my $xeb = $opts{XERRORBAR}; my $yeb = $opts{YERRORBAR}; delete ($opts{XERRORBAR}); delete ($opts{YERRORBAR}); # apply options $self->setparm(%opts); unless (exists($self->{BOX})) { $self->{BOX} = [$x->minmax, $y->minmax]; } # set up viewport, subpage, world coordinates $self->_setwindow; # draw labels $self->_drawlabels; # plot box plcol0 (1); # set to frame color plbox ($self->{XTICK}, $self->{NXSUB}, $self->{YTICK}, $self->{NYSUB}, $self->{XBOX}, $self->{YBOX}); # !!! note out of order call # set the color according to the color specified in the object # (we don't do this as an option, because then the frame might # get the color requested for the line/points plcol0 ($self->{CURRENT_COLOR_IDX}); # set line style for plot only (not box) pllsty ($self->{LINESTYLE}); # set line width for plot only (not box) plwid ($self->{LINEWIDTH}); # Plot lines if requested if ($self->{PLOTTYPE} =~ /LINE/) { plline ($x, $y); } # set line width back plwid (0); # plot points if requested if ($self->{PLOTTYPE} =~ /POINTS/) { my $c = $self->{SYMBOL}; unless (defined($c)) { # the default for $c is a PDL of ones with shape # equal to $x with the first dimension removed my $z = PDL->zeroes($x->nelem); $c = PDL->ones($z->zcover) unless defined($c); } plssym (0, $self->{SYMBOLSIZE}) if (defined($self->{SYMBOLSIZE})); if (defined($z)) { # if a color range plot requested my ($min, $max) = exists ($self->{ZRANGE}) ? @{$self->{ZRANGE}} : $z->minmax; plcolorpoints ($x, $y, $z, $c, $min, $max); } else { plsym ($x, $y, $c); } } # Plot error bars, if requested if (defined($xeb)) { # horizontal (X) error bars plerrx ($x->nelem, $x - $xeb/2, $x + $xeb/2, $y); } if (defined($yeb)) { # vertical (Y) error bars plerry ($y->nelem, $x, $y - $yeb/2, $y + $yeb/2); } # Flush the PLplot stream. plflush(); } # Handle sets of 2D strip plots sharing one X axis. Input is # $self -- PLplot object with existing options # $xs -- Ref to list of 1D PDLs with X values # $ys -- Ref to list of 1D PDLs with Y values # or a 2D PDL # %opts -- Options values sub stripplots { my $self = shift; my $xs = shift; my $yargs = shift; my %opts = @_; # NYTICK => number of y axis ticks my $nytick = $opts{NYTICK} || 2; delete ($opts{NYTICK}); # only process COLORMAP entries once my $zs = $opts{COLORMAP}; delete ($opts{COLORMAP}); # handle XLAB, YLAB and TITLE options my $title = $opts{TITLE} || ''; my $xlab = $opts{XLAB} || ''; my @ylabs = defined($opts{YLAB}) && (ref($opts{YLAB}) =~ /ARRAY/) ? @{$opts{YLAB}} : (); delete @opts{qw(TITLE XLAB YLAB)}; # Ensure we're dealing with an array reference my $ys; if (ref ($yargs) eq 'ARRAY') { $ys = $yargs; } elsif (ref ($yargs) =~ /PDL/) { $ys = [dog $yargs]; } else { barf("stripplots requires that its second argument be either a 2D piddle or\na reference to a list of 1D piddles, but you provided neither."); } # This doesn't work because $xs can be an anonymous array, too # # Let's be sure the user sent us what we expected: # foreach (@$ys) { # barf ("stripplots needs to have piddles for its y arguments!") # unless (ref =~ /PDL/); # barf("stripplots requires that the x and y dimensions agree!") # unless ($_->nelem == $xs->nelem); # } my $nplots = @$ys; # Use list of colors, or single color. If COLOR not specified, default to BLACK for each graph my @colors = (defined ($opts{COLOR}) && ref($opts{COLOR}) =~ /ARRAY/) ? @{$opts{COLOR}} : defined ($opts{COLOR}) ? ($opts{COLOR}) x $nplots : ('BLACK') x $nplots; delete @opts{qw(COLOR)}; my $y_base = defined($opts{Y_BASE}) ? $opts{Y_BASE} : 0.1; # Y offset to start bottom plot my $y_gutter = defined($opts{Y_GUTTER}) ? $opts{Y_GUTTER} : 0.02; # Y gap between plots delete @opts{qw(Y_BASE Y_GUTTER)}; # apply options $self->setparm(%opts); my ($xmin, $xmax); if (ref ($xs) =~ /PDL/) { ($xmin, $xmax) = $xs->minmax; } else { $xmin = pdl(map { $_->min } @$xs)->min; $xmax = pdl(map { $_->max } @$xs)->max; } SUBPAGE: for (my $subpage=0;$subpage<$nplots;$subpage++) { my $y = $ys->[$subpage]; my $x = ref ($xs) =~ /PDL/ ? $xs : $xs->[$subpage]; my $mask = $y->isgood; $y = $y->where($mask); $x = $x->where($mask); my $z = $zs->slice(":,($subpage)")->where($mask) if (defined($zs)); my $yeb = $yebs->slice(":,($subpage)")->where($mask) if (defined($yebs)); my $ylab = $ylabs[$subpage]; my $bottomplot = ($subpage == 0); my $topplot = ($subpage == $nplots-1); my $xbox = 'bc'; $xbox = 'cstnb' if ($bottomplot); my $box = $opts{BOX}; my $yrange = defined($box) ? $$box[3] - $$box[2] : $y->max - $y->min; my $del = $yrange ? $yrange * 0.05 : 1; my @ybounds = ($y->min - $del, $y->max + $del); my $ytick = ($yrange/$nytick); my @COLORMAP = (COLORMAP => $z) if defined($z); $self->xyplot($x, $y, COLOR => $colors[$subpage], BOX => defined($box) ? $box : [$xmin, $xmax, @ybounds], XBOX => $xbox, YBOX => 'BCNT', YTICK => $ytick, MAJTICKSIZE => 0.6, CHARSIZE => 0.4, @COLORMAP, VIEWPORT => [ 0.15, 0.9, $y_base + ($subpage * (0.8/$nplots)), $y_base - $y_gutter + (($subpage+1) * (0.8/$nplots)), ], ); $self->text($ylab, TEXTPOSITION => ['L', 4, 0.5, 0.5], COLOR => 'BLACK', CHARSIZE => 0.6) if (defined($ylab)); $self->text($xlab, TEXTPOSITION => ['B', 3, 0.5, 0.5], COLOR => 'BLACK', CHARSIZE => 0.6) if ($xlab && $bottomplot); $self->text($title, TEXTPOSITION => ['T', 2, 0.5, 0.5], COLOR => 'BLACK', CHARSIZE => 1.3) if ($title && $topplot); } } # Draw a color key or wedge showing the scale of map1 colors sub colorkey { my $self = shift; my $var = shift; my $orientation = shift; # 'v' (for vertical) or 'h' (for horizontal) my %opts = @_; # Set PLplot to right output stream plsstrm($self->{STREAMNUMBER}); # apply options $self->setparm(%opts); # set up viewport, subpage, world coordinates $self->_setwindow; # draw labels $self->_drawlabels; # Allow user to set X, Y box type for color key scale. D. Hunt 1/7/2009 my $xbox = exists($self->{XBOX}) ? $self->{XBOX} : 'TM'; my $ybox = exists($self->{YBOX}) ? $self->{YBOX} : 'TM'; my @box; plcol0 (1); # set to frame color my ($min, $max) = exists ($self->{ZRANGE}) ? @{$self->{ZRANGE}} : $var->minmax; # plot box if ($orientation eq 'v') { # set world coordinates based on input variable @box = (0, 1, $min, $max); plwind (@box); plbox (0, 0, 0, 0, '', $ybox); # !!! note out of order call } elsif ($orientation eq 'h') { @box = ($min, $max, 0, 1); plwind (@box); plbox (0, 0, 0, 0, $xbox, ''); # !!! note out of order call } else { die "Illegal orientation value: $orientation. Should be 'v' (vertical) or 'h' (horizontal)"; } # restore color setting plcol0 ($self->{CURRENT_COLOR_IDX}); # This is the number of colors shown in the color wedge. Make # this smaller for gif images as these are limited to 256 colors total. # D. Hunt 8/9/2006 my $ncols = ($self->{DEV} =~ /gif/) ? 32 : 128; if ($orientation eq 'v') { my $yinc = ($box[3] - $box[2])/$ncols; my $y0 = $box[2]; for (my $i=0;$i<$ncols;$i++) { $y0 = $box[2] + ($i * $yinc); my $y1 = $y0 + $yinc; PDL::Graphics::PLplot::plcol1($i/$ncols); # Instead of using plfill (which is not supported on some devices) # use multiple calls to plline to color in the space. D. Hunt 8/9/2006 foreach my $inc (0..9) { my $frac = $inc * 0.1; my $y = $y0 + (($y1 - $y0) * $frac); PDL::Graphics::PLplot::plline (PDL->new(0,1), PDL->new($y,$y)); } } } else { my $xinc = ($box[1] - $box[0])/$ncols; my $x0 = $box[0]; for (my $i=0;$i<$ncols;$i++) { $x0 = $box[0] + ($i * $xinc); my $x1 = $x0 + $xinc; PDL::Graphics::PLplot::plcol1($i/$ncols); # Instead of using plfill (which is not supported on some devices) # use multiple calls to plline to color in the space. D. Hunt 8/9/2006 foreach my $inc (0..9) { my $frac = $inc * 0.1; my $x = $x0 + (($x1 - $x0) * $frac); PDL::Graphics::PLplot::plline (PDL->new($x,$x), PDL->new(0,1)); } } } # Flush the PLplot stream. plflush(); } # handle shade plots of gridded (2D) data sub shadeplot { my $self = shift; my $z = shift; my $nsteps = shift; my %opts = @_; # Set PLplot to right output stream plsstrm($self->{STREAMNUMBER}); # apply options $self->setparm(%opts); my ($nx, $ny) = $z->dims; unless (exists($self->{BOX})) { $self->{BOX} = [0, $nx, 0, $ny]; } # set up plotting box $self->_setwindow; # draw labels $self->_drawlabels; # plot box plcol0 (1); # set to frame color plbox ($self->{XTICK}, $self->{NXSUB}, $self->{YTICK}, $self->{NYSUB}, $self->{XBOX}, $self->{YBOX}); # !!! note out of order call my ($min, $max) = exists ($self->{ZRANGE}) ? @{$self->{ZRANGE}} : $z->minmax; my $clevel = ((PDL->sequence($nsteps)*(($max - $min)/($nsteps-1))) + $min); # may add as options later. Now use constants my $fill_width = 2; my $cont_color = 0; my $cont_width = 0; my $rectangular = 1; # only false for non-linear coord mapping (not done yet in perl) # map X coords linearly to X range, Y coords linearly to Y range my $xmap = ((PDL->sequence($nx)*(($self->{BOX}[1] - $self->{BOX}[0])/($nx - 1))) + $self->{BOX}[0]); my $ymap = ((PDL->sequence($ny)*(($self->{BOX}[3] - $self->{BOX}[2])/($ny - 1))) + $self->{BOX}[2]); my $grid = plAllocGrid ($xmap, $ymap); plshades($z, @{$self->{BOX}}, $clevel, $fill_width, $cont_color, $cont_width, $rectangular, 0, \&pltr1, $grid); plFreeGrid ($grid); # Flush the PLplot stream. plflush(); } # handle histograms sub histogram { my $self = shift; my $x = shift; my $nbins = shift; my %opts = @_; # Set PLplot to right output stream plsstrm($self->{STREAMNUMBER}); # apply options $self->setparm(%opts); my ($min, $max) = $x->minmax; unless (exists($self->{BOX})) { $self->{BOX} = [$min, $max, 0, $x->nelem]; # box probably too tall! } # set up plotting box $self->_setwindow; # draw labels $self->_drawlabels; # plot box plcol0 (1); # set to frame color plbox ($self->{XTICK}, $self->{NXSUB}, $self->{YTICK}, $self->{NYSUB}, $self->{XBOX}, $self->{YBOX}); # !!! note out of order call # set line style for plot only (not box) pllsty ($self->{LINESTYLE}); # set line width for plot only (not box) plwid ($self->{LINEWIDTH}); # set color for histograms plcol0 ($self->{CURRENT_COLOR_IDX}); plhist ($x, $min, $max, $nbins, 1); # '1' is oldbins parm: dont call plenv! # set line width back plwid (0); # Flush the PLplot stream. plflush(); } # Draw bar graphs sub bargraph { my $self = shift; my $labels = shift; # ref to perl list of labels for bars my $values = shift; # pdl of values for bars my %opts = @_; # max number of readable labels on x axis my $maxlab = defined($opts{MAXBARLABELS}) ? $opts{MAXBARLABELS} : 20; delete ($opts{MAXBARLABELS}); # Set PLplot to right output stream plsstrm($self->{STREAMNUMBER}); my $xmax = scalar(@$labels); # apply options $self->setparm(%opts); my ($ymin, $ymax) = $values->minmax; unless (exists($self->{BOX})) { $self->{BOX} = [0, $xmax, $ymin, $ymax]; # box probably too tall! } # set up plotting box $self->_setwindow; # draw labels $self->_drawlabels; # plot box plcol0 (1); # set to frame color plbox ($self->{XTICK}, $self->{NXSUB}, $self->{YTICK}, $self->{NYSUB}, 'bc', $self->{YBOX}); # !!! note out of order call # Now respect TEXTPOSITION setting if TEXTMODE eq 'border' # This allows the user to tweak the label placement. D. Hunt 9/4/2007 my ($side, $disp, $foo, $just) = ('BV', 0.2, 0, 1.0); if (defined($self->{TEXTMODE}) && $self->{TEXTMODE} eq 'border') { ($side, $disp, $foo, $just) = @{$self->{TEXTPOSITION}}; } # plot labels plschr (0, $self->{CHARSIZE} * 0.7); # use smaller characters my $pos = 0; my $skip = int($xmax/$maxlab) + 1; for (my $i=0;$i<$xmax;$i+=$skip) { $pos = ((0.5+$i)/$xmax); my $lab = $$labels[$i]; plmtex ($disp, $pos, $just, $side, $lab); # !!! out of order parms } plcol0 ($self->{CURRENT_COLOR_IDX}); # set back to line color # set line style for plot only (not box) pllsty ($self->{LINESTYLE}); # set line width for plot only (not box) plwid ($self->{LINEWIDTH}); # draw bars if ($self->{UNFILLED_BARS}) { plunfbox (PDL->sequence($xmax)+0.5, $values); } else { plfbox (PDL->sequence($xmax)+0.5, $values); } # set line width back plwid (0); # set char size back plschr (0, $self->{CHARSIZE}); # Flush the PLplot stream. plflush(); } # Add text to a plot sub text { my $self = shift; my $text = shift; # Set PLplot to right output stream plsstrm($self->{STREAMNUMBER}); # apply options $self->setparm(@_); # set the color according to the color specified in the object plcol0 ($self->{CURRENT_COLOR_IDX}); # plot either relative to border, or inside view port if ($self->{TEXTMODE} eq 'border') { my ($side, $disp, $pos, $just) = @{$self->{TEXTPOSITION}}; plmtex ($disp, $pos, $just, $side, $text); # !!! out of order parms } elsif ($self->{TEXTMODE} eq 'plot') { my ($x, $y, $dx, $dy, $just) = @{$self->{TEXTPOSITION}}; plptex ($x, $y, $dx, $dy, $just, $text); } # Flush the PLplot stream. plflush(); } # Clear the current page. This should only be used with interactive devices! sub clear { my $self = shift; # Set PLplot to right output stream plsstrm($self->{STREAMNUMBER}); plclear(); return; } # Get mouse click coordinates (OO version). This should only be used with interactive devices! sub cursor { my $self = shift; # Set PLplot to right output stream plsstrm($self->{STREAMNUMBER}); # Flush the stream, to make sure the plot is visible & current plflush(); # Get the cursor position my %gin = plGetCursor(); # Return an array with the coordinates of the mouse click return ($gin{"wX"}, $gin{"wY"}, $gin{"pX"}, $gin{"pY"}, $gin{"dX"}, $gin{"dY"}); } # Explicitly close a plot and free the object sub close { my $self = shift; # Set PLplot to right output stream plsstrm($self->{STREAMNUMBER}); plend1 (); # Return this stream number to the pool. push (@plplot_stream_pool, $self->{STREAMNUMBER}); delete $self->{STREAMNUMBER}; return; } EOD # Used throughout when generating documentation. my $doc; # Necessary includes for .xs file pp_addhdr(<<'EOH'); #include #include #include #include EOH # The create_low_level_constants function is used to make the #define'd # constants in plplot.h available in Perl in the form of functions. It # should be then possible to write code like this: # # plParseOpts (\@ARGV, PL_PARSE_SKIP | PL_PARSE_NOPROGRAM); sub create_low_level_constants { my $defn = shift; my @lines = split (/\n/, $defn); foreach my $line (@lines) { next if (($line =~ /^\#/) or ($line =~ /^\s*$/)); foreach my $const ($line =~ /([^\s]+)/g) { my $func = <<"EOC"; int $const() PROTOTYPE: CODE: RETVAL = $const; OUTPUT: RETVAL EOC pp_addxs ($func); pp_add_exported ($const); } } } create_low_level_constants (<<'EODEF'); # Definitions used in plParseOpts PL_PARSE_PARTIAL PL_PARSE_FULL PL_PARSE_QUIET PL_PARSE_NODELETE PL_PARSE_SHOWALL PL_PARSE_OVERRIDE PL_PARSE_NOPROGRAM PL_PARSE_NODASH PL_PARSE_SKIP # Definitions for plmesh and plsurf3d DRAW_LINEX DRAW_LINEY DRAW_LINEXY MAG_COLOR BASE_CONT TOP_CONT SURF_CONT DRAW_SIDES FACETED MESH # fonts PL_FCI_SANS PL_FCI_MONO # Input event (especially keyboard) definitions for use from plplot # event handlers. PLK_BackSpace PLK_Tab PLK_Linefeed PLK_Return PLK_Escape PLK_Delete PLK_Clear PLK_Pause PLK_Scroll_Lock PLK_Home PLK_Left PLK_Up PLK_Right PLK_Down PLK_Prior PLK_Next PLK_End PLK_Begin PLK_Select PLK_Print PLK_Execute PLK_Insert PLK_Undo PLK_Redo PLK_Menu PLK_Find PLK_Cancel PLK_Help PLK_Break PLK_Mode_switch PLK_script_switch PLK_Num_Lock PLK_KP_Space PLK_KP_Tab PLK_KP_Enter PLK_KP_F1 PLK_KP_F2 PLK_KP_F3 PLK_KP_F4 PLK_KP_Equal PLK_KP_Multiply PLK_KP_Add PLK_KP_Separator PLK_KP_Subtract PLK_KP_Decimal PLK_KP_Divide PLK_KP_0 PLK_KP_1 PLK_KP_2 PLK_KP_3 PLK_KP_4 PLK_KP_5 PLK_KP_6 PLK_KP_7 PLK_KP_8 PLK_KP_9 PLK_F1 PLK_F2 PLK_F3 PLK_F4 PLK_F5 PLK_F6 PLK_F7 PLK_F8 PLK_F9 PLK_F10 PLK_F11 PLK_L1 PLK_F12 PLK_L2 PLK_F13 PLK_L3 PLK_F14 PLK_L4 PLK_F15 PLK_L5 PLK_F16 PLK_L6 PLK_F17 PLK_L7 PLK_F18 PLK_L8 PLK_F19 PLK_L9 PLK_F20 PLK_L10 PLK_F21 PLK_R1 PLK_F22 PLK_R2 PLK_F23 PLK_R3 PLK_F24 PLK_R4 PLK_F25 PLK_R5 PLK_F26 PLK_R6 PLK_F27 PLK_R7 PLK_F28 PLK_R8 PLK_F29 PLK_R9 PLK_F30 PLK_R10 PLK_F31 PLK_R11 PLK_F32 PLK_R12 PLK_R13 PLK_F33 PLK_F34 PLK_R14 PLK_F35 PLK_R15 PLK_Shift_L PLK_Shift_R PLK_Control_L PLK_Control_R PLK_Caps_Lock PLK_Shift_Lock PLK_Meta_L PLK_Meta_R PLK_Alt_L PLK_Alt_R PLK_Super_L PLK_Super_R PLK_Hyper_L PLK_Hyper_R # Type of gridding algorithm for plgriddata () GRID_CSA GRID_DTLI GRID_NNI GRID_NNIDW GRID_NNLI GRID_NNAIDW EODEF create_low_level_constants (<<'EODEF') if ($plversion->{'c_pllegend'}); # Definitions for plslabelfunc PL_X_AXIS PL_Y_AXIS PL_Z_AXIS # Definitions for colorbar PL_COLORBAR_SHADE PL_COLORBAR_SHADE_LABEL PL_COLORBAR_IMAGE PL_COLORBAR_GRADIENT PL_COLORBAR_CAP_LOW PL_COLORBAR_CAP_HIGH PL_COLORBAR_LABEL_LEFT PL_COLORBAR_LABEL_RIGHT PL_COLORBAR_LABEL_TOP PL_COLORBAR_LABEL_BOTTOM # Definitions for pllegend PL_LEGEND_BACKGROUND PL_LEGEND_BOUNDING_BOX PL_LEGEND_COLOR_BOX PL_LEGEND_LINE PL_LEGEND_NONE PL_LEGEND_ROW_MAJOR PL_LEGEND_SYMBOL PL_LEGEND_TEXT_LEFT PL_POSITION_BOTTOM PL_POSITION_INSIDE PL_POSITION_LEFT PL_POSITION_OUTSIDE PL_POSITION_RIGHT PL_POSITION_SUBPAGE PL_POSITION_TOP PL_POSITION_VIEWPORT EODEF # This subroutine is used to reorder PP_DEF arguments into the # standard plplot order. This is necessary for some plplot functions # with a mixture of PDL and non-PDL arguments, since pp_def requires # all non-PDL (OtherPars) arguments to go at the end of the argument list pp_addpm (<<'EOPM'); sub reorder { my $name = shift; my %reorder = ( plaxes => [0,1,6,2,3,7,4,5], plbox => [4,0,1,5,2,3], # 4th arg -> 0th arg, 0th arg -> 1st arg, etc plbox3 => [6,7,0,1,8,9,2,3,10,11,4,5], plmtex => [3,0,1,2,4], plmtex3 => [3,0,1,2,4], plstart => [2,0,1], plstripc => [13,14,0,1,2,3,4,5,6,7,8,9,10,11,12,15,16,17,18], plmap => [4,5,0,1,2,3], plmeridians => [6,0,1,2,3,4,5], plshades => [0,10,1,2,3,4,5,6,7,8,9,11,12], plshade1 => [0,15,1,2,3,4,5,6,7,8,9,10,11,12,13,14,16,17], ); my $ordering = $reorder{$name}; die "Cannot find argument reordering for $name" if (!defined($ordering)); my @out; for (my $i=0;$i<@_;$i++) { $out[$$ordering[$i]] = $_[$i]; } return @out; } # Routine for users to set normal plplot argument order sub plplot_use_standard_argument_order { $PDL::Graphics::PLplot::standard_order = shift; } EOPM pp_add_exported('plplot_use_standard_argument_order'); =head1 LOW-LEVEL INTERFACE =cut pp_addpm (<<'EOPM'); =pod The PDL low-level interface to the PLplot library closely mimics the C API. Users are referred to the PLplot User's Manual, distributed with the source PLplot tarball. This manual is also available on-line at the PLplot web site (L). There are three differences in the way the functions are called. The first one is due to a limitation in the pp_def wrapper of PDL, which forces all the non-piddle arguments to be at the end of the arguments list. It is the case of strings (C) arguments in the C API. This affects the following functions: plaxes plbox plbox3 plmtex plmtex3 plstart plstripc plmap plmeridians plshades plshade1 This difference can be got around by a call to plplot_use_standard_argument_order(1); This re-arranges the string arguments to their proper/intuitive position compared with the C plplot interface. This can be restored to it's default by calling: plplot_use_standard_argument_order(0); The second notable different between the C and the PDL APIs is that many of the PDL calls do not need arguments to specify the size of the the vectors and/or matrices being passed. These size parameters are deduced from the size of the piddles, when possible and are just omitted from the C call when translating it to perl. The third difference has to do with output parameters. In C these are passed in with the input parameters. In the perl interface, they are omitted. For example: C: pllegend(&p_legend_width, &p_legend_height, opt, position, x, y, plot_width, bg_color, bb_color, bb_style, nrow, ncolumn, nlegend, opt_array, text_offset, text_scale, text_spacing, text_justification, text_colors, (const char **)text, box_colors, box_patterns, box_scales, box_line_widths, line_colors, line_styles, line_widths, symbol_colors, symbol_scales, symbol_numbers, (const char **)symbols); perl: my ($legend_width, $legend_height) = pllegend ($position, $opt, $x, $y, $plot_width, $bg_color, $nlegend, \@opt_array, $text_offset, $text_scale, $text_spacing, $test_justification, \@text_colors, \@text, \@box_colors, \@box_patterns, \@box_scales, \@line_colors, \@line_styles, \@line_widths, \@symbol_colors, \@symbol_scales, \@symbol_numbers, \@symbols); Some of the API functions implemented in PDL have other specificities in comparison with the C API and will be discussed below. =cut EOPM pp_def ('pladv', NoPthread => 1, GenericTypes => [D], Pars => 'int page()', OtherPars => '', Code => 'c_pladv($page());' ); sub plaxes { plaxes_pp ($standard_order ? reorder ('plaxes', @_) : @_) } pp_add_exported('plaxes'); pp_def ('plaxes_pp', NoPthread => 1, GenericTypes => [D], Pars => 'double xzero();double yzero();double xtick();int nxsub();double ytick();int nysub()', OtherPars => 'char *xopt;char *yopt', Code => 'c_plaxes($xzero(),$yzero(),$COMP(xopt),$xtick(),$nxsub(),$COMP(yopt),$ytick(),$nysub());' ); pp_def ('plbin', NoPthread => 1, GenericTypes => [D], Pars => 'int nbin();double x(dima);double y(dima);int center()', OtherPars => '', Code => 'c_plbin($nbin(),$P(x),$P(y),$center());' ); pp_addxs (<<"EOC"); void plbop() CODE: c_plbop(); EOC pp_add_exported('plbop'); pp_addpm ('sub plbox { plbox_pp ($standard_order ? reorder ("plbox", @_) : @_) }'); pp_add_exported('plbox'); pp_def ('plbox_pp', NoPthread => 1, GenericTypes => [D], Pars => 'double xtick();int nxsub();double ytick();int nysub()', OtherPars => 'char *xopt;char *yopt', Code => 'c_plbox($COMP(xopt),$xtick(),$nxsub(),$COMP(yopt),$ytick(),$nysub());' ); pp_addpm ('sub plbox3 { plbox3_pp ($standard_order ? reorder ("plbox3", @_) : @_) }'); pp_add_exported('plbox3'); pp_def ('plbox3_pp', NoPthread => 1, GenericTypes => [D], Pars => 'double xtick();int nsubx();double ytick();int nsuby();double ztick();int nsubz()', OtherPars => 'char *xopt;char *xlabel;char *yopt;char *ylabel;char *zopt;char *zlabel', Code => 'c_plbox3($COMP(xopt),$COMP(xlabel),$xtick(),$nsubx(),$COMP(yopt),$COMP(ylabel),$ytick(),$nsuby(),$COMP(zopt),$COMP(zlabel),$ztick(),$nsubz());' ); pp_addxs (<<"EOC"); void plclear() CODE: c_plclear(); EOC pp_add_exported('plclear'); pp_def ('plcol0', NoPthread => 1, GenericTypes => [D], Pars => 'int icolzero()', OtherPars => '', Code => 'c_plcol0($icolzero());' ); pp_def ('plcol1', NoPthread => 1, GenericTypes => [D], Pars => 'double colone()', OtherPars => '', Code => 'c_plcol1($colone());' ); pp_def ('plcpstrm', NoPthread => 1, GenericTypes => [D], Pars => 'int iplsr();int flags()', OtherPars => '', Code => 'c_plcpstrm($iplsr(),$flags());' ); pp_def ('pldid2pc', NoPthread => 1, GenericTypes => [D], Pars => 'double xmin(dima);double ymin(dima);double xmax(dima);double ymax(dima)', OtherPars => '', Code => 'pldid2pc($P(xmin),$P(ymin),$P(xmax),$P(ymax));' ); pp_def ('pldip2dc', NoPthread => 1, GenericTypes => [D], Pars => 'double xmin(dima);double ymin(dima);double xmax(dima);double ymax(dima)', OtherPars => '', Code => 'pldip2dc($P(xmin),$P(ymin),$P(xmax),$P(ymax));' ); pp_addxs (<<"EOC"); void plend() CODE: c_plend(); EOC pp_add_exported('plend'); pp_addxs (<<"EOC"); void plend1() CODE: c_plend1(); EOC pp_add_exported('plend1'); pp_def ('plenv', NoPthread => 1, GenericTypes => [D], Pars => 'double xmin();double xmax();double ymin();double ymax();int just();int axis()', OtherPars => '', Code => 'c_plenv($xmin(),$xmax(),$ymin(),$ymax(),$just(),$axis());' ); pp_def ('plenv0', NoPthread => 1, GenericTypes => [D], Pars => 'double xmin();double xmax();double ymin();double ymax();int just();int axis()', OtherPars => '', Code => 'c_plenv0($xmin(),$xmax(),$ymin(),$ymax(),$just(),$axis());' ); pp_addxs (<<"EOC"); void pleop() CODE: c_pleop(); EOC pp_add_exported('pleop'); pp_def ('plerrx', NoPthread => 1, GenericTypes => [D], Pars => 'int n();double xmin(dima);double xmax(dima);double y(dima)', OtherPars => '', Code => 'c_plerrx($n(),$P(xmin),$P(xmax),$P(y));' ); pp_def ('plerry', NoPthread => 1, GenericTypes => [D], Pars => 'int n();double x(dima);double ymin(dima);double ymax(dima)', OtherPars => '', Code => 'c_plerry($n(),$P(x),$P(ymin),$P(ymax));' ); pp_addxs (<<"EOC"); void plfamadv() CODE: c_plfamadv(); EOC pp_add_exported('plfamadv'); pp_def ('plfill3', NoPthread => 1, GenericTypes => [D], Pars => 'int n();double x(dima);double y(dima);double z(dima)', OtherPars => '', Code => 'c_plfill3($n(),$P(x),$P(y),$P(z));' ); pp_addxs (<<"EOC"); void plflush() CODE: c_plflush(); EOC pp_add_exported('plflush'); pp_def ('plfont', NoPthread => 1, GenericTypes => [D], Pars => 'int ifont()', OtherPars => '', Code => 'c_plfont($ifont());' ); pp_def ('plfontld', NoPthread => 1, GenericTypes => [D], Pars => 'int fnt()', OtherPars => '', Code => 'c_plfontld($fnt());' ); pp_def ('plgchr', NoPthread => 1, GenericTypes => [D], Pars => 'double [o]p_def();double [o]p_ht()', OtherPars => '', Code => 'c_plgchr($P(p_def),$P(p_ht));' ); pp_def ('plgcompression', NoPthread => 1, GenericTypes => [D], Pars => 'int [o]compression()', OtherPars => '', Code => 'c_plgcompression($P(compression));' ); pp_def ('plgdidev', NoPthread => 1, GenericTypes => [D], Pars => 'double [o]p_mar();double [o]p_aspect();double [o]p_jx();double [o]p_jy()', OtherPars => '', Code => 'c_plgdidev($P(p_mar),$P(p_aspect),$P(p_jx),$P(p_jy));' ); pp_def ('plgdiori', NoPthread => 1, GenericTypes => [D], Pars => 'double [o]p_rot()', OtherPars => '', Code => 'c_plgdiori($P(p_rot));' ); pp_def ('plgdiplt', NoPthread => 1, GenericTypes => [D], Pars => 'double [o]p_xmin();double [o]p_ymin();double [o]p_xmax();double [o]p_ymax()', OtherPars => '', Code => 'c_plgdiplt($P(p_xmin),$P(p_ymin),$P(p_xmax),$P(p_ymax));' ); pp_def ('plgfam', NoPthread => 1, GenericTypes => [D], Pars => 'int [o]p_fam();int [o]p_num();int [o]p_bmax()', OtherPars => '', Code => 'c_plgfam($P(p_fam),$P(p_num),$P(p_bmax));' ); pp_def ('plglevel', NoPthread => 1, GenericTypes => [D], Pars => 'int [o]p_level()', OtherPars => '', Code => 'c_plglevel($P(p_level));' ); pp_def ('plgpage', NoPthread => 1, GenericTypes => [D], Pars => 'double [o]p_xp();double [o]p_yp();int [o]p_xleng();int [o]p_yleng();int [o]p_xoff();int [o]p_yoff()', OtherPars => '', Code => 'c_plgpage($P(p_xp),$P(p_yp),$P(p_xleng),$P(p_yleng),$P(p_xoff),$P(p_yoff));' ); pp_addxs (<<"EOC"); void plgra() CODE: c_plgra(); EOC pp_add_exported('plgra'); pp_def ('plgspa', NoPthread => 1, GenericTypes => [D], Pars => 'double [o]xmin();double [o]xmax();double [o]ymin();double [o]ymax()', OtherPars => '', Code => 'c_plgspa($P(xmin),$P(xmax),$P(ymin),$P(ymax));' ); pp_def ('plgvpd', NoPthread => 1, GenericTypes => [D], Pars => 'double [o]p_xmin();double [o]p_xmax();double [o]p_ymin();double [o]p_ymax()', OtherPars => '', Code => 'c_plgvpd($P(p_xmin),$P(p_xmax),$P(p_ymin),$P(p_ymax));' ); pp_def ('plgvpw', NoPthread => 1, GenericTypes => [D], Pars => 'double [o]p_xmin();double [o]p_xmax();double [o]p_ymin();double [o]p_ymax()', OtherPars => '', Code => 'c_plgvpw($P(p_xmin),$P(p_xmax),$P(p_ymin),$P(p_ymax));' ); pp_def ('plgxax', NoPthread => 1, GenericTypes => [D], Pars => 'int [o]p_digmax();int [o]p_digits()', OtherPars => '', Code => 'c_plgxax($P(p_digmax),$P(p_digits));' ); pp_def ('plgyax', NoPthread => 1, GenericTypes => [D], Pars => 'int [o]p_digmax();int [o]p_digits()', OtherPars => '', Code => 'c_plgyax($P(p_digmax),$P(p_digits));' ); pp_def ('plgzax', NoPthread => 1, GenericTypes => [D], Pars => 'int [o]p_digmax();int [o]p_digits()', OtherPars => '', Code => 'c_plgzax($P(p_digmax),$P(p_digits));' ); pp_addxs (<<"EOC"); void plinit() CODE: c_plinit(); EOC pp_add_exported('plinit'); pp_def ('pljoin', NoPthread => 1, GenericTypes => [D], Pars => 'double xone();double yone();double xtwo();double ytwo()', OtherPars => '', Code => 'c_pljoin($xone(),$yone(),$xtwo(),$ytwo());' ); pp_addxs (<<"EOC"); void pllab(xlabel,ylabel,tlabel) char * xlabel char * ylabel char * tlabel CODE: c_pllab(xlabel,ylabel,tlabel); EOC pp_add_exported('pllab'); pp_def ('pllightsource', NoPthread => 1, GenericTypes => [D], Pars => 'double x();double y();double z()', OtherPars => '', Code => 'c_pllightsource($x(),$y(),$z());' ); pp_def ('pllsty', NoPthread => 1, GenericTypes => [D], Pars => 'int lin()', OtherPars => '', Code => 'c_pllsty($lin());' ); pp_addpm ('sub plmtex { plmtex_pp ($standard_order ? reorder ("plmtex", @_) : @_) }'); pp_add_exported('plmtex'); pp_def ('plmtex_pp', NoPthread => 1, GenericTypes => [D], Pars => 'double disp();double pos();double just()', OtherPars => 'char *side;char *text', Code => 'c_plmtex($COMP(side),$disp(),$pos(),$just(),$COMP(text));' ); pp_addpm ('sub plmtex3 { plmtex3_pp ($standard_order ? reorder ("plmtex3", @_) : @_) }'); pp_add_exported('plmtex3'); pp_def ('plmtex3_pp', NoPthread => 1, GenericTypes => [D], Pars => 'double disp();double pos();double just()', OtherPars => 'char *side;char *text', Code => 'c_plmtex3($COMP(side),$disp(),$pos(),$just(),$COMP(text));' ); pp_def ('plpat', NoPthread => 1, GenericTypes => [D], Pars => 'int nlin();int inc(dima);int del(dima)', OtherPars => '', Code => 'c_plpat($nlin(),$P(inc),$P(del));' ); pp_def ('plprec', NoPthread => 1, GenericTypes => [D], Pars => 'int setp();int prec()', OtherPars => '', Code => 'c_plprec($setp(),$prec());' ); pp_def ('plpsty', NoPthread => 1, GenericTypes => [D], Pars => 'int patt()', OtherPars => '', Code => 'c_plpsty($patt());' ); pp_def ('plptex', NoPthread => 1, GenericTypes => [D], Pars => 'double x();double y();double dx();double dy();double just()', OtherPars => 'char *text', Code => 'c_plptex($x(),$y(),$dx(),$dy(),$just(),$COMP(text));' ); pp_def ('plptex3', NoPthread => 1, GenericTypes => [D], Pars => 'double x();double y();double z();double dx();double dy();double dz();double sx();double sy();double sz();double just()', OtherPars => 'char *text', Code => 'c_plptex3($x(),$y(),$z(),$dx(),$dy(),$dz(),$sx(),$sy(),$sz(),$just(),$COMP(text));' ); pp_addxs (<<"EOC"); void plreplot() CODE: c_plreplot(); EOC pp_add_exported('plreplot'); pp_def ('plschr', NoPthread => 1, GenericTypes => [D], Pars => 'double def();double scale()', OtherPars => '', Code => 'c_plschr($def(),$scale());' ); pp_def ('plscmap0n', NoPthread => 1, GenericTypes => [D], Pars => 'int ncolzero()', OtherPars => '', Code => 'c_plscmap0n($ncolzero());' ); pp_def ('plscmap1n', NoPthread => 1, GenericTypes => [D], Pars => 'int ncolone()', OtherPars => '', Code => 'c_plscmap1n($ncolone());' ); pp_def ('plscol0', NoPthread => 1, GenericTypes => [D], Pars => 'int icolzero();int r();int g();int b()', OtherPars => '', Code => 'c_plscol0($icolzero(),$r(),$g(),$b());' ); pp_def ('plscolbg', NoPthread => 1, GenericTypes => [D], Pars => 'int r();int g();int b()', OtherPars => '', Code => 'c_plscolbg($r(),$g(),$b());' ); pp_def ('plscolor', NoPthread => 1, GenericTypes => [D], Pars => 'int color()', OtherPars => '', Code => 'c_plscolor($color());' ); pp_def ('plscompression', NoPthread => 1, GenericTypes => [D], Pars => 'int compression()', OtherPars => '', Code => 'c_plscompression($compression());' ); pp_addxs (<<"EOC"); void plsdev(devname) char * devname CODE: c_plsdev(devname); EOC pp_add_exported('plsdev'); pp_def ('plsdidev', NoPthread => 1, GenericTypes => [D], Pars => 'double mar();double aspect();double jx();double jy()', OtherPars => '', Code => 'c_plsdidev($mar(),$aspect(),$jx(),$jy());' ); pp_def ('plsdimap', NoPthread => 1, GenericTypes => [D], Pars => 'int dimxmin();int dimxmax();int dimymin();int dimymax();double dimxpmm();double dimypmm()', OtherPars => '', Code => 'c_plsdimap($dimxmin(),$dimxmax(),$dimymin(),$dimymax(),$dimxpmm(),$dimypmm());' ); pp_def ('plsdiori', NoPthread => 1, GenericTypes => [D], Pars => 'double rot()', OtherPars => '', Code => 'c_plsdiori($rot());' ); pp_def ('plsdiplt', NoPthread => 1, GenericTypes => [D], Pars => 'double xmin();double ymin();double xmax();double ymax()', OtherPars => '', Code => 'c_plsdiplt($xmin(),$ymin(),$xmax(),$ymax());' ); pp_def ('plsdiplz', NoPthread => 1, GenericTypes => [D], Pars => 'double xmin();double ymin();double xmax();double ymax()', OtherPars => '', Code => 'c_plsdiplz($xmin(),$ymin(),$xmax(),$ymax());' ); pp_def ('pl_setcontlabelparam', NoPthread => 1, GenericTypes => [D], Pars => 'double offset();double size();double spacing();int active()', OtherPars => '', Code => 'c_pl_setcontlabelparam($offset(),$size(),$spacing(),$active());' ); pp_def ('pl_setcontlabelformat', NoPthread => 1, GenericTypes => [D], Pars => 'int lexp();int sigdig()', OtherPars => '', Code => 'c_pl_setcontlabelformat($lexp(),$sigdig());' ); pp_def ('plsfam', NoPthread => 1, GenericTypes => [D], Pars => 'int fam();int num();int bmax()', OtherPars => '', Code => 'c_plsfam($fam(),$num(),$bmax());' ); pp_addxs (<<"EOC"); void plsfnam(fnam) char * fnam CODE: c_plsfnam(fnam); EOC pp_add_exported('plsfnam'); pp_def ('plsmaj', NoPthread => 1, GenericTypes => [D], Pars => 'double def();double scale()', OtherPars => '', Code => 'c_plsmaj($def(),$scale());' ); pp_def ('plsmin', NoPthread => 1, GenericTypes => [D], Pars => 'double def();double scale()', OtherPars => '', Code => 'c_plsmin($def(),$scale());' ); pp_def ('plsori', NoPthread => 1, GenericTypes => [D], Pars => 'int ori()', OtherPars => '', Code => 'c_plsori($ori());' ); pp_def ('plspage', NoPthread => 1, GenericTypes => [D], Pars => 'double xp();double yp();int xleng();int yleng();int xoff();int yoff()', OtherPars => '', Code => 'c_plspage($xp(),$yp(),$xleng(),$yleng(),$xoff(),$yoff());' ); pp_def ('plspause', NoPthread => 1, GenericTypes => [D], Pars => 'int pause()', OtherPars => '', Code => 'c_plspause($pause());' ); pp_def ('plsstrm', NoPthread => 1, GenericTypes => [D], Pars => 'int strm()', OtherPars => '', Code => 'c_plsstrm($strm());' ); pp_def ('plssub', NoPthread => 1, GenericTypes => [D], Pars => 'int nx();int ny()', OtherPars => '', Code => 'c_plssub($nx(),$ny());' ); pp_def ('plssym', NoPthread => 1, GenericTypes => [D], Pars => 'double def();double scale()', OtherPars => '', Code => 'c_plssym($def(),$scale());' ); pp_def ('plstar', NoPthread => 1, GenericTypes => [D], Pars => 'int nx();int ny()', OtherPars => '', Code => 'c_plstar($nx(),$ny());' ); pp_addpm ('sub plstart { plstart_pp ($standard_order ? reorder ("plstart", @_) : @_) }'); pp_add_exported('plstart'); pp_def ('plstart_pp', NoPthread => 1, GenericTypes => [D], Pars => 'int nx();int ny()', OtherPars => 'char *devname', Code => 'c_plstart($COMP(devname),$nx(),$ny());' ); pp_def ('plstripa', NoPthread => 1, GenericTypes => [D], Pars => 'int id();int pen();double x();double y()', OtherPars => '', Code => 'c_plstripa($id(),$pen(),$x(),$y());' ); pp_def ('plstripd', NoPthread => 1, GenericTypes => [D], Pars => 'int id()', OtherPars => '', Code => 'c_plstripd($id());' ); pp_def ('plsvpa', NoPthread => 1, GenericTypes => [D], Pars => 'double xmin();double xmax();double ymin();double ymax()', OtherPars => '', Code => 'c_plsvpa($xmin(),$xmax(),$ymin(),$ymax());' ); pp_def ('plsxax', NoPthread => 1, GenericTypes => [D], Pars => 'int digmax();int digits()', OtherPars => '', Code => 'c_plsxax($digmax(),$digits());' ); pp_def ('plsxwin', NoPthread => 1, GenericTypes => [D], Pars => 'int window_id()', OtherPars => '', Code => 'plsxwin($window_id());' ); pp_def ('plsyax', NoPthread => 1, GenericTypes => [D], Pars => 'int digmax();int digits()', OtherPars => '', Code => 'c_plsyax($digmax(),$digits());' ); pp_def ('plszax', NoPthread => 1, GenericTypes => [D], Pars => 'int digmax();int digits()', OtherPars => '', Code => 'c_plszax($digmax(),$digits());' ); pp_addxs (<<"EOC"); void pltext() CODE: c_pltext(); EOC pp_add_exported('pltext'); pp_def ('plvasp', NoPthread => 1, GenericTypes => [D], Pars => 'double aspect()', OtherPars => '', Code => 'c_plvasp($aspect());' ); pp_def ('plvpas', NoPthread => 1, GenericTypes => [D], Pars => 'double xmin();double xmax();double ymin();double ymax();double aspect()', OtherPars => '', Code => 'c_plvpas($xmin(),$xmax(),$ymin(),$ymax(),$aspect());' ); pp_def ('plvpor', NoPthread => 1, GenericTypes => [D], Pars => 'double xmin();double xmax();double ymin();double ymax()', OtherPars => '', Code => 'c_plvpor($xmin(),$xmax(),$ymin(),$ymax());' ); pp_addxs (<<"EOC"); void plvsta() CODE: c_plvsta(); EOC pp_add_exported('plvsta'); pp_def ('plw3d', NoPthread => 1, GenericTypes => [D], Pars => 'double basex();double basey();double height();double xminzero();double xmaxzero();double yminzero();double ymaxzero();double zminzero();double zmaxzero();double alt();double az()', OtherPars => '', Code => 'c_plw3d($basex(),$basey(),$height(),$xminzero(),$xmaxzero(),$yminzero(),$ymaxzero(),$zminzero(),$zmaxzero(),$alt(),$az());' ); pp_def ('plwid', NoPthread => 1, GenericTypes => [D], Pars => 'int width()', OtherPars => '', Code => 'c_plwid($width());' ); pp_def ('plwind', NoPthread => 1, GenericTypes => [D], Pars => 'double xmin();double xmax();double ymin();double ymax()', OtherPars => '', Code => 'c_plwind($xmin(),$xmax(),$ymin(),$ymax());' ); pp_addxs (<<"EOC"); void plsetopt(opt,optarg) char * opt char * optarg CODE: c_plsetopt(opt,optarg); EOC pp_add_exported('plsetopt'); pp_def ('plP_gpixmm', NoPthread => 1, GenericTypes => [D], Pars => 'double p_x(dima);double p_y(dima)', OtherPars => '', Code => 'plP_gpixmm($P(p_x),$P(p_y));' ); pp_def ('plscolbga', NoPthread => 1, GenericTypes => [D], Pars => 'int r();int g();int b();double a()', OtherPars => '', Code => 'c_plscolbga($r(),$g(),$b(),$a());' ); pp_def ('plscol0a', NoPthread => 1, GenericTypes => [D], Pars => 'int icolzero();int r();int g();int b();double a()', OtherPars => '', Code => 'c_plscol0a($icolzero(),$r(),$g(),$b(),$a());' ); # C routine to draw lines with gaps. This is useful for map continents and other things. =head2 plline =cut $doc = <<'EOD'; =for ref Draws line segments along (x1,y1)->(x2,y2)->(x3,y3)->... =for bad If the nth value of either x or y are bad, then it will be skipped, breaking the line. In this way, you can specify multiple line segments with a single pair of x and y piddles. The usage is straight-forward: =for usage plline($x, $y); For example: =for example # Draw a sine wave $x = sequence(100)/10; $y = sin($x); # Draws the sine wave: plline($x, $y); # Set values above 3/4 to 'bad', effectively drawing a bunch of detached, # capped waves $y->setbadif($y > 3/4); plline($x, $y); =cut EOD pp_def ('plline', NoPthread => 1, Pars => 'x(n); y(n)', GenericTypes => [D], HandleBad => 1, NoBadifNaN => 1, Code => 'c_plline($SIZE(n),$P(x),$P(y));', BadCode => 'int i; int j; for (i=1;i<$SIZE(n);i++) { j = i-1; /* PP does not like using i-1 in a PDL ref. Use j instead. */ if ($ISGOOD(x(n=>i)) && $ISGOOD(x(n=>j)) && $ISGOOD(y(n=>i)) && $ISGOOD(y(n=>j))) { c_pljoin ($x(n=>j), $y(n=>j), $x(n=>i), $y(n=>i)); } }', Doc => $doc, ); =head2 plcolorpoints =cut $doc = <<'EOD'; =for ref PDL-specific: Implements what amounts to a threaded version of plsym. =for bad Bad values for z are simply skipped; all other bad values are not processed. In the following usage, all of the piddles must have the same dimensions: =for usage plcolorpoints($x, $y, $z, $symbol_index, $minz, $maxz) For example: =for example # Generate a parabola some points my $x = sequence(30) / 3; # Regular sampling my $y = $x**2; # Parabolic y my $z = 30 - $x**3; # Cubic coloration my $symbols = floor($x); # Use different symbols for each 1/3 of the plot # These should be integers. plcolorpoints($x, $y, $z, $symbols, -5, 20); # Thread over everything plcolorpoints($x, $y, 1, 1, -1, 2); # same color and symbol for all =cut EOD # C routine to draw points with a color scale pp_def ('plcolorpoints', NoPthread => 1, Pars => 'x(n); y(n); z(n); int sym(); minz(); maxz()', GenericTypes => [D], HandleBad => 1, Code => 'int i; int j; int ns = $SIZE(n); PLFLT zrange, ci; zrange = $maxz() - $minz(); for (i=0;ii) - $minz()) / zrange; /* get color idx in 0-1 range */ if (ci < 0) ci = 0; /* enforce bounds */ if (ci > 1) ci = 1; c_plcol1 (ci); /* set current color */ c_plsym (1, &$x(n=>i), &$y(n=>i), $sym()); /* plot it */ }', BadCode => 'int i; int j; int ns = $SIZE(n); PLFLT zrange, ci; zrange = $maxz() - $minz(); for (i=0;ii))) continue; ci = (zrange == 0.0) ? 0.5 : ($z(n=>i) - $minz()) / zrange; /* get color idx in 0-1 range */ if (ci < 0) ci = 0; /* enforce bounds */ if (ci > 1) ci = 1; c_plcol1 (ci); /* set current color */ c_plsym (1, &$x(n=>i), &$y(n=>i), $sym()); /* plot it */ }', Doc => $doc, ); pp_def ('plsmem', NoPthread => 1, GenericTypes => [B], Pars => 'int maxx();int maxy();image(3,x,y)', Code => 'c_plsmem($maxx(),$maxy(),$P(image));' ) if ($plversion->{'plsmem'}); # ## Box drawing primitive, taken from PLPLOT bar graph example # pp_def ('plfbox', NoPthread => 1, Pars => 'xo(); yo()', GenericTypes => [D], Code => 'PLFLT x[4], y[4]; x[0] = $xo() - 0.5; y[0] = 0.; x[1] = $xo() - 0.5; y[1] = $yo(); x[2] = $xo() + 0.5; y[2] = $yo(); x[3] = $xo() + 0.5; y[3] = 0.; plfill(4, x, y);', ); # ## Similar box drawing primitive, but without fill (just draw outline of box) # pp_def ('plunfbox', NoPthread => 1, Pars => 'xo(); yo()', GenericTypes => [D], Code => 'PLFLT x[4], y[4]; x[0] = $xo() - 0.5; y[0] = 0.; x[1] = $xo() - 0.5; y[1] = $yo(); x[2] = $xo() + 0.5; y[2] = $yo(); x[3] = $xo() + 0.5; y[3] = 0.; plline(4, x, y);', ); # ## Parse PLplot options given in @ARGV-like arrays # pp_def ('plParseOpts', NoPthread => 1, GenericTypes => [D], Pars => 'int [o] retval()', OtherPars => 'SV* argv; int mode', Doc => 'FIXME: documentation here!', Code => ' SV* sv = $COMP (argv); SV* dummy; AV* arr; int argc, newargc, i, retval; char** args; if ( !(SvROK (sv) && SvTYPE (SvRV (sv)) == SVt_PVAV)) { barf("plParseOpts requires an array ref"); } arr = (AV*) SvRV (sv); newargc = argc = av_len (arr) + 1; if (argc > 0) { New(1, args, argc, char *); if(args == NULL) croak("Failed to allocate memory in plParseOpts"); for (i = 0; i < argc; i++) { STRLEN len; args[i] = SvPV (* av_fetch (arr, i, 0), len); } $retval() = c_plparseopts (&newargc, (const char **)args, $COMP (mode)); for (i = 0; i < newargc; i++) av_push (arr, newSVpv (args[i], 0)); for (i = 0; i < argc; i++) dummy = av_shift (arr); /* assign to dummy to suppress compile warning */ Safefree (args); } ', ); # Plots a character at the specified points pp_def ('plpoin', NoPthread => 1, Pars => 'x(n); y(n); int code()', GenericTypes => [D], Code => 'c_plpoin($SIZE(n),$P(x),$P(y),$code());' ); # Plots a character at the specified points in 3 space pp_def ('plpoin3', NoPthread => 1, Pars => 'x(n); y(n); z(n); int code()', GenericTypes => [D], Code => 'c_plpoin3($SIZE(n),$P(x),$P(y),$P(z),$code());' ); # Draw a line in 3 space pp_def ('plline3', NoPthread => 1, Pars => 'x(n); y(n); z(n)', GenericTypes => [D], Code => 'c_plline3($SIZE(n),$P(x),$P(y),$P(z));' ); # Draws a polygon in 3 space pp_def ('plpoly3', NoPthread => 1, Pars => 'x(n); y(n); z(n); int draw(m); int ifcc()', GenericTypes => [D], Code => 'c_plpoly3($SIZE(n),$P(x),$P(y),$P(z),$P(draw),$ifcc());' ); # Plot a histogram from unbinned data pp_def ('plhist', NoPthread => 1, Pars => 'data(n); datmin(); datmax(); int nbin(); int oldwin()', GenericTypes => [D], Code => 'c_plhist($SIZE(n),$P(data),$datmin(),$datmax(),$nbin(),$oldwin());' ); # Area fill pp_def ('plfill', NoPthread => 1, Pars => 'x(n); y(n)', GenericTypes => [D], Code => 'c_plfill($SIZE(n),$P(x),$P(y));' ); # Area fill with color gradient pp_def ('plgradient', NoPthread => 1, Pars => 'x(n); y(n); angle();', GenericTypes => [D], Code => 'c_plgradient($SIZE(n),$P(x),$P(y),$angle());' ) if ($plversion->{'c_pllegend'}); # Plots a symbol at the specified points pp_def ('plsym', NoPthread => 1, Pars => 'x(n); y(n); int code()', GenericTypes => [D], Code => 'c_plsym($SIZE(n),$P(x),$P(y),$code());' ); # Plot shaded 3-d surface plot pp_def ('plsurf3d', NoPthread => 1, Pars => 'x(nx); y(ny); z(nx,ny); int opt(); clevel(nlevel);', GenericTypes => [D], Code => ' int i, j, size_x, size_y; PLFLT** zz; size_x = $SIZE(nx); size_y = $SIZE(ny); plAlloc2dGrid (&zz, size_x, size_y); for (i = 0; i < size_x; i++) for (j = 0; j < size_y; j++) zz[i][j] = $z(nx => i, ny => j); c_plsurf3d ($P(x), $P(y), (const PLFLT **)zz, size_x, size_y, $opt(), $P(clevel), $SIZE(nlevel)); plFree2dGrid (zz, size_x, size_y);' ); # Set line style pp_def ('plstyl', NoPthread => 1, Pars => 'int mark(nms); int space(nms)', GenericTypes => [D], Code => 'c_plstyl ($SIZE(nms), $P(mark), $P(space));' ); # PLplot standard random number generation. Using this # helps to keep the demo plots identical. if ($plversion->{'c_plseed'}) { pp_def ('plseed', NoPthread => 1, Pars => 'int seed()', Code => 'unsigned int useed = (unsigned int)$seed(); c_plseed(useed);' ); pp_def ('plrandd', NoPthread => 1, Pars => 'double [o]rand()', Code => '$rand() = c_plrandd();' ); } # $plversion->{'c_plseed'} # Plot contours # pltr0: Identity transformation # pltr1: Linear interpolation from singly dimensioned coord arrays # Linear interpolation from doubly dimensioned coord arrays for my $func ('pltr0', 'pltr1', 'pltr2') { pp_addxs (<<"EOC"); void $func (x, y, grid) double x double y long grid PPCODE: PLFLT tx, ty; $func (x, y, &tx, &ty, (PLPointer) grid); EXTEND (SP, 2); PUSHs (sv_2mortal (newSVnv ((double) tx))); PUSHs (sv_2mortal (newSVnv ((double) ty))); EOC pp_add_exported ($func); } # Allocates a PLcGrid object for use in pltr1 pp_def ('plAllocGrid', NoPthread => 1, Pars => "double xg(nx); double yg(ny); $pp_ptr_type [o] grid()", GenericTypes => [D], Doc => 'FIXME: documentation here!', Code => ' PLcGrid *grid; int i, nx, ny; nx = $SIZE(nx); ny = $SIZE(ny); New(1, grid, 1, PLcGrid); if(grid == NULL) croak("Failed to allocate memory for grid"); Newz(2, grid->xg, nx, PLFLT); if(grid->xg == NULL) croak("Failed to allocate memory for grid->xg"); Newz(3, grid->yg, ny, PLFLT); if(grid->yg == NULL) croak("Failed to allocate memory for grid->yg"); grid->nx = nx; grid->ny = ny; for (i = 0; i < nx; i++) grid->xg[i] = $xg(nx => i); for (i = 0; i < ny; i++) grid->yg[i] = $yg(ny => i); $grid() = ('.$int_ptr_type.') grid;' ); # Free a PLcGrid object pp_addxs (<<"EOC"); void plFreeGrid (pg) long pg PPCODE: PLcGrid* grid = (PLcGrid*) pg; free (grid->xg); free (grid->yg); free (grid); EOC pp_add_exported (plFreeGrid); # Allocates a PLcGrid2 object for use in pltr2 pp_def ('plAlloc2dGrid', NoPthread => 1, Pars => "double xg(nx,ny); double yg(nx,ny); $pp_ptr_type [o] grid()", GenericTypes => [D], Doc => 'FIXME: documentation here!', Code => ' PLcGrid2 *grid; int i, j, nx, ny; nx = $SIZE(nx); ny = $SIZE(ny); grid = (PLcGrid2*) malloc (sizeof (PLcGrid2)); plAlloc2dGrid (&(grid->xg), nx, ny); plAlloc2dGrid (&(grid->yg), nx, ny); for (i = 0; i < nx; i++) for (j = 0; j < ny; j++) { grid->xg[i][j] = $xg(nx => i, ny => j); grid->yg[i][j] = $yg(nx => i, ny => j); } grid->nx = nx; grid->ny = ny; $grid() = ('.$int_ptr_type.') grid;' ); # Free a PLcGrid2 object pp_addxs (<<"EOC"); void plFree2dGrid (pg) long pg PPCODE: PLcGrid2* grid = (PLcGrid2*) pg; plFree2dGrid (grid->xg, grid->nx, grid->ny); plFree2dGrid (grid->yg, grid->nx, grid->ny); free (grid); EOC pp_add_exported (plFree2dGrid); pp_addhdr (<<'EOH'); #define check_sub_pointer(subptr, errmsg) \ if (SvTRUE (subptr) \ && (! SvROK (subptr) || SvTYPE (SvRV (subptr)) != SVt_PVCV)) \ croak (errmsg); static SV* pltr_subroutine; static IV pltr0_iv; static IV pltr1_iv; static IV pltr2_iv; static void pltr_callback (PLFLT x, PLFLT y, PLFLT* tx, PLFLT* ty, PLPointer pltr_data) { I32 count; dSP; ENTER; SAVETMPS; PUSHMARK (SP); XPUSHs (sv_2mortal (newSVnv ((double) x))); XPUSHs (sv_2mortal (newSVnv ((double) y))); XPUSHs ((SV*) pltr_data); PUTBACK; count = call_sv (pltr_subroutine, G_ARRAY); SPAGAIN; if (count != 2) croak ("pltr: must return two scalars"); *ty = (PLFLT) POPn; *tx = (PLFLT) POPn; PUTBACK; FREETMPS; LEAVE; } static void* get_standard_pltrcb (SV* cb) { if ( !SvROK(cb) ) return NULL; /* Added to prevent bug in plshades for 0 input. D. Hunt 12/18/2008 */ IV sub = (IV) SvRV (cb); if (sub == pltr0_iv) return (void*) pltr0; else if (sub == pltr1_iv) return (void*) pltr1; else if (sub == pltr2_iv) return (void*) pltr2; else return SvTRUE (cb) ? (void*) pltr_callback : NULL; } static SV* defined_subroutine; static PLINT defined_callback (PLFLT x, PLFLT y) { I32 count, retval; dSP; ENTER; SAVETMPS; PUSHMARK (SP); XPUSHs (sv_2mortal (newSVnv ((double) x))); XPUSHs (sv_2mortal (newSVnv ((double) y))); PUTBACK; count = call_sv (defined_subroutine, G_SCALAR); SPAGAIN; if (count != 1) croak ("defined: must return one scalar"); retval = POPi; PUTBACK; FREETMPS; LEAVE; return retval; } static SV* mapform_subroutine; static void default_magic (pdl *p, int pa) { p->data = 0; } static void mapform_callback (PLINT n, PLFLT* x, PLFLT* y) { pdl *x_pdl, *y_pdl; PLFLT *tx, *ty; SV *x_sv, *y_sv; int dims, i; I32 count, ax; dSP; ENTER; SAVETMPS; dims = n; x_pdl = PDL->pdlnew (); PDL->add_deletedata_magic(x_pdl, default_magic, 0); PDL->setdims (x_pdl, &dims, 1); x_pdl->datatype = PDL_D; x_pdl->data = x; x_pdl->state |= PDL_DONTTOUCHDATA | PDL_ALLOCATED; x_sv = sv_newmortal (); PDL->SetSV_PDL (x_sv, x_pdl); y_pdl = PDL->pdlnew (); PDL->add_deletedata_magic(y_pdl, default_magic, 0); PDL->setdims (y_pdl, &dims, 1); y_pdl->datatype = PDL_D; y_pdl->data = y; y_pdl->state |= PDL_DONTTOUCHDATA | PDL_ALLOCATED; y_sv = sv_newmortal (); PDL->SetSV_PDL (y_sv, y_pdl); PUSHMARK (SP); XPUSHs (x_sv); XPUSHs (y_sv); PUTBACK; count = call_sv (mapform_subroutine, G_ARRAY); SPAGAIN; SP -= count ; ax = (SP - PL_stack_base) + 1; if (count != 2) croak ("mapform: must return two piddles"); tx = (PLFLT*) ((PDL->SvPDLV(ST(0)))->data); ty = (PLFLT*) ((PDL->SvPDLV(ST(1)))->data); for (i = 0; i < n; i++) { *(x + i) = *(tx + i); *(y + i) = *(ty + i); } PUTBACK; FREETMPS; LEAVE; } EOH pp_addhdr (<<'EOH') if ($plversion->{'c_pllegend'}); // Subroutines for adding transforms via plstransform static SV* xform_subroutine; static void xform_callback ( PLFLT x, PLFLT y, PLFLT *xt, PLFLT *yt, PLPointer data ) { SV *x_sv, *y_sv; // Perl scalars for the input x and y I32 count, ax; dSP; ENTER; SAVETMPS; x_sv = newSVnv((double)x); y_sv = newSVnv((double)y); PUSHMARK (SP); XPUSHs (x_sv); XPUSHs (y_sv); PUTBACK; count = call_sv (xform_subroutine, G_ARRAY); SPAGAIN; SP -= count ; ax = (SP - PL_stack_base) + 1; if (count != 2) croak ("xform: must return two perl scalars"); *xt = (PLFLT) SvNV(ST(0)); *yt = (PLFLT) SvNV(ST(1)); PUTBACK; FREETMPS; LEAVE; } // Subroutines for adding label formatting via plslabelfunc static SV* labelfunc_subroutine; void labelfunc_callback (PLINT axis, PLFLT value, char *label_text, PLINT length, void *data ) { SV *axis_sv, *value_sv, *length_sv; // Perl scalars for inputs I32 count, ax; dSP; ENTER; SAVETMPS; axis_sv = newSViv((IV)axis); value_sv = newSVnv((double)value); length_sv = newSViv((IV)length); PUSHMARK (SP); XPUSHs (axis_sv); XPUSHs (value_sv); XPUSHs (length_sv); PUTBACK; count = call_sv (labelfunc_subroutine, G_ARRAY); SPAGAIN; SP -= count ; ax = (SP - PL_stack_base) + 1; if (count != 1) croak ("labelfunc: must return one perl scalar"); // Copy label into output string snprintf( label_text, length, "%s", (char *)SvPV_nolen(ST(0)) ); PUTBACK; FREETMPS; LEAVE; } EOH # The init_pltr is used internally by the PLD::Graphics::PLplot # module to set the variables pltr{0,1,2}_iv to the "pointers" # of the Perl subroutines pltr{1,2,3}. These variables are later used by # get_standard_pltrcb to provide the pointers to the C function pltr{0,1,2}. # This accelerates functions like plcont and plshades when those standard # transformation functions are used. pp_def ('init_pltr', NoPthread => 1, GenericTypes => [D], Pars => '', OtherPars => 'SV* p0; SV* p1; SV* p2;', Doc => 'FIXME: documentation here!', Code => ' pltr0_iv = (IV) SvRV ($COMP(p0)); pltr1_iv = (IV) SvRV ($COMP(p1)); pltr2_iv = (IV) SvRV ($COMP(p2));'); pp_addpm (<<'EOPM'); init_pltr (\&pltr0, \&pltr1, \&pltr2); EOPM # plot continental outline in world coordinates pp_addpm ('sub plmap { plmap_pp ($standard_order ? reorder ("plmap", @_) : @_) }'); pp_add_exported('plmap'); pp_def ('plmap_pp', NoPthread => 1, Pars => 'minlong(); maxlong(); minlat(); maxlat();', # 0-3 OtherPars => 'SV* mapform; char* type;', # 4,5 GenericTypes => [D], Code => ' int use_xform; mapform_subroutine = $COMP(mapform); check_sub_pointer (mapform_subroutine, "plmap: mapform must be either 0 or a subroutine pointer"); use_xform = SvTRUE ($COMP(mapform)); plmap (use_xform ? mapform_callback : NULL, $COMP(type), $minlong(), $maxlong(), $minlat(), $maxlat());', ); # plot a string along a line pp_def ('plstring', NoPthread => 1, Pars => 'x(na); y(na);', OtherPars => 'char* string;', GenericTypes => [D], Code => 'c_plstring($SIZE(na), $P(x), $P(y), $COMP(string));', ) if ($plversion->{'c_pllegend'}); # plot a string along a 3D line pp_def ('plstring3', NoPthread => 1, Pars => 'x(na); y(na); z(na)', OtherPars => 'char* string;', GenericTypes => [D], Code => 'c_plstring3($SIZE(na), $P(x), $P(y), $P(z), $COMP(string));', ) if ($plversion->{'c_pllegend'}); # Plot the latitudes and longitudes on the background pp_addpm ('sub plmeridians { plmeridians_pp ($standard_order ? reorder ("plmeridians", @_) : @_) }'); pp_add_exported('plmeridians'); pp_def ('plmeridians_pp', NoPthread => 1, Pars => 'dlong(); dlat(); minlong(); maxlong(); minlat(); maxlat();', # 0-5 OtherPars => 'SV* mapform;', # 6 GenericTypes => [D], Code => ' mapform_subroutine = $COMP(mapform); check_sub_pointer (mapform_subroutine, "plmeridians: mapform must be either 0 or a subroutine pointer"); plmeridians (SvTRUE ($COMP(mapform)) ? mapform_callback : NULL, $dlong(), $dlat(), $minlong(), $maxlong(), $minlat(), $maxlat());' ); # Set the global coordinate transform # Shade regions on the basis of value pp_addpm ('sub plshades { plshades_pp ($standard_order ? reorder ("plshades", @_) : @_) }'); pp_add_exported('plshades'); pp_def ('plshades_pp', NoPthread => 1, Pars => 'z(x,y); xmin(); xmax(); ymin(); ymax(); clevel(l); int fill_width(); int cont_color(); int cont_width(); int rectangular()', # 0-9 OtherPars => 'SV* defined; SV* pltr; SV* pltr_data;', # 10-12 GenericTypes => [D], Code => ' int nx = $SIZE(x); int ny = $SIZE(y); int nlvl = $SIZE(l); int i, j; PLFLT **z; void (*pltrcb) (); PLPointer pltrdt; plAlloc2dGrid (&z, nx, ny); for (i = 0; i < nx; i++) for (j = 0; j < ny; j++) z[i][j] = (PLFLT) $z(x => i, y => j); defined_subroutine = $COMP(defined); check_sub_pointer (defined_subroutine, "plshades: defined must be either 0 or a subroutine pointer"); pltr_subroutine = $COMP(pltr); check_sub_pointer (pltr_subroutine, "plshades: pltr must be either 0 or a subroutine pointer"); pltrcb = get_standard_pltrcb ($COMP(pltr)); if (pltrcb != pltr_callback) pltrdt = (PLPointer) SvIV ($COMP(pltr_data)); else pltrdt = $COMP(pltr_data); c_plshades ((const PLFLT **)z, nx, ny, SvTRUE ($COMP(defined)) ? defined_callback : NULL, $xmin(), $xmax(), $ymin(), $ymax(), $P(clevel), nlvl, $fill_width(), $cont_color(), $cont_width(), plfill, $rectangular(), pltrcb, pltrdt); plFree2dGrid(z, nx, ny);', ); pp_def ('plcont', NoPthread => 1, GenericTypes => [D], Pars => 'f(nx,ny); int kx(); int lx(); int ky(); int ly(); ' . 'clevel(nlevel)', # 0-5 OtherPars => 'SV* pltr; SV* pltr_data;', # 6,7 Doc => 'FIXME: documentation here!', Code => ' int i, j, size_x, size_y; PLFLT** ff; void (*pltrcb) (); PLPointer pltrdt; size_x = $SIZE(nx); size_y = $SIZE(ny); plAlloc2dGrid (&ff, size_x, size_y); for (i = 0; i < size_x; i++) for (j = 0; j < size_y; j++) ff[i][j] = $f(nx => i, ny => j); pltr_subroutine = $COMP(pltr); check_sub_pointer (pltr_subroutine, "plcont: pltr must be either 0 or a subroutine pointer"); pltrcb = get_standard_pltrcb ($COMP(pltr)); if (pltrcb != pltr_callback) pltrdt = (PLPointer) SvIV ($COMP(pltr_data)); else pltrdt = $COMP(pltr_data); c_plcont ((const PLFLT **)ff, size_x, size_y, $kx(), $lx(), $ky(), $ly(), $P(clevel), $SIZE(nlevel), pltrcb, pltrdt); plFree2dGrid (ff, size_x, size_y);' ); # Surface mesh pp_def ('plmesh', NoPthread => 1, Pars => 'x(nx); y(ny); z(nx,ny); int opt()', GenericTypes => [D], Doc => 'FIXME: documentation here!', Code => ' int i, j, size_x, size_y; PLFLT** zz; size_x = $SIZE(nx); size_y = $SIZE(ny); plAlloc2dGrid (&zz, size_x, size_y); for (i = 0; i < size_x; i++) for (j = 0; j < size_y; j++) zz[i][j] = $z(nx => i, ny => j); c_plmesh ($P(x), $P(y), (const PLFLT **)zz, size_x, size_y, $opt()); plFree2dGrid (zz, size_x, size_y);' ); # Magnitude colored plot surface mesh with contour pp_def ('plmeshc', NoPthread => 1, Pars => 'x(nx); y(ny); z(nx,ny); int opt(); clevel(nlevel)', GenericTypes => [D], Doc => 'FIXME: documentation here!', Code => ' int i, j, size_x, size_y; PLFLT** zz; size_x = $SIZE(nx); size_y = $SIZE(ny); plAlloc2dGrid (&zz, size_x, size_y); for (i = 0; i < size_x; i++) for (j = 0; j < size_y; j++) zz[i][j] = $z(nx => i, ny => j); c_plmeshc ($P(x), $P(y), (const PLFLT **)zz, size_x, size_y, $opt(), $P(clevel), $SIZE(nlevel)); plFree2dGrid (zz, size_x, size_y);' ); # 3-d surface plot pp_def ('plot3d', NoPthread => 1, Pars => 'x(nx); y(ny); z(nx,ny); int opt(); int side()', GenericTypes => [D], Doc => 'FIXME: documentation here!', Code => ' int i, j, size_x, size_y; PLFLT** zz; size_x = $SIZE(nx); size_y = $SIZE(ny); plAlloc2dGrid (&zz, size_x, size_y); for (i = 0; i < size_x; i++) for (j = 0; j < size_y; j++) zz[i][j] = $z(nx => i, ny => j); c_plot3d ($P(x), $P(y), (const PLFLT **)zz, size_x, size_y, $opt(), $side()); plFree2dGrid (zz, size_x, size_y);' ); # Plots a 3-d representation of the function z[x][y] with contour pp_def ('plot3dc', NoPthread => 1, Pars => 'x(nx); y(ny); z(nx,ny); int opt(); clevel(nlevel)', GenericTypes => [D], Doc => 'FIXME: documentation here!', Code => ' int i, j, size_x, size_y; PLFLT** zz; size_x = $SIZE(nx); size_y = $SIZE(ny); plAlloc2dGrid (&zz, size_x, size_y); for (i = 0; i < size_x; i++) for (j = 0; j < size_y; j++) zz[i][j] = $z(nx => i, ny => j); c_plot3dc ($P(x), $P(y), (const PLFLT **)zz, size_x, size_y, $opt(), $P(clevel), $SIZE(nlevel)); plFree2dGrid (zz, size_x, size_y);' ); # Set color map1 colors using a piece-wise linear relationship pp_def ('plscmap1l', NoPthread => 1, Pars => 'int itype(); isty(n); coord1(n); coord2(n); coord3(n);' . ' int rev(nrev)', GenericTypes => [D], Doc => 'FIXME: documentation here!', Code => ' PLINT* rev; if ($SIZE(nrev) == 0) rev = NULL; else if ($SIZE(nrev) == $SIZE(n)) rev = $P(rev); else croak ("plscmap1l: rev must have either lenght == 0 or have the same length of the other input arguments"); c_plscmap1l ($itype(), $SIZE(n), $P(isty), $P(coord1), $P(coord2), $P(coord3), rev);' ); # Shade individual region on the basis of value pp_addpm ('sub plshade1 { plshade1_pp ($standard_order ? reorder ("plshade1", @_) : @_) }'); pp_add_exported('plshade1'); pp_def ('plshade1_pp', NoPthread => 1, GenericTypes => [D], Pars => 'a(nx,ny); left(); right(); bottom(); top(); shade_min();' . 'shade_max(); sh_cmap(); sh_color(); sh_width();' . 'min_color(); min_width(); max_color(); max_width();' . 'rectangular()', # 0-14 OtherPars => 'SV* defined; SV* pltr; SV* pltr_data;',# 15-17 Doc => 'FIXME: documentation here!', Code => ' int i, j, size_x, size_y; PLFLT* a; void (*pltrcb) (); PLPointer pltrdt; size_x = $SIZE(nx); size_y = $SIZE(ny); New(3, a, size_x * size_y, PLFLT); if(a == NULL) croak("Failed to allocate memory in plshade1_pp"); for (i = 0; i < size_x; i++) for (j = 0; j < size_y; j++) a[i * size_y + j] = (PLFLT) $a(nx => i, ny => j); defined_subroutine = $COMP(defined); check_sub_pointer (defined_subroutine, "plshade1: defined must be either 0 or a subroutine pointer"); pltr_subroutine = $COMP(pltr); check_sub_pointer (pltr_subroutine, "plshade1: pltr must be either 0 or a subroutine pointer"); pltrcb = get_standard_pltrcb ($COMP(pltr)); if (pltrcb != pltr_callback) pltrdt = (PLPointer) SvIV ($COMP(pltr_data)); else pltrdt = $COMP(pltr_data); c_plshade1 (a, size_x, size_y, SvTRUE ($COMP(defined)) ? defined_callback : NULL, $left(), $right(), $bottom(), $top(), $shade_min(), $shade_max(), $sh_cmap(), $sh_color(), $sh_width(), $min_color(), $min_width(), $max_color(), $max_width(), plfill, $rectangular(), pltrcb, pltrdt); Safefree (a);' ); # Plot gray-level image pp_def ('plimage', NoPthread => 1, GenericTypes => [D], Pars => 'idata(nx,ny); xmin(); xmax(); ymin(); ymax();' . 'zmin(); zmax(); Dxmin(); Dxmax(); Dymin(); Dymax();', Code => ' int i, j, size_x, size_y; PLFLT** idata; size_x = $SIZE(nx); size_y = $SIZE(ny); plAlloc2dGrid (&idata, size_x, size_y); for (i = 0; i < size_x; i++) for (j = 0; j < size_y; j++) idata[i][j] = $idata(nx => i, ny => j); plimage ((const PLFLT **)idata, size_x, size_y, $xmin(), $xmax(), $ymin(), $ymax(), $zmin(), $zmax(), $Dxmin(), $Dxmax(), $Dymin(), $Dymax()); plFree2dGrid (idata, size_x, size_y);' ); # Plot image with transformation if ($plversion->{'c_plimagefr'}) { pp_def ('plimagefr', NoPthread => 1, GenericTypes => [D], # plimagefr (idata, nx, ny, xmin, xmax, ymin, ymax, zmin, zmax, valuemin, valuemax, pltr, pltr_data); # plimagefr ($img, 0, $width, 0, $height, 0, 0, $img_min, $img_max, \&pltr2, $grid); Pars => 'idata(nx,ny); xmin(); xmax(); ymin(); ymax();' . 'zmin(); zmax(); valuemin(); valuemax();', # 0-8 OtherPars => 'SV* pltr; SV* pltr_data;', # 9,10 Code => ' int i, j, size_x, size_y; PLFLT** idata; void (*pltrcb) (); PLPointer pltrdt; size_x = $SIZE(nx); size_y = $SIZE(ny); pltr_subroutine = $COMP(pltr); check_sub_pointer (pltr_subroutine, "plimagefr: pltr must be either 0 or a subroutine pointer"); pltrcb = get_standard_pltrcb ($COMP(pltr)); if (pltrcb != pltr_callback) pltrdt = (PLPointer) SvIV ($COMP(pltr_data)); else pltrdt = $COMP(pltr_data); plAlloc2dGrid (&idata, size_x, size_y); for (i = 0; i < size_x; i++) for (j = 0; j < size_y; j++) idata[i][j] = $idata(nx => i, ny => j); c_plimagefr ((const PLFLT **)idata, size_x, size_y, $xmin(), $xmax(), $ymin(), $ymax(), $zmin(), $zmax(), $valuemin(), $valuemax(), (SvTRUE ($COMP(pltr)) ? pltrcb : NULL), (SvTRUE ($COMP(pltr)) ? pltrdt : NULL)); plFree2dGrid (idata, size_x, size_y);' ); } # $plversion->{'c_plimagefr'} # Set xor mode: # mode = 1-enter, 0-leave, status = 0 if not interactive device pp_addpm (<<'EOPM'); =head2 plxormod =for sig $status = plxormod ($mode) =for ref See the PLplot manual for reference. =cut EOPM pp_addxs (<<"EOC"); int plxormod (mode) int mode CODE: PLINT status; c_plxormod (mode, &status); RETVAL = status; OUTPUT: RETVAL EOC pp_add_exported ('plxormod'); # Wait for graphics input event and translate to world coordinates pp_addpm (<<'EOPM'); =head2 plGetCursor =for sig %gin = plGetCursor () =for ref plGetCursor waits for graphics input event and translate to world coordinates and returns a hash with the following keys: type: of event (CURRENTLY UNUSED) state: key or button mask keysym: key selected button: mouse button selected subwindow: subwindow (alias subpage, alias subplot) number string: translated string pX, pY: absolute device coordinates of pointer dX, dY: relative device coordinates of pointer wX, wY: world coordinates of pointer Returns an empty hash if no translation to world coordinates is possible. =cut EOPM pp_addxs (<<"EOC"); void plGetCursor () PPCODE: PLGraphicsIn gin; if (plGetCursor (&gin)) { EXTEND (SP, 24); PUSHs (sv_2mortal (newSVpv ("type", 0))); PUSHs (sv_2mortal (newSViv ((IV) gin.type))); PUSHs (sv_2mortal (newSVpv ("state", 0))); PUSHs (sv_2mortal (newSVuv ((UV) gin.state))); PUSHs (sv_2mortal (newSVpv ("keysym", 0))); PUSHs (sv_2mortal (newSVuv ((UV) gin.keysym))); PUSHs (sv_2mortal (newSVpv ("button", 0))); PUSHs (sv_2mortal (newSVuv ((UV) gin.button))); PUSHs (sv_2mortal (newSVpv ("subwindow", 0))); PUSHs (sv_2mortal (newSViv ((IV) gin.subwindow))); PUSHs (sv_2mortal (newSVpv ("string", 0))); PUSHs (sv_2mortal (newSVpv (gin.string, 0))); PUSHs (sv_2mortal (newSVpv ("pX", 0))); PUSHs (sv_2mortal (newSViv ((IV) gin.pX))); PUSHs (sv_2mortal (newSVpv ("pY", 0))); PUSHs (sv_2mortal (newSViv ((IV) gin.pY))); PUSHs (sv_2mortal (newSVpv ("dX", 0))); PUSHs (sv_2mortal (newSVnv ((double) gin.dX))); PUSHs (sv_2mortal (newSVpv ("dY", 0))); PUSHs (sv_2mortal (newSVnv ((double) gin.dY))); PUSHs (sv_2mortal (newSVpv ("wX", 0))); PUSHs (sv_2mortal (newSVnv ((double) gin.wX))); PUSHs (sv_2mortal (newSVpv ("wY", 0))); PUSHs (sv_2mortal (newSVnv ((double) gin.wY))); } EOC pp_add_exported ('plGetCursor'); pp_addpm (<<'EOPM'); =head2 plgstrm =for sig $strm = plgstrm () =for ref Returns the number of the current output stream. =cut EOPM pp_addxs (<<"EOC"); int plgstrm () CODE: PLINT strm; c_plgstrm (&strm); RETVAL = strm; OUTPUT: RETVAL EOC pp_add_exported ('plgstrm'); pp_addpm (<<'EOPM'); =head2 plgsdev =for sig $driver = plgdev () =for ref Returns the current driver name. =cut EOPM pp_addxs (<<"EOC"); char* plgdev () CODE: char driver[80]; c_plgdev (driver); RETVAL = driver; OUTPUT: RETVAL EOC pp_add_exported ('plgdev'); pp_addxs (<<"EOC"); char* plgfnam () CODE: char driver[80]; c_plgfnam (driver); RETVAL = driver; OUTPUT: RETVAL EOC pp_add_exported ('plgfnam'); pp_addpm (<<'EOPM'); =head2 plmkstrm =for sig $strm = plmkstrm () =for ref Creates a new stream and makes it the default. Returns the number of the created stream. =cut EOPM pp_addxs (<<"EOC"); int plmkstrm () CODE: PLINT strm; c_plmkstrm (&strm); RETVAL = strm; OUTPUT: RETVAL EOC pp_add_exported ('plmkstrm'); # Get the current library version number pp_addpm (<<'EOPM'); =head2 plgver =for sig $version = plgver () =for ref See the PLplot manual for reference. =cut EOPM pp_addxs (<<"EOC"); char* plgver () CODE: char ver[80]; c_plgver (ver); RETVAL = ver; OUTPUT: RETVAL EOC pp_add_exported ('plgver'); #---------------------------------------------------------------------------- pp_addpm ('sub plstripc { plstripc_pp ($standard_order ? reorder ("plstripc", @_) : @_) }'); pp_add_exported('plstripc'); pp_def ('plstripc_pp', NoPthread => 1, GenericTypes => [D], Pars => 'xmin(); xmax(); xjump(); ymin(); ymax();' . 'xlpos(); ylpos(); int y_ascl(); int acc();' . 'int colbox(); int collab();' . 'int colline(n); int styline(n); int [o] id()', # 0-12 OtherPars => 'char* xspec; char* yspec; SV* legline;' # 13-18 . 'char* labx; char* laby; char* labtop', Doc => 'FIXME: documentation here!', Code => ' I32 i; PLINT id; char* legline[4]; SV* sv_legline = $COMP(legline); AV* av_legline; if (! SvROK (sv_legline) || SvTYPE (SvRV (sv_legline)) != SVt_PVAV) croak ("plstripc: legline must be a reference to an array"); av_legline = (AV*) SvRV (sv_legline); if (av_len (av_legline) != 3) croak ("plstripc: legline must have four elements"); if ($SIZE(n) != 4) croak ("plstripc: colline and styline must have four elements"); for (i = 0; i < 4; i++) { SV** elem = av_fetch (av_legline, i, 0); legline[i] = (char *) SvPV_nolen (*elem); } c_plstripc (&id, $COMP(xspec), $COMP(yspec), $xmin(), $xmax(), $xjump(), $ymin(), $ymax(), $xlpos(), $ylpos(),$y_ascl(), $acc(), $colbox(), $collab(), $P(colline), $P(styline), (const char **)legline, $COMP(labx), $COMP(laby), $COMP(labtop)); $id() = (int) id;' ); #---------------------------------------------------------------------------- pp_def ('plgriddata', NoPthread => 1, GenericTypes => [D], Pars => 'x(npts); y(npts); z(npts); xg(nptsx); yg(nptsy);' . 'int type(); data(); [o] zg(nptsx,nptsy)', Doc => 'FIXME: documentation here!', Code => ' int i, j, size_x, size_y; PLFLT** zg; size_x = $SIZE(nptsx); size_y = $SIZE(nptsy); plAlloc2dGrid (&zg, size_x, size_y); c_plgriddata ($P(x), $P(y), $P(z), $SIZE(npts), $P(xg), size_x, $P(yg), size_y, zg, $type(), $data()); for (i = 0; i < size_x; i++) for (j = 0; j < size_y; j++) $zg(nptsx => i, nptsy => j) = zg[i][j]; plFree2dGrid (zg, size_x, size_y); ' ); #----------------------------------------------------------------------------- if ($plversion->{'c_pllegend'}) { # Draw an arc pp_addpm (<<'EOPM'); =head2 plarc =for sig plarc ($x, $y, $a, $b, $angle1, $angle2, $rotate, $fill); =for ref Draw a (possibly) filled arc centered at x, y with semimajor axis a and semiminor axis b, starting at angle1 and ending at angle2. See the PLplot manual for reference. =cut EOPM pp_addxs (<<"EOC"); int plarc (x, y, a, b, angle1, angle2, rotate, fill) double x double y double a double b double angle1 double angle2 double rotate int fill CODE: plarc (x, y, a, b, angle1, angle2, rotate, fill); EOC pp_add_exported ('plarc'); #---------------------------------------------------------------------------- pp_addpm (<<'EOPM'); =head2 plstransform =for sig plstransform ($subroutine_reference); =for ref Sets the default transformation routine for plotting. sub mapform { my ($x, $y) = @_; my $radius = 90.0 - $y; my $xp = $radius * cos ($x * pi / 180); my $yp = $radius * sin ($x * pi / 180); return ($xp, $yp); } plstransform (\&mapform); See the PLplot manual for more details. =cut EOPM pp_addxs (<<"EOC"); int plstransform (xform) SV* xform CODE: check_sub_pointer (xform_subroutine, "plstransform: xform must be either 0 or a subroutine pointer"); if (SvTRUE(xform)) xform_subroutine = SvRV(xform); plstransform (SvTRUE(xform) ? xform_callback : NULL, NULL); EOC pp_add_exported ('plstransform'); #---------------------------------------------------------------------------- pp_addpm (<<'EOPM'); =head2 plslabelfunc =for sig plslabelfunc ($subroutine_reference); =for ref # A custom axis labeling function for longitudes and latitudes. sub geolocation_labeler { my ($axis, $value, $length) = @_; my ($direction_label, $label_val); if (($axis == PL_Y_AXIS) && $value == 0) { return "Eq"; } elsif ($axis == PL_Y_AXIS) { $label_val = $value; $direction_label = ($label_val > 0) ? " N" : " S"; } elsif ($axis == PL_X_AXIS) { my $times = floor((abs($value) + 180.0 ) / 360.0); $label_val = ($value < 0) ? $value + 360.0 * $times : $value - 360.0 * $times; $direction_label = ($label_val > 0) ? " E" : ($label_val < 0) ? " W" : ""; } return substr (sprintf ("%.0f%s", abs($label_val), $direction_label), 0, $length); } plslabelfunc(\&geolocation_labeler); See the PLplot manual for more details. =cut EOPM # The PDL version of plslabelfunc only has one argument--the perl subroutine # to do the label translation: my $labeltext = perl_labelfunc($axis, $value, $length); # No 'data' argument is used, it is assumed that global data or a closure containing # necessary data can be used in 'perl_labelfunc'. pp_addxs (<<"EOC"); int plslabelfunc (labelfunc) SV* labelfunc CODE: check_sub_pointer (labelfunc_subroutine, "plslabelfunc: labelfunc must be either 0 or a subroutine pointer"); if (SvTRUE(labelfunc)) labelfunc_subroutine = SvRV(labelfunc); plslabelfunc (SvTRUE(labelfunc) ? labelfunc_callback : NULL, NULL); EOC pp_add_exported ('plslabelfunc'); #---------------------------------------------------------------------------- pp_addpm (<<'EOPM'); =head2 pllegend =for sig my ($legend_width, $legend_height) = pllegend ($position, $opt, $x, $y, $plot_width, $bg_color, $nlegend, \@opt_array, $text_offset, $text_scale, $text_spacing, $test_justification, \@text_colors, \@text, \@box_colors, \@box_patterns, \@box_scales, \@line_colors, \@line_styles, \@line_widths, \@symbol_colors, \@symbol_scales, \@symbol_numbers, \@symbols); =for ref See the PLplot manual for more details. =cut EOPM pp_addxs (<<"EOC"); void pllegend(opt, position, x, y, plot_width, bg_color, bb_color, bb_style, nrow, ncolumn, nlegend, opt_array_rv, text_offset, text_scale, text_spacing, text_justification, text_colors_rv, text_rv, box_colors_rv, box_patterns_rv, box_scales_rv, box_line_widths_rv, line_colors_rv, line_styles_rv, line_widths_rv, symbol_colors_rv, symbol_scales_rv, symbol_numbers_rv, symbols_rv) int opt int position double x double y double plot_width int bg_color int bb_color int bb_style int nrow int ncolumn int nlegend SV* opt_array_rv double text_offset double text_scale double text_spacing double text_justification SV* text_colors_rv SV* text_rv SV* box_colors_rv SV* box_patterns_rv SV* box_scales_rv SV* box_line_widths_rv SV* line_colors_rv SV* line_styles_rv SV* line_widths_rv SV* symbol_colors_rv SV* symbol_scales_rv SV* symbol_numbers_rv SV* symbols_rv PPCODE: int i; double p_legend_width; double p_legend_height; int opt_array[nlegend]; int text_colors[nlegend]; char *text[nlegend]; int box_colors[nlegend]; int box_patterns[nlegend]; double box_scales[nlegend]; int box_line_widths[nlegend]; int line_colors[nlegend]; int line_styles[nlegend]; int line_widths[nlegend]; int symbol_colors[nlegend]; double symbol_scales[nlegend]; int symbol_numbers[nlegend]; char *symbols[nlegend]; SV **elem; for (i = 0; i < nlegend; i++) { elem = av_fetch((AV *)SvRV(opt_array_rv), i, 0); opt_array[i] = (int)SvIV(*elem); elem = av_fetch((AV *)SvRV(text_colors_rv), i, 0); text_colors[i] = (int)SvIV(*elem); elem = av_fetch((AV *)SvRV(text_rv), i, 0); text[i] = (char *)SvPV_nolen(*elem); box_colors[i] = 0; if (SvROK(box_colors_rv)) { elem = av_fetch((AV *)SvRV(box_colors_rv), i, 0); if (elem && SvOK(*elem)) { box_colors[i] = (int)SvIV(*elem); } } box_patterns[i] = 0; if (SvROK(box_patterns_rv)) { elem = av_fetch((AV *)SvRV(box_patterns_rv), i, 0); if (elem && SvOK(*elem)) { box_patterns[i] = (int)SvIV(*elem); } } box_scales[i] = 0.0; if (SvROK(box_scales_rv)) { elem = av_fetch((AV *)SvRV(box_scales_rv), i, 0); if (elem && SvOK(*elem)) { box_scales[i] = (double)SvNV(*elem); } } box_line_widths[i] = 0; if (SvROK(box_line_widths_rv)) { elem = av_fetch((AV *)SvRV(box_line_widths_rv), i, 0); if (elem && SvOK(*elem)) { box_line_widths[i] = (int)SvIV(*elem); } } line_colors[i] = 0; if (SvROK(line_colors_rv)) { elem = av_fetch((AV *)SvRV(line_colors_rv), i, 0); if (elem && SvOK(*elem)) { line_colors[i] = (int)SvIV(*elem); } } line_styles[i] = 0; if (SvROK(line_styles_rv)) { elem = av_fetch((AV *)SvRV(line_styles_rv), i, 0); if (elem && SvOK(*elem)) { line_styles[i] = (int)SvIV(*elem); } } line_widths[i] = 0; if (SvROK(line_widths_rv)) { elem = av_fetch((AV *)SvRV(line_widths_rv), i, 0); if (elem && SvOK(*elem)) { line_widths[i] = (int)SvIV(*elem); } } symbol_colors[i] = 0; if (SvROK(symbol_colors_rv)) { elem = av_fetch((AV *)SvRV(symbol_colors_rv), i, 0); if (elem && SvOK(*elem)) { symbol_colors[i] = (int)SvIV(*elem); } } symbol_scales[i] = 0.0; if (SvROK(symbol_scales_rv)) { elem = av_fetch((AV *)SvRV(symbol_scales_rv), i, 0); if (elem && SvOK(*elem)) { symbol_scales[i] = (double)SvNV(*elem); } } symbol_numbers[i] = 0; if (SvROK(symbol_numbers_rv)) { elem = av_fetch((AV *)SvRV(symbol_numbers_rv), i, 0); if (elem && SvOK(*elem)) { symbol_numbers[i] = (int)SvIV(*elem); } } symbols[i] = "0"; if (SvROK(symbols_rv)) { elem = av_fetch((AV *)SvRV(symbols_rv), i, 0); if (elem && SvOK(*elem)) { symbols[i] = (char *)SvPV_nolen(*elem); } } } pllegend(&p_legend_width, &p_legend_height, opt, position, x, y, plot_width, bg_color, bb_color, bb_style, nrow, ncolumn, nlegend, opt_array, text_offset, text_scale, text_spacing, text_justification, text_colors, (const char **)text, box_colors, box_patterns, box_scales, box_line_widths, line_colors, line_styles, line_widths, symbol_colors, symbol_scales, symbol_numbers, (const char **)symbols); EXTEND (SP, 2); PUSHs (sv_2mortal (newSVnv (p_legend_width))); PUSHs (sv_2mortal (newSVnv (p_legend_height))); EOC pp_add_exported ('pllegend'); #---------------------------------------------------------------------------- pp_addpm (<<'EOPM'); =head2 plspal0 =for sig plspal0($filename); =for ref Set color palette 0 from the input .pal file. See the PLplot manual for more details. =cut EOPM # The PDL version of plspal0 pp_addxs (<<"EOC"); int plspal0 (filename) char* filename PPCODE: plspal0((const char *)filename); EOC pp_add_exported ('plspal0'); #---------------------------------------------------------------------------- # The PDL version of plspal1 pp_addpm (<<'EOPM'); =head2 plspal1 =for sig plspal1($filename); =for ref Set color palette 1 from the input .pal file. See the PLplot manual for more details. =cut EOPM pp_addxs (<<"EOC"); int plspal1 (filename, interpolate) char* filename int interpolate PPCODE: plspal1((const char *)filename, (PLBOOL)interpolate); EOC pp_add_exported ('plspal1'); } # if $plversion->{'c_pllegend'} #----------------------------------------------------------------------------------------- pp_addpm (<<'EOPM'); =head2 plbtime =for sig my ($year, $month, $day, $hour, $min, $sec) = plbtime($ctime); =for ref Calculate broken-down time from continuous time for current stream. =cut EOPM pp_addxs (<<"EOC"); void plbtime (ctime) double ctime PPCODE: PLINT year; PLINT month; PLINT day; PLINT hour; PLINT min; PLFLT sec; c_plbtime(&year, &month, &day, &hour, &min, &sec, ctime); EXTEND (SP, 6); PUSHs (sv_2mortal (newSViv (year))); PUSHs (sv_2mortal (newSViv (month))); PUSHs (sv_2mortal (newSViv (day))); PUSHs (sv_2mortal (newSViv (hour))); PUSHs (sv_2mortal (newSViv (min))); PUSHs (sv_2mortal (newSVnv (sec))); EOC pp_add_exported ('plbtime'); #----------------------------------------------------------------------------------------- pp_addpm (<<'EOPM'); =head2 plconfigtime =for sig plconfigtime($scale, $offset1, $offset2, $ccontrol, $ifbtime_offset, $year, $month, $day, $hour, $min, $sec); =for ref Configure transformation between continuous and broken-down time (and vice versa) for current stream. =cut EOPM pp_addxs (<<"EOC"); void plconfigtime(scale, offset1, offset2, ccontrol, ifbtime_offset, year, month, day, hour, min, sec) double scale double offset1 double offset2 int ccontrol int ifbtime_offset int year int month int day int hour int min double sec PPCODE: c_plconfigtime((PLFLT) scale, (PLFLT) offset1, (PLFLT) offset2, (PLINT) ccontrol, (PLBOOL) ifbtime_offset, (PLINT) year, (PLINT) month, (PLINT) day, (PLINT) hour, (PLINT) min, (PLFLT) sec); EOC pp_add_exported ('plconfigtime'); #----------------------------------------------------------------------------------------- pp_addpm (<<'EOPM'); =head2 plctime =for sig my $ctime = plctime($year, $month, $day, $hour, $min, $sec); =for ref Calculate continuous time from broken-down time for current stream. =cut EOPM pp_addxs (<<"EOC"); void plctime(year, month, day, hour, min, sec) int year int month int day int hour int min double sec PPCODE: PLFLT ctime; c_plctime(year, month, day, hour, min, sec, &ctime); EXTEND (SP, 1); PUSHs (sv_2mortal (newSVnv (ctime))); EOC pp_add_exported ('plctime'); #----------------------------------------------------------------------------------------- pp_addpm (<<'EOPM'); =head2 pltimefmt =for sig pltimefmt($fmt); =for ref Set format for date / time labels. See the PLplot manual for more details. =cut EOPM pp_addxs (<<"EOC"); void pltimefmt(fmt) char *fmt PPCODE: c_pltimefmt((const char *)fmt); EOC pp_add_exported ('pltimefmt'); #----------------------------------------------------------------------------------------- pp_addpm (<<'EOPM'); =head2 plsesc =for sig plsesc($esc); =for ref Set the escape character for text strings. See the PLplot manual for more details. =cut EOPM pp_addxs (<<"EOC"); void plsesc (esc) SV* esc PPCODE: char *esc_c; esc_c = (char *)SvPV_nolen(esc); c_plsesc((char)*esc_c); EOC pp_add_exported ('plsesc'); if ($plversion->{'plsvect'}) { # Vector field plots pp_def ('plvect', NoPthread => 1, GenericTypes => [D], Pars => 'u(nx,ny); v(nx,ny); scale();', OtherPars => 'SV* pltr; SV* pltr_data;', Doc => 'FIXME: documentation here!', Code => ' int i, j, size_x, size_y; PLFLT** u; PLFLT** v; void (*pltrcb) (); PLPointer pltrdt; size_x = $SIZE(nx); size_y = $SIZE(ny); plAlloc2dGrid (&u, size_x, size_y); plAlloc2dGrid (&v, size_x, size_y); for (i = 0; i < size_x; i++) for (j = 0; j < size_y; j++) { u[i][j] = $u(nx => i, ny => j); v[i][j] = $v(nx => i, ny => j); } pltr_subroutine = $COMP(pltr); check_sub_pointer (pltr_subroutine, "plvect: pltr must be either 0 or a subroutine pointer"); pltrcb = get_standard_pltrcb ($COMP(pltr)); if (pltrcb != pltr_callback) pltrdt = (PLPointer) SvIV ($COMP(pltr_data)); else pltrdt = $COMP(pltr_data); plvect ((const PLFLT **)u, (const PLFLT **)v, size_x, size_y, $scale(), pltrcb, pltrdt); plFree2dGrid (u, size_x, size_y); plFree2dGrid (v, size_x, size_y);' ); pp_def ('plsvect', NoPthread => 1, Pars => 'arrowx(npts); arrowy(npts); int fill()', GenericTypes => [D], Code => 'c_plsvect ($P(arrowx), $P(arrowy), $SIZE(npts), $fill());' ); } pp_def ('plhlsrgb', NoPthread => 1, GenericTypes => [D], Pars => 'double h();double l();double s();double [o]p_r();double [o]p_g();double [o]p_b()', Code => 'c_plhlsrgb($h(),$l(),$s(),$P(p_r),$P(p_g),$P(p_b));' ); # void c_plgcol0(PLINT icol0, PLINT *r, PLINT *g, PLINT *b); pp_def ('plgcol0', NoPthread => 1, Pars => 'int icolzero(); int [o]r(); int [o]g(); int [o]b()', Code => 'c_plgcol0($icolzero(),$P(r),$P(g),$P(b));' ); # void c_plgcolbg(PLINT *r, PLINT *g, PLINT *b); pp_def ('plgcolbg', NoPthread => 1, Pars => 'int [o]r(); int [o]g(); int [o]b()', Code => 'c_plgcolbg($P(r),$P(g),$P(b));' ); # void c_plscmap0(PLINT *r, PLINT *g, PLINT *b, PLINT ncol0); pp_def ('plscmap0', NoPthread => 1, Pars => 'int r(n); int g(n); int b(n)', Code => 'c_plscmap0($P(r),$P(g),$P(b), $SIZE(n));' ); # void c_plscmap1(PLINT *r, PLINT *g, PLINT *b, PLINT ncol1); pp_def ('plscmap1', NoPthread => 1, Pars => 'int r(n); int g(n); int b(n)', Code => 'c_plscmap1($P(r),$P(g),$P(b), $SIZE(n));' ); if (!$noalpha) { # void c_plgcol0a(PLINT icol0, PLINT *r, PLINT *g, PLINT *b, PLFLT *a); pp_def ('plgcol0a', NoPthread => 1, Pars => 'int icolzero(); int [o]r(); int [o]g(); int [o]b(); double [o]a()', Code => 'c_plgcol0a($icolzero(),$P(r),$P(g),$P(b),$P(a));' ); # void c_plgcolbga(PLINT *r, PLINT *g, PLINT *b, PLFLT *a); pp_def ('plgcolbga', NoPthread => 1, Pars => 'int [o]r(); int [o]g(); int [o]b(); double [o]a()', Code => 'c_plgcolbga($P(r),$P(g),$P(b),$P(a));' ); # void c_plscmap0a(PLINT *r, PLINT *g, PLINT *b, PLFLT *a, PLINT ncol0); pp_def ('plscmap0a', NoPthread => 1, Pars => 'int r(n); int g(n); int b(n); double a(n)', Code => 'c_plscmap0a($P(r),$P(g),$P(b),$P(a),$SIZE(n));' ); # void c_plscmap1a(PLINT *r, PLINT *g, PLINT *b, PLFLT *a, PLINT ncol1); pp_def ('plscmap1a', NoPthread => 1, Pars => 'int r(n); int g(n); int b(n); double a(n)', Code => 'c_plscmap1a($P(r),$P(g),$P(b),$P(a),$SIZE(n));' ); # Set color map1 colors using a piece-wise linear relationship, include alpha channel pp_def ('plscmap1la', NoPthread => 1, Pars => 'int itype(); isty(n); coord1(n); coord2(n); coord3(n); coord4(n);' . ' int rev(nrev)', GenericTypes => [D], Doc => 'FIXME: documentation here!', Code => ' PLINT* rev; if ($SIZE(nrev) == 0) rev = NULL; else if ($SIZE(nrev) == $SIZE(n)) rev = $P(rev); else croak ("plscmap1la: rev must have either length == 0 or have the same length of the other input arguments"); c_plscmap1la ($itype(), $SIZE(n), $P(isty), $P(coord1), $P(coord2), $P(coord3), $P(coord4), rev);' ); } # ## UNICODE font manipulation # if ($plversion->{'c_plsfont'}) { # plgfont(PLINT *p_family, PLINT *p_style, PLINT *p_weight); pp_def ('plgfont', NoPthread => 1, Pars => 'int [o]p_family(); int [o]p_style(); int [o]p_weight();', Code => 'c_plgfont($P(p_family),$P(p_style),$P(p_weight));' ); # plsfont (PLINT family, PLINT style, PLINT weight); pp_def ('plsfont', NoPthread => 1, Pars => 'int family(); int style(); int weight();', Code => 'c_plsfont($family(),$style(),$weight());' ); } # $plversion->{'c_plsfont'} # plcalc_world (PLFLT rx, PLFLT ry, PLFLT *wx, PLFLT *wy, PLINT *window); pp_def ('plcalc_world', NoPthread => 1, Pars => 'double rx(); double ry(); double [o]wx(); double [o]wy(); int [o]window()', Code => 'c_plcalc_world($rx(), $ry(), $P(wx), $P(wy), $P(window));' ); pp_addxs (<<"EOC"); unsigned int plgfci () CODE: { unsigned int RETVAL; unsigned int fci; c_plgfci(&fci); RETVAL = fci; XSprePUSH; PUSHu((UV)RETVAL); } XSRETURN(1); EOC pp_add_exported('', 'plgfci'); pp_addxs (<<'EOC'); void plsfci(fci) unsigned int fci CODE: c_plsfci(fci); EOC pp_add_exported('', 'plsfci'); pp_addpm (<<'EOPM'); =pod =head1 WARNINGS AND ERRORS PLplot gives many errors and warnings. Some of these are given by the PDL interface while others are internal PLplot messages. Below are some of these messages, and what you need to do to address them: =over =item * Box must be a ref to a four element array When specifying a box, you must pass a reference to a four-element array, or use an anonymous four-element array. # Gives trouble: $pl->xyplot($x, $y, BOX => (0, 0, 100, 200) ); # What you meant to say was: $pl->xyplot($x, $y, BOX => [0, 0, 100, 200] ); =item * Too many colors used! (max 15) =back =head1 AUTHORS Doug Hunt Rafael Laboissiere David Mertens =head1 SEE ALSO perl(1), PDL(1), L The other common graphics packages include L and L. =cut EOPM pp_done(); # Local Variables: # mode: cperl # End: