Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_pic_viewer.F90
Go to the documentation of this file.
1 !**************************************************************
2 ! Copyright INRIA
3 ! Authors :
4 ! CALVI project team
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 
22 
24 
25 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
26 #include "sll_assert.h"
27 #include "sll_memory.h"
28 #include "sll_working_precision.h"
29 
33 
34  use sll_m_utilities, only: &
37 
38  use sll_m_xdmf, only: &
40 
41  use sll_m_cartesian_meshes, only: &
43 
44  use sll_m_pic_visu, only: &
46 
47 #ifndef NOHDF5
48 
49  use sll_m_hdf5_io_serial, only: &
54 
55 #endif
56 
57  use sll_m_xml_io, only: &
61 
62  implicit none
63 
64  public :: &
67 
68  private
69 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
70 
71  type, public :: sll_t_pic_viewer_2d
72 
73  type(sll_t_cartesian_mesh_2d), pointer :: mesh
74  character(len=72) :: label
75 
76  end type
77 
79 
80  module procedure write_2d_field
81  module procedure write_2d_particles
82  module procedure write_2d_field_and_particles
83 
84  end interface sll_o_pic_viewer_write
85 
86 contains
87 
88  function sll_f_new_pic_viewer_2d(mesh, &
89  label &
90  ) result(viewer)
91 
92  type(sll_t_cartesian_mesh_2d), pointer :: mesh
93  type(sll_t_pic_viewer_2d), pointer :: viewer
94  character(len=*) :: label
95 
96  sll_int32 :: ierr
97 
98  sll_allocate(viewer, ierr)
99 
100  viewer%label = label
102  viewer, &
103  mesh)
104 
105  end function sll_f_new_pic_viewer_2d
106 
108  viewer, &
109  mesh)
110 
111  type(sll_t_pic_viewer_2d), pointer :: viewer
112  type(sll_t_cartesian_mesh_2d), pointer :: mesh
113 
114  viewer%mesh => mesh
115 
116  end subroutine initialize_pic_viewer_2d
117 
118  subroutine write_2d_field(viewer, field, iplot)
119 
120  type(sll_t_pic_viewer_2d) :: viewer
121  sll_real64, intent(in) :: field(:, :)
122  sll_int32, intent(in) :: iplot
123 
124  sll_int32 :: nx
125  sll_int32 :: ny
126  sll_int32 :: error
127  sll_int32 :: file_id
128  character(len=4) :: cplot
129 
130 #ifndef NOHDF5
131  type(sll_t_hdf5_ser_handle) :: hfile_id
132 #endif
133 
134  call sll_s_int2string(iplot, cplot)
135  nx = viewer%mesh%num_cells1
136  ny = viewer%mesh%num_cells2
137 
138  !Write the light data
139  call sll_s_xml_file_create(trim(viewer%label)//cplot//".xmf", file_id, error)
140  call sll_o_xml_field(file_id, 'omega', trim(viewer%label)//cplot//".h5:/fp", &
141  nx, ny, "HDF", "Node")
142  call sll_s_xml_file_close(file_id, error)
143 
144  !Write the heavy data
145 #ifndef NOHDF5
146  call sll_s_hdf5_ser_file_create(trim(viewer%label)//cplot//".h5", hfile_id, error)
147  call sll_o_hdf5_ser_write_array(hfile_id, field, "/fp", error)
148  call sll_s_hdf5_ser_file_close(hfile_id, error)
149 #endif
150 
151  end subroutine write_2d_field
152 
153  subroutine write_2d_particles(viewer, xp, yp, op, iplot, time)
154 
155  type(sll_t_pic_viewer_2d) :: viewer
156  sll_real64, intent(in) :: xp(:)
157  sll_real64, intent(in) :: yp(:)
158  sll_real64, intent(in) :: op(:)
159  sll_int32, intent(in) :: iplot
160  sll_real64, optional :: time
161 
162  sll_real64 :: xmin
163  sll_real64 :: ymin
164  sll_real64 :: xmax
165  sll_real64 :: ymax
166  sll_int32 :: error
167  sll_int32 :: file_id
168  character(len=4) :: cplot
169 
170 #ifndef NOHDF5
171  type(sll_t_hdf5_ser_handle) :: hfile_id
172 #endif
173 
174  call sll_s_int2string(iplot, cplot)
175 
176  xmin = viewer%mesh%eta1_min
177  xmax = viewer%mesh%eta1_max
178  ymin = viewer%mesh%eta2_min
179  ymax = viewer%mesh%eta2_max
180 
181  call sll_s_xml_file_create(trim(viewer%label)//cplot//".xmf", file_id, error)
182  call write_particles_2d_xml(file_id, trim(viewer%label)//cplot, xp, yp)
183  call sll_o_xml_field(file_id, 'omega', trim(viewer%label)//cplot//".h5:/op", &
184  size(op), "HDF", "Node")
185  call sll_s_xml_file_close(file_id, error)
186 
187 #ifndef NOHDF5
188  call sll_s_hdf5_ser_file_create(trim(viewer%label)//cplot//".h5", hfile_id, error)
189  call sll_o_hdf5_ser_write_array(hfile_id, xp, "/xp", error)
190  call sll_o_hdf5_ser_write_array(hfile_id, yp, "/yp", error)
191  call sll_o_hdf5_ser_write_array(hfile_id, op, "/op", error)
192  call sll_s_hdf5_ser_file_close(hfile_id, error)
193 #endif
194 
195  end subroutine write_2d_particles
196 
197  subroutine write_2d_field_and_particles(viewer, xp, yp, op, fp, iplot, time)
198 
199  type(sll_t_pic_viewer_2d) :: viewer
200  sll_real64, intent(in) :: xp(:)
201  sll_real64, intent(in) :: yp(:)
202  sll_real64, intent(in) :: op(:)
203  sll_real64, intent(in) :: fp(:, :)
204  sll_int32, intent(in) :: iplot
205  sll_real64, optional :: time
206 
207  sll_int32 :: error
208  sll_int32 :: nx
209  sll_int32 :: ny
210  sll_int32 :: file_id
211  character(len=4) :: cplot
212 
213 #ifndef NOHDF5
214  type(sll_t_hdf5_ser_handle) :: hfile_id
215 #endif
216 
217  nx = viewer%mesh%num_cells1
218  ny = viewer%mesh%num_cells2
219 
220  call sll_s_int2string(iplot, cplot)
221 
222  call sll_s_xml_file_create(trim(viewer%label)//cplot//".xmf", file_id, error)
223  call write_particles_2d_xml(file_id, trim(viewer%label)//cplot, xp, yp)
224  call sll_o_xml_field(file_id, 'omega', trim(viewer%label)//cplot//".h5:/op", &
225  size(op), "HDF", "Node")
226  write (file_id, "(a)") "</Grid>"
227  call write_grid_2d_xml(file_id, viewer)
228  call sll_o_xml_field(file_id, 'omega', trim(viewer%label)//cplot//".h5:/fp", &
229  nx, ny, "HDF", "Node")
230  call sll_s_xml_file_close(file_id, error)
231 
232 #ifndef NOHDF5
233  call sll_s_hdf5_ser_file_create(trim(viewer%label)//cplot//".h5", hfile_id, error)
234  call sll_o_hdf5_ser_write_array(hfile_id, xp, "/xp", error)
235  call sll_o_hdf5_ser_write_array(hfile_id, yp, "/yp", error)
236  call sll_o_hdf5_ser_write_array(hfile_id, op, "/op", error)
237  call sll_o_hdf5_ser_write_array(hfile_id, fp, "/fp", error)
238  call sll_s_hdf5_ser_file_close(hfile_id, error)
239 #endif
240 
241  end subroutine write_2d_field_and_particles
242 
243  subroutine write_grid_2d_xml(file_id, viewer)
244 
245  type(sll_t_pic_viewer_2d) :: viewer
246  sll_int32, intent(in) :: file_id
247  sll_real64 :: xmin
248  sll_real64 :: ymin
249  sll_real64 :: dx
250  sll_real64 :: dy
251  sll_int32 :: nx
252  sll_int32 :: ny
253 
254  xmin = viewer%mesh%eta1_min
255  ymin = viewer%mesh%eta2_min
256  dx = viewer%mesh%delta_eta1
257  dy = viewer%mesh%delta_eta2
258 
259  nx = viewer%mesh%num_cells1
260  ny = viewer%mesh%num_cells2
261 
262  write (file_id, "(a)") "<Grid Name='Grid' GridType='Uniform'>"
263  write (file_id, "(a,2i5,a)") "<Topology TopologyType='2DCoRectMesh' NumberOfElements='", ny, nx, "'/>"
264  write (file_id, "(a)") "<Geometry GeometryType='ORIGIN_DXDY'>"
265  write (file_id, "(a)") "<DataItem Dimensions='2' NumberType='Float' Format='XML'>"
266  write (file_id, "(2f10.5)") xmin, ymin
267  write (file_id, "(a)") "</DataItem>"
268  write (file_id, "(a)") "<DataItem Dimensions='2' NumberType='Float' Format='XML'>"
269  write (file_id, "(2f10.5)") dx, dy
270  write (file_id, "(a)") "</DataItem>"
271  write (file_id, "(a)") "</Geometry>"
272 
273  end subroutine write_grid_2d_xml
274 
275  subroutine write_particles_2d_xml(file_id, prefix, xp, yp)
276 
277  sll_int32, intent(in) :: file_id
278  character(len=*), intent(in) :: prefix
279  sll_real64, intent(in) :: xp(:)
280  sll_real64, intent(in) :: yp(:)
281 
282  sll_int32 :: n
283 
284  n = size(xp)
285  sll_assert(size(yp) == n)
286 
287  write (file_id, "(a)") "<Grid Name=""Particles"" Type=""Uniform"">"
288  write (file_id, "(a,i6,a)") "<Topology TopologyType=""Polyvertex"" NumberOfElements=""", n, """/>"
289  write (file_id, "(a)") "<Geometry Type=""X_Y"">"
290  write (file_id, "(a,i6,a)") "<DataItem Format=""HDF"" Dimensions=""", n, """>"
291  write (file_id, "(a)") prefix//".h5:/xp"
292  write (file_id, "(a)") "</DataItem>"
293  write (file_id, "(a,i6,a)") "<DataItem Format=""HDF"" Dimensions=""", n, """>"
294  write (file_id, "(a)") prefix//".h5:/yp"
295  write (file_id, "(a)") "</DataItem>"
296  write (file_id, "(a)") "</Geometry>"
297 
298  end subroutine write_particles_2d_xml
299 
300 end module sll_m_pic_viewer
301 
Write nD array of double precision floats or integers into HDF5 file.
plot particles centers with gnuplot
write a data attribute in the xml file
Cartesian mesh basic types.
Abstract class for coordinate transformations.
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.
This module provides some routines for plotting fields and particles during PIC simulations.
subroutine write_2d_field_and_particles(viewer, xp, yp, op, fp, iplot, time)
subroutine initialize_pic_viewer_2d(viewer, mesh)
subroutine write_grid_2d_xml(file_id, viewer)
subroutine write_particles_2d_xml(file_id, prefix, xp, yp)
type(sll_t_pic_viewer_2d) function, pointer, public sll_f_new_pic_viewer_2d(mesh, label)
subroutine write_2d_field(viewer, field, iplot)
subroutine write_2d_particles(viewer, xp, yp, op, iplot, time)
This module provides some routines for plotting during PIC simulations.
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.
Implements the functions to write xdmf file plotable by VisIt.
Definition: sll_m_xdmf.F90:24
subroutine, public sll_s_xdmf_corect2d_nodes(file_name, array, array_name, eta1_min, delta_eta1, eta2_min, delta_eta2, file_format, iplot, time)
Subroutine to write a 2D array in xdmf format The field is describe on a cartesian mesh Axis are perp...
Definition: sll_m_xdmf.F90:265
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