Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_xdmf_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 
31 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
32 #include "sll_working_precision.h"
33 
34  use sll_m_collective, only: &
37 
38  use sll_m_xml_io, only: &
43 
44 #ifndef NOHDF5
45  use sll_m_hdf5_io_parallel, only: &
50 
51 #endif
52  implicit none
53 
54  public :: &
58 
59  private
60 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
61 
63  interface sll_o_xdmf_open
64  module procedure sll_xdmf_open_2d_parallel
65  module procedure sll_xdmf_open_3d_parallel
66  end interface
67 
70  module procedure sll_xdmf_array_2d_parallel
71  module procedure sll_xdmf_array_3d_parallel
72  end interface
73 
75  interface sll_o_xdmf_close
76  module procedure sll_xdmf_close_parallel
77  end interface
78 contains
79 
81  subroutine sll_xdmf_open_2d_parallel(rank, &
82  file_name, &
83  mesh_name, &
84  nnodes_x1, &
85  nnodes_x2, &
86  file_id, &
87  error)
88 
89  sll_int32, intent(in) :: rank
90  character(len=*), intent(in) :: file_name
91  character(len=*), intent(in) :: mesh_name
92  sll_int32, intent(in) :: nnodes_x1
93  sll_int32, intent(in) :: nnodes_x2
94  sll_int32, intent(out) :: file_id
95  sll_int32, intent(out) :: error
96 
97  if (rank == 0) then
98  call sll_s_xml_file_create(trim(file_name), file_id, error)
99  call sll_o_xml_grid_geometry(file_id, trim(mesh_name), &
100  nnodes_x1, nnodes_x2)
101  end if
102 
103  end subroutine sll_xdmf_open_2d_parallel
104 
107  rank, &
108  file_name, &
109  mesh_name, &
110  nnodes_x1, &
111  nnodes_x2, &
112  nnodes_x3, &
113  file_id, &
114  error)
115 
116  sll_int32, intent(in) :: rank
117  character(len=*), intent(in) :: file_name
118  character(len=*), intent(in) :: mesh_name
119  sll_int32, intent(in) :: nnodes_x1
120  sll_int32, intent(in) :: nnodes_x2
121  sll_int32, intent(in) :: nnodes_x3
122  sll_int32, intent(out) :: file_id
123  sll_int32, intent(out) :: error
124 
125  if (rank == 0) then
126  call sll_s_xml_file_create(trim(file_name), file_id, error)
127  call sll_o_xml_grid_geometry(file_id, trim(mesh_name), &
128  nnodes_x1, nnodes_x2, nnodes_x3)
129  end if
130 
131  end subroutine sll_xdmf_open_3d_parallel
132 
134  subroutine sll_xdmf_array_2d_parallel(mesh_name, global_dims, offset, &
135  array, array_name, error, &
136  xmffile_id, center)
137 
138  character(len=*), intent(in) :: mesh_name
139  integer(i64), intent(in) :: global_dims(2)
140  integer(i64), intent(in) :: offset(2)
141  sll_real64, intent(in) :: array(:, :)
142  character(len=*), intent(in) :: array_name
143  sll_int32, intent(out) :: error
144  sll_int32, intent(in), optional :: xmffile_id
145  character(len=4), intent(in), optional :: center
146 
147  type(sll_t_hdf5_par_handle) :: handle
148  sll_int32 :: npoints_x1
149  sll_int32 :: npoints_x2
150  sll_int32 :: prank
151  sll_int32 :: comm
152 
153 #ifndef NOHDF5
154  comm = sll_v_world_collective%comm
155  call sll_s_hdf5_par_file_create(trim(mesh_name)//"-"//trim(array_name)//".h5", &
156  comm, handle, error)
157  call sll_o_hdf5_par_write_array(handle, global_dims, offset, &
158  array, "/"//trim(array_name), error)
159  call sll_s_hdf5_par_file_close(handle, error)
160 #endif
161 
163 
164  if (present(xmffile_id) .and. present(center) .and. prank == 0) then
165  npoints_x1 = int(global_dims(1), 4)
166  npoints_x2 = int(global_dims(2), 4)
167 
168 #ifndef NOHDF5
169  call sll_o_xml_field( &
170  xmffile_id, &
171  trim(array_name), &
172  trim(mesh_name)//"-"//trim(array_name)//".h5:/"//trim(array_name), &
173  npoints_x1, &
174  npoints_x2, &
175  'HDF', &
176  center)
177 #else
178  call sll_o_xml_field(xmffile_id, trim(array_name), &
179  trim(mesh_name)//"-"//trim(array_name)//".bin", &
180  npoints_x1, npoints_x2, 'Binary', center)
181 #endif
182  end if
183 
184  end subroutine sll_xdmf_array_2d_parallel
185 
187  subroutine sll_xdmf_array_3d_parallel(mesh_name, global_dims, offset, &
188  array, array_name, error, xmffile_id, center)
189 
190  character(len=*), intent(in) :: mesh_name
191  integer(i64), intent(in) :: global_dims(3)
192  integer(i64), intent(in) :: offset(3)
193  sll_real64, intent(in) :: array(:, :, :)
194  character(len=*), intent(in) :: array_name
195  sll_int32, intent(out) :: error
196  sll_int32, intent(in), optional :: xmffile_id
197  character(len=4), intent(in), optional :: center
198 
199  type(sll_t_hdf5_par_handle) :: handle
200  sll_int32 :: npoints_x1
201  sll_int32 :: npoints_x2
202  sll_int32 :: npoints_x3
203  sll_int32 :: prank
204  sll_int32 :: comm
205 
206  comm = sll_v_world_collective%comm
207 #ifndef NOHDF5
208  call sll_s_hdf5_par_file_create(trim(mesh_name)//"-"//trim(array_name)//".h5", &
209  comm, handle, error)
210  call sll_o_hdf5_par_write_array(handle, global_dims, offset, array, &
211  "/"//trim(array_name), error)
212  call sll_s_hdf5_par_file_close(handle, error)
213 #endif
214 
216  if (present(xmffile_id) .and. present(center) .and. prank == 0) then
217  npoints_x1 = int(global_dims(1), 4)
218  npoints_x2 = int(global_dims(2), 4)
219  npoints_x3 = int(global_dims(3), 4)
220 
221 #ifndef NOHDF5
222  call sll_o_xml_field( &
223  xmffile_id, &
224  trim(array_name), &
225  trim(mesh_name)//"-"//trim(array_name)//".h5:/"//trim(array_name), &
226  npoints_x1, &
227  npoints_x2, &
228  npoints_x3, &
229  'HDF', &
230  center)
231 #else
232  call sll_o_xml_field(xmffile_id, trim(array_name), &
233  trim(mesh_name)//"-"//trim(array_name)//".bin", &
234  npoints_x1, npoints_x2, npoints_x3, 'Binary', center)
235 #endif
236  end if
237 
238  end subroutine sll_xdmf_array_3d_parallel
239 
241  subroutine sll_xdmf_close_parallel(file_id, error)
242  sll_int32, intent(in) :: file_id
243  sll_int32, intent(out) :: error
244 
245  sll_int32 :: prank
247  if (prank == 0) then
248  call sll_s_xml_file_close(file_id, error)
249  end if
250  end subroutine sll_xdmf_close_parallel
251 
252 end module sll_m_xdmf_parallel
253 
Collectively write distributed nD array into HDF5 file.
write a data attribute in the xml file
write grid description in the xml 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.
Parallel version of sll_hdf5_io.
subroutine, public sll_s_hdf5_par_file_close(handle, error)
Close existing HDF5 file.
subroutine, public sll_s_hdf5_par_file_create(filename, comm, handle, error)
Create new HDF5 file.
Implements the functions to write xdmf file plotable by VisIt.
subroutine sll_xdmf_array_2d_parallel(mesh_name, global_dims, offset, array, array_name, error, xmffile_id, center)
Write 2d array in parallel hdf5 file and the matching line in XDMF file.
subroutine sll_xdmf_close_parallel(file_id, error)
Close the XML file and finish to write last lines.
subroutine sll_xdmf_array_3d_parallel(mesh_name, global_dims, offset, array, array_name, error, xmffile_id, center)
Write 3d array in binary or hdf5 file and the matching line in XDMF file.
subroutine sll_xdmf_open_3d_parallel(rank, file_name, mesh_name, nnodes_x1, nnodes_x2, nnodes_x3, file_id, error)
Open a XDMF format file for a 3d plot.
subroutine sll_xdmf_open_2d_parallel(rank, file_name, mesh_name, nnodes_x1, nnodes_x2, file_id, error)
Open a XDMF format file for a 2d plot.
Implements the functions to write xml file to store light data.
subroutine, public sll_s_xml_file_close(file_id, error)
Close the XML file and finish to write last lines. You give the file unit number.
subroutine, public sll_s_xml_file_create(filename, file_id, error)
Create the XML file and begin to write first lines. You get the file unit number.
    Report Typos and Errors