30 #include "sll_assert.h"
31 #include "sll_errors.h"
32 #include "sll_working_precision.h"
70 sll_real64,
dimension(:),
intent(in) :: array
71 character(len=*),
intent(in) :: array_name
72 sll_int32,
intent(in) :: iplot
78 character(len=4) :: cplot
79 character(len=8) :: gnu_status
84 sll_assert_always(iplot > 0)
92 gnu_status =
'replace'
99 open (file=trim(array_name)//
'.gnu', &
107 write (file_id,
"(a)")
"plot '"//trim(array_name)//cplot//
".dat' with linesp"
111 open (file=trim(array_name)//cplot//
'.dat', &
118 do ipoints = 1, npoints
119 write (file_id, *) sngl(array(ipoints))
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
137 character(len=4) :: cplot
138 character(len=8) :: gnu_status
141 sll_assert_always(
size(array2) == n)
144 sll_assert_always(iplot > 0)
152 gnu_status =
'replace'
159 open (file=trim(array_name)//
'.gnu', &
167 write (file_id,
"(a)")
"plot '"//trim(array_name)//cplot//
".dat' with linesp, &
168 & '"//trim(array_name)//cplot//
".dat' u 1:3 w l"
172 open (file=trim(array_name)//cplot//
'.dat', &
180 write (file_id, *) i, sngl(array1(i)), sngl(array2(i))
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
202 character(len=4) :: cplot
203 character(len=8) :: gnu_status
205 character(len=*),
parameter :: sub_name =
"sll_o_gnuplot_1d"
206 character(len=200) :: message
208 npoints =
size(x_array)
210 if (
present(iplot))
then
212 sll_assert_always(iplot > 0)
219 gnu_status =
'replace'
221 inquire (file=trim(array_name)//
'.gnu', exist=l_exist)
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)
232 open (file=trim(array_name)//
'.gnu', &
240 write (file_id,
"(a)")
"plot '"//trim(array_name)//cplot//
".dat' with linesp"
244 open (file=trim(array_name)//cplot//
'.dat', &
251 open (file=trim(array_name)//
'.dat', &
259 do ipoints = 1, npoints
260 write (file_id, *) sngl(x_array(ipoints)), sngl(y_array(ipoints))
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
297 sll_real64 :: dx, dy, x, y
298 character(len=4) :: fin
299 character(len=8) :: gnu_status
300 logical :: file_exists
303 sll_assert_always(iplot > 0)
311 gnu_status =
'replace'
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")
322 open (file=array_name//
'.gnu', &
330 write (file_id, *)
"splot '"//array_name//
"_"//fin//
".dat' w l"
334 open (file=array_name//
'_'//fin//
'.dat', &
341 dx = (xmax - xmin)/real(nx - 1, f64)
342 dy = (ymax - ymin)/real(ny - 1, f64)
347 write (file_id, *) sngl(x), sngl(y), sngl(array(i, j))
368 array, array_name, iplot, error)
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
381 character(len=4) :: fin
382 character(len=8) :: gnu_status
385 sll_assert_always(iplot > 0)
393 gnu_status =
'replace'
400 open (file=array_name//
'.gnu', &
408 write (file_id, *)
"splot '"//array_name//
"_"//fin//
".dat' w l"
412 open (file=array_name//
'_'//fin//
'.dat', &
421 write (file_id, *) sngl(xvec(i)), &
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
451 open (file=array_name//
'.gnu', &
458 write (file_id, *)
"set view 0,0"
459 write (file_id, *)
"splot '"//array_name//
"_mesh.dat' with lines"
463 open (file=array_name//
'_mesh.dat', &
472 write (file_id, *) sngl(xcoord(i, j)), &
473 sngl(ycoord(i, j)), 0.0_f32
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
506 character(len=4) :: fin
507 character(len=8) :: gnu_status
510 sll_assert_always(iplot > 0)
518 gnu_status =
'replace'
525 open (file=array_name//
'.gnu', &
533 write (file_id, *)
"set output '"//array_name//
"_"//fin//
".png'"
534 write (file_id, *)
"splot '"//array_name//
"_"//fin//
".dat' w l"
538 open (file=array_name//
'_'//fin//
'.dat', &
547 write (file_id, *) sngl(x(i, j)), &
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
575 character(len=4) :: cplot
576 character(len=8) :: gnu_status
578 sll_int32 :: num_cells
581 sll_real32 :: xs1, xs2, xs3
582 sll_real32 :: ys1, ys2, ys3
585 num_cells =
size(nodes, 2)
586 sll_assert_always(
size(nodes, 1) == 3)
589 sll_assert_always(plot_number > 0)
595 write (*,
"(/10x, 'Output file GNUplot ',a/)") field_name//
'_'//cplot//
'.dat'
598 if (plot_number == 1)
then
600 gnu_status =
'replace'
607 open (file=field_name//
'.gnu', &
615 write (file_id, *)
"set title 'Field "//field_name//
"'"
616 write (file_id, *)
"splot '"//field_name//
'_'//cplot//
".dat' w l"
620 open (file=field_name//
'_'//cplot//
'.dat', &
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)))
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))
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.