MFC
Exascale flow solver
Loading...
Searching...
No Matches
m_data_output.fpp.f90
Go to the documentation of this file.
1# 1 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
2!>
3!! @file
4!! @brief Contains module m_data_output
5
6!> @brief Writes post-processed grid and flow-variable data to Silo-HDF5 or binary database files
8
12 use m_mpi_proxy
14 use m_helper
16
17 implicit none
18
25
26 ! Include Silo-HDF5 interface library
27 include 'silo_f9x.inc'
28
29 ! Flow variable storage; q_root_sf gathers to rank 0 in 1D parallel runs
30 real(wp), allocatable, dimension(:,:,:), public :: q_sf
31 real(wp), allocatable, dimension(:,:,:) :: q_root_sf
32 real(wp), allocatable, dimension(:,:,:) :: cyl_q_sf
33
34 ! Single precision storage for flow variables
35 real(sp), allocatable, dimension(:,:,:), public :: q_sf_s
36 real(sp), allocatable, dimension(:,:,:) :: q_root_sf_s
37 real(sp), allocatable, dimension(:,:,:) :: cyl_q_sf_s
38
39 ! Spatial and data extents for VisIt visualization
40 real(wp), allocatable, dimension(:,:) :: spatial_extents
41 real(wp), allocatable, dimension(:,:) :: data_extents
42
43 ! Ghost zone layer sizes (lo/hi) for subdomain connectivity in VisIt
44 integer, allocatable, dimension(:) :: lo_offset
45 integer, allocatable, dimension(:) :: hi_offset
46
47 ! Track cell-boundary count per active coordinate direction
48 integer, allocatable, dimension(:) :: dims
49
50 ! Locations of various folders in the case's directory tree, associated with the choice of the formatted database format. These
51 ! include, in order, the location of the folder named after the selected formatted database format, and the locations of two
52 ! sub-directories of the latter, the first of which is named after the local processor rank, while the second is named 'root'.
53 ! The folder associated with the local processor rank contains only the data pertaining to the part of the domain taken care of
54 ! by the local processor. The root directory, on the other hand, will contain either the information about the connectivity
55 ! required to put the entire domain back together, or the actual data associated with the entire computational domain. This all
56 ! depends on dimensionality and the choice of the formatted database format.
57 character(LEN=path_len + name_len) :: dbdir
58 character(LEN=path_len + 2*name_len) :: proc_rank_dir
59 character(LEN=path_len + 2*name_len) :: rootdir
60
61 ! Handles of the formatted database master/root file, slave/local processor file and options list. The list of options is
62 ! explicitly used in the Silo- HDF5 database format to provide additional details about the contents of a formatted database
63 ! file, such as the previously described spatial and data extents.
64 integer :: dbroot
65 integer :: dbfile
66 integer :: optlist
67
68 ! The total number of flow variable(s) to be stored in a formatted database file. Note that this is only needed when using the
69 ! Binary format.
70 integer :: dbvars
71
72 ! Generic error flags utilized in the handling, checking and the reporting of the input and output operations errors with a
73 ! formatted database file
74 integer, private :: err
75
76contains
77
78 !> Allocate storage arrays, configure output directories, and count flow variables for formatted database output.
80
81 character(LEN=len_trim(case_dir) + 2*name_len) :: file_loc
82 logical :: dir_check
83 integer :: i
84
85 allocate (q_sf(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end))
86 if (grid_geometry == 3) then
87 allocate (cyl_q_sf(-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end,-offset_x%beg:m + offset_x%end))
88 end if
89
90 if (precision == 1) then
91 allocate (q_sf_s(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end))
92 if (grid_geometry == 3) then
93 allocate (cyl_q_sf_s(-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end,-offset_x%beg:m + offset_x%end))
94 end if
95 end if
96
97 if (n == 0) then
98 allocate (q_root_sf(0:m_root,0:0,0:0))
99 if (precision == 1) then
100 allocate (q_root_sf_s(0:m_root,0:0,0:0))
101 end if
102 end if
103
104 ! Allocating the spatial and data extents and also the variables for the offsets and the one bookkeeping the number of
105 ! cell-boundaries in each active coordinate direction. Note that all these variables are only needed by the Silo-HDF5 format
106 ! for multidimensional data.
107 if (format == 1) then
108 allocate (data_extents(1:2,0:num_procs - 1))
109
110 if (p > 0) then
111 allocate (spatial_extents(1:6,0:num_procs - 1))
112 allocate (lo_offset(1:3))
113 allocate (hi_offset(1:3))
114 allocate (dims(1:3))
115 else if (n > 0) then
116 allocate (spatial_extents(1:4,0:num_procs - 1))
117 allocate (lo_offset(1:2))
118 allocate (hi_offset(1:2))
119 allocate (dims(1:2))
120 else
121 allocate (spatial_extents(1:2,0:num_procs - 1))
122 allocate (lo_offset(1:1))
123 allocate (hi_offset(1:1))
124 allocate (dims(1:1))
125 end if
126 end if
127
128 ! The size of the ghost zone layer in each of the active coordinate directions was set in the module m_mpi_proxy.f90. The
129 ! results are now transferred to the local variables of this module when they are required by the Silo-HDF5 format, for
130 ! multidimensional data sets. With the same, latter, requirements, the variables bookkeeping the number of cell-boundaries
131 ! in each active coordinate direction are also set here.
132 if (format == 1) then
133 if (p > 0) then
134 if (grid_geometry == 3) then
135 lo_offset(:) = (/offset_y%beg, offset_z%beg, offset_x%beg/)
136 hi_offset(:) = (/offset_y%end, offset_z%end, offset_x%end/)
137 else
138 lo_offset(:) = (/offset_x%beg, offset_y%beg, offset_z%beg/)
139 hi_offset(:) = (/offset_x%end, offset_y%end, offset_z%end/)
140 end if
141
142 if (grid_geometry == 3) then
143 dims(:) = (/n + offset_y%beg + offset_y%end + 2, p + offset_z%beg + offset_z%end + 2, &
144 & m + offset_x%beg + offset_x%end + 2/)
145 else
146 dims(:) = (/m + offset_x%beg + offset_x%end + 2, n + offset_y%beg + offset_y%end + 2, &
147 & p + offset_z%beg + offset_z%end + 2/)
148 end if
149 else if (n > 0) then
150 lo_offset(:) = (/offset_x%beg, offset_y%beg/)
151 hi_offset(:) = (/offset_x%end, offset_y%end/)
152
153 dims(:) = (/m + offset_x%beg + offset_x%end + 2, n + offset_y%beg + offset_y%end + 2/)
154 else
155 lo_offset(:) = (/offset_x%beg/)
156 hi_offset(:) = (/offset_x%end/)
157 dims(:) = (/m + offset_x%beg + offset_x%end + 2/)
158 end if
159 end if
160
161 if (format == 1) then
162 dbdir = trim(case_dir) // '/silo_hdf5'
163
164 write (proc_rank_dir, '(A,I0)') '/p', proc_rank
165
166 proc_rank_dir = trim(dbdir) // trim(proc_rank_dir)
167
168 file_loc = trim(proc_rank_dir) // '/.'
169
170 call my_inquire(file_loc, dir_check)
171 if (dir_check .neqv. .true.) then
173 end if
174
175 if (proc_rank == 0) then
176 rootdir = trim(dbdir) // '/root'
177
178 file_loc = trim(rootdir) // '/.'
179
180 call my_inquire(file_loc, dir_check)
181 if (dir_check .neqv. .true.) then
182 call s_create_directory(trim(rootdir))
183 end if
184 end if
185 else
186 dbdir = trim(case_dir) // '/binary'
187
188 write (proc_rank_dir, '(A,I0)') '/p', proc_rank
189
190 proc_rank_dir = trim(dbdir) // trim(proc_rank_dir)
191
192 file_loc = trim(proc_rank_dir) // '/.'
193
194 call my_inquire(file_loc, dir_check)
195
196 if (dir_check .neqv. .true.) then
198 end if
199
200 if (n == 0 .and. proc_rank == 0) then
201 rootdir = trim(dbdir) // '/root'
202
203 file_loc = trim(rootdir) // '/.'
204
205 call my_inquire(file_loc, dir_check)
206
207 if (dir_check .neqv. .true.) then
208 call s_create_directory(trim(rootdir))
209 end if
210 end if
211 end if
212
213 if (bubbles_lagrange) then ! Lagrangian solver
214 if (lag_txt_wrt) then
215 dbdir = trim(case_dir) // '/lag_bubbles_post_process'
216 file_loc = trim(dbdir) // '/.'
217 call my_inquire(file_loc, dir_check)
218
219 if (dir_check .neqv. .true.) then
220 call s_create_directory(trim(dbdir))
221 end if
222 end if
223 end if
224
225 ! Contrary to the Silo-HDF5 database format, handles of the Binary database master/root and slave/local process files are
226 ! perfectly static throughout post-process. Hence, they are set here so that they do not have to be repetitively computed in
227 ! later procedures.
228 if (format == 2) then
229 if (n == 0 .and. proc_rank == 0) dbroot = 2
230 dbfile = 1
231 end if
232
233 if (format == 2) then
234 dbvars = 0
235
236 if ((model_eqns == 2) .or. (model_eqns == 3)) then
237 do i = 1, num_fluids
238 if (alpha_rho_wrt(i) .or. (cons_vars_wrt .or. prim_vars_wrt)) then
239 dbvars = dbvars + 1
240 end if
241 end do
242 end if
243
244 if ((rho_wrt .or. (model_eqns == 1 .and. (cons_vars_wrt .or. prim_vars_wrt))) .and. (.not. relativity)) then
245 dbvars = dbvars + 1
246 end if
247
248 if (relativity .and. (rho_wrt .or. prim_vars_wrt)) dbvars = dbvars + 1
249 if (relativity .and. (rho_wrt .or. cons_vars_wrt)) dbvars = dbvars + 1
250
251 do i = 1, eqn_idx%E - eqn_idx%mom%beg
252 if (mom_wrt(i) .or. cons_vars_wrt) dbvars = dbvars + 1
253 end do
254
255 do i = 1, eqn_idx%E - eqn_idx%mom%beg
256 if (vel_wrt(i) .or. prim_vars_wrt) dbvars = dbvars + 1
257 end do
258
259 do i = 1, eqn_idx%E - eqn_idx%mom%beg
260 if (flux_wrt(i)) dbvars = dbvars + 1
261 end do
262
263 if (e_wrt .or. cons_vars_wrt) dbvars = dbvars + 1
264 if (pres_wrt .or. prim_vars_wrt) dbvars = dbvars + 1
265 if (hypoelasticity) dbvars = dbvars + (num_dims*(num_dims + 1))/2
266 if (cont_damage) dbvars = dbvars + 1
267 if (hyper_cleaning) dbvars = dbvars + 1
268
269 if (mhd) then
270 if (n == 0) then
271 dbvars = dbvars + 2
272 else
273 dbvars = dbvars + 3
274 end if
275 end if
276
277 if ((model_eqns == 2) .or. (model_eqns == 3)) then
278 do i = 1, num_fluids - 1
279 if (alpha_wrt(i) .or. (cons_vars_wrt .or. prim_vars_wrt)) then
280 dbvars = dbvars + 1
281 end if
282 end do
283
284 if (alpha_wrt(num_fluids) .or. (cons_vars_wrt .or. prim_vars_wrt)) then
285 dbvars = dbvars + 1
286 end if
287 end if
288
289 if (gamma_wrt .or. (model_eqns == 1 .and. (cons_vars_wrt .or. prim_vars_wrt))) then
290 dbvars = dbvars + 1
291 end if
292
293 if (heat_ratio_wrt) dbvars = dbvars + 1
294
295 if (pi_inf_wrt .or. (model_eqns == 1 .and. (cons_vars_wrt .or. prim_vars_wrt))) then
296 dbvars = dbvars + 1
297 end if
298
299 if (pres_inf_wrt) dbvars = dbvars + 1
300 if (c_wrt) dbvars = dbvars + 1
301
302 if (p > 0) then
303 do i = 1, num_vels
304 if (omega_wrt(i)) dbvars = dbvars + 1
305 end do
306 else if (n > 0) then
307 do i = 1, num_vels
308 if (omega_wrt(i)) dbvars = dbvars + 1
309 end do
310 end if
311
312 if (schlieren_wrt) dbvars = dbvars + 1
313 end if
314
316
317 !> Compute the cell-index bounds for the user-specified partial output domain in each coordinate direction.
318 impure subroutine s_define_output_region
319
320 integer :: i
321 integer :: lower_bound, upper_bound
322
323# 323 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
324 if (m == 0) return ! Early return for y or z if simulation is 1D or 2D
325
326 lower_bound = -offset_x%beg
327 upper_bound = m + offset_x%end
328
329 do i = lower_bound, upper_bound
330 if (x_cc(i) > x_output%beg) then
331 x_output_idx%beg = i + offset_x%beg
332 exit
333 end if
334 end do
335
336 do i = upper_bound, lower_bound, -1
337 if (x_cc(i) < x_output%end) then
338 x_output_idx%end = i + offset_x%beg
339 exit
340 end if
341 end do
342
343 ! If no grid points are within the output region
344 if ((x_cc(lower_bound) > x_output%end) .or. (x_cc(upper_bound) < x_output%beg)) then
345 x_output_idx%beg = 0
346 x_output_idx%end = 0
347 end if
348# 323 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
349 if (n == 0) return ! Early return for y or z if simulation is 1D or 2D
350
351 lower_bound = -offset_y%beg
352 upper_bound = n + offset_y%end
353
354 do i = lower_bound, upper_bound
355 if (y_cc(i) > y_output%beg) then
356 y_output_idx%beg = i + offset_y%beg
357 exit
358 end if
359 end do
360
361 do i = upper_bound, lower_bound, -1
362 if (y_cc(i) < y_output%end) then
363 y_output_idx%end = i + offset_y%beg
364 exit
365 end if
366 end do
367
368 ! If no grid points are within the output region
369 if ((y_cc(lower_bound) > y_output%end) .or. (y_cc(upper_bound) < y_output%beg)) then
370 y_output_idx%beg = 0
371 y_output_idx%end = 0
372 end if
373# 323 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
374 if (p == 0) return ! Early return for y or z if simulation is 1D or 2D
375
376 lower_bound = -offset_z%beg
377 upper_bound = p + offset_z%end
378
379 do i = lower_bound, upper_bound
380 if (z_cc(i) > z_output%beg) then
381 z_output_idx%beg = i + offset_z%beg
382 exit
383 end if
384 end do
385
386 do i = upper_bound, lower_bound, -1
387 if (z_cc(i) < z_output%end) then
388 z_output_idx%end = i + offset_z%beg
389 exit
390 end if
391 end do
392
393 ! If no grid points are within the output region
394 if ((z_cc(lower_bound) > z_output%end) .or. (z_cc(upper_bound) < z_output%beg)) then
395 z_output_idx%beg = 0
396 z_output_idx%end = 0
397 end if
398# 348 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
399
400 end subroutine s_define_output_region
401
402 !> Open (or create) the Silo-HDF5 or Binary formatted database slave and master files for a given time step.
403 impure subroutine s_open_formatted_database_file(t_step)
404
405 integer, intent(in) :: t_step
406 character(LEN=len_trim(case_dir) + 3*name_len) :: file_loc
407 integer :: ierr
408
409 if (format == 1) then
410 write (file_loc, '(A,I0,A)') '/', t_step, '.silo'
411 file_loc = trim(proc_rank_dir) // trim(file_loc)
412
413 ierr = dbcreate(trim(file_loc), len_trim(file_loc), db_clobber, db_local, 'MFC v3.0', 8, db_hdf5, dbfile)
414
415 if (dbfile == -1) then
416 call s_mpi_abort('Unable to create Silo-HDF5 database ' // 'slave file ' // trim(file_loc) // '. ' // 'Exiting.')
417 end if
418
419 if (proc_rank == 0) then
420 write (file_loc, '(A,I0,A)') '/collection_', t_step, '.silo'
421 file_loc = trim(rootdir) // trim(file_loc)
422
423 ierr = dbcreate(trim(file_loc), len_trim(file_loc), db_clobber, db_local, 'MFC v3.0', 8, db_hdf5, dbroot)
424
425 if (dbroot == -1) then
426 call s_mpi_abort('Unable to create Silo-HDF5 database ' // 'master file ' // trim(file_loc) // '. ' &
427 & // 'Exiting.')
428 end if
429 end if
430 else
431 write (file_loc, '(A,I0,A)') '/', t_step, '.dat'
432 file_loc = trim(proc_rank_dir) // trim(file_loc)
433
434 open (dbfile, iostat=err, file=trim(file_loc), form='unformatted', status='replace')
435
436 if (err /= 0) then
437 call s_mpi_abort('Unable to create Binary database slave ' // 'file ' // trim(file_loc) // '. Exiting.')
438 end if
439
440 if (output_partial_domain) then
441 write (dbfile) x_output_idx%end - x_output_idx%beg, y_output_idx%end - y_output_idx%beg, &
442 & z_output_idx%end - z_output_idx%beg, dbvars
443 else
444 write (dbfile) m, n, p, dbvars
445 end if
446
447 if (n == 0 .and. proc_rank == 0) then
448 write (file_loc, '(A,I0,A)') '/', t_step, '.dat'
449 file_loc = trim(rootdir) // trim(file_loc)
450
451 open (dbroot, iostat=err, file=trim(file_loc), form='unformatted', status='replace')
452
453 if (err /= 0) then
454 call s_mpi_abort('Unable to create Binary database ' // 'master file ' // trim(file_loc) // '. Exiting.')
455 end if
456
457 if (output_partial_domain) then
458 write (dbroot) x_output_idx%end - x_output_idx%beg, 0, 0, dbvars
459 else
460 write (dbroot) m_root, 0, 0, dbvars
461 end if
462 end if
463 end if
464
465 end subroutine s_open_formatted_database_file
466
467 !> Open the interface data file for appending extracted interface coordinates.
468 impure subroutine s_open_intf_data_file()
469
470 character(LEN=path_len + 3*name_len) :: file_path
471
472 write (file_path, '(A)') '/intf_data.dat'
473 file_path = trim(case_dir) // trim(file_path)
474
475 open (211, file=trim(file_path), form='formatted', position='append', status='unknown')
476
477 end subroutine s_open_intf_data_file
478
479 !> Open the energy data file for appending volume-integrated energy budget quantities.
480 impure subroutine s_open_energy_data_file()
481
482 character(LEN=path_len + 3*name_len) :: file_path
483
484 write (file_path, '(A)') '/eng_data.dat'
485 file_path = trim(case_dir) // trim(file_path)
486
487 open (251, file=trim(file_path), form='formatted', position='append', status='unknown')
488
489 end subroutine s_open_energy_data_file
490
491 !> Write the computational grid (cell-boundary coordinates) to the formatted database slave and master files.
493
494 integer, intent(in) :: t_step
495
496 ! NAG compiler requires these to be statically sized
497 character(LEN=4*name_len), dimension(num_procs) :: meshnames
498 integer, dimension(num_procs) :: meshtypes
499 integer :: i
500 integer :: ierr
501
502 if (format == 1) then
503 ! For multidimensional data sets, the spatial extents of all of the grid(s) handled by the local processor(s) are
504 ! recorded so that they may be written, by root processor, to the formatted database master file.
505 if (num_procs > 1) then
507 else if (p > 0) then
508 if (grid_geometry == 3) then
509 spatial_extents(:,0) = (/minval(y_cb), minval(z_cb), minval(x_cb), maxval(y_cb), maxval(z_cb), maxval(x_cb)/)
510 else
511 spatial_extents(:,0) = (/minval(x_cb), minval(y_cb), minval(z_cb), maxval(x_cb), maxval(y_cb), maxval(z_cb)/)
512 end if
513 else if (n > 0) then
514 spatial_extents(:,0) = (/minval(x_cb), minval(y_cb), maxval(x_cb), maxval(y_cb)/)
515 else
516 spatial_extents(:,0) = (/minval(x_cb), maxval(x_cb)/)
517 end if
518
519 ! Next, the root processor proceeds to record all of the spatial extents in the formatted database master file. In
520 ! addition, it also records a sub-domain connectivity map so that the entire grid may be reassembled by looking at the
521 ! master file.
522 if (proc_rank == 0) then
523 do i = 1, num_procs
524 write (meshnames(i), '(A,I0,A,I0,A)') '../p', i - 1, '/', t_step, '.silo:rectilinear_grid'
525 end do
526
527 meshtypes = db_quad_rect
528
529 err = dbset2dstrlen(len(meshnames(1)))
530 err = dbmkoptlist(2, optlist)
531 err = dbaddiopt(optlist, dbopt_extents_size, size(spatial_extents, 1))
532 err = dbadddopt(optlist, dbopt_extents, spatial_extents)
533 err = dbputmmesh(dbroot, 'rectilinear_grid', 16, num_procs, meshnames, len_trim(meshnames), meshtypes, optlist, &
534 & ierr)
535 err = dbfreeoptlist(optlist)
536 end if
537
538 ! Finally, the local quadrilateral mesh, either 2D or 3D, along with its offsets that indicate the presence and size of
539 ! ghost zone layer(s), are put in the formatted database slave file.
540
541 if (p > 0) then
542 err = dbmkoptlist(2, optlist)
543 err = dbaddiaopt(optlist, dbopt_lo_offset, size(lo_offset), lo_offset)
544 err = dbaddiaopt(optlist, dbopt_hi_offset, size(hi_offset), hi_offset)
545 if (grid_geometry == 3) then
546 err = dbputqm(dbfile, 'rectilinear_grid', 16, 'x', 1, 'y', 1, 'z', 1, y_cb, z_cb, x_cb, dims, 3, db_double, &
547 & db_collinear, optlist, ierr)
548 else
549 err = dbputqm(dbfile, 'rectilinear_grid', 16, 'x', 1, 'y', 1, 'z', 1, x_cb, y_cb, z_cb, dims, 3, db_double, &
550 & db_collinear, optlist, ierr)
551 end if
552 err = dbfreeoptlist(optlist)
553 else if (n > 0) then
554 err = dbmkoptlist(2, optlist)
555 err = dbaddiaopt(optlist, dbopt_lo_offset, size(lo_offset), lo_offset)
556 err = dbaddiaopt(optlist, dbopt_hi_offset, size(hi_offset), hi_offset)
557 err = dbputqm(dbfile, 'rectilinear_grid', 16, 'x', 1, 'y', 1, 'z', 1, x_cb, y_cb, db_f77null, dims, 2, db_double, &
558 & db_collinear, optlist, ierr)
559 err = dbfreeoptlist(optlist)
560 else
561 err = dbmkoptlist(2, optlist)
562 err = dbaddiaopt(optlist, dbopt_lo_offset, size(lo_offset), lo_offset)
563 err = dbaddiaopt(optlist, dbopt_hi_offset, size(hi_offset), hi_offset)
564 err = dbputqm(dbfile, 'rectilinear_grid', 16, 'x', 1, 'y', 1, 'z', 1, x_cb, db_f77null, db_f77null, dims, 1, &
565 & db_double, db_collinear, optlist, ierr)
566 err = dbfreeoptlist(optlist)
567 end if
568 else if (format == 2) then
569 ! Multidimensional local grid data is written to the formatted database slave file. Recall that no master file to
570 ! maintained in multidimensions.
571 if (p > 0) then
572 if (precision == 1) then
573 write (dbfile) real(x_cb, sp), real(y_cb, sp), real(z_cb, sp)
574 else
575 if (output_partial_domain) then
576 write (dbfile) x_cb(x_output_idx%beg - 1:x_output_idx%end), y_cb(y_output_idx%beg - 1:y_output_idx%end), &
577 & z_cb(z_output_idx%beg - 1:z_output_idx%end)
578 else
579 write (dbfile) x_cb, y_cb, z_cb
580 end if
581 end if
582 else if (n > 0) then
583 if (precision == 1) then
584 write (dbfile) real(x_cb, sp), real(y_cb, sp)
585 else
586 if (output_partial_domain) then
587 write (dbfile) x_cb(x_output_idx%beg - 1:x_output_idx%end), y_cb(y_output_idx%beg - 1:y_output_idx%end)
588 else
589 write (dbfile) x_cb, y_cb
590 end if
591 end if
592
593 ! One-dimensional local grid data is written to the formatted database slave file. In addition, the local grid data
594 ! is put together by the root process and written to the master file.
595 else
596 if (precision == 1) then
597 write (dbfile) real(x_cb, sp)
598 else
599 if (output_partial_domain) then
600 write (dbfile) x_cb(x_output_idx%beg - 1:x_output_idx%end)
601 else
602 write (dbfile) x_cb
603 end if
604 end if
605
606 if (num_procs > 1) then
608 else
609 x_root_cb(:) = x_cb(:)
610 end if
611
612 if (proc_rank == 0) then
613 if (precision == 1) then
614 write (dbroot) real(x_root_cb, wp)
615 else
616 if (output_partial_domain) then
617 write (dbroot) x_root_cb(x_output_idx%beg - 1:x_output_idx%end)
618 else
619 write (dbroot) x_root_cb
620 end if
621 end if
622 end if
623 end if
624 end if
625
627
628 !> Write a single flow variable field to the formatted database slave and master files for a given time step.
629 impure subroutine s_write_variable_to_formatted_database_file(varname, t_step)
630
631 character(LEN=*), intent(in) :: varname
632 integer, intent(in) :: t_step
633
634 ! NAG compiler requires these to be statically sized
635 character(LEN=4*name_len), dimension(num_procs) :: varnames
636 integer, dimension(num_procs) :: vartypes
637 integer :: i, j, k
638 integer :: ierr
639
640 if (format == 1) then
641 ! Determining the extents of the flow variable on each local process and gathering all this information on root process
642 if (num_procs > 1) then
644 else
645 data_extents(:,0) = (/minval(q_sf), maxval(q_sf)/)
646 end if
647
648 if (proc_rank == 0) then
649 do i = 1, num_procs
650 write (varnames(i), '(A,I0,A,I0,A)') '../p', i - 1, '/', t_step, '.silo:' // trim(varname)
651 end do
652
653 vartypes = db_quadvar
654
655 err = dbset2dstrlen(len(varnames(1)))
656 err = dbmkoptlist(2, optlist)
657 err = dbaddiopt(optlist, dbopt_extents_size, 2)
658 err = dbadddopt(optlist, dbopt_extents, data_extents)
659 err = dbputmvar(dbroot, trim(varname), len_trim(varname), num_procs, varnames, len_trim(varnames), vartypes, &
660 & optlist, ierr)
661 err = dbfreeoptlist(optlist)
662 end if
663
664 if (wp == dp) then
665 if (precision == 1) then
666 do i = -offset_x%beg, m + offset_x%end
667 do j = -offset_y%beg, n + offset_y%end
668 do k = -offset_z%beg, p + offset_z%end
669 q_sf_s(i, j, k) = real(q_sf(i, j, k), sp)
670 end do
671 end do
672 end do
673 if (grid_geometry == 3) then
674 do i = -offset_x%beg, m + offset_x%end
675 do j = -offset_y%beg, n + offset_y%end
676 do k = -offset_z%beg, p + offset_z%end
677 cyl_q_sf_s(j, k, i) = q_sf_s(i, j, k)
678 end do
679 end do
680 end do
681 end if
682 else
683 if (grid_geometry == 3) then
684 do i = -offset_x%beg, m + offset_x%end
685 do j = -offset_y%beg, n + offset_y%end
686 do k = -offset_z%beg, p + offset_z%end
687 cyl_q_sf(j, k, i) = q_sf(i, j, k)
688 end do
689 end do
690 end do
691 end if
692 end if
693 else if (wp == sp) then
694 do i = -offset_x%beg, m + offset_x%end
695 do j = -offset_y%beg, n + offset_y%end
696 do k = -offset_z%beg, p + offset_z%end
697 q_sf_s(i, j, k) = q_sf(i, j, k)
698 end do
699 end do
700 end do
701 if (grid_geometry == 3) then
702 do i = -offset_x%beg, m + offset_x%end
703 do j = -offset_y%beg, n + offset_y%end
704 do k = -offset_z%beg, p + offset_z%end
705 cyl_q_sf_s(j, k, i) = q_sf_s(i, j, k)
706 end do
707 end do
708 end do
709 end if
710 end if
711
712# 662 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
713 if (precision == 1) then
714 if (p > 0) then
715 if (grid_geometry == 3) then
716 err = dbputqv1(dbfile, trim(varname), len_trim(varname), 'rectilinear_grid', 16, cyl_q_sf_s, &
717 & dims - 1, 3, db_f77null, 0, db_float, db_zonecent, db_f77null, ierr)
718 else
719 err = dbputqv1(dbfile, trim(varname), len_trim(varname), 'rectilinear_grid', 16, q_sf_s, &
720 & dims - 1, 3, db_f77null, 0, db_float, db_zonecent, db_f77null, ierr)
721 end if
722 else if (n > 0) then
723 err = dbputqv1(dbfile, trim(varname), len_trim(varname), 'rectilinear_grid', 16, q_sf_s, dims - 1, &
724 & 2, db_f77null, 0, db_float, db_zonecent, db_f77null, ierr)
725 else
726 err = dbputqv1(dbfile, trim(varname), len_trim(varname), 'rectilinear_grid', 16, q_sf_s, dims - 1, &
727 & 1, db_f77null, 0, db_float, db_zonecent, db_f77null, ierr)
728 end if
729 end if
730# 662 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
731 if (precision == 2) then
732 if (p > 0) then
733 if (grid_geometry == 3) then
734 err = dbputqv1(dbfile, trim(varname), len_trim(varname), 'rectilinear_grid', 16, cyl_q_sf, &
735 & dims - 1, 3, db_f77null, 0, db_double, db_zonecent, db_f77null, ierr)
736 else
737 err = dbputqv1(dbfile, trim(varname), len_trim(varname), 'rectilinear_grid', 16, q_sf, &
738 & dims - 1, 3, db_f77null, 0, db_double, db_zonecent, db_f77null, ierr)
739 end if
740 else if (n > 0) then
741 err = dbputqv1(dbfile, trim(varname), len_trim(varname), 'rectilinear_grid', 16, q_sf, dims - 1, &
742 & 2, db_f77null, 0, db_double, db_zonecent, db_f77null, ierr)
743 else
744 err = dbputqv1(dbfile, trim(varname), len_trim(varname), 'rectilinear_grid', 16, q_sf, dims - 1, &
745 & 1, db_f77null, 0, db_double, db_zonecent, db_f77null, ierr)
746 end if
747 end if
748# 680 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
749 else
750 ! Writing the name of the flow variable and its data, associated with the local processor, to the formatted database
751 ! slave file
752 if (precision == 1) then
753 write (dbfile) varname, real(q_sf, wp)
754 else
755 write (dbfile) varname, q_sf
756 end if
757
758 ! In 1D, the root process also takes care of gathering the flow variable data from all of the local processor(s) and
759 ! writes it to the formatted database master file.
760 if (n == 0) then
761 if (num_procs > 1) then
763 else
764 q_root_sf(:,:,:) = q_sf(:,:,:)
765 end if
766
767 if (proc_rank == 0) then
768 if (precision == 1) then
769 write (dbroot) varname, real(q_root_sf, wp)
770 else
771 write (dbroot) varname, q_root_sf
772 end if
773 end if
774 end if
775 end if
776
778
779 !> Write the post-processed results in the folder 'lag_bubbles_data'
780 impure subroutine s_write_lag_bubbles_results_to_text(t_step)
781
782 integer, intent(in) :: t_step
783 character(len=len_trim(case_dir) + 3*name_len) :: file_loc
784 integer :: id
785
786#ifdef MFC_MPI
787 real(wp), dimension(20) :: inputvals
788 real(wp) :: time_real
789 integer, dimension(MPI_STATUS_SIZE) :: status
790 integer(KIND=MPI_OFFSET_KIND) :: disp
791 integer :: view
792 logical :: file_exist
793 integer, dimension(2) :: gsizes, lsizes, start_idx_part
794 integer :: ifile
795 integer :: ierr
796 real(wp) :: file_time, file_dt
797 integer :: file_num_procs, file_tot_part
798 integer :: i
799 integer, dimension(:), allocatable :: proc_bubble_counts
800 real(wp), dimension(1:1,1:lag_io_vars) :: lag_io_null
801
802 lag_io_null = 0._wp
803
804 ! Construct file path
805 write (file_loc, '(A,I0,A)') 'lag_bubbles_', t_step, '.dat'
806 file_loc = trim(case_dir) // '/restart_data' // trim(mpiiofs) // trim(file_loc)
807
808 ! Check if file exists
809 inquire (file=trim(file_loc), exist=file_exist)
810 if (.not. file_exist) then
811 call s_mpi_abort('Restart file ' // trim(file_loc) // ' does not exist!')
812 end if
813
814 if (.not. parallel_io) return
815
816 if (proc_rank == 0) then
817 call mpi_file_open(mpi_comm_self, file_loc, mpi_mode_rdonly, mpi_info_int, ifile, ierr)
818
819 call mpi_file_read(ifile, file_tot_part, 1, mpi_integer, status, ierr)
820 call mpi_file_read(ifile, file_time, 1, mpi_p, status, ierr)
821 call mpi_file_read(ifile, file_dt, 1, mpi_p, status, ierr)
822 call mpi_file_read(ifile, file_num_procs, 1, mpi_integer, status, ierr)
823
824 call mpi_file_close(ifile, ierr)
825 end if
826
827 call mpi_bcast(file_tot_part, 1, mpi_integer, 0, mpi_comm_world, ierr)
828 call mpi_bcast(file_time, 1, mpi_p, 0, mpi_comm_world, ierr)
829 call mpi_bcast(file_dt, 1, mpi_p, 0, mpi_comm_world, ierr)
830 call mpi_bcast(file_num_procs, 1, mpi_integer, 0, mpi_comm_world, ierr)
831 time_real = file_time
832
833 allocate (proc_bubble_counts(file_num_procs))
834
835 if (proc_rank == 0) then
836 call mpi_file_open(mpi_comm_self, file_loc, mpi_mode_rdonly, mpi_info_int, ifile, ierr)
837
838 ! Skip to processor counts position
839 disp = int(sizeof(file_tot_part) + 2*sizeof(file_time) + sizeof(file_num_procs), mpi_offset_kind)
840 call mpi_file_seek(ifile, disp, mpi_seek_set, ierr)
841 call mpi_file_read(ifile, proc_bubble_counts, file_num_procs, mpi_integer, status, ierr)
842
843 call mpi_file_close(ifile, ierr)
844 end if
845
846 call mpi_bcast(proc_bubble_counts, file_num_procs, mpi_integer, 0, mpi_comm_world, ierr)
847
848 gsizes(1) = file_tot_part
849 gsizes(2) = lag_io_vars
850 lsizes(1) = file_tot_part
851 lsizes(2) = lag_io_vars
852 start_idx_part(1) = 0
853 start_idx_part(2) = 0
854
855 call mpi_type_create_subarray(2, gsizes, lsizes, start_idx_part, mpi_order_fortran, mpi_p, view, ierr)
856 call mpi_type_commit(view, ierr)
857
858 call mpi_file_open(mpi_comm_world, file_loc, mpi_mode_rdonly, mpi_info_int, ifile, ierr)
859
860 disp = int(sizeof(file_tot_part) + 2*sizeof(file_time) + sizeof(file_num_procs) &
861 & + file_num_procs*sizeof(proc_bubble_counts(1)), mpi_offset_kind)
862 call mpi_file_set_view(ifile, disp, mpi_p, view, 'native', mpi_info_null, ierr)
863
864 allocate (mpi_io_data_lg_bubbles(file_tot_part,1:lag_io_vars))
865
866 call mpi_file_read_all(ifile, mpi_io_data_lg_bubbles, lag_io_vars*file_tot_part, mpi_p, status, ierr)
867
868 write (file_loc, '(A,I0,A)') 'lag_bubbles_post_process_', t_step, '.dat'
869 file_loc = trim(case_dir) // '/lag_bubbles_post_process/' // trim(file_loc)
870
871 if (proc_rank == 0) then
872 open (unit=29, file=file_loc, form='formatted', position='rewind')
873
874 if (lag_header) then
875 write (29, '(A)', advance='no')
876 if (lag_id_wrt) write (29, '(A8)', advance='no') 'id, '
877 if (lag_pos_wrt) write (29, '(3(A17))', advance='no') 'px, ', 'py, ', 'pz, '
878 if (lag_pos_prev_wrt) write (29, '(3(A17))', advance='no') 'pvx, ', 'pvy, ', 'pvz, '
879 if (lag_vel_wrt) write (29, '(3(A17))', advance='no') 'vx, ', 'vy, ', 'vz, '
880 if (lag_rad_wrt) write (29, '(A17)', advance='no') 'radius, '
881 if (lag_rvel_wrt) write (29, '(A17)', advance='no') 'rvel, '
882 if (lag_r0_wrt) write (29, '(A17)', advance='no') 'r0, '
883 if (lag_rmax_wrt) write (29, '(A17)', advance='no') 'rmax, '
884 if (lag_rmin_wrt) write (29, '(A17)', advance='no') 'rmin, '
885 if (lag_dphidt_wrt) write (29, '(A17)', advance='no') 'dphidt, '
886 if (lag_pres_wrt) write (29, '(A17)', advance='no') 'pressure, '
887 if (lag_mv_wrt) write (29, '(A17)', advance='no') 'mv, '
888 if (lag_mg_wrt) write (29, '(A17)', advance='no') 'mg, '
889 if (lag_betat_wrt) write (29, '(A17)', advance='no') 'betaT, '
890 if (lag_betac_wrt) write (29, '(A17)', advance='no') 'betaC, '
891 write (29, '(A15)') 'time'
892 end if
893
894 do i = 1, file_tot_part
895 id = int(mpi_io_data_lg_bubbles(i, 1))
896 inputvals(1:20) = mpi_io_data_lg_bubbles(i,2:21)
897 if (id > 0) then
898 write (29, '(100(A))', advance='no') ''
899 if (lag_id_wrt) write (29, '(I6, A)', advance='no') id, ', '
900 if (lag_pos_wrt) write (29, '(3(E15.7, A))', advance='no') inputvals(1), ', ', inputvals(2), ', ', &
901 & inputvals(3), ', '
902 if (lag_pos_prev_wrt) write (29, '(3(E15.7, A))', advance='no') inputvals(4), ', ', inputvals(5), ', ', &
903 & inputvals(6), ', '
904 if (lag_vel_wrt) write (29, '(3(E15.7, A))', advance='no') inputvals(7), ', ', inputvals(8), ', ', &
905 & inputvals(9), ', '
906 if (lag_rad_wrt) write (29, '(E15.7, A)', advance='no') inputvals(10), ', '
907 if (lag_rvel_wrt) write (29, '(E15.7, A)', advance='no') inputvals(11), ', '
908 if (lag_r0_wrt) write (29, '(E15.7, A)', advance='no') inputvals(12), ', '
909 if (lag_rmax_wrt) write (29, '(E15.7, A)', advance='no') inputvals(13), ', '
910 if (lag_rmin_wrt) write (29, '(E15.7, A)', advance='no') inputvals(14), ', '
911 if (lag_dphidt_wrt) write (29, '(E15.7, A)', advance='no') inputvals(15), ', '
912 if (lag_pres_wrt) write (29, '(E15.7, A)', advance='no') inputvals(16), ', '
913 if (lag_mv_wrt) write (29, '(E15.7, A)', advance='no') inputvals(17), ', '
914 if (lag_mg_wrt) write (29, '(E15.7, A)', advance='no') inputvals(18), ', '
915 if (lag_betat_wrt) write (29, '(E15.7, A)', advance='no') inputvals(19), ', '
916 if (lag_betac_wrt) write (29, '(E15.7, A)', advance='no') inputvals(20), ', '
917 write (29, '(E15.7)') time_real
918 end if
919 end do
920 close (29)
921 end if
922
923 deallocate (mpi_io_data_lg_bubbles)
924
925 call s_mpi_barrier()
926
927 call mpi_file_close(ifile, ierr)
928#endif
929
931
932 !> Read Lagrangian bubble restart data and write bubble positions and scalar fields to the Silo database.
934
935 integer, intent(in) :: t_step
936 character(len=len_trim(case_dir) + 3*name_len) :: file_loc
937 integer :: id
938
939#ifdef MFC_MPI
940 real(wp) :: time_real
941 integer, dimension(MPI_STATUS_SIZE) :: status
942 integer(KIND=MPI_OFFSET_KIND) :: disp
943 integer :: view
944 logical :: file_exist
945 integer, dimension(2) :: gsizes, lsizes, start_idx_part
946 integer :: ifile, ierr, nbub
947 real(wp) :: file_time, file_dt
948 integer :: file_num_procs, file_tot_part
949 integer, dimension(:), allocatable :: proc_bubble_counts
950 real(wp), dimension(1:1,1:lag_io_vars) :: dummy
951 character(LEN=4*name_len), dimension(num_procs) :: meshnames
952 integer, dimension(num_procs) :: meshtypes
953 real(wp) :: dummy_data
954 integer :: i
955 real(wp), dimension(:), allocatable :: bub_id
956 real(wp), dimension(:), allocatable :: px, py, pz, ppx, ppy, ppz, vx, vy, vz
957 real(wp), dimension(:), allocatable :: radius, rvel, rnot, rmax, rmin, dphidt
958 real(wp), dimension(:), allocatable :: pressure, mv, mg, betat, betac
959
960 dummy = 0._wp
961 dummy_data = 0._wp
962
963 ! Construct file path
964 write (file_loc, '(A,I0,A)') 'lag_bubbles_', t_step, '.dat'
965 file_loc = trim(case_dir) // '/restart_data' // trim(mpiiofs) // trim(file_loc)
966
967 ! Check if file exists
968 inquire (file=trim(file_loc), exist=file_exist)
969 if (.not. file_exist) then
970 call s_mpi_abort('Restart file ' // trim(file_loc) // ' does not exist!')
971 end if
972
973 if (.not. parallel_io) return
974
975 if (proc_rank == 0) then
976 call mpi_file_open(mpi_comm_self, file_loc, mpi_mode_rdonly, mpi_info_int, ifile, ierr)
977
978 call mpi_file_read(ifile, file_tot_part, 1, mpi_integer, status, ierr)
979 call mpi_file_read(ifile, file_time, 1, mpi_p, status, ierr)
980 call mpi_file_read(ifile, file_dt, 1, mpi_p, status, ierr)
981 call mpi_file_read(ifile, file_num_procs, 1, mpi_integer, status, ierr)
982
983 call mpi_file_close(ifile, ierr)
984 end if
985
986 call mpi_bcast(file_tot_part, 1, mpi_integer, 0, mpi_comm_world, ierr)
987 call mpi_bcast(file_time, 1, mpi_p, 0, mpi_comm_world, ierr)
988 call mpi_bcast(file_dt, 1, mpi_p, 0, mpi_comm_world, ierr)
989 call mpi_bcast(file_num_procs, 1, mpi_integer, 0, mpi_comm_world, ierr)
990 time_real = file_time
991
992 allocate (proc_bubble_counts(file_num_procs))
993
994 if (proc_rank == 0) then
995 call mpi_file_open(mpi_comm_self, file_loc, mpi_mode_rdonly, mpi_info_int, ifile, ierr)
996
997 ! Skip to processor counts position
998 disp = int(sizeof(file_tot_part) + 2*sizeof(file_time) + sizeof(file_num_procs), mpi_offset_kind)
999 call mpi_file_seek(ifile, disp, mpi_seek_set, ierr)
1000 call mpi_file_read(ifile, proc_bubble_counts, file_num_procs, mpi_integer, status, ierr)
1001
1002 call mpi_file_close(ifile, ierr)
1003 end if
1004
1005 call mpi_bcast(proc_bubble_counts, file_num_procs, mpi_integer, 0, mpi_comm_world, ierr)
1006
1007 ! Set time variables from file
1008
1009 nbub = proc_bubble_counts(proc_rank + 1)
1010
1011 start_idx_part(1) = 0
1012 do i = 1, proc_rank
1013 start_idx_part(1) = start_idx_part(1) + proc_bubble_counts(i)
1014 end do
1015
1016 start_idx_part(2) = 0
1017 lsizes(1) = nbub
1018 lsizes(2) = lag_io_vars
1019
1020 gsizes(1) = file_tot_part
1021 gsizes(2) = lag_io_vars
1022
1023 if (nbub > 0) then
1024# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1025 allocate (bub_id(nbub))
1026# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1027 allocate (px(nbub))
1028# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1029 allocate (py(nbub))
1030# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1031 allocate (pz(nbub))
1032# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1033 allocate (ppx(nbub))
1034# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1035 allocate (ppy(nbub))
1036# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1037 allocate (ppz(nbub))
1038# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1039 allocate (vx(nbub))
1040# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1041 allocate (vy(nbub))
1042# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1043 allocate (vz(nbub))
1044# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1045 allocate (radius(nbub))
1046# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1047 allocate (rvel(nbub))
1048# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1049 allocate (rnot(nbub))
1050# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1051 allocate (rmax(nbub))
1052# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1053 allocate (rmin(nbub))
1054# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1055 allocate (dphidt(nbub))
1056# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1057 allocate (pressure(nbub))
1058# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1059 allocate (mv(nbub))
1060# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1061 allocate (mg(nbub))
1062# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1063 allocate (betat(nbub))
1064# 958 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1065 allocate (betac(nbub))
1066# 960 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1067 allocate (mpi_io_data_lg_bubbles(nbub,1:lag_io_vars))
1068
1069 call mpi_type_create_subarray(2, gsizes, lsizes, start_idx_part, mpi_order_fortran, mpi_p, view, ierr)
1070 call mpi_type_commit(view, ierr)
1071
1072 call mpi_file_open(mpi_comm_world, file_loc, mpi_mode_rdonly, mpi_info_int, ifile, ierr)
1073
1074 ! Skip extended header
1075 disp = int(sizeof(file_tot_part) + 2*sizeof(file_time) + sizeof(file_num_procs) &
1076 & + file_num_procs*sizeof(proc_bubble_counts(1)), mpi_offset_kind)
1077 call mpi_file_set_view(ifile, disp, mpi_p, view, 'native', mpi_info_int, ierr)
1078
1079 call mpi_file_read_all(ifile, mpi_io_data_lg_bubbles, lag_io_vars*nbub, mpi_p, status, ierr)
1080
1081 call mpi_file_close(ifile, ierr)
1082 call mpi_type_free(view, ierr)
1083
1084 ! Extract data from MPI_IO_DATA_lg_bubbles array Adjust these indices based on your actual data layout
1085# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1086 bub_id(:) = mpi_io_data_lg_bubbles(:,1)
1087# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1088 px(:) = mpi_io_data_lg_bubbles(:,2)
1089# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1090 py(:) = mpi_io_data_lg_bubbles(:,3)
1091# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1092 pz(:) = mpi_io_data_lg_bubbles(:,4)
1093# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1094 ppx(:) = mpi_io_data_lg_bubbles(:,5)
1095# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1096 ppy(:) = mpi_io_data_lg_bubbles(:,6)
1097# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1098 ppz(:) = mpi_io_data_lg_bubbles(:,7)
1099# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1100 vx(:) = mpi_io_data_lg_bubbles(:,8)
1101# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1102 vy(:) = mpi_io_data_lg_bubbles(:,9)
1103# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1104 vz(:) = mpi_io_data_lg_bubbles(:,10)
1105# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1106 radius(:) = mpi_io_data_lg_bubbles(:,11)
1107# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1108 rvel(:) = mpi_io_data_lg_bubbles(:,12)
1109# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1110 rnot(:) = mpi_io_data_lg_bubbles(:,13)
1111# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1112 rmax(:) = mpi_io_data_lg_bubbles(:,14)
1113# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1114 rmin(:) = mpi_io_data_lg_bubbles(:,15)
1115# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1116 dphidt(:) = mpi_io_data_lg_bubbles(:,16)
1117# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1118 pressure(:) = mpi_io_data_lg_bubbles(:,17)
1119# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1120 mv(:) = mpi_io_data_lg_bubbles(:,18)
1121# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1122 mg(:) = mpi_io_data_lg_bubbles(:,19)
1123# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1124 betat(:) = mpi_io_data_lg_bubbles(:,20)
1125# 982 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1126 betac(:) = mpi_io_data_lg_bubbles(:,21)
1127# 984 "/home/runner/work/MFC/MFC/src/post_process/m_data_output.fpp"
1128
1129 ! Next, the root processor proceeds to record all of the spatial extents in the formatted database master file. In
1130 ! addition, it also records a sub-domain connectivity map so that the entire grid may be reassembled by looking at the
1131 ! master file.
1132 if (proc_rank == 0) then
1133 do i = 1, num_procs
1134 write (meshnames(i), '(A,I0,A,I0,A)') '../p', i - 1, '/', t_step, '.silo:lag_bubbles'
1135 meshtypes(i) = db_pointmesh
1136 end do
1137 err = dbset2dstrlen(len(meshnames(1)))
1138 err = dbputmmesh(dbroot, 'lag_bubbles', 16, num_procs, meshnames, len_trim(meshnames), meshtypes, db_f77null, ierr)
1139 end if
1140
1141 err = dbputpm(dbfile, 'lag_bubbles', 11, 3, px, py, pz, nbub, db_double, db_f77null, ierr)
1142
1143 if (lag_id_wrt) call s_write_lag_variable_to_formatted_database_file('part_id', t_step, bub_id, nbub)
1144 if (lag_vel_wrt) then
1145 call s_write_lag_variable_to_formatted_database_file('part_vel1', t_step, vx, nbub)
1146 call s_write_lag_variable_to_formatted_database_file('part_vel2', t_step, vy, nbub)
1147 if (p > 0) call s_write_lag_variable_to_formatted_database_file('part_vel3', t_step, vz, nbub)
1148 end if
1149 if (lag_rad_wrt) call s_write_lag_variable_to_formatted_database_file('part_radius', t_step, radius, nbub)
1150 if (lag_rvel_wrt) call s_write_lag_variable_to_formatted_database_file('part_rdot', t_step, rvel, nbub)
1151 if (lag_r0_wrt) call s_write_lag_variable_to_formatted_database_file('part_r0', t_step, rnot, nbub)
1152 if (lag_rmax_wrt) call s_write_lag_variable_to_formatted_database_file('part_rmax', t_step, rmax, nbub)
1153 if (lag_rmin_wrt) call s_write_lag_variable_to_formatted_database_file('part_rmin', t_step, rmin, nbub)
1154 if (lag_dphidt_wrt) call s_write_lag_variable_to_formatted_database_file('part_dphidt', t_step, dphidt, nbub)
1155 if (lag_pres_wrt) call s_write_lag_variable_to_formatted_database_file('part_pressure', t_step, pressure, nbub)
1156 if (lag_mv_wrt) call s_write_lag_variable_to_formatted_database_file('part_mv', t_step, mv, nbub)
1157 if (lag_mg_wrt) call s_write_lag_variable_to_formatted_database_file('part_mg', t_step, mg, nbub)
1158 if (lag_betat_wrt) call s_write_lag_variable_to_formatted_database_file('part_betaT', t_step, betat, nbub)
1159 if (lag_betac_wrt) call s_write_lag_variable_to_formatted_database_file('part_betaC', t_step, betac, nbub)
1160
1161 deallocate (bub_id, px, py, pz, ppx, ppy, ppz, vx, vy, vz, radius, rvel, rnot, rmax, rmin, dphidt, pressure, mv, mg, &
1162 & betat, betac)
1163 deallocate (mpi_io_data_lg_bubbles)
1164 else
1165 call mpi_type_contiguous(0, mpi_p, view, ierr)
1166 call mpi_type_commit(view, ierr)
1167
1168 call mpi_file_open(mpi_comm_world, file_loc, mpi_mode_rdonly, mpi_info_int, ifile, ierr)
1169
1170 ! Skip extended header
1171 disp = int(sizeof(file_tot_part) + 2*sizeof(file_time) + sizeof(file_num_procs) &
1172 & + file_num_procs*sizeof(proc_bubble_counts(1)), mpi_offset_kind)
1173 call mpi_file_set_view(ifile, disp, mpi_p, view, 'native', mpi_info_int, ierr)
1174
1175 call mpi_file_read_all(ifile, dummy, 0, mpi_p, status, ierr)
1176
1177 call mpi_file_close(ifile, ierr)
1178 call mpi_type_free(view, ierr)
1179
1180 if (proc_rank == 0) then
1181 do i = 1, num_procs
1182 write (meshnames(i), '(A,I0,A,I0,A)') '../p', i - 1, '/', t_step, '.silo:lag_bubbles'
1183 meshtypes(i) = db_pointmesh
1184 end do
1185 err = dbset2dstrlen(len(meshnames(1)))
1186 err = dbputmmesh(dbroot, 'lag_bubbles', 16, num_procs, meshnames, len_trim(meshnames), meshtypes, db_f77null, ierr)
1187 end if
1188
1189 err = dbsetemptyok(1)
1190 err = dbputpm(dbfile, 'lag_bubbles', 11, 3, dummy_data, dummy_data, dummy_data, 0, db_double, db_f77null, ierr)
1191
1193 if (lag_vel_wrt) then
1194 call s_write_lag_variable_to_formatted_database_file('part_vel1', t_step)
1195 call s_write_lag_variable_to_formatted_database_file('part_vel2', t_step)
1196 if (p > 0) call s_write_lag_variable_to_formatted_database_file('part_vel3', t_step)
1197 end if
1198 if (lag_rad_wrt) call s_write_lag_variable_to_formatted_database_file('part_radius', t_step)
1204 if (lag_pres_wrt) call s_write_lag_variable_to_formatted_database_file('part_pressure', t_step)
1209 end if
1210#endif
1211
1213
1214 !> Write a single Lagrangian bubble point-variable to the Silo database slave and master files.
1215 subroutine s_write_lag_variable_to_formatted_database_file(varname, t_step, data, nBubs)
1216
1217 character(len=*), intent(in) :: varname
1218 integer, intent(in) :: t_step
1219 real(wp), dimension(1:), intent(in), optional :: data
1220 integer, intent(in), optional :: nBubs
1221 character(len=64), dimension(num_procs) :: var_names
1222 integer, dimension(num_procs) :: var_types
1223 real(wp) :: dummy_data
1224 integer :: ierr
1225 integer :: i
1226
1227 dummy_data = 0._wp
1228
1229 if (present(nbubs) .and. present(data)) then
1230 if (proc_rank == 0) then
1231 do i = 1, num_procs
1232 write (var_names(i), '(A,I0,A,I0,A)') '../p', i - 1, '/', t_step, '.silo:' // trim(varname)
1233 var_types(i) = db_pointvar
1234 end do
1235 err = dbset2dstrlen(len(var_names(1)))
1236 err = dbputmvar(dbroot, trim(varname), len_trim(varname), num_procs, var_names, len_trim(var_names), var_types, &
1237 & db_f77null, ierr)
1238 end if
1239
1240 err = dbputpv1(dbfile, trim(varname), len_trim(varname), 'lag_bubbles', 11, data, nbubs, db_double, db_f77null, ierr)
1241 else
1242 if (proc_rank == 0) then
1243 do i = 1, num_procs
1244 write (var_names(i), '(A,I0,A,I0,A)') '../p', i - 1, '/', t_step, '.silo:' // trim(varname)
1245 var_types(i) = db_pointvar
1246 end do
1247 err = dbset2dstrlen(len(var_names(1)))
1248 err = dbsetemptyok(1)
1249 err = dbputmvar(dbroot, trim(varname), len_trim(varname), num_procs, var_names, len_trim(var_names), var_types, &
1250 & db_f77null, ierr)
1251 end if
1252
1253 err = dbsetemptyok(1)
1254 err = dbputpv1(dbfile, trim(varname), len_trim(varname), 'lag_bubbles', 11, dummy_data, 0, db_double, db_f77null, ierr)
1255 end if
1256
1258
1259 !> Convert the binary immersed-boundary state file to per-body formatted text files
1260 impure subroutine s_write_ib_state_files()
1261
1262 character(len=len_trim(case_dir) + 4*name_len) :: in_file, out_file, file_loc
1263 integer :: iu_in, ios, i, rec_id
1264 integer, allocatable, dimension(:) :: iu_out
1265 real(wp) :: rec_time
1266 real(wp), dimension(3) :: rec_force, rec_torque
1267 real(wp), dimension(3) :: rec_vel, rec_angular_vel
1268 real(wp), dimension(3) :: rec_angles, rec_centroid
1269
1270 file_loc = trim(case_dir) // '/D'
1271
1272 in_file = trim(file_loc) // '/ib_state.dat'
1273 open (newunit=iu_in, file=trim(in_file), form='unformatted', access='stream', status='old', action='read', iostat=ios)
1274 if (ios /= 0) then
1275 call s_mpi_abort('Cannot open IB state input file: ' // trim(in_file))
1276 end if
1277
1278 allocate (iu_out(num_ibs))
1279 do i = 1, num_ibs
1280 write (out_file, '(A,I0,A)') trim(file_loc) // '/ib_', i, '.txt'
1281 open (newunit=iu_out(i), file=trim(out_file), form='formatted', status='replace', action='write', iostat=ios)
1282 if (ios /= 0) then
1283 call s_mpi_abort('Cannot open IB state output file: ' // trim(out_file))
1284 end if
1285 write (iu_out(i), &
1286 & '(A)') 'mytime fx fy fz Tau_x Tau_y Tau_z vx vy vz omega_x omega_y omega_z angle_x angle_y angle_z x_c y_c z_c'
1287 end do
1288
1289 do
1290 read (iu_in, iostat=ios) rec_time, rec_id, rec_force, rec_torque, rec_vel, rec_angular_vel, rec_angles, &
1291 & rec_centroid(1), rec_centroid(2), rec_centroid(3)
1292 if (ios /= 0) exit
1293
1294 if (rec_id >= 1 .and. rec_id <= num_ibs) then
1295 write (iu_out(rec_id), '(19(ES24.16E3,1X))') rec_time, rec_force(1), rec_force(2), rec_force(3), rec_torque(1), &
1296 & rec_torque(2), rec_torque(3), rec_vel(1), rec_vel(2), rec_vel(3), rec_angular_vel(1), &
1297 & rec_angular_vel(2), rec_angular_vel(3), rec_angles(1), rec_angles(2), rec_angles(3), rec_centroid(1), &
1298 & rec_centroid(2), rec_centroid(3)
1299 end if
1300 end do
1301
1302 close (iu_in)
1303 do i = 1, num_ibs
1304 close (iu_out(i))
1305 end do
1306 deallocate (iu_out)
1307
1308 end subroutine s_write_ib_state_files
1309
1310 !> Extract the volume-fraction interface contour from primitive fields and write the coordinates to the interface data file.
1311 impure subroutine s_write_intf_data_file(q_prim_vf)
1312
1313 type(scalar_field), dimension(sys_size), intent(in) :: q_prim_vf
1314 integer :: i, j, k, l, cent
1315 integer :: counter, root !< number of data points extracted to fit shape to SH perturbations
1316 real(wp), allocatable :: x_td(:), y_td(:), x_d1(:), y_d1(:), y_d(:), x_d(:)
1317 real(wp) :: axp, axm, ayp, aym, tgp, euc_d, thres, maxalph_loc, maxalph_glb
1318
1319 allocate (x_d1(m*n))
1320 allocate (y_d1(m*n))
1321 counter = 0
1322 maxalph_loc = 0._wp
1323 do k = 0, p
1324 do j = 0, n
1325 do i = 0, m
1326 if (q_prim_vf(eqn_idx%E + 2)%sf(i, j, k) > maxalph_loc) then
1327 maxalph_loc = q_prim_vf(eqn_idx%E + 2)%sf(i, j, k)
1328 end if
1329 end do
1330 end do
1331 end do
1332
1333 call s_mpi_allreduce_max(maxalph_loc, maxalph_glb)
1334 if (p > 0) then
1335 do l = 0, p
1336 if (z_cc(l) < dz(l) .and. z_cc(l) > 0) then
1337 cent = l
1338 end if
1339 end do
1340 else
1341 cent = 0
1342 end if
1343
1344 thres = 0.9_wp*maxalph_glb
1345 do k = 0, n
1346 do j = 0, m
1347 axp = q_prim_vf(eqn_idx%E + 2)%sf(j + 1, k, cent)
1348 axm = q_prim_vf(eqn_idx%E + 2)%sf(j, k, cent)
1349 ayp = q_prim_vf(eqn_idx%E + 2)%sf(j, k + 1, cent)
1350 aym = q_prim_vf(eqn_idx%E + 2)%sf(j, k, cent)
1351 if ((axp > thres .and. axm < thres) .or. (axp < thres .and. axm > thres) .or. (ayp > thres .and. aym < thres) &
1352 & .or. (ayp < thres .and. aym > thres)) then
1353 if (counter == 0) then
1354 counter = counter + 1
1355 x_d1(counter) = x_cc(j)
1356 y_d1(counter) = y_cc(k)
1357 else
1358 tgp = sqrt(dx(j)**2 + dy(k)**2)
1359 do i = 1, counter
1360 euc_d = sqrt((x_cc(j) - x_d1(i))**2 + (y_cc(k) - y_d1(i))**2)
1361 if (euc_d < tgp) then
1362 exit
1363 else if (i == counter) then
1364 counter = counter + 1
1365 x_d1(counter) = x_cc(j)
1366 y_d1(counter) = y_cc(k)
1367 end if
1368 end do
1369 end if
1370 end if
1371 end do
1372 end do
1373
1374 allocate (x_d(counter), y_d(counter))
1375
1376 do i = 1, counter
1377 y_d(i) = y_d1(i)
1378 x_d(i) = x_d1(i)
1379 end do
1380 root = 0
1381
1382 call s_mpi_gather_data(x_d, counter, x_td, root)
1383 call s_mpi_gather_data(y_d, counter, y_td, root)
1384 if (proc_rank == 0) then
1385 do i = 1, size(x_td)
1386 if (i == size(x_td)) then
1387 write (211, '(F12.9,1X,F12.9,1X,I4)') x_td(i), y_td(i), size(x_td)
1388 else
1389 write (211, '(F12.9,1X,F12.9,1X,F3.1)') x_td(i), y_td(i), 0._wp
1390 end if
1391 end do
1392 end if
1393
1394 end subroutine s_write_intf_data_file
1395
1396 !> Compute volume-integrated kinetic, potential, and internal energies and write the energy budget to the energy data file.
1397 impure subroutine s_write_energy_data_file(q_prim_vf, q_cons_vf)
1398
1399 type(scalar_field), dimension(sys_size), intent(in) :: q_prim_vf, q_cons_vf
1400 real(wp) :: elk, egk, elp, egint, vb, vl, pres_av, et
1401 real(wp) :: rho, pres, dv, tmp, gamma, pi_inf, maxma, maxma_glb, maxvel, c, ma, h, qv
1402 real(wp), dimension(num_vels) :: vel
1403 real(wp), dimension(num_fluids) :: adv
1404 integer :: i, j, k, l, s !< looping indices
1405
1406 egk = 0._wp
1407 elp = 0._wp
1408 egint = 0._wp
1409 vb = 0._wp
1410 maxvel = 0._wp
1411 maxma = 0._wp
1412 vl = 0._wp
1413 elk = 0._wp
1414 et = 0._wp
1415 vb = 0._wp
1416 dv = 0._wp
1417 pres_av = 0._wp
1418 pres = 0._wp
1419 c = 0._wp
1420
1421 do k = 0, p
1422 do j = 0, n
1423 do i = 0, m
1424 pres = 0._wp
1425 dv = dx(i)*dy(j)*dz(k)
1426 rho = 0._wp
1427 gamma = 0._wp
1428 pi_inf = 0._wp
1429 qv = 0._wp
1430 pres = q_prim_vf(eqn_idx%E)%sf(i, j, k)
1431 egint = egint + q_prim_vf(eqn_idx%E + 2)%sf(i, j, k)*(gammas(2)*pres)*dv
1432 do s = 1, num_vels
1433 vel(s) = q_prim_vf(num_fluids + s)%sf(i, j, k)
1434 egk = egk + 0.5_wp*q_prim_vf(eqn_idx%E + 2)%sf(i, j, k)*q_prim_vf(2)%sf(i, j, k)*vel(s)*vel(s)*dv
1435 elk = elk + 0.5_wp*q_prim_vf(eqn_idx%E + 1)%sf(i, j, k)*q_prim_vf(1)%sf(i, j, k)*vel(s)*vel(s)*dv
1436 if (abs(vel(s)) > maxvel) then
1437 maxvel = abs(vel(s))
1438 end if
1439 end do
1440 do l = 1, eqn_idx%adv%end - eqn_idx%E
1441 adv(l) = q_prim_vf(eqn_idx%E + l)%sf(i, j, k)
1442 gamma = gamma + adv(l)*gammas(l)
1443 pi_inf = pi_inf + adv(l)*pi_infs(l)
1444 rho = rho + adv(l)*q_prim_vf(l)%sf(i, j, k)
1445 qv = qv + adv(l)*q_prim_vf(l)%sf(i, j, k)*qvs(l)
1446 end do
1447
1448 h = ((gamma + 1._wp)*pres + pi_inf + qv)/rho
1449
1450 call s_compute_speed_of_sound(pres, rho, gamma, pi_inf, h, adv, 0._wp, 0._wp, c, qv)
1451
1452 ma = maxvel/c
1453 if (ma > maxma .and. (adv(1) > (1.0_wp - 1.0e-10_wp))) then
1454 maxma = ma
1455 end if
1456 vl = vl + adv(1)*dv
1457 vb = vb + adv(2)*dv
1458 pres_av = pres_av + adv(1)*pres*dv
1459 et = et + q_cons_vf(eqn_idx%E)%sf(i, j, k)*dv
1460 end do
1461 end do
1462 end do
1463
1464 tmp = pres_av
1465 call s_mpi_allreduce_sum(tmp, pres_av)
1466 tmp = vl
1467 call s_mpi_allreduce_sum(tmp, vl)
1468
1469 call s_mpi_allreduce_max(maxma, maxma_glb)
1470 tmp = elk
1471 call s_mpi_allreduce_sum(tmp, elk)
1472 tmp = egint
1473 call s_mpi_allreduce_sum(tmp, egint)
1474 tmp = egk
1475 call s_mpi_allreduce_sum(tmp, egk)
1476 tmp = vb
1477 call s_mpi_allreduce_sum(tmp, vb)
1478 tmp = et
1479 call s_mpi_allreduce_sum(tmp, et)
1480
1481 elp = pres_av/vl*vb
1482 if (proc_rank == 0) then
1483 write (251, '(10X, 8F24.8)') elp, egint, elk, egk, et, vb, vl, maxma_glb
1484 end if
1485
1486 end subroutine s_write_energy_data_file
1487
1488 !> Read IB state and write a Silo point mesh with per-body scalar fields.
1490
1491 integer, intent(in) :: t_step
1492 character(len=len_trim(case_dir) + 3*name_len) :: file_loc
1493
1494#ifdef MFC_MPI
1495 integer, parameter :: nfields_per_ib = 20
1496 real(wp) :: ib_buf(nfields_per_ib)
1497 real(wp), dimension(:,:), allocatable :: ib_data
1498 logical :: file_exist
1499 character(LEN=4*name_len), dimension(num_procs) :: meshnames
1500 integer, dimension(num_procs) :: meshtypes
1501 integer :: i, ios, file_unit
1502 integer :: ierr, nbodies
1503 real(wp), dimension(:), allocatable :: px, py, pz
1504 real(wp), dimension(:), allocatable :: force_x, force_y, force_z
1505 real(wp), dimension(:), allocatable :: torque_x, torque_y, torque_z
1506 real(wp), dimension(:), allocatable :: vel_x, vel_y, vel_z
1507 real(wp), dimension(:), allocatable :: omega_x, omega_y, omega_z
1508 real(wp), dimension(:), allocatable :: angle_x, angle_y, angle_z
1509 real(wp), dimension(:), allocatable :: ib_diameter
1510
1511 ! Build path to per-timestep IB state file
1512 write (file_loc, '(A,I0,A)') '/restart_data/ib_state_', t_step, '.dat'
1513 file_loc = trim(case_dir) // trim(file_loc)
1514
1515 inquire (file=trim(file_loc), exist=file_exist)
1516 if (.not. file_exist) then
1517 call s_mpi_abort('Restart file ' // trim(file_loc) // ' does not exist!')
1518 end if
1519
1520 nbodies = num_ibs
1521
1522 if (nbodies > 0) then
1523 allocate (ib_data(nbodies, nfields_per_ib))
1524 allocate (px(nbodies), py(nbodies), pz(nbodies))
1525 allocate (force_x(nbodies), force_y(nbodies), force_z(nbodies))
1526 allocate (torque_x(nbodies), torque_y(nbodies), torque_z(nbodies))
1527 allocate (vel_x(nbodies), vel_y(nbodies), vel_z(nbodies))
1528 allocate (omega_x(nbodies), omega_y(nbodies), omega_z(nbodies))
1529 allocate (angle_x(nbodies), angle_y(nbodies), angle_z(nbodies))
1530 allocate (ib_diameter(nbodies))
1531
1532 if (proc_rank == 0) then
1533 open (newunit=file_unit, file=trim(file_loc), form='unformatted', access='stream', status='old', iostat=ios)
1534 if (ios /= 0) call s_mpi_abort('Cannot open IB state file: ' // trim(file_loc))
1535
1536 do i = 1, nbodies
1537 read (file_unit, iostat=ios) ib_buf
1538 if (ios /= 0) call s_mpi_abort('Error reading IB state file')
1539 ib_data(i,:) = ib_buf(:)
1540 end do
1541
1542 close (file_unit)
1543 end if
1544
1545 call mpi_bcast(ib_data, nbodies*nfields_per_ib, mpi_p, 0, mpi_comm_world, ierr)
1546
1547 do i = 1, nbodies
1548 force_x(i) = ib_data(i, 2); force_y(i) = ib_data(i, 3); force_z(i) = ib_data(i, 4)
1549 torque_x(i) = ib_data(i, 5); torque_y(i) = ib_data(i, 6); torque_z(i) = ib_data(i, 7)
1550 vel_x(i) = ib_data(i, 8); vel_y(i) = ib_data(i, 9); vel_z(i) = ib_data(i, 10)
1551 omega_x(i) = ib_data(i, 11); omega_y(i) = ib_data(i, 12); omega_z(i) = ib_data(i, 13)
1552 angle_x(i) = ib_data(i, 14); angle_y(i) = ib_data(i, 15); angle_z(i) = ib_data(i, 16)
1553 px(i) = ib_data(i, 17); py(i) = ib_data(i, 18); pz(i) = ib_data(i, 19)
1554 ib_diameter(i) = ib_data(i, 20)*2.0_wp
1555 end do
1556
1557 if (proc_rank == 0) then
1558 do i = 1, num_procs
1559 write (meshnames(i), '(A,I0,A,I0,A)') '../p', i - 1, '/', t_step, '.silo:ib_bodies'
1560 meshtypes(i) = db_pointmesh
1561 end do
1562 err = dbset2dstrlen(len(meshnames(1)))
1563 err = dbputmmesh(dbroot, 'ib_bodies', 16, num_procs, meshnames, len_trim(meshnames), meshtypes, db_f77null, ierr)
1564 end if
1565
1566 err = dbputpm(dbfile, 'ib_bodies', 9, 3, px, py, pz, nbodies, db_double, db_f77null, ierr)
1567
1568 call s_write_ib_variable('ib_force_x', t_step, force_x, nbodies)
1569 call s_write_ib_variable('ib_force_y', t_step, force_y, nbodies)
1570 call s_write_ib_variable('ib_force_z', t_step, force_z, nbodies)
1571 call s_write_ib_variable('ib_torque_x', t_step, torque_x, nbodies)
1572 call s_write_ib_variable('ib_torque_y', t_step, torque_y, nbodies)
1573 call s_write_ib_variable('ib_torque_z', t_step, torque_z, nbodies)
1574 call s_write_ib_variable('ib_vel_x', t_step, vel_x, nbodies)
1575 call s_write_ib_variable('ib_vel_y', t_step, vel_y, nbodies)
1576 call s_write_ib_variable('ib_vel_z', t_step, vel_z, nbodies)
1577 call s_write_ib_variable('ib_omega_x', t_step, omega_x, nbodies)
1578 call s_write_ib_variable('ib_omega_y', t_step, omega_y, nbodies)
1579 call s_write_ib_variable('ib_omega_z', t_step, omega_z, nbodies)
1580 call s_write_ib_variable('ib_angle_x', t_step, angle_x, nbodies)
1581 call s_write_ib_variable('ib_angle_y', t_step, angle_y, nbodies)
1582 call s_write_ib_variable('ib_angle_z', t_step, angle_z, nbodies)
1583 call s_write_ib_variable('ib_diameter', t_step, ib_diameter, nbodies)
1584
1585 deallocate (ib_data, px, py, pz, force_x, force_y, force_z)
1586 deallocate (torque_x, torque_y, torque_z, vel_x, vel_y, vel_z)
1587 deallocate (omega_x, omega_y, omega_z, angle_x, angle_y, angle_z)
1588 deallocate (ib_diameter)
1589 end if
1590#endif
1591
1593
1594 !> Write a single IB point-variable to the Silo database slave and master files.
1595 subroutine s_write_ib_variable(varname, t_step, data, nBodies)
1596
1597 character(len=*), intent(in) :: varname
1598 integer, intent(in) :: t_step
1599 real(wp), dimension(:), intent(in) :: data
1600 integer, intent(in) :: nBodies
1601 character(len=4*name_len), dimension(num_procs) :: var_names
1602 integer, dimension(num_procs) :: var_types
1603 integer :: ierr, i
1604
1605 if (proc_rank == 0) then
1606 do i = 1, num_procs
1607 write (var_names(i), '(A,I0,A,I0,A)') '../p', i - 1, '/', t_step, '.silo:' // trim(varname)
1608 var_types(i) = db_pointvar
1609 end do
1610 err = dbset2dstrlen(len(var_names(1)))
1611 err = dbputmvar(dbroot, trim(varname), len_trim(varname), num_procs, var_names, len_trim(var_names), var_types, &
1612 & db_f77null, ierr)
1613 end if
1614
1615 err = dbputpv1(dbfile, trim(varname), len_trim(varname), 'ib_bodies', 9, data, nbodies, db_double, db_f77null, ierr)
1616
1617 end subroutine s_write_ib_variable
1618
1619 !> Close the formatted database slave file and, for the root process, the master file.
1621
1622 integer :: ierr
1623
1624 if (format == 1) then
1625 ierr = dbclose(dbfile)
1626 if (proc_rank == 0) ierr = dbclose(dbroot)
1627 else
1628 close (dbfile)
1629 if (n == 0 .and. proc_rank == 0) close (dbroot)
1630 end if
1631
1632 end subroutine s_close_formatted_database_file
1633
1634 !> Close the interface data file.
1635 impure subroutine s_close_intf_data_file()
1636
1637 close (211)
1638
1639 end subroutine s_close_intf_data_file
1640
1641 !> Close the energy data file.
1642 impure subroutine s_close_energy_data_file()
1643
1644 close (251)
1645
1646 end subroutine s_close_energy_data_file
1647
1648 !> Deallocate module arrays and release all data-output resources.
1650
1651 deallocate (q_sf)
1652 if (n == 0) deallocate (q_root_sf)
1653 if (grid_geometry == 3) then
1654 deallocate (cyl_q_sf)
1655 end if
1656
1657 ! Deallocating spatial and data extents and also the variables for the offsets and the one bookkeeping the number of
1658 ! cell-boundaries in each active coordinate direction. Note that all these variables were only needed by Silo-HDF5 format
1659 ! for multidimensional data.
1660 if (format == 1) then
1661 deallocate (spatial_extents)
1662 deallocate (data_extents)
1663 deallocate (lo_offset)
1664 deallocate (hi_offset)
1665 deallocate (dims)
1666 end if
1667
1668 end subroutine s_finalize_data_output_module
1669
1670end module m_data_output
type(scalar_field), dimension(sys_size), intent(inout) q_cons_vf
integer, intent(in) k
integer, intent(in) j
integer, intent(in) l
Platform-specific file and directory operations: create, delete, inquire, getcwd, and basename.
impure subroutine my_inquire(fileloc, dircheck)
Inquires on the existence of a directory.
impure subroutine s_create_directory(dir_name)
Create a directory and all its parents if it does not exist.
Writes post-processed grid and flow-variable data to Silo-HDF5 or binary database files.
impure subroutine, public s_write_grid_to_formatted_database_file(t_step)
Write the computational grid (cell-boundary coordinates) to the formatted database slave and master f...
real(wp), dimension(:,:,:), allocatable cyl_q_sf
integer, dimension(:), allocatable hi_offset
real(wp), dimension(:,:), allocatable data_extents
impure subroutine, public s_write_variable_to_formatted_database_file(varname, t_step)
Write a single flow variable field to the formatted database slave and master files for a given time ...
impure subroutine, public s_open_energy_data_file()
Open the energy data file for appending volume-integrated energy budget quantities.
impure subroutine, public s_open_intf_data_file()
Open the interface data file for appending extracted interface coordinates.
character(len=path_len+name_len) dbdir
impure subroutine, public s_write_ib_state_files()
Convert the binary immersed-boundary state file to per-body formatted text files.
impure subroutine, public s_write_energy_data_file(q_prim_vf, q_cons_vf)
Compute volume-integrated kinetic, potential, and internal energies and write the energy budget to th...
impure subroutine, public s_write_lag_bubbles_to_formatted_database_file(t_step)
Read Lagrangian bubble restart data and write bubble positions and scalar fields to the Silo database...
integer, dimension(:), allocatable lo_offset
impure subroutine, public s_write_intf_data_file(q_prim_vf)
Extract the volume-fraction interface contour from primitive fields and write the coordinates to the ...
real(wp), dimension(:,:,:), allocatable q_root_sf
subroutine s_write_ib_variable(varname, t_step, data, nbodies)
Write a single IB point-variable to the Silo database slave and master files.
impure subroutine, public s_initialize_data_output_module()
Allocate storage arrays, configure output directories, and count flow variables for formatted databas...
real(wp), dimension(:,:,:), allocatable, public q_sf
impure subroutine, public s_close_energy_data_file()
Close the energy data file.
impure subroutine, public s_close_formatted_database_file()
Close the formatted database slave file and, for the root process, the master file.
real(sp), dimension(:,:,:), allocatable cyl_q_sf_s
impure subroutine, public s_open_formatted_database_file(t_step)
Open (or create) the Silo-HDF5 or Binary formatted database slave and master files for a given time s...
impure subroutine, public s_close_intf_data_file()
Close the interface data file.
impure subroutine, public s_write_lag_bubbles_results_to_text(t_step)
Write the post-processed results in the folder 'lag_bubbles_data'.
impure subroutine, public s_finalize_data_output_module()
Deallocate module arrays and release all data-output resources.
real(sp), dimension(:,:,:), allocatable q_root_sf_s
integer, private err
integer, dimension(:), allocatable dims
real(wp), dimension(:,:), allocatable spatial_extents
impure subroutine, public s_define_output_region
Compute the cell-index bounds for the user-specified partial output domain in each coordinate directi...
impure subroutine, public s_write_ib_bodies_to_formatted_database_file(t_step)
Read IB state and write a Silo point mesh with per-body scalar fields.
character(len=path_len+2 *name_len) proc_rank_dir
subroutine s_write_lag_variable_to_formatted_database_file(varname, t_step, data, nbubs)
Write a single Lagrangian bubble point-variable to the Silo database slave and master files.
real(sp), dimension(:,:,:), allocatable, public q_sf_s
character(len=path_len+2 *name_len) rootdir
Shared derived types for field data, patch geometry, bubble dynamics, and MPI I/O structures.
Computes derived flow quantities (sound speed, vorticity, Schlieren, etc.) from conservative and prim...
Global parameters for the post-process: domain geometry, equation of state, and output database setti...
logical cont_damage
Continuum damage modeling.
logical hypoelasticity
Turn hypoelasticity on.
type(int_bounds_info) offset_y
integer num_fluids
Number of different fluids present in the flow.
real(wp), dimension(:,:), allocatable, public mpi_io_data_lg_bubbles
logical, dimension(3) flux_wrt
real(wp), dimension(:), allocatable y_cc
integer proc_rank
Rank of the local processor.
logical output_partial_domain
Specify portion of domain to output for post-processing.
real(wp), dimension(:), allocatable adv
Advection variables.
type(int_bounds_info) z_output_idx
Indices of domain to output for post-processing.
logical dummy
AMDFlang workaround for case-optimization + GPU-kernel bug.
real(wp), dimension(:), allocatable y_cb
character(len=name_len) mpiiofs
real(wp), dimension(:), allocatable dz
logical, dimension(3) mom_wrt
real(wp), dimension(:), allocatable x_root_cb
logical, dimension(num_fluids_max) alpha_wrt
logical, dimension(num_fluids_max) alpha_rho_wrt
integer model_eqns
Multicomponent flow model.
integer precision
Floating point precision of the database file(s).
real(wp), dimension(:), allocatable z_cb
type(bounds_info) z_output
Portion of domain to output for post-processing.
integer num_dims
Number of spatial dimensions.
type(int_bounds_info) x_output_idx
real(wp), dimension(:), allocatable x_cc
integer num_vels
Number of velocity components (different from num_dims for mhd).
real(wp), dimension(:), allocatable x_cb
real(wp), dimension(:), allocatable dy
type(int_bounds_info) offset_x
logical hyper_cleaning
Hyperbolic cleaning for MHD.
real(wp), dimension(:), allocatable z_cc
logical, dimension(3) omega_wrt
integer num_procs
Number of processors.
character(len=path_len) case_dir
Case folder location.
type(int_bounds_info) y_output_idx
type(int_bounds_info) offset_z
logical mhd
Magnetohydrodynamics.
logical parallel_io
Format of the data files.
logical, dimension(3) vel_wrt
logical relativity
Relativity for RMHD.
real(wp), dimension(:), allocatable dx
Cell-width distributions in the x-, y- and z-coordinate directions.
type(eqn_idx_info) eqn_idx
All conserved-variable equation index ranges and scalars.
integer num_ibs
Number of immersed boundaries.
Utility routines for bubble model setup, coordinate transforms, array sampling, and special functions...
MPI gather and scatter operations for distributing post-process grid and flow-variable data.
impure subroutine s_mpi_defragment_1d_grid_variable
Collect the sub-domain cell-boundary or cell-center location data from all processors and put back to...
impure subroutine s_mpi_defragment_1d_flow_variable(q_sf, q_root_sf)
Gather the sub-domain flow variable data from all processors and reassemble it for the entire computa...
impure subroutine s_mpi_gather_spatial_extents(spatial_extents)
Gather spatial extents from all ranks for Silo database metadata.
impure subroutine s_mpi_gather_data_extents(q_sf, data_extents)
Gather the Silo database metadata for the flow variable's extents to boost performance of the multidi...
Conservative-to-primitive variable conversion, mixture property evaluation, and pressure computation.
real(wp), dimension(:), allocatable, public gammas
subroutine s_compute_speed_of_sound(pres, rho, gamma, pi_inf, h, adv, vel_sum, c_c, c, qv)
Compute the speed of sound from thermodynamic state variables, supporting multiple equation-of-state ...
real(wp), dimension(:), allocatable, public qvs
real(wp), dimension(:), allocatable, public pi_infs
Derived type annexing a scalar field (SF).