Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_xdmf_serial_blocks.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 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
28 #include "sll_assert.h"
29 #include "sll_working_precision.h"
30 
31  use sll_m_utilities, only: &
33 
34  use sll_m_xml_io, only: &
37 
38  use mpi, only: &
39  mpi_comm_rank, &
40  mpi_comm_size, &
41  mpi_comm_world, &
42  mpi_integer, &
43  mpi_recv, &
44  mpi_send, &
45  mpi_status_size
46 
47 #ifndef NOHDF5
48  use sll_m_hdf5_io_serial, only: &
53 
54 #endif
55  implicit none
56 
57  public :: &
61 
62  private
63 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
64 
66  interface sll_xdmf_open
67  module procedure sll_s_xdmf_open_serial_blocks
68  end interface
69 
72  module procedure sll_s_xdmf_array_2d_serial_blocks
73  end interface
74 
76  interface sll_xdmf_close
77  module procedure sll_s_xdmf_close_serial_blocks
78  end interface
79 
80 contains
81 
83  subroutine sll_s_xdmf_open_serial_blocks(prefix, &
84  x1, &
85  x2, &
86  xmf, &
87  error)
88 
89  character(len=*), intent(in) :: prefix
90  sll_real64, intent(in) :: x1(:, :)
91  sll_real64, intent(in) :: x2(:, :)
92  sll_int32, intent(out) :: xmf
93  sll_int32 :: error
94  sll_int32 :: prank
95  character(len=4) :: crank
96 #ifndef NOHDF5
97  type(sll_t_hdf5_ser_handle) :: h5file_id
98 #endif
99 
100  call mpi_comm_rank(mpi_comm_world, prank, error)
101  call sll_s_int2string(prank, crank)
102  if (prank == 0) then
103  call sll_s_xml_file_create(trim(prefix)//".xmf", xmf, error)
104  end if
105 #ifndef NOHDF5
106  call sll_s_hdf5_ser_file_create(trim(prefix)//"-mesh-"//crank//".h5", h5file_id, error)
107  call sll_o_hdf5_ser_write_array(h5file_id, x1, "/x1", error)
108  call sll_o_hdf5_ser_write_array(h5file_id, x2, "/x2", error)
109  call sll_s_hdf5_ser_file_close(h5file_id, error)
110 #endif
111 
112  end subroutine sll_s_xdmf_open_serial_blocks
113 
115  subroutine sll_s_xdmf_array_2d_serial_blocks(xmf, prefix, &
116  array, array_name, error, &
117  center)
118 
119  sll_int32, intent(in) :: xmf
120  character(len=*), intent(in) :: prefix
121  sll_real64, intent(in) :: array(:, :)
122  character(len=*), intent(in) :: array_name
123 #ifndef NOHDF5
124  type(sll_t_hdf5_ser_handle) :: h5file_id
125 #endif
126  sll_int32 :: npts_x1
127  sll_int32 :: npts_x2
128  sll_int32, intent(out) :: error
129  character(len=4), optional :: center
130  sll_int32 :: prank
131  sll_int32 :: psize
132  character(len=4) :: crank
133  sll_int32 :: iproc
134  character(len=4) :: cproc
135  sll_int32, parameter :: tag = 1111
136  sll_int32 :: statut(mpi_status_size)
137 
138  call mpi_comm_rank(mpi_comm_world, prank, error)
139  call mpi_comm_size(mpi_comm_world, psize, error)
140  call sll_s_int2string(prank, crank)
141 
142 #ifndef NOHDF5
143  call sll_s_hdf5_ser_file_create(trim(prefix)//"-"//crank//".h5", h5file_id, error)
144  call sll_o_hdf5_ser_write_array(h5file_id, array, "/"//trim(array_name), error)
145  call sll_s_hdf5_ser_file_close(h5file_id, error)
146 #endif
147 
148  npts_x1 = size(array, 1)
149  npts_x2 = size(array, 2)
150 
151  if (prank == 0) then
152  write (xmf, '(a)') &
153  "<Grid Name=""AllDomain"" GridType=""Collection"" CollectionType=""Spatial"">"
154  do iproc = 0, psize - 1
155  call sll_s_int2string(iproc, cproc)
156  if (iproc > 0) then
157  call mpi_recv(npts_x1, 1, mpi_integer, iproc, tag, mpi_comm_world, statut, error)
158  call mpi_recv(npts_x2, 1, mpi_integer, iproc, tag, mpi_comm_world, statut, error)
159  end if
160 
161  write (xmf, '(a)') "<Grid Name=""SubDomain"" GridType=""Uniform"">"
162  write (xmf, '(a,2i6,a)') &
163  "<Topology TopologyType=""2DSMesh"" NumberOfElements='", npts_x2, npts_x1, "'/>"
164  write (xmf, '(a)') "<Geometry GeometryType=""X_Y"">"
165  write (xmf, '(a,2i6,a)') "<DataItem Dimensions='", npts_x2, npts_x1, &
166  "' NumberType=""Float"" Precision=""4"" Format=""HDF"">"
167  write (xmf, '(a)') trim(prefix)//"-mesh-"//cproc//".h5:/x1"
168  write (xmf, '(a)') "</DataItem>"
169  write (xmf, '(a,2i6,a)') "<DataItem Dimensions='", npts_x2, npts_x1, &
170  "' NumberType=""Float"" Precision=""4"" Format=""HDF"">"
171  write (xmf, '(a)') trim(prefix)//"-mesh-"//cproc//".h5:/x2"
172  write (xmf, '(a)') "</DataItem>"
173  write (xmf, '(a)') "</Geometry>"
174  write (xmf, '(a)') &
175  "<Attribute Name='"//trim(array_name)//"' AttributeType=""Scalar"" Center=""Node"">"
176  write (xmf, '(a,2i6,a)') "<DataItem Dimensions='", npts_x2, npts_x1, &
177  "' NumberType=""Float"" Precision=""8"" Format=""HDF"">"
178  write (xmf, '(a)') trim(prefix)//"-"//cproc//".h5:/"//trim(array_name)
179  write (xmf, '(a)') "</DataItem>"
180  write (xmf, '(a)') "</Attribute>"
181  write (xmf, '(a)') "</Grid>"
182  end do
183  else
184  call mpi_send(npts_x1, 1, mpi_integer, 0, tag, mpi_comm_world, error)
185  call mpi_send(npts_x2, 1, mpi_integer, 0, tag, mpi_comm_world, error)
186  end if
187 
188  return
189  sll_assert_always(present(center))
190 
191  end subroutine sll_s_xdmf_array_2d_serial_blocks
192 
194  subroutine sll_s_xdmf_close_serial_blocks(file_id, error)
195  sll_int32, intent(in) :: file_id
196  sll_int32, intent(out) :: error
197  sll_int32 :: prank
198  logical :: i_opened, i_exist
199  character(len=255) :: i_name
200 
201  call mpi_comm_rank(mpi_comm_world, prank, error)
202  if (prank == 0) then
203  INQUIRE (file_id, opened=i_opened, name=i_name, exist=i_exist)
204  call sll_s_xml_file_close(file_id, error)
205  end if
206  end subroutine sll_s_xdmf_close_serial_blocks
207 
208 end module sll_m_xdmf_serial_blocks
Write nD array of double precision floats or integers into HDF5 file.
Create a new xmdf file to plot parallel array using hdf5 serial library.
Implements the functions to write hdf5 file to store heavy data.
subroutine, public sll_s_hdf5_ser_file_create(filename, handle, error)
Create new HDF5 file.
subroutine, public sll_s_hdf5_ser_file_close(handle, error)
Close existing HDF5 file.
Some common numerical utilities.
subroutine, public sll_s_int2string(istep, cstep)
Convert an integer < 9999 to a 4 characters string.
Implements the functions to write xdmf file plotable by VisIt.
subroutine, public sll_s_xdmf_array_2d_serial_blocks(xmf, prefix, array, array_name, error, center)
Write 2d array in parallel hdf5 file and the matching line in XDMF file.
subroutine, public sll_s_xdmf_open_serial_blocks(prefix, x1, x2, xmf, error)
Open a XDMF format file for a 2d plot.
subroutine, public sll_s_xdmf_close_serial_blocks(file_id, error)
Close the XML file and finish to write last lines.
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