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 
905  subroutine plsvect2
906  end subroutine plsvect2
907 
908  end interface
909 
910  interface
911  subroutine plsvpa( xmin, xmax, ymin, ymax )
912  use plplot_flt
913  real(kind=plflt) :: xmin, xmax, ymin, ymax
914  end subroutine plsvpa
915  end interface
916 
917  interface
918  subroutine plsxax( digmax, digits )
919  integer :: digmax, digits
920  end subroutine plsxax
921  end interface
922 
923  interface
924  subroutine plsyax( digmax, digits )
925  integer :: digmax, digits
926  end subroutine plsyax
927  end interface
928 
929  interface plsym
930  module procedure plsym
931  end interface
932 
933  interface
934  subroutine plszax( digmax, digits )
935  integer :: digmax, digits
936  end subroutine plszax
937  end interface
938 
939  interface
940  subroutine pltext
941  end subroutine pltext
942  end interface
943 
944  interface
945  subroutine plvasp( aspect )
946  use plplot_flt
947  real(kind=plflt) :: aspect
948  end subroutine plvasp
949  end interface
950 
951  interface
952  subroutine plvpas( xmin, xmax, ymin, ymax, aspect )
953  use plplot_flt
954  real(kind=plflt) :: xmin, xmax, ymin, ymax, aspect
955  end subroutine plvpas
956  end interface
957 
958  interface
959  subroutine plvpor( xmin, xmax, ymin, ymax )
960  use plplot_flt
961  real(kind=plflt) :: xmin, xmax, ymin, ymax
962  end subroutine plvpor
963  end interface
964 
965  interface
966  subroutine plvsta
967  end subroutine plvsta
968  end interface
969 
970  interface
971  subroutine plw3d( basex, basey, height, xmin, xmax, ymin, ymax, zmin, zmax, alt, az )
972  use plplot_flt
973  real(kind=plflt) :: basex, basey, height, xmin, xmax, ymin, ymax, zmin, zmax, alt, az
974  end subroutine plw3d
975  end interface
976 
977  interface
978  subroutine plwidth( width )
979  use plplot_flt
980  real(kind=plflt) :: width
981  end subroutine plwidth
982  end interface
983 
984  interface
985  subroutine plwind( xmin, xmax, ymin, ymax )
986  use plplot_flt
987  real(kind=plflt) :: xmin, xmax, ymin, ymax
988  end subroutine plwind
989  end interface
990 
991  interface plxormod
992  module procedure plxormod
993  end interface
994 
995 
996  private :: convert_to_int
997  private :: convert_to_log
998 
999 ! -------------------------------------------------------------------
1000 contains
1001 ! -------------------------------------------------------------------
1002  integer function convert_to_int( logvalue )
1003  logical :: logvalue
1004  if ( logvalue ) then
1005  convert_to_int = 1
1006  else
1007  convert_to_int = 0
1008  endif
1009  end function convert_to_int
1010 
1011  logical function convert_to_log( intvalue )
1012  integer :: intvalue
1013  convert_to_log = intvalue.ne.0
1014  end function convert_to_log
1015 
1016  subroutine plbin( x, y, center )
1017  real(kind=plflt), dimension(:) :: x, y
1018  integer :: center
1019 
1020  call plbinf95( size(x), x, y, center )
1021  end subroutine plbin
1022 
1023  subroutine plcolorbar_1( p_colorbar_width, p_colorbar_height, &
1024  opt, position, x, y, &
1025  x_length, y_length, &
1026  bg_color, bb_color, bb_style, &
1027  low_cap_color, high_cap_color, &
1028  cont_color, cont_width, &
1029  n_labels, label_opts, labels, &
1030  n_axes, axis_opts, ticks, sub_ticks, &
1031  n_values, values)
1032  real (kind=plflt) :: p_colorbar_width, p_colorbar_height
1033  integer :: opt, position, bg_color, bb_color, bb_style, cont_color
1034  integer :: n_labels, n_axes
1035  real (kind=plflt) :: x, y, x_length, y_length, low_cap_color, high_cap_color, cont_width
1036  integer, dimension(:) :: label_opts, sub_ticks, n_values
1037  real (kind=plflt), dimension(:) :: ticks
1038  real (kind=plflt), dimension(:,:) :: values
1039  character(len=*), dimension(:) :: labels, axis_opts
1040 
1041  !
1042  ! Convert the text arrays and store the results in a convenient
1043  ! albeit global location. This way we avoid all manner of complications.
1044  ! (Though introducing a potentially nasty one: non-threadsafety)
1045  !
1046  call pllegend07_cnv_text( 3, n_labels, labels )
1047  call pllegend07_cnv_text( 4, n_axes, axis_opts )
1048 
1049  call plcolorbar07(p_colorbar_width, p_colorbar_height, &
1050  opt, position, x, y, &
1051  x_length, y_length, &
1052  bg_color, bb_color, bb_style, &
1053  low_cap_color, high_cap_color, &
1054  cont_color, cont_width, &
1055  n_labels, label_opts, n_axes, ticks, sub_ticks, &
1056  n_values, values)
1057  end subroutine plcolorbar_1
1058 
1059  subroutine plcolorbar_2( p_colorbar_width, p_colorbar_height, &
1060  opt, position, x, y, &
1061  x_length, y_length, &
1062  bg_color, bb_color, bb_style, &
1063  low_cap_color, high_cap_color, &
1064  cont_color, cont_width, &
1065  label_opts, labels, axis_opts, ticks, sub_ticks, &
1066  n_values, values)
1067  real (kind=plflt) :: p_colorbar_width, p_colorbar_height
1068  integer :: opt, position, bg_color, bb_color, bb_style, cont_color
1069  real (kind=plflt) :: x, y, x_length, y_length, low_cap_color, high_cap_color, cont_width
1070  integer, dimension(:) :: label_opts, sub_ticks, n_values
1071  real (kind=plflt), dimension(:) :: ticks
1072  real (kind=plflt), dimension(:,:) :: values
1073  character(len=*), dimension(:) :: labels, axis_opts
1074 
1075  integer :: n_labels, n_axes
1076 
1077  n_labels = size(label_opts,1)
1078  n_axes = size(axis_opts,1)
1079  !
1080  ! Convert the text arrays and store the results in a convenient
1081  ! albeit global location. This way we avoid all manner of complications.
1082  ! (Though introducing a potentially nasty one: non-threadsafety)
1083  !
1084  call pllegend07_cnv_text( 3, n_labels, labels )
1085  call pllegend07_cnv_text( 4, n_axes, axis_opts )
1086 
1087  call plcolorbar07(p_colorbar_width, p_colorbar_height, &
1088  opt, position, x, y, &
1089  x_length, y_length, &
1090  bg_color, bb_color, bb_style, &
1091  low_cap_color, high_cap_color, &
1092  cont_color, cont_width, &
1093  n_labels, label_opts, n_axes, ticks, sub_ticks, &
1094  n_values, values)
1095  end subroutine plcolorbar_2
1096 
1097  subroutine plcpstrm( iplsr, flags )
1098  integer :: iplsr
1099  logical :: flags
1100 
1101  integer :: iflags
1102 
1103  iflags = convert_to_int( flags )
1104  call plcpstrmf95( iplsr, iflags )
1105  end subroutine plcpstrm
1106 
1107  subroutine plerrx( xmin, xmax, y )
1108  real(kind=plflt), dimension(:) :: xmin, xmax, y
1109 
1110  call plerrxf95( size(xmin), xmin, xmax, y )
1111  end subroutine plerrx
1112 
1113  subroutine plerry( x, ymin, ymax )
1114  real(kind=plflt), dimension(:) :: x, ymin, ymax
1115 
1116  call plerryf95( size(x), x, ymin, ymax )
1117  end subroutine plerry
1118 
1119  subroutine plfill( x, y )
1120  real(kind=plflt), dimension(:) :: x, y
1121 
1122  call plfillf95( size(x), x, y )
1123  end subroutine plfill
1124 
1125  subroutine plfill3( x, y, z )
1126  real(kind=plflt), dimension(:) :: x, y, z
1127 
1128  call plfill3f95( size(x), x, y, z )
1129  end subroutine plfill3
1130 
1131  subroutine plgradient( x, y, angle )
1132  real(kind=plflt), dimension(:) :: x, y
1133  real(kind=plflt) :: angle
1134 
1135  call plgradientf95( size(x), x, y, angle )
1136  end subroutine plgradient
1137 
1138  subroutine plgriddata( x, y, z, xg, yg, zg, type, data )
1139  real(kind=plflt), dimension(:) :: x, y, z, xg, yg
1140  real(kind=plflt), dimension(:,:) :: zg
1141  real(kind=plflt) :: data
1142  integer :: type
1143 
1144  call plgriddataf95( x, y, z, size(x), xg, size(xg), yg, size(yg), zg, &
1145  type, data )
1146 
1147  return
1148  end subroutine plgriddata
1149 
1150  subroutine plhist( data, datmin, datmax, nbin, oldwin )
1151  real(kind=plflt), dimension(:) :: data
1152  real(kind=plflt) :: datmin, datmax
1153  integer :: nbin, oldwin
1154 
1155  call plhistf95( size(data), data, datmin, datmax, nbin, oldwin )
1156  end subroutine plhist
1157 
1158 ! subroutine plimagefr( idata, xmin, xmax, ymin, ymax, zmin, zmax, &
1159 ! dxmin, dxmax, dymin, dymax, valuemin, valuemax )
1160 ! real(kind=plflt), dimension(:,:) :: idata
1161 ! real(kind=plflt) :: xmin, xmax, ymin, ymax, zmin, zmax
1162 ! real(kind=plflt) :: dxmin, dxmax, dymin, dymax, &
1163 ! valuemin, valuemax
1164 !
1165 ! integer :: nx, ny
1166 !
1167 ! nx = size(idata,1)
1168 ! ny = size(idata,2)
1169 ! call plimagefrf95( idata, nx, ny, xmin, xmax, ymin, ymax, zmin, zmax, &
1170 ! dxmin, dxmax, dymin, dymax, valuemin, valuemax )
1171 ! end subroutine plimagefr
1172 
1173  subroutine plimage( idata, xmin, xmax, ymin, ymax, zmin, zmax, &
1174  dxmin, dxmax, dymin, dymax )
1175  real(kind=plflt), dimension(:,:) :: idata
1176  real(kind=plflt) :: xmin, xmax, ymin, ymax, zmin, zmax
1177  real(kind=plflt) :: dxmin, dxmax, dymin, dymax
1178 
1179  integer :: nx, ny
1180 
1181  nx = size(idata,1)
1182  ny = size(idata,2)
1183  call plimagef95( idata, nx, ny, xmin, xmax, ymin, ymax, zmin, zmax, &
1184  dxmin, dxmax, dymin, dymax )
1185  end subroutine plimage
1186 
1187  subroutine pllegend_1( legend_width, legend_height, &
1188  opt, position, x, y, &
1189  plot_width, bg_color, bb_color, bb_style, &
1190  nrow, ncolumn, nlegend, opt_array, &
1191  text_offset, text_scale, text_spacing, &
1192  text_justification, text_colors, text, &
1193  box_colors, box_patterns, box_scales, &
1194  box_line_widths, &
1195  line_colors, line_styles, line_widths, &
1196  symbol_colors, symbol_scales, &
1197  symbol_numbers, symbols )
1198 
1199  real(kind=plflt) :: legend_width, legend_height, plot_width, x, y
1200  real(kind=plflt) :: text_offset, text_scale, text_spacing, text_justification
1201  integer :: position, opt, bg_color, bb_color, bb_style
1202  integer :: nrow, ncolumn, nlegend
1203 
1204  character(len=*), dimension(:) :: text, symbols
1205 
1206  integer, dimension(:) :: opt_array, text_colors, box_colors
1207  integer, dimension(:) :: box_patterns
1208  real(kind=plflt), dimension(:) :: box_line_widths
1209  integer, dimension(:) :: line_colors, line_styles
1210  real(kind=plflt), dimension(:) :: line_widths
1211  integer, dimension(:) :: symbol_colors, symbol_numbers
1212  real(kind=plflt), dimension(:) :: box_scales, symbol_scales
1213 
1214  !
1215  ! Convert the text arrays and store the results in a convenient
1216  ! albeit global location. This way we avoid all manner of complications.
1217  ! (Though introducing a potentially nasty one: non-threadsafety)
1218  !
1219  call pllegend07_cnv_text( 1, nlegend, text )
1220  call pllegend07_cnv_text( 2, nlegend, symbols )
1221 
1222  call pllegend07( legend_width, legend_height, opt, position, x, y, &
1223  plot_width, bg_color, bb_color, bb_style, &
1224  nrow, ncolumn, nlegend, opt_array, &
1225  text_offset, text_scale, text_spacing, &
1226  text_justification, text_colors, &
1227  box_colors, box_patterns, box_scales, &
1228  box_line_widths, &
1229  line_colors, line_styles, line_widths, &
1230  symbol_colors, symbol_scales, &
1231  symbol_numbers )
1232 
1233  end subroutine pllegend_1
1234 
1235  subroutine pllegend_2( legend_width, legend_height, &
1236  opt, position, x, y, &
1237  plot_width, bg_color, bb_color, bb_style, &
1238  nrow, ncolumn, opt_array, &
1239  text_offset, text_scale, text_spacing, &
1240  text_justification, text_colors, text, &
1241  box_colors, box_patterns, box_scales, &
1242  box_line_widths, &
1243  line_colors, line_styles, line_widths, &
1244  symbol_colors, symbol_scales, &
1245  symbol_numbers, symbols )
1246 
1247  real(kind=plflt) :: legend_width, legend_height, plot_width, x, y
1248  real(kind=plflt) :: text_offset, text_scale, text_spacing, text_justification
1249  integer :: position, opt, bg_color, bb_color, bb_style
1250  integer :: nrow, ncolumn
1251 
1252  character(len=*), dimension(:) :: text, symbols
1253 
1254  integer, dimension(:) :: opt_array, text_colors, box_colors
1255  integer, dimension(:) :: box_patterns
1256  real(kind=plflt), dimension(:) :: box_line_widths
1257  integer, dimension(:) :: line_colors, line_styles
1258  real(kind=plflt), dimension(:) :: line_widths
1259  integer, dimension(:) :: symbol_colors, symbol_numbers
1260  real(kind=plflt), dimension(:) :: box_scales, symbol_scales
1261 
1262  integer :: nlegend
1263 
1264  !
1265  ! Determine number of legend entries
1266  !
1267  nlegend = min( size(opt_array), size(text) )
1268 
1269  call pllegend_1( legend_width, legend_height, &
1270  opt, position, x, y, &
1271  plot_width, bg_color, bb_color, bb_style, &
1272  nrow, ncolumn, nlegend, opt_array, &
1273  text_offset, text_scale, text_spacing, &
1274  text_justification, text_colors, text, &
1275  box_colors, box_patterns, box_scales, &
1276  box_line_widths, &
1277  line_colors, line_styles, line_widths, &
1278  symbol_colors, symbol_scales, &
1279  symbol_numbers, symbols )
1280 
1281  end subroutine pllegend_2
1282 
1283  subroutine plline( x, y )
1284  real(kind=plflt), dimension(:) :: x, y
1285 
1286  call pllinef95( size(x), x, y )
1287  end subroutine plline
1288 
1289  subroutine plline3( x, y, z )
1290  real(kind=plflt), dimension(:) :: x, y, z
1291 
1292  call plline3f95( size(x), x, y, z )
1293  end subroutine plline3
1294 
1295  subroutine plmap1(mapform,mapname,minx,maxx,miny,maxy)
1296  use plplot_flt
1297  use plplot_str
1298  implicit none
1299  real(kind=plflt) minx, maxx, miny, maxy
1300  character*(*) mapname
1301  external mapform
1302 
1303  call plstrf2c(mapname, string1)
1304 
1305  call plsetmapformc(mapform)
1306  s1 = transfer( string1, s1 )
1307  call plmap7(s1,minx,maxx,miny,maxy)
1308 
1309  end subroutine plmap1
1310 
1311  subroutine plmap2(mapname,minx,maxx,miny,maxy)
1312  use plplot_flt
1313  use plplot_str
1314  implicit none
1315  real(kind=plflt) minx, maxx, miny, maxy
1316  character*(*) mapname
1317 
1318  call plstrf2c(mapname, string1)
1319 
1320  call plclearmapformc()
1321  s1 = transfer( string1, s1 )
1322  call plmap7(s1,minx,maxx,miny,maxy)
1323 
1324  end subroutine plmap2
1325 
1326  subroutine plmeridians1(mapform,dlong,dlat,minlong,maxlong, &
1327  minlat,maxlat)
1328 
1329  implicit none
1330  real(kind=plflt) dlong, dlat, minlong, maxlong, minlat, maxlat
1331  external mapform
1332 
1333  call plsetmapformc(mapform)
1334  call plmeridians7(dlong,dlat,minlong,maxlong,minlat,maxlat)
1335 
1336  end subroutine plmeridians1
1337 
1338  subroutine plmeridians2(dlong,dlat,minlong,maxlong, &
1339  minlat,maxlat)
1340 
1341  implicit none
1342  real(kind=plflt) dlong, dlat, minlong, maxlong, minlat, maxlat
1343 
1344  call plclearmapformc
1345  call plmeridians7(dlong,dlat,minlong,maxlong,minlat,maxlat)
1346 
1347  end subroutine plmeridians2
1348 
1349  subroutine plmesh( x, y, z, opt )
1350  integer :: opt
1351  real(kind=plflt), dimension(:) :: x, y
1352  real(kind=plflt), dimension(:,:) :: z
1353 
1354  call plmeshf95( x, y, z, size(x), size(y), opt, size(x))
1355 
1356  end subroutine plmesh
1357 
1358  subroutine plmeshc( x, y, z, opt, clevel )
1359  integer :: opt
1360  real(kind=plflt), dimension(:) :: x, y, clevel
1361  real(kind=plflt), dimension(:,:) :: z
1362 
1363  call plmeshcf95( x, y, z, size(x), size(y), opt, &
1364  clevel, size(clevel), size(x))
1365 
1366  end subroutine plmeshc
1367 
1368  subroutine plot3d( x, y, z, opt, side )
1369  integer :: opt
1370  logical :: side
1371  real(kind=plflt), dimension(:) :: x, y
1372  real(kind=plflt), dimension(:,:) :: z
1373  integer :: iside
1374 
1375  iside = convert_to_int(side)
1376  call plot3df95( x, y, z, size(x), size(y), opt, iside, size(x))
1377 
1378  end subroutine plot3d
1379 
1380  subroutine plot3dc( x, y, z, opt, clevel )
1381  integer :: opt
1382  real(kind=plflt), dimension(:) :: x, y, clevel
1383  real(kind=plflt), dimension(:,:) :: z
1384 
1385  call plot3dcf95( x, y, z, size(x), size(y), opt, clevel, &
1386  size(clevel), size(x))
1387 
1388  end subroutine plot3dc
1389 
1390  subroutine plspause( lpause )
1391  logical :: lpause
1392 
1393  integer :: ipause
1394 
1395  ipause = convert_to_int( lpause )
1396  call plspausef95( ipause )
1397  end subroutine plspause
1398 
1399  subroutine plsurf3d( x, y, z, opt, clevel )
1400  integer :: opt
1401  real(kind=plflt), dimension(:) :: x, y, clevel
1402  real(kind=plflt), dimension(:,:) :: z
1403 
1404  call plsurf3df95( x, y, z, size(x), size(y), opt, clevel, &
1405  size(clevel), size(x))
1406 
1407  end subroutine plsurf3d
1408 
1409  subroutine plpoin( x, y, code )
1410  integer :: code
1411  real(kind=plflt), dimension(:) :: x, y
1412 
1413  call plpoinf95( size(x), x, y, code )
1414  end subroutine plpoin
1415 
1416  subroutine plpoin3( x, y, z, code )
1417  integer :: code
1418  real(kind=plflt), dimension(:) :: x, y, z
1419 
1420  call plpoin3f95( size(x), x, y, z, code )
1421  end subroutine plpoin3
1422 
1423  subroutine plpoly3( x, y, z, draw, ifcc )
1424  logical :: ifcc
1425  logical, dimension(:) :: draw
1426  real(kind=plflt), dimension(:) :: x, y, z
1427 
1428  integer, dimension(size(draw)) :: idraw
1429  integer :: i
1430  integer :: iifcc
1431 
1432  iifcc = convert_to_int( ifcc )
1433  do i = 1,size(draw)
1434  idraw(i) = convert_to_int( draw(i) )
1435  enddo
1436  call plpoly3f95( size(x), x, y, z, idraw, iifcc )
1437  end subroutine plpoly3
1438 
1439  real (kind=plflt) function plrandd()
1440  external plranddf95
1441  real(kind=plflt) :: plranddf95
1442 
1443  plrandd = plranddf95()
1444  end function plrandd
1445 
1446  subroutine plscmap0( r, g, b )
1447  integer, dimension(:) :: r, g, b
1448 
1449  call plscmap0f95( r, g, b, size(r) )
1450  end subroutine plscmap0
1451 
1452  subroutine plscmap0a( r, g, b, a )
1453  integer, dimension(:) :: r, g, b
1454  real(kind=plflt), dimension(:) :: a
1455 
1456  call plscmap0af95( r, g, b, a, size(r) )
1457  end subroutine plscmap0a
1458 
1459  subroutine plscmap1( r, g, b )
1460  integer, dimension(:) :: r, g, b
1461 
1462  call plscmap1f95( r, g, b, size(r) )
1463  end subroutine plscmap1
1464 
1465  subroutine plscmap1a( r, g, b, a )
1466  integer, dimension(:) :: r, g, b
1467  real(kind=plflt), dimension(:) :: a
1468 
1469  call plscmap1af95( r, g, b, a, size(r) )
1470  end subroutine plscmap1a
1471 
1472  subroutine plscmap1l( rgbtype, intensity, coord1, coord2, coord3, alt_hue_path)
1473  logical :: rgbtype
1474  real(kind=plflt), dimension(:) :: intensity, coord1, coord2, coord3
1475  logical, dimension(:) :: alt_hue_path
1476 
1477  integer, dimension(size(alt_hue_path)) :: ialt_hue_path
1478  integer :: i
1479  integer :: type
1480 
1481  type = convert_to_int( rgbtype )
1482  do i = 1,size(alt_hue_path)
1483  ialt_hue_path(i) = convert_to_int( alt_hue_path(i) )
1484  enddo
1485  call plscmap1lf95( type, size(intensity), intensity, coord1, coord2, coord3, ialt_hue_path )
1486  end subroutine plscmap1l
1487 
1488  subroutine plscmap1l2( rgbtype, intensity, coord1, coord2, coord3)
1489  logical :: rgbtype
1490  real(kind=plflt), dimension(:) :: intensity, coord1, coord2, coord3
1491 
1492  integer :: type
1493 
1494  type = convert_to_int( rgbtype )
1495  call plscmap1l2f95( type, size(intensity), intensity, coord1, coord2, coord3)
1496  end subroutine plscmap1l2
1497 
1498  subroutine plscmap1la( rgbtype, intensity, coord1, coord2, coord3, a, alt_hue_path)
1499  logical :: rgbtype
1500  real(kind=plflt), dimension(:) :: intensity, coord1, coord2, coord3, a
1501  logical, dimension(:) :: alt_hue_path
1502 
1503  integer, dimension(size(alt_hue_path)) :: ialt_hue_path
1504  integer :: i
1505  integer :: type
1506 
1507  type = convert_to_int( rgbtype )
1508  do i = 1,size(alt_hue_path)
1509  ialt_hue_path(i) = convert_to_int( alt_hue_path(i) )
1510  enddo
1511  call plscmap1laf95( type, size(intensity), intensity, coord1, coord2, coord3, a, ialt_hue_path )
1512  end subroutine plscmap1la
1513 
1514  subroutine plscmap1la2( rgbtype, intensity, coord1, coord2, coord3, a)
1515  logical :: rgbtype
1516  real(kind=plflt), dimension(:) :: intensity, coord1, coord2, coord3, a
1517 
1518  integer :: type
1519 
1520  type = convert_to_int( rgbtype )
1521  call plscmap1la2f95( type, size(intensity), intensity, coord1, coord2, coord3, a)
1522  end subroutine plscmap1la2
1523 
1524  subroutine plstripc(id, xspec, yspec, xmin, xmax, xjump, &
1525  ymin, ymax, xlpos, ylpos, y_ascl, acc, &
1526  colbox, collab, colline, styline, legline, &
1527  labx, laby, labtop)
1528 
1529 ! use plplot_str
1530  implicit none
1531  integer id, colbox, collab, colline(4), styline(4)
1532  character*(*) xspec, yspec, legline(4), labx, laby, labtop
1533  real(kind=plflt) xmin, xmax, xjump, ymin, ymax, xlpos, ylpos
1534  logical y_ascl, acc
1535  integer iy_ascl, iacc
1536 
1537 
1538  call plstrf2c(xspec, string1)
1539  call plstrf2c(yspec, string2)
1540  call plstrf2c(legline(1), string3)
1541  call plstrf2c(legline(2), string4)
1542  call plstrf2c(legline(3), string5)
1543  call plstrf2c(legline(4), string6)
1544  call plstrf2c(labx, string7)
1545  call plstrf2c(laby, string8)
1546  call plstrf2c(labtop, string9)
1547 
1548  iy_ascl = convert_to_int( y_ascl )
1549  iacc = convert_to_int( acc )
1550 
1551  s1 = transfer( string1, s1 )
1552  s2 = transfer( string2, s2 )
1553  s3 = transfer( string3, s3 )
1554  s4 = transfer( string4, s4 )
1555  s5 = transfer( string5, s5 )
1556  s6 = transfer( string6, s6 )
1557  s7 = transfer( string7, s7 )
1558  s8 = transfer( string8, s8 )
1559  s9 = transfer( string9, s9 )
1560  call plstripcf95(id, s1, s2, xmin, xmax, xjump, &
1561  ymin, ymax, xlpos, ylpos, iy_ascl, iacc, &
1562  colbox, collab, colline, styline, &
1563  s3, s4, s5, s6, &
1564  s7, s8, s9)
1565 
1566  end subroutine plstripc
1567 
1568  subroutine plsvect1( arrowx, arrowy, fill )
1569  logical :: fill
1570  real(kind=plflt), dimension(:) :: arrowx, arrowy
1571  integer ifill
1572  ifill = convert_to_int(fill)
1573 
1574  call plsvect1f95( arrowx, arrowy, size(arrowx), ifill )
1575  end subroutine plsvect1
1576 
1577  subroutine plsym( x, y, code )
1578  integer :: code
1579  real(kind=plflt), dimension(:) :: x, y
1580 
1581  call plsymf95( size(x), x, y, code )
1582  end subroutine plsym
1583 
1584  subroutine plxormod( mode, status )
1585  logical :: mode, status
1586  integer :: imode, istatus
1587  imode = convert_to_int(mode)
1588  call plxormodf95( imode, istatus)
1589  status = convert_to_log(istatus)
1590  end subroutine plxormod
1591 end module plplot