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