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/pre_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/pre_process/m_global_parameters.fpp" 2
17
18!> @brief Defines global parameters for the computational domain, simulation algorithm, and initial conditions
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 use m_helper_basic ! Functions to compare floating point numbers
27 use m_thermochem, only: num_species
28
29 implicit none
30
31 ! Logistics
32 integer :: num_procs !< Number of processors
33 character(LEN=path_len) :: case_dir !< Case folder location
34 logical :: old_grid !< Use existing grid data
35 logical :: old_ic, non_axis_sym !< Use existing IC data
36 integer :: t_step_old, t_step_start !< Existing IC/grid folder
38 integer :: n_start, n_start_old
39
40 ! Computational Domain Parameters
41
42 integer :: proc_rank !< Rank of the local processor Number of cells in the x-, y- and z-coordinate directions
43 integer :: m
44 integer :: n
45 integer :: p
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 !< Global number of cells in the domain
50 integer :: m_glb, n_glb, p_glb !< Global number of cells in each direction
51 integer :: num_dims !< Number of spatial dimensions
52 integer :: num_vels !< Number of velocity components (different from num_dims for mhd)
53 logical :: cyl_coord
54 integer :: grid_geometry !< Cylindrical coordinates (either axisymmetric or full 3D)
55 !> Locations of cell-centers (cc) in x-, y- and z-directions, respectively
56 real(wp), allocatable, dimension(:) :: x_cc, y_cc, z_cc
57 !> Locations of cell-boundaries (cb) in x-, y- and z-directions, respectively
58 real(wp), allocatable, dimension(:) :: x_cb, y_cb, z_cb
59 real(wp) :: dx, dy, dz !< Minimum cell-widths in the x-, y- and z-coordinate directions
60 type(bounds_info) :: x_domain, y_domain, z_domain !< Locations of the domain bounds in the x-, y- and z-coordinate directions
61 logical :: stretch_x, stretch_y, stretch_z !< Grid stretching flags for the x-, y- and z-coordinate directions
62 ! Grid stretching: a_x/a_y/a_z = rate, x_a/y_a/z_a = location
63 real(wp) :: a_x, a_y, a_z
64 integer :: loops_x, loops_y, loops_z
65 real(wp) :: x_a, y_a, z_a
66 real(wp) :: x_b, y_b, z_b
67
68 ! Simulation Algorithm Parameters
69 integer :: model_eqns !< Multicomponent flow model
70 logical :: relax !< activate phase change
71 integer :: relax_model !< Relax Model
72 real(wp) :: palpha_eps !< trigger parameter for the p relaxation procedure, phase change model
73 real(wp) :: ptgalpha_eps !< trigger parameter for the pTg relaxation procedure, phase change model
74 integer :: num_fluids !< Number of different fluids present in the flow
75 logical :: mpp_lim !< Alpha limiter
76 integer :: sys_size !< Number of unknowns in the system of equations
77 integer :: recon_type !< Reconstruction Type
78 integer :: weno_polyn !< Degree of the WENO polynomials (polyn)
79 integer :: muscl_polyn !< Degree of the MUSCL polynomials (polyn)
80 integer :: weno_order !< Order of accuracy for the WENO reconstruction
81 integer :: muscl_order !< Order of accuracy for the MUSCL reconstruction
82 logical :: hypoelasticity !< activate hypoelasticity
83 logical :: hyperelasticity !< activate hyperelasticity
84 logical :: elasticity !< elasticity modeling, true for hyper or hypo
85 logical :: mhd !< Magnetohydrodynamics
86 logical :: relativity !< Relativity for RMHD
87 integer :: b_size !< Number of components in the b tensor
88 integer :: tensor_size !< Number of components in the nonsymmetric tensor
89 logical :: pre_stress !< activate pre_stressed domain
90 logical :: cont_damage !< continuum damage modeling
91 logical :: hyper_cleaning !< Hyperbolic cleaning for MHD
92 logical :: igr !< Use information geometric regularization
93 integer :: igr_order !< IGR reconstruction order
94 logical, parameter :: chemistry = .false. !< Chemistry modeling
95 ! Annotations of the structure, i.e. the organization, of the state vectors
96 type(int_bounds_info) :: cont_idx !< Indexes of first & last continuity eqns.
97 type(int_bounds_info) :: mom_idx !< Indexes of first & last momentum eqns.
98 integer :: e_idx !< Index of total energy equation
99 integer :: alf_idx !< Index of void fraction
100 integer :: n_idx !< Index of number density
101 type(int_bounds_info) :: adv_idx !< Indexes of first & last advection eqns.
102 type(int_bounds_info) :: internalenergies_idx !< Indexes of first & last internal energy eqns.
103 type(bub_bounds_info) :: bub_idx !< Indexes of first & last bubble variable eqns.
104 integer :: gamma_idx !< Index of specific heat ratio func. eqn.
105 integer :: pi_inf_idx !< Index of liquid stiffness func. eqn.
106 type(int_bounds_info) :: b_idx !< Indexes of first and last magnetic field eqns.
107 type(int_bounds_info) :: stress_idx !< Indexes of elastic shear stress eqns.
108 type(int_bounds_info) :: xi_idx !< Indexes of first and last reference map eqns.
109 integer :: c_idx !< Index of the color function
110 type(int_bounds_info) :: species_idx !< Indexes of first & last concentration eqns.
111 integer :: damage_idx !< Index of damage state variable (D) for continuum damage model
112 integer :: psi_idx !< Index of hyperbolic cleaning state variable for MHD
113 ! Cell Indices for the (local) interior points (O-m, O-n, 0-p). Stands for "InDices With BUFFer".
115
116 ! Cell indices (InDices With BUFFer): includes buffer except in pre_process
118 type(int_bounds_info) :: bc_x, bc_y, bc_z !< Boundary conditions in the x-, y- and z-coordinate directions
119 integer :: shear_num !< Number of shear stress components
120 integer, dimension(3) :: shear_indices !< Indices of the stress components that represent shear stress
121 integer :: shear_bc_flip_num !< Number of shear stress components to reflect for boundary conditions
122 integer, dimension(3, 2) :: shear_bc_flip_indices !< Shear stress BC reflection indices (1:3, 1:shear_BC_flip_num)
123 logical :: parallel_io !< Format of the data files
124 logical :: file_per_process !< type of data output
125 integer :: precision !< Precision of output files
126 logical :: down_sample !< Down-sample the output data
127 logical :: mixlayer_vel_profile !< Set hyperbolic tangent streamwise velocity profile
128 real(wp) :: mixlayer_vel_coef !< Coefficient for the hyperbolic tangent streamwise velocity profile
129 logical :: mixlayer_perturb !< Superimpose instability waves to surrounding fluid flow
130 integer :: mixlayer_perturb_nk !< Number of Fourier modes for perturbation with mixlayer_perturb flag
131 real(wp) :: mixlayer_perturb_k0 !< Peak wavenumber for mixlayer perturbation (default: most unstable mode)
134 real(wp) :: pi_fac !< Factor for artificial pi_inf
135 logical :: viscous
137
138 ! Perturb density of surrounding air so as to break symmetry of grid
139 logical :: perturb_flow
140 integer :: perturb_flow_fluid !< Fluid to be perturbed with perturb_flow flag
141 real(wp) :: perturb_flow_mag !< Magnitude of perturbation with perturb_flow flag
142 logical :: perturb_sph
143 integer :: perturb_sph_fluid !< Fluid to be perturbed with perturb_sph flag
144 real(wp), dimension(num_fluids_max) :: fluid_rho
147 integer, allocatable, dimension(:) :: proc_coords !< Processor coordinates in MPI_CART_COMM
148 integer, allocatable, dimension(:) :: start_idx !< Starting cell-center index of local processor in global grid
149#ifdef MFC_MPI
150 type(mpi_io_var), public :: mpi_io_data
151 character(LEN=name_len) :: mpiiofs
152 integer :: mpi_info_int !< MPI info for parallel IO with Lustre file systems
153#endif
154
155 ! Initial Condition Parameters
156 integer :: num_patches !< Number of patches composing initial condition
157 type(ic_patch_parameters), dimension(num_patches_max) :: patch_icpp !< IC patch parameters (max: num_patches_max)
158 integer :: num_bc_patches !< Number of boundary condition patches
159 logical :: bc_io !< whether or not to save BC data
160 type(bc_patch_parameters), dimension(num_bc_patches_max) :: patch_bc !< Boundary condition patch parameters
161
162 ! Fluids Physical Parameters
163 type(physical_parameters), dimension(num_fluids_max) :: fluid_pp !< Stiffened gas EOS parameters and Reynolds numbers per fluid
164 ! Subgrid Bubble Parameters
166 real(wp) :: rhoref, pref !< Reference parameters for Tait EOS
168 !> @name Bubble modeling
169 !> @{
170 integer :: nb
171 real(wp) :: ca, web, re_inv, eu
172 real(wp), dimension(:), allocatable :: weight, r0
173 logical :: bubbles_euler
174 logical :: qbmm !< Quadrature moment method
175 integer :: nmom !< Number of carried moments
176 real(wp) :: sigr, sigv, rhorv !< standard deviations in R/V
177 logical :: adv_n !< Solve the number density equation and compute alpha from number density
178 !> @}
179
180 !> @name Immersed Boundaries
181 !> @{
182 logical :: ib !< Turn immersed boundaries on
183 integer :: num_ibs !< Number of immersed boundaries
184 integer :: np
185 type(ib_patch_parameters), dimension(num_patches_max) :: patch_ib !< Immersed boundary patch parameters
186 type(vec3_dt), allocatable, dimension(:) :: airfoil_grid_u, airfoil_grid_l
187 !> @}
188
189 !> @name Non-polytropic bubble gas compression
190 !> @{
191 logical :: polytropic
192 logical :: polydisperse
193 real(wp) :: poly_sigma
194 integer :: dist_type !< 1 = binormal, 2 = lognormal-normal
195 integer :: thermal !< 1 = adiabatic, 2 = isotherm, 3 = transfer
196 real(wp) :: phi_vg, phi_gv, pe_c, tw, k_vl, k_gl
197 real(wp) :: gam_m
198 real(wp), dimension(:), allocatable :: pb0, mass_g0, mass_v0, pe_t, k_v, k_g
199 real(wp), dimension(:), allocatable :: re_trans_t, re_trans_c, im_trans_t, im_trans_c, omegan
201 !> @}
202
203 !> @name Surface Tension Modeling
204 !> @{
205 real(wp) :: sigma
207 !> @}
208
209 !> @name Index variables used for m_variables_conversion
210 !> @{
211 integer :: momxb, momxe
212 integer :: advxb, advxe
213 integer :: contxb, contxe
214 integer :: intxb, intxe
215 integer :: bubxb, bubxe
216 integer :: strxb, strxe
217 integer :: xibeg, xiend
218 integer :: chemxb, chemxe
219 !> @}
220
221 integer, allocatable, dimension(:,:,:) :: logic_grid
222 type(pres_field) :: pb
223 type(pres_field) :: mv
224 real(wp) :: bx0 !< Constant magnetic field in the x-direction (1D)
225 integer :: buff_size !< Number of ghost cells for boundary condition storage
226 logical :: fft_wrt
227 logical :: dummy !< AMDFlang workaround for case-optimization + GPU-kernel bug
228
229contains
230
231 !> Assigns default values to user inputs prior to reading them in. This allows for an easier consistency check of these
232 !! parameters once they are read from the input file.
234
235 integer :: i !< Generic loop operator
236 ! Logistics
237
238 case_dir = '.'
239 old_grid = .false.
240 old_ic = .false.
241 t_step_old = dflt_int
242 t_step_start = dflt_int
243
244 cfl_adap_dt = .false.
245 cfl_const_dt = .false.
246 cfl_dt = .false.
247 n_start = dflt_int
248
249 ! Computational domain parameters
250 m = dflt_int; n = 0; p = 0
251
253
254 cyl_coord = .false.
255
256 x_domain%beg = dflt_real
257 x_domain%end = dflt_real
258 y_domain%beg = dflt_real
259 y_domain%end = dflt_real
260 z_domain%beg = dflt_real
261 z_domain%end = dflt_real
262
263 stretch_x = .false.
264 stretch_y = .false.
265 stretch_z = .false.
266
267 a_x = dflt_real
268 a_y = dflt_real
269 a_z = dflt_real
270 loops_x = 1
271 loops_y = 1
272 loops_z = 1
273 x_a = dflt_real
274 x_b = dflt_real
275 y_a = dflt_real
276 y_b = dflt_real
277 z_a = dflt_real
278 z_b = dflt_real
279
280 ! Simulation algorithm parameters
281 model_eqns = dflt_int
282 relax = .false.
283 relax_model = dflt_int
284 palpha_eps = dflt_real
285 ptgalpha_eps = dflt_real
286 num_fluids = dflt_int
287 recon_type = weno_type
288 weno_order = dflt_int
289 igr = .false.
290 igr_order = dflt_int
291 muscl_order = dflt_int
292
293 hypoelasticity = .false.
294 hyperelasticity = .false.
295 elasticity = .false.
296 pre_stress = .false.
297 b_size = dflt_int
298 tensor_size = dflt_int
299 cont_damage = .false.
300 hyper_cleaning = .false.
301
302 mhd = .false.
303 relativity = .false.
304
305 bc_x%beg = dflt_int; bc_x%end = dflt_int
306 bc_y%beg = dflt_int; bc_y%end = dflt_int
307 bc_z%beg = dflt_int; bc_z%end = dflt_int
308
309# 299 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
310# 300 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
311 bc_x%vb1 = 0._wp
312 bc_x%ve1 = 0._wp
313# 300 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
314 bc_x%vb2 = 0._wp
315 bc_x%ve2 = 0._wp
316# 300 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
317 bc_x%vb3 = 0._wp
318 bc_x%ve3 = 0._wp
319# 303 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
320# 299 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
321# 300 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
322 bc_y%vb1 = 0._wp
323 bc_y%ve1 = 0._wp
324# 300 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
325 bc_y%vb2 = 0._wp
326 bc_y%ve2 = 0._wp
327# 300 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
328 bc_y%vb3 = 0._wp
329 bc_y%ve3 = 0._wp
330# 303 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
331# 299 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
332# 300 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
333 bc_z%vb1 = 0._wp
334 bc_z%ve1 = 0._wp
335# 300 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
336 bc_z%vb2 = 0._wp
337 bc_z%ve2 = 0._wp
338# 300 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
339 bc_z%vb3 = 0._wp
340 bc_z%ve3 = 0._wp
341# 303 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
342# 304 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
343
344 parallel_io = .false.
345 file_per_process = .false.
346 precision = 2
347 down_sample = .false.
348 viscous = .false.
349 bubbles_lagrange = .false.
350 mixlayer_vel_profile = .false.
351 mixlayer_vel_coef = 1._wp
352 mixlayer_perturb = .false.
354 mixlayer_perturb_k0 = 0.4446_wp
355 perturb_flow = .false.
356 perturb_flow_fluid = dflt_int
357 perturb_flow_mag = dflt_real
358 perturb_sph = .false.
359 perturb_sph_fluid = dflt_int
360 fluid_rho = dflt_real
361 elliptic_smoothing_iters = dflt_int
362 elliptic_smoothing = .false.
363
364 fft_wrt = .false.
365 dummy = .false.
366
367 simplex_perturb = .false.
368 simplex_params%perturb_vel(:) = .false.
369 simplex_params%perturb_vel_freq(:) = dflt_real
370 simplex_params%perturb_vel_scale(:) = dflt_real
371 simplex_params%perturb_vel_offset(:,:) = dflt_real
372 simplex_params%perturb_dens(:) = .false.
373 simplex_params%perturb_dens_freq(:) = dflt_real
374 simplex_params%perturb_dens_scale(:) = dflt_real
375 simplex_params%perturb_dens_offset(:,:) = dflt_real
376
377 ! Initial condition parameters
378 num_patches = dflt_int
379
380 do i = 1, num_patches_max
381 patch_icpp(i)%geometry = dflt_int
382 patch_icpp(i)%model_scale(:) = 1._wp
383 patch_icpp(i)%model_translate(:) = 0._wp
384 patch_icpp(i)%model_filepath(:) = dflt_char
385 patch_icpp(i)%model_spc = num_ray
386 patch_icpp(i)%model_threshold = ray_tracing_threshold
387 patch_icpp(i)%x_centroid = dflt_real
388 patch_icpp(i)%y_centroid = dflt_real
389 patch_icpp(i)%z_centroid = dflt_real
390 patch_icpp(i)%length_x = dflt_real
391 patch_icpp(i)%length_y = dflt_real
392 patch_icpp(i)%length_z = dflt_real
393 patch_icpp(i)%radius = dflt_real
394 patch_icpp(i)%epsilon = dflt_real
395 patch_icpp(i)%beta = dflt_real
396 patch_icpp(i)%normal = dflt_real
397 patch_icpp(i)%radii = dflt_real
398 patch_icpp(i)%alter_patch = .false.
399 patch_icpp(i)%alter_patch(0) = .true.
400 patch_icpp(i)%smoothen = .false.
401 patch_icpp(i)%smooth_patch_id = i
402 patch_icpp(i)%smooth_coeff = dflt_real
403 patch_icpp(i)%alpha_rho = dflt_real
404 patch_icpp(i)%rho = dflt_real
405 patch_icpp(i)%vel = dflt_real
406 patch_icpp(i)%pres = dflt_real
407 patch_icpp(i)%alpha = dflt_real
408 patch_icpp(i)%gamma = dflt_real
409 patch_icpp(i)%pi_inf = dflt_real
410 patch_icpp(i)%cv = 0._wp
411 patch_icpp(i)%qv = 0._wp
412 patch_icpp(i)%qvp = 0._wp
413 patch_icpp(i)%tau_e = 0._wp
414 patch_icpp(i)%Bx = dflt_real
415 patch_icpp(i)%By = dflt_real
416 patch_icpp(i)%Bz = dflt_real
417 patch_icpp(i)%a(2) = dflt_real
418 patch_icpp(i)%a(3) = dflt_real
419 patch_icpp(i)%a(4) = dflt_real
420 patch_icpp(i)%a(5) = dflt_real
421 patch_icpp(i)%a(6) = dflt_real
422 patch_icpp(i)%a(7) = dflt_real
423 patch_icpp(i)%a(8) = dflt_real
424 patch_icpp(i)%a(9) = dflt_real
425 patch_icpp(i)%non_axis_sym = .false.
426 patch_icpp(i)%fourier_cos(:) = 0._wp
427 patch_icpp(i)%fourier_sin(:) = 0._wp
428 patch_icpp(i)%modal_clip_r_to_min = .false.
429 patch_icpp(i)%modal_r_min = 1.e-12_wp
430 patch_icpp(i)%modal_use_exp_form = .false.
431 patch_icpp(i)%sph_har_coeff(:,:) = 0._wp
432
433 ! should get all of r0's and v0's
434 patch_icpp(i)%r0 = dflt_real
435 patch_icpp(i)%v0 = dflt_real
436
437 patch_icpp(i)%p0 = dflt_real
438 patch_icpp(i)%m0 = dflt_real
439
440 patch_icpp(i)%hcid = dflt_int
441
442 if (chemistry) then
443 patch_icpp(i)%Y(:) = 0._wp
444 end if
445 end do
446
448 bc_io = .false.
449
450 do i = 1, num_bc_patches_max
451 patch_bc(i)%geometry = dflt_int
452 patch_bc(i)%type = dflt_int
453 patch_bc(i)%dir = dflt_int
454 patch_bc(i)%loc = dflt_int
455 patch_bc(i)%centroid(:) = dflt_real
456 patch_bc(i)%length(:) = dflt_real
457 patch_bc(i)%radius = dflt_real
458 end do
459
460 ! Tait EOS
461 rhoref = dflt_real
462 pref = dflt_real
463
464 ! Bubble modeling
465 bubbles_euler = .false.
466 polytropic = .true.
467 polydisperse = .false.
468
469 thermal = dflt_int
470 r0ref = dflt_real
471 nb = dflt_int
472
473 eu = dflt_real
474 ca = dflt_real
475 re_inv = dflt_real
476 web = dflt_real
477 poly_sigma = dflt_real
478 surface_tension = .false.
479
480 adv_n = .false.
481
482 qbmm = .false.
483 nmom = 1
484 sigr = dflt_real
485 sigv = dflt_real
486 rhorv = 0._wp
487 dist_type = dflt_int
488
489 r_g = dflt_real
490 r_v = dflt_real
491 phi_vg = dflt_real
492 phi_gv = dflt_real
493 pe_c = dflt_real
494 tw = dflt_real
495
496 ! surface tension modeling
497 sigma = dflt_real
498 pi_fac = 1._wp
499
500 ! Immersed Boundaries
501 ib = .false.
502 num_ibs = dflt_int
503
504 do i = 1, num_patches_max
505 patch_ib(i)%geometry = dflt_int
506 patch_ib(i)%x_centroid = dflt_real
507 patch_ib(i)%y_centroid = dflt_real
508 patch_ib(i)%z_centroid = dflt_real
509 patch_ib(i)%length_x = dflt_real
510 patch_ib(i)%length_y = dflt_real
511 patch_ib(i)%length_z = dflt_real
512 patch_ib(i)%radius = dflt_real
513 patch_ib(i)%theta = dflt_real
514 patch_ib(i)%c = dflt_real
515 patch_ib(i)%t = dflt_real
516 patch_ib(i)%m = dflt_real
517 patch_ib(i)%p = dflt_real
518 patch_ib(i)%slip = .false.
519
520 ! Proper default values for translating STL models
521 patch_ib(i)%model_scale(:) = 1._wp
522 patch_ib(i)%model_translate(:) = 0._wp
523 patch_ib(i)%model_rotate(:) = 0._wp
524 patch_ib(i)%model_filepath(:) = dflt_char
525 patch_ib(i)%model_spc = num_ray
526 patch_ib(i)%model_threshold = ray_tracing_threshold
527
528 ! Variables to handle moving immersed boundaries, defaulting to no movement
529 patch_ib(i)%moving_ibm = 0
530 patch_ib(i)%vel(:) = 0._wp
531 patch_ib(i)%angles(:) = 0._wp
532 patch_ib(i)%angular_vel(:) = 0._wp
533 patch_ib(i)%mass = dflt_real
534 patch_ib(i)%moment = dflt_real
535 patch_ib(i)%centroid_offset(:) = 0._wp
536
537 ! sets values of a rotation matrix which can be used when calculating rotations
538 patch_ib(i)%rotation_matrix = 0._wp
539 patch_ib(i)%rotation_matrix(1, 1) = 1._wp
540 patch_ib(i)%rotation_matrix(2, 2) = 1._wp
541 patch_ib(i)%rotation_matrix(3, 3) = 1._wp
542 patch_ib(i)%rotation_matrix_inverse = patch_ib(i)%rotation_matrix
543 end do
544
545 chem_params%gamma_method = 1
546 chem_params%transport_model = 1
547
548 ! Fluids physical parameters
549 do i = 1, num_fluids_max
550 fluid_pp(i)%gamma = dflt_real
551 fluid_pp(i)%pi_inf = dflt_real
552 fluid_pp(i)%cv = 0._wp
553 fluid_pp(i)%qv = 0._wp
554 fluid_pp(i)%qvp = 0._wp
555 fluid_pp(i)%G = 0._wp
556 end do
557
558 bx0 = dflt_real
559
560 ! Subgrid bubble parameters
561 bub_pp%R0ref = dflt_real; r0ref = dflt_real
562 bub_pp%p0ref = dflt_real; p0ref = dflt_real
563 bub_pp%rho0ref = dflt_real; rho0ref = dflt_real
564 bub_pp%T0ref = dflt_real; t0ref = dflt_real
565 bub_pp%ss = dflt_real; ss = dflt_real
566 bub_pp%pv = dflt_real; pv = dflt_real
567 bub_pp%vd = dflt_real; vd = dflt_real
568 bub_pp%mu_l = dflt_real; mu_l = dflt_real
569 bub_pp%mu_v = dflt_real; mu_v = dflt_real
570 bub_pp%mu_g = dflt_real; mu_g = dflt_real
571 bub_pp%gam_v = dflt_real; gam_v = dflt_real
572 bub_pp%gam_g = dflt_real; gam_g = dflt_real
573 bub_pp%M_v = dflt_real; m_v = dflt_real
574 bub_pp%M_g = dflt_real; m_g = dflt_real
575 bub_pp%k_v = dflt_real
576 bub_pp%k_g = dflt_real
577 bub_pp%cp_v = dflt_real; cp_v = dflt_real
578 bub_pp%cp_g = dflt_real; cp_g = dflt_real
579 bub_pp%R_v = dflt_real; r_v = dflt_real
580 bub_pp%R_g = dflt_real; r_g = dflt_real
581
583
584 !> Computation of parameters, allocation procedures, and/or any other tasks needed to properly setup the module
586
587 integer :: i, j, fac
588
589 if (recon_type == weno_type) then
590 weno_polyn = (weno_order - 1)/2
591 else if (recon_type == muscl_type) then
593 end if
594
595 ! Determining the layout of the state vectors and overall size of the system of equations, given the dimensionality and
596 ! choice of the equations of motion
597
598 ! Gamma/Pi_inf Model
599 if (model_eqns == 1) then
600 ! Setting number of fluids
601 num_fluids = 1
602
603 ! Annotating structure of the state and flux vectors belonging to the system of equations defined by the selected number
604 ! of spatial dimensions and the gamma/pi_inf model
605 cont_idx%beg = 1
606 cont_idx%end = cont_idx%beg
607 mom_idx%beg = cont_idx%end + 1
608 mom_idx%end = cont_idx%end + num_vels
609 e_idx = mom_idx%end + 1
610 adv_idx%beg = e_idx + 1
611 adv_idx%end = adv_idx%beg + 1
612 gamma_idx = adv_idx%beg
613 pi_inf_idx = adv_idx%end
614 sys_size = adv_idx%end
615
616 ! Volume Fraction Model (5-equation model)
617 else if (model_eqns == 2) then
618 ! Annotating structure of the state and flux vectors belonging to the system of equations defined by the selected number
619 ! of spatial dimensions and the volume fraction model
620 cont_idx%beg = 1
621 cont_idx%end = num_fluids
622 mom_idx%beg = cont_idx%end + 1
623 mom_idx%end = cont_idx%end + num_vels
624 e_idx = mom_idx%end + 1
625
626 if (igr) then
627 ! Volume fractions are stored in the indices immediately following the energy equation. IGR tracks a total of (N-1)
628 ! volume fractions for N fluids, hence the "-1" in adv_idx%end. If num_fluids = 1 then adv_idx%end < adv_idx%beg,
629 ! which skips all loops over the volume fractions since there is no volume fraction to track
630 adv_idx%beg = e_idx + 1
631 adv_idx%end = e_idx + num_fluids - 1
632 else
633 ! Volume fractions are stored in the indices immediately following the energy equation. WENO/MUSCL + Riemann tracks
634 ! a total of (N) volume fractions for N fluids, hence the lack of "-1" in adv_idx%end
635 adv_idx%beg = e_idx + 1
636 adv_idx%end = e_idx + num_fluids
637 end if
638
639 sys_size = adv_idx%end
640
641 if (bubbles_euler) then
642 alf_idx = adv_idx%end
643 else
644 alf_idx = 1
645 end if
646
647 if (bubbles_euler) then
648 bub_idx%beg = sys_size + 1
649 if (qbmm) then
650 if (nnode == 4) then
651 nmom = 6 !< Already set as a parameter
652 end if
653 bub_idx%end = adv_idx%end + nb*nmom
654 else
655 if (.not. polytropic) then
656 bub_idx%end = sys_size + 4*nb
657 else
658 bub_idx%end = sys_size + 2*nb
659 end if
660 end if
661 sys_size = bub_idx%end
662
663 if (adv_n) then
664 n_idx = bub_idx%end + 1
666 end if
667
668 allocate (bub_idx%rs(nb), bub_idx%vs(nb))
669 allocate (bub_idx%ps(nb), bub_idx%ms(nb))
670
671 if (qbmm) then
672 allocate (bub_idx%moms(nb, nmom))
673 allocate (bub_idx%fullmom(nb,0:nmom,0:nmom))
674
675 do i = 1, nb
676 do j = 1, nmom
677 bub_idx%moms(i, j) = bub_idx%beg + (j - 1) + (i - 1)*nmom
678 end do
679 bub_idx%fullmom(i, 0, 0) = bub_idx%moms(i, 1)
680 bub_idx%fullmom(i, 1, 0) = bub_idx%moms(i, 2)
681 bub_idx%fullmom(i, 0, 1) = bub_idx%moms(i, 3)
682 bub_idx%fullmom(i, 2, 0) = bub_idx%moms(i, 4)
683 bub_idx%fullmom(i, 1, 1) = bub_idx%moms(i, 5)
684 bub_idx%fullmom(i, 0, 2) = bub_idx%moms(i, 6)
685 bub_idx%rs(i) = bub_idx%fullmom(i, 1, 0)
686 end do
687 else
688 do i = 1, nb
689 if (.not. polytropic) then
690 fac = 4
691 else
692 fac = 2
693 end if
694
695 bub_idx%rs(i) = bub_idx%beg + (i - 1)*fac
696 bub_idx%vs(i) = bub_idx%rs(i) + 1
697
698 if (.not. polytropic) then
699 bub_idx%ps(i) = bub_idx%vs(i) + 1
700 bub_idx%ms(i) = bub_idx%ps(i) + 1
701 end if
702 end do
703 end if
704 end if
705
706 if (mhd) then
707 b_idx%beg = sys_size + 1
708 if (n == 0) then
709 b_idx%end = sys_size + 2 ! 1D: By, Bz
710 else
711 b_idx%end = sys_size + 3 ! 2D/3D: Bx, By, Bz
712 end if
713 sys_size = b_idx%end
714 end if
715
716 ! Volume Fraction Model (6-equation model)
717 else if (model_eqns == 3) then
718 ! Annotating structure of the state and flux vectors belonging to the system of equations defined by the selected number
719 ! of spatial dimensions and the volume fraction model
720 cont_idx%beg = 1
721 cont_idx%end = num_fluids
722 mom_idx%beg = cont_idx%end + 1
723 mom_idx%end = cont_idx%end + num_vels
724 e_idx = mom_idx%end + 1
725 adv_idx%beg = e_idx + 1
726 adv_idx%end = e_idx + num_fluids
727 internalenergies_idx%beg = adv_idx%end + 1
730 else if (model_eqns == 4) then
731 ! 4 equation model with subgrid bubbles_euler
732 cont_idx%beg = 1 ! one continuity equation
733 cont_idx%end = 1 ! num_fluids
734 mom_idx%beg = cont_idx%end + 1 ! one momentum equation in each direction
735 mom_idx%end = cont_idx%end + num_vels
736 e_idx = mom_idx%end + 1 ! one energy equation
737 adv_idx%beg = e_idx + 1
738 adv_idx%end = adv_idx%beg ! one volume advection equation
739 alf_idx = adv_idx%end
740 sys_size = alf_idx ! adv_idx%end
741
742 if (bubbles_euler) then
743 bub_idx%beg = sys_size + 1
744 bub_idx%end = sys_size + 2*nb
745 if (.not. polytropic) then
746 bub_idx%end = sys_size + 4*nb
747 end if
748 sys_size = bub_idx%end
749
750 allocate (bub_idx%rs(nb), bub_idx%vs(nb))
751 allocate (bub_idx%ps(nb), bub_idx%ms(nb))
752 allocate (weight(nb), r0(nb))
753
754 do i = 1, nb
755 if (.not. polytropic) then
756 fac = 4
757 else
758 fac = 2
759 end if
760
761 bub_idx%rs(i) = bub_idx%beg + (i - 1)*fac
762 bub_idx%vs(i) = bub_idx%rs(i) + 1
763
764 if (.not. polytropic) then
765 bub_idx%ps(i) = bub_idx%vs(i) + 1
766 bub_idx%ms(i) = bub_idx%ps(i) + 1
767 end if
768 end do
769
770 if (nb == 1) then
771 weight(:) = 1._wp
772 r0(:) = 1._wp
773 else if (nb < 1) then
774 stop 'Invalid value of nb'
775 end if
776
777 if (polytropic) then
778 rhoref = 1._wp
779 pref = 1._wp
780 end if
781 end if
782 end if
783
784 if (model_eqns == 2 .or. model_eqns == 3) then
785 if (hypoelasticity .or. hyperelasticity) then
786 elasticity = .true.
787 stress_idx%beg = sys_size + 1
788 stress_idx%end = sys_size + (num_dims*(num_dims + 1))/2
789 if (cyl_coord) stress_idx%end = stress_idx%end + 1
790 ! number of stresses is 1 in 1D, 3 in 2D, 4 in 2D-Axisym, 6 in 3D
791 sys_size = stress_idx%end
792
793 ! shear stress index is 2 for 2D and 2,4,5 for 3D
794 if (num_dims == 1) then
795 shear_num = 0
796 else if (num_dims == 2) then
797 shear_num = 1
798 shear_indices(1) = stress_idx%beg - 1 + 2
801 ! Both x-dir and y-dir: flip tau_xy only
802 else if (num_dims == 3) then
803 shear_num = 3
804 shear_indices(1:3) = stress_idx%beg - 1 + (/2, 4, 5/)
806 shear_bc_flip_indices(1,1:2) = shear_indices((/1, 2/))
807 shear_bc_flip_indices(2,1:2) = shear_indices((/1, 3/))
808 shear_bc_flip_indices(3,1:2) = shear_indices((/2, 3/))
809 ! x-dir: flip tau_xy and tau_xz y-dir: flip tau_xy and tau_yz z-dir: flip tau_xz and tau_yz
810 end if
811 end if
812
813 if (hyperelasticity) then
814 ! number of entries in the symmetric btensor plus the jacobian
815 b_size = (num_dims*(num_dims + 1))/2 + 1
816 tensor_size = num_dims**2 + 1
817 xi_idx%beg = sys_size + 1
818 xi_idx%end = sys_size + num_dims
819 ! adding three more equations for the \xi field and the elastic energy
820 sys_size = xi_idx%end + 1
821 end if
822
823 if (surface_tension) then
824 c_idx = sys_size + 1
826 end if
827
828 if (cont_damage) then
829 damage_idx = sys_size + 1
831 end if
832
833 if (hyper_cleaning) then
834 psi_idx = sys_size + 1
836 end if
837 end if
838
839 if (chemistry) then
840 species_idx%beg = sys_size + 1
841 species_idx%end = sys_size + num_species
843 end if
844
845 momxb = mom_idx%beg
846 momxe = mom_idx%end
847 advxb = adv_idx%beg
848 advxe = adv_idx%end
849 contxb = cont_idx%beg
850 contxe = cont_idx%end
851 bubxb = bub_idx%beg
852 bubxe = bub_idx%end
853 strxb = stress_idx%beg
854 strxe = stress_idx%end
857 xibeg = xi_idx%beg
858 xiend = xi_idx%end
859 chemxb = species_idx%beg
860 chemxe = species_idx%end
861
864
865#ifdef MFC_MPI
866 if (qbmm .and. .not. polytropic) then
867 allocate (mpi_io_data%view(1:sys_size + 2*nb*nnode))
868 allocate (mpi_io_data%var(1:sys_size + 2*nb*nnode))
869 else
870 allocate (mpi_io_data%view(1:sys_size))
871 allocate (mpi_io_data%var(1:sys_size))
872 end if
873
874 if (.not. down_sample) then
875 do i = 1, sys_size
876 allocate (mpi_io_data%var(i)%sf(0:m,0:n,0:p))
877 mpi_io_data%var(i)%sf => null()
878 end do
879 end if
880 if (qbmm .and. .not. polytropic) then
881 do i = sys_size + 1, sys_size + 2*nb*nnode
882 allocate (mpi_io_data%var(i)%sf(0:m,0:n,0:p))
883 mpi_io_data%var(i)%sf => null()
884 end do
885 end if
886#endif
887
888 ! Allocating grid variables for the x-direction
889 allocate (x_cc(0:m), x_cb(-1:m))
890 ! Allocating grid variables for the y- and z-directions
891 if (n > 0) then
892 allocate (y_cc(0:n), y_cb(-1:n))
893 if (p > 0) then
894 allocate (z_cc(0:p), z_cb(-1:p))
895 end if
896 end if
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
906 if (.not. igr) then
907 allocate (logic_grid(0:m,0:n,0:p))
908 end if
909
911
912 !> Configure MPI parallel I/O settings and allocate processor coordinate arrays.
913 impure subroutine s_initialize_parallel_io
914
915#ifdef MFC_MPI
916 integer :: ierr !< Generic flag used to identify and report MPI errors
917#endif
918
919 num_dims = 1 + min(1, n) + min(1, p)
920
921 if (mhd) then
922 num_vels = 3
923 else
925 end if
926
927 allocate (proc_coords(1:num_dims))
928
929 if (parallel_io .neqv. .true.) return
930
931#ifdef MFC_MPI
932 ! Option for Lustre file system (Darter/Comet/Stampede)
933 write (mpiiofs, '(A)') '/lustre_'
934 mpiiofs = trim(mpiiofs)
935 call mpi_info_create(mpi_info_int, ierr)
936 call mpi_info_set(mpi_info_int, 'romio_ds_write', 'disable', ierr)
937
938 ! Option for UNIX file system (Hooke/Thomson) WRITE(mpiiofs, '(A)') '/ufs_' mpiiofs = TRIM(mpiiofs) mpi_info_int =
939 ! MPI_INFO_NULL
940
941 allocate (start_idx(1:num_dims))
942#endif
943
944 end subroutine s_initialize_parallel_io
945
946 !> Deallocate all global grid, index, and equation-of-state parameter arrays.
948
949 integer :: i
950
951 ! Deallocating grid variables for the x-direction
952
953 deallocate (x_cc, x_cb)
954 ! Deallocating grid variables for the y- and z-directions
955 if (n > 0) then
956 deallocate (y_cc, y_cb)
957 if (p > 0) then
958 deallocate (z_cc, z_cb)
959 end if
960 end if
961
962 deallocate (proc_coords)
963
964#ifdef MFC_MPI
965 if (parallel_io) then
966 deallocate (start_idx)
967 do i = 1, sys_size
968 mpi_io_data%var(i)%sf => null()
969 end do
970
971 deallocate (mpi_io_data%var)
972 deallocate (mpi_io_data%view)
973 end if
974#endif
975
977
978end module m_global_parameters
integer, intent(in) j
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 for mixlayer perturbation (default: most unstable mode).
integer grid_geometry
Cylindrical coordinates (either axisymmetric or full 3D).
type(int_bounds_info) xi_idx
Indexes of first and last reference map eqns.
logical cont_damage
continuum damage modeling
integer p_glb
Global number of cells in each direction.
logical igr
Use information geometric regularization.
integer, dimension(3, 2) shear_bc_flip_indices
Shear stress BC reflection indices (1:3, 1:shear_BC_flip_num).
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
Deallocate all global grid, index, and equation-of-state parameter arrays.
integer thermal
1 = adiabatic, 2 = isotherm, 3 = transfer
real(wp), dimension(:), allocatable im_trans_c
integer perturb_flow_fluid
Fluid to be perturbed with perturb_flow flag.
type(int_bounds_info), dimension(1:3) idwint
integer recon_type
Reconstruction Type.
logical, parameter chemistry
Chemistry modeling.
integer mpi_info_int
MPI info for parallel IO with Lustre file systems.
type(ib_patch_parameters), dimension(num_patches_max) patch_ib
Immersed boundary patch parameters.
type(int_bounds_info) mom_idx
Indexes of first & last momentum eqns.
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.
real(wp), dimension(:), allocatable pb0
integer num_fluids
Number of different fluids present in the flow.
logical pre_stress
activate pre_stressed domain
real(wp), dimension(:), allocatable im_trans_t
integer dist_type
1 = binormal, 2 = lognormal-normal
integer weno_polyn
Degree of the WENO polynomials (polyn).
real(wp), dimension(:), allocatable y_cc
logical non_axis_sym
Use existing IC data.
type(int_bounds_info) stress_idx
Indexes of elastic shear stress eqns.
integer proc_rank
Rank of the local processor Number of cells in the x-, y- and z-coordinate directions.
real(wp), dimension(:), allocatable re_trans_t
integer n_idx
Index of number density.
logical dummy
AMDFlang workaround for case-optimization + GPU-kernel bug.
integer, dimension(:), allocatable proc_coords
Processor coordinates in MPI_CART_COMM.
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.
real(wp), dimension(:), allocatable weight
integer gamma_idx
Index of specific heat ratio func. eqn.
real(wp), dimension(:), allocatable k_v
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.
type(int_bounds_info) cont_idx
Indexes of first & last continuity eqns.
integer num_patches
Number of patches composing initial condition.
logical ib
Turn immersed boundaries on.
type(int_bounds_info) b_idx
Indexes of first and last magnetic field eqns.
integer num_bc_patches
Number of boundary condition patches.
integer tensor_size
Number of components in the nonsymmetric tensor.
type(bc_patch_parameters), dimension(num_bc_patches_max) patch_bc
Boundary condition patch parameters.
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
Precision of output files.
logical hyperelasticity
activate hyperelasticity
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(vec3_dt), dimension(:), allocatable airfoil_grid_u
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
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
Locations of cell-centers (cc) in x-, y- and z-directions, respectively.
real(wp), dimension(:), allocatable k_g
integer mixlayer_perturb_nk
Number of Fourier modes for perturbation with mixlayer_perturb flag.
type(chemistry_parameters) chem_params
integer num_vels
Number of velocity components (different from num_dims for mhd).
integer perturb_sph_fluid
Fluid to be perturbed with perturb_sph flag.
real(wp), dimension(:), allocatable x_cb
Locations of cell-boundaries (cb) in x-, y- and z-directions, respectively.
logical relax
activate phase change
logical qbmm
Quadrature moment method.
logical old_grid
Use existing grid data.
type(bub_bounds_info) bub_idx
Indexes of first & last bubble variable eqns.
integer damage_idx
Index of damage state variable (D) for continuum damage model.
real(wp) pi_fac
Factor for artificial pi_inf.
logical hyper_cleaning
Hyperbolic cleaning for MHD.
real(wp), dimension(num_fluids_max) fluid_rho
type(int_bounds_info) internalenergies_idx
Indexes of first & last internal energy eqns.
real(wp), dimension(:), allocatable z_cc
real(wp) pref
Reference parameters for Tait EOS.
real(wp) bx0
Constant magnetic field in the x-direction (1D).
integer b_size
Number of components in the b tensor.
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(vec3_dt), dimension(:), allocatable airfoil_grid_l
type(ic_patch_parameters), dimension(num_patches_max) patch_icpp
IC patch parameters (max: num_patches_max).
type(int_bounds_info) adv_idx
Indexes of first & last advection eqns.
integer weno_order
Order of accuracy for the WENO reconstruction.
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.
integer nmom
Number of carried moments.
integer e_idx
Index of total energy equation.
type(cell_num_bounds) cells_bounds
logical down_sample
Down-sample the output data.
real(wp), dimension(:), allocatable mass_g0
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.
logical elasticity
elasticity modeling, true for hyper or hypo
type(mpi_io_var), public mpi_io_data
integer c_idx
Index of the color function.
real(wp) mixlayer_vel_coef
Coefficient for the hyperbolic tangent streamwise velocity profile.
impure subroutine s_initialize_parallel_io
Configure MPI parallel I/O settings and allocate processor coordinate arrays.
logical mpp_lim
Alpha limiter.
real(wp), dimension(:), allocatable omegan
integer, dimension(:,:,:), allocatable logic_grid
integer igr_order
IGR reconstruction order.
integer psi_idx
Index of hyperbolic cleaning state variable for MHD.
integer muscl_polyn
Degree of the MUSCL polynomials (polyn).
real(wp), dimension(:), allocatable re_trans_c
type(subgrid_bubble_physical_parameters) bub_pp
type(int_bounds_info) species_idx
Indexes of first & last concentration eqns.
real(wp) rhorv
standard deviations in R/V
logical relativity
Relativity for RMHD.
real(wp), dimension(:), allocatable pe_t
integer alf_idx
Index of void fraction.
real(wp), dimension(:), allocatable mass_v0
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.
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.
subroutine, public s_configure_coordinate_bounds(recon_type, weno_polyn, muscl_polyn, igr_order, buff_size, idwint, idwbuff, viscous, bubbles_lagrange, m, n, p, num_dims, igr, ib)
Compute ghost-cell buffer size and set interior/buffered coordinate index bounds.
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-.
Derived type adding initial condition (ic) patch parameters as attributes NOTE: The requirements for ...
Integer bounds for variables.
Derived type annexing the physical parameters (PP) of the fluids. These include the specific heat rat...
Derived type for bubble variables pb and mv at quadrature nodes (qbmm).
Derived type annexing the physical parameters required for sub-grid bubble models.
Generic 3-component vector (e.g., spatial coordinates or field components) Named _dt (derived types: ...