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 ! NOTE: Abstract interface enables dynamic dispatch without repeated model_eqns checks
23 type(scalar_field), allocatable, dimension(:) :: q_prim_vf !< primitive variables
24 type(scalar_field), allocatable, dimension(:) :: q_cons_vf !< conservative variables
25 type(scalar_field) :: q_t_sf !< Temperature field
26 type(integer_field), dimension(:,:), allocatable :: bc_type !< bc_type fields
27 !> @cond
28#ifdef MFC_MIXED_PRECISION
29 integer(kind=1), allocatable, dimension(:,:,:) :: patch_id_fp
30#else
31 !> @endcond
32 integer, allocatable, dimension(:,:,:) :: patch_id_fp
33 !> @cond
34#endif
35 !> @endcond
36
37contains
38
39 !> Computation of parameters, allocation procedures, and/or any other tasks needed to properly setup the module
41
42 integer :: i, j, k, l
43
44 allocate (q_prim_vf(1:sys_size))
45 allocate (q_cons_vf(1:sys_size))
46
47 do i = 1, sys_size
48 allocate (q_prim_vf(i)%sf(idwbuff(1)%beg:idwbuff(1)%end,idwbuff(2)%beg:idwbuff(2)%end,idwbuff(3)%beg:idwbuff(3)%end))
49 allocate (q_cons_vf(i)%sf(idwbuff(1)%beg:idwbuff(1)%end,idwbuff(2)%beg:idwbuff(2)%end,idwbuff(3)%beg:idwbuff(3)%end))
50 end do
51
52 if (chemistry) then
53 allocate (q_t_sf%sf(0:m,0:n,0:p))
54 end if
55
56 allocate (patch_id_fp(0:m,0:n,0:p))
57
58 if (qbmm .and. .not. polytropic) then
59 allocate (pb%sf(0:m,0:n,0:p,1:nnode,1:nb))
60 allocate (mv%sf(0:m,0:n,0:p,1:nnode,1:nb))
61 end if
62
63 do i = 1, sys_size
64 q_cons_vf(i)%sf = -1.e-6_stp ! real(dflt_real, kind=stp) ! TODO :: remove this magic number
65 q_prim_vf(i)%sf = -1.e-6_stp ! real(dflt_real, kind=stp)
66 end do
67
68 allocate (bc_type(1:num_dims,1:2))
69
70 allocate (bc_type(1, 1)%sf(0:0,0:n,0:p))
71 allocate (bc_type(1, 2)%sf(0:0,0:n,0:p))
72
73 do l = 0, p
74 do k = 0, n
75 bc_type(1, 1)%sf(0, k, l) = int(min(bc_x%beg, 0), kind=1)
76 bc_type(1, 2)%sf(0, k, l) = int(min(bc_x%end, 0), kind=1)
77 end do
78 end do
79
80 if (n > 0) then
81 allocate (bc_type(2, 1)%sf(-buff_size:m + buff_size,0:0,0:p))
82 allocate (bc_type(2, 2)%sf(-buff_size:m + buff_size,0:0,0:p))
83
84 do l = 0, p
85 do j = -buff_size, m + buff_size
86 bc_type(2, 1)%sf(j, 0, l) = int(min(bc_y%beg, 0), kind=1)
87 bc_type(2, 2)%sf(j, 0, l) = int(min(bc_y%end, 0), kind=1)
88 end do
89 end do
90
91 if (p > 0) then
92 allocate (bc_type(3, 1)%sf(-buff_size:m + buff_size,-buff_size:n + buff_size,0:0))
93 allocate (bc_type(3, 2)%sf(-buff_size:m + buff_size,-buff_size:n + buff_size,0:0))
94
95 do k = -buff_size, n + buff_size
96 do j = -buff_size, m + buff_size
97 bc_type(3, 1)%sf(j, k, 0) = int(min(bc_z%beg, 0), kind=1)
98 bc_type(3, 2)%sf(j, k, 0) = int(min(bc_z%end, 0), kind=1)
99 end do
100 end do
101 end if
102 end if
103
104 ! Initial damage state is always zero
105 if (cont_damage) then
106 q_cons_vf(damage_idx)%sf = 0._wp
107 q_prim_vf(damage_idx)%sf = 0._wp
108 end if
109
110 ! Initial hyper_cleaning state is always zero TODO more general
111 if (hyper_cleaning) then
112 q_cons_vf(psi_idx)%sf = 0._wp
113 q_prim_vf(psi_idx)%sf = 0._wp
114 end if
115
116 ! Setting default values for patch identities bookkeeping variable. This is necessary to avoid any confusion in the
117 ! assessment of the extent of application that the overwrite permissions give a patch when it is being applied in the
118 ! domain.
119 patch_id_fp = 0
120
122
123 !> Iterate over patches and, depending on the geometry type, call the related subroutine to setup the said geometry on the grid
124 !! using the primitive variables included with the patch parameters. The subroutine is complete once the primitive variables are
125 !! converted to conservative ones.
153 end subroutine s_generate_initial_condition
154
155 !> Deallocation procedures for the module
157
158 integer :: i
159
160 do i = 1, sys_size
161 deallocate (q_prim_vf(i)%sf)
162 deallocate (q_cons_vf(i)%sf)
163 end do
164
165 deallocate (q_prim_vf)
166 deallocate (q_cons_vf)
167
168 if (chemistry) then
169 deallocate (q_t_sf%sf)
170 end if
171
172 deallocate (patch_id_fp)
173
174 deallocate (bc_type(1, 1)%sf)
175 deallocate (bc_type(1, 2)%sf)
176
177 if (n > 0) then
178 deallocate (bc_type(2, 1)%sf)
179 deallocate (bc_type(2, 2)%sf)
180 end if
181
182 if (p > 0) then
183 deallocate (bc_type(3, 1)%sf)
184 deallocate (bc_type(3, 2)%sf)
185 end if
186
187 deallocate (bc_type)
188
190
191end 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.
logical cont_damage
continuum damage modeling
type(int_bounds_info), dimension(1:3) idwint
logical, parameter chemistry
Chemistry modeling.
type(int_bounds_info) bc_z
Boundary conditions in the x-, y- and z-coordinate directions.
integer sys_size
Number of unknowns in the system of equations.
integer num_bc_patches
Number of boundary condition patches.
type(int_bounds_info), dimension(1:3) idwbuff
integer buff_size
Number of ghost cells for boundary condition storage.
integer num_dims
Number of spatial dimensions.
logical qbmm
Quadrature moment method.
integer damage_idx
Index of damage state variable (D) for continuum damage model.
logical hyper_cleaning
Hyperbolic cleaning for MHD.
integer psi_idx
Index of hyperbolic cleaning state variable for MHD.
logical mixlayer_perturb
Superimpose instability waves to surrounding fluid flow.
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.
integer, dimension(:,:,:), allocatable patch_id_fp
type(scalar_field), dimension(:), allocatable q_cons_vf
conservative variables
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.
type(integer_field), dimension(:,:), allocatable bc_type
bc_type fields
type(scalar_field) q_t_sf
Temperature field.
type(scalar_field), dimension(:), allocatable q_prim_vf
primitive variables
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_elliptic_smoothing(q_prim_vf, bc_type)
Iteratively smooth all primitive variable fields using a discrete elliptic (Laplacian) filter.
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.
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,...
Derived type annexing an integer scalar field (SF).
Derived type annexing a scalar field (SF).