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
9 use m_derived_types ! Definitions of the derived types
10
11 use m_global_parameters ! Global parameters for the code
12
13 use m_mpi_proxy !< message passing interface (mpi) module proxy
14
15 use m_helper
16
17 use m_variables_conversion ! Subroutines to change the state variables from
18 ! one form to another
19
21
23
24 use m_perturbation ! Subroutines to perturb initial flow fields
25
26 use m_chemistry
27
29
30 implicit none
31
32 ! NOTE: The abstract interface allows for the declaration of a pointer to
33 ! a procedure such that the choice of the model equations does not have to
34 ! be queried every time the patch primitive variables are to be assigned in
35 ! a cell in the computational domain.
36 type(scalar_field), allocatable, dimension(:) :: q_prim_vf !< primitive variables
37
38 type(scalar_field), allocatable, dimension(:) :: q_cons_vf !< conservative variables
39
40 type(scalar_field) :: q_t_sf !< Temperature field
41
42 type(integer_field), dimension(:, :), allocatable :: bc_type !< bc_type fields
43
44!> @cond
45#ifdef MFC_MIXED_PRECISION
46 integer(kind=1), allocatable, dimension(:, :, :) :: patch_id_fp
47#else
48!> @endcond
49 integer, allocatable, dimension(:, :, :) :: patch_id_fp
50!> @cond
51#endif
52!> @endcond
53
54contains
55
56 !> Computation of parameters, allocation procedures, and/or
57 !! any other tasks needed to properly setup the module
59
60 integer :: i, j, k, l !< generic loop iterators
61
62 ! Allocating the primitive and conservative variables
63 allocate (q_prim_vf(1:sys_size))
64 allocate (q_cons_vf(1:sys_size))
65
66 do i = 1, sys_size
67 allocate (q_prim_vf(i)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
68 idwbuff(2)%beg:idwbuff(2)%end, &
69 idwbuff(3)%beg:idwbuff(3)%end))
70 allocate (q_cons_vf(i)%sf(idwbuff(1)%beg:idwbuff(1)%end, &
71 idwbuff(2)%beg:idwbuff(2)%end, &
72 idwbuff(3)%beg:idwbuff(3)%end))
73 end do
74
75 if (chemistry) then
76 allocate (q_t_sf%sf(0:m, 0:n, 0:p))
77 end if
78
79 ! Allocating the patch identities bookkeeping variable
80 allocate (patch_id_fp(0:m, 0:n, 0:p))
81
82 if (qbmm .and. .not. polytropic) then
83 !Allocate bubble pressure pb and vapor mass mv for non-polytropic qbmm at all quad nodes and R0 bins
84 allocate (pb%sf(0:m, &
85 0:n, &
86 0:p, 1:nnode, 1:nb))
87 allocate (mv%sf(0:m, &
88 0:n, &
89 0:p, 1:nnode, 1:nb))
90 end if
91
92 ! Setting default values for conservative and primitive variables so
93 ! that in the case that the initial condition is wrongly laid out on
94 ! the grid the simulation component will catch the problem on start-
95 ! up. The conservative variables do not need to be similarly treated
96 ! since they are computed directly from the primitive variables.
97 do i = 1, sys_size
98 q_cons_vf(i)%sf = -1.e-6_stp ! real(dflt_real, kind=stp) ! TODO :: remove this magic number
99 q_prim_vf(i)%sf = -1.e-6_stp ! real(dflt_real, kind=stp)
100 end do
101
102 ! Allocating arrays to store the bc types
103 allocate (bc_type(1:num_dims, 1:2))
104
105 allocate (bc_type(1, 1)%sf(0:0, 0:n, 0:p))
106 allocate (bc_type(1, 2)%sf(0:0, 0:n, 0:p))
107
108 do l = 0, p
109 do k = 0, n
110 bc_type(1, 1)%sf(0, k, l) = int(min(bc_x%beg, 0), kind=1)
111 bc_type(1, 2)%sf(0, k, l) = int(min(bc_x%end, 0), kind=1)
112 end do
113 end do
114
115 if (n > 0) then
116 allocate (bc_type(2, 1)%sf(-buff_size:m + buff_size, 0:0, 0:p))
117 allocate (bc_type(2, 2)%sf(-buff_size:m + buff_size, 0:0, 0:p))
118
119 do l = 0, p
120 do j = -buff_size, m + buff_size
121 bc_type(2, 1)%sf(j, 0, l) = int(min(bc_y%beg, 0), kind=1)
122 bc_type(2, 2)%sf(j, 0, l) = int(min(bc_y%end, 0), kind=1)
123 end do
124 end do
125
126 if (p > 0) then
127 allocate (bc_type(3, 1)%sf(-buff_size:m + buff_size, -buff_size:n + buff_size, 0:0))
128 allocate (bc_type(3, 2)%sf(-buff_size:m + buff_size, -buff_size:n + buff_size, 0:0))
129
130 do k = -buff_size, n + buff_size
131 do j = -buff_size, m + buff_size
132 bc_type(3, 1)%sf(j, k, 0) = int(min(bc_z%beg, 0), kind=1)
133 bc_type(3, 2)%sf(j, k, 0) = int(min(bc_z%end, 0), kind=1)
134 end do
135 end do
136 end if
137 end if
138
139 ! Initial damage state is always zero
140 if (cont_damage) then
141 q_cons_vf(damage_idx)%sf = 0._wp
142 q_prim_vf(damage_idx)%sf = 0._wp
143 end if
144
145 ! Initial hyper_cleaning state is always zero
146 ! TODO more general
147 if (hyper_cleaning) then
148 q_cons_vf(psi_idx)%sf = 0._wp
149 q_prim_vf(psi_idx)%sf = 0._wp
150 end if
151
152 ! Setting default values for patch identities bookkeeping variable.
153 ! This is necessary to avoid any confusion in the assessment of the
154 ! extent of application that the overwrite permissions give a patch
155 ! when it is being applied in the domain.
156 patch_id_fp = 0
157
159
160 !> This subroutine peruses the patches and depending on the
161 !! type of geometry associated with a particular patch, it
162 !! calls the related subroutine to setup the said geometry
163 !! on the grid using the primitive variables included with
164 !! the patch parameters. The subroutine is complete once the
165 !! primitive variables are converted to conservative ones.
167
168 integer :: i
169
170 ! Converting the conservative variables to the primitive ones given
171 ! preexisting initial condition data files were read in on start-up
172 if (old_ic) then
174 q_t_sf, &
175 q_prim_vf, &
176 idwbuff)
177 end if
178
180
182
188
189 ! Converting the primitive variables to the conservative ones
191
193
194 if (qbmm .and. .not. polytropic) then
195 !Initialize pb and mv
197 call s_initialize_pb(q_cons_vf, mv%sf, pb%sf)
198 end if
199
200 end subroutine s_generate_initial_condition
201
202 !> Deallocation procedures for the module
204
205 integer :: i !< Generic loop iterator
206
207 ! Dellocating the primitive and conservative variables
208 do i = 1, sys_size
209 deallocate (q_prim_vf(i)%sf)
210 deallocate (q_cons_vf(i)%sf)
211 end do
212
213 deallocate (q_prim_vf)
214 deallocate (q_cons_vf)
215
216 if (chemistry) then
217 deallocate (q_t_sf%sf)
218 end if
219
220 ! Deallocating the patch identities bookkeeping variable
221 deallocate (patch_id_fp)
222
223 deallocate (bc_type(1, 1)%sf)
224 deallocate (bc_type(1, 2)%sf)
225
226 if (n > 0) then
227 deallocate (bc_type(2, 1)%sf)
228 deallocate (bc_type(2, 2)%sf)
229 end if
230
231 if (p > 0) then
232 deallocate (bc_type(3, 1)%sf)
233 deallocate (bc_type(3, 2)%sf)
234 end if
235
236 deallocate (bc_type)
237
239
240end 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)
Iterates over all boundary condition patches and dispatches 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)
Computes 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
The number of cells that are necessary to be able to store enough boundary conditions data to march t...
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)
Dispatches each initial condition patch to its geometry-specific initialization routine.
Assembles initial conditions by layering prioritized patches via constructive solid geometry.
type(scalar_field), dimension(:), allocatable q_cons_vf
conservative variables
type(integer_field), dimension(:, :), allocatable bc_type
bc_type fields
integer, dimension(:, :, :), allocatable patch_id_fp
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(scalar_field) q_t_sf
Temperature field.
type(scalar_field), dimension(:), allocatable q_prim_vf
primitive variables
impure subroutine s_generate_initial_condition
This subroutine peruses the patches and depending on the type of geometry associated with a particula...
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 smooths all primitive variable fields using a discrete elliptic (Laplacian) filter.
impure subroutine s_perturb_sphere(q_prim_vf)
Randomly perturbs partial density fields at the interface of a spherical volume fraction region.
impure subroutine s_perturb_surrounding_flow(q_prim_vf)
Adds random noise to the velocity and void fraction of the surrounding flow field.
subroutine s_perturb_simplex(q_prim_vf)
Perturbs velocity and volume fraction fields using multi-octave simplex noise.
subroutine s_perturb_mixlayer(q_prim_vf)
This subroutine computes velocity perturbations for a temporal mixing layer with a hyperbolic tangent...
Conservative-to-primitive variable conversion, mixture property evaluation, and pressure computation.
subroutine, public s_initialize_mv(qk_cons_vf, mv)
Initializes bubble mass-vapor values at quadrature nodes from the conserved moment statistics.
subroutine, public s_initialize_pb(qk_cons_vf, mv, pb)
Initializes 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)
The following procedure handles the conversion between the primitive variables and the conservative v...
subroutine, public s_convert_conservative_to_primitive_variables(qk_cons_vf, q_t_sf, qk_prim_vf, ibounds)
The following procedure handles the conversion between the conservative variables and the primitive v...
Derived type annexing an integer scalar field (SF).
Derived type annexing a scalar field (SF).