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