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