1# 1 "/home/runner/work/MFC/MFC/src/pre_process/m_data_output.fpp"
15 use m_mpi_proxy !< message passing interface (mpi) module proxy
18 use mpi !< message passing interface (mpi) module
33 use m_thermochem,
only: species_names
59 dimension(sys_size), &
63 dimension(1:num_dims, -1:1), &
69 character(LEN=path_len + 2*name_len),
private ::
t_step_dir
86 dimension(sys_size), &
91 dimension(1:num_dims, -1:1), &
96 character(LEN=15) :: fmt
97 character(LEN=3) :: status
100 int(floor(log10(real(
sys_size, wp)))) + 1) :: file_num
104 character(LEN=len_trim(t_step_dir) + name_len) :: file_loc
107 integer :: i,
j,
k,
l, r, c
110 real(wp),
dimension(nb) :: nrtmp
112 real(wp) :: gamma, lit_gamma, pi_inf, qv
116 real(wp) :: rhoyks(1:num_species)
144 open (1, file=trim(file_loc), form=
'unformatted', status=status)
152 open (1, file=trim(file_loc), form=
'unformatted', &
160 open (1, file=trim(file_loc), form=
'unformatted', &
169 write (file_num,
'(I0)') i
170 file_loc = trim(
t_step_dir)//
'/q_cons_vf'//trim(file_num) &
172 open (1, file=trim(file_loc), form=
'unformatted', &
182 write (file_num,
'(I0)') r + (i - 1)*nnode +
sys_size
183 file_loc = trim(
t_step_dir)//
'/pb'//trim(file_num) &
185 open (1, file=trim(file_loc), form=
'unformatted', &
187 write (1)
pb%sf(:, :, :, r, i)
194 write (file_num,
'(I0)') r + (i - 1)*nnode +
sys_size
195 file_loc = trim(
t_step_dir)//
'/mv'//trim(file_num) &
197 open (1, file=trim(file_loc), form=
'unformatted', &
199 write (1)
mv%sf(:, :, :, r, i)
219 inquire (file=trim(file_loc), exist=file_exist)
226 if (
n == 0 .and.
p == 0)
then
229 write (file_loc,
'(A,I0,A,I2.2,A,I6.6,A)') trim(
t_step_dir)//
'/prim.', i,
'.',
proc_rank,
'.', t_step,
'.dat'
231 open (2, file=trim(file_loc))
235 do c = 1, num_species
242 lit_gamma = 1._wp/gamma + 1._wp
253 else if (i ==
mom_idx%beg)
then
257 else if (i ==
e_idx)
then
266 pi_inf, gamma, rho, qv, rhoyks, pres, t, pres_mag=pres_mag)
267 write (2, fmt)
x_cb(
j), pres
271 else if (i ==
mom_idx%beg + 2)
then
273 else if (i ==
b_idx%beg)
then
275 else if (i ==
b_idx%beg + 1)
then
305 write (file_loc,
'(A,I0,A,I2.2,A,I6.6,A)') trim(
t_step_dir)//
'/cons.', i,
'.',
proc_rank,
'.', t_step,
'.dat'
307 open (2, file=trim(file_loc))
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'
319 open (2, file=trim(file_loc))
321 write (2, fmt)
x_cb(
j),
pb%sf(
j, 0, 0, r, i)
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'
330 open (2, file=trim(file_loc))
332 write (2, fmt)
x_cb(
j),
mv%sf(
j, 0, 0, r, i)
347 if ((
n > 0) .and. (
p == 0))
then
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))
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'
365 open (2, file=trim(file_loc))
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'
378 open (2, file=trim(file_loc))
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))
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'
418 open (2, file=trim(file_loc))
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'
433 open (2, file=trim(file_loc))
458 dimension(sys_size), &
462 dimension(1:num_dims, -1:1), &
463 intent(in) :: bc_type
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
475 character(LEN=path_len + 2*name_len) :: file_loc
476 logical :: file_exist, dir_check
479 integer :: i,
j,
k,
l
480 real(wp) :: loc_violations, glb_violations
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
488 if ((mod(
m + 1, 3) > 0) .or. (mod(
n + 1, 3) > 0) .or. (mod(
p + 1, 3) > 0))
then
489 loc_violations = 1._wp
491 call s_mpi_allreduce_sum(loc_violations, glb_violations)
492 if (
proc_rank == 0 .and. nint(glb_violations) > 0)
then
493 print *,
"WARNING: Attempting to downsample data but there are"// &
494 "processors with local problem sizes that are not divisible by 3."
498 m_ds, n_ds, p_ds, m_glb_ds, n_glb_ds, p_glb_ds)
503 file_loc = trim(
case_dir)//
'/restart_data/lustre_0'
505 if (dir_check .neqv. .true.)
then
527 inquire (file=trim(file_loc), exist=file_exist)
528 if (file_exist .and.
proc_rank == 0)
then
531 if (file_exist)
call mpi_file_delete(file_loc,
mpi_info_int, ierr)
532 call mpi_file_open(mpi_comm_self, file_loc, ior(mpi_mode_wronly, mpi_mode_create), &
537 data_size = (m_ds + 3)*(n_ds + 3)*(p_ds + 3)
538 m_glb_save = m_glb_ds + 3
539 n_glb_save = n_glb_ds + 3
540 p_glb_save = p_glb_ds + 3
543 data_size = (
m + 1)*(
n + 1)*(
p + 1)
544 m_glb_save =
m_glb + 1
545 n_glb_save =
n_glb + 1
546 p_glb_save =
p_glb + 1
550 m_mok = int(m_glb_save, mpi_offset_kind)
551 n_mok = int(n_glb_save, mpi_offset_kind)
552 p_mok = int(p_glb_save, mpi_offset_kind)
553 wp_mok = int(8._wp, mpi_offset_kind)
554 mok = int(1._wp, mpi_offset_kind)
555 str_mok = int(name_len, mpi_offset_kind)
556 nvars_mok = int(
sys_size, mpi_offset_kind)
561 var_mok = int(i, mpi_offset_kind)
563 call mpi_file_write_all(ifile,
mpi_io_data%var(i)%sf, data_size*mpi_io_type, &
564 mpi_io_p, status, ierr)
569 var_mok = int(i, mpi_offset_kind)
571 call mpi_file_write_all(ifile,
mpi_io_data%var(i)%sf, data_size*mpi_io_type, &
572 mpi_io_p, status, ierr)
578 var_mok = int(i, mpi_offset_kind)
580 call mpi_file_write_all(ifile,
q_cons_temp(i)%sf, data_size*mpi_io_type, &
581 mpi_io_p, status, ierr)
585 var_mok = int(i, mpi_offset_kind)
587 call mpi_file_write_all(ifile,
mpi_io_data%var(i)%sf, data_size*mpi_io_type, &
588 mpi_io_p, status, ierr)
593 call mpi_file_close(ifile, ierr)
600 write (file_loc,
'(I0,A)')
n_start,
'.dat'
605 inquire (file=trim(file_loc), exist=file_exist)
606 if (file_exist .and.
proc_rank == 0)
then
609 call mpi_file_open(mpi_comm_world, file_loc, ior(mpi_mode_wronly, mpi_mode_create), &
613 data_size = (
m + 1)*(
n + 1)*(
p + 1)
616 m_mok = int(
m_glb + 1, mpi_offset_kind)
617 n_mok = int(
n_glb + 1, mpi_offset_kind)
618 p_mok = int(
p_glb + 1, mpi_offset_kind)
619 wp_mok = int(8._wp, mpi_offset_kind)
620 mok = int(1._wp, mpi_offset_kind)
621 str_mok = int(name_len, mpi_offset_kind)
622 nvars_mok = int(
sys_size, mpi_offset_kind)
627 var_mok = int(i, mpi_offset_kind)
630 disp = m_mok*max(mok, n_mok)*max(mok, p_mok)*wp_mok*(var_mok - 1)
632 call mpi_file_set_view(ifile, disp, mpi_io_p,
mpi_io_data%view(i), &
634 call mpi_file_write_all(ifile,
mpi_io_data%var(i)%sf, data_size*mpi_io_type, &
635 mpi_io_p, status, ierr)
640 var_mok = int(i, mpi_offset_kind)
643 disp = m_mok*max(mok, n_mok)*max(mok, p_mok)*wp_mok*(var_mok - 1)
645 call mpi_file_set_view(ifile, disp, mpi_io_p,
mpi_io_data%view(i), &
647 call mpi_file_write_all(ifile,
mpi_io_data%var(i)%sf, data_size*mpi_io_type, &
648 mpi_io_p, status, ierr)
654 var_mok = int(i, mpi_offset_kind)
657 disp = m_mok*max(mok, n_mok)*max(mok, p_mok)*wp_mok*(var_mok - 1)
659 call mpi_file_set_view(ifile, disp, mpi_io_p,
mpi_io_data%view(i), &
661 call mpi_file_write_all(ifile,
mpi_io_data%var(i)%sf, data_size*mpi_io_type, &
662 mpi_io_p, status, ierr)
667 call mpi_file_close(ifile, ierr)
685 character(LEN=len_trim(case_dir) + 2*name_len) :: file_loc
686 character(len=15) :: temp
687 character(LEN=1),
dimension(3),
parameter :: coord = (/
'x',
'y',
'z'/)
693 integer :: m_ds, n_ds, p_ds
737 open (1, file=
'indices.dat', status=
'unknown')
739 write (1,
'(A)')
"Warning: The creation of file is currently experimental."
740 write (1,
'(A)')
"This file may contain errors and not support all features."
742 write (1,
'(A3,A20,A20)')
"#",
"Conservative",
"Primitive"
745 write (temp,
'(I0)') i -
contxb + 1
746 write (1,
'(I3,A20,A20)') i,
"\alpha_{"//trim(temp)//
"} \rho_{"//trim(temp)//
"}",
"\alpha_{"//trim(temp)//
"} \rho"
749 write (1,
'(I3,A20,A20)') i,
"\rho u_"//coord(i -
momxb + 1),
"u_"//coord(i -
momxb + 1)
752 write (1,
'(I3,A20,A20)') i,
"\rho U",
"p"
755 write (temp,
'(I0)') i -
contxb + 1
756 write (1,
'(I3,A20,A20)') i,
"\alpha_{"//trim(temp)//
"}",
"\alpha_{"//trim(temp)//
"}"
759 do i = 1, num_species
760 write (1,
'(I3,A20,A20)')
chemxb + i - 1,
"Y_{"//trim(species_names(i))//
"} \rho",
"Y_{"//trim(species_names(i))//
"}"
765 if (
momxb /= 0)
write (1,
'("[",I2,",",I2,"]",A)')
momxb,
momxe,
" Momentum"
766 if (
e_idx /= 0)
write (1,
'("[",I2,",",I2,"]",A)')
e_idx,
e_idx,
" Energy/Pressure"
767 if (
advxb /= 0)
write (1,
'("[",I2,",",I2,"]",A)')
advxb,
advxe,
" Advection"
769 if (
bubxb /= 0)
write (1,
'("[",I2,",",I2,"]",A)')
bubxb,
bubxe,
" Bubbles_euler"
770 if (
strxb /= 0)
write (1,
'("[",I2,",",I2,"]",A)')
strxb,
strxe,
" Stress"
771 if (
intxb /= 0)
write (1,
'("[",I2,",",I2,"]",A)')
intxb,
intxe,
" Internal Energies"
777 m_ds = int((
m + 1)/3) - 1
778 n_ds = int((
n + 1)/3) - 1
779 p_ds = int((
p + 1)/3) - 1
783 allocate (
q_cons_temp(i)%sf(-1:m_ds + 1, -1:n_ds + 1, -1:p_ds + 1))
Interface for the conservative data.
type(scalar_field), dimension(sys_size), intent(inout) q_cons_vf
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).