Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_gnuplot.F90
Go to the documentation of this file.
1 !**************************************************************
2 ! Copyright INRIA
3 !
4 ! This code SeLaLib (for Semi-Lagrangian-Library)
5 ! is a parallel library for simulating the plasma turbulence
6 ! in a tokamak.
7 !
8 ! This software is governed by the CeCILL-B license
9 ! under French law and abiding by the rules of distribution
10 ! of free software. You can use, modify and redistribute
11 ! the software under the terms of the CeCILL-B license as
12 ! circulated by CEA, CNRS and INRIA at the following URL
13 ! "http://www.cecill.info".
14 !**************************************************************
15 
25 !
28 
29 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
30 #include "sll_assert.h"
31 #include "sll_errors.h"
32 #include "sll_working_precision.h"
33 
34  use sll_m_utilities, only: &
36 
37  implicit none
38 
39  public :: &
43 
44  private
45 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
46 
48  interface sll_o_gnuplot_1d
49  module procedure sll_gnuplot_write_1d
50  module procedure sll_s_gnuplot_write
51  module procedure sll_gnuplot_write_two_arrays_1d
52  end interface
53 
55  interface sll_o_gnuplot_2d
56  module procedure sll_gnuplot_corect_2d
57  module procedure sll_gnuplot_rect_2d
58  module procedure sll_gnuplot_curv_2d
59  module procedure sll_gnuplot_mesh_2d
60  module procedure write_unstructured_field
61  end interface
62 
63 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
64 contains
65 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
66 
68  subroutine sll_s_gnuplot_write(array, array_name, iplot)
69 
70  sll_real64, dimension(:), intent(in) :: array
71  character(len=*), intent(in) :: array_name
72  sll_int32, intent(in) :: iplot
73 
74  sll_int32 :: file_id ! file unit number
75  sll_int32 :: error ! error code (not used here..)
76  sll_int32 :: npoints
77  sll_int32 :: ipoints
78  character(len=4) :: cplot
79  character(len=8) :: gnu_status
80 
81  npoints = size(array)
82 
83  ! Check that plot index is strictly positive
84  sll_assert_always(iplot > 0)
85 
86  ! Convert plot index to string
87  call sll_s_int2string(iplot, cplot)
88 
89  ! Determine Gnuplot file status
90  if (iplot == 1) then
91  ! A new ASCII file will be created (replaced if already existing)
92  gnu_status = 'replace'
93  else
94  ! A pre-existing ASCII file will be appended
95  gnu_status = 'old'
96  end if
97 
98  ! Open Gnuplot file
99  open (file=trim(array_name)//'.gnu', &
100  status=gnu_status, &
101  form='formatted', &
102  position='append', &
103  newunit=file_id, &
104  iostat=error)
105 
106  ! Write Gnuplot instructions, than close file
107  write (file_id, "(a)") "plot '"//trim(array_name)//cplot//".dat' with linesp"
108  close (file_id)
109 
110  ! Create new data file
111  open (file=trim(array_name)//cplot//'.dat', &
112  status='replace', &
113  form='formatted', &
114  newunit=file_id, &
115  iostat=error)
116 
117  ! Write arrays, then close file
118  do ipoints = 1, npoints
119  write (file_id, *) sngl(array(ipoints))
120  end do
121  close (file_id)
122 
123  end subroutine sll_s_gnuplot_write
124 
126  subroutine sll_gnuplot_write_two_arrays_1d(array_name, array1, array2, iplot)
127 
128  character(len=*), intent(in) :: array_name
129  sll_real64, dimension(:), intent(in) :: array1
130  sll_real64, dimension(:), intent(in) :: array2
131  sll_int32, intent(in) :: iplot
132 
133  sll_int32 :: file_id ! file unit number
134  sll_int32 :: error ! error code (not used here..)
135  sll_int32 :: n
136  sll_int32 :: i
137  character(len=4) :: cplot
138  character(len=8) :: gnu_status
139 
140  n = size(array1)
141  sll_assert_always(size(array2) == n)
142 
143  ! Check that plot index is strictly positive
144  sll_assert_always(iplot > 0)
145 
146  ! Convert plot index to string
147  call sll_s_int2string(iplot, cplot)
148 
149  ! Determine Gnuplot file status
150  if (iplot == 1) then
151  ! A new ASCII file will be created (replaced if already existing)
152  gnu_status = 'replace'
153  else
154  ! A pre-existing ASCII file will be appended
155  gnu_status = 'old'
156  end if
157 
158  ! Open Gnuplot file
159  open (file=trim(array_name)//'.gnu', &
160  status=gnu_status, &
161  form='formatted', &
162  position='append', &
163  newunit=file_id, &
164  iostat=error)
165 
166  ! Write Gnuplot instructions, than close file
167  write (file_id, "(a)") "plot '"//trim(array_name)//cplot//".dat' with linesp, &
168  & '"//trim(array_name)//cplot//".dat' u 1:3 w l"
169  close (file_id)
170 
171  ! Create new data file
172  open (file=trim(array_name)//cplot//'.dat', &
173  status='replace', &
174  form='formatted', &
175  newunit=file_id, &
176  iostat=error)
177 
178  ! Write arrays, then close file
179  do i = 1, n
180  write (file_id, *) i, sngl(array1(i)), sngl(array2(i))
181  end do
182  close (file_id)
183 
184  end subroutine sll_gnuplot_write_two_arrays_1d
185 
191  subroutine sll_gnuplot_write_1d(y_array, x_array, array_name, iplot)
192 
193  sll_real64, dimension(:), intent(in) :: y_array
194  sll_real64, dimension(:), intent(in) :: x_array
195  character(len=*), intent(in) :: array_name
196  sll_int32, optional, intent(in) :: iplot
197 
198  sll_int32 :: error ! error code
199  sll_int32 :: file_id ! file unit number
200  sll_int32 :: npoints
201  sll_int32 :: ipoints
202  character(len=4) :: cplot
203  character(len=8) :: gnu_status
204  logical :: l_exist
205  character(len=*), parameter :: sub_name = "sll_o_gnuplot_1d"
206  character(len=200) :: message
207 
208  npoints = size(x_array)
209 
210  if (present(iplot)) then
211  ! Check that plot index is strictly positive
212  sll_assert_always(iplot > 0)
213 
214  ! Convert plot index to string
215  call sll_s_int2string(iplot, cplot)
216 
217  ! Determine Gnuplot file status
218  if (iplot == 1) then
219  gnu_status = 'replace'
220  else
221  inquire (file=trim(array_name)//'.gnu', exist=l_exist)
222  if (l_exist) then
223  gnu_status = 'old'
224  else
225  message = "The file "//trim(array_name)//'.gnu' &
226  & //" does not exist and iplot>1. A first call with iplot=1 is needed"
227  sll_error(sub_name, message)
228  end if
229  end if
230 
231  ! Open Gnuplot file
232  open (file=trim(array_name)//'.gnu', &
233  status=gnu_status, &
234  form='formatted', &
235  position='append', &
236  newunit=file_id, &
237  iostat=error)
238 
239  ! Write Gnuplot instructions, than close file
240  write (file_id, "(a)") "plot '"//trim(array_name)//cplot//".dat' with linesp"
241  close (file_id)
242 
243  ! Create new data file, numbered
244  open (file=trim(array_name)//cplot//'.dat', &
245  status='replace', &
246  form='formatted', &
247  newunit=file_id, &
248  iostat=error)
249  else
250  ! Create new data file, not numbered
251  open (file=trim(array_name)//'.dat', &
252  status='replace', &
253  form='formatted', &
254  newunit=file_id, &
255  iostat=error)
256  end if
257 
258  ! Write array, then close file
259  do ipoints = 1, npoints
260  write (file_id, *) sngl(x_array(ipoints)), sngl(y_array(ipoints))
261  end do
262  close (file_id)
263 
264  end subroutine sll_gnuplot_write_1d
265 
280  subroutine sll_gnuplot_corect_2d(xmin, xmax, nx, &
281  ymin, ymax, ny, &
282  array, array_name, &
283  iplot, error)
284  sll_real64, intent(in) :: xmin
285  sll_real64, intent(in) :: xmax
286  sll_real64, intent(in) :: ymin
287  sll_real64, intent(in) :: ymax
288  sll_int32, intent(in) :: nx
289  sll_int32, intent(in) :: ny
290  sll_real64, intent(in) :: array(:, :)
291  character(len=*), intent(in) :: array_name
292  sll_int32, intent(in) :: iplot
293  sll_int32, intent(out) :: error
294 
295  sll_int32 :: file_id
296  sll_int32 :: i, j
297  sll_real64 :: dx, dy, x, y
298  character(len=4) :: fin
299  character(len=8) :: gnu_status
300  logical :: file_exists
301 
302  ! Check that plot index is strictly positive
303  sll_assert_always(iplot > 0)
304 
305  ! Convert plot index to string
306  call sll_s_int2string(iplot, fin)
307 
308  ! Determine Gnuplot file status
309  if (iplot == 1) then
310  ! A new ASCII file will be created (replaced if already existing)
311  gnu_status = 'replace'
312  else
313  inquire (file=array_name//'.gnu', exist=file_exists)
314  if (.not. file_exists) then
315  sll_error("sll_gnuplot_corect_2d", "The gnu file does not exist, perharps you did not call this suboutine with iplot=1 for the first plot")
316  end if
317  ! A pre-existing ASCII file will be appended
318  gnu_status = 'old'
319  end if
320 
321  ! Open Gnuplot file
322  open (file=array_name//'.gnu', &
323  status=gnu_status, &
324  form='formatted', &
325  position='append', &
326  newunit=file_id, &
327  iostat=error)
328 
329  ! Write Gnuplot instructions, than close file
330  write (file_id, *) "splot '"//array_name//"_"//fin//".dat' w l"
331  close (file_id)
332 
333  ! Create new data file
334  open (file=array_name//'_'//fin//'.dat', &
335  status='replace', &
336  form='formatted', &
337  newunit=file_id, &
338  iostat=error)
339 
340  ! Write array, then close file
341  dx = (xmax - xmin)/real(nx - 1, f64)
342  dy = (ymax - ymin)/real(ny - 1, f64)
343  x = xmin
344  do i = 1, nx
345  y = ymin
346  do j = 1, ny
347  write (file_id, *) sngl(x), sngl(y), sngl(array(i, j))
348  y = y + dy
349  end do
350  x = x + dx
351  write (file_id, *)
352  end do
353  close (file_id)
354 
355  end subroutine sll_gnuplot_corect_2d
356 
367  subroutine sll_gnuplot_rect_2d(nx, xvec, ny, yvec, &
368  array, array_name, iplot, error)
369 
370  sll_int32, intent(in) :: nx
371  sll_int32, intent(in) :: ny
372  sll_real64, intent(in) :: xvec(nx)
373  sll_real64, intent(in) :: yvec(ny)
374  sll_real64, intent(in) :: array(nx, ny)
375  character(len=*), intent(in) :: array_name
376  sll_int32, intent(in) :: iplot
377  sll_int32, intent(out) :: error
378 
379  sll_int32 :: file_id
380  sll_int32 :: i, j
381  character(len=4) :: fin
382  character(len=8) :: gnu_status
383 
384  ! Check that plot index is strictly positive
385  sll_assert_always(iplot > 0)
386 
387  ! Convert plot index to string
388  call sll_s_int2string(iplot, fin)
389 
390  ! Determine Gnuplot file status
391  if (iplot == 1) then
392  ! A new ASCII file will be created (replaced if already existing)
393  gnu_status = 'replace'
394  else
395  ! A pre-existing ASCII file will be appended
396  gnu_status = 'old'
397  end if
398 
399  ! Open Gnuplot file
400  open (file=array_name//'.gnu', &
401  status=gnu_status, &
402  form='formatted', &
403  position='append', &
404  newunit=file_id, &
405  iostat=error)
406 
407  ! Write Gnuplot instructions, than close file
408  write (file_id, *) "splot '"//array_name//"_"//fin//".dat' w l"
409  close (file_id)
410 
411  ! Create new data file
412  open (file=array_name//'_'//fin//'.dat', &
413  status='replace', &
414  form='formatted', &
415  newunit=file_id, &
416  iostat=error)
417 
418  ! Write array, then close file
419  do i = 1, nx
420  do j = 1, ny
421  write (file_id, *) sngl(xvec(i)), &
422  sngl(yvec(j)), &
423  sngl(array(i, j))
424  end do
425  write (file_id, *)
426  end do
427  close (file_id)
428 
429  end subroutine sll_gnuplot_rect_2d
430 
438  subroutine sll_gnuplot_mesh_2d(nx, ny, xcoord, ycoord, array_name, error)
439 
440  sll_int32, intent(in) :: nx
441  sll_int32, intent(in) :: ny
442  sll_real64, dimension(nx, ny), intent(in) :: xcoord
443  sll_real64, dimension(nx, ny), intent(in) :: ycoord
444  character(len=*), intent(in) :: array_name
445  sll_int32, intent(out) :: error
446 
447  sll_int32 :: file_id
448  sll_int32 :: i, j
449 
450  ! Create new Gnuplot file (replace if existing)
451  open (file=array_name//'.gnu', &
452  status='replace', &
453  form='formatted', &
454  newunit=file_id, &
455  iostat=error)
456 
457  ! Write Gnuplot instructions, than close file
458  write (file_id, *) "set view 0,0"
459  write (file_id, *) "splot '"//array_name//"_mesh.dat' with lines"
460  close (file_id)
461 
462  ! Create new data file (replace if existing)
463  open (file=array_name//'_mesh.dat', &
464  status='replace', &
465  form='formatted', &
466  newunit=file_id, &
467  iostat=error)
468 
469  ! Write array, then close file
470  do i = 1, nx
471  do j = 1, ny
472  write (file_id, *) sngl(xcoord(i, j)), &
473  sngl(ycoord(i, j)), 0.0_f32
474  end do
475  write (file_id, *)
476  end do
477  close (file_id)
478 
479  end subroutine sll_gnuplot_mesh_2d
480 
493  subroutine sll_gnuplot_curv_2d(nx, ny, x, y, array, array_name, iplot, error)
494 
495  sll_int32, intent(in) :: nx
496  sll_int32, intent(in) :: ny
497  sll_real64, intent(in) :: x(nx, ny)
498  sll_real64, intent(in) :: y(nx, ny)
499  sll_real64, intent(in) :: array(nx, ny)
500  character(len=*), intent(in) :: array_name
501  sll_int32, intent(in) :: iplot
502  sll_int32, intent(out) :: error
503 
504  sll_int32 :: file_id
505  sll_int32 :: i, j
506  character(len=4) :: fin
507  character(len=8) :: gnu_status
508 
509  ! Check that plot index is strictly positive
510  sll_assert_always(iplot > 0)
511 
512  ! Convert plot index to string
513  call sll_s_int2string(iplot, fin)
514 
515  ! Determine Gnuplot file status
516  if (iplot == 1) then
517  ! A new ASCII file will be created (replaced if already existing)
518  gnu_status = 'replace'
519  else
520  ! A pre-existing ASCII file will be appended
521  gnu_status = 'old'
522  end if
523 
524  ! Open Gnuplot file
525  open (file=array_name//'.gnu', &
526  status=gnu_status, &
527  form='formatted', &
528  position='append', &
529  newunit=file_id, &
530  iostat=error)
531 
532  ! Write Gnuplot instructions, then close file
533  write (file_id, *) "set output '"//array_name//"_"//fin//".png'"
534  write (file_id, *) "splot '"//array_name//"_"//fin//".dat' w l"
535  close (file_id)
536 
537  ! Create new data file
538  open (file=array_name//'_'//fin//'.dat', &
539  status='replace', &
540  form='formatted', &
541  newunit=file_id, &
542  iostat=error)
543 
544  ! Write array, then close file
545  do i = 1, nx
546  do j = 1, ny
547  write (file_id, *) sngl(x(i, j)), &
548  sngl(y(i, j)), &
549  sngl(array(i, j))
550  end do
551  write (file_id, *)
552  end do
553  close (file_id)
554 
555  end subroutine sll_gnuplot_curv_2d
556 
563  subroutine write_unstructured_field(field_at_node, &
564  field_name, &
565  coord, &
566  nodes, &
567  plot_number)
568 
569  sll_real64, dimension(:), intent(in) :: field_at_node
570  character(len=*), intent(in) :: field_name
571  sll_real64, dimension(:, :), intent(in) :: coord
572  sll_int32, dimension(:, :), intent(in) :: nodes
573  sll_int32, intent(in) :: plot_number
574 
575  character(len=4) :: cplot
576  character(len=8) :: gnu_status
577  sll_int32 :: file_id
578  sll_int32 :: num_cells
579  sll_int32 :: ierr
580  sll_int32 :: i
581  sll_real32 :: xs1, xs2, xs3
582  sll_real32 :: ys1, ys2, ys3
583 
584  ! Get number of nodes, and verify that three coordinates are given
585  num_cells = size(nodes, 2)
586  sll_assert_always(size(nodes, 1) == 3)
587 
588  ! Check that plot index is strictly positive
589  sll_assert_always(plot_number > 0)
590 
591  ! Convert plot index to string
592  call sll_s_int2string(plot_number, cplot)
593 
594  ! Print message
595  write (*, "(/10x, 'Output file GNUplot ',a/)") field_name//'_'//cplot//'.dat'
596 
597  ! Determine Gnuplot file status
598  if (plot_number == 1) then
599  ! A new ASCII file will be created (replaced if already existing)
600  gnu_status = 'replace'
601  else
602  ! A pre-existing ASCII file will be appended
603  gnu_status = 'old'
604  end if
605 
606  ! Open Gnuplot file
607  open (file=field_name//'.gnu', &
608  status=gnu_status, &
609  form='formatted', &
610  position='append', &
611  newunit=file_id, &
612  iostat=ierr)
613 
614  ! Write Gnuplot instructions, then close file
615  write (file_id, *) "set title 'Field "//field_name//"'"
616  write (file_id, *) "splot '"//field_name//'_'//cplot//".dat' w l"
617  close (file_id)
618 
619  ! Create new data file
620  open (file=field_name//'_'//cplot//'.dat', &
621  status='replace', &
622  form='formatted', &
623  newunit=file_id, &
624  iostat=ierr)
625 
626  ! Write field, then close file
627  do i = 1, num_cells
628 
629  xs1 = sngl(coord(1, nodes(1, i))); ys1 = sngl(coord(2, nodes(1, i)))
630  xs2 = sngl(coord(1, nodes(2, i))); ys2 = sngl(coord(2, nodes(2, i)))
631  xs3 = sngl(coord(1, nodes(3, i))); ys3 = sngl(coord(2, nodes(3, i)))
632 
633  write (file_id, "(3e12.3)") xs1, ys1, field_at_node(nodes(1, i))
634  write (file_id, "(3e12.3)") xs2, ys2, field_at_node(nodes(2, i))
635  write (file_id, "(3e12.3)") xs3, ys3, field_at_node(nodes(3, i))
636  write (file_id, "(3e12.3)") xs1, ys1, field_at_node(nodes(1, i)) !CLOSE THE SQUARE
637  write (file_id, *)
638  write (file_id, *)
639 
640  end do
641  close (file_id)
642 
643  end subroutine write_unstructured_field
644 
645 end module sll_m_gnuplot
646 
write file plotable by gnuplot to visualize 2d field
Write file for gnuplot to display 2d field.
Implements the functions to write data file plotable by GNUplot.
subroutine sll_gnuplot_rect_2d(nx, xvec, ny, yvec, array, array_name, iplot, error)
Write a data file plotable by gnuplot to visualize a 2d field on structured rectangular mesh where sp...
subroutine sll_gnuplot_write_two_arrays_1d(array_name, array1, array2, iplot)
Write two arrays to display with gnuplot.
subroutine sll_gnuplot_curv_2d(nx, ny, x, y, array, array_name, iplot, error)
write a data file plotable by gnuplot.
subroutine sll_gnuplot_corect_2d(xmin, xmax, nx, ymin, ymax, ny, array, array_name, iplot, error)
Write a data file plotable by gnuplot to visualize a 2d field.
subroutine, public sll_s_gnuplot_write(array, array_name, iplot)
Write an array to display with gnuplot.
subroutine sll_gnuplot_write_1d(y_array, x_array, array_name, iplot)
This subroutine write a data file to plot a 1d curve.
subroutine write_unstructured_field(field_at_node, field_name, coord, nodes, plot_number)
Write a field on unstructures mesh of triangles.
subroutine sll_gnuplot_mesh_2d(nx, ny, xcoord, ycoord, array_name, error)
Write a data file plotable by gnuplot to visualize a 2d curvilinear mesh.
Some common numerical utilities.
subroutine, public sll_s_int2string(istep, cstep)
Convert an integer < 9999 to a 4 characters string.
    Report Typos and Errors