Afivo 0.3
All Classes Namespaces Functions Variables Pages
m_write_silo.f90
1!> This module contains wrapper functions to simplify writing Silo files
2!> \todo Document the functions in this module
4
5 implicit none
6 private
7
8 include 'silo_f9x.inc'
9
10 integer, parameter :: dp = kind(0.0d0)
11 integer, parameter :: line_len = 200
12 integer, parameter :: db_type = db_pdb
13
14 public :: silo_mkdir
15 public :: silo_create_file
16 public :: silo_open_file
17 public :: silo_close_file
18 public :: silo_set_time_varying
19 public :: silo_add_curve
20 public :: silo_add_grid
21 public :: silo_add_var
22 public :: silo_set_mmesh_grid
23 public :: silo_set_mmesh_var
24
25contains
26
27 subroutine silo_create_file(filename, dbix)
28 character(len=*), intent(in) :: filename
29 integer, intent(out) :: dbix
30 integer :: ierr
31 character(len=line_len) :: fileinfo
32
33 fileinfo = "A silo file"
34 ierr = dbcreate(trim(filename), len_trim(filename), db_clobber, &
35 db_local, fileinfo, len_trim(fileinfo), db_type, dbix)
36 if (ierr /= 0) then
37 print *, "Error creating file ", trim(filename)
38 error stop
39 end if
40 end subroutine silo_create_file
41
42 subroutine silo_open_file(filename, dbix)
43 character(len=*), intent(in) :: filename
44 integer :: dbix, ierr
45
46 ierr = dbopen(trim(filename), len_trim(filename), db_type, &
47 db_append, dbix)
48 if (ierr /= 0) then
49 print *, "Error opening file ", trim(filename)
50 error stop
51 end if
52 end subroutine silo_open_file
53
54 subroutine silo_close_file(dbix)
55 integer, intent(in) :: dbix
56 integer :: ierr
57
58 ierr = dbclose(dbix)
59 if (ierr /= 0) then
60 print *, "Error closing file with index", dbix
61 error stop
62 end if
63 end subroutine silo_close_file
64
65 subroutine silo_mkdir(dbix, dirname)
66 character(len=*), intent(in) :: dirname
67 integer, intent(in) :: dbix
68 integer :: ierr, iostat
69
70 ierr = dbmkdir(dbix, trim(dirname), len_trim(dirname), iostat)
71 if (ierr /= 0) then
72 print *, "Error creating directory ", dirname
73 error stop
74 end if
75 end subroutine silo_mkdir
76
77 !> Write two entries to the Silo file so that Visit treats it as a
78 !> time-varying database
79 subroutine silo_set_time_varying(dbix)
80 integer, intent(in) :: dbix
81 integer :: ierr
82 integer, parameter :: int_bool(1) = 1
83 integer, parameter :: dims(1) = 1
84 character(len=*), parameter :: name1 = "/ConnectivityIsTimeVarying"
85 character(len=*), parameter :: name2 = "/MetadataIsTimeVarying"
86
87 interface
88 integer function dbwrite(dbid, varname, lvarname, var, dims, &
89 ndims, datatype)
90 use, intrinsic :: iso_c_binding
91 integer(c_int) :: dbid, var(*), lvarname, dims(*), ndims, datatype
92 character(kind=c_char) :: varname(*)
93 end function dbwrite
94 end interface
95
96 ierr = dbwrite(dbix, name1, len(name1), int_bool, dims, 1, db_int);
97 if (ierr /= 0) print *, "Error writing ", trim(name1)
98 ierr = dbwrite(dbix, name2, len(name2), int_bool, dims, 1, db_int);
99 if (ierr /= 0) print *, "Error writing ", trim(name2)
100 end subroutine silo_set_time_varying
101
102 !> Add a curve-object (pairs of x-y values) to the Silo file
103 subroutine silo_add_curve(dbix, curvename, xvals, yvals)
104 character(len=*), intent(in) :: curvename
105 integer, intent(in) :: dbix
106 real(dp), intent(in) :: xvals(:), yvals(:)
107 integer :: iostat, ierr, dboptix
108
109 interface
110 integer (c_int) function dbputcurve(dbid, curvename, lcurvename, &
111 xvals, yvals, datatype, npoints, optlist_id, status)
112 use, intrinsic :: iso_c_binding
113 integer(c_int) :: dbid, lcurvename, datatype, npoints, status, optlist_id
114 real(c_double) :: xvals(*), yvals(*)
115 character(kind=c_char) :: curvename(*)
116 end function dbputcurve
117
118 integer (c_int) function dbaddiopt(optlist_id, option, ivalue)
119 use, intrinsic :: iso_c_binding
120 integer(c_int), intent(in) :: optlist_id, option, ivalue(*)
121 end function dbaddiopt
122 end interface
123
124 if (size(shape(xvals)) /= 1 .or. size(shape(yvals)) /= 1) then
125 print *, "Input is not a one-dimensional array. "
126 end if
127
128 if (size(xvals) /= size(yvals)) then
129 print *, "x and y arrays not of the same dimensions"
130 print *, size(xvals), " and ", size(yvals)
131 end if
132
133 ierr = dbmkoptlist(20, dboptix)
134 if (ierr /= 0) print *, &
135 "Error creating options list in SILO_add_curve ", dboptix
136
137 ierr = dbputcurve(dbix, trim(curvename), len_trim(curvename), &
138 xvals, yvals, db_double, size(xvals), dboptix, iostat)
139
140 if (ierr /= 0) print *, &
141 "ERROR: curve object not added to SILO file"
142
143 ierr = dbfreeoptlist(dboptix)
144 if (ierr /= 0) print *, &
145 "Error dbfreeoptlist in SILO_add_curve", ierr
146 end subroutine silo_add_curve
147
148 subroutine silo_add_grid(dbix, gridname, n_dim, N_r, r_min, dr, &
149 lo_offset, hi_offset, n_cycle)
150 character(len=*), intent(in) :: gridname
151 integer, intent(in) :: dbix, n_dim, n_r(:)
152 integer, intent(in) :: lo_offset(n_dim), hi_offset(n_dim)
153 real(dp), intent(in) :: r_min(:), dr(:)
154 real(dp), allocatable :: x_coords(:), y_coords(:), z_coords(:)
155 integer :: i, ierr, iostat, dboptix
156 integer, optional, intent(in) :: n_cycle
157
158 interface
159 function dbputqm(dbid, name, lname, xname, lxname, yname, &
160 lyname, zname, lzname, x, y, z, dims, ndims, &
161 datatype, coordtype, optlist_id, status)
162 use, intrinsic :: iso_c_binding
163 integer(c_int) :: dbid, lname, lxname, lyname, lzname, dims(*), ndims
164 integer(c_int) :: datatype, coordtype, optlist_id, status, dbputqm
165 real(c_double) :: x(*), y(*), z(*)
166 character(kind=c_char) :: name(*), xname(*), yname(*), zname(*)
167 end function dbputqm
168
169 integer (c_int) function dbaddiopt(optlist_id, option, ivalue)
170 use, intrinsic :: iso_c_binding
171 integer(c_int), intent(in) :: optlist_id, option, ivalue(*)
172 end function dbaddiopt
173 end interface
174
175 if (n_dim < 1 .or. n_dim > 3) then
176 print *, "Cannot add grid for which n_dim < 1 or n_dim > 3"
177 return
178 end if
179
180 allocate(x_coords(n_r(1)))
181 do i = 1, n_r(1)
182 x_coords(i) = r_min(1) + (i-1) * dr(1)
183 end do
184
185 if (n_dim > 1) then
186 allocate(y_coords(n_r(2)))
187 do i = 1, n_r(2)
188 y_coords(i) = r_min(2) + (i-1) * dr(2)
189 end do
190 else
191 allocate(y_coords(0))
192 end if
193
194 if (n_dim > 2) then
195 allocate(z_coords(n_r(3)))
196 do i = 1, n_r(3)
197 z_coords(i) = r_min(3) + (i-1) * dr(3)
198 end do
199 else
200 allocate(z_coords(0))
201 end if
202
203 ! Make option list
204 ierr = dbmkoptlist(20, dboptix)
205 if (ierr /= 0) print *, &
206 "Error creating options list in SILO_add_grid ", dboptix
207
208 ! Make iteration number explicit
209 ierr = dbaddiopt(dboptix, dbopt_cycle, [n_cycle])
210 if (ierr /= 0) print *, &
211 "Error passing iteration number SILO_add_grid ", dboptix
212
213 ! Set integer options
214 ierr = dbaddiopt(dboptix, dbopt_nspace, [n_dim])
215 if (ierr /= 0) print *, &
216 "Error dbaddiopt in SILO_add_grid: DBOPT_NSPACE", ierr
217
218 ierr = dbaddiopt(dboptix, dbopt_lo_offset, lo_offset)
219 if (ierr /= 0) print *, &
220 "Error dbaddiopt in SILO_add_grid: DBOPT_LO_OFFSET", ierr
221
222 ierr = dbaddiopt(dboptix, dbopt_hi_offset, hi_offset)
223 if (ierr /= 0) print *, &
224 "Error dbaddiopt in SILO_add_grid: DBOPT_HI_OFFSET", ierr
225
226 ! Write the grid structure
227 ierr = dbputqm(dbix, trim(gridname), len_trim(gridname), &
228 'x', 1, 'y', 1, 'z', 1, x_coords, y_coords, z_coords, &
229 n_r, n_dim, db_double, db_collinear, dboptix, iostat)
230 if (ierr /= 0) print *, &
231 "Error dbputqm is SILO_add_grid", ierr
232
233 ierr = dbfreeoptlist(dboptix)
234 if (ierr /= 0) print *, &
235 "Error dbfreeoptlist is SILO_add_grid", ierr
236 end subroutine silo_add_grid
237
238 subroutine silo_add_var(dbix, dataname, gridname, &
239 d_packed, d_shape, n_cycle)
240 character(len=*), intent(in) :: gridname, dataname
241 real(dp), intent(in) :: d_packed(:)
242 integer, intent(in) :: dbix, d_shape(:)
243
244 integer :: dboptix, ierr, iostat
245 real(dp) :: dummy(1)
246 integer, optional, intent(in) :: n_cycle
247
248 interface
249 function dbputqv1(dbid, name, lname, meshname, lmeshname, &
250 var, dims, ndims, mixvar, mixlen, datatype, &
251 centering, optlist_id, status)
252 use, intrinsic :: iso_c_binding
253 integer(c_int) :: dbid, lname, lmeshname, dims(*), ndims, mixlen
254 integer(c_int) :: centering, optlist_id, status, datatype, dbputqv1
255 real(c_double) :: var(*), mixvar(*)
256 character(kind=c_char) :: name(*), meshname(*)
257 end function dbputqv1
258
259 integer (c_int) function dbaddiopt(optlist_id, option, ivalue)
260 use, intrinsic :: iso_c_binding
261 integer(c_int), intent(in) :: optlist_id, option, ivalue(*)
262 end function dbaddiopt
263 end interface
264
265 if (size(d_packed) /= product(d_shape)) then
266 error stop "Error: d_packed does not correspond to d_shape"
267 end if
268
269 if (size(d_shape) < 1 .or. size(d_shape) > 3) then
270 error stop "Error: size(d_shape) < 1 or size(d_shape) > 3"
271 end if
272
273 ierr = dbmkoptlist(10, dboptix)
274 if (ierr /= 0) then
275 error stop "Error creating options list in SILO_add_var"
276 end if
277
278 ! Make iteration number explicit
279 if (present(n_cycle)) then
280 ierr = dbaddiopt(dboptix, dbopt_cycle, [n_cycle])
281 if (ierr /= 0) print *, &
282 "Error passing iteration number SILO_add_var ", dboptix
283 end if
284
285 ierr = dbaddiopt(dboptix, dbopt_hide_from_gui, [1])
286
287 ! Write the data to the grid
288 ierr = dbputqv1(dbix, trim(dataname), len_trim(dataname), &
289 trim(gridname), len_trim(gridname), d_packed, d_shape, &
290 size(d_shape), dummy, 0, db_double, db_zonecent, dboptix, iostat)
291 if (ierr /= 0) then
292 print *, "Error dbputqv1 in SILO_add_var", ierr
293 error stop
294 end if
295
296 ierr = dbfreeoptlist(dboptix)
297 end subroutine silo_add_var
298
299 subroutine silo_set_mmesh_grid(dbix, mmname, gridnames, n_cycle, time)
300 character(len=*), intent(in) :: mmname, gridnames(:)
301 integer, intent(in) :: dbix
302 integer, intent(in), optional :: n_cycle
303 real(dp), intent(in), optional :: time
304
305 integer :: i, ierr,length
306 integer :: dboptix, iostat, old_str_len
307 integer :: n_grids, name_len, total_len
308 integer, allocatable :: m_types(:), name_lengths(:)
309 character(len=:), allocatable :: mnames
310
311 interface
312 function dbputmmesh(dbid, name, lname, nmesh, meshnames, lmeshnames, &
313 meshtypes, optlist_id, status)
314 use, intrinsic :: iso_c_binding
315 integer(c_int) :: dbputmmesh, lmeshnames(*)
316 integer(c_int) :: dbid, lname, nmesh, meshtypes(*), optlist_id, status
317 character(kind=c_char) :: name(*), meshnames(*)
318 end function dbputmmesh
319 end interface
320
321 n_grids = size(gridnames)
322 if (n_grids < 1) then
323 error stop "SILO_set_mmesh_grid: error too few grids (<1)"
324 end if
325
326 name_len = len(gridnames(1))
327 total_len = name_len * n_grids
328 allocate(character(total_len) :: mnames)
329 allocate(name_lengths(n_grids))
330 allocate(m_types(n_grids))
331
332 do i = 1, n_grids
333 mnames((i-1)*name_len+1:i*name_len) = trim(gridnames(i)) // char(0)
334 end do
335
336 old_str_len = dbset2dstrlen(name_len)
337 m_types = db_quad_rect
338 name_lengths = name_len
339
340 ierr = dbmkoptlist(10, dboptix)
341
342 if (present(n_cycle)) then
343 ierr = dbaddiopt(dboptix, dbopt_cycle, n_cycle)
344 end if
345
346 if (present(time)) then
347 ierr = dbadddopt(dboptix, dbopt_dtime, time)
348 end if
349
350 ierr = dbputmmesh(dbix, trim(mmname), len_trim(mmname), n_grids, &
351 mnames(1:total_len), name_lengths, m_types, dboptix, iostat)
352 if (ierr /= 0) then
353 error stop "Error calling dbputmmesh"
354 end if
355
356 ierr = dbfreeoptlist(dboptix)
357 length = dbset2dstrlen(old_str_len)
358 end subroutine silo_set_mmesh_grid
359
360 subroutine silo_set_mmesh_var(dbix, mvname, mmname, datanames, n_cycle, time)
361 character(len=*), intent(in) :: mvname, mmname, datanames(:)
362 integer, intent(in) :: dbix
363 integer, intent(in), optional :: n_cycle
364 real(dp), intent(in), optional :: time
365
366
367 integer :: i, ierr, dboptix, iostat,length
368 integer :: old_str_len, n_grids, name_len, total_len
369 integer, allocatable :: m_types(:), name_lengths(:)
370 character(:), allocatable :: dnames
371
372 interface
373 function dbputmvar(dbid, name, lname, nlevels, meshnames, &
374 lmnames, meshtypes, optlist_id, status)
375 use, intrinsic :: iso_c_binding
376 integer(c_int) :: dbputmvar, lmnames(*)
377 integer(c_int) :: dbid, lname, nlevels, meshtypes(*)
378 integer(c_int) :: optlist_id, status
379 character(kind=c_char) :: name(*), meshnames(*)
380 end function dbputmvar
381 end interface
382
383 n_grids = size(datanames)
384 if (n_grids < 1) then
385 error stop "SILO_set_mmesh_var: error too few grids (<1)"
386 end if
387
388 name_len = len(datanames(1))
389 total_len = name_len * n_grids
390 allocate(character(total_len) :: dnames)
391 allocate(name_lengths(n_grids))
392 allocate(m_types(n_grids))
393
394 do i = 1, n_grids
395 dnames((i-1)*name_len+1:i*name_len) = trim(datanames(i)) // char(0)
396 end do
397 old_str_len = dbset2dstrlen(name_len)
398 m_types = db_quadvar
399 name_lengths = name_len
400
401 ierr = dbmkoptlist(10, dboptix)
402 if (ierr /= 0) then
403 error stop "Error creating options list in SILO_set_mmesh_var"
404 end if
405
406 if (present(n_cycle)) then
407 ierr = dbaddiopt(dboptix, dbopt_cycle, n_cycle)
408 end if
409
410 if (present(time)) then
411 ierr = dbadddopt(dboptix, dbopt_dtime, time)
412 end if
413
414 ! Jannis: valgrind complains about a memory leak due to this call, but if
415 ! that leak is there, it seems to be Silo's fault.
416 ierr = dbaddcopt(dboptix, dbopt_mmesh_name, &
417 trim(mmname), len_trim(mmname))
418 if (ierr /= 0) print *, &
419 "Error dbaddiopt in SILO_set_mmesh_var: DBOPT_MMESH_NAME", ierr
420
421 ierr = dbputmvar(dbix, trim(mvname), len_trim(mvname), n_grids, &
422 dnames(1:total_len), name_lengths, m_types, dboptix, iostat)
423 if (ierr /= 0) then
424 error stop "Error calling dbputmvar"
425 end if
426
427 ierr = dbfreeoptlist(dboptix)
428 length = dbset2dstrlen(old_str_len)
429 end subroutine silo_set_mmesh_var
430
431 ! !> Add a curve-object (pairs of x-y values) to the Silo file
432 ! subroutine SILO_add_curve(dbix, curvename, xvals, yvals, xname, yname)
433 ! character(len=*), intent(in) :: curvename
434 ! integer, intent(in) :: dbix
435 ! real(dp), intent(in) :: xvals(:), yvals(:)
436 ! character(len=*), intent(in) :: xname, yname
437 ! integer :: iostat, ierr, dboptix
438 !
439 ! interface
440 ! integer (c_int) function dbputcurve(dbid, curvename, lcurvename, &
441 ! xvals, yvals, datatype, npoints, optlist_id, status)
442 ! use, intrinsic :: iso_c_binding
443 ! integer(c_int) :: dbid, lcurvename, datatype, npoints, status, optlist_id
444 ! real(c_double) :: xvals(*), yvals(*)
445 ! character(kind=c_char) :: curvename(*)
446 ! end function dbputcurve
447 !
448 ! integer (c_int) function dbaddiopt(optlist_id, option, ivalue)
449 ! use, intrinsic :: iso_c_binding
450 ! integer(c_int), intent(in) :: optlist_id, option, ivalue(*)
451 ! end function dbaddiopt
452 ! end interface
453 !
454 ! if (size(xvals) /= size(yvals)) &
455 ! error stop "SILO_add_curve: x and y arrays unequal size"
456 !
457 ! ierr = dbmkoptlist(20, dboptix)
458 ! if (ierr /= 0) error stop "creating options list in SILO_add_curve "
459 !
460 ! ierr = dbaddcopt(dboptix, DBOPT_XLABEL, trim(xname), len_trim(xname))
461 ! if (ierr /= 0) error stop "adding option in SILO_add_curve"
462 ! ierr = dbaddcopt(dboptix, DBOPT_YLABEL, trim(yname), len_trim(yname))
463 ! if (ierr /= 0) error stop "adding option in SILO_add_curve"
464 !
465 ! ierr = dbputcurve(dbix, trim(curvename), len_trim(curvename), &
466 ! xvals, yvals, DB_DOUBLE, size(xvals), dboptix, iostat)
467 ! if (ierr /= 0) error stop "curve object not added to SILO file"
468 !
469 ! ierr = dbfreeoptlist(dboptix)
470 ! if (ierr /= 0) error stop "dbfreeoptlist in SILO_add_curve"
471 ! end subroutine SILO_add_curve
472
473end module m_write_silo
This module contains wrapper functions to simplify writing Silo files.
Definition: m_write_silo.f90:3