MFC
Exascale flow solver
Loading...
Searching...
No Matches
m_global_parameters.fpp.f90
Go to the documentation of this file.
1# 1 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
2!>
3!! @file
4!! @brief Contains module m_global_parameters
5
6# 1 "/home/runner/work/MFC/MFC/src/common/include/case.fpp" 1
7! This file exists so that Fypp can be run without generating case.fpp files for
8! each target. This is useful when generating documentation, for example. This
9! should also let MFC be built with CMake directly, without invoking mfc.sh.
10
11! For pre-process.
12# 9 "/home/runner/work/MFC/MFC/src/common/include/case.fpp"
13
14! For moving immersed boundaries in simulation
15# 14 "/home/runner/work/MFC/MFC/src/common/include/case.fpp"
16# 6 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp" 2
17
18!> @brief Global parameters for the post-process: domain geometry, equation of state, and output database settings
20
21#ifdef MFC_MPI
22 use mpi !< message passing interface (mpi) module
23#endif
24
25 use m_derived_types !< definitions of the derived types
26
27 use m_helper_basic !< functions to compare floating point numbers
28
29 use m_thermochem, only: num_species, species_names
30
31 implicit none
32
33 !> @name Logistics
34 !> @{
35 integer :: num_procs !< Number of processors
36 character(LEN=path_len) :: case_dir !< Case folder location
37 !> @}
38
39 ! Computational Domain Parameters
40
41 integer :: proc_rank !< Rank of the local processor
42
43 !> @name Number of cells in the x-, y- and z-coordinate directions
44 !> @{
45 integer :: m, m_root
46 integer :: n
47 integer :: p
48 !> @}
49
50 !> @name Max and min number of cells in a direction of each combination of x-,y-, and z-
52
53 integer(kind=8) :: nglobal ! Total number of cells in global domain
54
55 !> @name Cylindrical coordinates (either axisymmetric or full 3D)
56 !> @{
57 logical :: cyl_coord
58 integer :: grid_geometry
59 !> @}
60
61 !> @name Global number of cells in each direction
62 !> @{
63 integer :: m_glb, n_glb, p_glb
64 !> @}
65
66 integer :: num_dims !< Number of spatial dimensions
67 integer :: num_vels !< Number of velocity components (different from num_dims for mhd)
68
69 !> @name Cell-boundary locations in the x-, y- and z-coordinate directions
70 !> @{
71 real(wp), allocatable, dimension(:) :: x_cb, x_root_cb, y_cb, z_cb
72 !> @}
73
74 !> @name Cell-center locations in the x-, y- and z-coordinate directions
75 !> @{
76 real(wp), allocatable, dimension(:) :: x_cc, x_root_cc, y_cc, z_cc
77 real(sp), allocatable, dimension(:) :: x_root_cc_s, x_cc_s
78 !> @}
79
80 !> Cell-width distributions in the x-, y- and z-coordinate directions
81 !> @{
82 real(wp), allocatable, dimension(:) :: dx, dy, dz
83 !> @}
84
85 integer :: buff_size !<
86 !! Number of cells in buffer region. For the variables which feature a buffer
87 !! region, this region is used to store information outside the computational
88 !! domain based on the boundary conditions.
89
90 integer :: t_step_start !< First time-step directory
91 integer :: t_step_stop !< Last time-step directory
92 integer :: t_step_save !< Interval between consecutive time-step directory
93
94 !> @name IO options for adaptive time-stepping
95 !> @{
97 real(wp) :: t_save
98 real(wp) :: t_stop
99 real(wp) :: cfl_target
100 integer :: n_save
101 integer :: n_start
102 !> @}
103
104 ! NOTE: The variables m_root, x_root_cb and x_root_cc contain the grid data
105 ! of the defragmented computational domain. They are only used in 1D. For
106 ! serial simulations, they are equal to m, x_cb and x_cc, respectively.
107
108 !> @name Simulation Algorithm Parameters
109 !> @{
110 integer :: model_eqns !< Multicomponent flow model
111 integer :: num_fluids !< Number of different fluids present in the flow
112 logical :: relax !< phase change
113 integer :: relax_model !< Phase change relaxation model
114 logical :: mpp_lim !< Maximum volume fraction limiter
115 integer :: sys_size !< Number of unknowns in the system of equations
116 integer :: recon_type !< Which type of reconstruction to use
117 integer :: weno_order !< Order of accuracy for the WENO reconstruction
118 integer :: muscl_order !< Order of accuracy for the MUSCL reconstruction
119 logical :: mixture_err !< Mixture error limiter
120 logical :: alt_soundspeed !< Alternate sound speed
121 logical :: mhd !< Magnetohydrodynamics
122 logical :: relativity !< Relativity for RMHD
123 logical :: hypoelasticity !< Turn hypoelasticity on
124 logical :: hyperelasticity !< Turn hyperelasticity on
125 logical :: elasticity !< elasticity modeling, true for hyper or hypo
126 integer :: b_size !< Number of components in the b tensor
127 integer :: tensor_size !< Number of components in the nonsymmetric tensor
128 logical :: cont_damage !< Continuum damage modeling
129 logical :: hyper_cleaning !< Hyperbolic cleaning for MHD
130 logical :: igr !< enable IGR
131 integer :: igr_order !< IGR reconstruction order
132 logical, parameter :: chemistry = .false. !< Chemistry modeling
133 !> @}
134
135 integer :: avg_state !< Average state evaluation method
136
137 !> @name Annotations of the structure, i.e. the organization, of the state vectors
138 !> @{
139 type(int_bounds_info) :: cont_idx !< Indexes of first & last continuity eqns.
140 type(int_bounds_info) :: mom_idx !< Indexes of first & last momentum eqns.
141 integer :: e_idx !< Index of energy equation
142 integer :: n_idx !< Index of number density
143 integer :: beta_idx !< Index of lagrange bubbles beta
144 type(int_bounds_info) :: adv_idx !< Indexes of first & last advection eqns.
145 type(int_bounds_info) :: internalenergies_idx !< Indexes of first & last internal energy eqns.
146 type(bub_bounds_info) :: bub_idx !< Indexes of first & last bubble variable eqns.
147 integer :: gamma_idx !< Index of specific heat ratio func. eqn.
148 integer :: alf_idx !< Index of specific heat ratio func. eqn.
149 integer :: pi_inf_idx !< Index of liquid stiffness func. eqn.
150 type(int_bounds_info) :: b_idx !< Indexes of first and last magnetic field eqns.
151 type(int_bounds_info) :: stress_idx !< Indices of elastic stresses
152 type(int_bounds_info) :: xi_idx !< Indexes of first and last reference map eqns.
153 integer :: c_idx !< Index of color function
154 type(int_bounds_info) :: species_idx !< Indexes of first & last concentration eqns.
155 integer :: damage_idx !< Index of damage state variable (D) for continuum damage model
156 integer :: psi_idx !< Index of hyperbolic cleaning state variable for MHD
157 !> @}
158
159 ! Cell Indices for the (local) interior points (O-m, O-n, 0-p).
160 ! Stands for "InDices With BUFFer".
162
163 ! Cell Indices for the entire (local) domain. In simulation, this includes
164 ! the buffer region. idwbuff and idwint are the same otherwise.
165 ! Stands for "InDices With BUFFer".
167
168 integer :: num_bc_patches
169 logical :: bc_io
170 !> @name Boundary conditions in the x-, y- and z-coordinate directions
171 !> @{
173 !> @}
174
175 integer :: shear_num !! Number of shear stress components
176 integer, dimension(3) :: shear_indices !<
177 !! Indices of the stress components that represent shear stress
178 integer :: shear_bc_flip_num !<
179 !! Number of shear stress components to reflect for boundary conditions
180 integer, dimension(3, 2) :: shear_bc_flip_indices !<
181 !! Indices of shear stress components to reflect for boundary conditions.
182 !! Size: (1:3, 1:shear_BC_flip_num) for (x/y/z, [indices])
183
184 logical :: parallel_io !< Format of the data files
185 logical :: sim_data
186 logical :: file_per_process !< output format
187
188 integer, allocatable, dimension(:) :: proc_coords !<
189 !! Processor coordinates in MPI_CART_COMM
190
191 integer, allocatable, dimension(:) :: start_idx !<
192 !! Starting cell-center index of local processor in global grid
193
194 integer :: num_ibs !< Number of immersed boundaries
195
196#ifdef MFC_MPI
197
198 type(mpi_io_var), public :: mpi_io_data
202 real(wp), allocatable, dimension(:, :), public :: mpi_io_data_lg_bubbles
203
204#endif
205
206 !> @name MPI info for parallel IO with Lustre file systems
207 !> @{
208 character(LEN=name_len) :: mpiiofs
209 integer :: mpi_info_int
210 !> @}
211
212 type(physical_parameters), dimension(num_fluids_max) :: fluid_pp !<
213 !! Database of the physical parameters of each of the fluids that is present
214 !! in the flow. These include the stiffened gas equation of state parameters,
215 !! and the Reynolds numbers.
216
217 ! Subgrid Bubble Parameters
219
220 real(wp), allocatable, dimension(:) :: adv !< Advection variables
221
222 ! Formatted Database File(s) Structure Parameters
223
224 integer :: format !< Format of the database file(s)
225
226 integer :: precision !< Floating point precision of the database file(s)
227 logical :: down_sample !< down sampling of the database file(s)
228
229 logical :: output_partial_domain !< Specify portion of domain to output for post-processing
230
231 type(bounds_info) :: x_output, y_output, z_output !< Portion of domain to output for post-processing
232 type(int_bounds_info) :: x_output_idx, y_output_idx, z_output_idx !< Indices of domain to output for post-processing
233
234 !> @name Size of the ghost zone layer in the x-, y- and z-coordinate directions.
235 !! The definition of the ghost zone layers is only necessary when using the
236 !! Silo database file format in multidimensions. These zones provide VisIt
237 !! with the subdomain connectivity information that it requires in order to
238 !! produce smooth plots.
239 !> @{
241 !> @}
242
243 !> @name The list of all possible flow variables that may be written to a database
244 !! file. It includes partial densities, density, momentum, velocity, energy,
245 !! pressure, volume fraction(s), specific heat ratio function, specific heat
246 !! ratio, liquid stiffness function, liquid stiffness, primitive variables,
247 !! conservative variables, speed of sound, the vorticity,
248 !! and the numerical Schlieren function.
249 !> @{
250 logical, dimension(num_fluids_max) :: alpha_rho_wrt
251 logical :: rho_wrt
252 logical, dimension(3) :: mom_wrt
253 logical, dimension(3) :: vel_wrt
254 integer :: flux_lim
255 logical, dimension(3) :: flux_wrt
256 logical :: e_wrt
257 logical, dimension(num_fluids_max) :: alpha_rho_e_wrt
258 logical :: fft_wrt
259 logical :: dummy !< AMDFlang workaround: keep a dummy logical to avoid a compiler case-optimization bug when a parameter+GPU-kernel conditional is false
260 logical :: pres_wrt
261 logical, dimension(num_fluids_max) :: alpha_wrt
262 logical :: gamma_wrt
263 logical :: heat_ratio_wrt
264 logical :: pi_inf_wrt
265 logical :: pres_inf_wrt
266 logical :: prim_vars_wrt
267 logical :: cons_vars_wrt
268 logical :: c_wrt
269 logical, dimension(3) :: omega_wrt
270 logical :: qm_wrt
271 logical :: liutex_wrt
272 logical :: schlieren_wrt
273 logical :: cf_wrt
274 logical :: ib
275 logical :: chem_wrt_y(1:num_species)
276 logical :: chem_wrt_t
277 logical :: lag_header
278 logical :: lag_txt_wrt
279 logical :: lag_db_wrt
280 logical :: lag_id_wrt
281 logical :: lag_pos_wrt
283 logical :: lag_vel_wrt
284 logical :: lag_rad_wrt
285 logical :: lag_rvel_wrt
286 logical :: lag_r0_wrt
287 logical :: lag_rmax_wrt
288 logical :: lag_rmin_wrt
289 logical :: lag_dphidt_wrt
290 logical :: lag_pres_wrt
291 logical :: lag_mv_wrt
292 logical :: lag_mg_wrt
293 logical :: lag_betat_wrt
294 logical :: lag_betac_wrt
295 !> @}
296
297 real(wp), dimension(num_fluids_max) :: schlieren_alpha !<
298 !! Amplitude coefficients of the numerical Schlieren function that are used
299 !! to adjust the intensity of numerical Schlieren renderings for individual
300 !! fluids. This enables waves and interfaces of varying strengths and in all
301 !! of the fluids to be made simultaneously visible on a single plot.
302
303 integer :: fd_order !<
304 !! The order of the finite-difference (fd) approximations of the first-order
305 !! derivatives that need to be evaluated when vorticity and/or the numerical
306 !! Schlieren function are to be outputted to the formatted database file(s).
307
308 integer :: fd_number !<
309 !! The finite-difference number is given by MAX(1, fd_order/2). Essentially,
310 !! it is a measure of the half-size of the finite-difference stencil for the
311 !! selected order of accuracy.
312
313 !> @name Reference parameters for Tait EOS
314 !> @{
315 real(wp) :: rhoref, pref
316 !> @}
317
319 !> @name Bubble modeling variables and parameters
320 !> @{
321 integer :: nb
322 real(wp) :: eu, ca, web, re_inv
323 real(wp), dimension(:), allocatable :: weight, r0
324 logical :: bubbles_euler
325 logical :: qbmm
326 logical :: polytropic
327 logical :: polydisperse
328 logical :: adv_n
329 integer :: thermal !< 1 = adiabatic, 2 = isotherm, 3 = transfer
330 real(wp) :: phi_vg, phi_gv, pe_c, tw, k_vl, k_gl
331 real(wp) :: gam_m
332 real(wp), dimension(:), allocatable :: pb0, mass_g0, mass_v0, pe_t, k_v, k_g
333 real(wp), dimension(:), allocatable :: re_trans_t, re_trans_c, im_trans_t, im_trans_c, omegan
334 real(wp) :: r0ref, p0ref, rho0ref, t0ref, ss, pv, vd, mu_l, mu_v, mu_g, &
336 real(wp) :: g
337 real(wp) :: poly_sigma
338 real(wp) :: sigr
339 integer :: nmom
340 !> @}
341
342 !> @name surface tension coefficient
343 !> @{
344
345 real(wp) :: sigma
347 !> @}
348
349 !> @name Index variables used for m_variables_conversion
350 !> @{
351 integer :: momxb, momxe
352 integer :: advxb, advxe
353 integer :: contxb, contxe
354 integer :: intxb, intxe
355 integer :: bubxb, bubxe
356 integer :: strxb, strxe
357 integer :: xibeg, xiend
358 integer :: chemxb, chemxe
359 !> @}
360
361 !> @name Lagrangian bubbles
362 !> @{
364 !> @}
365
366 real(wp) :: bx0 !< Constant magnetic field in the x-direction (1D)
367
368 real(wp) :: wall_time, wall_time_avg !< Wall time measurements
369
370contains
371
372 !> Assigns default values to user inputs prior to reading
373 !! them in. This allows for an easier consistency check of
374 !! these parameters once they are read from the input file.
376
377 integer :: i !< Generic loop iterator
378
379 ! Logistics
380 case_dir = '.'
381
382 ! Computational domain parameters
383 m = dflt_int; n = 0; p = 0
385
386 m_root = dflt_int
387 cyl_coord = .false.
388
389 t_step_start = dflt_int
390 t_step_stop = dflt_int
391 t_step_save = dflt_int
392
393 cfl_adap_dt = .false.
394 cfl_const_dt = .false.
395 cfl_dt = .false.
396 cfl_target = dflt_real
397 t_save = dflt_real
398 n_start = dflt_int
399 t_stop = dflt_real
400
401 ! Simulation algorithm parameters
402 model_eqns = dflt_int
403 num_fluids = dflt_int
404 recon_type = weno_type
405 weno_order = dflt_int
406 muscl_order = dflt_int
407 mixture_err = .false.
408 alt_soundspeed = .false.
409 relax = .false.
410 relax_model = dflt_int
411
412 mhd = .false.
413 relativity = .false.
414
415 hypoelasticity = .false.
416 hyperelasticity = .false.
417 elasticity = .false.
418 b_size = dflt_int
419 tensor_size = dflt_int
420 cont_damage = .false.
421 hyper_cleaning = .false.
422 igr = .false.
423
424 bc_x%beg = dflt_int; bc_x%end = dflt_int
425 bc_y%beg = dflt_int; bc_y%end = dflt_int
426 bc_z%beg = dflt_int; bc_z%end = dflt_int
427 bc_io = .false.
428 num_bc_patches = dflt_int
429
430# 420 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
431# 421 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
432 bc_x%vb1 = 0._wp
433 bc_x%ve1 = 0._wp
434# 421 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
435 bc_x%vb2 = 0._wp
436 bc_x%ve2 = 0._wp
437# 421 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
438 bc_x%vb3 = 0._wp
439 bc_x%ve3 = 0._wp
440# 424 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
441# 420 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
442# 421 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
443 bc_y%vb1 = 0._wp
444 bc_y%ve1 = 0._wp
445# 421 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
446 bc_y%vb2 = 0._wp
447 bc_y%ve2 = 0._wp
448# 421 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
449 bc_y%vb3 = 0._wp
450 bc_y%ve3 = 0._wp
451# 424 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
452# 420 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
453# 421 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
454 bc_z%vb1 = 0._wp
455 bc_z%ve1 = 0._wp
456# 421 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
457 bc_z%vb2 = 0._wp
458 bc_z%ve2 = 0._wp
459# 421 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
460 bc_z%vb3 = 0._wp
461 bc_z%ve3 = 0._wp
462# 424 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
463# 425 "/home/runner/work/MFC/MFC/src/post_process/m_global_parameters.fpp"
464
465 chem_params%gamma_method = 1
466 chem_params%transport_model = 1
467
468 ! Fluids physical parameters
469 do i = 1, num_fluids_max
470 fluid_pp(i)%gamma = dflt_real
471 fluid_pp(i)%pi_inf = dflt_real
472 fluid_pp(i)%cv = 0._wp
473 fluid_pp(i)%qv = 0._wp
474 fluid_pp(i)%qvp = 0._wp
475 fluid_pp(i)%G = dflt_real
476 end do
477
478 ! Subgrid bubble parameters
479 bub_pp%R0ref = dflt_real; r0ref = dflt_real
480 bub_pp%p0ref = dflt_real; p0ref = dflt_real
481 bub_pp%rho0ref = dflt_real; rho0ref = dflt_real
482 bub_pp%T0ref = dflt_real; t0ref = dflt_real
483 bub_pp%ss = dflt_real; ss = dflt_real
484 bub_pp%pv = dflt_real; pv = dflt_real
485 bub_pp%vd = dflt_real; vd = dflt_real
486 bub_pp%mu_l = dflt_real; mu_l = dflt_real
487 bub_pp%mu_v = dflt_real; mu_v = dflt_real
488 bub_pp%mu_g = dflt_real; mu_g = dflt_real
489 bub_pp%gam_v = dflt_real; gam_v = dflt_real
490 bub_pp%gam_g = dflt_real; gam_g = dflt_real
491 bub_pp%M_v = dflt_real; m_v = dflt_real
492 bub_pp%M_g = dflt_real; m_g = dflt_real
493 bub_pp%k_v = dflt_real;
494 bub_pp%k_g = dflt_real;
495 bub_pp%cp_v = dflt_real; cp_v = dflt_real
496 bub_pp%cp_g = dflt_real; cp_g = dflt_real
497 bub_pp%R_v = dflt_real; r_v = dflt_real
498 bub_pp%R_g = dflt_real; r_g = dflt_real
499
500 ! Formatted database file(s) structure parameters
501 format = dflt_int
502
503 precision = dflt_int
504 down_sample = .false.
505
506 alpha_rho_wrt = .false.
507 alpha_rho_e_wrt = .false.
508 rho_wrt = .false.
509 mom_wrt = .false.
510 vel_wrt = .false.
511 chem_wrt_y = .false.
512 chem_wrt_t = .false.
513 flux_lim = dflt_int
514 flux_wrt = .false.
515 parallel_io = .false.
516 file_per_process = .false.
517 e_wrt = .false.
518 fft_wrt = .false.
519 dummy = .false.
520 pres_wrt = .false.
521 alpha_wrt = .false.
522 gamma_wrt = .false.
523 heat_ratio_wrt = .false.
524 pi_inf_wrt = .false.
525 pres_inf_wrt = .false.
526 prim_vars_wrt = .false.
527 cons_vars_wrt = .false.
528 c_wrt = .false.
529 omega_wrt = .false.
530 qm_wrt = .false.
531 liutex_wrt = .false.
532 schlieren_wrt = .false.
533 sim_data = .false.
534 cf_wrt = .false.
535 ib = .false.
536 lag_txt_wrt = .false.
537 lag_header = .true.
538 lag_db_wrt = .false.
539 lag_id_wrt = .true.
540 lag_pos_wrt = .true.
541 lag_pos_prev_wrt = .false.
542 lag_vel_wrt = .true.
543 lag_rad_wrt = .true.
544 lag_rvel_wrt = .false.
545 lag_r0_wrt = .false.
546 lag_rmax_wrt = .false.
547 lag_rmin_wrt = .false.
548 lag_dphidt_wrt = .false.
549 lag_pres_wrt = .false.
550 lag_mv_wrt = .false.
551 lag_mg_wrt = .false.
552 lag_betat_wrt = .false.
553 lag_betac_wrt = .false.
554
555 schlieren_alpha = dflt_real
556
557 fd_order = dflt_int
558 avg_state = dflt_int
559
560 ! Tait EOS
561 rhoref = dflt_real
562 pref = dflt_real
563
564 ! Bubble modeling
565 bubbles_euler = .false.
566 qbmm = .false.
567 r0ref = dflt_real
568 nb = dflt_int
569 polydisperse = .false.
570 poly_sigma = dflt_real
571 sigr = dflt_real
572 sigma = dflt_real
573 surface_tension = .false.
574 adv_n = .false.
575
576 ! Lagrangian bubbles modeling
577 bubbles_lagrange = .false.
578
579 ! IBM
580 num_ibs = dflt_int
581
582 ! Output partial domain
583 output_partial_domain = .false.
584 x_output%beg = dflt_real
585 x_output%end = dflt_real
586 y_output%beg = dflt_real
587 y_output%end = dflt_real
588 z_output%beg = dflt_real
589 z_output%end = dflt_real
590
591 ! MHD
592 bx0 = dflt_real
593
595
596 !> Computation of parameters, allocation procedures, and/or
597 !! any other tasks needed to properly setup the module
599
600 integer :: i, j, fac
601
602 ! Setting m_root equal to m in the case of a 1D serial simulation
603 if (n == 0) m_root = m_glb
604
605 ! Gamma/Pi_inf Model
606 if (model_eqns == 1) then
607
608 ! Setting number of fluids
609 num_fluids = 1
610
611 ! Annotating structure of the state and flux vectors belonging
612 ! to the system of equations defined by the selected number of
613 ! spatial dimensions and the gamma/pi_inf model
614 cont_idx%beg = 1
615 cont_idx%end = cont_idx%beg
616 mom_idx%beg = cont_idx%end + 1
617 mom_idx%end = cont_idx%end + num_vels
618 e_idx = mom_idx%end + 1
619 adv_idx%beg = e_idx + 1
620 adv_idx%end = adv_idx%beg + 1
621 gamma_idx = adv_idx%beg
622 pi_inf_idx = adv_idx%end
623 sys_size = adv_idx%end
624
625 ! Volume Fraction Model (5-equation model)
626 else if (model_eqns == 2) then
627
628 ! Annotating structure of the state and flux vectors belonging
629 ! to the system of equations defined by the selected number of
630 ! spatial dimensions and the volume fraction model
631 cont_idx%beg = 1
632 cont_idx%end = num_fluids
633 mom_idx%beg = cont_idx%end + 1
634 mom_idx%end = cont_idx%end + num_vels
635 e_idx = mom_idx%end + 1
636
637 if (igr) then
638 ! Volume fractions are stored in the indices immediately following
639 ! the energy equation. IGR tracks a total of (N-1) volume fractions
640 ! for N fluids, hence the "-1" in adv_idx%end. If num_fluids = 1
641 ! then adv_idx%end < adv_idx%beg, which skips all loops over the
642 ! volume fractions since there is no volume fraction to track
643 adv_idx%beg = e_idx + 1 ! Alpha for fluid 1
644 adv_idx%end = e_idx + num_fluids - 1
645 else
646 ! Volume fractions are stored in the indices immediately following
647 ! the energy equation. WENO/MUSCL + Riemann tracks a total of (N)
648 ! volume fractions for N fluids, hence the lack of "-1" in adv_idx%end
649 adv_idx%beg = e_idx + 1
650 adv_idx%end = e_idx + num_fluids
651 end if
652
653 sys_size = adv_idx%end
654
655 if (bubbles_euler) then
656 alf_idx = adv_idx%end
657 else
658 alf_idx = 1
659 end if
660
661 if (qbmm) then
662 nmom = 6
663 end if
664
665 if (bubbles_euler) then
666
667 bub_idx%beg = sys_size + 1
668 if (qbmm) then
669 bub_idx%end = adv_idx%end + nb*nmom
670 else
671 if (.not. polytropic) then
672 bub_idx%end = sys_size + 4*nb
673 else
674 bub_idx%end = sys_size + 2*nb
675 end if
676 end if
677 sys_size = bub_idx%end
678
679 if (adv_n) then
680 n_idx = bub_idx%end + 1
682 end if
683
684 allocate (bub_idx%rs(nb), bub_idx%vs(nb))
685 allocate (bub_idx%ps(nb), bub_idx%ms(nb))
686
687 if (qbmm) then
688 allocate (bub_idx%moms(nb, nmom))
689 do i = 1, nb
690 do j = 1, nmom
691 bub_idx%moms(i, j) = bub_idx%beg + (j - 1) + (i - 1)*nmom
692 end do
693 bub_idx%rs(i) = bub_idx%moms(i, 2)
694 bub_idx%vs(i) = bub_idx%moms(i, 3)
695 end do
696 else
697 do i = 1, nb
698 if (polytropic .neqv. .true.) then
699 fac = 4
700 else
701 fac = 2
702 end if
703
704 bub_idx%rs(i) = bub_idx%beg + (i - 1)*fac
705 bub_idx%vs(i) = bub_idx%rs(i) + 1
706
707 if (polytropic .neqv. .true.) then
708 bub_idx%ps(i) = bub_idx%vs(i) + 1
709 bub_idx%ms(i) = bub_idx%ps(i) + 1
710 end if
711 end do
712 end if
713 end if
714
715 if (bubbles_lagrange) then
716 beta_idx = sys_size + 1
718 end if
719
720 if (mhd) then
721 b_idx%beg = sys_size + 1
722 if (n == 0) then
723 b_idx%end = sys_size + 2 ! 1D: By, Bz
724 else
725 b_idx%end = sys_size + 3 ! 2D/3D: Bx, By, Bz
726 end if
727 sys_size = b_idx%end
728 end if
729
730 ! Volume Fraction Model (6-equation model)
731 else if (model_eqns == 3) then
732
733 ! Annotating structure of the state and flux vectors belonging
734 ! to the system of equations defined by the selected number of
735 ! spatial dimensions and the volume fraction model
736 cont_idx%beg = 1
737 cont_idx%end = num_fluids
738 mom_idx%beg = cont_idx%end + 1
739 mom_idx%end = cont_idx%end + num_vels
740 e_idx = mom_idx%end + 1
741 adv_idx%beg = e_idx + 1
742 adv_idx%end = e_idx + num_fluids
743 internalenergies_idx%beg = adv_idx%end + 1
746 alf_idx = 1 ! dummy, cannot actually have a void fraction
747
748 else if (model_eqns == 4) then
749 cont_idx%beg = 1 ! one continuity equation
750 cont_idx%end = 1 !num_fluids
751 mom_idx%beg = cont_idx%end + 1 ! one momentum equation in each
752 mom_idx%end = cont_idx%end + num_vels
753 e_idx = mom_idx%end + 1 ! one energy equation
754 adv_idx%beg = e_idx + 1
755 adv_idx%end = adv_idx%beg !one volume advection equation
756 alf_idx = adv_idx%end
757 sys_size = alf_idx !adv_idx%end
758
759 if (bubbles_euler) then
760 bub_idx%beg = sys_size + 1
761 bub_idx%end = sys_size + 2*nb
762 if (polytropic .neqv. .true.) then
763 bub_idx%end = sys_size + 4*nb
764 end if
765 sys_size = bub_idx%end
766
767 allocate (bub_idx%rs(nb), bub_idx%vs(nb))
768 allocate (bub_idx%ps(nb), bub_idx%ms(nb))
769 allocate (weight(nb), r0(nb))
770
771 do i = 1, nb
772 if (polytropic .neqv. .true.) then
773 fac = 4
774 else
775 fac = 2
776 end if
777
778 bub_idx%rs(i) = bub_idx%beg + (i - 1)*fac
779 bub_idx%vs(i) = bub_idx%rs(i) + 1
780
781 if (polytropic .neqv. .true.) then
782 bub_idx%ps(i) = bub_idx%vs(i) + 1
783 bub_idx%ms(i) = bub_idx%ps(i) + 1
784 end if
785 end do
786
787 if (nb == 1) then
788 weight(:) = 1._wp
789 r0(:) = 1._wp
790 else if (nb < 1) then
791 stop 'Invalid value of nb'
792 end if
793
794 if (polytropic) then
795 rhoref = 1._wp
796 pref = 1._wp
797 end if
798 end if
799 end if
800
801 if (model_eqns == 2 .or. model_eqns == 3) then
802
803 if (hypoelasticity .or. hyperelasticity) then
804 elasticity = .true.
805 stress_idx%beg = sys_size + 1
806 stress_idx%end = sys_size + (num_dims*(num_dims + 1))/2
807 if (cyl_coord) stress_idx%end = stress_idx%end + 1
808 ! number of stresses is 1 in 1D, 3 in 2D, 4 in 2D-Axisym, 6 in 3D
809 sys_size = stress_idx%end
810
811 ! shear stress index is 2 for 2D and 2,4,5 for 3D
812 if (num_dims == 1) then
813 shear_num = 0
814 else if (num_dims == 2) then
815 shear_num = 1
816 shear_indices(1) = stress_idx%beg - 1 + 2
819 ! Both x-dir and y-dir: flip tau_xy only
820 else if (num_dims == 3) then
821 shear_num = 3
822 shear_indices(1:3) = stress_idx%beg - 1 + (/2, 4, 5/)
824 shear_bc_flip_indices(1, 1:2) = shear_indices((/1, 2/))
825 shear_bc_flip_indices(2, 1:2) = shear_indices((/1, 3/))
826 shear_bc_flip_indices(3, 1:2) = shear_indices((/2, 3/))
827 ! x-dir: flip tau_xy and tau_xz
828 ! y-dir: flip tau_xy and tau_yz
829 ! z-dir: flip tau_xz and tau_yz
830 end if
831 end if
832
833 if (hyperelasticity) then
834 xi_idx%beg = sys_size + 1
835 xi_idx%end = sys_size + num_dims
836 ! adding three more equations for the \xi field and the elastic energy
837 sys_size = xi_idx%end + 1
838 ! number of entries in the symmetric btensor plus the jacobian
839 b_size = (num_dims*(num_dims + 1))/2 + 1
840 tensor_size = num_dims**2 + 1
841 end if
842
843 if (surface_tension) then
844 c_idx = sys_size + 1
846 end if
847
848 if (cont_damage) then
849 damage_idx = sys_size + 1
851 else
852 damage_idx = dflt_int
853 end if
854
855 if (hyper_cleaning) then
856 psi_idx = sys_size + 1
858 else
859 psi_idx = dflt_int
860 end if
861
862 end if
863
864 if (chemistry) then
865 species_idx%beg = sys_size + 1
866 species_idx%end = sys_size + num_species
868 else
869 species_idx%beg = 1
870 species_idx%end = 1
871 end if
872
873 if (output_partial_domain) then
874 x_output_idx%beg = 0
875 x_output_idx%end = 0
876 y_output_idx%beg = 0
877 y_output_idx%end = 0
878 z_output_idx%beg = 0
879 z_output_idx%end = 0
880 end if
881
882 momxb = mom_idx%beg
883 momxe = mom_idx%end
884 advxb = adv_idx%beg
885 advxe = adv_idx%end
886 contxb = cont_idx%beg
887 contxe = cont_idx%end
888 bubxb = bub_idx%beg
889 bubxe = bub_idx%end
890 strxb = stress_idx%beg
891 strxe = stress_idx%end
894 xibeg = xi_idx%beg
895 xiend = xi_idx%end
896 chemxb = species_idx%beg
897 chemxe = species_idx%end
898
899#ifdef MFC_MPI
900 allocate (mpi_io_data%view(1:sys_size))
901 allocate (mpi_io_data%var(1:sys_size))
902 do i = 1, sys_size
903 if (down_sample) then
904 allocate (mpi_io_data%var(i)%sf(-1:m + 1, -1:n + 1, -1:p + 1))
905 else
906 allocate (mpi_io_data%var(i)%sf(0:m, 0:n, 0:p))
907 end if
908 mpi_io_data%var(i)%sf => null()
909 end do
910
911 if (ib) allocate (mpi_io_ib_data%var%sf(0:m, 0:n, 0:p))
912#endif
913
914 ! Size of the ghost zone layer is non-zero only when post-processing
915 ! the raw simulation data of a parallel multidimensional computation
916 ! in the Silo-HDF5 format. If this is the case, one must also verify
917 ! whether the raw simulation data is 2D or 3D. In the 2D case, size
918 ! of the z-coordinate direction ghost zone layer must be zeroed out.
919 if (num_procs == 1 .or. format /= 1) then
920
921 offset_x%beg = 0
922 offset_x%end = 0
923 offset_y%beg = 0
924 offset_y%end = 0
925 offset_z%beg = 0
926 offset_z%end = 0
927
928 elseif (n == 0) then
929
930 offset_y%beg = 0
931 offset_y%end = 0
932 offset_z%beg = 0
933 offset_z%end = 0
934
935 elseif (p == 0) then
936
937 offset_z%beg = 0
938 offset_z%end = 0
939
940 end if
941
942 ! Determining the finite-difference number and the buffer size. Note
943 ! that the size of the buffer is unrelated to the order of the WENO
944 ! scheme. Rather, it is directly dependent on maximum size of ghost
945 ! zone layers and possibly the order of the finite difference scheme
946 ! used for the computation of vorticity and/or numerical Schlieren
947 ! function.
948 buff_size = max(offset_x%beg, offset_x%end, offset_y%beg, &
949 offset_y%end, offset_z%beg, offset_z%end)
950
951 if (any(omega_wrt) .or. schlieren_wrt .or. qm_wrt .or. liutex_wrt) then
952 fd_number = max(1, fd_order/2)
954 end if
955
956 ! Configuring Coordinate Direction Indexes
957 idwint(1)%beg = 0; idwint(2)%beg = 0; idwint(3)%beg = 0
958 idwint(1)%end = m; idwint(2)%end = n; idwint(3)%end = p
959
960 idwbuff(1)%beg = -buff_size
961 if (num_dims > 1) then; idwbuff(2)%beg = -buff_size; else; idwbuff(2)%beg = 0; end if
962 if (num_dims > 2) then; idwbuff(3)%beg = -buff_size; else; idwbuff(3)%beg = 0; end if
963
964 idwbuff(1)%end = idwint(1)%end - idwbuff(1)%beg
965 idwbuff(2)%end = idwint(2)%end - idwbuff(2)%beg
966 idwbuff(3)%end = idwint(3)%end - idwbuff(3)%beg
967
968 ! Allocating single precision grid variables if needed
969 allocate (x_cc_s(-buff_size:m + buff_size))
970
971 ! Allocating the grid variables in the x-coordinate direction
972 allocate (x_cb(-1 - offset_x%beg:m + offset_x%end))
973 allocate (x_cc(-buff_size:m + buff_size))
974 allocate (dx(-buff_size:m + buff_size))
975
976 ! Allocating grid variables in the y- and z-coordinate directions
977 if (n > 0) then
978
979 allocate (y_cb(-1 - offset_y%beg:n + offset_y%end))
980 allocate (y_cc(-buff_size:n + buff_size))
981 allocate (dy(-buff_size:n + buff_size))
982
983 if (p > 0) then
984 allocate (z_cb(-1 - offset_z%beg:p + offset_z%end))
985 allocate (z_cc(-buff_size:p + buff_size))
986 allocate (dz(-buff_size:p + buff_size))
987 end if
988
989 ! Allocating the grid variables, only used for the 1D simulations,
990 ! and containing the defragmented computational domain grid data
991 else
992
993 allocate (x_root_cb(-1:m_root))
994 allocate (x_root_cc(0:m_root))
995
996 if (precision == 1) then
997 allocate (x_root_cc_s(0:m_root))
998 end if
999
1000 end if
1001
1002 allocate (adv(num_fluids))
1003
1004 if (cyl_coord .neqv. .true.) then ! Cartesian grid
1005 grid_geometry = 1
1006 elseif (cyl_coord .and. p == 0) then ! Axisymmetric cylindrical grid
1007 grid_geometry = 2
1008 else ! Fully 3D cylindrical grid
1009 grid_geometry = 3
1010 end if
1011
1013
1014 !> Subroutine to initialize parallel infrastructure
1016
1017#ifdef MFC_MPI
1018 integer :: ierr !< Generic flag used to identify and report MPI errors
1019#endif
1020
1021 num_dims = 1 + min(1, n) + min(1, p)
1022
1023 if (mhd) then
1024 num_vels = 3
1025 else
1027 end if
1028
1029 allocate (proc_coords(1:num_dims))
1030
1031 if (parallel_io .neqv. .true.) return
1032
1033#ifdef MFC_MPI
1034
1035 ! Option for Lustre file system (Darter/Comet/Stampede)
1036 write (mpiiofs, '(A)') '/lustre_'
1037 mpiiofs = trim(mpiiofs)
1038 call mpi_info_create(mpi_info_int, ierr)
1039 call mpi_info_set(mpi_info_int, 'romio_ds_write', 'disable', ierr)
1040
1041 ! Option for UNIX file system (Hooke/Thomson)
1042 ! WRITE(mpiiofs, '(A)') '/ufs_'
1043 ! mpiiofs = TRIM(mpiiofs)
1044 ! mpi_info_int = MPI_INFO_NULL
1045
1046 allocate (start_idx(1:num_dims))
1047
1048#endif
1049
1050 end subroutine s_initialize_parallel_io
1051
1052 !> Deallocation procedures for the module
1054
1055 integer :: i
1056
1057 ! Deallocating the grid variables for the x-coordinate direction
1058 deallocate (x_cc, x_cb, dx)
1059
1060 ! Deallocating grid variables for the y- and z-coordinate directions
1061 if (n > 0) then
1062 deallocate (y_cc, y_cb, dy)
1063 if (p > 0) then
1064 deallocate (z_cc, z_cb, dz)
1065 end if
1066 else
1067 ! Deallocating the grid variables, only used for the 1D simulations,
1068 ! and containing the defragmented computational domain grid data
1069 deallocate (x_root_cb, x_root_cc)
1070 end if
1071
1072 deallocate (proc_coords)
1073
1074 deallocate (adv)
1075
1076#ifdef MFC_MPI
1077
1078 if (parallel_io) then
1079 deallocate (start_idx)
1080 do i = 1, sys_size
1081 mpi_io_data%var(i)%sf => null()
1082 end do
1083
1084 deallocate (mpi_io_data%var)
1085 deallocate (mpi_io_data%view)
1086 end if
1087
1088 if (ib) mpi_io_ib_data%var%sf => null()
1089#endif
1090
1092
1093end module m_global_parameters
integer, intent(in) j
Shared derived types for field data, patch geometry, bubble dynamics, and MPI I/O structures.
Global parameters for the post-process: domain geometry, equation of state, and output database setti...
real(sp), dimension(:), allocatable x_cc_s
type(mpi_io_levelset_norm_var), public mpi_io_levelsetnorm_data
type(int_bounds_info) xi_idx
Indexes of first and last reference map eqns.
logical cont_damage
Continuum damage modeling.
integer, dimension(3, 2) shear_bc_flip_indices
Indices of shear stress components to reflect for boundary conditions. Size: (1:3,...
logical hypoelasticity
Turn hypoelasticity on.
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
Deallocation procedures for the module.
integer thermal
1 = adiabatic, 2 = isotherm, 3 = transfer
integer avg_state
Average state evaluation method.
real(wp), dimension(:), allocatable im_trans_c
type(int_bounds_info), dimension(1:3) idwint
integer recon_type
Which type of reconstruction to use.
logical, parameter chemistry
Chemistry modeling.
integer beta_idx
Index of lagrange bubbles beta.
type(int_bounds_info) offset_y
type(int_bounds_info) mom_idx
Indexes of first & last momentum eqns.
real(wp), dimension(num_fluids_max) schlieren_alpha
Amplitude coefficients of the numerical Schlieren function that are used to adjust the intensity of n...
real(wp), dimension(:), allocatable pb0
integer num_fluids
Number of different fluids present in the flow.
real(wp), dimension(:), allocatable im_trans_t
logical, dimension(3) flux_wrt
real(wp), dimension(:), allocatable y_cc
integer format
Format of the database file(s).
type(int_bounds_info) stress_idx
Indices of elastic stresses.
integer proc_rank
Rank of the local processor.
real(wp), dimension(:), allocatable re_trans_t
logical mixture_err
Mixture error limiter.
logical output_partial_domain
Specify portion of domain to output for post-processing.
real(wp), dimension(:), allocatable adv
Advection variables.
integer n_idx
Index of number density.
real(wp), dimension(:), allocatable x_root_cc
type(int_bounds_info) z_output_idx
Indices of domain to output for post-processing.
type(mpi_io_ib_var), public mpi_io_ib_data
logical dummy
AMDFlang workaround: keep a dummy logical to avoid a compiler case-optimization bug when a parameter+...
integer, dimension(:), allocatable proc_coords
Processor coordinates in MPI_CART_COMM.
real(wp), dimension(:), allocatable y_cb
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.
real(wp), dimension(:), allocatable dz
real(wp), dimension(:, :), allocatable, public mpi_io_data_lg_bubbles
real(wp), dimension(:), allocatable weight
integer gamma_idx
Index of specific heat ratio func. eqn.
real(wp), dimension(:), allocatable k_v
integer muscl_order
Order of accuracy for the MUSCL reconstruction.
logical alt_soundspeed
Alternate sound speed.
integer relax_model
Phase change relaxation model.
logical, dimension(3) mom_wrt
type(int_bounds_info) cont_idx
Indexes of first & last continuity eqns.
real(wp), dimension(:), allocatable x_root_cb
integer fd_number
The finite-difference number is given by MAX(1, fd_order/2). Essentially, it is a measure of the half...
logical, dimension(num_fluids_max) alpha_wrt
logical, dimension(num_fluids_max) alpha_rho_wrt
type(mpi_io_levelset_var), public mpi_io_levelset_data
type(int_bounds_info) b_idx
Indexes of first and last magnetic field eqns.
logical, dimension(num_fluids_max) alpha_rho_e_wrt
integer tensor_size
Number of components in the nonsymmetric tensor.
type(int_bounds_info), dimension(1:3) idwbuff
integer model_eqns
Multicomponent flow model.
integer buff_size
Number of cells in buffer region. For the variables which feature a buffer region,...
integer precision
Floating point precision of the database file(s).
logical hyperelasticity
Turn hyperelasticity on.
real(wp), dimension(:), allocatable z_cb
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....
type(bounds_info) z_output
Portion of domain to output for post-processing.
integer num_dims
Number of spatial dimensions.
integer shear_bc_flip_num
Number of shear stress components to reflect for boundary conditions.
real(wp), dimension(:), allocatable r0
type(int_bounds_info) x_output_idx
impure subroutine s_initialize_global_parameters_module
Computation of parameters, allocation procedures, and/or any other tasks needed to properly setup the...
integer pi_inf_idx
Index of liquid stiffness func. eqn.
real(wp), dimension(:), allocatable x_cc
real(wp), dimension(:), allocatable k_g
type(chemistry_parameters) chem_params
integer fd_order
The order of the finite-difference (fd) approximations of the first-order derivatives that need to be...
integer num_vels
Number of velocity components (different from num_dims for mhd).
real(wp), dimension(:), allocatable x_cb
real(wp), dimension(:), allocatable dy
integer t_step_save
Interval between consecutive time-step directory.
type(int_bounds_info) offset_x
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.
logical hyper_cleaning
Hyperbolic cleaning for MHD.
type(int_bounds_info) internalenergies_idx
Indexes of first & last internal energy eqns.
real(wp), dimension(:), allocatable z_cc
real(wp) bx0
Constant magnetic field in the x-direction (1D).
integer b_size
Number of components in the b tensor.
logical, dimension(3) omega_wrt
integer num_procs
Number of processors.
character(len=path_len) case_dir
Case folder location.
type(int_bounds_info) y_output_idx
type(int_bounds_info) adv_idx
Indexes of first & last advection eqns.
integer weno_order
Order of accuracy for the WENO reconstruction.
type(int_bounds_info) offset_z
logical mhd
Magnetohydrodynamics.
integer, dimension(3) shear_indices
Indices of the stress components that represent shear stress.
logical parallel_io
Format of the data files.
integer e_idx
Index of energy equation.
type(cell_num_bounds) cells_bounds
logical down_sample
down sampling of the database file(s)
real(wp), dimension(:), allocatable mass_g0
logical file_per_process
output format
integer t_step_start
First time-step directory.
real(wp) wall_time_avg
Wall time measurements.
logical elasticity
elasticity modeling, true for hyper or hypo
type(mpi_io_var), public mpi_io_data
integer c_idx
Index of color function.
impure subroutine s_initialize_parallel_io
Subroutine to initialize parallel infrastructure.
logical mpp_lim
Maximum volume fraction limiter.
real(wp), dimension(:), allocatable omegan
integer igr_order
IGR reconstruction order.
integer psi_idx
Index of hyperbolic cleaning state variable for MHD.
real(wp), dimension(:), allocatable re_trans_c
logical, dimension(3) vel_wrt
type(subgrid_bubble_physical_parameters) bub_pp
type(int_bounds_info) species_idx
Indexes of first & last concentration eqns.
logical, dimension(1:num_species) chem_wrt_y
logical relativity
Relativity for RMHD.
real(sp), dimension(:), allocatable x_root_cc_s
real(wp), dimension(:), allocatable dx
Cell-width distributions in the x-, y- and z-coordinate directions.
real(wp), dimension(:), allocatable pe_t
integer alf_idx
Index of specific heat ratio func. eqn.
real(wp), dimension(:), allocatable mass_v0
integer num_ibs
Number of immersed boundaries.
integer t_step_stop
Last time-step directory.
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.
Derived type adding beginning (beg) and end bounds info as attributes.
bounds for the bubble dynamic variables
Max and min number of cells in a direction of each combination of x-,y-, and z-.
Integer bounds for variables.
Derived type annexing the physical parameters (PP) of the fluids. These include the specific heat rat...
Derived type annexing the physical parameters required for sub-grid bubble models.