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
9 use m_derived_types !< definitions of the derived types
10
11 use m_global_parameters !< global parameters for the code
12
13 use m_helper
14
15 use m_mpi_proxy !< message passing interface (mpi) module proxy
16
17#ifdef MFC_MPI
18 use mpi !< message passing interface (mpi) module
19#endif
20
22
24
25 use m_helper
26
28
30
32
33 use m_thermochem, only: species_names
34
35 use m_helper
36
37 implicit none
38
39 private;
40 public :: s_write_serial_data_files, &
45
46 type(scalar_field), allocatable, dimension(:) :: q_cons_temp
47
48 abstract interface
49
50 !> Interface for the conservative data
51 !! @param q_cons_vf Conservative variables
52 impure subroutine s_write_abstract_data_files(q_cons_vf, q_prim_vf, bc_type)
53
54 import :: scalar_field, integer_field, sys_size, m, n, p, &
56
57 ! Conservative variables
58 type(scalar_field), &
59 dimension(sys_size), &
60 intent(inout) :: q_cons_vf, q_prim_vf
61
62 type(integer_field), &
63 dimension(1:num_dims, -1:1), &
64 intent(in) :: bc_type
65
66 end subroutine s_write_abstract_data_files
67 end interface
68
69 character(LEN=path_len + 2*name_len), private :: t_step_dir !<
70 !! Time-step folder into which grid and initial condition data will be placed
71
72 character(LEN=path_len + 2*name_len), public :: restart_dir !<
73 !! Restart data folder
74
75 procedure(s_write_abstract_data_files), pointer :: s_write_data_files => null()
76
77contains
78
79 !> Writes grid and initial condition data files to the "0"
80 !! time-step directory in the local processor rank folder
81 !! @param q_cons_vf Conservative variables
82 !! @param q_prim_vf Primitive variables
83 !! @param bc_type Boundary condition types
84 impure subroutine s_write_serial_data_files(q_cons_vf, q_prim_vf, bc_type)
85 type(scalar_field), &
86 dimension(sys_size), &
87 intent(inout) :: q_cons_vf, q_prim_vf
88
89 ! BC types
90 type(integer_field), &
91 dimension(1:num_dims, -1:1), &
92 intent(in) :: bc_type
93
94 logical :: file_exist !< checks if file exists
95
96 character(LEN=15) :: fmt
97 character(LEN=3) :: status
98
99 character(len= &
100 int(floor(log10(real(sys_size, wp)))) + 1) :: file_num !< Used to store
101 !! the number, in character form, of the currently
102 !! manipulated conservative variable data file
103
104 character(LEN=len_trim(t_step_dir) + name_len) :: file_loc !<
105 !! Generic string used to store the address of a particular file
106
107 integer :: i, j, k, l, r, c !< Generic loop iterator
108 integer :: t_step
109
110 real(wp), dimension(nb) :: nrtmp !< Temporary bubble concentration
111 real(wp) :: nbub !< Temporary bubble number density
112 real(wp) :: gamma, lit_gamma, pi_inf, qv !< Temporary EOS params
113 real(wp) :: rho !< Temporary density
114 real(wp) :: pres, t !< Temporary pressure
115
116 real(wp) :: rhoyks(1:num_species) !< Temporary species mass fractions
117
118 real(wp) :: pres_mag
119
120 pres_mag = 0._wp
121
122 t = dflt_t_guess
123
124 t_step = 0
125
126 ! Outputting the Locations of the Cell-boundaries
127
128 if (old_grid) then
129 status = 'old'
130 else
131 status = 'new'
132 end if
133
134 if (bc_io) then
135 if (igr) then
137 else
139 end if
140 end if
141
142 ! x-coordinate direction
143 file_loc = trim(t_step_dir)//'/x_cb.dat'
144 open (1, file=trim(file_loc), form='unformatted', status=status)
145 write (1) x_cb(-1:m)
146 close (1)
147
148 ! y- and z-coordinate directions
149 if (n > 0) then
150 ! y-coordinate direction
151 file_loc = trim(t_step_dir)//'/y_cb.dat'
152 open (1, file=trim(file_loc), form='unformatted', &
153 status=status)
154 write (1) y_cb(-1:n)
155 close (1)
156
157 ! z-coordinate direction
158 if (p > 0) then
159 file_loc = trim(t_step_dir)//'/z_cb.dat'
160 open (1, file=trim(file_loc), form='unformatted', &
161 status=status)
162 write (1) z_cb(-1:p)
163 close (1)
164 end if
165 end if
166
167 ! Outputting Conservative Variables
168 do i = 1, sys_size
169 write (file_num, '(I0)') i
170 file_loc = trim(t_step_dir)//'/q_cons_vf'//trim(file_num) &
171 //'.dat'
172 open (1, file=trim(file_loc), form='unformatted', &
173 status=status)
174 write (1) q_cons_vf(i)%sf(0:m, 0:n, 0:p)
175 close (1)
176 end do
177
178 !Outputting pb and mv for non-polytropic qbmm
179 if (qbmm .and. .not. polytropic) then
180 do i = 1, nb
181 do r = 1, nnode
182 write (file_num, '(I0)') r + (i - 1)*nnode + sys_size
183 file_loc = trim(t_step_dir)//'/pb'//trim(file_num) &
184 //'.dat'
185 open (1, file=trim(file_loc), form='unformatted', &
186 status=status)
187 write (1) pb%sf(:, :, :, r, i)
188 close (1)
189 end do
190 end do
191
192 do i = 1, nb
193 do r = 1, nnode
194 write (file_num, '(I0)') r + (i - 1)*nnode + sys_size
195 file_loc = trim(t_step_dir)//'/mv'//trim(file_num) &
196 //'.dat'
197 open (1, file=trim(file_loc), form='unformatted', &
198 status=status)
199 write (1) mv%sf(:, :, :, r, i)
200 close (1)
201 end do
202 end do
203 end if
204
205 gamma = gammas(1)
206 lit_gamma = gs_min(1)
207 pi_inf = pi_infs(1)
208 qv = qvs(1)
209
210 if (precision == 1) then
211 fmt = "(2F30.3)"
212 else
213 fmt = "(2F40.14)"
214 end if
215
216 write (t_step_dir, '(A,I0,A,I0)') trim(case_dir)//'/D'
217 file_loc = trim(t_step_dir)//'/.'
218
219 inquire (file=trim(file_loc), exist=file_exist)
220
221 if (.not. file_exist) call s_create_directory(trim(t_step_dir))
222
223 if (cfl_dt) t_step = n_start
224
225 !1D
226 if (n == 0 .and. p == 0) then
227 if (model_eqns == 2) then
228 do i = 1, sys_size
229 write (file_loc, '(A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir)//'/prim.', i, '.', proc_rank, '.', t_step, '.dat'
230
231 open (2, file=trim(file_loc))
232 do j = 0, m
233
234 if (chemistry) then
235 do c = 1, num_species
236 rhoyks(c) = q_cons_vf(chemxb + c - 1)%sf(j, 0, 0)
237 end do
238 end if
239
240 call s_convert_to_mixture_variables(q_cons_vf, j, 0, 0, rho, gamma, pi_inf, qv)
241
242 lit_gamma = 1._wp/gamma + 1._wp
243
244 if ((i >= chemxb) .and. (i <= chemxe)) then
245 write (2, fmt) x_cb(j), q_cons_vf(i)%sf(j, 0, 0)/rho
246 else if (((i >= cont_idx%beg) .and. (i <= cont_idx%end)) &
247 .or. &
248 ((i >= adv_idx%beg) .and. (i <= adv_idx%end)) &
249 .or. &
250 ((i >= chemxb) .and. (i <= chemxe)) &
251 ) then
252 write (2, fmt) x_cb(j), q_cons_vf(i)%sf(j, 0, 0)
253 else if (i == mom_idx%beg) then !u
254 write (2, fmt) x_cb(j), q_cons_vf(mom_idx%beg)%sf(j, 0, 0)/rho
255 else if (i == stress_idx%beg) then !tau_e
256 write (2, fmt) x_cb(j), q_cons_vf(stress_idx%beg)%sf(j, 0, 0)/rho
257 else if (i == e_idx) then !p
258 if (mhd) then
259 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, 0, 0)**2)
260 end if
261
262 call s_compute_pressure( &
263 q_cons_vf(e_idx)%sf(j, 0, 0), &
264 q_cons_vf(alf_idx)%sf(j, 0, 0), &
265 0.5_wp*(q_cons_vf(mom_idx%beg)%sf(j, 0, 0)**2._wp)/rho, &
266 pi_inf, gamma, rho, qv, rhoyks, pres, t, pres_mag=pres_mag)
267 write (2, fmt) x_cb(j), pres
268 else if (mhd) then
269 if (i == mom_idx%beg + 1) then ! v
270 write (2, fmt) x_cb(j), q_cons_vf(mom_idx%beg + 1)%sf(j, 0, 0)/rho
271 else if (i == mom_idx%beg + 2) then ! w
272 write (2, fmt) x_cb(j), q_cons_vf(mom_idx%beg + 2)%sf(j, 0, 0)/rho
273 else if (i == b_idx%beg) then ! By
274 write (2, fmt) x_cb(j), q_cons_vf(b_idx%beg)%sf(j, 0, 0)/rho
275 else if (i == b_idx%beg + 1) then ! Bz
276 write (2, fmt) x_cb(j), q_cons_vf(b_idx%beg + 1)%sf(j, 0, 0)/rho
277 end if
278 else if ((i >= bub_idx%beg) .and. (i <= bub_idx%end) .and. bubbles_euler) then
279
280 if (qbmm) then
281 nbub = q_cons_vf(bubxb)%sf(j, 0, 0)
282 else
283 if (adv_n) then
284 nbub = q_cons_vf(n_idx)%sf(j, 0, 0)
285 else
286 do k = 1, nb
287 nrtmp(k) = q_cons_vf(bub_idx%rs(k))%sf(j, 0, 0)
288 end do
289
290 call s_comp_n_from_cons(real(q_cons_vf(alf_idx)%sf(j, 0, 0), kind=wp), nrtmp, nbub, weight)
291 end if
292 end if
293 write (2, fmt) x_cb(j), q_cons_vf(i)%sf(j, 0, 0)/nbub
294 else if (i == n_idx .and. adv_n .and. bubbles_euler) then
295 write (2, fmt) x_cb(j), q_cons_vf(i)%sf(j, 0, 0)
296 else if (i == damage_idx) then
297 write (2, fmt) x_cb(j), q_cons_vf(i)%sf(j, 0, 0)
298 end if
299 end do
300 close (2)
301 end do
302 end if
303
304 do i = 1, sys_size
305 write (file_loc, '(A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir)//'/cons.', i, '.', proc_rank, '.', t_step, '.dat'
306
307 open (2, file=trim(file_loc))
308 do j = 0, m
309 write (2, fmt) x_cb(j), q_cons_vf(i)%sf(j, 0, 0)
310 end do
311 close (2)
312 end do
313
314 if (qbmm .and. .not. polytropic) then
315 do i = 1, nb
316 do r = 1, nnode
317 write (file_loc, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir)//'/pres.', i, '.', r, '.', proc_rank, '.', t_step, '.dat'
318
319 open (2, file=trim(file_loc))
320 do j = 0, m
321 write (2, fmt) x_cb(j), pb%sf(j, 0, 0, r, i)
322 end do
323 close (2)
324 end do
325 end do
326 do i = 1, nb
327 do r = 1, nnode
328 write (file_loc, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir)//'/mv.', i, '.', r, '.', proc_rank, '.', t_step, '.dat'
329
330 open (2, file=trim(file_loc))
331 do j = 0, m
332 write (2, fmt) x_cb(j), mv%sf(j, 0, 0, r, i)
333 end do
334 close (2)
335 end do
336 end do
337 end if
338 end if
339
340 if (precision == 1) then
341 fmt = "(3F30.7)"
342 else
343 fmt = "(3F40.14)"
344 end if
345
346 ! 2D
347 if ((n > 0) .and. (p == 0)) then
348 do i = 1, sys_size
349 write (file_loc, '(A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir)//'/cons.', i, '.', proc_rank, '.', t_step, '.dat'
350 open (2, file=trim(file_loc))
351 do j = 0, m
352 do k = 0, n
353 write (2, fmt) x_cb(j), y_cb(k), q_cons_vf(i)%sf(j, k, 0)
354 end do
355 write (2, *)
356 end do
357 close (2)
358 end do
359
360 if (qbmm .and. .not. polytropic) then
361 do i = 1, nb
362 do r = 1, nnode
363 write (file_loc, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir)//'/pres.', i, '.', r, '.', proc_rank, '.', t_step, '.dat'
364
365 open (2, file=trim(file_loc))
366 do j = 0, m
367 do k = 0, n
368 write (2, fmt) x_cb(j), y_cb(k), pb%sf(j, k, 0, r, i)
369 end do
370 end do
371 close (2)
372 end do
373 end do
374 do i = 1, nb
375 do r = 1, nnode
376 write (file_loc, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir)//'/mv.', i, '.', r, '.', proc_rank, '.', t_step, '.dat'
377
378 open (2, file=trim(file_loc))
379 do j = 0, m
380 do k = 0, n
381 write (2, fmt) x_cb(j), y_cb(k), mv%sf(j, k, 0, r, i)
382 end do
383 end do
384 close (2)
385 end do
386 end do
387 end if
388 end if
389
390 if (precision == 1) then
391 fmt = "(4F30.7)"
392 else
393 fmt = "(4F40.14)"
394 end if
395
396 ! 3D
397 if (p > 0) then
398 do i = 1, sys_size
399 write (file_loc, '(A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir)//'/cons.', i, '.', proc_rank, '.', t_step, '.dat'
400 open (2, file=trim(file_loc))
401 do j = 0, m
402 do k = 0, n
403 do l = 0, p
404 write (2, fmt) x_cb(j), y_cb(k), z_cb(l), q_cons_vf(i)%sf(j, k, l)
405 end do
406 write (2, *)
407 end do
408 write (2, *)
409 end do
410 close (2)
411 end do
412
413 if (qbmm .and. .not. polytropic) then
414 do i = 1, nb
415 do r = 1, nnode
416 write (file_loc, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir)//'/pres.', i, '.', r, '.', proc_rank, '.', t_step, '.dat'
417
418 open (2, file=trim(file_loc))
419 do j = 0, m
420 do k = 0, n
421 do l = 0, p
422 write (2, fmt) x_cb(j), y_cb(k), z_cb(l), pb%sf(j, k, l, r, i)
423 end do
424 end do
425 end do
426 close (2)
427 end do
428 end do
429 do i = 1, nb
430 do r = 1, nnode
431 write (file_loc, '(A,I0,A,I0,A,I2.2,A,I6.6,A)') trim(t_step_dir)//'/mv.', i, '.', r, '.', proc_rank, '.', t_step, '.dat'
432
433 open (2, file=trim(file_loc))
434 do j = 0, m
435 do k = 0, n
436 do l = 0, p
437 write (2, fmt) x_cb(j), y_cb(k), z_cb(l), mv%sf(j, k, l, r, i)
438 end do
439 end do
440 end do
441 close (2)
442 end do
443 end do
444 end if
445 end if
446
447 end subroutine s_write_serial_data_files
448
449 !> Writes grid and initial condition data files in parallel to the "0"
450 !! time-step directory in the local processor rank folder
451 !! @param q_cons_vf Conservative variables
452 !! @param q_prim_vf Primitive variables
453 !! @param bc_type Boundary condition types
454 impure subroutine s_write_parallel_data_files(q_cons_vf, q_prim_vf, bc_type)
455
456 ! Conservative variables
457 type(scalar_field), &
458 dimension(sys_size), &
459 intent(inout) :: q_cons_vf, q_prim_vf
460
461 type(integer_field), &
462 dimension(1:num_dims, -1:1), &
463 intent(in) :: bc_type
464
465#ifdef MFC_MPI
466
467 integer :: ifile, ierr, data_size
468 integer, dimension(MPI_STATUS_SIZE) :: status
469 integer(KIND=MPI_OFFSET_KIND) :: disp
470 integer(KIND=MPI_OFFSET_KIND) :: m_mok, n_mok, p_mok
471 integer(KIND=MPI_OFFSET_KIND) :: wp_mok, var_mok, str_mok
472 integer(KIND=MPI_OFFSET_KIND) :: nvars_mok
473 integer(KIND=MPI_OFFSET_KIND) :: mok
474
475 character(LEN=path_len + 2*name_len) :: file_loc
476 logical :: file_exist, dir_check
477
478 ! Generic loop iterators
479 integer :: i, j, k, l
480 real(wp) :: loc_violations, glb_violations
481
482 ! Downsample variables
483 integer :: m_ds, n_ds, p_ds
484 integer :: m_glb_ds, n_glb_ds, p_glb_ds
485 integer :: m_glb_save, n_glb_save, p_glb_save ! Size of array being saved
486
487 loc_violations = 0._wp
488
489 if (down_sample) then
490 if ((mod(m + 1, 3) > 0) .or. (mod(n + 1, 3) > 0) .or. (mod(p + 1, 3) > 0)) then
491 loc_violations = 1._wp
492 end if
493 call s_mpi_allreduce_sum(loc_violations, glb_violations)
494 if (proc_rank == 0 .and. nint(glb_violations) > 0) then
495 print *, "WARNING: Attempting to downsample data but there are"// &
496 "processors with local problem sizes that are not divisible by 3."
497 end if
500 m_ds, n_ds, p_ds, m_glb_ds, n_glb_ds, p_glb_ds)
501 end if
502
503 if (file_per_process) then
504 if (proc_rank == 0) then
505 file_loc = trim(case_dir)//'/restart_data/lustre_0'
506 call my_inquire(file_loc, dir_check)
507 if (dir_check .neqv. .true.) then
508 call s_create_directory(trim(file_loc))
509 end if
510 call s_create_directory(trim(file_loc))
511 end if
512 call s_mpi_barrier()
514
515 ! Initialize MPI data I/O
516 if (down_sample) then
517 call s_initialize_mpi_data_ds(q_cons_temp)
518 else
519 call s_initialize_mpi_data(q_cons_vf)
520 end if
521
522 ! Open the file to write all flow variables
523 if (cfl_dt) then
524 write (file_loc, '(I0,A,i7.7,A)') n_start, '_', proc_rank, '.dat'
525 else
526 write (file_loc, '(I0,A,i7.7,A)') t_step_start, '_', proc_rank, '.dat'
527 end if
528 file_loc = trim(restart_dir)//'/lustre_0'//trim(mpiiofs)//trim(file_loc)
529 inquire (file=trim(file_loc), exist=file_exist)
530 if (file_exist .and. proc_rank == 0) then
531 call mpi_file_delete(file_loc, mpi_info_int, ierr)
532 end if
533 if (file_exist) call mpi_file_delete(file_loc, mpi_info_int, ierr)
534 call mpi_file_open(mpi_comm_self, file_loc, ior(mpi_mode_wronly, mpi_mode_create), &
535 mpi_info_int, ifile, ierr)
536
537 if (down_sample) then
538 ! Size of local arrays
539 data_size = (m_ds + 3)*(n_ds + 3)*(p_ds + 3)
540 m_glb_save = m_glb_ds + 3
541 n_glb_save = n_glb_ds + 3
542 p_glb_save = p_glb_ds + 3
543 else
544 ! Size of local arrays
545 data_size = (m + 1)*(n + 1)*(p + 1)
546 m_glb_save = m_glb + 1
547 n_glb_save = n_glb + 1
548 p_glb_save = p_glb + 1
549 end if
550
551 ! Resize some integers so MPI can write even the biggest files
552 m_mok = int(m_glb_save, mpi_offset_kind)
553 n_mok = int(n_glb_save, mpi_offset_kind)
554 p_mok = int(p_glb_save, mpi_offset_kind)
555 wp_mok = int(8._wp, mpi_offset_kind)
556 mok = int(1._wp, mpi_offset_kind)
557 str_mok = int(name_len, mpi_offset_kind)
558 nvars_mok = int(sys_size, mpi_offset_kind)
559
560 ! Write the data for each variable
561 if (bubbles_euler) then
562 do i = 1, sys_size! adv_idx%end
563 var_mok = int(i, mpi_offset_kind)
564
565 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, &
566 mpi_io_p, status, ierr)
567 end do
568 !Additional variables pb and mv for non-polytropic qbmm
569 if (qbmm .and. .not. polytropic) then
570 do i = sys_size + 1, sys_size + 2*nb*nnode
571 var_mok = int(i, mpi_offset_kind)
572
573 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, &
574 mpi_io_p, status, ierr)
575 end do
576 end if
577 else
578 if (down_sample) then
579 do i = 1, sys_size
580 var_mok = int(i, mpi_offset_kind)
581
582 call mpi_file_write_all(ifile, q_cons_temp(i)%sf, data_size*mpi_io_type, &
583 mpi_io_p, status, ierr)
584 end do
585 else
586 do i = 1, sys_size
587 var_mok = int(i, mpi_offset_kind)
588
589 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, &
590 mpi_io_p, status, ierr)
591 end do
592 end if
593 end if
594
595 call mpi_file_close(ifile, ierr)
596
597 else
598 call s_initialize_mpi_data(q_cons_vf)
599
600 ! Open the file to write all flow variables
601 if (cfl_dt) then
602 write (file_loc, '(I0,A)') n_start, '.dat'
603 else
604 write (file_loc, '(I0,A)') t_step_start, '.dat'
605 end if
606 file_loc = trim(restart_dir)//trim(mpiiofs)//trim(file_loc)
607 inquire (file=trim(file_loc), exist=file_exist)
608 if (file_exist .and. proc_rank == 0) then
609 call mpi_file_delete(file_loc, mpi_info_int, ierr)
610 end if
611 call mpi_file_open(mpi_comm_world, file_loc, ior(mpi_mode_wronly, mpi_mode_create), &
612 mpi_info_int, ifile, ierr)
613
614 ! Size of local arrays
615 data_size = (m + 1)*(n + 1)*(p + 1)
616
617 ! Resize some integers so MPI can write even the biggest files
618 m_mok = int(m_glb + 1, mpi_offset_kind)
619 n_mok = int(n_glb + 1, mpi_offset_kind)
620 p_mok = int(p_glb + 1, mpi_offset_kind)
621 wp_mok = int(8._wp, mpi_offset_kind)
622 mok = int(1._wp, mpi_offset_kind)
623 str_mok = int(name_len, mpi_offset_kind)
624 nvars_mok = int(sys_size, mpi_offset_kind)
625
626 ! Write the data for each variable
627 if (bubbles_euler) then
628 do i = 1, sys_size! adv_idx%end
629 var_mok = int(i, mpi_offset_kind)
630
631 ! Initial displacement to skip at beginning of file
632 disp = m_mok*max(mok, n_mok)*max(mok, p_mok)*wp_mok*(var_mok - 1)
633
634 call mpi_file_set_view(ifile, disp, mpi_io_p, mpi_io_data%view(i), &
635 'native', mpi_info_int, ierr)
636 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, &
637 mpi_io_p, status, ierr)
638 end do
639 !Additional variables pb and mv for non-polytropic qbmm
640 if (qbmm .and. .not. polytropic) then
641 do i = sys_size + 1, sys_size + 2*nb*nnode
642 var_mok = int(i, mpi_offset_kind)
643
644 ! Initial displacement to skip at beginning of file
645 disp = m_mok*max(mok, n_mok)*max(mok, p_mok)*wp_mok*(var_mok - 1)
646
647 call mpi_file_set_view(ifile, disp, mpi_io_p, mpi_io_data%view(i), &
648 'native', mpi_info_int, ierr)
649 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, &
650 mpi_io_p, status, ierr)
651 end do
652 end if
653 else
654 do i = 1, sys_size !TODO: check if this is right
655 ! do i = 1, adv_idx%end
656 var_mok = int(i, mpi_offset_kind)
657
658 ! Initial displacement to skip at beginning of file
659 disp = m_mok*max(mok, n_mok)*max(mok, p_mok)*wp_mok*(var_mok - 1)
660
661 call mpi_file_set_view(ifile, disp, mpi_io_p, mpi_io_data%view(i), &
662 'native', mpi_info_int, ierr)
663 call mpi_file_write_all(ifile, mpi_io_data%var(i)%sf, data_size*mpi_io_type, &
664 mpi_io_p, status, ierr)
665 end do
666
667 end if
668
669 call mpi_file_close(ifile, ierr)
670 end if
671#endif
672
673 if (bc_io) then
674 if (igr) then
676 else
677 call s_write_parallel_boundary_condition_files(q_prim_vf, bc_type)
678 end if
679 end if
680
681 end subroutine s_write_parallel_data_files
682
683 !> Computation of parameters, allocation procedures, and/or
684 !! any other tasks needed to properly setup the module
686 ! Generic string used to store the address of a particular file
687 character(LEN=len_trim(case_dir) + 2*name_len) :: file_loc
688 character(len=15) :: temp
689 character(LEN=1), dimension(3), parameter :: coord = (/'x', 'y', 'z'/)
690
691 ! Generic logical used to check the existence of directories
692 logical :: dir_check
693 integer :: i
694
695 integer :: m_ds, n_ds, p_ds !< down sample dimensions
696
697 if (parallel_io .neqv. .true.) then
698 ! Setting the address of the time-step directory
699 write (t_step_dir, '(A,I0,A)') '/p_all/p', proc_rank, '/0'
700 t_step_dir = trim(case_dir)//trim(t_step_dir)
701
702 ! Checking the existence of the time-step directory, removing it, if
703 ! it exists, and creating a new copy. Note that if preexisting grid
704 ! and/or initial condition data are to be read in from the very same
705 ! location, then the above described steps are not executed here but
706 ! rather in the module m_start_up.f90.
707 if (old_grid .neqv. .true.) then
708
709 file_loc = trim(t_step_dir)//'/'
710
711 call my_inquire(file_loc, dir_check)
712
713 if (dir_check) call s_delete_directory(trim(t_step_dir))
714
716
717 end if
718
720 else
721 write (restart_dir, '(A)') '/restart_data'
722 restart_dir = trim(case_dir)//trim(restart_dir)
723
724 if ((old_grid .neqv. .true.) .and. (proc_rank == 0)) then
725
726 file_loc = trim(restart_dir)//'/'
727 call my_inquire(file_loc, dir_check)
728
729 if (dir_check) call s_delete_directory(trim(restart_dir))
731 end if
732
733 call s_mpi_barrier()
734
736
737 end if
738
739 open (1, file='indices.dat', status='unknown')
740
741 write (1, '(A)') "Warning: The creation of file is currently experimental."
742 write (1, '(A)') "This file may contain errors and not support all features."
743
744 write (1, '(A3,A20,A20)') "#", "Conservative", "Primitive"
745 write (1, '(A)') " "
746 do i = contxb, contxe
747 write (temp, '(I0)') i - contxb + 1
748 write (1, '(I3,A20,A20)') i, "\alpha_{"//trim(temp)//"} \rho_{"//trim(temp)//"}", "\alpha_{"//trim(temp)//"} \rho"
749 end do
750 do i = momxb, momxe
751 write (1, '(I3,A20,A20)') i, "\rho u_"//coord(i - momxb + 1), "u_"//coord(i - momxb + 1)
752 end do
753 do i = e_idx, e_idx
754 write (1, '(I3,A20,A20)') i, "\rho U", "p"
755 end do
756 do i = advxb, advxe
757 write (temp, '(I0)') i - contxb + 1
758 write (1, '(I3,A20,A20)') i, "\alpha_{"//trim(temp)//"}", "\alpha_{"//trim(temp)//"}"
759 end do
760 if (chemistry) then
761 do i = 1, num_species
762 write (1, '(I3,A20,A20)') chemxb + i - 1, "Y_{"//trim(species_names(i))//"} \rho", "Y_{"//trim(species_names(i))//"}"
763 end do
764 end if
765
766 write (1, '(A)') ""
767 if (momxb /= 0) write (1, '("[",I2,",",I2,"]",A)') momxb, momxe, " Momentum"
768 if (e_idx /= 0) write (1, '("[",I2,",",I2,"]",A)') e_idx, e_idx, " Energy/Pressure"
769 if (advxb /= 0) write (1, '("[",I2,",",I2,"]",A)') advxb, advxe, " Advection"
770 if (contxb /= 0) write (1, '("[",I2,",",I2,"]",A)') contxb, contxe, " Continuity"
771 if (bubxb /= 0) write (1, '("[",I2,",",I2,"]",A)') bubxb, bubxe, " Bubbles_euler"
772 if (strxb /= 0) write (1, '("[",I2,",",I2,"]",A)') strxb, strxe, " Stress"
773 if (intxb /= 0) write (1, '("[",I2,",",I2,"]",A)') intxb, intxe, " Internal Energies"
774 if (chemxb /= 0) write (1, '("[",I2,",",I2,"]",A)') chemxb, chemxe, " Chemistry"
775
776 close (1)
777
778 if (down_sample) then
779 m_ds = int((m + 1)/3) - 1
780 n_ds = int((n + 1)/3) - 1
781 p_ds = int((p + 1)/3) - 1
782
783 allocate (q_cons_temp(1:sys_size))
784 do i = 1, sys_size
785 allocate (q_cons_temp(i)%sf(-1:m_ds + 1, -1:n_ds + 1, -1:p_ds + 1))
786 end do
787 end if
788
790
791 !> Resets s_write_data_files pointer
793
794 integer :: i
795
796 s_write_data_files => null()
797
798 if (down_sample) then
799 do i = 1, sys_size
800 deallocate (q_cons_temp(i)%sf)
801 end do
802 deallocate (q_cons_temp)
803 end if
804
805 end subroutine s_finalize_data_output_module
806
807end 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)
Writes boundary condition type and buffer data to serial (unformatted) restart files.
subroutine, public s_write_parallel_boundary_condition_files(q_prim_vf, bc_type)
Writes 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)
The purpose of this procedure is to populate the buffers of the primitive variables,...
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 deletes 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)
Creates 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)
Introduces 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.
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
Locations of cell-boundaries (cb) in x-, y- and z-directions, respectively.
integer num_dims
Number of spatial dimensions.
real(wp), dimension(:), allocatable x_cb
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)
Downsamples 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)
Computes 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)
This procedure conditionally calculates the appropriate pressure.
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).