Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_gnuplot_parallel.F90
Go to the documentation of this file.
1 !**************************************************************
2 ! Copyright INRIA
3 ! Authors :
4 ! Pierre Navaro
5 !
6 ! This code SeLaLib (for Semi-Lagrangian-Library)
7 ! is a parallel library for simulating the plasma turbulence
8 ! in a tokamak.
9 !
10 ! This software is governed by the CeCILL-B license
11 ! under French law and abiding by the rules of distribution
12 ! of free software. You can use, modify and redistribute
13 ! the software under the terms of the CeCILL-B license as
14 ! circulated by CEA, CNRS and INRIA at the following URL
15 ! "http://www.cecill.info".
16 !**************************************************************
17 
27 #define MPI_MASTER 0
28 
30 
31 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
32 #include "sll_assert.h"
33 #include "sll_working_precision.h"
34 
35 #ifdef __INTEL_COMPILER
36  use ifport
37 #endif
38 
39  use sll_m_ascii_io, only: &
41 
42  use sll_m_utilities, only: &
45 
46  use sll_m_collective, only: &
50 
51  implicit none
52 
53  public :: &
57 
58  private
59 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
60 
63  module procedure sll_s_gnuplot_curv_2d_parallel
64  module procedure sll_s_gnuplot_rect_2d_parallel
65  end interface sll_o_gnuplot_2d_parallel
66 
67 contains
68 
70  subroutine sll_s_gnuplot_curv_2d_parallel(array_x, array_y, array, &
71  array_name, iplot, error)
72 
73  sll_real64, dimension(:, :), intent(in) :: array_x
74  sll_real64, dimension(:, :), intent(in) :: array_y
75  sll_real64, dimension(:, :), intent(in) :: array
76  character(len=*), intent(in) :: array_name
77  sll_int32, intent(in) :: iplot
78  sll_int32, intent(out):: error
79 
80  character(len=4) :: fin
81  sll_int32, save :: gnu_id
82  sll_int32 :: file_id
83  sll_int32 :: i, j
84  character(len=4) :: cproc
85  sll_int32 :: iproc, nproc
86  logical :: dir_e
87  logical, save :: first_call = .true.
88 
91 
92  call sll_s_int2string(iproc, cproc)
93  call sll_s_int2string(iplot, fin)
94 
95 #ifdef __INTEL_COMPILER
96  if (makedirqq(cproc)) then
97  print *, ' Make directory '//cproc
98  end if
99 #else
100  inquire (file=cproc//"/"".", exist=dir_e)
101  if (.not. dir_e) then
102  call execute_command_line("mkdir -p "//cproc)
103  end if
104 #endif
105 
106  sll_assert_always(size(array_x, 1) == size(array_y, 1))
107  sll_assert_always(size(array_x, 2) == size(array_y, 2))
108  sll_assert_always(size(array, 1) == size(array_x, 1))
109  sll_assert_always(size(array, 2) == size(array_x, 2))
110  sll_assert_always(size(array, 1) == size(array_y, 1))
111  sll_assert_always(size(array, 2) == size(array_y, 2))
112 
113  call sll_s_new_file_id(file_id, error)
114  call sll_s_ascii_file_create(cproc//"/"//array_name//'_'//fin//'.dat', file_id, error)
115 
116  ! we need to get rid of size() calls for anything that is not an argument
117  ! check...
118  do i = 1, size(array, 1)
119  do j = 1, size(array, 2)
120  write (file_id, *) sngl(array_x(i, j)), sngl(array_y(i, j)), sngl(array(i, j))
121  end do
122  write (file_id, *)
123  end do
124  close (file_id)
125 
126  if (iproc == mpi_master) then
127 
128  if (first_call) then
129  call sll_s_new_file_id(gnu_id, error)
130  open (gnu_id, file=array_name//".gnu")
131  rewind(gnu_id)
132  first_call = .false.
133  else
134  open (gnu_id, file=array_name//".gnu", position="append")
135  end if
136 
137  write (gnu_id, *) "set title 'Time = ", iplot, "'"
138  write (gnu_id, "(a)", advance='no') &
139  "splot '"//"0000/"//array_name//'_'//fin//".dat' w l"
140  do iproc = 1, nproc - 1
141  call sll_s_int2string(iproc, cproc)
142  write (gnu_id, "(a)", advance='no') &
143  ",'"//cproc//"/"//array_name//'_'//fin//".dat' w l "
144  end do
145  write (gnu_id, *)
146  close (gnu_id)
147 
148  end if
149 
150  end subroutine sll_s_gnuplot_curv_2d_parallel
151 
164  subroutine sll_s_gnuplot_rect_2d_parallel(x_min, delta_x, &
165  y_min, delta_y, &
166  npts_x, npts_y, &
167  array, array_name, &
168  iplot, error)
169 
170  sll_real64, intent(in) :: x_min
171  sll_real64, intent(in) :: delta_x
172  sll_real64, intent(in) :: y_min
173  sll_real64, intent(in) :: delta_y
174  sll_int32, intent(in) :: npts_x
175  sll_int32, intent(in) :: npts_y
176  sll_real64, intent(in) :: array(:, :)
177  character(len=*), intent(in) :: array_name
178  sll_int32, intent(in) :: iplot
179  sll_int32, intent(out) :: error
180 
181  character(len=4) :: fin
182  sll_int32, save :: gnu_id
183  sll_int32 :: file_id
184  sll_int32 :: i, j
185  sll_real64 :: x, y
186  character(len=4) :: cproc
187  sll_int32 :: iproc, nproc
188  logical :: dir_e
189  logical, save :: first_call = .true.
190 
193 
194  call sll_s_int2string(iproc, cproc)
195  call sll_s_int2string(iplot, fin)
196 
197 #ifdef __INTEL_COMPILER
198  if (makedirqq(cproc)) then
199  print *, ' Make directory '//cproc
200  end if
201 #else
202  inquire (file=cproc//"/"".", exist=dir_e)
203  if (.not. dir_e) then
204  call execute_command_line("mkdir -p "//cproc)
205  end if
206 #endif
207 
208  call sll_s_new_file_id(file_id, error)
209  call sll_s_ascii_file_create(cproc//"/"//array_name//'_'//fin//'.dat', file_id, error)
210  do j = 1, npts_y
211  y = y_min + real(j - 1, f64)*delta_y
212  do i = 1, npts_x
213  x = x_min + real(i - 1, f64)*delta_x
214  write (file_id, *) sngl(x), sngl(y), sngl(array(i, j))
215  end do
216  write (file_id, *)
217  end do
218  close (file_id)
219 
220  if (iproc == mpi_master) then
221 
222  if (first_call) then
223  call sll_s_new_file_id(gnu_id, error)
224  open (gnu_id, file=array_name//".gnu")
225  rewind(gnu_id)
226  first_call = .false.
227  else
228  open (gnu_id, file=array_name//".gnu", position="append")
229  end if
230 
231  write (gnu_id, *) "set title 'Time = ", iplot, "'"
232  write (gnu_id, *) "set output '"//array_name//'_'//fin//".png'"
233  write (gnu_id, "(a)", advance='no') &
234  "splot '"//"0000/"//array_name//'_'//fin//".dat' w l"
235  do iproc = 1, nproc - 1
236  call sll_s_int2string(iproc, cproc)
237  write (gnu_id, "(a)", advance='no') &
238  ",'"//cproc//"/"//array_name//'_'//fin//".dat' w l "
239  end do
240  write (gnu_id, *)
241  close (gnu_id)
242 
243  end if
244 
245  end subroutine sll_s_gnuplot_rect_2d_parallel
246 
247 ! it would be nice to have a function like the one commented below,
248 ! (unfinished) which would work with layouts. But this means that the level
249 ! of abstraction of this module would be elevated above remap. This is very
250 ! tricky since this module is used as a reasonably low-level module. Yet some
251 ! solution to this would be desirable.
252 
253 !!$!> write a data file plotable by gnuplot to visualize a 2d field
254 !!$subroutine sll_gnuplot_draw_parallel_array_2d( &
255 !!$ sll_t_layout_2d, &
256 !!$ x_min_global, &
257 !!$ delta_x, &
258 !!$ y_min_global, &
259 !!$ delta_y, &
260 !!$ array, &
261 !!$ array_name, &
262 !!$ iplot, &
263 !!$ error)
264 !!$
265 !!$ type(sll_t_layout_2d), pointer :: sll_t_layout_2d
266 !!$ sll_real64, intent(in) :: x_min_global !< Box corners
267 !!$ sll_real64, intent(in) :: delta_x !< step size
268 !!$ sll_real64, intent(in) :: y_min_global !< Box corners
269 !!$ sll_real64, intent(in) :: delta_y !< step size
270 !!$ sll_real64, dimension(:,:), intent(in) :: array !< data
271 !!$ character(len=*), intent(in) :: array_name !< field name
272 !!$ sll_int32, intent(in) :: iplot !< plot counter
273 !!$ sll_int32, intent(in) :: error !< error code
274 !!$
275 !!$ type(sll_collective_t), pointer :: col
276 !!$ sll_int32, dimension(2) :: gi ! for storing global indices of corner
277 !!$ character(len=4) :: fin
278 !!$ sll_int32, save :: gnu_id
279 !!$ sll_int32 :: file_id
280 !!$ sll_int32 :: i, j
281 !!$ sll_real64 :: x, y
282 !!$ character(len=4) :: cproc
283 !!$ sll_int32 :: comm, iproc, nproc
284 !!$ logical :: dir_e
285 !!$
286 !!$ col =>get_layout_collective(sll_t_layout_2d)
287 !!$ nproc = sll_f_get_collective_size(col)
288 !!$ iproc = sll_f_get_collective_rank(col)
289 !!$ call sll_s_int2string(iproc, cproc)
290 !!$ comm = sll_v_world_collective%comm
291 !!$ call sll_s_int2string(iplot, fin)
292 !!$
293 !!$ gi(1:2) = sll_o_local_to_global(sll_t_layout_2d, (/ 1, 1 /) )
294 !!$
295 !!$ !shouldn't we allow this same instruction when iplot is zero??
296 !!$ if (iplot == 1) then
297 !!$ inquire(file=cproc//"/"".", exist=dir_e)
298 !!$ if (dir_e) then
299 !!$ write(*,*) "directory "//cproc//" exists!"
300 !!$ else
301 !!$ call execute_command_line("mkdir -p "//cproc)
302 !!$ end if
303 !!$ end if
304 !!$
305 !!$ call sll_s_new_file_id(file_id, error)
306 !!$ call sll_s_ascii_file_create(cproc//"/"//array_name//'_'//fin//'.dat', &
307 !!$ file_id, error )
308 !!$ do i = 1, size(array,1)
309 !!$ x = x_min+(i-1)*delta_x
310 !!$ do j = 1, size(array,2)
311 !!$ y = y_min + (j-1)*delta_y
312 !!$ write(file_id,*) x, y, sngl(array(i,j))
313 !!$ end do
314 !!$ write(file_id,*)
315 !!$enddo
316 !!$close(file_id)
317 !!$
318 !!$if (iproc == MPI_MASTER) then
319 !!$
320 !!$ if (iplot == 1) call sll_s_new_file_id(gnu_id, error)
321 !!$ open(gnu_id,file=array_name//".gnu", position="append")
322 !!$ if (iplot == 1) rewind(gnu_id)
323 !!$
324 !!$ write(gnu_id,*)"set title 'Time = ",iplot,"'"
325 !!$ write(gnu_id,"(a)",advance='no') &
326 !!$ "splot '"//"0000/"//array_name//'_'//fin//".dat' w l"
327 !!$ do iproc = 1, nproc - 1
328 !!$ call sll_s_int2string(iproc, cproc)
329 !!$ write(gnu_id,"(a)",advance='no') &
330 !!$ ",'"//cproc//"/"//array_name//'_'//fin//".dat' w l "
331 !!$ end do
332 !!$ write(gnu_id,*)
333 !!$ close(gnu_id)
334 !!$
335 !!$end if
336 !!$
337 !!$end subroutine sll_s_gnuplot_rect_2d_parallel
338 
339 end module sll_m_gnuplot_parallel
Create a gnuplot file for a 2d mesh (cartesian or curvilinear)
Module that contains routines to write data in ASCII format file.
subroutine, public sll_s_ascii_file_create(file_name, file_id, error)
Create ASCII file.
Parallelizing facility.
integer(kind=i32) function, public sll_f_get_collective_rank(col)
Determines the rank of the calling process in the communicator.
type(sll_t_collective_t), pointer, public sll_v_world_collective
Control of subset assignment. Processes with the same color are in the same new communicator.
integer(kind=i32) function, public sll_f_get_collective_size(col)
Determines the size of the group associated with a communicator.
parallel version of sll_m_gnuplot
subroutine, public sll_s_gnuplot_curv_2d_parallel(array_x, array_y, array, array_name, iplot, error)
write a data file plotable by gnuplot to visualize a 2d field
subroutine, public sll_s_gnuplot_rect_2d_parallel(x_min, delta_x, y_min, delta_y, npts_x, npts_y, array, array_name, iplot, error)
write a data file plotable by gnuplot to visualize a 2d field
Some common numerical utilities.
subroutine, public sll_s_int2string(istep, cstep)
Convert an integer < 9999 to a 4 characters string.
subroutine, public sll_s_new_file_id(file_id, error)
Get a file unit number free before creating a file.
    Report Typos and Errors