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# 8 "/home/runner/work/MFC/MFC/src/common/include/case.fpp"
13
14! For moving immersed boundaries in simulation
15# 12 "/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 ! Shared state: generated_decls, sys_size, eqn_idx, b_size, tensor_size, chemistry, elasticity, shear_*
29
30 implicit none
31
32 ! Logistics
33 integer :: num_procs !< Number of processors
34 logical :: non_axis_sym !< Use existing IC data
35 logical :: cfl_dt
36
37 ! Computational Domain Parameters
38
39 integer :: proc_rank !< Rank of the local processor Number of cells in the x-, y- and z-coordinate directions
40
41 !> @name Max and min number of cells in a direction of each combination of x-,y-, and z-
43 integer(kind=8) :: nglobal !< Global number of cells in the domain
44 integer :: m_glb, n_glb, p_glb !< Global number of cells in each direction
45 integer :: grid_geometry !< Cylindrical coordinates (either axisymmetric or full 3D)
46 !> Locations of cell-centers (cc) in x-, y- and z-directions, respectively
47 real(wp), allocatable, dimension(:) :: x_cc, y_cc, z_cc
48 !> Locations of cell-boundaries (cb) in x-, y- and z-directions, respectively
49 real(wp), allocatable, dimension(:) :: x_cb, y_cb, z_cb
50 real(wp) :: dx, dy, dz !< Minimum cell-widths in the x-, y- and z-coordinate directions
51 type(bounds_info) :: x_domain, y_domain, z_domain !< Locations of the domain bounds in the x-, y- and z-coordinate directions
52
53 ! Simulation Algorithm Parameters
54 ! sys_size, eqn_idx, b_size, tensor_size, chemistry, elasticity, shear_*: in m_global_parameters_common
55 ! weno_polyn, muscl_polyn, num_dims, num_vels: in m_global_parameters_common
56 ! Annotations of the structure, i.e. the organization, of the state vectors
57 type(qbmm_idx_info) :: qbmm_idx !< QBMM moment index mappings.
58 ! Cell Indices for the (local) interior points (O-m, O-n, 0-p). Stands for "InDices With BUFFer".
59 type(int_bounds_info) :: idwint(1:3)
60
61 ! Cell indices (InDices With BUFFer): includes buffer except in pre_process
63 type(int_bounds_info) :: bc_x, bc_y, bc_z !< Boundary conditions in the x-, y- and z-coordinate directions
64 ! simplex_params: auto-generated in generated_decls.fpp
65
66 ! fluid_rho (perturbs surrounding-air density to break grid symmetry): auto-generated in generated_decls.fpp
67 ! proc_coords, start_idx, mpiiofs, mpi_info_int: in m_global_parameters_common
68#ifdef MFC_MPI
69 type(mpi_io_var), public :: mpi_io_data
70#endif
71
72 ! Initial Condition Parameters patch_icpp, patch_bc: auto-generated in generated_decls.fpp
73 logical :: bc_io !< whether or not to save BC data
74
75 ! Fluids Physical Parameters fluid_pp, bub_pp: auto-generated in generated_decls.fpp
77 !> @name Bubble modeling
78 !> @{
79 real(wp) :: eu
80 real(wp), dimension(:), allocatable :: weight, r0
81 integer :: nmom !< Number of carried moments
82 !> @}
83
84 ! patch_ib, ib_airfoil, stl_models: auto-generated in generated_decls.fpp
85
86 !> @name Non-polytropic bubble gas compression
87 !> @{
88 real(wp) :: phi_vg, phi_gv, pe_c, tw, k_vl, k_gl
89 real(wp) :: gam_m
90 real(wp), dimension(:), allocatable :: pb0, mass_g0, mass_v0, pe_t, k_v, k_g
91 real(wp), dimension(:), allocatable :: re_trans_t, re_trans_c, im_trans_t, im_trans_c, omegan
93 !> @}
94
95 integer, allocatable, dimension(:,:,:) :: logic_grid
96 type(pres_field) :: pb
97 type(pres_field) :: mv
98 integer :: buff_size !< Number of ghost cells for boundary condition storage
99
100contains
101
102 !> Assigns default values to user inputs prior to reading them in. This allows for an easier consistency check of these
103 !! parameters once they are read from the input file.
105
106 integer :: i !< Generic loop operator
107
108 ! Shared defaults (case_dir, m/n/p, cyl_coord, cfl flags, model_eqns, elasticity, BC blocks,
109 ! recon/weno/muscl/num_fluids/igr/mhd/relativity under case-opt guard, Tait EOS, bubble flags,
110 ! IB flags, parallel I/O flags, fft_wrt)
111
113
114 ! Boundary conditions (bc_x/y/z are per-target declarations, not visible in common)
115 bc_x%beg = dflt_int; bc_x%end = dflt_int
116 bc_y%beg = dflt_int; bc_y%end = dflt_int
117 bc_z%beg = dflt_int; bc_z%end = dflt_int
118
119# 109 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
120# 110 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
121 bc_x%vb1 = 0._wp
122 bc_x%ve1 = 0._wp
123# 110 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
124 bc_x%vb2 = 0._wp
125 bc_x%ve2 = 0._wp
126# 110 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
127 bc_x%vb3 = 0._wp
128 bc_x%ve3 = 0._wp
129# 113 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
130# 109 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
131# 110 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
132 bc_y%vb1 = 0._wp
133 bc_y%ve1 = 0._wp
134# 110 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
135 bc_y%vb2 = 0._wp
136 bc_y%ve2 = 0._wp
137# 110 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
138 bc_y%vb3 = 0._wp
139 bc_y%ve3 = 0._wp
140# 113 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
141# 109 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
142# 110 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
143 bc_z%vb1 = 0._wp
144 bc_z%ve1 = 0._wp
145# 110 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
146 bc_z%vb2 = 0._wp
147 bc_z%ve2 = 0._wp
148# 110 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
149 bc_z%vb3 = 0._wp
150 bc_z%ve3 = 0._wp
151# 113 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
152# 114 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
153
154# 116 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
155 bc_x%isothermal_in = .false.
156 bc_x%isothermal_out = .false.
157 bc_x%Twall_in = dflt_real
158 bc_x%Twall_out = dflt_real
159# 116 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
160 bc_y%isothermal_in = .false.
161 bc_y%isothermal_out = .false.
162 bc_y%Twall_in = dflt_real
163 bc_y%Twall_out = dflt_real
164# 116 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
165 bc_z%isothermal_in = .false.
166 bc_z%isothermal_out = .false.
167 bc_z%Twall_in = dflt_real
168 bc_z%Twall_out = dflt_real
169# 121 "/home/runner/work/MFC/MFC/src/pre_process/m_global_parameters.fpp"
170
172
173 ! Logistics (pre-specific)
174 old_grid = .false.
175 old_ic = .false.
176 t_step_old = dflt_int
177 cfl_dt = .false.
178
179 ! Computational domain parameters (pre-specific)
180 x_domain%beg = dflt_real
181 x_domain%end = dflt_real
182 y_domain%beg = dflt_real
183 y_domain%end = dflt_real
184 z_domain%beg = dflt_real
185 z_domain%end = dflt_real
186
187 stretch_x = .false.
188 stretch_y = .false.
189 stretch_z = .false.
190
191 a_x = dflt_real
192 a_y = dflt_real
193 a_z = dflt_real
194 loops_x = 1
195 loops_y = 1
196 loops_z = 1
197 x_a = dflt_real
198 x_b = dflt_real
199 y_a = dflt_real
200 y_b = dflt_real
201 z_a = dflt_real
202 z_b = dflt_real
203
204 ! Simulation algorithm parameters (pre-specific)
205 palpha_eps = dflt_real
206 ptgalpha_eps = dflt_real
207 igr_order = dflt_int
208 pre_stress = .false.
209
210 precision = 2
211 viscous = .false.
212 mixlayer_vel_profile = .false.
213 mixlayer_vel_coef = 1._wp
214 mixlayer_perturb = .false.
216 mixlayer_perturb_k0 = 0.4446_wp
217 perturb_flow = .false.
218 perturb_flow_fluid = dflt_int
219 perturb_flow_mag = dflt_real
220 perturb_sph = .false.
221 perturb_sph_fluid = dflt_int
222 fluid_rho = dflt_real
223 elliptic_smoothing_iters = dflt_int
224 elliptic_smoothing = .false.
225
226 simplex_perturb = .false.
227 simplex_params%perturb_vel(:) = .false.
228 simplex_params%perturb_vel_freq(:) = dflt_real
229 simplex_params%perturb_vel_scale(:) = dflt_real
230 simplex_params%perturb_vel_offset(:,:) = dflt_real
231 simplex_params%perturb_dens(:) = .false.
232 simplex_params%perturb_dens_freq(:) = dflt_real
233 simplex_params%perturb_dens_scale(:) = dflt_real
234 simplex_params%perturb_dens_offset(:,:) = dflt_real
235
236 ! Initial condition parameters
237 num_patches = dflt_int
238
239 do i = 1, num_patches_max
240 patch_icpp(i)%geometry = dflt_int
241 patch_icpp(i)%model_id = 0
242 patch_icpp(i)%x_centroid = dflt_real
243 patch_icpp(i)%y_centroid = dflt_real
244 patch_icpp(i)%z_centroid = dflt_real
245 patch_icpp(i)%length_x = dflt_real
246 patch_icpp(i)%length_y = dflt_real
247 patch_icpp(i)%length_z = dflt_real
248 patch_icpp(i)%radius = dflt_real
249 patch_icpp(i)%epsilon = dflt_real
250 patch_icpp(i)%beta = dflt_real
251 patch_icpp(i)%normal = dflt_real
252 patch_icpp(i)%radii = dflt_real
253 patch_icpp(i)%alter_patch = .false.
254 patch_icpp(i)%alter_patch(0) = .true.
255 patch_icpp(i)%smoothen = .false.
256 patch_icpp(i)%smooth_patch_id = i
257 patch_icpp(i)%smooth_coeff = dflt_real
258 patch_icpp(i)%alpha_rho = dflt_real
259 patch_icpp(i)%rho = dflt_real
260 patch_icpp(i)%vel = dflt_real
261 patch_icpp(i)%pres = dflt_real
262 patch_icpp(i)%alpha = dflt_real
263 patch_icpp(i)%gamma = dflt_real
264 patch_icpp(i)%pi_inf = dflt_real
265 patch_icpp(i)%cv = 0._wp
266 patch_icpp(i)%qv = 0._wp
267 patch_icpp(i)%qvp = 0._wp
268 patch_icpp(i)%tau_e = 0._wp
269 patch_icpp(i)%Bx = dflt_real
270 patch_icpp(i)%By = dflt_real
271 patch_icpp(i)%Bz = dflt_real
272 patch_icpp(i)%a(2) = dflt_real
273 patch_icpp(i)%a(3) = dflt_real
274 patch_icpp(i)%a(4) = dflt_real
275 patch_icpp(i)%a(5) = dflt_real
276 patch_icpp(i)%a(6) = dflt_real
277 patch_icpp(i)%a(7) = dflt_real
278 patch_icpp(i)%a(8) = dflt_real
279 patch_icpp(i)%a(9) = dflt_real
280 patch_icpp(i)%non_axis_sym = .false.
281 patch_icpp(i)%fourier_cos(:) = 0._wp
282 patch_icpp(i)%fourier_sin(:) = 0._wp
283 patch_icpp(i)%modal_clip_r_to_min = .false.
284 patch_icpp(i)%modal_r_min = 1.e-12_wp
285 patch_icpp(i)%modal_use_exp_form = .false.
286 patch_icpp(i)%sph_har_coeff(:,:) = 0._wp
287
288 ! should get all of r0's and v0's
289 patch_icpp(i)%r0 = dflt_real
290 patch_icpp(i)%v0 = dflt_real
291
292 patch_icpp(i)%p0 = dflt_real
293 patch_icpp(i)%m0 = dflt_real
294
295 patch_icpp(i)%hcid = dflt_int
296
297 if (chemistry) then
298 patch_icpp(i)%Y(:) = 0._wp
299 end if
300 end do
301
303 bc_io = .false.
304
305 do i = 1, num_bc_patches_max
306 patch_bc(i)%geometry = dflt_int
307 patch_bc(i)%type = dflt_int
308 patch_bc(i)%dir = dflt_int
309 patch_bc(i)%loc = dflt_int
310 patch_bc(i)%centroid(:) = dflt_real
311 patch_bc(i)%length(:) = dflt_real
312 patch_bc(i)%radius = dflt_real
313 end do
314
315 ! Bubble modeling (pre-specific)
316 polytropic = .true.
317 thermal = dflt_int
318 nb = dflt_int
319
320 eu = dflt_real
321 ca = dflt_real
322 re_inv = dflt_real
323 web = dflt_real
324
325 nmom = 1
326 sigr = dflt_real
327 sigv = dflt_real
328 rhorv = 0._wp
329 dist_type = dflt_int
330
331 r_g = dflt_real
332 r_v = dflt_real
333 phi_vg = dflt_real
334 phi_gv = dflt_real
335 pe_c = dflt_real
336 tw = dflt_real
337
338 pi_fac = 1._wp
339
340 do i = 1, num_ib_patches_max_namelist
341 patch_ib(i)%geometry = dflt_int
342 patch_ib(i)%x_centroid = dflt_real
343 patch_ib(i)%y_centroid = dflt_real
344 patch_ib(i)%z_centroid = dflt_real
345 patch_ib(i)%length_x = dflt_real
346 patch_ib(i)%length_y = dflt_real
347 patch_ib(i)%length_z = dflt_real
348 patch_ib(i)%radius = dflt_real
349 patch_ib(i)%airfoil_id = 0
350 patch_ib(i)%model_id = 0
351 patch_ib(i)%slip = .false.
352
353 ! Variables to handle moving immersed boundaries, defaulting to no movement
354 patch_ib(i)%moving_ibm = 0
355 patch_ib(i)%vel(:) = 0._wp
356 patch_ib(i)%angles(:) = 0._wp
357 patch_ib(i)%angular_vel(:) = 0._wp
358 patch_ib(i)%mass = dflt_real
359 patch_ib(i)%moment = dflt_real
360 patch_ib(i)%centroid_offset(:) = 0._wp
361
362 ! sets values of a rotation matrix which can be used when calculating rotations
363 patch_ib(i)%rotation_matrix = 0._wp
364 patch_ib(i)%rotation_matrix(1, 1) = 1._wp
365 patch_ib(i)%rotation_matrix(2, 2) = 1._wp
366 patch_ib(i)%rotation_matrix(3, 3) = 1._wp
367 patch_ib(i)%rotation_matrix_inverse = patch_ib(i)%rotation_matrix
368 end do
369
370 do i = 1, num_ib_airfoils_max
371 ib_airfoil(i)%c = dflt_real
372 ib_airfoil(i)%p = dflt_real
373 ib_airfoil(i)%t = dflt_real
374 ib_airfoil(i)%m = dflt_real
375 end do
376
378
379 do i = 1, num_stl_models_max
380 stl_models(i)%model_filepath(:) = dflt_char
381 stl_models(i)%model_translate(:) = 0._wp
382 stl_models(i)%model_scale(:) = 1._wp
383 stl_models(i)%model_threshold = ray_tracing_threshold
384 end do
385
386 chem_params%gamma_method = 1
387 chem_params%transport_model = 1
388
389 ! Fluids physical parameters
390 do i = 1, num_fluids_max
391 fluid_pp(i)%gamma = dflt_real
392 fluid_pp(i)%pi_inf = dflt_real
393 fluid_pp(i)%cv = 0._wp
394 fluid_pp(i)%qv = 0._wp
395 fluid_pp(i)%qvp = 0._wp
396 fluid_pp(i)%G = 0._wp
397 fluid_pp(i)%non_newtonian = .false.
398 fluid_pp(i)%K = dflt_real
399 fluid_pp(i)%nn = dflt_real
400 fluid_pp(i)%tau0 = 0._wp
401 fluid_pp(i)%hb_m = dflt_real
402 fluid_pp(i)%mu_min = dflt_real
403 fluid_pp(i)%mu_max = dflt_real
404 fluid_pp(i)%mu_bulk = 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
430
431 !> Computation of parameters, allocation procedures, and/or any other tasks needed to properly setup the module
433
434 integer :: i, j, fac
435
436 if (recon_type == recon_type_weno) then
437 weno_polyn = (weno_order - 1)/2
438 else if (recon_type == recon_type_muscl) then
440 end if
441
442 ! Gamma/Pi_inf: force num_fluids=1 (pre_process-specific side effect of the gamma-law model)
443 if (model_eqns == model_eqns_gamma_law) num_fluids = 1
444
445 ! Pre-process sets nmom to 6 for qbmm before the shared eqn_idx setup
446 ! (guards match the original site: 5-equation bubbles with 4-node qbmm)
447 if (model_eqns == model_eqns_5eq .and. bubbles_euler .and. qbmm .and. nnode == 4) nmom = 6
448
449 ! Populate eqn_idx, sys_size, b_size, tensor_size, elasticity, shear_* (shared logic)
451
452 ! Per-target (pre_process): qbmm_idx allocations and fills
453 if (model_eqns == model_eqns_5eq .and. bubbles_euler) then
454 allocate (qbmm_idx%rs(nb), qbmm_idx%vs(nb))
455 allocate (qbmm_idx%ps(nb), qbmm_idx%ms(nb))
456
457 if (qbmm) then
458 allocate (qbmm_idx%moms(nb, nmom))
459 allocate (qbmm_idx%fullmom(nb,0:nmom,0:nmom))
460
461 do i = 1, nb
462 do j = 1, nmom
463 qbmm_idx%moms(i, j) = eqn_idx%bub%beg + (j - 1) + (i - 1)*nmom
464 end do
465 qbmm_idx%fullmom(i, 0, 0) = qbmm_idx%moms(i, 1)
466 qbmm_idx%fullmom(i, 1, 0) = qbmm_idx%moms(i, 2)
467 qbmm_idx%fullmom(i, 0, 1) = qbmm_idx%moms(i, 3)
468 qbmm_idx%fullmom(i, 2, 0) = qbmm_idx%moms(i, 4)
469 qbmm_idx%fullmom(i, 1, 1) = qbmm_idx%moms(i, 5)
470 qbmm_idx%fullmom(i, 0, 2) = qbmm_idx%moms(i, 6)
471 qbmm_idx%rs(i) = qbmm_idx%fullmom(i, 1, 0)
472 end do
473 else
474 do i = 1, nb
475 if (.not. polytropic) then
476 fac = 4
477 else
478 fac = 2
479 end if
480
481 qbmm_idx%rs(i) = eqn_idx%bub%beg + (i - 1)*fac
482 qbmm_idx%vs(i) = qbmm_idx%rs(i) + 1
483
484 if (.not. polytropic) then
485 qbmm_idx%ps(i) = qbmm_idx%vs(i) + 1
486 qbmm_idx%ms(i) = qbmm_idx%ps(i) + 1
487 end if
488 end do
489 end if
490 end if
491
492 if (model_eqns == model_eqns_4eq .and. bubbles_euler) then
493 allocate (qbmm_idx%rs(nb), qbmm_idx%vs(nb))
494 allocate (qbmm_idx%ps(nb), qbmm_idx%ms(nb))
495 allocate (weight(nb), r0(nb))
496
497 do i = 1, nb
498 if (.not. polytropic) then
499 fac = 4
500 else
501 fac = 2
502 end if
503
504 qbmm_idx%rs(i) = eqn_idx%bub%beg + (i - 1)*fac
505 qbmm_idx%vs(i) = qbmm_idx%rs(i) + 1
506
507 if (.not. polytropic) then
508 qbmm_idx%ps(i) = qbmm_idx%vs(i) + 1
509 qbmm_idx%ms(i) = qbmm_idx%ps(i) + 1
510 end if
511 end do
512
513 if (nb == 1) then
514 weight(:) = 1._wp
515 r0(:) = 1._wp
516 else if (nb < 1) then
517 stop 'Invalid value of nb'
518 end if
519
520 if (polytropic) then
521 rhoref = 1._wp
522 pref = 1._wp
523 end if
524 end if
525
528
529#ifdef MFC_MPI
530 if (qbmm .and. .not. polytropic) then
531 allocate (mpi_io_data%view(1:sys_size + 2*nb*nnode))
532 allocate (mpi_io_data%var(1:sys_size + 2*nb*nnode))
533 else
534 allocate (mpi_io_data%view(1:sys_size))
535 allocate (mpi_io_data%var(1:sys_size))
536 end if
537
538 if (.not. down_sample) then
539 do i = 1, sys_size
540 allocate (mpi_io_data%var(i)%sf(0:m,0:n,0:p))
541 mpi_io_data%var(i)%sf => null()
542 end do
543 end if
544 if (qbmm .and. .not. polytropic) then
545 do i = sys_size + 1, sys_size + 2*nb*nnode
546 allocate (mpi_io_data%var(i)%sf(0:m,0:n,0:p))
547 mpi_io_data%var(i)%sf => null()
548 end do
549 end if
550#endif
551
552 ! Allocating grid variables for the x-direction
553 allocate (x_cc(0:m), x_cb(-1:m))
554 ! Allocating grid variables for the y- and z-directions
555 if (n > 0) then
556 allocate (y_cc(0:n), y_cb(-1:n))
557 if (p > 0) then
558 allocate (z_cc(0:p), z_cb(-1:p))
559 end if
560 end if
561
562 if (cyl_coord .neqv. .true.) then ! Cartesian grid
563 grid_geometry = 1
564 else if (cyl_coord .and. p == 0) then ! Axisymmetric cylindrical grid
565 grid_geometry = 2
566 else ! Fully 3D cylindrical grid
567 grid_geometry = 3
568 end if
569
570 if (.not. igr) then
571 allocate (logic_grid(0:m,0:n,0:p))
572 end if
573
575
576 !> Configure MPI parallel I/O settings and allocate processor coordinate arrays.
581 end subroutine s_initialize_parallel_io
582
583 !> Deallocate all global grid, index, and equation-of-state parameter arrays.
585
586 integer :: i
587
588 if (bubbles_euler) then
589 deallocate (qbmm_idx%rs, qbmm_idx%vs, qbmm_idx%ps, qbmm_idx%ms)
590 if (qbmm) deallocate (qbmm_idx%moms, qbmm_idx%fullmom)
591 end if
592
593 ! Deallocating grid variables for the x-direction
594 deallocate (x_cc, x_cb)
595 ! Deallocating grid variables for the y- and z-directions
596 if (n > 0) then
597 deallocate (y_cc, y_cb)
598 if (p > 0) then
599 deallocate (z_cc, z_cb)
600 end if
601 end if
602
603 ! Shared: deallocate proc_coords and start_idx
605
606#ifdef MFC_MPI
607 if (parallel_io) then
608 do i = 1, sys_size
609 mpi_io_data%var(i)%sf => null()
610 end do
611
612 deallocate (mpi_io_data%var)
613 deallocate (mpi_io_data%view)
614 end if
615#endif
616
618
619end module m_global_parameters
integer, intent(in) j
Shared derived types for field data, patch geometry, bubble dynamics, and MPI I/O structures.
Shared global parameters and equation-index setup for all three executables. Each per-target m_global...
integer muscl_polyn
Degree of the MUSCL polynomials.
type(ib_stl_parameters), dimension(num_stl_models_max) stl_models
Per-STL model parameters.
impure subroutine s_finalize_global_parameters_common
Shared finalize core: deallocate proc_coords and start_idx. Per-target finalize routines call this fi...
type(ic_patch_parameters), dimension(num_patches_max) patch_icpp
IC patch parameters.
type(subgrid_bubble_physical_parameters) bub_pp
Subgrid bubble physical parameters.
integer sys_size
Number of unknowns in system of equations.
type(physical_parameters), dimension(num_fluids_max) fluid_pp
Per-fluid stiffened-gas EOS parameters, Reynolds numbers, and shear modulus.
type(bc_patch_parameters), dimension(num_bc_patches_max) patch_bc
Boundary condition patch parameters.
integer num_dims
Number of spatial dimensions.
type(eqn_idx_info) eqn_idx
All conserved-variable equation index ranges and scalars.
real(wp), dimension(num_fluids_max) fluid_rho
impure subroutine s_initialize_eqn_idx(nmom_in, nb_in)
Initialize equation-index state (eqn_idx, sys_size, b_size, tensor_size) from the namelist parameters...
type(ib_patch_parameters), dimension(num_ib_patches_max_namelist) patch_ib
Immersed boundary patch parameters.
impure subroutine s_initialize_parallel_io_common
Configure MPI parallel I/O settings and allocate processor coordinate arrays. Shared across all three...
integer weno_polyn
Degree of the WENO polynomials.
type(ib_airfoil_parameters), dimension(num_ib_airfoils_max) ib_airfoil
Per-airfoil NACA user inputs.
impure subroutine s_assign_common_defaults
Assign default values to the user-input parameters that are shared across all three executables (pre_...
Defines global parameters for the computational domain, simulation algorithm, and initial conditions.
integer grid_geometry
Cylindrical coordinates (either axisymmetric or full 3D).
integer p_glb
Global number of cells in each direction.
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.
real(wp), dimension(:), allocatable im_trans_c
type(int_bounds_info), dimension(1:3) idwint
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
real(wp), dimension(:), allocatable im_trans_t
type(qbmm_idx_info) qbmm_idx
QBMM moment index mappings.
real(wp), dimension(:), allocatable y_cc
logical non_axis_sym
Use existing IC data.
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
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.
real(wp), dimension(:), allocatable weight
real(wp), dimension(:), allocatable k_v
type(int_bounds_info), dimension(1:3) idwbuff
integer buff_size
Number of ghost cells for boundary condition storage.
real(wp), dimension(:), allocatable z_cb
type(int_bounds_info) bc_y
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...
real(wp), dimension(:), allocatable x_cc
Locations of cell-centers (cc) in x-, y- and z-directions, respectively.
real(wp), dimension(:), allocatable k_g
type(chemistry_parameters) chem_params
type(int_bounds_info) bc_x
real(wp), dimension(:), allocatable x_cb
Locations of cell-boundaries (cb) in x-, y- and z-directions, respectively.
real(wp), dimension(:), allocatable z_cc
integer num_procs
Number of processors.
integer nmom
Number of carried moments.
type(cell_num_bounds) cells_bounds
real(wp), dimension(:), allocatable mass_g0
type(mpi_io_var), public mpi_io_data
impure subroutine s_initialize_parallel_io
Configure MPI parallel I/O settings and allocate processor coordinate arrays.
real(wp), dimension(:), allocatable omegan
integer, dimension(:,:,:), allocatable logic_grid
real(wp), dimension(:), allocatable re_trans_c
real(wp), dimension(:), allocatable pe_t
real(wp), dimension(:), allocatable mass_v0
integer(kind=8) nglobal
Global number of cells in the domain.
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.
Max and min number of cells in a direction of each combination of x-,y-, and z-.
Integer bounds for variables.
Derived type for bubble variables pb and mv at quadrature nodes (qbmm).
QBMM moment index mappings - separate from bub beg/end so eqn_idx contains no allocatables.