MFC
Exascale flow solver
Loading...
Searching...
No Matches
m_start_up.fpp.f90
Go to the documentation of this file.
1# 1 "/home/runner/work/MFC/MFC/src/pre_process/m_start_up.fpp"
2!>
3!! @file
4!! @brief Contains module m_start_up
5
6!> @brief Reads and validates user inputs, loads existing grid/IC data, and initializes pre-process modules
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_mpi_proxy !< message passing interface (mpi) module proxy
14
15 use m_mpi_common
16
17 use m_variables_conversion !< subroutines to change the state variables from
18 !! one form to another
19
20 use m_grid !< procedures to generate (non-)uniform grids
21
22 use m_initial_condition !< procedures to generate initial condition
23
24 use m_data_output !< procedures to write the grid data and the
25 !! conservative variables to files
26
27 use m_compile_specific !< compile-specific procedures
28
30
32
33 use m_phase_change !< phase-change module
34
35 use m_helper_basic !< functions to compare floating point numbers
36
37 use m_helper
38
39#ifdef MFC_MPI
40 use mpi !< message passing interface (mpi) module
41#endif
42
44
46
47 use m_helper
48
50
51 use m_checker
52
54
56
57 implicit none
58
59 private;
60 public :: s_read_input_file, &
74
75 abstract interface
76
77 !> @brief Abstract interface for reading grid data files in serial or parallel.
79
81
82 !> @brief Abstract interface for reading initial condition data files in serial or parallel.
83 !! @param q_cons_vf Conservative variables
84 impure subroutine s_read_abstract_ic_data_files(q_cons_vf_in)
85
87
88 type(scalar_field), &
89 dimension(sys_size), &
90 intent(inout) :: q_cons_vf_in
91
93
94 end interface
95
96 character(LEN=path_len + name_len) :: proc_rank_dir !<
97 !! Location of the folder associated with the rank of the local processor
98
99 character(LEN=path_len + 2*name_len), private :: t_step_dir !<
100 !! Possible location of time-step folder containing preexisting grid and/or
101 !! conservative variables data to be used as starting point for pre-process
102
105
106contains
107
108 !> Reads the configuration file pre_process.inp, in order to
109 !! populate the parameters in module m_global_parameters.f90
110 !! with the user provided inputs
111 impure subroutine s_read_input_file
112
113 character(LEN=name_len) :: file_loc !<
114 !! Generic string used to store the address of a particular file
115
116 logical :: file_check !<
117 !! Generic logical used for the purpose of asserting whether a file
118 !! is or is not present in the designated location
119
120 integer :: iostatus
121 !! Integer to check iostat of file read
122
123 character(len=1000) :: line
124
125 ! Namelist for all of the parameters to be inputted by the user
126 namelist /user_inputs/ case_dir, old_grid, old_ic, &
129 a_z, x_a, y_a, z_a, x_b, y_b, z_b, &
151
152 ! Inquiring the status of the pre_process.inp file
153 file_loc = 'pre_process.inp'
154 inquire (file=trim(file_loc), exist=file_check)
155
156 ! Checking whether the input file is there. If it is, the input file
157 ! is read. If not, the program is terminated.
158 if (file_check) then
159 open (1, file=trim(file_loc), form='formatted', &
160 status='old', action='read')
161 read (1, nml=user_inputs, iostat=iostatus)
162 if (iostatus /= 0) then
163 backspace(1)
164 read (1, fmt='(A)') line
165 print *, 'Invalid line in namelist: '//trim(line)
166 call s_mpi_abort('Invalid line in pre_process.inp. It is '// &
167 'likely due to a datatype mismatch. Exiting.')
168 end if
169 close (1)
170
172
173 ! Store m,n,p into global m,n,p
174 m_glb = m
175 n_glb = n
176 p_glb = p
177
178 nglobal = int(m_glb + 1, kind=8)*int(n_glb + 1, kind=8)*int(p_glb + 1, kind=8)
179
180 if (cfl_adap_dt .or. cfl_const_dt) cfl_dt = .true.
181
182 if (any((/bc_x%beg, bc_x%end, bc_y%beg, bc_y%end, bc_z%beg, bc_z%end/) == bc_dirichlet) .or. &
183 num_bc_patches > 0) then
184 bc_io = .true.
185 end if
186
187 else
188 call s_mpi_abort('File pre_process.inp is missing. Exiting.')
189 end if
190
191 end subroutine s_read_input_file
192
193 !> Checking that the user inputs make sense, i.e. that the
194 !! individual choices are compatible with the code's options
195 !! and that the combination of these choices results into a
196 !! valid configuration for the pre-process
197 impure subroutine s_check_input_file
198
199 character(LEN=len_trim(case_dir)) :: file_loc !<
200 !! Generic string used to store the address of a particular file
201
202 logical :: dir_check !<
203 !! Logical variable used to test the existence of folders
204
205 ! Checking the existence of the case folder
206 case_dir = adjustl(case_dir)
207
208 file_loc = trim(case_dir)//'/.'
209
210 call my_inquire(file_loc, dir_check)
211
212 if (dir_check .neqv. .true.) then
213 print '(A)', 'WARNING: Ensure that compiler flags/choices in Makefiles match your compiler! '
214 print '(A)', 'WARNING: Ensure that preprocessor flags are enabled! '
215 call s_mpi_abort('Unsupported choice for the value of case_dir.'// &
216 'Exiting.')
217 end if
218
220 call s_check_inputs()
221
222 ! Check all the patch properties
223 call s_check_patches()
224
225 if (ib) call s_check_ib_patches()
226
227 end subroutine s_check_input_file
228
229 !> The goal of this subroutine is to read in any preexisting
230 !! grid data as well as based on the imported grid, complete
231 !! the necessary global computational domain parameters.
233
234 ! Generic string used to store the address of a particular file
235 character(LEN=len_trim(case_dir) + 3*name_len) :: file_loc
236
237 ! Logical variable used to test the existence of folders
238 logical :: dir_check
239
240 ! Generic logical used for the purpose of asserting whether a file
241 ! is or is not present in the designated location
242 logical :: file_check
243
244 ! Setting address of the local processor rank and time-step directory
245 write (proc_rank_dir, '(A,I0)') '/p_all/p', proc_rank
247
248 write (t_step_dir, '(A,I0)') '/', t_step_start
249 t_step_dir = trim(proc_rank_dir)//trim(t_step_dir)
250
251 ! Inquiring as to the existence of the time-step directory
252 file_loc = trim(t_step_dir)//'/.'
253 call my_inquire(file_loc, dir_check)
254
255 ! If the time-step directory is missing, the pre-process exits
256 if (dir_check .neqv. .true.) then
257 call s_mpi_abort('Time-step folder '//trim(t_step_dir)// &
258 ' is missing. Exiting.')
259 end if
260
261 ! Reading the Grid Data File for the x-direction
262
263 ! Checking whether x_cb.dat exists
264 file_loc = trim(t_step_dir)//'/x_cb.dat'
265 inquire (file=trim(file_loc), exist=file_check)
266
267 ! If it exists, x_cb.dat is read
268 if (file_check) then
269 open (1, file=trim(file_loc), form='unformatted', &
270 status='old', action='read')
271 read (1) x_cb(-1:m)
272 close (1)
273 else
274 call s_mpi_abort('File x_cb.dat is missing in '// &
275 trim(t_step_dir)//'. Exiting.')
276 end if
277
278 ! Computing cell-center locations
279 x_cc(0:m) = (x_cb(0:m) + x_cb(-1:(m - 1)))/2._wp
280
281 ! Computing minimum cell-width
282 dx = minval(x_cb(0:m) - x_cb(-1:m - 1))
283 if (num_procs > 1) call s_mpi_reduce_min(dx)
284
285 ! Setting locations of domain bounds
286 x_domain%beg = x_cb(-1)
287 x_domain%end = x_cb(m)
288
289 ! Reading the Grid Data File for the y-direction
290
291 if (n > 0) then
292
293 ! Checking whether y_cb.dat exists
294 file_loc = trim(t_step_dir)//'/y_cb.dat'
295 inquire (file=trim(file_loc), exist=file_check)
296
297 ! If it exists, y_cb.dat is read
298 if (file_check) then
299 open (1, file=trim(file_loc), form='unformatted', &
300 status='old', action='read')
301 read (1) y_cb(-1:n)
302 close (1)
303 else
304 call s_mpi_abort('File y_cb.dat is missing in '// &
305 trim(t_step_dir)//'. Exiting.')
306 end if
307
308 ! Computing cell-center locations
309 y_cc(0:n) = (y_cb(0:n) + y_cb(-1:(n - 1)))/2._wp
310
311 ! Computing minimum cell-width
312 dy = minval(y_cb(0:n) - y_cb(-1:n - 1))
313 if (num_procs > 1) call s_mpi_reduce_min(dy)
314
315 ! Setting locations of domain bounds
316 y_domain%beg = y_cb(-1)
317 y_domain%end = y_cb(n)
318
319 ! Reading the Grid Data File for the z-direction
320 if (p > 0) then
321
322 ! Checking whether z_cb.dat exists
323 file_loc = trim(t_step_dir)//'/z_cb.dat'
324 inquire (file=trim(file_loc), exist=file_check)
325
326 ! If it exists, z_cb.dat is read
327 if (file_check) then
328 open (1, file=trim(file_loc), form='unformatted', &
329 status='old', action='read')
330 read (1) z_cb(-1:p)
331 close (1)
332 else
333 call s_mpi_abort('File z_cb.dat is missing in '// &
334 trim(t_step_dir)//'. Exiting.')
335 end if
336
337 ! Computing cell-center locations
338 z_cc(0:p) = (z_cb(0:p) + z_cb(-1:(p - 1)))/2._wp
339
340 ! Computing minimum cell-width
341 dz = minval(z_cb(0:p) - z_cb(-1:p - 1))
342 if (num_procs > 1) call s_mpi_reduce_min(dz)
343
344 ! Setting locations of domain bounds
345 z_domain%beg = z_cb(-1)
346 z_domain%end = z_cb(p)
347
348 end if
349
350 end if
351
352 ! If only the preexisting grid data files are read in and there will
353 ! not be any preexisting initial condition data files imported, then
354 ! the directory associated with the rank of the local processor may
355 ! be cleaned to make room for the new pre-process data. In addition,
356 ! the time-step directory that will contain the new grid and initial
357 ! condition data are also generated.
358 if (old_ic .neqv. .true.) then
359 call s_delete_directory(trim(proc_rank_dir)//'/*')
360 call s_create_directory(trim(proc_rank_dir)//'/0')
361 end if
362
363 end subroutine s_read_serial_grid_data_files
364
365 !> Cell-boundary data are checked for consistency by looking
366 !! at the (non-)uniform cell-width distributions for all the
367 !! active coordinate directions and making sure that all of
368 !! the cell-widths are positively valued
369 impure subroutine s_check_grid_data_files
370
371 ! Cell-boundary Data Consistency Check in x-direction
372
373 if (any(x_cb(0:m) - x_cb(-1:m - 1) <= 0._wp)) then
374 call s_mpi_abort('x_cb.dat in '//trim(t_step_dir)// &
375 ' contains non-positive cell-spacings. Exiting.')
376 end if
377
378 ! Cell-boundary Data Consistency Check in y-direction
379
380 if (n > 0) then
381
382 if (any(y_cb(0:n) - y_cb(-1:n - 1) <= 0._wp)) then
383 call s_mpi_abort('y_cb.dat in '//trim(t_step_dir)// &
384 ' contains non-positive cell-spacings. '// &
385 'Exiting.')
386 end if
387
388 ! Cell-boundary Data Consistency Check in z-direction
389
390 if (p > 0) then
391
392 if (any(z_cb(0:p) - z_cb(-1:p - 1) <= 0._wp)) then
393 call s_mpi_abort('z_cb.dat in '//trim(t_step_dir)// &
394 ' contains non-positive cell-spacings'// &
395 ' .Exiting.')
396 end if
397
398 end if
399
400 end if
401
402 end subroutine s_check_grid_data_files
403
404 !> The goal of this subroutine is to read in any preexisting
405 !! initial condition data files so that they may be used by
406 !! the pre-process as a starting point in the creation of an
407 !! all new initial condition.
408 !! @param q_cons_vf_in Conservative variables
409 impure subroutine s_read_serial_ic_data_files(q_cons_vf_in)
410
411 type(scalar_field), &
412 dimension(sys_size), &
413 intent(inout) :: q_cons_vf_in
414
415 character(LEN=len_trim(case_dir) + 3*name_len) :: file_loc !<
416 ! Generic string used to store the address of a particular file
417
418 character(len= &
419 int(floor(log10(real(sys_size, wp)))) + 1) :: file_num !<
420 !! Used to store the variable position, in character form, of the
421 !! currently manipulated conservative variable file
422
423 logical :: file_check !<
424 !! Generic logical used for the purpose of asserting whether a file
425 !! is or is not present in the designated location
426
427 integer :: i, r !< Generic loop iterator
428
429 ! Reading the Conservative Variables Data Files
430 do i = 1, sys_size
431
432 ! Checking whether data file associated with variable position
433 ! of the currently manipulated conservative variable exists
434 write (file_num, '(I0)') i
435 file_loc = trim(t_step_dir)//'/q_cons_vf'// &
436 trim(file_num)//'.dat'
437 inquire (file=trim(file_loc), exist=file_check)
438
439 ! If it exists, the data file is read
440 if (file_check) then
441 open (1, file=trim(file_loc), form='unformatted', &
442 status='old', action='read')
443 read (1) q_cons_vf_in(i)%sf
444 close (1)
445 else
446 call s_mpi_abort('File q_cons_vf'//trim(file_num)// &
447 '.dat is missing in '//trim(t_step_dir)// &
448 '. Exiting.')
449 end if
450
451 end do
452
453 !Read bubble variables pb and mv for non-polytropic qbmm
454 if (qbmm .and. .not. polytropic) then
455 do i = 1, nb
456 do r = 1, nnode
457 ! Checking whether data file associated with variable position
458 ! of the currently manipulated bubble variable exists
459 write (file_num, '(I0)') sys_size + r + (i - 1)*nnode
460 file_loc = trim(t_step_dir)//'/pb'// &
461 trim(file_num)//'.dat'
462 inquire (file=trim(file_loc), exist=file_check)
463
464 ! If it exists, the data file is read
465 if (file_check) then
466 open (1, file=trim(file_loc), form='unformatted', &
467 status='old', action='read')
468 read (1) pb%sf(:, :, :, r, i)
469 close (1)
470 else
471 call s_mpi_abort('File pb'//trim(file_num)// &
472 '.dat is missing in '//trim(t_step_dir)// &
473 '. Exiting.')
474 end if
475 end do
476
477 end do
478
479 do i = 1, nb
480 do r = 1, 4
481 ! Checking whether data file associated with variable position
482 ! of the currently manipulated bubble variable exists
483 write (file_num, '(I0)') sys_size + r + (i - 1)*4
484 file_loc = trim(t_step_dir)//'/mv'// &
485 trim(file_num)//'.dat'
486 inquire (file=trim(file_loc), exist=file_check)
487
488 ! If it exists, the data file is read
489 if (file_check) then
490 open (1, file=trim(file_loc), form='unformatted', &
491 status='old', action='read')
492 read (1) mv%sf(:, :, :, r, i)
493 close (1)
494 else
495 call s_mpi_abort('File mv'//trim(file_num)// &
496 '.dat is missing in '//trim(t_step_dir)// &
497 '. Exiting.')
498 end if
499 end do
500
501 end do
502 end if
503
504 ! Since the preexisting grid and initial condition data files have
505 ! been read in, the directory associated with the rank of the local
506 ! process may be cleaned out to make room for new pre-process data.
507 ! In addition, the time-step folder that will contain the new grid
508 ! and initial condition data are also generated.
509 call s_create_directory(trim(proc_rank_dir)//'/*')
510 call s_create_directory(trim(proc_rank_dir)//'/0')
511
512 end subroutine s_read_serial_ic_data_files
513
514 !> Cell-boundary data are checked for consistency by looking
515 !! at the (non-)uniform cell-width distributions for all the
516 !! active coordinate directions and making sure that all of
517 !! the cell-widths are positively valued
519
520#ifdef MFC_MPI
521
522 real(wp), allocatable, dimension(:) :: x_cb_glb, y_cb_glb, z_cb_glb
523
524 integer :: ifile, ierr, data_size
525 integer, dimension(MPI_STATUS_SIZE) :: status
526
527 character(LEN=path_len + 2*name_len) :: file_loc
528 logical :: file_exist
529
530 allocate (x_cb_glb(-1:m_glb))
531 allocate (y_cb_glb(-1:n_glb))
532 allocate (z_cb_glb(-1:p_glb))
533
534 ! Read in cell boundary locations in x-direction
535 file_loc = trim(case_dir)//'/restart_data'//trim(mpiiofs)//'x_cb.dat'
536 inquire (file=trim(file_loc), exist=file_exist)
537
538 if (file_exist) then
539 data_size = m_glb + 2
540 call mpi_file_open(mpi_comm_world, file_loc, mpi_mode_rdonly, mpi_info_int, ifile, ierr)
541 call mpi_file_read_all(ifile, x_cb_glb, data_size, mpi_p, status, ierr)
542 call mpi_file_close(ifile, ierr)
543 else
544 call s_mpi_abort('File '//trim(file_loc)//' is missing. Exiting. ')
545 end if
546
547 ! Assigning local cell boundary locations
548 x_cb(-1:m) = x_cb_glb((start_idx(1) - 1):(start_idx(1) + m))
549 ! Computing cell center locations
550 x_cc(0:m) = (x_cb(0:m) + x_cb(-1:(m - 1)))/2._wp
551 ! Computing minimum cell width
552 dx = minval(x_cb(0:m) - x_cb(-1:(m - 1)))
553 if (num_procs > 1) call s_mpi_reduce_min(dx)
554 ! Setting locations of domain bounds
555 x_domain%beg = x_cb(-1)
556 x_domain%end = x_cb(m)
557
558 if (n > 0) then
559 ! Read in cell boundary locations in y-direction
560 file_loc = trim(case_dir)//'/restart_data'//trim(mpiiofs)//'y_cb.dat'
561 inquire (file=trim(file_loc), exist=file_exist)
562
563 if (file_exist) then
564 data_size = n_glb + 2
565 call mpi_file_open(mpi_comm_world, file_loc, mpi_mode_rdonly, mpi_info_int, ifile, ierr)
566 call mpi_file_read_all(ifile, y_cb_glb, data_size, mpi_p, status, ierr)
567 call mpi_file_close(ifile, ierr)
568 else
569 call s_mpi_abort('File '//trim(file_loc)//' is missing. Exiting. ')
570 end if
571
572 ! Assigning local cell boundary locations
573 y_cb(-1:n) = y_cb_glb((start_idx(2) - 1):(start_idx(2) + n))
574 ! Computing cell center locations
575 y_cc(0:n) = (y_cb(0:n) + y_cb(-1:(n - 1)))/2._wp
576 ! Computing minimum cell width
577 dy = minval(y_cb(0:n) - y_cb(-1:(n - 1)))
578 if (num_procs > 1) call s_mpi_reduce_min(dy)
579 ! Setting locations of domain bounds
580 y_domain%beg = y_cb(-1)
581 y_domain%end = y_cb(n)
582
583 if (p > 0) then
584 ! Read in cell boundary locations in z-direction
585 file_loc = trim(case_dir)//'/restart_data'//trim(mpiiofs)//'z_cb.dat'
586 inquire (file=trim(file_loc), exist=file_exist)
587
588 if (file_exist) then
589 data_size = p_glb + 2
590 call mpi_file_open(mpi_comm_world, file_loc, mpi_mode_rdonly, mpi_info_int, ifile, ierr)
591 call mpi_file_read_all(ifile, z_cb_glb, data_size, mpi_p, status, ierr)
592 call mpi_file_close(ifile, ierr)
593 else
594 call s_mpi_abort('File '//trim(file_loc)//' is missing. Exiting. ')
595 end if
596
597 ! Assigning local cell boundary locations
598 z_cb(-1:p) = z_cb_glb((start_idx(3) - 1):(start_idx(3) + p))
599 ! Computing cell center locations
600 z_cc(0:p) = (z_cb(0:p) + z_cb(-1:(p - 1)))/2._wp
601 ! Computing minimum cell width
602 dz = minval(z_cb(0:p) - z_cb(-1:(p - 1)))
603 if (num_procs > 1) call s_mpi_reduce_min(dz)
604 ! Setting locations of domain bounds
605 z_domain%beg = z_cb(-1)
606 z_domain%end = z_cb(p)
607
608 end if
609 end if
610
611 deallocate (x_cb_glb, y_cb_glb, z_cb_glb)
612
613#endif
614
616
617 !> The goal of this subroutine is to read in any preexisting
618 !! initial condition data files so that they may be used by
619 !! the pre-process as a starting point in the creation of an
620 !! all new initial condition.
621 !! @param q_cons_vf_in Conservative variables
622 impure subroutine s_read_parallel_ic_data_files(q_cons_vf_in)
623
624 type(scalar_field), &
625 dimension(sys_size), &
626 intent(inout) :: q_cons_vf_in
627
628#ifdef MFC_MPI
629
630 integer :: ifile, ierr, data_size
631 integer, dimension(MPI_STATUS_SIZE) :: status
632 integer(KIND=MPI_OFFSET_KIND) :: disp
633 integer(KIND=MPI_OFFSET_KIND) :: m_mok, n_mok, p_mok
634 integer(KIND=MPI_OFFSET_KIND) :: wp_mok, var_mok, str_mok
635 integer(KIND=MPI_OFFSET_KIND) :: nvars_mok
636 integer(KIND=MPI_OFFSET_KIND) :: mok
637
638 character(LEN=path_len + 2*name_len) :: file_loc
639 logical :: file_exist
640
641 integer :: i
642
643 ! Open the file to read
644 if (cfl_adap_dt) then
645 write (file_loc, '(I0,A)') n_start, '.dat'
646 else
647 write (file_loc, '(I0,A)') t_step_start, '.dat'
648 end if
649 file_loc = trim(restart_dir)//trim(mpiiofs)//trim(file_loc)
650 inquire (file=trim(file_loc), exist=file_exist)
651
652 if (file_exist) then
653 call mpi_file_open(mpi_comm_world, file_loc, mpi_mode_rdonly, mpi_info_int, ifile, ierr)
654
655 call s_initialize_mpi_data(q_cons_vf_in)
656
657 ! Size of local arrays
658 data_size = (m + 1)*(n + 1)*(p + 1)
659
660 ! Resize some integers so MPI can read even the biggest files
661 m_mok = int(m_glb + 1, mpi_offset_kind)
662 n_mok = int(n_glb + 1, mpi_offset_kind)
663 p_mok = int(p_glb + 1, mpi_offset_kind)
664 wp_mok = int(8._wp, mpi_offset_kind)
665 mok = int(1._wp, mpi_offset_kind)
666 str_mok = int(name_len, mpi_offset_kind)
667 nvars_mok = int(sys_size, mpi_offset_kind)
668
669 ! Read the data for each variable
670 do i = 1, sys_size
671 var_mok = int(i, mpi_offset_kind)
672
673 ! Initial displacement to skip at beginning of file
674 disp = m_mok*max(mok, n_mok)*max(mok, p_mok)*wp_mok*(var_mok - 1)
675
676 call mpi_file_set_view(ifile, disp, mpi_p, mpi_io_data%view(i), &
677 'native', mpi_info_int, ierr)
678 call mpi_file_read(ifile, mpi_io_data%var(i)%sf, data_size, &
679 mpi_p, status, ierr)
680 end do
681
682 if (qbmm .and. .not. polytropic) then
683 do i = sys_size + 1, sys_size + 2*nb*4
684 var_mok = int(i, mpi_offset_kind)
685
686 ! Initial displacement to skip at beginning of file
687 disp = m_mok*max(mok, n_mok)*max(mok, p_mok)*wp_mok*(var_mok - 1)
688
689 call mpi_file_set_view(ifile, disp, mpi_p, mpi_io_data%view(i), &
690 'native', mpi_info_int, ierr)
691 call mpi_file_read(ifile, mpi_io_data%var(i)%sf, data_size, &
692 mpi_p, status, ierr)
693 end do
694 end if
695
696 call s_mpi_barrier()
697
698 call mpi_file_close(ifile, ierr)
699
700 else
701 call s_mpi_abort('File '//trim(file_loc)//' is missing. Exiting. ')
702 end if
703
704 call s_mpi_barrier()
705
706#endif
707
708 end subroutine s_read_parallel_ic_data_files
709
710 !> @brief Initializes all pre-process modules, allocates data structures, and sets I/O procedure pointers.
711 impure subroutine s_initialize_modules
712 ! Computation of parameters, allocation procedures, and/or any other tasks
713 ! needed to properly setup the modules
715 if (bubbles_euler .or. bubbles_lagrange) then
717 end if
723 call s_initialize_perturbation_module()
727
728 ! Create the D directory if it doesn't exit, to store
729 ! the serial data files
730 call s_create_directory('D')
731
732 ! Associate pointers for serial or parallel I/O
733 if (parallel_io .neqv. .true.) then
738 else
743 end if
744
745 end subroutine s_initialize_modules
746
747 !> @brief Reads an existing grid from data files or generates a new grid from user inputs.
748 impure subroutine s_read_grid()
749
750 if (old_grid) then
753 else
754 if (parallel_io .neqv. .true.) then
755 call s_generate_grid()
756 else
757 if (proc_rank == 0) call s_generate_grid()
758 call s_mpi_barrier()
761 end if
762 end if
763
764 end subroutine s_read_grid
765
766 !> @brief Generates or reads the initial condition, applies relaxation if needed, and writes output data files.
767 impure subroutine s_apply_initial_condition(start, finish)
768
769 real(wp), intent(inout) :: start, finish
770
771 integer :: j, k
772
773 ! Setting up the grid and the initial condition. If the grid is read in from
774 ! preexisting grid data files, it is checked for consistency. If the grid is
775 ! not read in, it is generated from scratch according to the inputs provided
776 ! by the user. The initial condition may also be read in. It in turn is not
777 ! checked for consistency since it WILL further be edited by the pre-process
778 ! and also because it may be incomplete at the time it is read in. Finally,
779 ! when the grid and initial condition are completely setup, they are written
780 ! to their respective data files.
781
782 ! Setting up grid and initial condition
783 call cpu_time(start)
784
786
788
789 ! hard-coded psi
790 if (hyper_cleaning) then
791 do j = 0, m
792 do k = 0, n
793 q_cons_vf(psi_idx)%sf(j, k, 0) = 1d-2*exp(-(x_cc(j)**2 + y_cc(k)**2)/(2.0*0.05**2))
794 q_prim_vf(psi_idx)%sf(j, k, 0) = q_cons_vf(psi_idx)%sf(j, k, 0)
795 end do
796 end do
797 end if
798
799 if (relax) then
800 if (proc_rank == 0) then
801 print *, 'initial condition might have been altered due to enforcement of &
802& pTg-equilirium (relax = "T" activated)'
803 end if
804
806 end if
807
809
810 call cpu_time(finish)
811 end subroutine s_apply_initial_condition
812
813 !> @brief Gathers processor timing data and writes elapsed wall-clock time to a summary file.
814 impure subroutine s_save_data(proc_time, time_avg, time_final, file_exists)
815
816 real(wp), dimension(:), intent(inout) :: proc_time
817 real(wp), intent(inout) :: time_avg, time_final
818 logical, intent(inout) :: file_exists
819
820 call s_mpi_barrier()
821
822 if (num_procs > 1) then
823 call mpi_bcast_time_step_values(proc_time, time_avg)
824 end if
825
826 if (proc_rank == 0) then
827 time_final = 0._wp
828 if (num_procs == 1) then
829 time_final = time_avg
830 print *, "Elapsed Time", time_final
831 else
832 time_final = maxval(proc_time)
833 print *, "Elapsed Time", time_final
834 end if
835 inquire (file='pre_time_data.dat', exist=file_exists)
836 if (file_exists) then
837 open (1, file='pre_time_data.dat', position='append', status='old')
838 write (1, *) num_procs, time_final
839 close (1)
840 else
841 open (1, file='pre_time_data.dat', status='new')
842 write (1, *) num_procs, time_final
843 close (1)
844 end if
845 end if
846 end subroutine s_save_data
847
848 !> @brief Initializes MPI, reads and validates user inputs on rank 0, and decomposes the computational domain.
849 impure subroutine s_initialize_mpi_domain
850 ! Initialization of the MPI environment
851
852 call s_mpi_initialize()
853
854 ! Rank 0 processor assigns default values to user inputs prior to reading
855 ! those in from the input file. Next, the user inputs are read in and their
856 ! consistency is checked. The detection of any inconsistencies automatically
857 ! leads to the termination of the pre-process.
858
859 if (proc_rank == 0) then
861 call s_read_input_file()
862 call s_check_input_file()
863
864 print '(" Pre-processing a ", I0, "x", I0, "x", I0, " case on ", I0, " rank(s)")', m, n, p, num_procs
865 end if
866
867 ! Broadcasting the user inputs to all of the processors and performing the
868 ! parallel computational domain decomposition. Neither procedure has to be
869 ! carried out if pre-process is in fact not truly executed in parallel.
873 end subroutine s_initialize_mpi_domain
874
875 !> @brief Finalizes all pre-process modules, deallocates resources, and shuts down MPI.
876 impure subroutine s_finalize_modules
877 ! Disassociate pointers for serial and parallel I/O
878 s_generate_grid => null()
879 s_read_grid_data_files => null()
880 s_read_ic_data_files => null()
881 s_write_data_files => null()
882
883 ! Deallocation procedures for the modules
890 call s_finalize_perturbation_module()
894 ! Finalization of the MPI environment
895 call s_mpi_finalize()
896 end subroutine s_finalize_modules
897
898end module m_start_up
Abstract interface for reading grid data files in serial or parallel.
Abstract interface for reading initial condition data files in serial or parallel.
integer, intent(in) k
integer, intent(in) j
Assigns initial primitive variables to computational cells based on patch geometry.
impure subroutine, public s_initialize_assign_variables_module
Allocates volume fraction sum and sets the patch primitive variable assignment procedure pointer.
impure subroutine, public s_finalize_assign_variables_module
Nullifies the patch primitive variable assignment procedure pointer.
Noncharacteristic and processor boundary condition application for ghost cells and buffer regions.
impure subroutine, public s_initialize_boundary_common_module()
Allocates and sets up boundary condition buffer arrays for all coordinate directions.
subroutine, public s_finalize_boundary_common_module()
Deallocates boundary condition buffer arrays allocated during module initialization.
Applies spatially varying boundary condition patches along domain edges and faces.
Validates geometry parameters and constraints for immersed boundary patches.
impure subroutine, public s_check_ib_patches
Validates the geometry parameters of all active and inactive immersed boundary patches.
Validates geometry parameters and constraints for initial condition patches.
impure subroutine, public s_check_patches
Validates the geometry parameters of all active and inactive initial condition patches.
Shared input validation checks for grid dimensions and AMD GPU compiler limits.
impure subroutine, public s_check_inputs_common
Checks compatibility of parameters in the input file. Used by all three stages.
Checks pre-process input file parameters for compatibility and correctness.
impure subroutine, public s_check_inputs
Checks compatibility of parameters in the input file. Used by the pre_process stage.
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.
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.
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.
real(wp) perturb_flow_mag
Magnitude of perturbation with perturb_flow flag.
real(wp) mixlayer_perturb_k0
Peak wavenumber of prescribed energy spectra with mixlayer_perturb flag Default value (k0 = 0....
logical cont_damage
continuum damage modeling
integer p_glb
Global number of cells in each direction.
logical igr
Use information geometric regularization.
logical hypoelasticity
activate hypoelasticity
impure subroutine s_assign_default_values_to_user_inputs
Assigns default values to user inputs prior to reading them in. This allows for an easier consistency...
impure subroutine s_finalize_global_parameters_module
Deallocates all global grid, index, and equation-of-state parameter arrays.
integer perturb_flow_fluid
Fluid to be perturbed with perturb_flow flag.
integer recon_type
Reconstruction Type.
integer mpi_info_int
MPI info for parallel IO with Lustre file systems.
type(ib_patch_parameters), dimension(num_patches_max) patch_ib
real(wp) dz
Minimum cell-widths in the x-, y- and z-coordinate directions.
type(int_bounds_info) bc_z
Boundary conditions in the x-, y- and z-coordinate directions.
integer num_fluids
Number of different fluids present in the flow.
logical pre_stress
activate pre_stressed domain
real(wp), dimension(:), allocatable y_cc
integer proc_rank
Rank of the local processor.
logical bc_io
whether or not to save BC data
real(wp), dimension(:), allocatable y_cb
type(bounds_info) z_domain
Locations of the domain bounds in the x-, y- and z-coordinate directions.
character(len=name_len) mpiiofs
integer, dimension(:), allocatable start_idx
Starting cell-center index of local processor in global grid.
integer sys_size
Number of unknowns in the system of equations.
type(simplex_noise_params) simplex_params
integer muscl_order
Order of accuracy for the MUSCL reconstruction.
real(wp) ptgalpha_eps
trigger parameter for the pTg relaxation procedure, phase change model
integer relax_model
Relax Model.
integer num_patches
Number of patches composing initial condition.
logical ib
Turn immersed boundaries on.
integer num_bc_patches
Number of boundary condition patches.
type(bc_patch_parameters), dimension(num_bc_patches_max) patch_bc
integer model_eqns
Multicomponent flow model.
integer precision
Precision of output files.
logical hyperelasticity
activate hyperelasticity
real(wp), dimension(:), allocatable z_cb
Locations of cell-boundaries (cb) in x-, y- and z-directions, respectively.
type(physical_parameters), dimension(num_fluids_max) fluid_pp
Database of the physical parameters of each of the fluids that is present in the flow....
impure subroutine s_initialize_global_parameters_module
Computation of parameters, allocation procedures, and/or any other tasks needed to properly setup the...
real(wp), dimension(:), allocatable x_cc
integer mixlayer_perturb_nk
Number of Fourier modes for perturbation with mixlayer_perturb flag.
integer perturb_sph_fluid
Fluid to be perturbed with perturb_sph flag.
real(wp), dimension(:), allocatable x_cb
logical relax
activate phase change
logical qbmm
Quadrature moment method.
logical old_grid
Use existing grid data.
real(wp) pi_fac
Factor for artificial pi_inf.
logical hyper_cleaning
Hyperbolic cleaning for MHD.
real(wp), dimension(num_fluids_max) fluid_rho
real(wp), dimension(:), allocatable z_cc
Locations of cell-centers (cc) in x-, y- and z-directions, respectively.
real(wp) pref
Reference parameters for Tait EOS.
real(wp) bx0
Constant magnetic field in the x-direction (1D).
logical stretch_z
Grid stretching flags for the x-, y- and z-coordinate directions.
logical adv_n
Solve the number density equation and compute alpha from number density.
integer num_procs
Number of processors.
character(len=path_len) case_dir
Case folder location.
type(ic_patch_parameters), dimension(num_patches_max) patch_icpp
Database of the initial condition patch parameters (icpp) for each of the patches employed in the con...
integer weno_order
Order of accuracy for the WENO reconstruction.
logical mhd
Magnetohydrodynamics.
logical parallel_io
Format of the data files.
type(cell_num_bounds) cells_bounds
logical down_sample
Down-sample the output data.
logical file_per_process
type of data output
real(wp) palpha_eps
trigger parameter for the p relaxation procedure, phase change model
integer t_step_start
Existing IC/grid folder.
type(mpi_io_var), public mpi_io_data
real(wp) mixlayer_vel_coef
Coefficient for the hyperbolic tangent streamwise velocity profile.
impure subroutine s_initialize_parallel_io
Configures MPI parallel I/O settings and allocates processor coordinate arrays.
logical mpp_lim
Alpha limiter.
integer igr_order
IGR reconstruction order.
integer psi_idx
Index of hyperbolic cleaning state variable for MHD.
type(subgrid_bubble_physical_parameters) bub_pp
real(wp) rhorv
standard deviations in R/V
logical relativity
Relativity for RMHD.
integer num_ibs
Number of immersed boundaries.
logical mixlayer_vel_profile
Set hyperbolic tangent streamwise velocity profile.
integer(kind=8) nglobal
Global number of cells in the domain.
logical mixlayer_perturb
Superimpose instability waves to surrounding fluid flow.
Generates uniform or stretched rectilinear grids with hyperbolic-tangent spacing.
Definition m_grid.f90:6
impure subroutine, public s_generate_serial_grid
The following subroutine generates either a uniform or non-uniform rectilinear grid in serial,...
Definition m_grid.f90:48
impure subroutine, public s_initialize_grid_module
Computation of parameters, allocation procedures, and/or any other tasks needed to properly setup the...
Definition m_grid.f90:336
procedure(s_generate_abstract_grid), pointer, public s_generate_grid
Definition m_grid.f90:38
impure subroutine, public s_generate_parallel_grid
The following subroutine generates either a uniform or non-uniform rectilinear grid in parallel,...
Definition m_grid.f90:187
impure subroutine, public s_finalize_grid_module
Deallocation procedures for the module.
Definition m_grid.f90:347
Basic floating-point utilities: approximate equality, default detection, and coordinate bounds.
elemental subroutine, public s_update_cell_bounds(bounds, m, n, p)
Updates the min and max number of cells in each set of axes.
Utility routines for bubble model setup, coordinate transforms, array sampling, and special functions...
impure subroutine, public s_initialize_bubbles_model()
bubbles_euler + polytropic bubbles_euler + non-polytropic bubbles_lagrange + non-polytropic
Allocate memory and read initial condition data for IC extrusion.
Assembles initial conditions by layering prioritized patches via constructive solid geometry.
type(scalar_field), dimension(:), allocatable q_cons_vf
conservative variables
type(integer_field), dimension(:, :), allocatable bc_type
bc_type fields
impure subroutine s_initialize_initial_condition_module
Computation of parameters, allocation procedures, and/or any other tasks needed to properly setup the...
impure subroutine s_finalize_initial_condition_module
Deallocation procedures for the module.
type(scalar_field), dimension(:), allocatable q_prim_vf
primitive variables
impure subroutine s_generate_initial_condition
This subroutine peruses the patches and depending on the type of geometry associated with a particula...
MPI communication layer: domain decomposition, halo exchange, reductions, and parallel I/O setup.
impure subroutine s_mpi_abort(prnt, code)
The subroutine terminates the MPI execution environment.
impure subroutine s_initialize_mpi_common_module
The computation of parameters, the allocation of memory, the association of pointers and/or the execu...
impure subroutine s_mpi_barrier
Halts all processes until all have reached barrier.
impure subroutine s_mpi_initialize
The subroutine initializes the MPI execution environment and queries both the number of processors wh...
impure subroutine s_initialize_mpi_data(q_cons_vf, ib_markers, beta)
impure subroutine s_mpi_finalize
The subroutine finalizes the MPI execution environment.
impure subroutine mpi_bcast_time_step_values(proc_time, time_avg)
Gathers per-rank time step wall-clock times onto rank 0 for performance reporting.
impure subroutine s_mpi_reduce_min(var_loc)
The following subroutine takes the inputted variable and determines its minimum value on the entire c...
impure subroutine s_finalize_mpi_common_module
Module deallocation and/or disassociation procedures.
subroutine s_mpi_decompose_computational_domain
The purpose of this procedure is to optimally decompose the computational domain among the available ...
Broadcasts user inputs and decomposes the domain across MPI ranks for pre-processing.
impure subroutine s_mpi_bcast_user_inputs
Since only processor with rank 0 is in charge of reading and checking the consistency of the user pro...
Phase transition relaxation solvers for liquid-vapor flows with cavitation and boiling.
impure subroutine, public s_finalize_relaxation_solver_module
This subroutine finalizes the phase change module.
subroutine, public s_infinite_relaxation_k(q_cons_vf)
This subroutine is created to activate either the pT- (N fluids) or the pTg-equilibrium (2 fluids for...
impure subroutine, public s_initialize_phasechange_module
The purpose of this subroutine is to initialize the phase change module by setting the parameters nee...
Reads and validates user inputs, loads existing grid/IC data, and initializes pre-process modules.
impure subroutine, public s_read_serial_ic_data_files(q_cons_vf_in)
The goal of this subroutine is to read in any preexisting initial condition data files so that they m...
impure subroutine, public s_initialize_modules
Initializes all pre-process modules, allocates data structures, and sets I/O procedure pointers.
impure subroutine, public s_save_data(proc_time, time_avg, time_final, file_exists)
Gathers processor timing data and writes elapsed wall-clock time to a summary file.
impure subroutine, public s_apply_initial_condition(start, finish)
Generates or reads the initial condition, applies relaxation if needed, and writes output data files.
character(len=path_len+name_len) proc_rank_dir
Location of the folder associated with the rank of the local processor.
impure subroutine, public s_read_serial_grid_data_files
The goal of this subroutine is to read in any preexisting grid data as well as based on the imported ...
impure subroutine, public s_read_parallel_ic_data_files(q_cons_vf_in)
The goal of this subroutine is to read in any preexisting initial condition data files so that they m...
procedure(s_read_abstract_ic_data_files), pointer, public s_read_ic_data_files
impure subroutine, public s_read_grid()
Reads an existing grid from data files or generates a new grid from user inputs.
impure subroutine, public s_check_grid_data_files
Cell-boundary data are checked for consistency by looking at the (non-)uniform cell-width distributio...
impure subroutine, public s_initialize_mpi_domain
Initializes MPI, reads and validates user inputs on rank 0, and decomposes the computational domain.
impure subroutine, public s_finalize_modules
Finalizes all pre-process modules, deallocates resources, and shuts down MPI.
impure subroutine, public s_read_input_file
Reads the configuration file pre_process.inp, in order to populate the parameters in module m_global_...
impure subroutine, public s_check_input_file
Checking that the user inputs make sense, i.e. that the individual choices are compatible with the co...
procedure(s_read_abstract_grid_data_files), pointer, public s_read_grid_data_files
impure subroutine, public s_read_parallel_grid_data_files
Cell-boundary data are checked for consistency by looking at the (non-)uniform cell-width distributio...
Conservative-to-primitive variable conversion, mixture property evaluation, and pressure computation.
impure subroutine, public s_initialize_variables_conversion_module
The computation of parameters, the allocation of memory, the association of pointers and/or the execu...
impure subroutine s_finalize_variables_conversion_module()
Deallocates fluid property arrays and post-processing fields allocated during module initialization.
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).