MFC
Exascale flow solver
Loading...
Searching...
No Matches
m_boundary_conditions.fpp.f90
Go to the documentation of this file.
1# 1 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
2!>
3!! @file
4!! @brief Contains module m_boundary_conditions
5
6!> @brief Applies spatially varying boundary condition patches along domain edges and faces
8
11#ifdef MFC_MPI
12 use mpi
13#endif
16
17 implicit none
18
21 real(wp) :: radius
23 private; public :: s_apply_boundary_patches
24
25contains
26 !> Apply a line-segment boundary condition patch along a domain edge in 2D.
27 impure subroutine s_line_segment_bc(patch_id, bc_type)
28
29 type(integer_field), dimension(1:num_dims,1:2), intent(inout) :: bc_type
30 integer, intent(in) :: patch_id
31 integer :: j
32
33 ! Patch is a line segment along y on the x-boundary face
34
35 if (patch_bc(patch_id)%dir == 1) then
36 y_centroid = patch_bc(patch_id)%centroid(2)
37 length_y = patch_bc(patch_id)%length(2)
38
39 y_boundary%beg = y_centroid - 0.5_wp*length_y
40 y_boundary%end = y_centroid + 0.5_wp*length_y
41
42 ! Apply patch if x boundary is a domain boundary
43# 43 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
44 if (patch_bc(patch_id)%loc == -1 .and. bc_x%beg < 0) then
45 do j = 0, n
46 if (y_cc(j) > y_boundary%beg .and. y_cc(j) < y_boundary%end) then
47 bc_type(1, 1)%sf(0, j, 0) = patch_bc(patch_id)%type
48 end if
49 end do
50 end if
51# 43 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
52 if (patch_bc(patch_id)%loc == 1 .and. bc_x%end < 0) then
53 do j = 0, n
54 if (y_cc(j) > y_boundary%beg .and. y_cc(j) < y_boundary%end) then
55 bc_type(1, 2)%sf(0, j, 0) = patch_bc(patch_id)%type
56 end if
57 end do
58 end if
59# 51 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
60 end if
61
62 ! Patch is a line segment along x on the y-boundary face
63 if (patch_bc(patch_id)%dir == 2) then
64 x_centroid = patch_bc(patch_id)%centroid(1)
65 length_x = patch_bc(patch_id)%length(1)
66
67 x_boundary%beg = x_centroid - 0.5_wp*length_x
68 x_boundary%end = x_centroid + 0.5_wp*length_x
69
70 ! Apply patch if y boundary is a domain boundary
71# 63 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
72 if (patch_bc(patch_id)%loc == -1 .and. bc_y%beg < 0) then
73 do j = 0, m
74 if (x_cc(j) > x_boundary%beg .and. x_cc(j) < x_boundary%end) then
75 bc_type(2, 1)%sf(j, 0, 0) = patch_bc(patch_id)%type
76 end if
77 end do
78 end if
79# 63 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
80 if (patch_bc(patch_id)%loc == 1 .and. bc_y%end < 0) then
81 do j = 0, m
82 if (x_cc(j) > x_boundary%beg .and. x_cc(j) < x_boundary%end) then
83 bc_type(2, 2)%sf(j, 0, 0) = patch_bc(patch_id)%type
84 end if
85 end do
86 end if
87# 71 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
88 end if
89
90 end subroutine s_line_segment_bc
91
92 !> Apply a circular boundary condition patch on a domain face in 3D.
93 impure subroutine s_circle_bc(patch_id, bc_type)
94
95 type(integer_field), dimension(1:num_dims,1:2), intent(inout) :: bc_type
96 integer, intent(in) :: patch_id
97 integer :: j, k
98
99 if (patch_bc(patch_id)%dir == 1) then
100 y_centroid = patch_bc(patch_id)%centroid(2)
101 z_centroid = patch_bc(patch_id)%centroid(3)
102 radius = patch_bc(patch_id)%radius
103 ! Patch is a circle at x_beg and x_beg is a domain boundary
104# 88 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
105 if (patch_bc(patch_id)%loc == -1 .and. bc_x%beg < 0) then
106 do k = 0, p
107 do j = 0, n
108 if ((z_cc(k) - z_centroid)**2._wp + (y_cc(j) - y_centroid)**2._wp <= radius**2._wp) then
109 bc_type(1, 1)%sf(0, j, k) = patch_bc(patch_id)%type
110 end if
111 end do
112 end do
113 end if
114# 88 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
115 if (patch_bc(patch_id)%loc == 1 .and. bc_x%end < 0) then
116 do k = 0, p
117 do j = 0, n
118 if ((z_cc(k) - z_centroid)**2._wp + (y_cc(j) - y_centroid)**2._wp <= radius**2._wp) then
119 bc_type(1, 2)%sf(0, j, k) = patch_bc(patch_id)%type
120 end if
121 end do
122 end do
123 end if
124# 98 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
125 end if
126 if (patch_bc(patch_id)%dir == 2) then
127 x_centroid = patch_bc(patch_id)%centroid(1)
128 z_centroid = patch_bc(patch_id)%centroid(3)
129 radius = patch_bc(patch_id)%radius
130 ! Patch is a circle at y_beg and y_beg is a domain boundary
131# 105 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
132 if (patch_bc(patch_id)%loc == -1 .and. bc_y%beg < 0) then
133 do k = 0, p
134 do j = 0, m
135 if ((z_cc(k) - z_centroid)**2._wp + (x_cc(j) - x_centroid)**2._wp <= radius**2._wp) then
136 bc_type(2, 1)%sf(j, 0, k) = patch_bc(patch_id)%type
137 end if
138 end do
139 end do
140 end if
141# 105 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
142 if (patch_bc(patch_id)%loc == 1 .and. bc_y%end < 0) then
143 do k = 0, p
144 do j = 0, m
145 if ((z_cc(k) - z_centroid)**2._wp + (x_cc(j) - x_centroid)**2._wp <= radius**2._wp) then
146 bc_type(2, 2)%sf(j, 0, k) = patch_bc(patch_id)%type
147 end if
148 end do
149 end do
150 end if
151# 115 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
152 end if
153 if (patch_bc(patch_id)%dir == 3) then
154 x_centroid = patch_bc(patch_id)%centroid(1)
155 y_centroid = patch_bc(patch_id)%centroid(2)
156 radius = patch_bc(patch_id)%radius
157# 121 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
158 if (patch_bc(patch_id)%loc == -1 .and. bc_z%beg < 0) then
159 do k = 0, n
160 do j = 0, m
161 if ((y_cc(k) - y_centroid)**2._wp + (x_cc(j) - x_centroid)**2._wp <= radius**2._wp) then
162 bc_type(3, 1)%sf(j, k, 0) = patch_bc(patch_id)%type
163 end if
164 end do
165 end do
166 end if
167# 121 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
168 if (patch_bc(patch_id)%loc == 1 .and. bc_z%end < 0) then
169 do k = 0, n
170 do j = 0, m
171 if ((y_cc(k) - y_centroid)**2._wp + (x_cc(j) - x_centroid)**2._wp <= radius**2._wp) then
172 bc_type(3, 2)%sf(j, k, 0) = patch_bc(patch_id)%type
173 end if
174 end do
175 end do
176 end if
177# 131 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
178 end if
179
180 end subroutine s_circle_bc
181
182 !> Apply a rectangular boundary condition patch on a domain face in 3D.
183 impure subroutine s_rectangle_bc(patch_id, bc_type)
184
185 type(integer_field), dimension(1:num_dims,1:2), intent(inout) :: bc_type
186 integer, intent(in) :: patch_id
187 integer :: j, k
188
189 if (patch_bc(patch_id)%dir == 1) then
190 y_centroid = patch_bc(patch_id)%centroid(2)
191 z_centroid = patch_bc(patch_id)%centroid(3)
192 length_y = patch_bc(patch_id)%length(2)
193 length_z = patch_bc(patch_id)%length(3)
194
195 y_boundary%beg = y_centroid - 0.5_wp*length_y
196 y_boundary%end = y_centroid + 0.5_wp*length_y
197
198 z_boundary%beg = z_centroid - 0.5_wp*length_z
199 z_boundary%end = z_centroid + 0.5_wp*length_z
200 ! Patch is a rectangle on the x-boundary face
201# 155 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
202 if (patch_bc(patch_id)%loc == -1 .and. bc_x%beg < 0) then
203 do k = 0, p
204 do j = 0, n
205 if (y_boundary%beg <= y_cc(j) .and. y_boundary%end >= y_cc(j) .and. z_boundary%beg <= z_cc(k) &
206 & .and. z_boundary%end >= z_cc(k)) then
207 bc_type(1, 1)%sf(0, j, k) = patch_bc(patch_id)%type
208 end if
209 end do
210 end do
211 end if
212# 155 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
213 if (patch_bc(patch_id)%loc == 1 .and. bc_x%end < 0) then
214 do k = 0, p
215 do j = 0, n
216 if (y_boundary%beg <= y_cc(j) .and. y_boundary%end >= y_cc(j) .and. z_boundary%beg <= z_cc(k) &
217 & .and. z_boundary%end >= z_cc(k)) then
218 bc_type(1, 2)%sf(0, j, k) = patch_bc(patch_id)%type
219 end if
220 end do
221 end do
222 end if
223# 166 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
224 end if
225 if (patch_bc(patch_id)%dir == 2) then
226 x_centroid = patch_bc(patch_id)%centroid(1)
227 z_centroid = patch_bc(patch_id)%centroid(3)
228 length_x = patch_bc(patch_id)%length(1)
229 length_z = patch_bc(patch_id)%length(3)
230
231 x_boundary%beg = x_centroid - 0.5_wp*length_x
232 x_boundary%end = x_centroid + 0.5_wp*length_x
233
234 z_boundary%beg = z_centroid - 0.5_wp*length_z
235 z_boundary%end = z_centroid + 0.5_wp*length_z
236 ! Patch is a rectangle on the y-boundary face
237# 180 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
238 if (patch_bc(patch_id)%loc == -1 .and. bc_y%beg < 0) then
239 do k = 0, p
240 do j = 0, m
241 if (x_boundary%beg <= x_cc(j) .and. x_boundary%end >= x_cc(j) .and. z_boundary%beg <= z_cc(k) &
242 & .and. z_boundary%end >= z_cc(k)) then
243 bc_type(2, 1)%sf(j, 0, k) = patch_bc(patch_id)%type
244 end if
245 end do
246 end do
247 end if
248# 180 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
249 if (patch_bc(patch_id)%loc == 1 .and. bc_y%end < 0) then
250 do k = 0, p
251 do j = 0, m
252 if (x_boundary%beg <= x_cc(j) .and. x_boundary%end >= x_cc(j) .and. z_boundary%beg <= z_cc(k) &
253 & .and. z_boundary%end >= z_cc(k)) then
254 bc_type(2, 2)%sf(j, 0, k) = patch_bc(patch_id)%type
255 end if
256 end do
257 end do
258 end if
259# 191 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
260 end if
261 if (patch_bc(patch_id)%dir == 3) then
262 x_centroid = patch_bc(patch_id)%centroid(1)
263 y_centroid = patch_bc(patch_id)%centroid(2)
264 length_x = patch_bc(patch_id)%length(1)
265 length_y = patch_bc(patch_id)%length(2)
266
267 x_boundary%beg = x_centroid - 0.5_wp*length_x
268 x_boundary%end = x_centroid + 0.5_wp*length_x
269
270 y_boundary%beg = y_centroid - 0.5_wp*length_y
271 y_boundary%end = y_centroid + 0.5_wp*length_y
272# 204 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
273 if (patch_bc(patch_id)%loc == -1 .and. bc_z%beg < 0) then
274 do k = 0, n
275 do j = 0, m
276 if (x_boundary%beg <= x_cc(j) .and. x_boundary%end >= x_cc(j) .and. y_boundary%beg <= y_cc(k) &
277 & .and. y_boundary%end >= y_cc(k)) then
278 bc_type(3, 1)%sf(j, k, 0) = patch_bc(patch_id)%type
279 end if
280 end do
281 end do
282 end if
283# 204 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
284 if (patch_bc(patch_id)%loc == 1 .and. bc_z%end < 0) then
285 do k = 0, n
286 do j = 0, m
287 if (x_boundary%beg <= x_cc(j) .and. x_boundary%end >= x_cc(j) .and. y_boundary%beg <= y_cc(k) &
288 & .and. y_boundary%end >= y_cc(k)) then
289 bc_type(3, 2)%sf(j, k, 0) = patch_bc(patch_id)%type
290 end if
291 end do
292 end do
293 end if
294# 215 "/home/runner/work/MFC/MFC/src/pre_process/m_boundary_conditions.fpp"
295 end if
296
297 end subroutine s_rectangle_bc
298
299 !> Iterate over all boundary condition patches and dispatch them by geometry type.
300 impure subroutine s_apply_boundary_patches(q_prim_vf, bc_type)
301
302 type(scalar_field), dimension(sys_size) :: q_prim_vf
303 type(integer_field), dimension(1:num_dims,1:2) :: bc_type
304 integer :: i
305
306 !> Apply 2D patches to 3D domain
307
308 if (p > 0) then
309 do i = 1, num_bc_patches
310 if (proc_rank == 0) then
311 print *, 'Processing boundary condition patch', i
312 end if
313
314 if (patch_bc(i)%geometry == 2) then
315 call s_circle_bc(i, bc_type)
316 else if (patch_bc(i)%geometry == 3) then
317 call s_rectangle_bc(i, bc_type)
318 end if
319 end do
320 !> Apply 1D patches to 2D domain
321 else if (n > 0) then
322 do i = 1, num_bc_patches
323 if (proc_rank == 0) then
324 print *, 'Processing boundary condition patch', i
325 end if
326
327 if (patch_bc(i)%geometry == 1) then
328 call s_line_segment_bc(i, bc_type)
329 end if
330 end do
331 end if
332
333 end subroutine s_apply_boundary_patches
334
335end module m_boundary_conditions
integer, intent(in) k
integer, intent(in) j
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.
impure subroutine s_circle_bc(patch_id, bc_type)
Apply a circular boundary condition patch on a domain face in 3D.
impure subroutine s_line_segment_bc(patch_id, bc_type)
Apply a line-segment boundary condition patch along a domain edge in 2D.
impure subroutine s_rectangle_bc(patch_id, bc_type)
Apply a rectangular boundary condition patch on a domain face in 3D.
Platform-specific file and directory operations: create, delete, inquire, getcwd, and basename.
Rank-staggered file access delays to prevent I/O contention on parallel file systems.
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) bc_z
Boundary conditions in the x-, y- and z-coordinate directions.
real(wp), dimension(:), allocatable y_cc
integer proc_rank
Rank of the local processor Number of cells in the x-, y- and z-coordinate directions.
type(int_bounds_info) bc_y
real(wp), dimension(:), allocatable x_cc
Locations of cell-centers (cc) in x-, y- and z-directions, respectively.
type(int_bounds_info) bc_x
real(wp), dimension(:), allocatable z_cc
Derived type adding beginning (beg) and end bounds info as attributes.
Derived type annexing an integer scalar field (SF).
Derived type annexing a scalar field (SF).