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