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 ! $Id: sfstubsf95.f90 12836 2013-12-09 19:33:53Z arjenmarkus $
3 ! sfstubsf95.f
4 !
5 ! Copyright (C) 2005, 2006 Arjen Markus
6 ! Copyright (C) 2006 Alan W. Irwin
7 !
8 ! This file is part of PLplot.
9 !
10 ! PLplot is free software; you can redistribute it and/or modify
11 ! it under the terms of the GNU Library General Public License as published
12 ! by the Free Software Foundation; either version 2 of the License, or
13 ! (at your option) any later version.
14 !
15 ! PLplot is distributed in the hope that it will be useful,
16 ! but WITHOUT ANY WARRANTY; without even the implied warranty of
17 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 ! GNU Library General Public License for more details.
19 !
20 ! You should have received a copy of the GNU Library General Public License
21 ! along with PLplot; if not, write to the Free Software
22 ! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
23 !
24 !
25 ! This file contains the interfaces for Fortran 95:
26 ! - it includes the actual FORTRAN routines from the FORTRAN 95 bindings
27 ! - it includes interfaces to the C routines from these bindings
28 ! - it defines a few Fortran 95 specific items and interfaces
29 !
30 ! NB
31 ! This module is written in fixed form.
32 ! To enable a redefinition of certain interfaces, we actually have
33 ! two modules.
34 !
35 ! NB
36 ! The INTENT attribute is not currently used. This is a matter to
37 ! be looked into.
38 !
39 ! NB
40 ! It is assumed in the current implementation that all arrays are
41 ! passed with correct dimensions. It would be wiser, perhaps, to
42 ! use the minimum dimension instead of just the dimension of the
43 ! first array.
44 !
45 ! NOTE:
46 ! Some of the C routines will have to be renamed (via macros)
47 !
48 !***********************************************************************
49 
50 !
51 ! Parameters for identifying the kind of PLplot's real
52 ! numbers (a configuration parameter)
53 ! Use whatever name suits you better.
54 !
55 module plplot_flt
56  include 'plflt.inc'
57 end module
58 
59 !
60 ! Parameters and variables for strings / arrays for
61 ! string conversion
62 !
63 module plplot_str
64  integer :: maxleni
65  parameter(maxlen = 320)
66  parameter(maxleni = 80)
67  character (len = maxlen) :: string1, string2, string3
68  character (len = maxlen) :: string4, string5, string6
69  character (len = maxlen) :: string7, string8, string9
70  integer, dimension(maxleni) :: s1, s2, s3, s4, s5, s6, s7, s8, s9
71 
72  character(len=1), parameter :: PL_END_OF_STRING = achar(0)
73 end module
74 
75 module plplotp
76  use plplot_flt
77  use plplot_str
78  use plplot_strutils
79  implicit none
80 
81  interface plcont
82  module procedure plcontour_0
83  module procedure plcontour_1
84  module procedure plcontour_2
85  module procedure plcontour_tr
86  module procedure plcontour_0_all
87  module procedure plcontour_1_all
88  module procedure plcontour_2_all
89  module procedure plcontour_tr_all
90  end interface
93 
94  interface plvect
95  module procedure plvectors_0
96  module procedure plvectors_1
97  module procedure plvectors_2
98  module procedure plvectors_tr
99  end interface
101 
102  interface plshade
103  module procedure plshade_single_0
104  module procedure plshade_single_1
105  module procedure plshade_single_2
106  module procedure plshade_single_tr
107  end interface
109 
110  interface plshades
111  module procedure plshades_multiple_0
112  module procedure plshades_multiple_1
113  module procedure plshades_multiple_1r
114  module procedure plshades_multiple_2
115  module procedure plshades_multiple_2r
116  module procedure plshades_multiple_tr
117  module procedure plshades_multiple_trr
118  end interface
122 
123  interface plimagefr
124  module procedure plimagefr_0
125  module procedure plimagefr_1
126  module procedure plimagefr_2
127  module procedure plimagefr_tr
128  end interface
130 
131 contains
132  include 'sfstubs.f90'
133 end module plplotp
134 
136  use plplot_flt
137  type :: plgraphicsin
138  integer type ! of event (CURRENTLY UNUSED)
139  integer state ! key or button mask
140  integer keysym ! key selected
141  integer button ! mouse button selected
142  integer subwindow ! subwindow (alias subpage, alias subplot) number
143  character(len=16) string ! translated string
144  integer pX, pY ! absolute device coordinates of pointer
145  real(kind=plflt) dX, dY ! relative device coordinates of pointer
146  real(kind=plflt) wX, wY ! world coordinates of pointer
147  end type plgraphicsin
148 end module plplot_types
149 
150 module plplot
151  use plplotp
152  use plplot_flt
153  use plplot_types
154  use plplot_strutils
155  !
156  ! To be added: renaming list
157  !
158 
159  implicit none
160  include 'plplot_parameters.h'
161 
162  !
163  ! To be added: alternative interfaces
164  !
165  interface
166  subroutine pladv( sub )
167  integer :: sub
168  end subroutine pladv
169  end interface
170 
171  interface plbin
172  module procedure plbin
173  end interface
174 
175  interface
176  subroutine plbop
177  end subroutine plbop
178  end interface
179 
180  interface
181  subroutine plcalc_world( rx, ry, wx, wy, window )
182  use plplot_flt
183  real(kind=plflt) :: rx, ry, wx, wy
184  integer :: window
185  end subroutine plcalc_world
186  end interface
187 
188  interface
189  subroutine plclear
190  end subroutine plclear
191  end interface
192 
193  interface
194  subroutine plcol0( icol )
195  integer :: icol
196  end subroutine plcol0
197  end interface
198 
199  interface
200  subroutine plcol1( col )
201  use plplot_flt
202  real(kind=plflt) :: col
203  end subroutine plcol1
204  end interface
205 
206  interface plcolorbar
207  module procedure plcolorbar_1
208  module procedure plcolorbar_2
209  end interface
210 
211  interface plcpstrm
212  module procedure plcpstrm
213  end interface
214 
215  interface
216  subroutine plend
217  end subroutine plend
218  end interface
219 
220  interface
221  subroutine plend1
222  end subroutine plend1
223  end interface
224 
225  interface
226  subroutine plenv( xmin, xmax, ymin, ymax, just, axis )
227  use plplot_flt
228  real(kind=plflt) :: xmin, xmax, ymin, ymax
229  integer :: just, axis
230  end subroutine plenv
231  end interface
232 
233  interface
234  subroutine pleop
235  end subroutine pleop
236  end interface
237 
238  interface plerrx
239  module procedure plerrx
240  end interface
241 
242  interface plerry
243  module procedure plerry
244  end interface
245 
246  interface plfamadv
247  subroutine plfamadv
248  end subroutine plfamadv
249  end interface
250 
251  interface plfill
252  module procedure plfill
253  end interface
254 
255  interface plfill3
256  module procedure plfill3
257  end interface
258 
259  interface
260  subroutine plflush
261  end subroutine plflush
262  end interface
263 
264  interface
265  subroutine plfont( font )
266  integer :: font
267  end subroutine plfont
268  end interface
269 
270  interface
271  subroutine plfontld( charset )
272  integer :: charset
273  end subroutine plfontld
274  end interface
275 
276  interface
277  subroutine plgchr( chrdef, chrht )
278  use plplot_flt
279  real(kind=plflt) :: chrdef, chrht
280  end subroutine plgchr
281  end interface
282 
283  interface
284  subroutine plgcmap1_range( min_color, max_color )
285  use plplot_flt
286  real(kind=plflt) :: min_color, max_color
287  end subroutine plgcmap1_range
288  end interface
289 
290  interface
291  subroutine plgcol0( icol, r, g, b )
292  integer :: icol, r, g, b
293  end subroutine plgcol0
294  end interface
295 
296  interface
297  subroutine plgcol0a( icol, r, g, b, a )
298  use plplot_flt
299  integer :: icol, r, g, b
300  real(kind=plflt) :: a
301  end subroutine plgcol0a
302  end interface
303 
304  interface
305  subroutine plgcolbg( r, g, b )
306  integer :: r, g, b
307  end subroutine plgcolbg
308  end interface
309 
310  interface
311  subroutine plgcolbga( r, g, b, a )
312  use plplot_flt
313  integer :: r, g, b
314  real(kind=plflt) :: a
315  end subroutine plgcolbga
316  end interface
317 
318  interface
319  subroutine plgcompression( compression )
320  integer :: compression
321  end subroutine plgcompression
322  end interface
323 
324  interface
325  subroutine plgdidev( mar, aspect, jx, jy )
326  use plplot_flt
327  real(kind=plflt) :: mar, aspect, jx, jy
328  end subroutine plgdidev
329  end interface
330 
331  interface
332  subroutine plgdiori( rot )
333  use plplot_flt
334  real(kind=plflt) :: rot
335  end subroutine plgdiori
336  end interface
337 
338  interface
339  subroutine plgdiplt( xmin, xmax, ymin, ymax )
340  use plplot_flt
341  real(kind=plflt) :: xmin, xmax, ymin, ymax
342  end subroutine plgdiplt
343  end interface
344 
345  interface
346  subroutine plgetcursor( gin )
347  use plplot_flt
348  use plplot_types
349  type(plgraphicsin) :: gin
350  end subroutine plgetcursor
351  end interface
352 
353  interface
354  subroutine plgfam( fam, num, bmax )
355  integer :: fam, num, bmax
356  end subroutine plgfam
357  end interface
358 
359  interface
360  subroutine plgfci( fci )
361  use plplot_flt
362  integer(kind=plunicode) :: fci
363  end subroutine plgfci
364  end interface
365 
366  interface
367  subroutine plgfont( family, style, weight )
368  integer :: family, style, weight
369  end subroutine plgfont
370  end interface
371 
372  interface
373  subroutine plglevel( level )
374  integer :: level
375  end subroutine plglevel
376  end interface
377 
378  interface
379  subroutine plgpage( xpmm, ypmm, xwid, ywid, xoff, yoff )
380  use plplot_flt
381  real(kind=plflt) :: xpmm, ypmm
382  integer :: xwid, ywid, xoff, yoff
383  end subroutine plgpage
384  end interface
385 
386  interface
387  subroutine plgra
388  end subroutine plgra
389  end interface
390 
391  interface plgradient
392  module procedure plgradient
393  end interface
394 
395  interface plgriddata
396  module procedure plgriddata
397  end interface
398 
399  interface
400  subroutine plgspa( xmin, xmax, ymin, ymax )
401  use plplot_flt
402  real(kind=plflt) :: xmin, xmax, ymin, ymax
403  end subroutine plgspa
404  end interface
405 
406  interface
407  subroutine plgstrm( strm )
408  integer :: strm
409  end subroutine plgstrm
410  end interface
411 
412  interface
413  subroutine plgvpd( xmin, xmax, ymin, ymax )
414  use plplot_flt
415  real(kind=plflt) :: xmin, xmax, ymin, ymax
416  end subroutine plgvpd
417  end interface
418 
419  interface
420  subroutine plgvpw( xmin, xmax, ymin, ymax )
421  use plplot_flt
422  real(kind=plflt) :: xmin, xmax, ymin, ymax
423  end subroutine plgvpw
424  end interface
425 
426  interface
427  subroutine plgxax( digmax, digits )
428  integer :: digmax, digits
429  end subroutine plgxax
430  end interface
431 
432  interface
433  subroutine plgyax( digmax, digits )
434  integer :: digmax, digits
435  end subroutine plgyax
436  end interface
437 
438  interface
439  subroutine plgzax( digmax, digits )
440  integer :: digmax, digits
441  end subroutine plgzax
442  end interface
443 
444  interface plhist
445  module procedure plhist
446  end interface
447 
448  interface
449  subroutine plhls( h, l, s )
450  use plplot_flt
451  real(kind=plflt) :: h, l, s
452  end subroutine plhls
453  end interface
454 
455  interface
456  subroutine plhlsrgb( h, l, s, r, g, b )
457  use plplot_flt
458  real(kind=plflt) :: h, l, s, r, g, b
459  end subroutine plhlsrgb
460  end interface
461 
462  interface
463  subroutine plinit
464  end subroutine plinit
465  end interface
466 
467  interface
468  subroutine pljoin( x1, y1, x2, y2 )
469  use plplot_flt
470  real(kind=plflt) :: x1, y1, x2, y2
471  end subroutine pljoin
472  end interface
473 
474  interface
475  subroutine pllightsource( x, y, z )
476  use plplot_flt
477  real(kind=plflt) :: x, y, z
478  end subroutine pllightsource
479  end interface
480 
481  interface pllegend
482  module procedure pllegend_1
483  module procedure pllegend_2
484  end interface
485 
486  interface plline
487  module procedure plline
488  end interface
489 
490  interface plline3
491  module procedure plline3
492  end interface
493 
494  interface pllsty
495  subroutine pllsty( lin )
496  integer :: lin
497  end subroutine pllsty
498  end interface
499 
500  interface plmap
501  module procedure plmap1, plmap2
502  end interface plmap
503 
504  interface plmeridians
505  module procedure plmeridians1, plmeridians2
506  end interface plmeridians
507 
508  interface plmesh
509  module procedure plmesh
510  end interface
511 
512  interface plmeshc
513  module procedure plmeshc
514  end interface
515 
516  interface
517  subroutine plmkstrm( strm )
518  integer :: strm
519  end subroutine plmkstrm
520  end interface
521 
522  interface
523  subroutine plpat( nlin, inc, del )
524  integer :: nlin, inc, del
525  end subroutine plpat
526  end interface
527 
528  interface
529  subroutine plpath( n, x1, y1, x2, y2 )
530  use plplot_flt
531  integer :: n
532  real(kind=plflt) :: x1, y1, x2, y2
533  end subroutine plpath
534  end interface
535 
536  interface plot3d
537  module procedure plot3d
538  end interface
539 
540  interface plot3dc
541  module procedure plot3dc
542  end interface
543 
544  interface plpoin
545  module procedure plpoin
546  end interface
547 
548  interface plpoin3
549  module procedure plpoin3
550  end interface
551 
552  interface plpoly3
553  module procedure plpoly3
554  end interface
555 
556  interface
557  subroutine plprec( setp, prec )
558  integer :: setp, prec
559  end subroutine plprec
560  end interface
561 
562  interface
563  subroutine plpsty( patt )
564  integer :: patt
565  end subroutine plpsty
566  end interface
567 
568  interface
569  subroutine plreplot
570  end subroutine plreplot
571  end interface
572 
573  !
574  ! Note: plrgb and plrgb1 can be merged
575  !
576  interface
577  subroutine plrgb( r, g, b )
578  use plplot_flt
579  real(kind=plflt) :: r, g, b
580  end subroutine plrgb
581  end interface
582 
583  interface
584  subroutine plrgb1( r, g, b )
585  integer :: r, g, b
586  end subroutine plrgb1
587  end interface
588 
589  interface
590  subroutine plrgbhls( r, g, b, h, l, s )
591  use plplot_flt
592  real(kind=plflt) :: r, g, b, h, l, s
593  end subroutine plrgbhls
594  end interface
595 
596  interface
597  subroutine plschr( chrdef, chrht )
598  use plplot_flt
599  real(kind=plflt) :: chrdef, chrht
600  end subroutine plschr
601  end interface
602 
603  interface plscmap0
604  module procedure plscmap0
605  end interface
606 
607  interface plscmap0a
608  module procedure plscmap0a
609  end interface
610 
611  interface
612  subroutine plscmap0n( n )
613  integer :: n
614  end subroutine plscmap0n
615  end interface
616 
617  interface plscmap1
618  module procedure plscmap1
619  end interface
620 
621  interface plscmap1a
622  module procedure plscmap1a
623  end interface
624 
625  interface plscmap1l
626  module procedure plscmap1l
627  module procedure plscmap1l2
628  end interface
629 
630  interface plscmap1la
631  module procedure plscmap1la
632  module procedure plscmap1la2
633  end interface
634 
635  interface
636  subroutine plscmap1n( n )
637  integer :: n
638  end subroutine plscmap1n
639  end interface
640 
641  interface
642  subroutine plscmap1_range( min_color, max_color )
643  use plplot_flt
644  real(kind=plflt) :: min_color, max_color
645  end subroutine plscmap1_range
646  end interface
647 
648  interface
649  subroutine plscol0( icol, r, g, b )
650  integer :: icol, r, g, b
651  end subroutine plscol0
652  end interface
653 
654  interface
655  subroutine plscol0a( icol, r, g, b, a )
656  use plplot_flt
657  integer :: icol, r, g, b
658  real(kind=plflt) :: a
659  end subroutine plscol0a
660  end interface
661 
662  interface
663  subroutine plscolbg( r, g, b )
664  integer :: r, g, b
665  end subroutine plscolbg
666  end interface
667 
668  interface
669  subroutine plscolbga( r, g, b, a )
670  use plplot_flt
671  integer :: r, g, b
672  real(kind=plflt) :: a
673  end subroutine plscolbga
674  end interface
675 
676  interface
677  subroutine plscolor( color )
678  integer :: color
679  end subroutine plscolor
680  end interface
681 
682  interface
683  subroutine plscompression( compression )
684  integer :: compression
685  end subroutine plscompression
686  end interface
687 
688  interface
689  subroutine plsdidev( mar, aspect, jx, jy )
690  use plplot_flt
691  real(kind=plflt) :: mar, aspect, jx, jy
692  end subroutine plsdidev
693  end interface
694 
695  interface
696  subroutine plsdimap( dimxmi, dimxmax, diymin, dimymax, dimxpmm, diypmm )
697  use plplot_flt
698  real(kind=plflt) :: dimxmi, dimxmax, diymin, dimymax, dimxpmm, diypmm
699  end subroutine plsdimap
700  end interface
701 
702  interface
703  subroutine plsdiori( rot )
704  use plplot_flt
705  real(kind=plflt) :: rot
706  end subroutine plsdiori
707  end interface
708 
709  interface
710  subroutine plsdiplt( xmin, xmax, ymin, ymax )
711  use plplot_flt
712  real(kind=plflt) :: xmin, xmax, ymin, ymax
713  end subroutine plsdiplt
714  end interface
715 
716  interface
717  subroutine plsdiplz( xmin, xmax, ymin, ymax )
718  use plplot_flt
719  real(kind=plflt) :: xmin, xmax, ymin, ymax
720  end subroutine plsdiplz
721  end interface
722 
723  interface
724  subroutine plseed( s )
725  integer :: s
726  end subroutine plseed
727  end interface
728 
729  ! TODO: character-version
730  interface
731  subroutine plsesc( esc )
732  integer :: esc
733  end subroutine plsesc
734  end interface
735 
736  !
737  ! TODO: F95-specific form for these routines
738  !
739  interface plsetmapformc
740  subroutine plsetmapformc( mapform )
741  use plplot_flt
742  interface
743  subroutine mapform( n, x, y )
744  use plplot_flt
745  integer :: n
746  real(kind=plflt), dimension(*) :: x, y
747  end subroutine mapform
748  end interface
749  end subroutine plsetmapformc
750  end interface
751 
752  interface
753  subroutine plsfam( fam, num, bmax )
754  integer :: fam, num, bmax
755  end subroutine plsfam
756  end interface
757 
758  interface
759  subroutine plsfci( fci )
760  use plplot_flt
761  integer(kind=plunicode) :: fci
762  end subroutine plsfci
763  end interface
764 
765  interface
766  subroutine plsfont( family, style, weight )
767  integer :: family, style, weight
768  end subroutine plsfont
769  end interface
770 
771  interface plslabelfunc
772  subroutine plslabelfunc_on( labelfunc )
773  interface
774  subroutine labelfunc(axis, value, label, length)
775  use plplot_flt
776  implicit none
777  integer :: axis, length
778  real(kind=plflt) :: value
779  character*(length) label
780  end subroutine labelfunc
781  end interface
782  end subroutine plslabelfunc_on
783 
784  subroutine plslabelfunc_off( dummy )
785  implicit none
786  integer :: dummy
787  end subroutine plslabelfunc_off
788 
789  subroutine plslabelfunc_none
790  end subroutine plslabelfunc_none
791 
792  end interface
793 
794  interface
795  subroutine plsmaj( def, scale )
796  use plplot_flt
797  real(kind=plflt) :: def, scale
798  end subroutine plsmaj
799  end interface
800 
801  ! plsmem: void * argument tricky - TODO
802  ! plsmema: void * argument tricky - TODO
803 
804  interface
805  subroutine plsmin( def, scale )
806  use plplot_flt
807  real(kind=plflt) :: def, scale
808  end subroutine plsmin
809  end interface
810 
811  interface
812  subroutine plsori( rot )
813  integer :: rot
814  end subroutine plsori
815  end interface
816 
817  interface
818  subroutine plspage( xpmm, ypmm, xwid, ywid, xoff, yoff )
819  use plplot_flt
820  real(kind=plflt) :: xpmm, ypmm
821  integer :: xwid, ywid, xoff, yoff
822  end subroutine plspage
823  end interface
824 
825  interface plspause
826  module procedure plspause
827  end interface
828 
829  interface
830  subroutine plsstrm( strm )
831  integer :: strm
832  end subroutine plsstrm
833  end interface
834 
835  interface
836  subroutine plssub( nx, ny )
837  integer :: nx, ny
838  end subroutine plssub
839  end interface
840 
841  interface
842  subroutine plssym( def, scale )
843  use plplot_flt
844  real(kind=plflt) :: def, scale
845  end subroutine plssym
846  end interface
847 
848  interface
849  subroutine plstar( nx, ny )
850  integer :: nx, ny
851  end subroutine plstar
852  end interface
853 
854  interface plstransform
855  subroutine plstransform1( transformfunc )
856  interface
857  subroutine transformfunc(x, y, xt, yt)
858  use plplot_flt
859  implicit none
860  real(kind=plflt) :: x, y, xt, yt
861  end subroutine transformfunc
862  end interface
863  end subroutine plstransform1
864 
865  subroutine plstransform2( dummy )
866  implicit none
867  integer :: dummy
868  end subroutine plstransform2
869 
870  subroutine plstransform3
871  end subroutine plstransform3
872 
873  end interface
874 
875  interface
876  subroutine plstripa( id, pen, x, y )
877  use plplot_flt
878  integer :: id, pen
879  real(kind=plflt) :: x, y
880  end subroutine plstripa
881  end interface
882 
883  interface
884  subroutine plstripd( id )
885  integer :: id
886  end subroutine plstripd
887  end interface
888 
889  interface
890  subroutine plstyl( n, mark, space )
891  integer :: n, mark, space
892  end subroutine plstyl
893  end interface
894 
895  interface plsurf3d
896  module procedure plsurf3d
897  end interface
898 
899  interface plstripc
900  module procedure plstripc
901  end interface
902 
903  interface plsvect
904  module procedure plsvect1
905 
906  subroutine plsvect2
907  end subroutine plsvect2
908 
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