MFC
Exascale flow solver
Loading...
Searching...
No Matches
m_compute_levelset.fpp.f90
Go to the documentation of this file.
1# 1 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
2!>
3!!@file
4!! @brief Contains module m_compute_levelset
5
6# 1 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp" 1
7# 1 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp" 1
8# 1 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp" 1
9# 2 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
10# 3 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
11# 4 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
12# 5 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
13# 6 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
14
15# 8 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
16# 9 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
17# 10 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
18
19# 17 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
20
21# 46 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
22
23# 58 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
24
25# 68 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
26
27# 98 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
28
29# 110 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
30
31# 120 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
32! New line at end of file is required for FYPP
33# 2 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp" 2
34# 1 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp" 1
35# 1 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp" 1
36# 2 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
37# 3 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
38# 4 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
39# 5 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
40# 6 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
41
42# 8 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
43# 9 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
44# 10 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
45
46# 17 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
47
48# 46 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
49
50# 58 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
51
52# 68 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
53
54# 98 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
55
56# 110 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
57
58# 120 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
59! New line at end of file is required for FYPP
60# 2 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp" 2
61
62# 4 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
63# 5 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
64# 6 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
65# 7 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
66# 8 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
67
68# 20 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
69
70# 43 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
71
72# 48 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
73
74# 53 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
75
76# 58 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
77
78# 63 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
79
80# 68 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
81
82# 76 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
83
84# 81 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
85
86# 86 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
87
88# 91 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
89
90# 96 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
91
92# 101 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
93
94# 106 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
95
96# 111 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
97
98# 116 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
99
100# 121 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
101
102# 151 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
103
104# 192 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
105
106# 207 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
107
108# 232 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
109
110# 243 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
111
112# 245 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
113# 255 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
114
115# 283 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
116
117# 293 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
118
119# 303 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
120
121# 312 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
122
123# 329 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
124
125# 339 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
126
127# 346 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
128
129# 352 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
130
131# 358 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
132
133# 364 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
134
135# 370 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
136
137# 376 "/home/runner/work/MFC/MFC/src/common/include/omp_macros.fpp"
138! New line at end of file is required for FYPP
139# 3 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp" 2
140# 1 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp" 1
141# 1 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp" 1
142# 2 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
143# 3 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
144# 4 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
145# 5 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
146# 6 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
147
148# 8 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
149# 9 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
150# 10 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
151
152# 17 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
153
154# 46 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
155
156# 58 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
157
158# 68 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
159
160# 98 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
161
162# 110 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
163
164# 120 "/home/runner/work/MFC/MFC/src/common/include/shared_parallel_macros.fpp"
165! New line at end of file is required for FYPP
166# 2 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp" 2
167
168# 7 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
169
170# 17 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
171
172# 22 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
173
174# 27 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
175
176# 32 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
177
178# 37 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
179
180# 42 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
181
182# 47 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
183
184# 52 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
185
186# 57 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
187
188# 62 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
189
190# 73 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
191
192# 78 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
193
194# 83 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
195
196# 88 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
197
198# 103 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
199
200# 131 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
201
202# 160 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
203
204# 175 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
205
206# 192 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
207
208# 213 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
209
210# 241 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
211
212# 256 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
213
214# 266 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
215
216# 275 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
217
218# 291 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
219
220# 301 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
221
222# 308 "/home/runner/work/MFC/MFC/src/common/include/acc_macros.fpp"
223! New line at end of file is required for FYPP
224# 4 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp" 2
225
226# 21 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
227
228# 37 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
229
230# 50 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
231
232# 76 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
233
234# 91 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
235
236# 102 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
237
238# 115 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
239
240# 143 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
241
242# 154 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
243
244# 165 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
245
246# 176 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
247
248# 187 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
249
250# 198 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
251
252# 208 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
253
254# 214 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
255
256# 220 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
257
258# 226 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
259
260# 232 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
261
262# 234 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
263# 235 "/home/runner/work/MFC/MFC/src/common/include/parallel_macros.fpp"
264! New line at end of file is required for FYPP
265# 2 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp" 2
266
267# 14 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
268
269! Caution:
270! This macro requires the use of a binding script to set CUDA_VISIBLE_DEVICES, such that we have one GPU device per MPI rank.
271! That's because for both cudaMemAdvise (preferred location) and cudaMemPrefetchAsync we use location = device_id = 0.
272! For an example see misc/nvidia_uvm/bind.sh.
273# 63 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
274
275# 81 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
276
277# 88 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
278
279# 111 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
280
281# 127 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
282
283# 153 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
284
285# 159 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
286
287# 167 "/home/runner/work/MFC/MFC/src/common/include/macros.fpp"
288! New line at end of file is required for FYPP
289# 6 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp" 2
290
291!> @brief Computes signed-distance level-set fields and surface normals for immersed-boundary patch geometries
293
294 use m_ib_patches !< the ib patch parameters
295
296 use m_model !< subroutine(s) related to stl files
297
298 use m_derived_types !< definitions of the derived types
299
300 use m_global_parameters !< definitions of the global parameters
301
302 use m_mpi_proxy !< message passing interface (mpi) module proxy
303
304 use m_helper_basic !< functions to compare floating point numbers
305
306 implicit none
307
308 private; public :: s_apply_levelset
309
310contains
311
312 !> @brief Dispatches level-set distance and normal computations for all ghost points based on their patch geometry type.
313 impure subroutine s_apply_levelset(gps, num_gps)
314
315 type(ghost_point), dimension(:), intent(inout) :: gps
316 integer, intent(in) :: num_gps
317
318 integer :: i, patch_id, patch_geometry
319
320 if (num_gps < 1) then
321 return
322 end if
323
324
325# 40 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
326#if defined(MFC_OpenACC)
327# 40 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
328!$acc update device(gps(1:num_gps))
329# 40 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
330#elif defined(MFC_OpenMP)
331# 40 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
332!$omp target update to(gps(1:num_gps))
333# 40 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
334#endif
335
336 ! 3D Patch Geometries
337 if (p > 0) then
338
339
340# 45 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
341
342# 45 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
343#if defined(MFC_OpenACC)
344# 45 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
345!$acc parallel loop gang vector default(present) private(i, patch_id, patch_geometry) copy(gps) copyin(patch_ib(1:num_ibs), Np)
346# 45 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
347#elif defined(MFC_OpenMP)
348# 45 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
349
350# 45 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
351
352# 45 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
353
354# 45 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
355!$omp target teams loop defaultmap(firstprivate:scalar) bind(teams,parallel) defaultmap(tofrom:aggregate) defaultmap(tofrom:allocatable) defaultmap(tofrom:pointer) private(i, patch_id, patch_geometry) map(tofrom:gps) map(to:patch_ib(1:num_ibs), Np)
356# 45 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
357#endif
358# 45 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
359
360 do i = 1, num_gps
361
362 patch_id = gps(i)%ib_patch_id
363 patch_geometry = patch_ib(patch_id)%geometry
364
365 if (patch_geometry == 8) then
366 call s_sphere_levelset(gps(i))
367 elseif (patch_geometry == 9) then
368 call s_cuboid_levelset(gps(i))
369 elseif (patch_geometry == 10) then
370 call s_cylinder_levelset(gps(i))
371 elseif (patch_geometry == 11) then
372 call s_3d_airfoil_levelset(gps(i))
373 end if
374 end do
375
376# 61 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
377
378# 61 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
379#if defined(MFC_OpenACC)
380# 61 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
381!$acc end parallel loop
382# 61 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
383#elif defined(MFC_OpenMP)
384# 61 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
385
386# 61 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
387
388# 61 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
389!$omp end target teams loop
390# 61 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
391#endif
392# 61 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
393
394
395 ! 2D Patch Geometries
396 elseif (n > 0) then
397
398
399# 66 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
400
401# 66 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
402#if defined(MFC_OpenACC)
403# 66 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
404!$acc parallel loop gang vector default(present) private(i, patch_id, patch_geometry) copy(gps) copyin(Np, patch_ib(1:num_ibs))
405# 66 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
406#elif defined(MFC_OpenMP)
407# 66 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
408
409# 66 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
410
411# 66 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
412
413# 66 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
414!$omp target teams loop defaultmap(firstprivate:scalar) bind(teams,parallel) defaultmap(tofrom:aggregate) defaultmap(tofrom:allocatable) defaultmap(tofrom:pointer) private(i, patch_id, patch_geometry) map(tofrom:gps) map(to:Np, patch_ib(1:num_ibs))
415# 66 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
416#endif
417# 66 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
418
419 do i = 1, num_gps
420
421 patch_id = gps(i)%ib_patch_id
422 patch_geometry = patch_ib(patch_id)%geometry
423
424 if (patch_geometry == 2) then
425 call s_circle_levelset(gps(i))
426 elseif (patch_geometry == 3) then
427 call s_rectangle_levelset(gps(i))
428 elseif (patch_geometry == 4) then
429 call s_airfoil_levelset(gps(i))
430 elseif (patch_geometry == 6) then
431 call s_ellipse_levelset(gps(i))
432 end if
433 end do
434
435# 82 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
436
437# 82 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
438#if defined(MFC_OpenACC)
439# 82 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
440!$acc end parallel loop
441# 82 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
442#elif defined(MFC_OpenMP)
443# 82 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
444
445# 82 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
446
447# 82 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
448!$omp end target teams loop
449# 82 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
450#endif
451# 82 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
452
453
454 end if
455
456 ! STL models computed on the CPU for now
457 do i = 1, num_gps
458 patch_id = gps(i)%ib_patch_id
459 patch_geometry = patch_ib(patch_id)%geometry
460
461 if (patch_geometry == 5 .or. patch_geometry == 12) then
462 call s_model_levelset(gps(i))
463
464# 93 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
465#if defined(MFC_OpenACC)
466# 93 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
467!$acc update device(gps(i))
468# 93 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
469#elif defined(MFC_OpenMP)
470# 93 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
471!$omp target update to(gps(i))
472# 93 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
473#endif
474 end if
475 end do
476
477
478# 97 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
479#if defined(MFC_OpenACC)
480# 97 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
481!$acc update host(gps(1:num_gps))
482# 97 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
483#elif defined(MFC_OpenMP)
484# 97 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
485!$omp target update from(gps(1:num_gps))
486# 97 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
487#endif
488
489 end subroutine s_apply_levelset
490
491 !> @brief Computes the signed distance and outward normal from a ghost point to a circular immersed boundary.
492 subroutine s_circle_levelset(gp)
493
494
495# 104 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
496#if MFC_OpenACC
497# 104 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
498!$acc routine seq
499# 104 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
500#elif MFC_OpenMP
501# 104 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
502
503# 104 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
504
505# 104 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
506!$omp declare target device_type(any)
507# 104 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
508#endif
509
510 type(ghost_point), intent(inout) :: gp
511
512 real(wp) :: radius, dist
513 real(wp), dimension(2) :: center
514 real(wp), dimension(3) :: dist_vec
515
516 integer :: i, j, ib_patch_id !< Loop index variables
517
518 ib_patch_id = gp%ib_patch_id
519 i = gp%loc(1)
520 j = gp%loc(2)
521
522 radius = patch_ib(ib_patch_id)%radius
523
524 dist_vec(1) = x_cc(i) - patch_ib(ib_patch_id)%x_centroid
525 dist_vec(2) = y_cc(j) - patch_ib(ib_patch_id)%y_centroid
526 dist_vec(3) = 0._wp
527 dist = sqrt(sum(dist_vec**2))
528
529 gp%levelset = dist - radius
530 if (f_approx_equal(dist, 0._wp)) then
531 gp%levelset_norm = 0._wp
532 else
533 gp%levelset_norm = dist_vec(:)/dist
534 end if
535
536 end subroutine s_circle_levelset
537
538 !> @brief Computes the signed distance and outward normal from a ghost point to a 2D NACA airfoil surface.
539 subroutine s_airfoil_levelset(gp)
540
541
542# 137 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
543#if MFC_OpenACC
544# 137 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
545!$acc routine seq
546# 137 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
547#elif MFC_OpenMP
548# 137 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
549
550# 137 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
551
552# 137 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
553!$omp declare target device_type(any)
554# 137 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
555#endif
556
557 type(ghost_point), intent(inout) :: gp
558
559 real(wp) :: dist, global_dist
560 integer :: global_id
561 real(wp), dimension(3) :: dist_vec
562
563 real(wp), dimension(1:3) :: xy_local, offset !< x and y coordinates in local IB frame
564 real(wp), dimension(1:2) :: center
565 real(wp), dimension(1:3, 1:3) :: rotation, inverse_rotation
566
567 integer :: i, j, k, ib_patch_id !< Loop index variables
568
569 ib_patch_id = gp%ib_patch_id
570 i = gp%loc(1)
571 j = gp%loc(2)
572
573 center(1) = patch_ib(ib_patch_id)%x_centroid
574 center(2) = patch_ib(ib_patch_id)%y_centroid
575 inverse_rotation(:, :) = patch_ib(ib_patch_id)%rotation_matrix_inverse(:, :)
576 rotation(:, :) = patch_ib(ib_patch_id)%rotation_matrix(:, :)
577 offset(:) = patch_ib(ib_patch_id)%centroid_offset(:)
578
579 xy_local = [x_cc(i) - center(1), y_cc(j) - center(2), 0._wp] ! get coordinate frame centered on IB
580 xy_local = matmul(inverse_rotation, xy_local) ! rotate the frame into the IB's coordinate
581 xy_local = xy_local - offset ! airfoils are a patch that require a centroid offset
582
583 if (xy_local(2) >= 0._wp) then
584 ! finds the location on the airfoil grid with the minimum distance (closest)
585 do k = 1, np
586 dist_vec(1) = xy_local(1) - airfoil_grid_u(k)%x
587 dist_vec(2) = xy_local(2) - airfoil_grid_u(k)%y
588 dist_vec(3) = 0._wp
589 dist = sqrt(sum(dist_vec**2))
590 if (k == 1) then
591 global_dist = dist
592 global_id = k
593 else
594 if (dist < global_dist) then
595 global_dist = dist
596 global_id = k
597 end if
598 end if
599 end do
600 dist_vec(1) = xy_local(1) - airfoil_grid_u(global_id)%x
601 dist_vec(2) = xy_local(2) - airfoil_grid_u(global_id)%y
602 dist_vec(3) = 0
603 dist = global_dist
604 else
605 do k = 1, np
606 dist_vec(1) = xy_local(1) - airfoil_grid_l(k)%x
607 dist_vec(2) = xy_local(2) - airfoil_grid_l(k)%y
608 dist_vec(3) = 0
609 dist = sqrt(sum(dist_vec**2))
610 if (k == 1) then
611 global_dist = dist
612 global_id = k
613 else
614 if (dist < global_dist) then
615 global_dist = dist
616 global_id = k
617 end if
618 end if
619 end do
620 dist_vec(1) = xy_local(1) - airfoil_grid_l(global_id)%x
621 dist_vec(2) = xy_local(2) - airfoil_grid_l(global_id)%y
622 dist_vec(3) = 0
623 dist = global_dist
624 end if
625
626 gp%levelset = dist
627 if (f_approx_equal(dist, 0._wp)) then
628 gp%levelset_norm = 0._wp
629 else
630 gp%levelset_norm = matmul(rotation, dist_vec(:))/dist ! convert the normal vector back to global grid coordinates
631 end if
632
633 end subroutine s_airfoil_levelset
634
635 !> @brief Computes the signed distance and outward normal from a ghost point to a 3D extruded airfoil surface including spanwise end caps.
637
638
639# 220 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
640#if MFC_OpenACC
641# 220 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
642!$acc routine seq
643# 220 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
644#elif MFC_OpenMP
645# 220 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
646
647# 220 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
648
649# 220 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
650!$omp declare target device_type(any)
651# 220 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
652#endif
653
654 type(ghost_point), intent(inout) :: gp
655
656 real(wp) :: dist, dist_surf, dist_side, global_dist
657 integer :: global_id
658 real(wp) :: lz, z_max, z_min
659 real(wp), dimension(3) :: dist_vec
660
661 real(wp), dimension(1:3) :: xyz_local, center, offset, normal !< x, y, z coordinates in local IB frame
662 real(wp), dimension(1:3, 1:3) :: rotation, inverse_rotation
663
664 real(wp) :: length_z
665
666 integer :: i, j, k, l, ib_patch_id !< Loop index variables
667
668 ib_patch_id = gp%ib_patch_id
669 i = gp%loc(1)
670 j = gp%loc(2)
671 l = gp%loc(3)
672
673 center(1) = patch_ib(ib_patch_id)%x_centroid
674 center(2) = patch_ib(ib_patch_id)%y_centroid
675 center(3) = patch_ib(ib_patch_id)%z_centroid
676 lz = patch_ib(ib_patch_id)%length_z
677 inverse_rotation(:, :) = patch_ib(ib_patch_id)%rotation_matrix_inverse(:, :)
678 rotation(:, :) = patch_ib(ib_patch_id)%rotation_matrix(:, :)
679 offset(:) = patch_ib(ib_patch_id)%centroid_offset(:)
680
681 z_max = lz/2
682 z_min = -lz/2
683
684 xyz_local = [x_cc(i), y_cc(j), z_cc(l)] - center
685 xyz_local = matmul(inverse_rotation, xyz_local) ! rotate the frame into the IB's coordinates
686 xyz_local = xyz_local - offset ! airfoils are a patch that require a centroid offset
687
688 if (xyz_local(2) >= 0._wp) then
689 do k = 1, np
690 dist_vec(1) = xyz_local(1) - airfoil_grid_u(k)%x
691 dist_vec(2) = xyz_local(2) - airfoil_grid_u(k)%y
692 dist_vec(3) = 0._wp
693 dist_surf = sqrt(sum(dist_vec**2))
694 if (k == 1) then
695 global_dist = dist_surf
696 global_id = k
697 else
698 if (dist_surf < global_dist) then
699 global_dist = dist_surf
700 global_id = k
701 end if
702 end if
703 end do
704 dist_vec(1) = xyz_local(1) - airfoil_grid_u(global_id)%x
705 dist_vec(2) = xyz_local(2) - airfoil_grid_u(global_id)%y
706 dist_vec(3) = 0._wp
707 dist_surf = global_dist
708 else
709 do k = 1, np
710 dist_vec(1) = xyz_local(1) - airfoil_grid_l(k)%x
711 dist_vec(2) = xyz_local(2) - airfoil_grid_l(k)%y
712 dist_vec(3) = 0
713 dist_surf = sqrt(sum(dist_vec**2))
714 if (k == 1) then
715 global_dist = dist_surf
716 global_id = k
717 else
718 if (dist_surf < global_dist) then
719 global_dist = dist_surf
720 global_id = k
721 end if
722 end if
723 end do
724 dist_vec(1) = xyz_local(1) - airfoil_grid_l(global_id)%x
725 dist_vec(2) = xyz_local(2) - airfoil_grid_l(global_id)%y
726 dist_vec(3) = 0._wp
727 dist_surf = global_dist
728 end if
729
730 dist_side = min(abs(xyz_local(3) - z_min), abs(z_max - xyz_local(3)))
731
732 if (dist_side < dist_surf) then
733 gp%levelset = dist_side
734 normal = 0._wp
735 if (f_approx_equal(dist_side, abs(xyz_local(3) - z_min))) then
736 normal(3) = -1._wp
737 else
738 normal(3) = 1._wp
739 end if
740 gp%levelset_norm = matmul(rotation, normal)
741 else
742 gp%levelset = dist_surf
743 if (f_approx_equal(dist_surf, 0._wp)) then
744 gp%levelset_norm = 0._wp
745 else
746 gp%levelset_norm = matmul(rotation, dist_vec(:)/dist_surf)
747 end if
748 end if
749
750 end subroutine s_3d_airfoil_levelset
751
752 !> Subroutine for computing the levelset values at a ghost point belonging to the rectangle IB
753 subroutine s_rectangle_levelset(gp)
754
755
756# 323 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
757#if MFC_OpenACC
758# 323 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
759!$acc routine seq
760# 323 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
761#elif MFC_OpenMP
762# 323 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
763
764# 323 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
765
766# 323 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
767!$omp declare target device_type(any)
768# 323 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
769#endif
770
771 type(ghost_point), intent(inout) :: gp
772
773 real(wp) :: top_right(2), bottom_left(2)
774 real(wp) :: min_dist
775 real(wp) :: side_dists(4)
776
777 real(wp) :: length_x, length_y
778 real(wp), dimension(1:3) :: xy_local, dist_vec !< x and y coordinates in local IB frame
779 real(wp), dimension(2) :: center !< x and y coordinates in local IB frame
780 real(wp), dimension(1:3, 1:3) :: rotation, inverse_rotation
781
782 integer :: i, j, k !< Loop index variables
783 integer :: idx !< Shortest path direction indicator
784 integer :: ib_patch_id !< patch ID
785
786 ib_patch_id = gp%ib_patch_id
787 i = gp%loc(1)
788 j = gp%loc(2)
789
790 length_x = patch_ib(ib_patch_id)%length_x
791 length_y = patch_ib(ib_patch_id)%length_y
792 center(1) = patch_ib(ib_patch_id)%x_centroid
793 center(2) = patch_ib(ib_patch_id)%y_centroid
794 inverse_rotation(:, :) = patch_ib(ib_patch_id)%rotation_matrix_inverse(:, :)
795 rotation(:, :) = patch_ib(ib_patch_id)%rotation_matrix(:, :)
796
797 top_right(1) = length_x/2
798 top_right(2) = length_y/2
799 bottom_left(1) = -length_x/2
800 bottom_left(2) = -length_y/2
801
802 ! convert grid to local coordinates
803 xy_local = [x_cc(i) - center(1), y_cc(j) - center(2), 0._wp]
804 xy_local = matmul(inverse_rotation, xy_local)
805
806 side_dists(1) = bottom_left(1) - xy_local(1)
807 side_dists(2) = top_right(1) - xy_local(1)
808 side_dists(3) = bottom_left(2) - xy_local(2)
809 side_dists(4) = top_right(2) - xy_local(2)
810 min_dist = side_dists(1)
811 idx = 1
812
813 do k = 2, 4
814 if (abs(side_dists(k)) < abs(min_dist)) then
815 idx = k
816 min_dist = side_dists(idx)
817 end if
818 end do
819
820 gp%levelset = side_dists(idx)
821 dist_vec = 0._wp
822 if (.not. f_approx_equal(side_dists(idx), 0._wp)) then
823 if (idx == 1 .or. idx == 2) then
824 ! vector points along the x axis
825 dist_vec(1) = side_dists(idx)/abs(side_dists(idx))
826 else
827 ! vector points along the y axis
828 dist_vec(2) = side_dists(idx)/abs(side_dists(idx))
829 end if
830 ! convert the normal vector back into the global coordinate system
831 gp%levelset_norm = matmul(rotation, dist_vec)
832 else
833 gp%levelset_norm = 0._wp
834 end if
835
836 end subroutine s_rectangle_levelset
837
838 !> @brief Computes the signed distance and outward normal from a ghost point to an elliptical immersed boundary via a quadratic projection.
839 subroutine s_ellipse_levelset(gp)
840
841
842# 395 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
843#if MFC_OpenACC
844# 395 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
845!$acc routine seq
846# 395 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
847#elif MFC_OpenMP
848# 395 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
849
850# 395 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
851
852# 395 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
853!$omp declare target device_type(any)
854# 395 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
855#endif
856
857 type(ghost_point), intent(inout) :: gp
858
859 real(wp) :: ellipse_coeffs(2) ! a and b in the ellipse equation
860 real(wp) :: quadratic_coeffs(3) ! A, B, C in the quadratic equation to compute levelset
861
862 real(wp) :: length_x, length_y
863 real(wp), dimension(1:3) :: xy_local, normal_vector !< x and y coordinates in local IB frame
864 real(wp), dimension(2) :: center !< x and y coordinates in local IB frame
865 real(wp), dimension(1:3, 1:3) :: rotation, inverse_rotation
866
867 integer :: i, j, k !< Loop index variables
868 integer :: idx !< Shortest path direction indicator
869 integer :: ib_patch_id !< patch ID
870
871 ib_patch_id = gp%ib_patch_id
872 i = gp%loc(1)
873 j = gp%loc(2)
874
875 length_x = patch_ib(ib_patch_id)%length_x
876 length_y = patch_ib(ib_patch_id)%length_y
877 center(1) = patch_ib(ib_patch_id)%x_centroid
878 center(2) = patch_ib(ib_patch_id)%y_centroid
879 inverse_rotation(:, :) = patch_ib(ib_patch_id)%rotation_matrix_inverse(:, :)
880 rotation(:, :) = patch_ib(ib_patch_id)%rotation_matrix(:, :)
881
882 ellipse_coeffs(1) = 0.5_wp*length_x
883 ellipse_coeffs(2) = 0.5_wp*length_y
884
885 xy_local = [x_cc(i) - center(1), y_cc(j) - center(2), 0._wp]
886 xy_local = matmul(inverse_rotation, xy_local)
887
888 normal_vector = xy_local
889 normal_vector(2) = normal_vector(2)*(ellipse_coeffs(1)/ellipse_coeffs(2))**2._wp ! get the normal direction via the coordinate transformation method
890 normal_vector = normal_vector/sqrt(dot_product(normal_vector, normal_vector)) ! normalize the vector
891 gp%levelset_norm = matmul(rotation, normal_vector) ! save after rotating the vector to the global frame
892
893 ! use the normal vector to set up the quadratic equation for the levelset, using A, B, and C in indices 1, 2, and 3
894 quadratic_coeffs(1) = (normal_vector(1)/ellipse_coeffs(1))**2 + (normal_vector(2)/ellipse_coeffs(2))**2
895 quadratic_coeffs(2) = 2._wp*((xy_local(1)*normal_vector(1)/(ellipse_coeffs(1)**2)) + (xy_local(2)*normal_vector(2)/(ellipse_coeffs(2)**2)))
896 quadratic_coeffs(3) = (xy_local(1)/ellipse_coeffs(1))**2._wp + (xy_local(2)/ellipse_coeffs(2))**2._wp - 1._wp
897
898 ! compute the levelset with the quadratic equation [ -B + sqrt(B^2 - 4AC) ] / 2A
899 gp%levelset = -0.5_wp*(-quadratic_coeffs(2) + sqrt(quadratic_coeffs(2)**2._wp - 4._wp*quadratic_coeffs(1)*quadratic_coeffs(3)))/quadratic_coeffs(1)
900
901 end subroutine s_ellipse_levelset
902
903 !> @brief Computes the signed distance and outward normal from a ghost point to the nearest face of a cuboid immersed boundary.
904 subroutine s_cuboid_levelset(gp)
905
906
907# 446 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
908#if MFC_OpenACC
909# 446 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
910!$acc routine seq
911# 446 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
912#elif MFC_OpenMP
913# 446 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
914
915# 446 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
916
917# 446 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
918!$omp declare target device_type(any)
919# 446 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
920#endif
921
922 type(ghost_point), intent(inout) :: gp
923
924 real(wp) :: Right, Left, Bottom, Top, Front, Back
925 real(wp) :: min_dist
926 real(wp) :: dist_left, dist_right, dist_bottom, dist_top, dist_back, dist_front
927
928 real(wp), dimension(3) :: center
929 real(wp) :: length_x, length_y, length_z
930 real(wp), dimension(1:3) :: xyz_local, dist_vec !< x and y coordinates in local IB frame
931 real(wp), dimension(1:3, 1:3) :: rotation, inverse_rotation
932
933 integer :: i, j, k !< Loop index variables
934 integer :: ib_patch_id !< patch ID
935
936 ib_patch_id = gp%ib_patch_id
937 i = gp%loc(1)
938 j = gp%loc(2)
939 k = gp%loc(3)
940
941 length_x = patch_ib(ib_patch_id)%length_x
942 length_y = patch_ib(ib_patch_id)%length_y
943 length_z = patch_ib(ib_patch_id)%length_z
944
945 center(1) = patch_ib(ib_patch_id)%x_centroid
946 center(2) = patch_ib(ib_patch_id)%y_centroid
947 center(3) = patch_ib(ib_patch_id)%z_centroid
948
949 inverse_rotation(:, :) = patch_ib(ib_patch_id)%rotation_matrix_inverse(:, :)
950 rotation(:, :) = patch_ib(ib_patch_id)%rotation_matrix(:, :)
951
952 right = length_x/2
953 left = -length_x/2
954 top = length_y/2
955 bottom = -length_y/2
956 front = length_z/2
957 back = -length_z/2
958
959 xyz_local = [x_cc(i), y_cc(j), z_cc(k)] - center ! get coordinate frame centered on IB
960 xyz_local = matmul(inverse_rotation, xyz_local) ! rotate the frame into the IB's coordinate
961
962 dist_left = left - xyz_local(1)
963 dist_right = xyz_local(1) - right
964 dist_bottom = bottom - xyz_local(2)
965 dist_top = xyz_local(2) - top
966 dist_back = back - xyz_local(3)
967 dist_front = xyz_local(3) - front
968
969 min_dist = min(abs(dist_left), abs(dist_right), abs(dist_bottom), &
970 abs(dist_top), abs(dist_back), abs(dist_front))
971
972 ! TODO :: The way that this is written, it looks like we will
973 ! trigger at the first size that is close to the minimum distance,
974 ! meaning corners where side_dists are the same will
975 ! trigger on what may not actually be the minimum,
976 ! leading to undesired behavior. This should be resolved
977 ! and this code should be cleaned up. It also means that
978 ! rotating the box 90 degrees will cause tests to fail.
979
980 dist_vec = 0._wp
981
982 if (f_approx_equal(min_dist, abs(dist_left))) then
983 gp%levelset = dist_left
984 if (.not. f_approx_equal(dist_left, 0._wp)) then
985 dist_vec(1) = dist_left/abs(dist_left)
986 end if
987 else if (f_approx_equal(min_dist, abs(dist_right))) then
988 gp%levelset = dist_right
989 if (.not. f_approx_equal(dist_right, 0._wp)) then
990 dist_vec(1) = -dist_right/abs(dist_right)
991 end if
992 else if (f_approx_equal(min_dist, abs(dist_bottom))) then
993 gp%levelset = dist_bottom
994 if (.not. f_approx_equal(dist_bottom, 0._wp)) then
995 dist_vec(2) = dist_bottom/abs(dist_bottom)
996 end if
997 else if (f_approx_equal(min_dist, abs(dist_top))) then
998 gp%levelset = dist_top
999 if (.not. f_approx_equal(dist_top, 0._wp)) then
1000 dist_vec(2) = -dist_top/abs(dist_top)
1001 end if
1002 else if (f_approx_equal(min_dist, abs(dist_back))) then
1003 gp%levelset = dist_back
1004 if (.not. f_approx_equal(dist_back, 0._wp)) then
1005 dist_vec(3) = dist_back/abs(dist_back)
1006 end if
1007 else if (f_approx_equal(min_dist, abs(dist_front))) then
1008 gp%levelset = dist_front
1009 if (.not. f_approx_equal(dist_front, 0._wp)) then
1010 dist_vec(3) = -dist_front/abs(dist_front)
1011 end if
1012 end if
1013
1014 gp%levelset_norm = matmul(rotation, dist_vec)
1015
1016 end subroutine s_cuboid_levelset
1017
1018 !> @brief Computes the signed distance and outward normal from a ghost point to a spherical immersed boundary.
1019 subroutine s_sphere_levelset(gp)
1020
1021
1022# 547 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1023#if MFC_OpenACC
1024# 547 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1025!$acc routine seq
1026# 547 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1027#elif MFC_OpenMP
1028# 547 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1029
1030# 547 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1031
1032# 547 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1033!$omp declare target device_type(any)
1034# 547 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1035#endif
1036
1037 type(ghost_point), intent(inout) :: gp
1038
1039 real(wp) :: radius, dist
1040 real(wp), dimension(3) :: dist_vec, center
1041
1042 integer :: i, j, k, ib_patch_id !< Loop index variables
1043
1044 ib_patch_id = gp%ib_patch_id
1045 i = gp%loc(1)
1046 j = gp%loc(2)
1047 k = gp%loc(3)
1048
1049 radius = patch_ib(ib_patch_id)%radius
1050 center(1) = patch_ib(ib_patch_id)%x_centroid
1051 center(2) = patch_ib(ib_patch_id)%y_centroid
1052 center(3) = patch_ib(ib_patch_id)%z_centroid
1053
1054 dist_vec(1) = x_cc(i) - center(1)
1055 dist_vec(2) = y_cc(j) - center(2)
1056 dist_vec(3) = z_cc(k) - center(3)
1057 dist = sqrt(sum(dist_vec**2))
1058 gp%levelset = dist - radius
1059 if (f_approx_equal(dist, 0._wp)) then
1060 gp%levelset_norm = (/1, 0, 0/)
1061 else
1062 gp%levelset_norm = dist_vec(:)/dist
1063 end if
1064
1065 end subroutine s_sphere_levelset
1066
1067 !> @brief Computes the signed distance and outward normal from a ghost point to a cylindrical immersed boundary surface and end caps.
1068 subroutine s_cylinder_levelset(gp)
1069
1070
1071# 582 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1072#if MFC_OpenACC
1073# 582 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1074!$acc routine seq
1075# 582 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1076#elif MFC_OpenMP
1077# 582 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1078
1079# 582 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1080
1081# 582 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1082!$omp declare target device_type(any)
1083# 582 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1084#endif
1085
1086 type(ghost_point), intent(inout) :: gp
1087
1088 real(wp) :: radius
1089 real(wp), dimension(3) :: dist_sides_vec, dist_surface_vec, length
1090 real(wp), dimension(2) :: boundary
1091 real(wp) :: dist_side, dist_surface, side_pos
1092 integer :: i, j, k !< Loop index variables
1093 integer :: ib_patch_id !< patch ID
1094
1095 real(wp), dimension(1:3) :: xyz_local, center !< x and y coordinates in local IB frame
1096 real(wp), dimension(1:3, 1:3) :: rotation, inverse_rotation
1097
1098 ib_patch_id = gp%ib_patch_id
1099 i = gp%loc(1)
1100 j = gp%loc(2)
1101 k = gp%loc(3)
1102
1103 radius = patch_ib(ib_patch_id)%radius
1104 center(1) = patch_ib(ib_patch_id)%x_centroid
1105 center(2) = patch_ib(ib_patch_id)%y_centroid
1106 center(3) = patch_ib(ib_patch_id)%z_centroid
1107 length(1) = patch_ib(ib_patch_id)%length_x
1108 length(2) = patch_ib(ib_patch_id)%length_y
1109 length(3) = patch_ib(ib_patch_id)%length_z
1110
1111 inverse_rotation(:, :) = patch_ib(ib_patch_id)%rotation_matrix_inverse(:, :)
1112 rotation(:, :) = patch_ib(ib_patch_id)%rotation_matrix(:, :)
1113
1114 if (.not. f_approx_equal(length(1), 0._wp)) then
1115 boundary(1) = -0.5_wp*length(1)
1116 boundary(2) = 0.5_wp*length(1)
1117 dist_sides_vec = (/1, 0, 0/)
1118 dist_surface_vec = (/0, 1, 1/)
1119 else if (.not. f_approx_equal(length(2), 0._wp)) then
1120 boundary(1) = -0.5_wp*length(2)
1121 boundary(2) = 0.5_wp*length(2)
1122 dist_sides_vec = (/0, 1, 0/)
1123 dist_surface_vec = (/1, 0, 1/)
1124 else if (.not. f_approx_equal(length(3), 0._wp)) then
1125 boundary(1) = -0.5_wp*length(3)
1126 boundary(2) = 0.5_wp*length(3)
1127 dist_sides_vec = (/0, 0, 1/)
1128 dist_surface_vec = (/1, 1, 0/)
1129 end if
1130
1131 xyz_local = [x_cc(i), y_cc(j), z_cc(k)] - center ! get coordinate frame centered on IB
1132 xyz_local = matmul(inverse_rotation, xyz_local) ! rotate the frame into the IB's coordinates
1133
1134 ! get distance to flat edge of cylinder
1135 side_pos = dot_product(xyz_local, dist_sides_vec)
1136 dist_side = min(abs(side_pos - boundary(1)), &
1137 abs(boundary(2) - side_pos))
1138 ! get distance to curved side of cylinder
1139 dist_surface = norm2(xyz_local*dist_surface_vec) &
1140 - radius
1141
1142 if (dist_side < abs(dist_surface)) then
1143 ! if the closest edge is flat
1144 gp%levelset = -dist_side
1145 if (f_approx_equal(dist_side, abs(side_pos - boundary(1)))) then
1146 gp%levelset_norm = matmul(rotation, -dist_sides_vec)
1147 else
1148 gp%levelset_norm = matmul(rotation, dist_sides_vec)
1149 end if
1150 else
1151 gp%levelset = dist_surface
1152 xyz_local = xyz_local*dist_surface_vec
1153 xyz_local = xyz_local/max(norm2(xyz_local), sgm_eps)
1154 gp%levelset_norm = matmul(rotation, xyz_local)
1155 end if
1156
1157 end subroutine s_cylinder_levelset
1158
1159 !> The STL patch is a 2/3D geometry that is imported from an STL file.
1160 !! @param gp Ghost point to compute levelset for
1161 subroutine s_model_levelset(gp)
1162
1163# 660 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1164#if MFC_OpenACC
1165# 660 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1166!$acc routine seq
1167# 660 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1168#elif MFC_OpenMP
1169# 660 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1170
1171# 660 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1172
1173# 660 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1174!$omp declare target device_type(any)
1175# 660 "/home/runner/work/MFC/MFC/src/simulation/m_compute_levelset.fpp"
1176#endif
1177
1178 type(ghost_point), intent(inout) :: gp
1179
1180 integer :: i, j, k, patch_id, boundary_edge_count, total_vertices
1181 logical :: interpolate
1182 real(wp), dimension(1:3) :: center, xyz_local
1183 real(wp) :: normals(1:3) !< Boundary normal buffer
1184 real(wp) :: distance
1185 real(wp), dimension(1:3, 1:3) :: inverse_rotation, rotation
1186
1187 patch_id = gp%ib_patch_id
1188 i = gp%loc(1)
1189 j = gp%loc(2)
1190 k = gp%loc(3)
1191
1192 ! load in model values
1193 interpolate = models(patch_id)%interpolate
1194 boundary_edge_count = models(patch_id)%boundary_edge_count
1195 total_vertices = models(patch_id)%total_vertices
1196
1197 center = 0._wp
1198 if (.not. f_is_default(patch_ib(patch_id)%x_centroid)) center(1) = patch_ib(patch_id)%x_centroid
1199 if (.not. f_is_default(patch_ib(patch_id)%y_centroid)) center(2) = patch_ib(patch_id)%y_centroid
1200 if (p > 0) then
1201 if (.not. f_is_default(patch_ib(patch_id)%z_centroid)) center(3) = patch_ib(patch_id)%z_centroid
1202 end if
1203 inverse_rotation(:, :) = patch_ib(patch_id)%rotation_matrix_inverse(:, :)
1204 rotation(:, :) = patch_ib(patch_id)%rotation_matrix(:, :)
1205
1206 ! determine where we are located in space
1207 xyz_local = (/x_cc(i) - center(1), y_cc(j) - center(2), 0._wp/)
1208 if (p > 0) then
1209 xyz_local(3) = z_cc(k) - center(3)
1210 end if
1211 xyz_local = matmul(inverse_rotation, xyz_local)
1212
1213 if (grid_geometry == 3) then
1214 xyz_local = f_convert_cyl_to_cart(xyz_local)
1215 end if
1216
1217 ! 3D models
1218 if (p > 0) then
1219
1220 ! Get the boundary normals and shortest distance between the cell center and the model boundary
1221 call f_distance_normals_3d(models(patch_id)%model, xyz_local, normals, distance)
1222
1223 ! Get the shortest distance between the cell center and the interpolated model boundary
1224 if (interpolate) then
1225 gp%levelset = f_interpolated_distance(models(patch_id)%interpolated_boundary_v, total_vertices, xyz_local)
1226 else
1227 gp%levelset = distance
1228 end if
1229
1230 ! Correct the sign of the levelset
1231 gp%levelset = -abs(gp%levelset)
1232
1233 ! Assign the levelset_norm
1234 gp%levelset_norm = matmul(rotation, normals(1:3))
1235 else
1236 ! 2D models
1237 if (interpolate) then
1238 ! Get the shortest distance between the cell center and the model boundary
1239 gp%levelset = f_interpolated_distance(models(patch_id)%interpolated_boundary_v, total_vertices, xyz_local)
1240 else
1241 ! Get the shortest distance between the cell center and the interpolated model boundary
1242 gp%levelset = f_distance(models(patch_id)%boundary_v, boundary_edge_count, xyz_local)
1243 end if
1244
1245 ! Correct the sign of the levelset
1246 gp%levelset = -abs(gp%levelset)
1247
1248 ! Get the boundary normals
1249 call f_normals(models(patch_id)%boundary_v, &
1250 boundary_edge_count, &
1251 xyz_local, &
1252 normals)
1253
1254 ! Assign the levelset_norm
1255 gp%levelset_norm = matmul(rotation, normals(1:3))
1256
1257 end if
1258
1259 end subroutine s_model_levelset
1260
1261end module m_compute_levelset
Computes signed-distance level-set fields and surface normals for immersed-boundary patch geometries.
impure subroutine, public s_apply_levelset(gps, num_gps)
Dispatches level-set distance and normal computations for all ghost points based on their patch geome...
subroutine s_airfoil_levelset(gp)
Computes the signed distance and outward normal from a ghost point to a 2D NACA airfoil surface.
subroutine s_sphere_levelset(gp)
Computes the signed distance and outward normal from a ghost point to a spherical immersed boundary.
subroutine s_circle_levelset(gp)
Computes the signed distance and outward normal from a ghost point to a circular immersed boundary.
subroutine s_ellipse_levelset(gp)
Computes the signed distance and outward normal from a ghost point to an elliptical immersed boundary...
subroutine s_cylinder_levelset(gp)
Computes the signed distance and outward normal from a ghost point to a cylindrical immersed boundary...
subroutine s_3d_airfoil_levelset(gp)
Computes the signed distance and outward normal from a ghost point to a 3D extruded airfoil surface i...
subroutine s_rectangle_levelset(gp)
Subroutine for computing the levelset values at a ghost point belonging to the rectangle IB.
subroutine s_model_levelset(gp)
The STL patch is a 2/3D geometry that is imported from an STL file.
subroutine s_cuboid_levelset(gp)
Computes the signed distance and outward normal from a ghost point to the nearest face of a cuboid im...
Shared derived types for field data, patch geometry, bubble dynamics, and MPI I/O structures.
Global parameters for the computational domain, fluid properties, and simulation algorithm configurat...
type(ib_patch_parameters), dimension(num_patches_max) patch_ib
type(vec3_dt), dimension(:), allocatable airfoil_grid_u
real(wp), dimension(:), allocatable, target y_cc
real(wp), dimension(:), allocatable, target z_cc
type(vec3_dt), dimension(:), allocatable airfoil_grid_l
real(wp), dimension(:), allocatable, target x_cc
Basic floating-point utilities: approximate equality, default detection, and coordinate bounds.
logical elemental function, public f_approx_equal(a, b, tol_input)
This procedure checks if two floating point numbers of wp are within tolerance.
logical elemental function, public f_is_default(var)
Checks if a real(wp) variable is of default value.
Allocate memory and read initial condition data for IC extrusion.
type(t_model_array), dimension(:), allocatable, target, public models
pure real(wp) function, dimension(1:3), public f_convert_cyl_to_cart(cyl)
Converts a 3D cylindrical coordinate vector (x, r, theta) to Cartesian (x, y, z).
Binary STL file reader and processor for immersed boundary geometry.
subroutine, public f_distance_normals_3d(model, point, normals, distance)
This procedure determines the levelset distance and normals of the 3D models without interpolation.
real(wp) function, public f_interpolated_distance(interpolated_boundary_v, total_vertices, point)
This procedure determines the levelset of interpolated 2D models.
subroutine, public f_normals(boundary_v, boundary_edge_count, point, normals)
This procedure determines the levelset normals of 2D models without interpolation.
real(wp) function, public f_distance(boundary_v, boundary_edge_count, point)
This procedure determines the levelset distance of 2D models without interpolation.
MPI halo exchange, domain decomposition, and buffer packing/unpacking for the simulation solver.
Ghost Point for Immersed Boundaries.