MFC
Exascale flow solver
Loading...
Searching...
No Matches
m_boundary_primitives.fpp.f90
Go to the documentation of this file.
1# 1 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
2!>
3!! @file
4!! @brief Contains module m_boundary_primitives
5
6!> @brief Per-cell noncharacteristic boundary condition primitives applied in the ghost cells
7# 1 "/home/runner/work/MFC/MFC/src/common/include/case.fpp" 1
8! This file exists so that Fypp can be run without generating case.fpp files for
9! each target. This is useful when generating documentation, for example. This
10! should also let MFC be built with CMake directly, without invoking mfc.sh.
11
12! For pre-process.
13# 8 "/home/runner/work/MFC/MFC/src/common/include/case.fpp"
14
15! For moving immersed boundaries in simulation
16# 12 "/home/runner/work/MFC/MFC/src/common/include/case.fpp"
17# 7 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp" 2
18# 1 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp" 1
19# 1 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp" 1
20# 1 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp" 1
21# 2 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
22# 3 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
23# 4 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
24# 5 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
25# 6 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
26
27# 8 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
28# 9 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
29# 10 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
30
31# 17 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
32
33# 46 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
34
35# 58 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
36
37# 68 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
38
39# 98 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
40
41# 110 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
42
43# 120 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
44! New line at end of file is required for FYPP
45# 2 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp" 2
46# 1 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp" 1
47# 1 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp" 1
48# 2 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
49# 3 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
50# 4 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
51# 5 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
52# 6 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
53
54# 8 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
55# 9 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
56# 10 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
57
58# 17 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
59
60# 46 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
61
62# 58 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
63
64# 68 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
65
66# 98 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
67
68# 110 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
69
70# 120 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
71! New line at end of file is required for FYPP
72# 2 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp" 2
73
74# 4 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
75# 5 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
76# 6 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
77# 7 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
78# 8 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
79
80# 20 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
81
82# 43 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
83
84# 48 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
85
86# 53 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
87
88# 58 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
89
90# 63 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
91
92# 68 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
93
94# 76 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
95
96# 81 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
97
98# 86 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
99
100# 91 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
101
102# 96 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
103
104# 101 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
105
106# 106 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
107
108# 111 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
109
110# 116 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
111
112# 121 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
113
114# 151 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
115
116# 192 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
117
118# 206 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
119
120# 231 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
121
122# 242 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
123
124# 244 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
125# 255 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
126
127# 284 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
128
129# 294 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
130
131# 304 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
132
133# 313 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
134
135# 330 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
136
137# 340 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
138
139# 347 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
140
141# 353 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
142
143# 359 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
144
145# 365 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
146
147# 371 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
148
149# 377 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
150! New line at end of file is required for FYPP
151# 3 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp" 2
152# 1 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp" 1
153# 1 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp" 1
154# 2 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
155# 3 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
156# 4 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
157# 5 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
158# 6 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
159
160# 8 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
161# 9 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
162# 10 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
163
164# 17 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
165
166# 46 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
167
168# 58 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
169
170# 68 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
171
172# 98 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
173
174# 110 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
175
176# 120 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
177! New line at end of file is required for FYPP
178# 2 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp" 2
179
180# 7 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
181
182# 17 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
183
184# 22 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
185
186# 27 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
187
188# 32 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
189
190# 37 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
191
192# 42 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
193
194# 47 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
195
196# 52 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
197
198# 57 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
199
200# 62 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
201
202# 73 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
203
204# 78 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
205
206# 83 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
207
208# 88 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
209
210# 103 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
211
212# 131 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
213
214# 160 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
215
216# 175 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
217
218# 193 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
219
220# 215 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
221
222# 244 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
223
224# 259 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
225
226# 269 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
227
228# 278 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
229
230# 294 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
231
232# 304 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
233
234# 311 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
235! New line at end of file is required for FYPP
236# 4 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp" 2
237
238! GPU parallel region (scalar reductions, maxval/minval)
239# 23 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
240
241! GPU parallel loop over threads (most common GPU macro)
242# 43 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
243
244! Required closing for GPU_PARALLEL_LOOP
245# 55 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
246
247! Mark routine for device compilation
248# 112 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
249
250! Declare device-resident data
251# 130 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
252
253! Inner loop within a GPU parallel region
254# 145 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
255
256! Scoped GPU data region
257# 164 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
258
259! Host code with device pointers (for MPI with GPU buffers)
260# 193 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
261
262! Allocate device memory (unscoped)
263# 207 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
264
265! Free device memory
266# 219 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
267
268! Atomic operation on device
269# 231 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
270
271! End atomic capture block
272# 242 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
273
274! Copy data between host and device
275# 254 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
276
277! Synchronization barrier
278# 266 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
279
280! Import GPU library module (openacc or omp_lib)
281# 275 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
282
283! Emit code only for AMD compiler
284# 282 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
285
286! Emit code for non-Cray compilers
287# 289 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
288
289! Emit code only for Cray compiler
290# 296 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
291
292! Emit code for non-NVIDIA compilers
293# 303 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
294
295# 305 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
296# 306 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
297! New line at end of file is required for FYPP
298# 2 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp" 2
299
300# 14 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
301
302! Caution: This macro requires the use of a binding script to set CUDA_VISIBLE_DEVICES, such that we have one GPU device per MPI
303! rank. That's because for both cudaMemAdvise (preferred location) and cudaMemPrefetchAsync we use location = device_id = 0. For an
304! example see misc/nvidia_uvm/bind.sh. NVIDIA unified memory page placement hint
305# 57 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
306
307! Allocate and create GPU device memory
308# 77 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
309
310! Free GPU device memory and deallocate
311# 85 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
312
313! Cray-specific GPU pointer setup for vector fields
314# 109 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
315
316! Cray-specific GPU pointer setup for scalar fields
317# 125 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
318
319! Cray-specific GPU pointer setup for acoustic source spatials
320# 150 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
321
322# 156 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
323
324# 163 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
325! New line at end of file is required for FYPP
326# 8 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp" 2
327
329
332 use m_constants
333
334 implicit none
335
336 type(scalar_field), dimension(:,:), allocatable :: bc_buffers
337
338# 18 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
339#if defined(MFC_OpenACC)
340# 18 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
341!$acc declare create(bc_buffers)
342# 18 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
343#elif defined(MFC_OpenMP)
344# 18 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
345!$omp declare target (bc_buffers)
346# 18 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
347#endif
348
349contains
350
351 !> Fill ghost cells by copying the nearest boundary cell value along the specified direction.
352 subroutine s_ghost_cell_extrapolation(q_prim_vf, bc_dir, bc_loc, k, l, q_T_sf)
353
354
355# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
356#ifdef _CRAYFTN
357# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
358#if MFC_OpenACC
359# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
360!$acc routine seq
361# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
362#elif MFC_OpenMP
363# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
364
365# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
366
367# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
368!$omp declare target device_type(any)
369# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
370#else
371# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
372!DIR$ INLINEALWAYS s_ghost_cell_extrapolation
373# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
374#endif
375# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
376#elif MFC_OpenACC
377# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
378!$acc routine seq
379# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
380#elif MFC_OpenMP
381# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
382
383# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
384
385# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
386!$omp declare target device_type(any)
387# 25 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
388#endif
389 type(scalar_field), dimension(sys_size), intent(inout) :: q_prim_vf
390 integer, intent(in) :: bc_dir, bc_loc
391 integer, intent(in) :: k, l
392 integer :: j, i
393 type(scalar_field), optional, intent(inout) :: q_T_sf
394
395 if (bc_dir == 1) then !< x-direction
396 if (bc_loc == -1) then ! bc_x%beg
397 do i = 1, sys_size
398 do j = 1, buff_size
399 q_prim_vf(i)%sf(-j, k, l) = q_prim_vf(i)%sf(0, k, l)
400 end do
401 end do
402 if (chemistry .and. present(q_t_sf)) then
403 do j = 1, buff_size
404 q_t_sf%sf(-j, k, l) = q_t_sf%sf(0, k, l)
405 end do
406 end if
407 else !< bc_x%end
408 do i = 1, sys_size
409 do j = 1, buff_size
410 q_prim_vf(i)%sf(m + j, k, l) = q_prim_vf(i)%sf(m, k, l)
411 end do
412 end do
413 if (chemistry .and. present(q_t_sf)) then
414 do j = 1, buff_size
415 q_t_sf%sf(m + j, k, l) = q_t_sf%sf(m, k, l)
416 end do
417 end if
418 end if
419 else if (bc_dir == 2) then !< y-direction
420 if (bc_loc == -1) then !< bc_y%beg
421 do i = 1, sys_size
422 do j = 1, buff_size
423 q_prim_vf(i)%sf(k, -j, l) = q_prim_vf(i)%sf(k, 0, l)
424 end do
425 end do
426
427 if (chemistry .and. present(q_t_sf)) then
428 do j = 1, buff_size
429 q_t_sf%sf(k, -j, l) = q_t_sf%sf(k, 0, l)
430 end do
431 end if
432 else !< bc_y%end
433 do i = 1, sys_size
434 do j = 1, buff_size
435 q_prim_vf(i)%sf(k, n + j, l) = q_prim_vf(i)%sf(k, n, l)
436 end do
437 end do
438 if (chemistry .and. present(q_t_sf)) then
439 do j = 1, buff_size
440 q_t_sf%sf(k, n + j, l) = q_t_sf%sf(k, n, l)
441 end do
442 end if
443 end if
444 else if (bc_dir == 3) then !< z-direction
445 if (bc_loc == -1) then !< bc_z%beg
446 do i = 1, sys_size
447 do j = 1, buff_size
448 q_prim_vf(i)%sf(k, l, -j) = q_prim_vf(i)%sf(k, l, 0)
449 end do
450 end do
451 if (chemistry .and. present(q_t_sf)) then
452 do j = 1, buff_size
453 q_t_sf%sf(k, l, -j) = q_t_sf%sf(k, l, 0)
454 end do
455 end if
456 else !< bc_z%end
457 do i = 1, sys_size
458 do j = 1, buff_size
459 q_prim_vf(i)%sf(k, l, p + j) = q_prim_vf(i)%sf(k, l, p)
460 end do
461 end do
462 if (chemistry .and. present(q_t_sf)) then
463 do j = 1, buff_size
464 q_t_sf%sf(k, l, p + j) = q_t_sf%sf(k, l, p)
465 end do
466 end if
467 end if
468 end if
469
470 end subroutine s_ghost_cell_extrapolation
471
472 !> Apply reflective (symmetry) boundary conditions by mirroring primitive variables and flipping the normal velocity component.
473 subroutine s_symmetry(q_prim_vf, bc_dir, bc_loc, k, l, pb_in, mv_in, q_T_sf)
474
475
476# 112 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
477#if MFC_OpenACC
478# 112 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
479!$acc routine seq
480# 112 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
481#elif MFC_OpenMP
482# 112 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
483
484# 112 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
485
486# 112 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
487!$omp declare target device_type(any)
488# 112 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
489#endif
490 type(scalar_field), dimension(sys_size), intent(inout) :: q_prim_vf
491 real(stp), optional, dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:,1:), intent(inout) :: pb_in, mv_in
492 integer, intent(in) :: bc_dir, bc_loc
493 integer, intent(in) :: k, l
494 integer :: j, q, i
495 type(scalar_field), optional, intent(inout) :: q_T_sf
496
497 if (bc_dir == 1) then !< x-direction
498 if (bc_loc == -1) then !< bc_x%beg
499 do j = 1, buff_size
500 do i = 1, eqn_idx%cont%end
501 q_prim_vf(i)%sf(-j, k, l) = q_prim_vf(i)%sf(j - 1, k, l)
502 end do
503
504 q_prim_vf(eqn_idx%mom%beg)%sf(-j, k, l) = -q_prim_vf(eqn_idx%mom%beg)%sf(j - 1, k, l)
505
506 do i = eqn_idx%mom%beg + 1, sys_size
507 q_prim_vf(i)%sf(-j, k, l) = q_prim_vf(i)%sf(j - 1, k, l)
508 end do
509
510 if (chemistry .and. present(q_t_sf)) then
511 q_t_sf%sf(-j, k, l) = q_t_sf%sf(j - 1, k, l)
512 end if
513
514 if (elasticity) then
515 do i = 1, shear_bc_flip_num
516 q_prim_vf(shear_bc_flip_indices(1, i))%sf(-j, k, l) = -q_prim_vf(shear_bc_flip_indices(1, &
517 & i))%sf(j - 1, k, l)
518 end do
519 end if
520
521 if (hyperelasticity) then
522 q_prim_vf(eqn_idx%xi%beg)%sf(-j, k, l) = -q_prim_vf(eqn_idx%xi%beg)%sf(j - 1, k, l)
523 end if
524 end do
525
526 if (qbmm .and. .not. polytropic .and. present(pb_in) .and. present(mv_in)) then
527 do i = 1, nb
528 do q = 1, nnode
529 do j = 1, buff_size
530 pb_in(-j, k, l, q, i) = pb_in(j - 1, k, l, q, i)
531 mv_in(-j, k, l, q, i) = mv_in(j - 1, k, l, q, i)
532 end do
533 end do
534 end do
535 end if
536 else !< bc_x%end
537 do j = 1, buff_size
538 do i = 1, eqn_idx%cont%end
539 q_prim_vf(i)%sf(m + j, k, l) = q_prim_vf(i)%sf(m - (j - 1), k, l)
540 end do
541
542 q_prim_vf(eqn_idx%mom%beg)%sf(m + j, k, l) = -q_prim_vf(eqn_idx%mom%beg)%sf(m - (j - 1), k, l)
543
544 do i = eqn_idx%mom%beg + 1, sys_size
545 q_prim_vf(i)%sf(m + j, k, l) = q_prim_vf(i)%sf(m - (j - 1), k, l)
546 end do
547
548 if (chemistry .and. present(q_t_sf)) then
549 q_t_sf%sf(m + j, k, l) = q_t_sf%sf(m - (j - 1), k, l)
550 end if
551
552 if (elasticity) then
553 do i = 1, shear_bc_flip_num
554 q_prim_vf(shear_bc_flip_indices(1, i))%sf(m + j, k, l) = -q_prim_vf(shear_bc_flip_indices(1, &
555 & i))%sf(m - (j - 1), k, l)
556 end do
557 end if
558
559 if (hyperelasticity) then
560 q_prim_vf(eqn_idx%xi%beg)%sf(m + j, k, l) = -q_prim_vf(eqn_idx%xi%beg)%sf(m - (j - 1), k, l)
561 end if
562 end do
563 if (qbmm .and. .not. polytropic .and. present(pb_in) .and. present(mv_in)) then
564 do i = 1, nb
565 do q = 1, nnode
566 do j = 1, buff_size
567 pb_in(m + j, k, l, q, i) = pb_in(m - (j - 1), k, l, q, i)
568 mv_in(m + j, k, l, q, i) = mv_in(m - (j - 1), k, l, q, i)
569 end do
570 end do
571 end do
572 end if
573 end if
574 else if (bc_dir == 2) then !< y-direction
575 if (bc_loc == -1) then !< bc_y%beg
576 do j = 1, buff_size
577 do i = 1, eqn_idx%mom%beg
578 q_prim_vf(i)%sf(k, -j, l) = q_prim_vf(i)%sf(k, j - 1, l)
579 end do
580
581 q_prim_vf(eqn_idx%mom%beg + 1)%sf(k, -j, l) = -q_prim_vf(eqn_idx%mom%beg + 1)%sf(k, j - 1, l)
582
583 do i = eqn_idx%mom%beg + 2, sys_size
584 q_prim_vf(i)%sf(k, -j, l) = q_prim_vf(i)%sf(k, j - 1, l)
585 end do
586
587 if (chemistry .and. present(q_t_sf)) then
588 q_t_sf%sf(k, -j, l) = q_t_sf%sf(k, j - 1, l)
589 end if
590
591 if (elasticity) then
592 do i = 1, shear_bc_flip_num
593 q_prim_vf(shear_bc_flip_indices(2, i))%sf(k, -j, l) = -q_prim_vf(shear_bc_flip_indices(2, i))%sf(k, &
594 & j - 1, l)
595 end do
596 end if
597
598 if (hyperelasticity) then
599 q_prim_vf(eqn_idx%xi%beg + 1)%sf(k, -j, l) = -q_prim_vf(eqn_idx%xi%beg + 1)%sf(k, j - 1, l)
600 end if
601 end do
602
603 if (qbmm .and. .not. polytropic .and. present(pb_in) .and. present(mv_in)) then
604 do i = 1, nb
605 do q = 1, nnode
606 do j = 1, buff_size
607 pb_in(k, -j, l, q, i) = pb_in(k, j - 1, l, q, i)
608 mv_in(k, -j, l, q, i) = mv_in(k, j - 1, l, q, i)
609 end do
610 end do
611 end do
612 end if
613 else !< bc_y%end
614 do j = 1, buff_size
615 do i = 1, eqn_idx%mom%beg
616 q_prim_vf(i)%sf(k, n + j, l) = q_prim_vf(i)%sf(k, n - (j - 1), l)
617 end do
618
619 q_prim_vf(eqn_idx%mom%beg + 1)%sf(k, n + j, l) = -q_prim_vf(eqn_idx%mom%beg + 1)%sf(k, n - (j - 1), l)
620
621 do i = eqn_idx%mom%beg + 2, sys_size
622 q_prim_vf(i)%sf(k, n + j, l) = q_prim_vf(i)%sf(k, n - (j - 1), l)
623 end do
624
625 if (chemistry .and. present(q_t_sf)) then
626 q_t_sf%sf(k, n + j, l) = q_t_sf%sf(k, n - (j - 1), l)
627 end if
628
629 if (elasticity) then
630 do i = 1, shear_bc_flip_num
631 q_prim_vf(shear_bc_flip_indices(2, i))%sf(k, n + j, l) = -q_prim_vf(shear_bc_flip_indices(2, &
632 & i))%sf(k, n - (j - 1), l)
633 end do
634 end if
635
636 if (hyperelasticity) then
637 q_prim_vf(eqn_idx%xi%beg + 1)%sf(k, n + j, l) = -q_prim_vf(eqn_idx%xi%beg + 1)%sf(k, n - (j - 1), l)
638 end if
639 end do
640
641 if (qbmm .and. .not. polytropic .and. present(pb_in) .and. present(mv_in)) then
642 do i = 1, nb
643 do q = 1, nnode
644 do j = 1, buff_size
645 pb_in(k, n + j, l, q, i) = pb_in(k, n - (j - 1), l, q, i)
646 mv_in(k, n + j, l, q, i) = mv_in(k, n - (j - 1), l, q, i)
647 end do
648 end do
649 end do
650 end if
651 end if
652 else if (bc_dir == 3) then !< z-direction
653 if (bc_loc == -1) then !< bc_z%beg
654 do j = 1, buff_size
655 do i = 1, eqn_idx%mom%beg + 1
656 q_prim_vf(i)%sf(k, l, -j) = q_prim_vf(i)%sf(k, l, j - 1)
657 end do
658
659 q_prim_vf(eqn_idx%mom%end)%sf(k, l, -j) = -q_prim_vf(eqn_idx%mom%end)%sf(k, l, j - 1)
660
661 do i = eqn_idx%E, sys_size
662 q_prim_vf(i)%sf(k, l, -j) = q_prim_vf(i)%sf(k, l, j - 1)
663 end do
664
665 if (chemistry .and. present(q_t_sf)) then
666 q_t_sf%sf(k, l, -j) = q_t_sf%sf(k, l, j - 1)
667 end if
668
669 if (elasticity) then
670 do i = 1, shear_bc_flip_num
671 q_prim_vf(shear_bc_flip_indices(3, i))%sf(k, l, -j) = -q_prim_vf(shear_bc_flip_indices(3, i))%sf(k, &
672 & l, j - 1)
673 end do
674 end if
675
676 if (hyperelasticity) then
677 q_prim_vf(eqn_idx%xi%end)%sf(k, l, -j) = -q_prim_vf(eqn_idx%xi%end)%sf(k, l, j - 1)
678 end if
679 end do
680
681 if (qbmm .and. .not. polytropic .and. present(pb_in) .and. present(mv_in)) then
682 do i = 1, nb
683 do q = 1, nnode
684 do j = 1, buff_size
685 pb_in(k, l, -j, q, i) = pb_in(k, l, j - 1, q, i)
686 mv_in(k, l, -j, q, i) = mv_in(k, l, j - 1, q, i)
687 end do
688 end do
689 end do
690 end if
691 else !< bc_z%end
692 do j = 1, buff_size
693 do i = 1, eqn_idx%mom%beg + 1
694 q_prim_vf(i)%sf(k, l, p + j) = q_prim_vf(i)%sf(k, l, p - (j - 1))
695 end do
696
697 q_prim_vf(eqn_idx%mom%end)%sf(k, l, p + j) = -q_prim_vf(eqn_idx%mom%end)%sf(k, l, p - (j - 1))
698
699 do i = eqn_idx%E, sys_size
700 q_prim_vf(i)%sf(k, l, p + j) = q_prim_vf(i)%sf(k, l, p - (j - 1))
701 end do
702
703 if (chemistry .and. present(q_t_sf)) then
704 q_t_sf%sf(k, l, p + j) = q_t_sf%sf(k, l, p - (j - 1))
705 end if
706
707 if (elasticity) then
708 do i = 1, shear_bc_flip_num
709 q_prim_vf(shear_bc_flip_indices(3, i))%sf(k, l, p + j) = -q_prim_vf(shear_bc_flip_indices(3, &
710 & i))%sf(k, l, p - (j - 1))
711 end do
712 end if
713
714 if (hyperelasticity) then
715 q_prim_vf(eqn_idx%xi%end)%sf(k, l, p + j) = -q_prim_vf(eqn_idx%xi%end)%sf(k, l, p - (j - 1))
716 end if
717 end do
718
719 if (qbmm .and. .not. polytropic .and. present(pb_in) .and. present(mv_in)) then
720 do i = 1, nb
721 do q = 1, nnode
722 do j = 1, buff_size
723 pb_in(k, l, p + j, q, i) = pb_in(k, l, p - (j - 1), q, i)
724 mv_in(k, l, p + j, q, i) = mv_in(k, l, p - (j - 1), q, i)
725 end do
726 end do
727 end do
728 end if
729 end if
730 end if
731
732 end subroutine s_symmetry
733
734 !> Apply periodic boundary conditions by copying values from the opposite domain boundary.
735 subroutine s_periodic(q_prim_vf, bc_dir, bc_loc, k, l, pb_in, mv_in, q_T_sf)
736
737
738# 360 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
739#if MFC_OpenACC
740# 360 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
741!$acc routine seq
742# 360 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
743#elif MFC_OpenMP
744# 360 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
745
746# 360 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
747
748# 360 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
749!$omp declare target device_type(any)
750# 360 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
751#endif
752 type(scalar_field), dimension(sys_size), intent(inout) :: q_prim_vf
753 real(stp), optional, dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:,1:), intent(inout) :: pb_in, mv_in
754 integer, intent(in) :: bc_dir, bc_loc
755 integer, intent(in) :: k, l
756 integer :: j, q, i
757 type(scalar_field), optional, intent(inout) :: q_T_sf
758
759 if (bc_dir == 1) then !< x-direction
760 if (bc_loc == -1) then !< bc_x%beg
761 do i = 1, sys_size
762 do j = 1, buff_size
763 q_prim_vf(i)%sf(-j, k, l) = q_prim_vf(i)%sf(m - (j - 1), k, l)
764 end do
765 end do
766
767 if (chemistry .and. present(q_t_sf)) then
768 do j = 1, buff_size
769 q_t_sf%sf(-j, k, l) = q_t_sf%sf(m - (j - 1), k, l)
770 end do
771 end if
772
773 if (qbmm .and. .not. polytropic .and. present(pb_in) .and. present(mv_in)) then
774 do i = 1, nb
775 do q = 1, nnode
776 do j = 1, buff_size
777 pb_in(-j, k, l, q, i) = pb_in(m - (j - 1), k, l, q, i)
778 mv_in(-j, k, l, q, i) = mv_in(m - (j - 1), k, l, q, i)
779 end do
780 end do
781 end do
782 end if
783 else !< bc_x%end
784 do i = 1, sys_size
785 do j = 1, buff_size
786 q_prim_vf(i)%sf(m + j, k, l) = q_prim_vf(i)%sf(j - 1, k, l)
787 end do
788 end do
789
790 if (chemistry .and. present(q_t_sf)) then
791 do j = 1, buff_size
792 q_t_sf%sf(m + j, k, l) = q_t_sf%sf(j - 1, k, l)
793 end do
794 end if
795
796 if (qbmm .and. .not. polytropic .and. present(pb_in) .and. present(mv_in)) then
797 do i = 1, nb
798 do q = 1, nnode
799 do j = 1, buff_size
800 pb_in(m + j, k, l, q, i) = pb_in(j - 1, k, l, q, i)
801 mv_in(m + j, k, l, q, i) = mv_in(j - 1, k, l, q, i)
802 end do
803 end do
804 end do
805 end if
806 end if
807 else if (bc_dir == 2) then !< y-direction
808 if (bc_loc == -1) then !< bc_y%beg
809 do i = 1, sys_size
810 do j = 1, buff_size
811 q_prim_vf(i)%sf(k, -j, l) = q_prim_vf(i)%sf(k, n - (j - 1), l)
812 end do
813 end do
814
815 if (chemistry .and. present(q_t_sf)) then
816 do j = 1, buff_size
817 q_t_sf%sf(k, -j, l) = q_t_sf%sf(k, n - (j - 1), l)
818 end do
819 end if
820
821 if (qbmm .and. .not. polytropic .and. present(pb_in) .and. present(mv_in)) then
822 do i = 1, nb
823 do q = 1, nnode
824 do j = 1, buff_size
825 pb_in(k, -j, l, q, i) = pb_in(k, n - (j - 1), l, q, i)
826 mv_in(k, -j, l, q, i) = mv_in(k, n - (j - 1), l, q, i)
827 end do
828 end do
829 end do
830 end if
831 else !< bc_y%end
832 do i = 1, sys_size
833 do j = 1, buff_size
834 q_prim_vf(i)%sf(k, n + j, l) = q_prim_vf(i)%sf(k, j - 1, l)
835 end do
836 end do
837
838 if (chemistry .and. present(q_t_sf)) then
839 do j = 1, buff_size
840 q_t_sf%sf(k, n + j, l) = q_t_sf%sf(k, j - 1, l)
841 end do
842 end if
843
844 if (qbmm .and. .not. polytropic .and. present(pb_in) .and. present(mv_in)) then
845 do i = 1, nb
846 do q = 1, nnode
847 do j = 1, buff_size
848 pb_in(k, n + j, l, q, i) = pb_in(k, (j - 1), l, q, i)
849 mv_in(k, n + j, l, q, i) = mv_in(k, (j - 1), l, q, i)
850 end do
851 end do
852 end do
853 end if
854 end if
855 else if (bc_dir == 3) then !< z-direction
856 if (bc_loc == -1) then !< bc_z%beg
857 do i = 1, sys_size
858 do j = 1, buff_size
859 q_prim_vf(i)%sf(k, l, -j) = q_prim_vf(i)%sf(k, l, p - (j - 1))
860 end do
861 end do
862
863 if (chemistry .and. present(q_t_sf)) then
864 do j = 1, buff_size
865 q_t_sf%sf(k, l, -j) = q_t_sf%sf(k, l, p - (j - 1))
866 end do
867 end if
868
869 if (qbmm .and. .not. polytropic .and. present(pb_in) .and. present(mv_in)) then
870 do i = 1, nb
871 do q = 1, nnode
872 do j = 1, buff_size
873 pb_in(k, l, -j, q, i) = pb_in(k, l, p - (j - 1), q, i)
874 mv_in(k, l, -j, q, i) = mv_in(k, l, p - (j - 1), q, i)
875 end do
876 end do
877 end do
878 end if
879 else !< bc_z%end
880 do i = 1, sys_size
881 do j = 1, buff_size
882 q_prim_vf(i)%sf(k, l, p + j) = q_prim_vf(i)%sf(k, l, j - 1)
883 end do
884 end do
885
886 if (chemistry .and. present(q_t_sf)) then
887 do j = 1, buff_size
888 q_t_sf%sf(k, l, p + j) = q_t_sf%sf(k, l, j - 1)
889 end do
890 end if
891
892 if (qbmm .and. .not. polytropic .and. present(pb_in) .and. present(mv_in)) then
893 do i = 1, nb
894 do q = 1, nnode
895 do j = 1, buff_size
896 pb_in(k, l, p + j, q, i) = pb_in(k, l, j - 1, q, i)
897 mv_in(k, l, p + j, q, i) = mv_in(k, l, j - 1, q, i)
898 end do
899 end do
900 end do
901 end if
902 end if
903 end if
904
905 end subroutine s_periodic
906
907 !> Apply axis boundary conditions for cylindrical coordinates by reflecting values across the axis with azimuthal phase shift.
908 subroutine s_axis(q_prim_vf, pb_in, mv_in, k, l)
909
910
911# 519 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
912#if MFC_OpenACC
913# 519 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
914!$acc routine seq
915# 519 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
916#elif MFC_OpenMP
917# 519 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
918
919# 519 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
920
921# 519 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
922!$omp declare target device_type(any)
923# 519 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
924#endif
925 type(scalar_field), dimension(sys_size), intent(inout) :: q_prim_vf
926 real(stp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:,1:), optional, intent(inout) :: pb_in, mv_in
927 integer, intent(in) :: k, l
928 integer :: j, q, i
929
930 do j = 1, buff_size
931 if (z_cc(l) < pi) then
932 do i = 1, eqn_idx%mom%beg
933 q_prim_vf(i)%sf(k, -j, l) = q_prim_vf(i)%sf(k, j - 1, l + ((p + 1)/2))
934 end do
935
936 q_prim_vf(eqn_idx%mom%beg + 1)%sf(k, -j, l) = -q_prim_vf(eqn_idx%mom%beg + 1)%sf(k, j - 1, l + ((p + 1)/2))
937
938 q_prim_vf(eqn_idx%mom%end)%sf(k, -j, l) = -q_prim_vf(eqn_idx%mom%end)%sf(k, j - 1, l + ((p + 1)/2))
939
940 do i = eqn_idx%E, sys_size
941 q_prim_vf(i)%sf(k, -j, l) = q_prim_vf(i)%sf(k, j - 1, l + ((p + 1)/2))
942 end do
943 else
944 do i = 1, eqn_idx%mom%beg
945 q_prim_vf(i)%sf(k, -j, l) = q_prim_vf(i)%sf(k, j - 1, l - ((p + 1)/2))
946 end do
947
948 q_prim_vf(eqn_idx%mom%beg + 1)%sf(k, -j, l) = -q_prim_vf(eqn_idx%mom%beg + 1)%sf(k, j - 1, l - ((p + 1)/2))
949
950 q_prim_vf(eqn_idx%mom%end)%sf(k, -j, l) = -q_prim_vf(eqn_idx%mom%end)%sf(k, j - 1, l - ((p + 1)/2))
951
952 do i = eqn_idx%E, sys_size
953 q_prim_vf(i)%sf(k, -j, l) = q_prim_vf(i)%sf(k, j - 1, l - ((p + 1)/2))
954 end do
955 end if
956 end do
957
958 if (qbmm .and. .not. polytropic .and. present(pb_in) .and. present(mv_in)) then
959 do i = 1, nb
960 do q = 1, nnode
961 do j = 1, buff_size
962 if (z_cc(l) < pi) then
963 pb_in(k, -j, l, q, i) = pb_in(k, j - 1, l + ((p + 1)/2), q, i)
964 mv_in(k, -j, l, q, i) = mv_in(k, j - 1, l + ((p + 1)/2), q, i)
965 else
966 pb_in(k, -j, l, q, i) = pb_in(k, j - 1, l - ((p + 1)/2), q, i)
967 mv_in(k, -j, l, q, i) = mv_in(k, j - 1, l - ((p + 1)/2), q, i)
968 end if
969 end do
970 end do
971 end do
972 end if
973
974 end subroutine s_axis
975
976 !> Apply slip wall boundary conditions by extrapolating scalars and reflecting the wall-normal velocity component.
977 subroutine s_slip_wall(q_prim_vf, bc_dir, bc_loc, k, l, q_T_sf)
978
979
980# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
981#ifdef _CRAYFTN
982# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
983#if MFC_OpenACC
984# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
985!$acc routine seq
986# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
987#elif MFC_OpenMP
988# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
989
990# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
991
992# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
993!$omp declare target device_type(any)
994# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
995#else
996# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
997!DIR$ INLINEALWAYS s_slip_wall
998# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
999#endif
1000# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1001#elif MFC_OpenACC
1002# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1003!$acc routine seq
1004# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1005#elif MFC_OpenMP
1006# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1007
1008# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1009
1010# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1011!$omp declare target device_type(any)
1012# 574 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1013#endif
1014 type(scalar_field), dimension(sys_size), intent(inout) :: q_prim_vf
1015 integer, intent(in) :: bc_dir, bc_loc
1016 integer, intent(in) :: k, l
1017 integer :: j, i
1018 type(scalar_field), optional, intent(inout) :: q_T_sf
1019
1020 if (bc_dir == 1) then !< x-direction
1021 if (bc_loc == -1) then !< bc_x%beg
1022 do i = 1, sys_size
1023 do j = 1, buff_size
1024 if (i == eqn_idx%mom%beg) then
1025 q_prim_vf(i)%sf(-j, k, l) = -q_prim_vf(i)%sf(j - 1, k, l) + 2._wp*bc_x%vb1
1026 else
1027 q_prim_vf(i)%sf(-j, k, l) = q_prim_vf(i)%sf(0, k, l)
1028 end if
1029 end do
1030 end do
1031
1032 if (chemistry .and. present(q_t_sf)) then
1033 if (bc_x%isothermal_in) then
1034 do j = 1, buff_size
1035 q_t_sf%sf(-j, k, l) = 2._wp*bc_x%Twall_in - q_t_sf%sf(j - 1, k, l)
1036 end do
1037 else
1038 do j = 1, buff_size
1039 q_t_sf%sf(-j, k, l) = q_t_sf%sf(0, k, l)
1040 end do
1041 end if
1042 end if
1043 else !< bc_x%end
1044 do i = 1, sys_size
1045 do j = 1, buff_size
1046 if (i == eqn_idx%mom%beg) then
1047 q_prim_vf(i)%sf(m + j, k, l) = -q_prim_vf(i)%sf(m - (j - 1), k, l) + 2._wp*bc_x%ve1
1048 else
1049 q_prim_vf(i)%sf(m + j, k, l) = q_prim_vf(i)%sf(m, k, l)
1050 end if
1051 end do
1052 end do
1053
1054 if (chemistry .and. present(q_t_sf)) then
1055 if (bc_x%isothermal_out) then
1056 do j = 1, buff_size
1057 q_t_sf%sf(m + j, k, l) = 2._wp*bc_x%Twall_out - q_t_sf%sf(m - (j - 1), k, l)
1058 end do
1059 else
1060 do j = 1, buff_size
1061 q_t_sf%sf(m + j, k, l) = q_t_sf%sf(m, k, l)
1062 end do
1063 end if
1064 end if
1065 end if
1066 else if (bc_dir == 2) then !< y-direction
1067 if (bc_loc == -1) then !< bc_y%beg
1068 do i = 1, sys_size
1069 do j = 1, buff_size
1070 if (i == eqn_idx%mom%beg + 1) then
1071 q_prim_vf(i)%sf(k, -j, l) = -q_prim_vf(i)%sf(k, j - 1, l) + 2._wp*bc_y%vb2
1072 else
1073 q_prim_vf(i)%sf(k, -j, l) = q_prim_vf(i)%sf(k, 0, l)
1074 end if
1075 end do
1076 end do
1077
1078 if (chemistry .and. present(q_t_sf)) then
1079 if (bc_y%isothermal_in) then
1080 do j = 1, buff_size
1081 q_t_sf%sf(k, -j, l) = 2._wp*bc_y%Twall_in - q_t_sf%sf(k, j - 1, l)
1082 end do
1083 else
1084 do j = 1, buff_size
1085 q_t_sf%sf(k, -j, l) = q_t_sf%sf(k, 0, l)
1086 end do
1087 end if
1088 end if
1089 else !< bc_y%end
1090 do i = 1, sys_size
1091 do j = 1, buff_size
1092 if (i == eqn_idx%mom%beg + 1) then
1093 q_prim_vf(i)%sf(k, n + j, l) = -q_prim_vf(i)%sf(k, n - (j - 1), l) + 2._wp*bc_y%ve2
1094 else
1095 q_prim_vf(i)%sf(k, n + j, l) = q_prim_vf(i)%sf(k, n, l)
1096 end if
1097 end do
1098 end do
1099
1100 if (chemistry .and. present(q_t_sf)) then
1101 if (bc_y%isothermal_out) then
1102 do j = 1, buff_size
1103 q_t_sf%sf(k, n + j, l) = 2._wp*bc_y%Twall_out - q_t_sf%sf(k, n - (j - 1), l)
1104 end do
1105 else
1106 do j = 1, buff_size
1107 q_t_sf%sf(k, n + j, l) = q_t_sf%sf(k, n, l)
1108 end do
1109 end if
1110 end if
1111 end if
1112 else if (bc_dir == 3) then !< z-direction
1113 if (bc_loc == -1) then !< bc_z%beg
1114 do i = 1, sys_size
1115 do j = 1, buff_size
1116 if (i == eqn_idx%mom%end) then
1117 q_prim_vf(i)%sf(k, l, -j) = -q_prim_vf(i)%sf(k, l, j - 1) + 2._wp*bc_z%vb3
1118 else
1119 q_prim_vf(i)%sf(k, l, -j) = q_prim_vf(i)%sf(k, l, 0)
1120 end if
1121 end do
1122 end do
1123
1124 if (chemistry .and. present(q_t_sf)) then
1125 if (bc_z%isothermal_in) then
1126 do j = 1, buff_size
1127 q_t_sf%sf(k, l, -j) = 2._wp*bc_z%Twall_in - q_t_sf%sf(k, l, j - 1)
1128 end do
1129 else
1130 do j = 1, buff_size
1131 q_t_sf%sf(k, l, -j) = q_t_sf%sf(k, l, 0)
1132 end do
1133 end if
1134 end if
1135 else !< bc_z%end
1136 do i = 1, sys_size
1137 do j = 1, buff_size
1138 if (i == eqn_idx%mom%end) then
1139 q_prim_vf(i)%sf(k, l, p + j) = -q_prim_vf(i)%sf(k, l, p - (j - 1)) + 2._wp*bc_z%ve3
1140 else
1141 q_prim_vf(i)%sf(k, l, p + j) = q_prim_vf(i)%sf(k, l, p)
1142 end if
1143 end do
1144 end do
1145
1146 if (chemistry .and. present(q_t_sf)) then
1147 if (bc_z%isothermal_out) then
1148 do j = 1, buff_size
1149 q_t_sf%sf(k, l, p + j) = 2._wp*bc_z%Twall_out - q_t_sf%sf(k, l, p - (j - 1))
1150 end do
1151 else
1152 do j = 1, buff_size
1153 q_t_sf%sf(k, l, p + j) = q_t_sf%sf(k, l, p)
1154 end do
1155 end if
1156 end if
1157 end if
1158 end if
1159
1160 end subroutine s_slip_wall
1161
1162 !> Apply no-slip wall boundary conditions by reflecting and negating all velocity components at the wall.
1163 subroutine s_no_slip_wall(q_prim_vf, bc_dir, bc_loc, k, l, q_T_sf)
1164
1165
1166# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1167#ifdef _CRAYFTN
1168# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1169#if MFC_OpenACC
1170# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1171!$acc routine seq
1172# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1173#elif MFC_OpenMP
1174# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1175
1176# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1177
1178# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1179!$omp declare target device_type(any)
1180# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1181#else
1182# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1183!DIR$ INLINEALWAYS s_no_slip_wall
1184# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1185#endif
1186# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1187#elif MFC_OpenACC
1188# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1189!$acc routine seq
1190# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1191#elif MFC_OpenMP
1192# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1193
1194# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1195
1196# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1197!$omp declare target device_type(any)
1198# 726 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1199#endif
1200
1201 type(scalar_field), dimension(sys_size), intent(inout) :: q_prim_vf
1202 integer, intent(in) :: bc_dir, bc_loc
1203 integer, intent(in) :: k, l
1204 integer :: j, i
1205 type(scalar_field), optional, intent(inout) :: q_T_sf
1206
1207 if (bc_dir == 1) then !< x-direction
1208 if (bc_loc == -1) then !< bc_x%beg
1209 do i = 1, sys_size
1210 do j = 1, buff_size
1211 if (i == eqn_idx%mom%beg) then
1212 q_prim_vf(i)%sf(-j, k, l) = -q_prim_vf(i)%sf(j - 1, k, l) + 2._wp*bc_x%vb1
1213 else if (i == eqn_idx%mom%beg + 1 .and. num_dims > 1) then
1214 q_prim_vf(i)%sf(-j, k, l) = -q_prim_vf(i)%sf(j - 1, k, l) + 2._wp*bc_x%vb2
1215 else if (i == eqn_idx%mom%beg + 2 .and. num_dims > 2) then
1216 q_prim_vf(i)%sf(-j, k, l) = -q_prim_vf(i)%sf(j - 1, k, l) + 2._wp*bc_x%vb3
1217 else
1218 q_prim_vf(i)%sf(-j, k, l) = q_prim_vf(i)%sf(0, k, l)
1219 end if
1220 end do
1221 end do
1222
1223 if (chemistry .and. present(q_t_sf)) then
1224 if (bc_x%isothermal_in) then
1225 do j = 1, buff_size
1226 q_t_sf%sf(-j, k, l) = 2._wp*bc_x%Twall_in - q_t_sf%sf(j - 1, k, l)
1227 end do
1228 else
1229 do j = 1, buff_size
1230 q_t_sf%sf(-j, k, l) = q_t_sf%sf(0, k, l)
1231 end do
1232 end if
1233 end if
1234 else !< bc_x%end
1235 do i = 1, sys_size
1236 do j = 1, buff_size
1237 if (i == eqn_idx%mom%beg) then
1238 q_prim_vf(i)%sf(m + j, k, l) = -q_prim_vf(i)%sf(m - (j - 1), k, l) + 2._wp*bc_x%ve1
1239 else if (i == eqn_idx%mom%beg + 1 .and. num_dims > 1) then
1240 q_prim_vf(i)%sf(m + j, k, l) = -q_prim_vf(i)%sf(m - (j - 1), k, l) + 2._wp*bc_x%ve2
1241 else if (i == eqn_idx%mom%beg + 2 .and. num_dims > 2) then
1242 q_prim_vf(i)%sf(m + j, k, l) = -q_prim_vf(i)%sf(m - (j - 1), k, l) + 2._wp*bc_x%ve3
1243 else
1244 q_prim_vf(i)%sf(m + j, k, l) = q_prim_vf(i)%sf(m, k, l)
1245 end if
1246 end do
1247 end do
1248
1249 if (chemistry .and. present(q_t_sf)) then
1250 if (bc_x%isothermal_out) then
1251 do j = 1, buff_size
1252 q_t_sf%sf(m + j, k, l) = 2._wp*bc_x%Twall_out - q_t_sf%sf(m - (j - 1), k, l)
1253 end do
1254 else
1255 do j = 1, buff_size
1256 q_t_sf%sf(m + j, k, l) = q_t_sf%sf(m, k, l)
1257 end do
1258 end if
1259 end if
1260 end if
1261 else if (bc_dir == 2) then !< y-direction
1262 if (bc_loc == -1) then !< bc_y%beg
1263 do i = 1, sys_size
1264 do j = 1, buff_size
1265 if (i == eqn_idx%mom%beg) then
1266 q_prim_vf(i)%sf(k, -j, l) = -q_prim_vf(i)%sf(k, j - 1, l) + 2._wp*bc_y%vb1
1267 else if (i == eqn_idx%mom%beg + 1 .and. num_dims > 1) then
1268 q_prim_vf(i)%sf(k, -j, l) = -q_prim_vf(i)%sf(k, j - 1, l) + 2._wp*bc_y%vb2
1269 else if (i == eqn_idx%mom%beg + 2 .and. num_dims > 2) then
1270 q_prim_vf(i)%sf(k, -j, l) = -q_prim_vf(i)%sf(k, j - 1, l) + 2._wp*bc_y%vb3
1271 else
1272 q_prim_vf(i)%sf(k, -j, l) = q_prim_vf(i)%sf(k, 0, l)
1273 end if
1274 end do
1275 end do
1276 if (chemistry .and. present(q_t_sf)) then
1277 if (bc_y%isothermal_in) then
1278 do j = 1, buff_size
1279 q_t_sf%sf(k, -j, l) = 2._wp*bc_y%Twall_in - q_t_sf%sf(k, j - 1, l)
1280 end do
1281 else
1282 do j = 1, buff_size
1283 q_t_sf%sf(k, -j, l) = q_t_sf%sf(k, 0, l)
1284 end do
1285 end if
1286 end if
1287 else !< bc_y%end
1288 do i = 1, sys_size
1289 do j = 1, buff_size
1290 if (i == eqn_idx%mom%beg) then
1291 q_prim_vf(i)%sf(k, n + j, l) = -q_prim_vf(i)%sf(k, n - (j - 1), l) + 2._wp*bc_y%ve1
1292 else if (i == eqn_idx%mom%beg + 1 .and. num_dims > 1) then
1293 q_prim_vf(i)%sf(k, n + j, l) = -q_prim_vf(i)%sf(k, n - (j - 1), l) + 2._wp*bc_y%ve2
1294 else if (i == eqn_idx%mom%beg + 2 .and. num_dims > 2) then
1295 q_prim_vf(i)%sf(k, n + j, l) = -q_prim_vf(i)%sf(k, n - (j - 1), l) + 2._wp*bc_y%ve3
1296 else
1297 q_prim_vf(i)%sf(k, n + j, l) = q_prim_vf(i)%sf(k, n, l)
1298 end if
1299 end do
1300 end do
1301 if (chemistry .and. present(q_t_sf)) then
1302 if (bc_y%isothermal_out) then
1303 do j = 1, buff_size
1304 q_t_sf%sf(k, n + j, l) = 2._wp*bc_y%Twall_out - q_t_sf%sf(k, n - (j - 1), l)
1305 end do
1306 else
1307 do j = 1, buff_size
1308 q_t_sf%sf(k, n + j, l) = q_t_sf%sf(k, n, l)
1309 end do
1310 end if
1311 end if
1312 end if
1313 else if (bc_dir == 3) then !< z-direction
1314 if (bc_loc == -1) then !< bc_z%beg
1315 do i = 1, sys_size
1316 do j = 1, buff_size
1317 if (i == eqn_idx%mom%beg) then
1318 q_prim_vf(i)%sf(k, l, -j) = -q_prim_vf(i)%sf(k, l, j - 1) + 2._wp*bc_z%vb1
1319 else if (i == eqn_idx%mom%beg + 1 .and. num_dims > 1) then
1320 q_prim_vf(i)%sf(k, l, -j) = -q_prim_vf(i)%sf(k, l, j - 1) + 2._wp*bc_z%vb2
1321 else if (i == eqn_idx%mom%beg + 2 .and. num_dims > 2) then
1322 q_prim_vf(i)%sf(k, l, -j) = -q_prim_vf(i)%sf(k, l, j - 1) + 2._wp*bc_z%vb3
1323 else
1324 q_prim_vf(i)%sf(k, l, -j) = q_prim_vf(i)%sf(k, l, 0)
1325 end if
1326 end do
1327 end do
1328 if (chemistry .and. present(q_t_sf)) then
1329 if (bc_z%isothermal_in) then
1330 do j = 1, buff_size
1331 q_t_sf%sf(k, l, -j) = 2._wp*bc_z%Twall_in - q_t_sf%sf(k, l, j - 1)
1332 end do
1333 else
1334 do j = 1, buff_size
1335 q_t_sf%sf(k, l, -j) = q_t_sf%sf(k, l, 0)
1336 end do
1337 end if
1338 end if
1339 else !< bc_z%end
1340 do i = 1, sys_size
1341 do j = 1, buff_size
1342 if (i == eqn_idx%mom%beg) then
1343 q_prim_vf(i)%sf(k, l, p + j) = -q_prim_vf(i)%sf(k, l, p - (j - 1)) + 2._wp*bc_z%ve1
1344 else if (i == eqn_idx%mom%beg + 1 .and. num_dims > 1) then
1345 q_prim_vf(i)%sf(k, l, p + j) = -q_prim_vf(i)%sf(k, l, p - (j - 1)) + 2._wp*bc_z%ve2
1346 else if (i == eqn_idx%mom%beg + 2 .and. num_dims > 2) then
1347 q_prim_vf(i)%sf(k, l, p + j) = -q_prim_vf(i)%sf(k, l, p - (j - 1)) + 2._wp*bc_z%ve3
1348 else
1349 q_prim_vf(i)%sf(k, l, p + j) = q_prim_vf(i)%sf(k, l, p)
1350 end if
1351 end do
1352 end do
1353 if (chemistry .and. present(q_t_sf)) then
1354 if (bc_z%isothermal_out) then
1355 do j = 1, buff_size
1356 q_t_sf%sf(k, l, p + j) = 2._wp*bc_z%Twall_out - q_t_sf%sf(k, l, p - (j - 1))
1357 end do
1358 else
1359 do j = 1, buff_size
1360 q_t_sf%sf(k, l, p + j) = q_t_sf%sf(k, l, p)
1361 end do
1362 end if
1363 end if
1364 end if
1365 end if
1366
1367 end subroutine s_no_slip_wall
1368
1369 !> Apply Dirichlet boundary conditions by prescribing ghost cell values from stored boundary buffers.
1370 subroutine s_dirichlet(q_prim_vf, bc_dir, bc_loc, k, l, q_T_sf)
1371
1372
1373# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1374#ifdef _CRAYFTN
1375# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1376#if MFC_OpenACC
1377# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1378!$acc routine seq
1379# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1380#elif MFC_OpenMP
1381# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1382
1383# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1384
1385# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1386!$omp declare target device_type(any)
1387# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1388#else
1389# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1390!DIR$ INLINEALWAYS s_dirichlet
1391# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1392#endif
1393# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1394#elif MFC_OpenACC
1395# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1396!$acc routine seq
1397# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1398#elif MFC_OpenMP
1399# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1400
1401# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1402
1403# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1404!$omp declare target device_type(any)
1405# 899 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1406#endif
1407 type(scalar_field), dimension(sys_size), intent(inout) :: q_prim_vf
1408 integer, intent(in) :: bc_dir, bc_loc
1409 integer, intent(in) :: k, l
1410 integer :: j, i
1411 type(scalar_field), optional, intent(inout) :: q_T_sf
1412
1413#ifdef MFC_SIMULATION
1414 if (bc_dir == 1) then !< x-direction
1415 if (bc_loc == -1) then ! bc_x%beg
1416 do i = 1, sys_size
1417 do j = 1, buff_size
1418 q_prim_vf(i)%sf(-j, k, l) = bc_buffers(1, 1)%sf(i, k, l)
1419 end do
1420 end do
1421 if (chemistry .and. present(q_t_sf)) then
1422 do j = 1, buff_size
1423 q_t_sf%sf(-j, k, l) = bc_buffers(1, 1)%sf(sys_size + 1, k, l)
1424 end do
1425 end if
1426 else !< bc_x%end
1427 do i = 1, sys_size
1428 do j = 1, buff_size
1429 q_prim_vf(i)%sf(m + j, k, l) = bc_buffers(1, 2)%sf(i, k, l)
1430 end do
1431 end do
1432 if (chemistry .and. present(q_t_sf)) then
1433 do j = 1, buff_size
1434 q_t_sf%sf(m + j, k, l) = bc_buffers(1, 2)%sf(sys_size + 1, k, l)
1435 end do
1436 end if
1437 end if
1438 else if (bc_dir == 2) then !< y-direction
1439# 933 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1440 if (bc_loc == -1) then !< bc_y%beg
1441 do i = 1, sys_size
1442 do j = 1, buff_size
1443 q_prim_vf(i)%sf(k, -j, l) = bc_buffers(2, 1)%sf(k, i, l)
1444 end do
1445 end do
1446 if (chemistry .and. present(q_t_sf)) then
1447 do j = 1, buff_size
1448 q_t_sf%sf(k, -j, l) = bc_buffers(2, 1)%sf(k, sys_size + 1, l)
1449 end do
1450 end if
1451 else !< bc_y%end
1452 do i = 1, sys_size
1453 do j = 1, buff_size
1454 q_prim_vf(i)%sf(k, n + j, l) = bc_buffers(2, 2)%sf(k, i, l)
1455 end do
1456 end do
1457 if (chemistry .and. present(q_t_sf)) then
1458 do j = 1, buff_size
1459 q_t_sf%sf(k, n + j, l) = bc_buffers(2, 2)%sf(k, sys_size + 1, l)
1460 end do
1461 end if
1462 end if
1463# 957 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1464 else if (bc_dir == 3) then !< z-direction
1465# 959 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1466 if (bc_loc == -1) then !< bc_z%beg
1467 do i = 1, sys_size
1468 do j = 1, buff_size
1469 q_prim_vf(i)%sf(k, l, -j) = bc_buffers(3, 1)%sf(k, l, i)
1470 end do
1471 end do
1472 if (chemistry .and. present(q_t_sf)) then
1473 do j = 1, buff_size
1474 q_t_sf%sf(k, l, -j) = bc_buffers(3, 1)%sf(k, l, sys_size + 1)
1475 end do
1476 end if
1477 else !< bc_z%end
1478 do i = 1, sys_size
1479 do j = 1, buff_size
1480 q_prim_vf(i)%sf(k, l, p + j) = bc_buffers(3, 2)%sf(k, l, i)
1481 end do
1482 end do
1483 if (chemistry .and. present(q_t_sf)) then
1484 do j = 1, buff_size
1485 q_t_sf%sf(k, l, p + j) = bc_buffers(3, 2)%sf(k, l, sys_size + 1)
1486 end do
1487 end if
1488 end if
1489# 983 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1490 end if
1491#else
1492 call s_ghost_cell_extrapolation(q_prim_vf, bc_dir, bc_loc, k, l, q_t_sf)
1493#endif
1494
1495 end subroutine s_dirichlet
1496
1497 !> Extrapolate QBMM bubble pressure and mass-vapor variables into ghost cells by copying boundary values.
1498 subroutine s_qbmm_extrapolation(bc_dir, bc_loc, k, l, pb_in, mv_in)
1499
1500
1501# 993 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1502#if MFC_OpenACC
1503# 993 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1504!$acc routine seq
1505# 993 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1506#elif MFC_OpenMP
1507# 993 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1508
1509# 993 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1510
1511# 993 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1512!$omp declare target device_type(any)
1513# 993 "/home/runner/work/MFC/MFC/src/common/m_boundary_primitives.fpp"
1514#endif
1515 real(stp), optional, dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:,1:), intent(inout) :: pb_in, mv_in
1516 integer, intent(in) :: bc_dir, bc_loc
1517 integer, intent(in) :: k, l
1518 integer :: j, q, i
1519
1520 if (bc_dir == 1) then !< x-direction
1521 if (bc_loc == -1) then ! bc_x%beg
1522 do i = 1, nb
1523 do q = 1, nnode
1524 do j = 1, buff_size
1525 pb_in(-j, k, l, q, i) = pb_in(0, k, l, q, i)
1526 mv_in(-j, k, l, q, i) = mv_in(0, k, l, q, i)
1527 end do
1528 end do
1529 end do
1530 else !< bc_x%end
1531 do i = 1, nb
1532 do q = 1, nnode
1533 do j = 1, buff_size
1534 pb_in(m + j, k, l, q, i) = pb_in(m, k, l, q, i)
1535 mv_in(m + j, k, l, q, i) = mv_in(m, k, l, q, i)
1536 end do
1537 end do
1538 end do
1539 end if
1540 else if (bc_dir == 2) then !< y-direction
1541 if (bc_loc == -1) then !< bc_y%beg
1542 do i = 1, nb
1543 do q = 1, nnode
1544 do j = 1, buff_size
1545 pb_in(k, -j, l, q, i) = pb_in(k, 0, l, q, i)
1546 mv_in(k, -j, l, q, i) = mv_in(k, 0, l, q, i)
1547 end do
1548 end do
1549 end do
1550 else !< bc_y%end
1551 do i = 1, nb
1552 do q = 1, nnode
1553 do j = 1, buff_size
1554 pb_in(k, n + j, l, q, i) = pb_in(k, n, l, q, i)
1555 mv_in(k, n + j, l, q, i) = mv_in(k, n, l, q, i)
1556 end do
1557 end do
1558 end do
1559 end if
1560 else if (bc_dir == 3) then !< z-direction
1561 if (bc_loc == -1) then !< bc_z%beg
1562 do i = 1, nb
1563 do q = 1, nnode
1564 do j = 1, buff_size
1565 pb_in(k, l, -j, q, i) = pb_in(k, l, 0, q, i)
1566 mv_in(k, l, -j, q, i) = mv_in(k, l, 0, q, i)
1567 end do
1568 end do
1569 end do
1570 else !< bc_z%end
1571 do i = 1, nb
1572 do q = 1, nnode
1573 do j = 1, buff_size
1574 pb_in(k, l, p + j, q, i) = pb_in(k, l, p, q, i)
1575 mv_in(k, l, p + j, q, i) = mv_in(k, l, p, q, i)
1576 end do
1577 end do
1578 end do
1579 end if
1580 end if
1581
1582 end subroutine s_qbmm_extrapolation
1583
1584end module m_boundary_primitives
Per-cell noncharacteristic boundary condition primitives applied in the ghost cells.
subroutine s_slip_wall(q_prim_vf, bc_dir, bc_loc, k, l, q_t_sf)
Apply slip wall boundary conditions by extrapolating scalars and reflecting the wall-normal velocity ...
type(scalar_field), dimension(:,:), allocatable bc_buffers
subroutine s_no_slip_wall(q_prim_vf, bc_dir, bc_loc, k, l, q_t_sf)
Apply no-slip wall boundary conditions by reflecting and negating all velocity components at the wall...
subroutine s_axis(q_prim_vf, pb_in, mv_in, k, l)
Apply axis boundary conditions for cylindrical coordinates by reflecting values across the axis with ...
subroutine s_qbmm_extrapolation(bc_dir, bc_loc, k, l, pb_in, mv_in)
Extrapolate QBMM bubble pressure and mass-vapor variables into ghost cells by copying boundary values...
subroutine s_periodic(q_prim_vf, bc_dir, bc_loc, k, l, pb_in, mv_in, q_t_sf)
Apply periodic boundary conditions by copying values from the opposite domain boundary.
subroutine s_ghost_cell_extrapolation(q_prim_vf, bc_dir, bc_loc, k, l, q_t_sf)
Fill ghost cells by copying the nearest boundary cell value along the specified direction.
subroutine s_symmetry(q_prim_vf, bc_dir, bc_loc, k, l, pb_in, mv_in, q_t_sf)
Apply reflective (symmetry) boundary conditions by mirroring primitive variables and flipping the nor...
subroutine s_dirichlet(q_prim_vf, bc_dir, bc_loc, k, l, q_t_sf)
Apply Dirichlet boundary conditions by prescribing ghost cell values from stored boundary buffers.
Compile-time constant parameters: default values, tolerances, and physical constants.
integer, parameter nnode
Number of QBMM nodes.
real(wp), parameter pi
Pi.
Shared derived types for field data, patch geometry, bubble dynamics, and MPI I/O structures.
Global parameters for the post-process: domain geometry, equation of state, and output database setti...
integer buff_size
Number of ghost cells for boundary condition storage.
real(wp), dimension(:), allocatable z_cc
Derived type annexing a scalar field (SF).