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/pre_process/m_data_output.fpp"
2!>
3!! @file
4!! @brief Contains module m_data_output
5
6!> @brief Writes grid and initial condition data to serial or parallel output files
8
11 use m_helper
12 use m_mpi_proxy
13#ifdef MFC_MPI
14 use mpi
15#endif
16
19 use m_helper
23 use m_thermochem, only: species_names
24 use m_helper
26
27 implicit none
28
29 private
32
33 type(scalar_field), allocatable, dimension(:) :: q_cons_temp
34
35 abstract interface
36
37 !> Interface for the conservative data
38 impure subroutine s_write_abstract_data_files(q_cons_vf, q_prim_vf, bc_type, q_T_sf)
39
40 import :: scalar_field, integer_field, sys_size, m, n, p, pres_field, num_dims
41
42 type(scalar_field), dimension(sys_size), intent(inout) :: q_cons_vf, q_prim_vf
43 type(integer_field), dimension(1:num_dims,-1:1), intent(in) :: bc_type
44 type(scalar_field), intent(inout), optional :: q_t_sf
45
46 end subroutine s_write_abstract_data_files
47 end interface
48
49 !> Time-step folder into which grid and initial condition data will be placed
50 character(LEN=path_len + 2*name_len), private :: t_step_dir
51 character(LEN=path_len + 2*name_len), public :: restart_dir !< Restart data folder
52 procedure(s_write_abstract_data_files), pointer :: s_write_data_files => null()
53
54contains
55
56 !> Writes grid and initial condition data files to the "0" time-step directory in the local processor rank folder
57 impure subroutine s_write_serial_data_files(q_cons_vf, q_prim_vf, bc_type, q_T_sf)
58
59 type(scalar_field), dimension(sys_size), intent(inout) :: q_cons_vf, q_prim_vf
60 type(integer_field), dimension(1:num_dims,-1:1), intent(in) :: bc_type
61 logical :: file_exist
62 character(LEN=15) :: fmt
63 character(LEN=3) :: status
64 character(len=int(floor(log10(real(sys_size, wp)))) + 1) :: file_num
65 character(LEN=len_trim(t_step_dir) + name_len) :: file_loc
66 integer :: i, j, k, l, r, c
67 integer :: t_step
68 real(wp), dimension(nb) :: nrtmp
69 real(wp) :: nbub
70 real(wp) :: gamma, lit_gamma, pi_inf, qv
71 real(wp) :: rho
72 real(wp) :: pres, t
73 real(wp) :: rhoyks(1:num_species)
74 real(wp) :: pres_mag
75 type(scalar_field), intent(inout), optional :: q_t_sf
76
77 pres_mag = 0._wp
78
79 t = dflt_t_guess
80
81 t_step = 0
82
83 if (old_grid) then
84 status = 'old'
85 else
86 status = 'new'
87 end if
88
89 if (bc_io) then
90 if (igr) then
91 call s_write_serial_boundary_condition_files(q_cons_vf, bc_type, t_step_dir, old_grid)
92 else
93 call s_write_serial_boundary_condition_files(q_prim_vf, bc_type, t_step_dir, old_grid, q_t_sf)
94 end if
95 end if
96
97 file_loc = trim(t_step_dir) // '/x_cb.dat'
98 open (1, file=trim(file_loc), form='unformatted', status=status)
99 write (1) x_cb(-1:m)
100 close (1)
101
102 if (n > 0) then
103 file_loc = trim(t_step_dir) // '/y_cb.dat'
104 open (1, file=trim(file_loc), form='unformatted', status=status)
105 write (1) y_cb(-1:n)
106 close (1)
107
108 if (p > 0) then
109 file_loc = trim(t_step_dir) // '/z_cb.dat'
110 open (1, file=trim(file_loc), form='unformatted', status=status)
111 write (1) z_cb(-1:p)
112 close (1)
113 end if
114 end if
115
116 do i = 1, sys_size
117 write (file_num, '(I0)') i
118 file_loc = trim(t_step_dir) // '/q_cons_vf' // trim(file_num) // '.dat'
119 open (1, file=trim(file_loc), form='unformatted', status=status)
120 write (1) q_cons_vf(i)%sf(0:m,0:n,0:p)
121 close (1)
122 end do
123
124 if (qbmm .and. .not. polytropic) then
125 do i = 1, nb
126 do r = 1, nnode
127 write (file_num, '(I0)') r + (i - 1)*nnode + sys_size
128 file_loc = trim(t_step_dir) // '/pb' // trim(file_num) // '.dat'
129 open (1, file=trim(file_loc), form='unformatted', status=status)
130 write (1) pb%sf(:,:,:,r, i)
131 close (1)
132 end do
133 end do
134
135 do i = 1, nb
136 do r = 1, nnode
137 write (file_num, '(I0)') r + (i - 1)*nnode + sys_size
138 file_loc = trim(t_step_dir) // '/mv' // trim(file_num) // '.dat'
139 open (1, file=trim(file_loc), form='unformatted', status=status)
140 write (1) mv%sf(:,:,:,r, i)
141 close (1)
142 end do
143 end do
144 end if
145
146 gamma = gammas(1)
147 lit_gamma = gs_min(1)
148 pi_inf = pi_infs(1)
149 qv = qvs(1)
150
151 if (precision == precision_single) then
152 fmt = "(2F30.3)"
153 else
154 fmt = "(2F40.14)"
155 end if
156
157 write (t_step_dir, '(A,I0,A,I0)') trim(case_dir) // '/D'
158 file_loc = trim(t_step_dir) // '/.'
159
160 inquire (file=trim(file_loc), exist=file_exist)
161
162 if (.not. file_exist) call s_create_directory(trim(t_step_dir))
163
164 if (cfl_dt) t_step = n_start
165
166 if (n == 0 .and. p == 0) then
167 if (model_eqns == model_eqns_5eq) then
168 do i = 1, sys_size
169 write (file_loc, '(A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/prim.', i, '.', proc_rank, '.', t_step, '.dat'
170
171 open (2, file=trim(file_loc))
172 do j = 0, m
173 if (chemistry) then
174 do c = 1, num_species
175 rhoyks(c) = q_cons_vf(eqn_idx%species%beg + c - 1)%sf(j, 0, 0)
176 end do
177 end if
178
179 call s_convert_to_mixture_variables(q_cons_vf, j, 0, 0, rho, gamma, pi_inf, qv)
180
181 lit_gamma = 1._wp/gamma + 1._wp
182
183 if ((i >= eqn_idx%species%beg) .and. (i <= eqn_idx%species%end)) then
184 write (2, fmt) x_cb(j), q_cons_vf(i)%sf(j, 0, 0)/rho
185 else if (((i >= eqn_idx%cont%beg) .and. (i <= eqn_idx%cont%end)) .or. ((i >= eqn_idx%adv%beg) &
186 & .and. (i <= eqn_idx%adv%end)) .or. ((i >= eqn_idx%species%beg) .and. (i <= eqn_idx%species%end) &
187 & )) then
188 write (2, fmt) x_cb(j), q_cons_vf(i)%sf(j, 0, 0)
189 else if (i == eqn_idx%mom%beg) then ! u
190 write (2, fmt) x_cb(j), q_cons_vf(eqn_idx%mom%beg)%sf(j, 0, 0)/rho
191 else if (i == eqn_idx%stress%beg) then ! tau_e
192 write (2, fmt) x_cb(j), q_cons_vf(eqn_idx%stress%beg)%sf(j, 0, 0)/rho
193 else if (i == eqn_idx%E) then ! p
194 if (mhd) then
195 pres_mag = 0.5_wp*(bx0**2 + q_cons_vf(eqn_idx%B%beg)%sf(j, 0, &
196 & 0)**2 + q_cons_vf(eqn_idx%B%beg + 1)%sf(j, 0, 0)**2)
197 end if
198
199 call s_compute_pressure(q_cons_vf(eqn_idx%E)%sf(j, 0, 0), q_cons_vf(eqn_idx%alf)%sf(j, 0, 0), &
200 & 0.5_wp*(q_cons_vf(eqn_idx%mom%beg)%sf(j, 0, 0)**2._wp)/rho, pi_inf, gamma, &
201 & rho, qv, rhoyks, pres, t, pres_mag=pres_mag)
202 write (2, fmt) x_cb(j), pres
203 else if (mhd) then
204 if (i == eqn_idx%mom%beg + 1) then ! v
205 write (2, fmt) x_cb(j), q_cons_vf(eqn_idx%mom%beg + 1)%sf(j, 0, 0)/rho
206 else if (i == eqn_idx%mom%beg + 2) then ! w
207 write (2, fmt) x_cb(j), q_cons_vf(eqn_idx%mom%beg + 2)%sf(j, 0, 0)/rho
208 else if (i == eqn_idx%B%beg) then ! By
209 write (2, fmt) x_cb(j), q_cons_vf(eqn_idx%B%beg)%sf(j, 0, 0)/rho
210 else if (i == eqn_idx%B%beg + 1) then ! Bz
211 write (2, fmt) x_cb(j), q_cons_vf(eqn_idx%B%beg + 1)%sf(j, 0, 0)/rho
212 end if
213 else if ((i >= eqn_idx%bub%beg) .and. (i <= eqn_idx%bub%end) .and. bubbles_euler) then
214 if (qbmm) then
215 nbub = q_cons_vf(eqn_idx%bub%beg)%sf(j, 0, 0)
216 else
217 if (adv_n) then
218 nbub = q_cons_vf(eqn_idx%n)%sf(j, 0, 0)
219 else
220 do k = 1, nb
221 nrtmp(k) = q_cons_vf(qbmm_idx%rs(k))%sf(j, 0, 0)
222 end do
223
224 call s_comp_n_from_cons(real(q_cons_vf(eqn_idx%alf)%sf(j, 0, 0), kind=wp), nrtmp, nbub, weight)
225 end if
226 end if
227 write (2, fmt) x_cb(j), q_cons_vf(i)%sf(j, 0, 0)/nbub
228 else if (i == eqn_idx%n .and. adv_n .and. bubbles_euler) then
229 write (2, fmt) x_cb(j), q_cons_vf(i)%sf(j, 0, 0)
230 else if (i == eqn_idx%damage) then
231 write (2, fmt) x_cb(j), q_cons_vf(i)%sf(j, 0, 0)
232 end if
233 end do
234 close (2)
235 end do
236 end if
237
238 do i = 1, sys_size
239 write (file_loc, '(A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/cons.', i, '.', proc_rank, '.', t_step, '.dat'
240
241 open (2, file=trim(file_loc))
242 do j = 0, m
243 write (2, fmt) x_cb(j), q_cons_vf(i)%sf(j, 0, 0)
244 end do
245 close (2)
246 end do
247
248 if (qbmm .and. .not. polytropic) then
249 do i = 1, nb
250 do r = 1, nnode
251 write (file_loc, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/pres.', i, '.', r, '.', proc_rank, &
252 & '.', t_step, '.dat'
253
254 open (2, file=trim(file_loc))
255 do j = 0, m
256 write (2, fmt) x_cb(j), pb%sf(j, 0, 0, r, i)
257 end do
258 close (2)
259 end do
260 end do
261 do i = 1, nb
262 do r = 1, nnode
263 write (file_loc, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/mv.', i, '.', r, '.', proc_rank, &
264 & '.', t_step, '.dat'
265
266 open (2, file=trim(file_loc))
267 do j = 0, m
268 write (2, fmt) x_cb(j), mv%sf(j, 0, 0, r, i)
269 end do
270 close (2)
271 end do
272 end do
273 end if
274 end if
275
276 if (precision == precision_single) then
277 fmt = "(3F30.7)"
278 else
279 fmt = "(3F40.14)"
280 end if
281
282 if ((n > 0) .and. (p == 0)) then
283 do i = 1, sys_size
284 write (file_loc, '(A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/cons.', i, '.', proc_rank, '.', t_step, '.dat'
285 open (2, file=trim(file_loc))
286 do j = 0, m
287 do k = 0, n
288 write (2, fmt) x_cb(j), y_cb(k), q_cons_vf(i)%sf(j, k, 0)
289 end do
290 write (2, *)
291 end do
292 close (2)
293 end do
294
295 if (qbmm .and. .not. polytropic) then
296 do i = 1, nb
297 do r = 1, nnode
298 write (file_loc, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/pres.', i, '.', r, '.', proc_rank, &
299 & '.', t_step, '.dat'
300
301 open (2, file=trim(file_loc))
302 do j = 0, m
303 do k = 0, n
304 write (2, fmt) x_cb(j), y_cb(k), pb%sf(j, k, 0, r, i)
305 end do
306 end do
307 close (2)
308 end do
309 end do
310 do i = 1, nb
311 do r = 1, nnode
312 write (file_loc, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/mv.', i, '.', r, '.', proc_rank, &
313 & '.', t_step, '.dat'
314
315 open (2, file=trim(file_loc))
316 do j = 0, m
317 do k = 0, n
318 write (2, fmt) x_cb(j), y_cb(k), mv%sf(j, k, 0, r, i)
319 end do
320 end do
321 close (2)
322 end do
323 end do
324 end if
325 end if
326
327 if (precision == precision_single) then
328 fmt = "(4F30.7)"
329 else
330 fmt = "(4F40.14)"
331 end if
332
333 if (p > 0) then
334 do i = 1, sys_size
335 write (file_loc, '(A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/cons.', i, '.', proc_rank, '.', t_step, '.dat'
336 open (2, file=trim(file_loc))
337 do j = 0, m
338 do k = 0, n
339 do l = 0, p
340 write (2, fmt) x_cb(j), y_cb(k), z_cb(l), q_cons_vf(i)%sf(j, k, l)
341 end do
342 write (2, *)
343 end do
344 write (2, *)
345 end do
346 close (2)
347 end do
348
349 if (qbmm .and. .not. polytropic) then
350 do i = 1, nb
351 do r = 1, nnode
352 write (file_loc, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/pres.', i, '.', r, '.', proc_rank, &
353 & '.', t_step, '.dat'
354
355 open (2, file=trim(file_loc))
356 do j = 0, m
357 do k = 0, n
358 do l = 0, p
359 write (2, fmt) x_cb(j), y_cb(k), z_cb(l), pb%sf(j, k, l, r, i)
360 end do
361 end do
362 end do
363 close (2)
364 end do
365 end do
366 do i = 1, nb
367 do r = 1, nnode
368 write (file_loc, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir) // '/mv.', i, '.', r, '.', proc_rank, &
369 & '.', t_step, '.dat'
370
371 open (2, file=trim(file_loc))
372 do j = 0, m
373 do k = 0, n
374 do l = 0, p
375 write (2, fmt) x_cb(j), y_cb(k), z_cb(l), mv%sf(j, k, l, r, i)
376 end do
377 end do
378 end do
379 close (2)
380 end do
381 end do
382 end if
383 end if
384
385 end subroutine s_write_serial_data_files
386
387 !> Writes grid and initial condition data files in parallel to the "0" time-step directory in the local processor rank folder
388 impure subroutine s_write_parallel_data_files(q_cons_vf, q_prim_vf, bc_type, q_T_sf)
389
390 type(scalar_field), dimension(sys_size), intent(inout) :: q_cons_vf, q_prim_vf
391 type(integer_field), dimension(1:num_dims,-1:1), intent(in) :: bc_type
392 type(scalar_field), optional, intent(inout) :: q_t_sf
393
394#ifdef MFC_MPI
395 integer :: ifile, ierr, data_size
396 integer, dimension(MPI_STATUS_SIZE) :: status
397 integer(KIND=MPI_OFFSET_KIND) :: disp
398 integer(KIND=MPI_OFFSET_KIND) :: m_mok, n_mok, p_mok
399 integer(KIND=MPI_OFFSET_KIND) :: wp_mok, var_mok, str_mok
400 integer(KIND=MPI_OFFSET_KIND) :: nvars_mok
401 integer(KIND=MPI_OFFSET_KIND) :: mok
402 character(LEN=path_len + 2*name_len) :: file_loc
403 logical :: file_exist, dir_check
404 integer :: i, j, k, l
405 real(wp) :: loc_violations, glb_violations
406 integer :: m_ds, n_ds, p_ds
407 integer :: m_glb_ds, n_glb_ds, p_glb_ds
408 integer :: m_glb_save, n_glb_save, p_glb_save !< Size of array being saved
409
410 loc_violations = 0._wp
411
412 if (down_sample) then
413 if ((mod(m + 1, 3) > 0) .or. (mod(n + 1, 3) > 0) .or. (mod(p + 1, 3) > 0)) then
414 loc_violations = 1._wp
415 end if
416 call s_mpi_allreduce_sum(loc_violations, glb_violations)
417 if (proc_rank == 0 .and. nint(glb_violations) > 0) then
418 print *, &
419 & "WARNING: Attempting to downsample data but there are" &
420 & // "processors with local problem sizes that are not divisible by 3."
421 end if
423 call s_downsample_data(q_cons_vf, q_cons_temp, m_ds, n_ds, p_ds, m_glb_ds, n_glb_ds, p_glb_ds)
424 end if
425
426 if (file_per_process) then
427 if (proc_rank == 0) then
428 file_loc = trim(case_dir) // '/restart_data/lustre_0'
429 call my_inquire(file_loc, dir_check)
430 if (dir_check .neqv. .true.) then
431 call s_create_directory(trim(file_loc))
432 end if
433 call s_create_directory(trim(file_loc))
434 end if
435 call s_mpi_barrier()
437
438 if (down_sample) then
439 call s_initialize_mpi_data_ds(q_cons_temp)
440 else
441 call s_initialize_mpi_data(q_cons_vf)
442 end if
443
444 if (cfl_dt) then
445 write (file_loc, '(I0,A,i7.7,A)') n_start, '_', proc_rank, '.dat'
446 else
447 write (file_loc, '(I0,A,i7.7,A)') t_step_start, '_', proc_rank, '.dat'
448 end if
449 file_loc = trim(restart_dir) // '/lustre_0' // trim(mpiiofs) // trim(file_loc)
450 inquire (file=trim(file_loc), exist=file_exist)
451 if (file_exist .and. proc_rank == 0) then
452 call mpi_file_delete(file_loc, mpi_info_int, ierr)
453 end if
454 if (file_exist) call mpi_file_delete(file_loc, mpi_info_int, ierr)
455 call mpi_file_open(mpi_comm_self, file_loc, ior(mpi_mode_wronly, mpi_mode_create), mpi_info_int, ifile, ierr)
456
457 if (down_sample) then
458 data_size = (m_ds + 3)*(n_ds + 3)*(p_ds + 3)
459 m_glb_save = m_glb_ds + 3
460 n_glb_save = n_glb_ds + 3
461 p_glb_save = p_glb_ds + 3
462 else
463 data_size = (m + 1)*(n + 1)*(p + 1)
464 m_glb_save = m_glb + 1
465 n_glb_save = n_glb + 1
466 p_glb_save = p_glb + 1
467 end if
468
469 ! Resize some integers so MPI can write even the biggest files
470 m_mok = int(m_glb_save, mpi_offset_kind)
471 n_mok = int(n_glb_save, mpi_offset_kind)
472 p_mok = int(p_glb_save, mpi_offset_kind)
473 wp_mok = int(storage_size(0._stp)/8, mpi_offset_kind)
474 mok = int(1._wp, mpi_offset_kind)
475 str_mok = int(name_len, mpi_offset_kind)
476 nvars_mok = int(sys_size, mpi_offset_kind)
477
478 if (bubbles_euler) then
479 do i = 1, sys_size
480 var_mok = int(i, mpi_offset_kind)
481
482 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, mpi_io_p, status, ierr)
483 end do
484 if (qbmm .and. .not. polytropic) then
485 do i = sys_size + 1, sys_size + 2*nb*nnode
486 var_mok = int(i, mpi_offset_kind)
487
488 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, mpi_io_p, status, ierr)
489 end do
490 end if
491 else
492 if (down_sample) then
493 do i = 1, sys_size
494 var_mok = int(i, mpi_offset_kind)
495
496 call mpi_file_write_all(ifile, q_cons_temp(i)%sf, data_size*mpi_io_type, mpi_io_p, status, ierr)
497 end do
498 else
499 do i = 1, sys_size
500 var_mok = int(i, mpi_offset_kind)
501
502 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, mpi_io_p, status, ierr)
503 end do
504 end if
505 end if
506
507 call mpi_file_close(ifile, ierr)
508 else
509 call s_initialize_mpi_data(q_cons_vf)
510
511 if (cfl_dt) then
512 write (file_loc, '(I0,A)') n_start, '.dat'
513 else
514 write (file_loc, '(I0,A)') t_step_start, '.dat'
515 end if
516 file_loc = trim(restart_dir) // trim(mpiiofs) // trim(file_loc)
517 inquire (file=trim(file_loc), exist=file_exist)
518 if (file_exist .and. proc_rank == 0) then
519 call mpi_file_delete(file_loc, mpi_info_int, ierr)
520 end if
521 call mpi_file_open(mpi_comm_world, file_loc, ior(mpi_mode_wronly, mpi_mode_create), mpi_info_int, ifile, ierr)
522
523 data_size = (m + 1)*(n + 1)*(p + 1)
524
525 ! Resize some integers so MPI can write even the biggest files
526 m_mok = int(m_glb + 1, mpi_offset_kind)
527 n_mok = int(n_glb + 1, mpi_offset_kind)
528 p_mok = int(p_glb + 1, mpi_offset_kind)
529 wp_mok = int(storage_size(0._stp)/8, mpi_offset_kind)
530 mok = int(1._wp, mpi_offset_kind)
531 str_mok = int(name_len, mpi_offset_kind)
532 nvars_mok = int(sys_size, mpi_offset_kind)
533
534 if (bubbles_euler) then
535 do i = 1, sys_size
536 var_mok = int(i, mpi_offset_kind)
537
538 disp = m_mok*max(mok, n_mok)*max(mok, p_mok)*wp_mok*(var_mok - 1)
539
540 call mpi_file_set_view(ifile, disp, mpi_io_p, mpi_io_data%view(i), 'native', mpi_info_int, ierr)
541 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, mpi_io_p, status, ierr)
542 end do
543 if (qbmm .and. .not. polytropic) then
544 do i = sys_size + 1, sys_size + 2*nb*nnode
545 var_mok = int(i, mpi_offset_kind)
546
547 disp = m_mok*max(mok, n_mok)*max(mok, p_mok)*wp_mok*(var_mok - 1)
548
549 call mpi_file_set_view(ifile, disp, mpi_io_p, mpi_io_data%view(i), 'native', mpi_info_int, ierr)
550 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, mpi_io_p, status, ierr)
551 end do
552 end if
553 else
554 do i = 1, sys_size
555 var_mok = int(i, mpi_offset_kind)
556
557 disp = m_mok*max(mok, n_mok)*max(mok, p_mok)*wp_mok*(var_mok - 1)
558
559 call mpi_file_set_view(ifile, disp, mpi_io_p, mpi_io_data%view(i), 'native', mpi_info_int, ierr)
560 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, mpi_io_p, status, ierr)
561 end do
562 end if
563
564 call mpi_file_close(ifile, ierr)
565 end if
566#endif
567
568 if (bc_io) then
569 if (igr) then
570 call s_write_parallel_boundary_condition_files(q_cons_vf, bc_type)
571 else
572 call s_write_parallel_boundary_condition_files(q_prim_vf, bc_type, q_t_sf)
573 end if
574 end if
575
576 end subroutine s_write_parallel_data_files
577
578 !> Computation of parameters, allocation procedures, and/or any other tasks needed to properly setup the module
580
581 character(LEN=len_trim(case_dir) + 2*name_len) :: file_loc
582 character(len=15) :: temp
583 character(LEN=1), dimension(3), parameter :: coord = (/'x', 'y', 'z'/)
584 logical :: dir_check
585 integer :: i, iu
586 integer :: m_ds, n_ds, p_ds
587
588 if (parallel_io .neqv. .true.) then
589 write (t_step_dir, '(A,I0,A)') '/p_all/p', proc_rank, '/0'
590 t_step_dir = trim(case_dir) // trim(t_step_dir)
591
592 if (old_grid .neqv. .true.) then
593 file_loc = trim(t_step_dir) // '/'
594
595 call my_inquire(file_loc, dir_check)
596
597 if (dir_check) call s_delete_directory(trim(t_step_dir))
598
600 end if
601
603 else
604 write (restart_dir, '(A)') '/restart_data'
605 restart_dir = trim(case_dir) // trim(restart_dir)
606
607 if ((old_grid .neqv. .true.) .and. (proc_rank == 0)) then
608 file_loc = trim(restart_dir) // '/'
609 call my_inquire(file_loc, dir_check)
610
611 if (dir_check) call s_delete_directory(trim(restart_dir))
613 end if
614
615 call s_mpi_barrier()
616
618 end if
619
620 open (newunit=iu, file='indices.dat', status='unknown')
621
622 write (iu, '(A)') "Warning: The creation of file is currently experimental."
623 write (iu, '(A)') "This file may contain errors and not support all features."
624
625 write (iu, '(A3,A20,A20)') "#", "Conservative", "Primitive"
626 write (iu, '(A)') " "
627 do i = eqn_idx%cont%beg, eqn_idx%cont%end
628 write (temp, '(I0)') i - eqn_idx%cont%beg + 1
629 write (iu, '(I3,A20,A20)') i, "\alpha_{" // trim(temp) // "} \rho_{" // trim(temp) // "}", &
630 & "\alpha_{" // trim(temp) // "} \rho"
631 end do
632 do i = eqn_idx%mom%beg, eqn_idx%mom%end
633 write (iu, '(I3,A20,A20)') i, "\rho u_" // coord(i - eqn_idx%mom%beg + 1), "u_" // coord(i - eqn_idx%mom%beg + 1)
634 end do
635 if (eqn_idx%E /= 0) write (iu, '(I3,A20,A20)') eqn_idx%E, "\rho U", "p"
636 do i = eqn_idx%adv%beg, eqn_idx%adv%end
637 write (temp, '(I0)') i - eqn_idx%cont%beg + 1
638 write (iu, '(I3,A20,A20)') i, "\alpha_{" // trim(temp) // "}", "\alpha_{" // trim(temp) // "}"
639 end do
640 if (chemistry) then
641 do i = 1, num_species
642 write (iu, '(I3,A20,A20)') eqn_idx%species%beg + i - 1, "Y_{" // trim(species_names(i)) // "} \rho", &
643 & "Y_{" // trim(species_names(i)) // "}"
644 end do
645 end if
646
647 write (iu, '(A)') ""
648 call write_range(eqn_idx%cont%beg, eqn_idx%cont%end, " Continuity")
649 call write_range(eqn_idx%mom%beg, eqn_idx%mom%end, " Momentum")
650 call write_range(eqn_idx%E, eqn_idx%E, " Energy/Pressure")
651 call write_range(eqn_idx%adv%beg, eqn_idx%adv%end, " Advection")
652 call write_range(eqn_idx%bub%beg, eqn_idx%bub%end, " Bubbles")
653 call write_range(eqn_idx%stress%beg, eqn_idx%stress%end, " Stress")
654 call write_range(eqn_idx%int_en%beg, eqn_idx%int_en%end, " Internal Energies")
655 call write_range(eqn_idx%xi%beg, eqn_idx%xi%end, " Reference Map")
656 call write_range(eqn_idx%B%beg, eqn_idx%B%end, " Magnetic Field")
657 call write_range(eqn_idx%c, eqn_idx%c, " Color Function")
658 call write_range(eqn_idx%species%beg, eqn_idx%species%end, " Chemistry")
659
660 close (iu)
661
662 if (down_sample) then
663 m_ds = int((m + 1)/3) - 1
664 n_ds = int((n + 1)/3) - 1
665 p_ds = int((p + 1)/3) - 1
666
667 allocate (q_cons_temp(1:sys_size))
668 do i = 1, sys_size
669 allocate (q_cons_temp(i)%sf(-1:m_ds + 1,-1:n_ds + 1,-1:p_ds + 1))
670 end do
671 end if
672
673 contains
674
675 subroutine write_range(beg, end, label)
676
677 integer, intent(in) :: beg, end
678 character(*), intent(in) :: label
679
680 if (beg /= 0) write (iu, '("[",I0,",",I0,"]",A)') beg, end, label
681
682 end subroutine write_range
683
685
686 !> Resets s_write_data_files pointer
688
689 integer :: i
690
691 s_write_data_files => null()
692
693 if (down_sample) then
694 do i = 1, sys_size
695 deallocate (q_cons_temp(i)%sf)
696 end do
697 deallocate (q_cons_temp)
698 end if
699
700 end subroutine s_finalize_data_output_module
701
702end module m_data_output
Interface for the conservative data.
subroutine write_range(beg, end, label)
type(scalar_field), dimension(sys_size), intent(inout) q_cons_vf
integer, intent(in) k
integer, intent(in) j
integer, intent(in) l
Noncharacteristic and processor boundary condition application for ghost cells and buffer regions.
impure subroutine, public s_populate_variables_buffers(bc_type, q_prim_vf, pb_in, mv_in, q_t_sf)
Populate the buffers of the primitive variables based on the selected boundary conditions.
Applies spatially varying boundary condition patches along domain edges and faces.
Platform-specific file and directory operations: create, delete, inquire, getcwd, and basename.
impure subroutine s_delete_directory(dir_name)
Recursively delete a directory using a platform-specific system command.
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.
Compile-time constant parameters: default values, tolerances, and physical constants.
integer, parameter model_eqns_5eq
integer, parameter name_len
Maximum name length.
real(wp), parameter dflt_t_guess
Default guess for temperature (when a previous value is not available).
integer, parameter nnode
Number of QBMM nodes.
integer, parameter precision_single
Writes grid and initial condition data to serial or parallel output files.
type(scalar_field), dimension(:), allocatable q_cons_temp
procedure(s_write_abstract_data_files), pointer, public s_write_data_files
impure subroutine, public s_write_serial_data_files(q_cons_vf, q_prim_vf, bc_type, q_t_sf)
Writes grid and initial condition data files to the "0" time-step directory in the local processor ra...
character(len=path_len+2 *name_len), private t_step_dir
Time-step folder into which grid and initial condition data will be placed.
impure subroutine, public s_initialize_data_output_module
Computation of parameters, allocation procedures, and/or any other tasks needed to properly setup the...
impure subroutine, public s_finalize_data_output_module
Resets s_write_data_files pointer.
impure subroutine, public s_write_parallel_data_files(q_cons_vf, q_prim_vf, bc_type, q_t_sf)
Writes grid and initial condition data files in parallel to the "0" time-step directory in the local ...
character(len=path_len+2 *name_len), public restart_dir
Restart data folder.
Rank-staggered file access delays to prevent I/O contention on parallel file systems.
impure subroutine, public delayfileaccess(processrank)
Introduce a rank-dependent busy-wait delay to stagger parallel file access and reduce I/O contention.
Shared derived types for field data, patch geometry, bubble dynamics, and MPI I/O structures.
Defines global parameters for the computational domain, simulation algorithm, and initial conditions.
integer p_glb
Global number of cells in each direction.
type(qbmm_idx_info) qbmm_idx
QBMM moment index mappings.
integer proc_rank
Rank of the local processor Number of cells in the x-, y- and z-coordinate directions.
logical bc_io
whether or not to save BC data
real(wp), dimension(:), allocatable y_cb
real(wp), dimension(:), allocatable weight
real(wp), dimension(:), allocatable z_cb
real(wp), dimension(:), allocatable x_cb
Locations of cell-boundaries (cb) in x-, y- and z-directions, respectively.
type(mpi_io_var), public mpi_io_data
Utility routines for bubble model setup, coordinate transforms, array sampling, and special functions...
subroutine, public s_downsample_data(q_cons_vf, q_cons_temp, m_ds, n_ds, p_ds, m_glb_ds, n_glb_ds, p_glb_ds)
Downsample conservative variable fields by a factor of 3 in each direction using volume averaging.
subroutine, public s_comp_n_from_cons(vftmp, nrtmp, ntmp, weights)
Compute the bubble number density from the conservative void fraction and weighted bubble radii.
Broadcasts user inputs and decomposes the domain across MPI ranks for pre-processing.
Conservative-to-primitive variable conversion, mixture property evaluation, and pressure computation.
real(wp), dimension(:), allocatable, public gammas
real(wp), dimension(:), allocatable, public gs_min
real(wp), dimension(:), allocatable, public qvs
subroutine, public s_compute_pressure(energy, alf, dyn_p, pi_inf, gamma, rho, qv, rhoyks, pres, t, stress, mom, g, pres_mag)
Compute the pressure from the appropriate equation of state.
real(wp), dimension(:), allocatable, public pi_infs
subroutine, public s_convert_to_mixture_variables(q_vf, i, j, k, rho, gamma, pi_inf, qv, re_k, g_k, g)
Dispatch to the s_convert_mixture_to_mixture_variables and s_convert_species_to_mixture_variables sub...
Derived type annexing an integer scalar field (SF).
Derived type for bubble variables pb and mv at quadrature nodes (qbmm).
Derived type annexing a scalar field (SF).