Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_xdmf_light_serial.F90
Go to the documentation of this file.
1 !------------------------------------------------------------------------------
2 ! SELALIB
3 !------------------------------------------------------------------------------
4 ! MODULE: sll_m_xdmf_light_serial
5 !
6 ! DESCRIPTION:
11 !------------------------------------------------------------------------------
13 
14 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
15 #include "sll_working_precision.h"
16 
17  use sll_m_io_utilities, only: &
19 
20  use sll_m_xml, only: &
23 
24  implicit none
25 
26  public :: &
28 
29  private
30 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
31 
32 !==============================================================================
33 
35  type :: t_xdmf_grid
36  type(sll_t_xml_element), pointer :: xml_grid
37  sll_int32, allocatable :: dims(:)
38  end type t_xdmf_grid
39 
40  !----------------------------------------------------------------------------
43 
44  sll_real64 :: time = 0.0_f64
45  type(sll_t_xml_document), allocatable :: xml_doc
46  type(sll_t_xml_element), pointer :: xml_domain => null()
47  type(t_xdmf_grid), allocatable :: grids(:)
48 
49  contains
50 
51  procedure :: init => t_xdmf__init
52  procedure :: write => t_xdmf__write
53  procedure :: delete => t_xdmf__delete
54  procedure :: add_grid => t_xdmf__add_grid
55  procedure :: add_field => t_xdmf__add_field
56 
57  end type sll_t_xdmf_file
58 
59 !==============================================================================
60 contains
61 !==============================================================================
62 
63  !----------------------------------------------------------------------------
65  subroutine t_xdmf__init(self, time)
66  class(sll_t_xdmf_file), intent(out) :: self
67  sll_real64, intent(in) :: time
68 
69  type(sll_t_xml_element), pointer :: root
70 
71  ! Fill in fields
72  self%time = time
73 
74  ! Allocate XML document
75  allocate (self%xml_doc)
76 
77  ! Add header
78  call self%xml_doc%add_header_line("<?xml version='1.0' ?>")
79  call self%xml_doc%add_header_line("<!DOCTYPE Xdmf SYSTEM 'Xdmf.dtd' []>")
80 
81  ! Create root element (only one may exist!)
82  root => self%xml_doc%new_element('Xdmf')
83  call root%add_attribute('Version', '2.0')
84 
85  ! Add Domain element to root (only one for our purposes)
86  self%xml_domain => root%new_element('Domain')
87 
88  end subroutine t_xdmf__init
89 
90  !----------------------------------------------------------------------------
92  subroutine t_xdmf__write(self, fname)
93  class(sll_t_xdmf_file), intent(in) :: self
94  character(len=*), intent(in) :: fname
95 
96  call self%xml_doc%write(fname)
97 
98  end subroutine t_xdmf__write
99 
100  !----------------------------------------------------------------------------
102  subroutine t_xdmf__delete(self)
103  class(sll_t_xdmf_file), intent(inout) :: self
104 
105  self%time = 0.0_f64
106  self%xml_domain => null()
107 
108  if (allocated(self%xml_doc)) then
109  call self%xml_doc%delete()
110  deallocate (self%xml_doc)
111  end if
112 
113  if (allocated(self%grids)) then
114  deallocate (self%grids)
115  end if
116 
117  end subroutine t_xdmf__delete
118 
119  !----------------------------------------------------------------------------
121  subroutine t_xdmf__add_grid(self, grid_name, x1_path, x2_path, dims, gid)
122  class(sll_t_xdmf_file), intent(inout) :: self
123  character(len=*), intent(in) :: grid_name
124  character(len=*), intent(in) :: x1_path
125  character(len=*), intent(in) :: x2_path
126  sll_int32, intent(in) :: dims(2)
127  sll_int32, intent(out) :: gid
128 
129  sll_int32 :: ng
130  character(len=64) :: time_str
131  character(len=:), allocatable :: dims_str
132  type(t_xdmf_grid), allocatable :: tmp(:)
133  type(sll_t_xml_element), pointer :: grid, geometry
134  type(sll_t_xml_element), pointer :: time, topology, dataitem
135 
136  ! Prepare strings with data
137  write (time_str, '(g20.12)') self%time ! TODO: investigate format options
138  call sll_s_ints_to_string((/dims(2), dims(1)/), dims_str)
139 
140  ! Add new grid to domain
141  grid => self%xml_domain%new_element('Grid')
142  call grid%add_attribute('Name', trim(grid_name))
143  call grid%add_attribute('GridType', 'Uniform') ! only option for now
144 
145  ! Add time to grid
146  time => grid%new_element('Time')
147  call time%add_attribute('Value', trim(adjustl(time_str)))
148 
149  ! Add topology to grid
150  topology => grid%new_element('Topology')
151  call topology%add_attribute('TopologyType', '2DSMesh') ! only option now
152  call topology%add_attribute('NumberOfElements', dims_str)
153 
154  ! Add geometry to grid
155  geometry => grid%new_element('Geometry')
156  call geometry%add_attribute('GeometryType', 'X_Y') ! only option for now
157 
158  ! Add X axis to geometry
159  dataitem => geometry%new_element('DataItem')
160  call dataitem%add_attribute('Dimensions', dims_str)
161  call dataitem%add_attribute('NumberType', 'Float')
162  call dataitem%add_attribute('Precision', '8')
163  call dataitem%add_attribute('Format', 'HDF') ! only option for now
164  call dataitem%add_chardata(trim(x1_path)) ! only option for now
165 
166  ! Add Y axis to geometry
167  dataitem => geometry%new_element('DataItem')
168  call dataitem%add_attribute('Dimensions', dims_str)
169  call dataitem%add_attribute('NumberType', 'Float')
170  call dataitem%add_attribute('Precision', '8')
171  call dataitem%add_attribute('Format', 'HDF') ! only option for now
172  call dataitem%add_chardata(trim(x2_path)) ! only option for now
173 
174  ! Determine size of 'self%grids' array
175  if (allocated(self%grids)) then
176  ng = size(self%grids)
177  allocate (tmp(ng + 1))
178  tmp(1:ng) = self%grids(1:ng)
179  else
180  ng = 0
181  allocate (tmp(1))
182  end if
183 
184  ! Extend 'self%grids' array and store pointer to new grid
185  allocate (tmp(ng + 1)%dims(size(dims)))
186  tmp(ng + 1)%xml_grid => grid
187  tmp(ng + 1)%dims = dims
188  call move_alloc(from=tmp, to=self%grids)
189 
190  ! Return the grid id, i.e. the last index in the 'self%grids' array
191  gid = ng + 1
192 
193  end subroutine t_xdmf__add_grid
194 
195  !----------------------------------------------------------------------------
197  subroutine t_xdmf__add_field(self, grid_id, field_name, field_path)
198  class(sll_t_xdmf_file), intent(inout) :: self
199  sll_int32, intent(in) :: grid_id
200  character(len=*), intent(in) :: field_name
201  character(len=*), intent(in) :: field_path
202 
203  character(len=:), allocatable :: dims_str
204  type(sll_t_xml_element), pointer :: field, dataitem
205 
206  ! Prepare strings with data
207  call sll_s_ints_to_string(self%grids(grid_id)%dims, dims_str)
208 
209  ! Create new field (scalar, nodal)
210  field => self%grids(grid_id)%xml_grid%new_element('Attribute')
211  call field%add_attribute('Name', trim(field_name))
212  call field%add_attribute('AttributeType', 'Scalar') ! only option for now
213  call field%add_attribute('Center', 'Node') ! only option for now
214 
215  ! Add new dataitem
216  dataitem => field%new_element('DataItem')
217  call dataitem%add_attribute('Dimensions', dims_str)
218  call dataitem%add_attribute('NumberType', 'Float')
219  call dataitem%add_attribute('Precision', '8')
220  call dataitem%add_attribute('Format', 'HDF') ! only option for now
221  call dataitem%add_chardata(trim(adjustl(field_path)))
222 
223  end subroutine t_xdmf__add_field
224 
225 end module sll_m_xdmf_light_serial
Collection of functions/subroutines operating on files and strings.
subroutine, public sll_s_ints_to_string(ints, str)
Write an array of integers to a single string: . Numbers are separated by a blank space; ....
Construct the XML component of an XDMF database (sequential).
subroutine t_xdmf__init(self, time)
Initialize XDMF file (set time, allocate storage, store pointer to domain)
subroutine t_xdmf__delete(self)
Delete XDMF file (deallocate storage, nullify pointers)
subroutine t_xdmf__write(self, fname)
Write XML file.
subroutine t_xdmf__add_grid(self, grid_name, x1_path, x2_path, dims, gid)
Add grid to file (new grid ID is returned)
subroutine t_xdmf__add_field(self, grid_id, field_name, field_path)
Add field to grid (selected by its grid ID)
Facilities for constructing an XML tree and printing it to file.
Definition: sll_m_xml.F90:13
XML document type.
Definition: sll_m_xml.F90:87
    Report Typos and Errors