PLplot  5.10.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros
sfstubsf95.f90
Go to the documentation of this file.
1 !***********************************************************************
2 ! sfstubsf95.f
3 !
4 ! Copyright (C) 2005, 2006 Arjen Markus
5 ! Copyright (C) 2006-2014 Alan W. Irwin
6 !
7 ! This file is part of PLplot.
8 !
9 ! PLplot is free software; you can redistribute it and/or modify
10 ! it under the terms of the GNU Library General Public License as published
11 ! by the Free Software Foundation; either version 2 of the License, or
12 ! (at your option) any later version.
13 !
14 ! PLplot is distributed in the hope that it will be useful,
15 ! but WITHOUT ANY WARRANTY; without even the implied warranty of
16 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 ! GNU Library General Public License for more details.
18 !
19 ! You should have received a copy of the GNU Library General Public License
20 ! along with PLplot; if not, write to the Free Software
21 ! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
22 !
23 !
24 ! This file contains the interfaces for Fortran 95:
25 ! - it includes the actual FORTRAN routines from the FORTRAN 95 bindings
26 ! - it includes interfaces to the C routines from these bindings
27 ! - it defines a few Fortran 95 specific items and interfaces
28 !
29 ! NB
30 ! This module is written in fixed form.
31 ! To enable a redefinition of certain interfaces, we actually have
32 ! two modules.
33 !
34 ! NB
35 ! The INTENT attribute is not currently used. This is a matter to
36 ! be looked into.
37 !
38 ! NB
39 ! It is assumed in the current implementation that all arrays are
40 ! passed with correct dimensions. It would be wiser, perhaps, to
41 ! use the minimum dimension instead of just the dimension of the
42 ! first array.
43 !
44 ! NOTE:
45 ! Some of the C routines will have to be renamed (via macros)
46 !
47 !***********************************************************************
48 
49 !
50 ! Parameters for identifying the kind of PLplot's real
51 ! numbers (a configuration parameter)
52 ! Use whatever name suits you better.
53 !
54 module plplot_flt
55  include 'plflt.inc'
56 end module
57 
58 !
59 ! Parameters and variables for strings / arrays for
60 ! string conversion
61 !
62 module plplot_str
63  integer :: maxleni
64  parameter(maxlen = 320)
65  parameter(maxleni = 80)
66  character (len = maxlen) :: string1, string2, string3
67  character (len = maxlen) :: string4, string5, string6
68  character (len = maxlen) :: string7, string8, string9
69  integer, dimension(maxleni) :: s1, s2, s3, s4, s5, s6, s7, s8, s9
70 
71  character(len=1), parameter :: PL_END_OF_STRING = achar(0)
72 end module
73 
74 module plplotp
75  use plplot_flt
76  use plplot_str
77  use plplot_strutils
78  implicit none
79 
80  interface plcont
81  module procedure plcontour_0
82  module procedure plcontour_1
83  module procedure plcontour_2
84  module procedure plcontour_tr
85  module procedure plcontour_0_all
86  module procedure plcontour_1_all
87  module procedure plcontour_2_all
88  module procedure plcontour_tr_all
89  end interface
92 
93  interface plvect
94  module procedure plvectors_0
95  module procedure plvectors_1
96  module procedure plvectors_2
97  module procedure plvectors_tr
98  end interface
100 
101  interface plshade
102  module procedure plshade_single_0
103  module procedure plshade_single_1
104  module procedure plshade_single_2
105  module procedure plshade_single_tr
106  end interface
108 
109  interface plshades
110  module procedure plshades_multiple_0
111  module procedure plshades_multiple_1
112  module procedure plshades_multiple_1r
113  module procedure plshades_multiple_2
114  module procedure plshades_multiple_2r
115  module procedure plshades_multiple_tr
116  module procedure plshades_multiple_trr
117  end interface
121 
122  interface plimagefr
123  module procedure plimagefr_0
124  module procedure plimagefr_1
125  module procedure plimagefr_2
126  module procedure plimagefr_tr
127  end interface
129 
130 contains
131  include 'sfstubs.f90'
132 end module plplotp
133 
135  use plplot_flt
136  type :: plgraphicsin
137  integer type ! of event (CURRENTLY UNUSED)
138  integer state ! key or button mask
139  integer keysym ! key selected
140  integer button ! mouse button selected
141  integer subwindow ! subwindow (alias subpage, alias subplot) number
142  character(len=16) string ! translated string
143  integer pX, pY ! absolute device coordinates of pointer
144  real(kind=plflt) dX, dY ! relative device coordinates of pointer
145  real(kind=plflt) wX, wY ! world coordinates of pointer
146  end type plgraphicsin
147 end module plplot_types
148 
149 module plplot
150  use plplotp
151  use plplot_flt
152  use plplot_types
153  use plplot_strutils
154  !
155  ! To be added: renaming list
156  !
157 
158  implicit none
159  include 'plplot_parameters.h'
160 
161  !
162  ! To be added: alternative interfaces
163  !
164  interface
165  subroutine pladv( sub )
166  integer :: sub
167  end subroutine pladv
168  end interface
169 
170  interface plbin
171  module procedure plbin
172  end interface
173 
174  interface
175  subroutine plbop
176  end subroutine plbop
177  end interface
178 
179  interface
180  subroutine plcalc_world( rx, ry, wx, wy, window )
181  use plplot_flt
182  real(kind=plflt) :: rx, ry, wx, wy
183  integer :: window
184  end subroutine plcalc_world
185  end interface
186 
187  interface
188  subroutine plclear
189  end subroutine plclear
190  end interface
191 
192  interface
193  subroutine plcol0( icol )
194  integer :: icol
195  end subroutine plcol0
196  end interface
197 
198  interface
199  subroutine plcol1( col )
200  use plplot_flt
201  real(kind=plflt) :: col
202  end subroutine plcol1
203  end interface
204 
205  interface plcolorbar
206  module procedure plcolorbar_1
207  module procedure plcolorbar_2
208  end interface
209 
210  interface plcpstrm
211  module procedure plcpstrm
212  end interface
213 
214  interface
215  subroutine plend
216  end subroutine plend
217  end interface
218 
219  interface
220  subroutine plend1
221  end subroutine plend1
222  end interface
223 
224  interface
225  subroutine plenv( xmin, xmax, ymin, ymax, just, axis )
226  use plplot_flt
227  real(kind=plflt) :: xmin, xmax, ymin, ymax
228  integer :: just, axis
229  end subroutine plenv
230  end interface
231 
232  interface
233  subroutine pleop
234  end subroutine pleop
235  end interface
236 
237  interface plerrx
238  module procedure plerrx
239  end interface
240 
241  interface plerry
242  module procedure plerry
243  end interface
244 
245  interface plfamadv
246  subroutine plfamadv
247  end subroutine plfamadv
248  end interface
249 
250  interface plfill
251  module procedure plfill
252  end interface
253 
254  interface plfill3
255  module procedure plfill3
256  end interface
257 
258  interface
259  subroutine plflush
260  end subroutine plflush
261  end interface
262 
263  interface
264  subroutine plfont( font )
265  integer :: font
266  end subroutine plfont
267  end interface
268 
269  interface
270  subroutine plfontld( charset )
271  integer :: charset
272  end subroutine plfontld
273  end interface
274 
275  interface
276  subroutine plgchr( chrdef, chrht )
277  use plplot_flt
278  real(kind=plflt) :: chrdef, chrht
279  end subroutine plgchr
280  end interface
281 
282  interface
283  subroutine plgcmap1_range( min_color, max_color )
284  use plplot_flt
285  real(kind=plflt) :: min_color, max_color
286  end subroutine plgcmap1_range
287  end interface
288 
289  interface
290  subroutine plgcol0( icol, r, g, b )
291  integer :: icol, r, g, b
292  end subroutine plgcol0
293  end interface
294 
295  interface
296  subroutine plgcol0a( icol, r, g, b, a )
297  use plplot_flt
298  integer :: icol, r, g, b
299  real(kind=plflt) :: a
300  end subroutine plgcol0a
301  end interface
302 
303  interface
304  subroutine plgcolbg( r, g, b )
305  integer :: r, g, b
306  end subroutine plgcolbg
307  end interface
308 
309  interface
310  subroutine plgcolbga( r, g, b, a )
311  use plplot_flt
312  integer :: r, g, b
313  real(kind=plflt) :: a
314  end subroutine plgcolbga
315  end interface
316 
317  interface
318  subroutine plgcompression( compression )
319  integer :: compression
320  end subroutine plgcompression
321  end interface
322 
323  interface
324  subroutine plgdidev( mar, aspect, jx, jy )
325  use plplot_flt
326  real(kind=plflt) :: mar, aspect, jx, jy
327  end subroutine plgdidev
328  end interface
329 
330  interface
331  subroutine plgdiori( rot )
332  use plplot_flt
333  real(kind=plflt) :: rot
334  end subroutine plgdiori
335  end interface
336 
337  interface
338  subroutine plgdiplt( xmin, xmax, ymin, ymax )
339  use plplot_flt
340  real(kind=plflt) :: xmin, xmax, ymin, ymax
341  end subroutine plgdiplt
342  end interface
343 
344  interface
345  subroutine plgetcursor( gin )
346  use plplot_flt
347  use plplot_types
348  type(plgraphicsin) :: gin
349  end subroutine plgetcursor
350  end interface
351 
352  interface
353  subroutine plgfam( fam, num, bmax )
354  integer :: fam, num, bmax
355  end subroutine plgfam
356  end interface
357 
358  interface
359  subroutine plgfci( fci )
360  use plplot_flt
361  integer(kind=plunicode) :: fci
362  end subroutine plgfci
363  end interface
364 
365  interface
366  subroutine plgfont( family, style, weight )
367  integer :: family, style, weight
368  end subroutine plgfont
369  end interface
370 
371  interface
372  subroutine plglevel( level )
373  integer :: level
374  end subroutine plglevel
375  end interface
376 
377  interface
378  subroutine plgpage( xpmm, ypmm, xwid, ywid, xoff, yoff )
379  use plplot_flt
380  real(kind=plflt) :: xpmm, ypmm
381  integer :: xwid, ywid, xoff, yoff
382  end subroutine plgpage
383  end interface
384 
385  interface
386  subroutine plgra
387  end subroutine plgra
388  end interface
389 
390  interface plgradient
391  module procedure plgradient
392  end interface
393 
394  interface plgriddata
395  module procedure plgriddata
396  end interface
397 
398  interface
399  subroutine plgspa( xmin, xmax, ymin, ymax )
400  use plplot_flt
401  real(kind=plflt) :: xmin, xmax, ymin, ymax
402  end subroutine plgspa
403  end interface
404 
405  interface
406  subroutine plgstrm( strm )
407  integer :: strm
408  end subroutine plgstrm
409  end interface
410 
411  interface
412  subroutine plgvpd( xmin, xmax, ymin, ymax )
413  use plplot_flt
414  real(kind=plflt) :: xmin, xmax, ymin, ymax
415  end subroutine plgvpd
416  end interface
417 
418  interface
419  subroutine plgvpw( xmin, xmax, ymin, ymax )
420  use plplot_flt
421  real(kind=plflt) :: xmin, xmax, ymin, ymax
422  end subroutine plgvpw
423  end interface
424 
425  interface
426  subroutine plgxax( digmax, digits )
427  integer :: digmax, digits
428  end subroutine plgxax
429  end interface
430 
431  interface
432  subroutine plgyax( digmax, digits )
433  integer :: digmax, digits
434  end subroutine plgyax
435  end interface
436 
437  interface
438  subroutine plgzax( digmax, digits )
439  integer :: digmax, digits
440  end subroutine plgzax
441  end interface
442 
443  interface plhist
444  module procedure plhist
445  end interface
446 
447  interface
448  subroutine plhls( h, l, s )
449  use plplot_flt
450  real(kind=plflt) :: h, l, s
451  end subroutine plhls
452  end interface
453 
454  interface
455  subroutine plhlsrgb( h, l, s, r, g, b )
456  use plplot_flt
457  real(kind=plflt) :: h, l, s, r, g, b
458  end subroutine plhlsrgb
459  end interface
460 
461  interface
462  subroutine plinit
463  end subroutine plinit
464  end interface
465 
466  interface
467  subroutine pljoin( x1, y1, x2, y2 )
468  use plplot_flt
469  real(kind=plflt) :: x1, y1, x2, y2
470  end subroutine pljoin
471  end interface
472 
473  interface
474  subroutine pllightsource( x, y, z )
475  use plplot_flt
476  real(kind=plflt) :: x, y, z
477  end subroutine pllightsource
478  end interface
479 
480  interface pllegend
481  module procedure pllegend_1
482  module procedure pllegend_2
483  end interface
484 
485  interface plline
486  module procedure plline
487  end interface
488 
489  interface plline3
490  module procedure plline3
491  end interface
492 
493  interface pllsty
494  subroutine pllsty( lin )
495  integer :: lin
496  end subroutine pllsty
497  end interface
498 
499  interface plmap
500  module procedure plmap1, plmap2
501  end interface plmap
502 
503  interface plmeridians
504  module procedure plmeridians1, plmeridians2
505  end interface plmeridians
506 
507  interface plmesh
508  module procedure plmesh
509  end interface
510 
511  interface plmeshc
512  module procedure plmeshc
513  end interface
514 
515  interface
516  subroutine plmkstrm( strm )
517  integer :: strm
518  end subroutine plmkstrm
519  end interface
520 
521  interface
522  subroutine plpat( nlin, inc, del )
523  integer :: nlin, inc, del
524  end subroutine plpat
525  end interface
526 
527  interface
528  subroutine plpath( n, x1, y1, x2, y2 )
529  use plplot_flt
530  integer :: n
531  real(kind=plflt) :: x1, y1, x2, y2
532  end subroutine plpath
533  end interface
534 
535  interface plot3d
536  module procedure plot3d
537  end interface
538 
539  interface plot3dc
540  module procedure plot3dc
541  end interface
542 
543  interface plpoin
544  module procedure plpoin
545  end interface
546 
547  interface plpoin3
548  module procedure plpoin3
549  end interface
550 
551  interface plpoly3
552  module procedure plpoly3
553  end interface
554 
555  interface
556  subroutine plprec( setp, prec )
557  integer :: setp, prec
558  end subroutine plprec
559  end interface
560 
561  interface
562  subroutine plpsty( patt )
563  integer :: patt
564  end subroutine plpsty
565  end interface
566 
567  interface
568  subroutine plreplot
569  end subroutine plreplot
570  end interface
571 
572  !
573  ! Note: plrgb and plrgb1 can be merged
574  !
575  interface
576  subroutine plrgb( r, g, b )
577  use plplot_flt
578  real(kind=plflt) :: r, g, b
579  end subroutine plrgb
580  end interface
581 
582  interface
583  subroutine plrgb1( r, g, b )
584  integer :: r, g, b
585  end subroutine plrgb1
586  end interface
587 
588  interface
589  subroutine plrgbhls( r, g, b, h, l, s )
590  use plplot_flt
591  real(kind=plflt) :: r, g, b, h, l, s
592  end subroutine plrgbhls
593  end interface
594 
595  interface
596  subroutine plschr( chrdef, chrht )
597  use plplot_flt
598  real(kind=plflt) :: chrdef, chrht
599  end subroutine plschr
600  end interface
601 
602  interface plscmap0
603  module procedure plscmap0
604  end interface
605 
606  interface plscmap0a
607  module procedure plscmap0a
608  end interface
609 
610  interface
611  subroutine plscmap0n( n )
612  integer :: n
613  end subroutine plscmap0n
614  end interface
615 
616  interface plscmap1
617  module procedure plscmap1
618  end interface
619 
620  interface plscmap1a
621  module procedure plscmap1a
622  end interface
623 
624  interface plscmap1l
625  module procedure plscmap1l
626  module procedure plscmap1l2
627  end interface
628 
629  interface plscmap1la
630  module procedure plscmap1la
631  module procedure plscmap1la2
632  end interface
633 
634  interface
635  subroutine plscmap1n( n )
636  integer :: n
637  end subroutine plscmap1n
638  end interface
639 
640  interface
641  subroutine plscmap1_range( min_color, max_color )
642  use plplot_flt
643  real(kind=plflt) :: min_color, max_color
644  end subroutine plscmap1_range
645  end interface
646 
647  interface
648  subroutine plscol0( icol, r, g, b )
649  integer :: icol, r, g, b
650  end subroutine plscol0
651  end interface
652 
653  interface
654  subroutine plscol0a( icol, r, g, b, a )
655  use plplot_flt
656  integer :: icol, r, g, b
657  real(kind=plflt) :: a
658  end subroutine plscol0a
659  end interface
660 
661  interface
662  subroutine plscolbg( r, g, b )
663  integer :: r, g, b
664  end subroutine plscolbg
665  end interface
666 
667  interface
668  subroutine plscolbga( r, g, b, a )
669  use plplot_flt
670  integer :: r, g, b
671  real(kind=plflt) :: a
672  end subroutine plscolbga
673  end interface
674 
675  interface
676  subroutine plscolor( color )
677  integer :: color
678  end subroutine plscolor
679  end interface
680 
681  interface
682  subroutine plscompression( compression )
683  integer :: compression
684  end subroutine plscompression
685  end interface
686 
687  interface
688  subroutine plsdidev( mar, aspect, jx, jy )
689  use plplot_flt
690  real(kind=plflt) :: mar, aspect, jx, jy
691  end subroutine plsdidev
692  end interface
693 
694  interface
695  subroutine plsdimap( dimxmi, dimxmax, diymin, dimymax, dimxpmm, diypmm )
696  use plplot_flt
697  real(kind=plflt) :: dimxmi, dimxmax, diymin, dimymax, dimxpmm, diypmm
698  end subroutine plsdimap
699  end interface
700 
701  interface
702  subroutine plsdiori( rot )
703  use plplot_flt
704  real(kind=plflt) :: rot
705  end subroutine plsdiori
706  end interface
707 
708  interface
709  subroutine plsdiplt( xmin, xmax, ymin, ymax )
710  use plplot_flt
711  real(kind=plflt) :: xmin, xmax, ymin, ymax
712  end subroutine plsdiplt
713  end interface
714 
715  interface
716  subroutine plsdiplz( xmin, xmax, ymin, ymax )
717  use plplot_flt
718  real(kind=plflt) :: xmin, xmax, ymin, ymax
719  end subroutine plsdiplz
720  end interface
721 
722  interface
723  subroutine plseed( s )
724  integer :: s
725  end subroutine plseed
726  end interface
727 
728  ! TODO: character-version
729  interface
730  subroutine plsesc( esc )
731  integer :: esc
732  end subroutine plsesc
733  end interface
734 
735  !
736  ! TODO: F95-specific form for these routines
737  !
738  interface plsetmapformc
739  subroutine plsetmapformc( mapform )
740  use plplot_flt
741  interface
742  subroutine mapform( n, x, y )
743  use plplot_flt
744  integer :: n
745  real(kind=plflt), dimension(*) :: x, y
746  end subroutine mapform
747  end interface
748  end subroutine plsetmapformc
749  end interface
750 
751  interface
752  subroutine plsfam( fam, num, bmax )
753  integer :: fam, num, bmax
754  end subroutine plsfam
755  end interface
756 
757  interface
758  subroutine plsfci( fci )
759  use plplot_flt
760  integer(kind=plunicode) :: fci
761  end subroutine plsfci
762  end interface
763 
764  interface
765  subroutine plsfont( family, style, weight )
766  integer :: family, style, weight
767  end subroutine plsfont
768  end interface
769 
770  interface plslabelfunc
771  subroutine plslabelfunc_on( labelfunc )
772  interface
773  subroutine labelfunc(axis, value, label, length)
774  use plplot_flt
775  implicit none
776  integer :: axis, length
777  real(kind=plflt) :: value
778  character*(length) label
779  end subroutine labelfunc
780  end interface
781  end subroutine plslabelfunc_on
782 
783  subroutine plslabelfunc_off( dummy )
784  implicit none
785  integer :: dummy
786  end subroutine plslabelfunc_off
787 
788  subroutine plslabelfunc_none
789  end subroutine plslabelfunc_none
790 
791  end interface
792 
793  interface
794  subroutine plsmaj( def, scale )
795  use plplot_flt
796  real(kind=plflt) :: def, scale
797  end subroutine plsmaj
798  end interface
799 
800  ! plsmem: void * argument tricky - TODO
801  ! plsmema: void * argument tricky - TODO
802 
803  interface
804  subroutine plsmin( def, scale )
805  use plplot_flt
806  real(kind=plflt) :: def, scale
807  end subroutine plsmin
808  end interface
809 
810  interface
811  subroutine plsori( rot )
812  integer :: rot
813  end subroutine plsori
814  end interface
815 
816  interface
817  subroutine plspage( xpmm, ypmm, xwid, ywid, xoff, yoff )
818  use plplot_flt
819  real(kind=plflt) :: xpmm, ypmm
820  integer :: xwid, ywid, xoff, yoff
821  end subroutine plspage
822  end interface
823 
824  interface plspause
825  module procedure plspause
826  end interface
827 
828  interface
829  subroutine plsstrm( strm )
830  integer :: strm
831  end subroutine plsstrm
832  end interface
833 
834  interface
835  subroutine plssub( nx, ny )
836  integer :: nx, ny
837  end subroutine plssub
838  end interface
839 
840  interface
841  subroutine plssym( def, scale )
842  use plplot_flt
843  real(kind=plflt) :: def, scale
844  end subroutine plssym
845  end interface
846 
847  interface
848  subroutine plstar( nx, ny )
849  integer :: nx, ny
850  end subroutine plstar
851  end interface
852 
853  interface plstransform
854  subroutine plstransform1( transformfunc )
855  interface
856  subroutine transformfunc(x, y, xt, yt)
857  use plplot_flt
858  implicit none
859  real(kind=plflt) :: x, y, xt, yt
860  end subroutine transformfunc
861  end interface
862  end subroutine plstransform1
863 
864  subroutine plstransform2( dummy )
865  implicit none
866  integer :: dummy
867  end subroutine plstransform2
868 
869  subroutine plstransform3
870  end subroutine plstransform3
871 
872  end interface
873 
874  interface
875  subroutine plstripa( id, pen, x, y )
876  use plplot_flt
877  integer :: id, pen
878  real(kind=plflt) :: x, y
879  end subroutine plstripa
880  end interface
881 
882  interface
883  subroutine plstripd( id )
884  integer :: id
885  end subroutine plstripd
886  end interface
887 
888  interface
889  subroutine plstyl( n, mark, space )
890  integer :: n, mark, space
891  end subroutine plstyl
892  end interface
893 
894  interface plsurf3d
895  module procedure plsurf3d
896  end interface
897 
898  interface plstripc
899  module procedure plstripc
900  end interface
901 
902  interface plsvect
903  module procedure plsvect1
904  end interface
905 
906  interface plsvect
907  subroutine plsvect2
908  end subroutine plsvect2
909  end interface
910 
911  interface
912  subroutine plsvpa( xmin, xmax, ymin, ymax )
913  use plplot_flt
914  real(kind=plflt) :: xmin, xmax, ymin, ymax
915  end subroutine plsvpa
916  end interface
917 
918  interface
919  subroutine plsxax( digmax, digits )
920  integer :: digmax, digits
921  end subroutine plsxax
922  end interface
923 
924  interface
925  subroutine plsyax( digmax, digits )
926  integer :: digmax, digits
927  end subroutine plsyax
928  end interface
929 
930  interface plsym
931  module procedure plsym
932  end interface
933 
934  interface
935  subroutine plszax( digmax, digits )
936  integer :: digmax, digits
937  end subroutine plszax
938  end interface
939 
940  interface
941  subroutine pltext
942  end subroutine pltext
943  end interface
944 
945  interface
946  subroutine plvasp( aspect )
947  use plplot_flt
948  real(kind=plflt) :: aspect
949  end subroutine plvasp
950  end interface
951 
952  interface
953  subroutine plvpas( xmin, xmax, ymin, ymax, aspect )
954  use plplot_flt
955  real(kind=plflt) :: xmin, xmax, ymin, ymax, aspect
956  end subroutine plvpas
957  end interface
958 
959  interface
960  subroutine plvpor( xmin, xmax, ymin, ymax )
961  use plplot_flt
962  real(kind=plflt) :: xmin, xmax, ymin, ymax
963  end subroutine plvpor
964  end interface
965 
966  interface
967  subroutine plvsta
968  end subroutine plvsta
969  end interface
970 
971  interface
972  subroutine plw3d( basex, basey, height, xmin, xmax, ymin, ymax, zmin, zmax, alt, az )
973  use plplot_flt
974  real(kind=plflt) :: basex, basey, height, xmin, xmax, ymin, ymax, zmin, zmax, alt, az
975  end subroutine plw3d
976  end interface
977 
978  interface
979  subroutine plwidth( width )
980  use plplot_flt
981  real(kind=plflt) :: width
982  end subroutine plwidth
983  end interface
984 
985  interface
986  subroutine plwind( xmin, xmax, ymin, ymax )
987  use plplot_flt
988  real(kind=plflt) :: xmin, xmax, ymin, ymax
989  end subroutine plwind
990  end interface
991 
992  interface plxormod
993  module procedure plxormod
994  end interface
995 
996 
997  private :: convert_to_int
998  private :: convert_to_log
999 
1000 ! -------------------------------------------------------------------
1001 contains
1002 ! -------------------------------------------------------------------
1003  integer function convert_to_int( logvalue )
1004  logical :: logvalue
1005  if ( logvalue ) then
1006  convert_to_int = 1
1007  else
1008  convert_to_int = 0
1009  endif
1010  end function convert_to_int
1011 
1012  logical function convert_to_log( intvalue )
1013  integer :: intvalue
1014  convert_to_log = intvalue.ne.0
1015  end function convert_to_log
1016 
1017  subroutine plbin( x, y, center )
1018  real(kind=plflt), dimension(:) :: x, y
1019  integer :: center
1020 
1021  call plbinf95( size(x), x, y, center )
1022  end subroutine plbin
1023 
1024  subroutine plcolorbar_1( p_colorbar_width, p_colorbar_height, &
1025  opt, position, x, y, &
1026  x_length, y_length, &
1027  bg_color, bb_color, bb_style, &
1028  low_cap_color, high_cap_color, &
1029  cont_color, cont_width, &
1030  n_labels, label_opts, labels, &
1031  n_axes, axis_opts, ticks, sub_ticks, &
1032  n_values, values)
1033  real (kind=plflt) :: p_colorbar_width, p_colorbar_height
1034  integer :: opt, position, bg_color, bb_color, bb_style, cont_color
1035  integer :: n_labels, n_axes
1036  real (kind=plflt) :: x, y, x_length, y_length, low_cap_color, high_cap_color, cont_width
1037  integer, dimension(:) :: label_opts, sub_ticks, n_values
1038  real (kind=plflt), dimension(:) :: ticks
1039  real (kind=plflt), dimension(:,:) :: values
1040  character(len=*), dimension(:) :: labels, axis_opts
1041 
1042  !
1043  ! Convert the text arrays and store the results in a convenient
1044  ! albeit global location. This way we avoid all manner of complications.
1045  ! (Though introducing a potentially nasty one: non-threadsafety)
1046  !
1047  call pllegend07_cnv_text( 3, n_labels, labels )
1048  call pllegend07_cnv_text( 4, n_axes, axis_opts )
1049 
1050  call plcolorbar07(p_colorbar_width, p_colorbar_height, &
1051  opt, position, x, y, &
1052  x_length, y_length, &
1053  bg_color, bb_color, bb_style, &
1054  low_cap_color, high_cap_color, &
1055  cont_color, cont_width, &
1056  n_labels, label_opts, n_axes, ticks, sub_ticks, &
1057  n_values, values)
1058  end subroutine plcolorbar_1
1059 
1060  subroutine plcolorbar_2( p_colorbar_width, p_colorbar_height, &
1061  opt, position, x, y, &
1062  x_length, y_length, &
1063  bg_color, bb_color, bb_style, &
1064  low_cap_color, high_cap_color, &
1065  cont_color, cont_width, &
1066  label_opts, labels, axis_opts, ticks, sub_ticks, &
1067  n_values, values)
1068  real (kind=plflt) :: p_colorbar_width, p_colorbar_height
1069  integer :: opt, position, bg_color, bb_color, bb_style, cont_color
1070  real (kind=plflt) :: x, y, x_length, y_length, low_cap_color, high_cap_color, cont_width
1071  integer, dimension(:) :: label_opts, sub_ticks, n_values
1072  real (kind=plflt), dimension(:) :: ticks
1073  real (kind=plflt), dimension(:,:) :: values
1074  character(len=*), dimension(:) :: labels, axis_opts
1075 
1076  integer :: n_labels, n_axes
1077 
1078  n_labels = size(label_opts,1)
1079  n_axes = size(axis_opts,1)
1080  !
1081  ! Convert the text arrays and store the results in a convenient
1082  ! albeit global location. This way we avoid all manner of complications.
1083  ! (Though introducing a potentially nasty one: non-threadsafety)
1084  !
1085  call pllegend07_cnv_text( 3, n_labels, labels )
1086  call pllegend07_cnv_text( 4, n_axes, axis_opts )
1087 
1088  call plcolorbar07(p_colorbar_width, p_colorbar_height, &
1089  opt, position, x, y, &
1090  x_length, y_length, &
1091  bg_color, bb_color, bb_style, &
1092  low_cap_color, high_cap_color, &
1093  cont_color, cont_width, &
1094  n_labels, label_opts, n_axes, ticks, sub_ticks, &
1095  n_values, values)
1096  end subroutine plcolorbar_2
1097 
1098  subroutine plcpstrm( iplsr, flags )
1099  integer :: iplsr
1100  logical :: flags
1101 
1102  integer :: iflags
1103 
1104  iflags = convert_to_int( flags )
1105  call plcpstrmf95( iplsr, iflags )
1106  end subroutine plcpstrm
1107 
1108  subroutine plerrx( xmin, xmax, y )
1109  real(kind=plflt), dimension(:) :: xmin, xmax, y
1110 
1111  call plerrxf95( size(xmin), xmin, xmax, y )
1112  end subroutine plerrx
1113 
1114  subroutine plerry( x, ymin, ymax )
1115  real(kind=plflt), dimension(:) :: x, ymin, ymax
1116 
1117  call plerryf95( size(x), x, ymin, ymax )
1118  end subroutine plerry
1119 
1120  subroutine plfill( x, y )
1121  real(kind=plflt), dimension(:) :: x, y
1122 
1123  call plfillf95( size(x), x, y )
1124  end subroutine plfill
1125 
1126  subroutine plfill3( x, y, z )
1127  real(kind=plflt), dimension(:) :: x, y, z
1128 
1129  call plfill3f95( size(x), x, y, z )
1130  end subroutine plfill3
1131 
1132  subroutine plgradient( x, y, angle )
1133  real(kind=plflt), dimension(:) :: x, y
1134  real(kind=plflt) :: angle
1135 
1136  call plgradientf95( size(x), x, y, angle )
1137  end subroutine plgradient
1138 
1139  subroutine plgriddata( x, y, z, xg, yg, zg, type, data )
1140  real(kind=plflt), dimension(:) :: x, y, z, xg, yg
1141  real(kind=plflt), dimension(:,:) :: zg
1142  real(kind=plflt) :: data
1143  integer :: type
1144 
1145  call plgriddataf95( x, y, z, size(x), xg, size(xg), yg, size(yg), zg, &
1146  type, data )
1147 
1148  return
1149  end subroutine plgriddata
1150 
1151  subroutine plhist( data, datmin, datmax, nbin, oldwin )
1152  real(kind=plflt), dimension(:) :: data
1153  real(kind=plflt) :: datmin, datmax
1154  integer :: nbin, oldwin
1155 
1156  call plhistf95( size(data), data, datmin, datmax, nbin, oldwin )
1157  end subroutine plhist
1158 
1159 ! subroutine plimagefr( idata, xmin, xmax, ymin, ymax, zmin, zmax, &
1160 ! dxmin, dxmax, dymin, dymax, valuemin, valuemax )
1161 ! real(kind=plflt), dimension(:,:) :: idata
1162 ! real(kind=plflt) :: xmin, xmax, ymin, ymax, zmin, zmax
1163 ! real(kind=plflt) :: dxmin, dxmax, dymin, dymax, &
1164 ! valuemin, valuemax
1165 !
1166 ! integer :: nx, ny
1167 !
1168 ! nx = size(idata,1)
1169 ! ny = size(idata,2)
1170 ! call plimagefrf95( idata, nx, ny, xmin, xmax, ymin, ymax, zmin, zmax, &
1171 ! dxmin, dxmax, dymin, dymax, valuemin, valuemax )
1172 ! end subroutine plimagefr
1173 
1174  subroutine plimage( idata, xmin, xmax, ymin, ymax, zmin, zmax, &
1175  dxmin, dxmax, dymin, dymax )
1176  real(kind=plflt), dimension(:,:) :: idata
1177  real(kind=plflt) :: xmin, xmax, ymin, ymax, zmin, zmax
1178  real(kind=plflt) :: dxmin, dxmax, dymin, dymax
1179 
1180  integer :: nx, ny
1181 
1182  nx = size(idata,1)
1183  ny = size(idata,2)
1184  call plimagef95( idata, nx, ny, xmin, xmax, ymin, ymax, zmin, zmax, &
1185  dxmin, dxmax, dymin, dymax )
1186  end subroutine plimage
1187 
1188  subroutine pllegend_1( legend_width, legend_height, &
1189  opt, position, x, y, &
1190  plot_width, bg_color, bb_color, bb_style, &
1191  nrow, ncolumn, nlegend, opt_array, &
1192  text_offset, text_scale, text_spacing, &
1193  text_justification, text_colors, text, &
1194  box_colors, box_patterns, box_scales, &
1195  box_line_widths, &
1196  line_colors, line_styles, line_widths, &
1197  symbol_colors, symbol_scales, &
1198  symbol_numbers, symbols )
1199 
1200  real(kind=plflt) :: legend_width, legend_height, plot_width, x, y
1201  real(kind=plflt) :: text_offset, text_scale, text_spacing, text_justification
1202  integer :: position, opt, bg_color, bb_color, bb_style
1203  integer :: nrow, ncolumn, nlegend
1204 
1205  character(len=*), dimension(:) :: text, symbols
1206 
1207  integer, dimension(:) :: opt_array, text_colors, box_colors
1208  integer, dimension(:) :: box_patterns
1209  real(kind=plflt), dimension(:) :: box_line_widths
1210  integer, dimension(:) :: line_colors, line_styles
1211  real(kind=plflt), dimension(:) :: line_widths
1212  integer, dimension(:) :: symbol_colors, symbol_numbers
1213  real(kind=plflt), dimension(:) :: box_scales, symbol_scales
1214 
1215  !
1216  ! Convert the text arrays and store the results in a convenient
1217  ! albeit global location. This way we avoid all manner of complications.
1218  ! (Though introducing a potentially nasty one: non-threadsafety)
1219  !
1220  call pllegend07_cnv_text( 1, nlegend, text )
1221  call pllegend07_cnv_text( 2, nlegend, symbols )
1222 
1223  call pllegend07( legend_width, legend_height, opt, position, x, y, &
1224  plot_width, bg_color, bb_color, bb_style, &
1225  nrow, ncolumn, nlegend, opt_array, &
1226  text_offset, text_scale, text_spacing, &
1227  text_justification, text_colors, &
1228  box_colors, box_patterns, box_scales, &
1229  box_line_widths, &
1230  line_colors, line_styles, line_widths, &
1231  symbol_colors, symbol_scales, &
1232  symbol_numbers )
1233 
1234  end subroutine pllegend_1
1235 
1236  subroutine pllegend_2( legend_width, legend_height, &
1237  opt, position, x, y, &
1238  plot_width, bg_color, bb_color, bb_style, &
1239  nrow, ncolumn, opt_array, &
1240  text_offset, text_scale, text_spacing, &
1241  text_justification, text_colors, text, &
1242  box_colors, box_patterns, box_scales, &
1243  box_line_widths, &
1244  line_colors, line_styles, line_widths, &
1245  symbol_colors, symbol_scales, &
1246  symbol_numbers, symbols )
1247 
1248  real(kind=plflt) :: legend_width, legend_height, plot_width, x, y
1249  real(kind=plflt) :: text_offset, text_scale, text_spacing, text_justification
1250  integer :: position, opt, bg_color, bb_color, bb_style
1251  integer :: nrow, ncolumn
1252 
1253  character(len=*), dimension(:) :: text, symbols
1254 
1255  integer, dimension(:) :: opt_array, text_colors, box_colors
1256  integer, dimension(:) :: box_patterns
1257  real(kind=plflt), dimension(:) :: box_line_widths
1258  integer, dimension(:) :: line_colors, line_styles
1259  real(kind=plflt), dimension(:) :: line_widths
1260  integer, dimension(:) :: symbol_colors, symbol_numbers
1261  real(kind=plflt), dimension(:) :: box_scales, symbol_scales
1262 
1263  integer :: nlegend
1264 
1265  !
1266  ! Determine number of legend entries
1267  !
1268  nlegend = min( size(opt_array), size(text) )
1269 
1270  call pllegend_1( legend_width, legend_height, &
1271  opt, position, x, y, &
1272  plot_width, bg_color, bb_color, bb_style, &
1273  nrow, ncolumn, nlegend, opt_array, &
1274  text_offset, text_scale, text_spacing, &
1275  text_justification, text_colors, text, &
1276  box_colors, box_patterns, box_scales, &
1277  box_line_widths, &
1278  line_colors, line_styles, line_widths, &
1279  symbol_colors, symbol_scales, &
1280  symbol_numbers, symbols )
1281 
1282  end subroutine pllegend_2
1283 
1284  subroutine plline( x, y )
1285  real(kind=plflt), dimension(:) :: x, y
1286 
1287  call pllinef95( size(x), x, y )
1288  end subroutine plline
1289 
1290  subroutine plline3( x, y, z )
1291  real(kind=plflt), dimension(:) :: x, y, z
1292 
1293  call plline3f95( size(x), x, y, z )
1294  end subroutine plline3
1295 
1296  subroutine plmap1(mapform,mapname,minx,maxx,miny,maxy)
1297  use plplot_flt
1298  use plplot_str
1299  implicit none
1300  real(kind=plflt) minx, maxx, miny, maxy
1301  character*(*) mapname
1302  external mapform
1303 
1304  call plstrf2c(mapname, string1)
1305 
1306  call plsetmapformc(mapform)
1307  s1 = transfer( string1, s1 )
1308  call plmap7(s1,minx,maxx,miny,maxy)
1309 
1310  end subroutine plmap1
1311 
1312  subroutine plmap2(mapname,minx,maxx,miny,maxy)
1313  use plplot_flt
1314  use plplot_str
1315  implicit none
1316  real(kind=plflt) minx, maxx, miny, maxy
1317  character*(*) mapname
1318 
1319  call plstrf2c(mapname, string1)
1320 
1321  call plclearmapformc()
1322  s1 = transfer( string1, s1 )
1323  call plmap7(s1,minx,maxx,miny,maxy)
1324 
1325  end subroutine plmap2
1326 
1327  subroutine plmeridians1(mapform,dlong,dlat,minlong,maxlong, &
1328  minlat,maxlat)
1329 
1330  implicit none
1331  real(kind=plflt) dlong, dlat, minlong, maxlong, minlat, maxlat
1332  external mapform
1333 
1334  call plsetmapformc(mapform)
1335  call plmeridians7(dlong,dlat,minlong,maxlong,minlat,maxlat)
1336 
1337  end subroutine plmeridians1
1338 
1339  subroutine plmeridians2(dlong,dlat,minlong,maxlong, &
1340  minlat,maxlat)
1341 
1342  implicit none
1343  real(kind=plflt) dlong, dlat, minlong, maxlong, minlat, maxlat
1344 
1345  call plclearmapformc
1346  call plmeridians7(dlong,dlat,minlong,maxlong,minlat,maxlat)
1347 
1348  end subroutine plmeridians2
1349 
1350  subroutine plmesh( x, y, z, opt )
1351  integer :: opt
1352  real(kind=plflt), dimension(:) :: x, y
1353  real(kind=plflt), dimension(:,:) :: z
1354 
1355  call plmeshf95( x, y, z, size(x), size(y), opt, size(x))
1356 
1357  end subroutine plmesh
1358 
1359  subroutine plmeshc( x, y, z, opt, clevel )
1360  integer :: opt
1361  real(kind=plflt), dimension(:) :: x, y, clevel
1362  real(kind=plflt), dimension(:,:) :: z
1363 
1364  call plmeshcf95( x, y, z, size(x), size(y), opt, &
1365  clevel, size(clevel), size(x))
1366 
1367  end subroutine plmeshc
1368 
1369  subroutine plot3d( x, y, z, opt, side )
1370  integer :: opt
1371  logical :: side
1372  real(kind=plflt), dimension(:) :: x, y
1373  real(kind=plflt), dimension(:,:) :: z
1374  integer :: iside
1375 
1376  iside = convert_to_int(side)
1377  call plot3df95( x, y, z, size(x), size(y), opt, iside, size(x))
1378 
1379  end subroutine plot3d
1380 
1381  subroutine plot3dc( x, y, z, opt, clevel )
1382  integer :: opt
1383  real(kind=plflt), dimension(:) :: x, y, clevel
1384  real(kind=plflt), dimension(:,:) :: z
1385 
1386  call plot3dcf95( x, y, z, size(x), size(y), opt, clevel, &
1387  size(clevel), size(x))
1388 
1389  end subroutine plot3dc
1390 
1391  subroutine plspause( lpause )
1392  logical :: lpause
1393 
1394  integer :: ipause
1395 
1396  ipause = convert_to_int( lpause )
1397  call plspausef95( ipause )
1398  end subroutine plspause
1399 
1400  subroutine plsurf3d( x, y, z, opt, clevel )
1401  integer :: opt
1402  real(kind=plflt), dimension(:) :: x, y, clevel
1403  real(kind=plflt), dimension(:,:) :: z
1404 
1405  call plsurf3df95( x, y, z, size(x), size(y), opt, clevel, &
1406  size(clevel), size(x))
1407 
1408  end subroutine plsurf3d
1409 
1410  subroutine plpoin( x, y, code )
1411  integer :: code
1412  real(kind=plflt), dimension(:) :: x, y
1413 
1414  call plpoinf95( size(x), x, y, code )
1415  end subroutine plpoin
1416 
1417  subroutine plpoin3( x, y, z, code )
1418  integer :: code
1419  real(kind=plflt), dimension(:) :: x, y, z
1420 
1421  call plpoin3f95( size(x), x, y, z, code )
1422  end subroutine plpoin3
1423 
1424  subroutine plpoly3( x, y, z, draw, ifcc )
1425  logical :: ifcc
1426  logical, dimension(:) :: draw
1427  real(kind=plflt), dimension(:) :: x, y, z
1428 
1429  integer, dimension(size(draw)) :: idraw
1430  integer :: i
1431  integer :: iifcc
1432 
1433  iifcc = convert_to_int( ifcc )
1434  do i = 1,size(draw)
1435  idraw(i) = convert_to_int( draw(i) )
1436  enddo
1437  call plpoly3f95( size(x), x, y, z, idraw, iifcc )
1438  end subroutine plpoly3
1439 
1440  real (kind=plflt) function plrandd()
1441  external plranddf95
1442  real(kind=plflt) :: plranddf95
1443 
1444  plrandd = plranddf95()
1445  end function plrandd
1446 
1447  subroutine plscmap0( r, g, b )
1448  integer, dimension(:) :: r, g, b
1449 
1450  call plscmap0f95( r, g, b, size(r) )
1451  end subroutine plscmap0
1452 
1453  subroutine plscmap0a( r, g, b, a )
1454  integer, dimension(:) :: r, g, b
1455  real(kind=plflt), dimension(:) :: a
1456 
1457  call plscmap0af95( r, g, b, a, size(r) )
1458  end subroutine plscmap0a
1459 
1460  subroutine plscmap1( r, g, b )
1461  integer, dimension(:) :: r, g, b
1462 
1463  call plscmap1f95( r, g, b, size(r) )
1464  end subroutine plscmap1
1465 
1466  subroutine plscmap1a( r, g, b, a )
1467  integer, dimension(:) :: r, g, b
1468  real(kind=plflt), dimension(:) :: a
1469 
1470  call plscmap1af95( r, g, b, a, size(r) )
1471  end subroutine plscmap1a
1472 
1473  subroutine plscmap1l( rgbtype, intensity, coord1, coord2, coord3, alt_hue_path)
1474  logical :: rgbtype
1475  real(kind=plflt), dimension(:) :: intensity, coord1, coord2, coord3
1476  logical, dimension(:) :: alt_hue_path
1477 
1478  integer, dimension(size(alt_hue_path)) :: ialt_hue_path
1479  integer :: i
1480  integer :: type
1481 
1482  type = convert_to_int( rgbtype )
1483  do i = 1,size(alt_hue_path)
1484  ialt_hue_path(i) = convert_to_int( alt_hue_path(i) )
1485  enddo
1486  call plscmap1lf95( type, size(intensity), intensity, coord1, coord2, coord3, ialt_hue_path )
1487  end subroutine plscmap1l
1488 
1489  subroutine plscmap1l2( rgbtype, intensity, coord1, coord2, coord3)
1490  logical :: rgbtype
1491  real(kind=plflt), dimension(:) :: intensity, coord1, coord2, coord3
1492 
1493  integer :: type
1494 
1495  type = convert_to_int( rgbtype )
1496  call plscmap1l2f95( type, size(intensity), intensity, coord1, coord2, coord3)
1497  end subroutine plscmap1l2
1498 
1499  subroutine plscmap1la( rgbtype, intensity, coord1, coord2, coord3, a, alt_hue_path)
1500  logical :: rgbtype
1501  real(kind=plflt), dimension(:) :: intensity, coord1, coord2, coord3, a
1502  logical, dimension(:) :: alt_hue_path
1503 
1504  integer, dimension(size(alt_hue_path)) :: ialt_hue_path
1505  integer :: i
1506  integer :: type
1507 
1508  type = convert_to_int( rgbtype )
1509  do i = 1,size(alt_hue_path)
1510  ialt_hue_path(i) = convert_to_int( alt_hue_path(i) )
1511  enddo
1512  call plscmap1laf95( type, size(intensity), intensity, coord1, coord2, coord3, a, ialt_hue_path )
1513  end subroutine plscmap1la
1514 
1515  subroutine plscmap1la2( rgbtype, intensity, coord1, coord2, coord3, a)
1516  logical :: rgbtype
1517  real(kind=plflt), dimension(:) :: intensity, coord1, coord2, coord3, a
1518 
1519  integer :: type
1520 
1521  type = convert_to_int( rgbtype )
1522  call plscmap1la2f95( type, size(intensity), intensity, coord1, coord2, coord3, a)
1523  end subroutine plscmap1la2
1524 
1525  subroutine plstripc(id, xspec, yspec, xmin, xmax, xjump, &
1526  ymin, ymax, xlpos, ylpos, y_ascl, acc, &
1527  colbox, collab, colline, styline, legline, &
1528  labx, laby, labtop)
1529 
1530 ! use plplot_str
1531  implicit none
1532  integer id, colbox, collab, colline(4), styline(4)
1533  character*(*) xspec, yspec, legline(4), labx, laby, labtop
1534  real(kind=plflt) xmin, xmax, xjump, ymin, ymax, xlpos, ylpos
1535  logical y_ascl, acc
1536  integer iy_ascl, iacc
1537 
1538 
1539  call plstrf2c(xspec, string1)
1540  call plstrf2c(yspec, string2)
1541  call plstrf2c(legline(1), string3)
1542  call plstrf2c(legline(2), string4)
1543  call plstrf2c(legline(3), string5)
1544  call plstrf2c(legline(4), string6)
1545  call plstrf2c(labx, string7)
1546  call plstrf2c(laby, string8)
1547  call plstrf2c(labtop, string9)
1548 
1549  iy_ascl = convert_to_int( y_ascl )
1550  iacc = convert_to_int( acc )
1551 
1552  s1 = transfer( string1, s1 )
1553  s2 = transfer( string2, s2 )
1554  s3 = transfer( string3, s3 )
1555  s4 = transfer( string4, s4 )
1556  s5 = transfer( string5, s5 )
1557  s6 = transfer( string6, s6 )
1558  s7 = transfer( string7, s7 )
1559  s8 = transfer( string8, s8 )
1560  s9 = transfer( string9, s9 )
1561  call plstripcf95(id, s1, s2, xmin, xmax, xjump, &
1562  ymin, ymax, xlpos, ylpos, iy_ascl, iacc, &
1563  colbox, collab, colline, styline, &
1564  s3, s4, s5, s6, &
1565  s7, s8, s9)
1566 
1567  end subroutine plstripc
1568 
1569  subroutine plsvect1( arrowx, arrowy, fill )
1570  logical :: fill
1571  real(kind=plflt), dimension(:) :: arrowx, arrowy
1572  integer ifill
1573  ifill = convert_to_int(fill)
1574 
1575  call plsvect1f95( arrowx, arrowy, size(arrowx), ifill )
1576  end subroutine plsvect1
1577 
1578  subroutine plsym( x, y, code )
1579  integer :: code
1580  real(kind=plflt), dimension(:) :: x, y
1581 
1582  call plsymf95( size(x), x, y, code )
1583  end subroutine plsym
1584 
1585  subroutine plxormod( mode, status )
1586  logical :: mode, status
1587  integer :: imode, istatus
1588  imode = convert_to_int(mode)
1589  call plxormodf95( imode, istatus)
1590  status = convert_to_log(istatus)
1591  end subroutine plxormod
1592 end module plplot