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