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