MFC
Exascale flow solver
Loading...
Searching...
No Matches
m_initial_condition.fpp.f90
Go to the documentation of this file.
1# 1 "/home/runner/work/MFC/MFC/src/pre_process/m_initial_condition.fpp"
2!>
3!! @file
4!! @brief Contains module m_initial_condition
5
6!> @brief Assembles initial conditions by layering prioritized patches via constructive solid geometry
8
11 use m_mpi_proxy
12 use m_helper
17 use m_chemistry
19
20 implicit none
21
22 type(ic_context) :: ic !< Initial-condition state (fields, bc types, patch ids)
23
24contains
25
26 !> Computation of parameters, allocation procedures, and/or any other tasks needed to properly setup the module
28
29 integer :: i, j, k, l
30
31 allocate (ic%q_prim_vf(1:sys_size))
32 allocate (ic%q_cons_vf(1:sys_size))
33
34 do i = 1, sys_size
35 allocate (ic%q_prim_vf(i)%sf(idwbuff(1)%beg:idwbuff(1)%end,idwbuff(2)%beg:idwbuff(2)%end,idwbuff(3)%beg:idwbuff(3)%end))
36 allocate (ic%q_cons_vf(i)%sf(idwbuff(1)%beg:idwbuff(1)%end,idwbuff(2)%beg:idwbuff(2)%end,idwbuff(3)%beg:idwbuff(3)%end))
37 end do
38
39 if (chemistry) then
40 allocate (ic%q_T_sf%sf(idwbuff(1)%beg:idwbuff(1)%end,idwbuff(2)%beg:idwbuff(2)%end,idwbuff(3)%beg:idwbuff(3)%end))
41 end if
42
43 allocate (ic%patch_id_fp(0:m,0:n,0:p))
44
45 if (qbmm .and. .not. polytropic) then
46 allocate (pb%sf(0:m,0:n,0:p,1:nnode,1:nb))
47 allocate (mv%sf(0:m,0:n,0:p,1:nnode,1:nb))
48 end if
49
50 do i = 1, sys_size
51 ic%q_cons_vf(i)%sf = -1.e-6_stp ! real(dflt_real, kind=stp) ! TODO :: remove this magic number
52 ic%q_prim_vf(i)%sf = -1.e-6_stp ! real(dflt_real, kind=stp)
53 end do
54
55 allocate (ic%bc_type(1:num_dims,1:2))
56
57 allocate (ic%bc_type(1, 1)%sf(0:0,0:n,0:p))
58 allocate (ic%bc_type(1, 2)%sf(0:0,0:n,0:p))
59
60 do l = 0, p
61 do k = 0, n
62 ic%bc_type(1, 1)%sf(0, k, l) = int(min(bc_x%beg, 0), kind=1)
63 ic%bc_type(1, 2)%sf(0, k, l) = int(min(bc_x%end, 0), kind=1)
64 end do
65 end do
66
67 if (n > 0) then
68 allocate (ic%bc_type(2, 1)%sf(-buff_size:m + buff_size,0:0,0:p))
69 allocate (ic%bc_type(2, 2)%sf(-buff_size:m + buff_size,0:0,0:p))
70
71 do l = 0, p
72 do j = -buff_size, m + buff_size
73 ic%bc_type(2, 1)%sf(j, 0, l) = int(min(bc_y%beg, 0), kind=1)
74 ic%bc_type(2, 2)%sf(j, 0, l) = int(min(bc_y%end, 0), kind=1)
75 end do
76 end do
77
78 if (p > 0) then
79 allocate (ic%bc_type(3, 1)%sf(-buff_size:m + buff_size,-buff_size:n + buff_size,0:0))
80 allocate (ic%bc_type(3, 2)%sf(-buff_size:m + buff_size,-buff_size:n + buff_size,0:0))
81
82 do k = -buff_size, n + buff_size
83 do j = -buff_size, m + buff_size
84 ic%bc_type(3, 1)%sf(j, k, 0) = int(min(bc_z%beg, 0), kind=1)
85 ic%bc_type(3, 2)%sf(j, k, 0) = int(min(bc_z%end, 0), kind=1)
86 end do
87 end do
88 end if
89 end if
90
91 ! Initial damage state is always zero
92 if (cont_damage) then
93 ic%q_cons_vf(eqn_idx%damage)%sf = 0._wp
94 ic%q_prim_vf(eqn_idx%damage)%sf = 0._wp
95 end if
96
97 ! Initial hyper_cleaning state is always zero TODO more general
98 if (hyper_cleaning) then
99 ic%q_cons_vf(eqn_idx%psi)%sf = 0._wp
100 ic%q_prim_vf(eqn_idx%psi)%sf = 0._wp
101 end if
102
103 ! Setting default values for patch identities bookkeeping variable. This is necessary to avoid any confusion in the
104 ! assessment of the extent of application that the overwrite permissions give a patch when it is being applied in the
105 ! domain.
106 ic%patch_id_fp = 0
107
109
110 !> Iterate over patches and, depending on the geometry type, call the related subroutine to setup the said geometry on the grid
111 !! using the primitive variables included with the patch parameters. The subroutine is complete once the primitive variables are
112 !! converted to conservative ones.
114
115 integer :: i
116
117 if (old_ic) then
118 call s_convert_conservative_to_primitive_variables(ic%q_cons_vf, ic%q_T_sf, ic%q_prim_vf, idwbuff)
119 end if
120
121 call s_apply_icpp_patches(ic%patch_id_fp, ic%q_prim_vf)
122
123 if (num_bc_patches > 0) call s_apply_boundary_patches(ic%q_prim_vf, ic%bc_type)
124
125 if (perturb_flow) call s_perturb_surrounding_flow(ic%q_prim_vf)
126 if (perturb_sph) call s_perturb_sphere(ic%q_prim_vf)
127 if (mixlayer_perturb) call s_perturb_mixlayer(ic%q_prim_vf)
128 if (simplex_perturb) call s_perturb_simplex(ic%q_prim_vf)
129 if (chemistry) call s_compute_t_from_primitives(ic%q_T_sf, ic%q_prim_vf, idwint)
130
131 if (elliptic_smoothing .and. chemistry) then
132 call s_elliptic_smoothing(ic%q_prim_vf, ic%bc_type, ic%q_T_sf)
133 call s_compute_t_from_primitives(ic%q_T_sf, ic%q_prim_vf, idwint)
134 else if (elliptic_smoothing) then
135 call s_elliptic_smoothing(ic%q_prim_vf, ic%bc_type)
136 end if
137
138 call s_convert_primitive_to_conservative_variables(ic%q_prim_vf, ic%q_cons_vf)
139
140 if (qbmm .and. .not. polytropic) then
141 call s_initialize_mv(ic%q_cons_vf, mv%sf)
142 call s_initialize_pb(ic%q_cons_vf, mv%sf, pb%sf)
143 end if
144
145 end subroutine s_generate_initial_condition
146
147 !> Deallocation procedures for the module
149
150 integer :: i
151
152 do i = 1, sys_size
153 deallocate (ic%q_prim_vf(i)%sf)
154 deallocate (ic%q_cons_vf(i)%sf)
155 end do
156
157 deallocate (ic%q_prim_vf)
158 deallocate (ic%q_cons_vf)
159
160 if (chemistry) then
161 deallocate (ic%q_T_sf%sf)
162 end if
163
164 deallocate (ic%patch_id_fp)
165
166 deallocate (ic%bc_type(1, 1)%sf)
167 deallocate (ic%bc_type(1, 2)%sf)
168
169 if (n > 0) then
170 deallocate (ic%bc_type(2, 1)%sf)
171 deallocate (ic%bc_type(2, 2)%sf)
172 end if
173
174 if (p > 0) then
175 deallocate (ic%bc_type(3, 1)%sf)
176 deallocate (ic%bc_type(3, 2)%sf)
177 end if
178
179 deallocate (ic%bc_type)
180
182
183end module m_initial_condition
integer, intent(in) k
integer, intent(in) j
integer, intent(in) l
Assigns initial primitive variables to computational cells based on patch geometry.
Applies spatially varying boundary condition patches along domain edges and faces.
impure subroutine, public s_apply_boundary_patches(q_prim_vf, bc_type)
Iterate over all boundary condition patches and dispatch them by geometry type.
Multi-species chemistry interface for thermodynamic properties, reaction rates, and transport coeffic...
subroutine s_compute_t_from_primitives(q_t_sf, q_prim_vf, bounds)
Compute the temperature field from primitive variables using the ideal gas law and mixture molecular ...
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.
type(int_bounds_info), dimension(1:3) idwint
type(int_bounds_info) bc_z
Boundary conditions in the x-, y- and z-coordinate directions.
type(int_bounds_info), dimension(1:3) idwbuff
integer buff_size
Number of ghost cells for boundary condition storage.
type(int_bounds_info) bc_y
type(int_bounds_info) bc_x
Utility routines for bubble model setup, coordinate transforms, array sampling, and special functions...
Allocate memory and read initial condition data for IC extrusion.
impure subroutine, public s_apply_icpp_patches(patch_id_fp, q_prim_vf)
Dispatch each initial condition patch to its geometry-specific initialization routine.
Assembles initial conditions by layering prioritized patches via constructive solid geometry.
type(ic_context) ic
Initial-condition state (fields, bc types, patch ids).
impure subroutine s_initialize_initial_condition_module
Computation of parameters, allocation procedures, and/or any other tasks needed to properly setup the...
impure subroutine s_finalize_initial_condition_module
Deallocation procedures for the module.
impure subroutine s_generate_initial_condition
Iterate over patches and, depending on the geometry type, call the related subroutine to setup the sa...
Broadcasts user inputs and decomposes the domain across MPI ranks for pre-processing.
Perturbs initial mean flow fields with random noise, mixing-layer instabilities, or simplex noise.
impure subroutine s_perturb_sphere(q_prim_vf)
Randomly perturb partial density fields at the interface of a spherical volume fraction region.
impure subroutine s_perturb_surrounding_flow(q_prim_vf)
Add random noise to the velocity and void fraction of the surrounding flow field.
impure subroutine s_elliptic_smoothing(q_prim_vf, bc_type, q_t_sf)
Iteratively smooth all primitive variable fields using a discrete elliptic (Laplacian) filter.
subroutine s_perturb_simplex(q_prim_vf)
Perturb velocity and volume fraction fields using multi-octave simplex noise.
subroutine s_perturb_mixlayer(q_prim_vf)
Compute velocity perturbations for a temporal mixing layer with a hyperbolic tangent mean streamwise ...
Conservative-to-primitive variable conversion, mixture property evaluation, and pressure computation.
subroutine, public s_initialize_mv(qk_cons_vf, mv)
Initialize bubble mass-vapor values at quadrature nodes from the conserved moment statistics.
subroutine, public s_initialize_pb(qk_cons_vf, mv, pb)
Initialize bubble internal pressures at quadrature nodes using isothermal relations from the Preston ...
impure subroutine, public s_convert_primitive_to_conservative_variables(q_prim_vf, q_cons_vf)
Convert primitives (rho, u, p, alpha) to conserved variables (rho*alpha, rho*u, E,...
subroutine, public s_convert_conservative_to_primitive_variables(qk_cons_vf, q_t_sf, qk_prim_vf, ibounds)
Convert conserved variables (rho*alpha, rho*u, E, alpha) to primitives (rho, u, p,...
Initial-condition state assembled by pre_process: working primitive and conservative fields,...