Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
signed_distance.f90
Go to the documentation of this file.
1! Copyright (c) 2024, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
35 use num_types, only : dp, rp
36 use field, only : field_t
37 use tri, only : tri_t
38 use tri_mesh, only : tri_mesh_t
39 use aabb_tree, only : aabb_tree_t
42 use utils, only : neko_error, neko_warning
43
44 implicit none
45 private
46 public :: signed_distance_field
47
48contains
49
60 subroutine signed_distance_field(field_data, object, max_distance)
61 type(field_t), intent(inout) :: field_data
62 class(*), intent(in) :: object
63 real(kind=dp), intent(in), optional :: max_distance
64
65 real(kind=dp) :: max_dist
66
67 if (present(max_distance)) then
68 max_dist = max_distance
69 else
70 max_dist = huge(0.0_dp)
71 end if
72
73 select type (object)
74 type is (tri_mesh_t)
75 call signed_distance_field_tri_mesh(field_data, object, max_dist)
76
77 class default
78 call neko_error("signed_distance_field: Object type not supported.")
79 end select
80
81 end subroutine signed_distance_field
82
93 subroutine signed_distance_field_tri_mesh(field_data, mesh, max_distance)
94 type(field_t), intent(inout) :: field_data
95 type(tri_mesh_t), intent(in) :: mesh
96 real(kind=dp), intent(in) :: max_distance
97
98 integer :: total_size
99 integer :: id
100 type(aabb_tree_t) :: search_tree
101 real(kind=dp), dimension(3) :: p
102 real(kind=dp) :: distance
103
104 ! Zero the field
105 field_data%x = 0.0_dp
106 total_size = field_data%dof%size()
107
108 call search_tree%init(mesh%nelv)
109 call search_tree%build(mesh%el)
110
111 if (search_tree%get_size() .ne. mesh%nelv) then
112 call neko_error("signed_distance_field_tri_mesh: &
113 & Error building the search tree.")
114 end if
115
116 do id = 1, total_size
117 p(1) = field_data%dof%x(id, 1, 1, 1)
118 p(2) = field_data%dof%y(id, 1, 1, 1)
119 p(3) = field_data%dof%z(id, 1, 1, 1)
120
121 distance = tri_mesh_aabb_tree(search_tree, mesh%el, p, max_distance)
122
123 field_data%x(id, 1, 1, 1) = real(distance, kind=rp)
124 end do
125
126 if (neko_bcknd_device .eq. 1) then
127 call neko_warning("signed_distance_field_tri_mesh:&
128 & Device version not implemented.")
129 call device_memcpy(field_data%x, field_data%x_d, field_data%size(), &
130 host_to_device, sync = .false.)
131 end if
132
133 end subroutine signed_distance_field_tri_mesh
134
145 function tri_mesh_brute_force(mesh, p, max_distance) result(distance)
146 use tri, only : tri_t
147 use point, only : point_t
148 use num_types, only : dp
149
150 implicit none
151
152 type(tri_mesh_t), intent(in) :: mesh
153 real(kind=dp), intent(in) :: p(3)
154 real(kind=dp), intent(in) :: max_distance
155
156 integer :: id
157 real(kind=dp) :: distance, weighted_sign
158 real(kind=dp) :: cd, cs
159 real(kind=dp) :: tol = 1e-6_dp
160
161 distance = 1e10_dp
162 weighted_sign = 0.0_dp
163
164 do id = 1, mesh%nelv
165 call element_distance(mesh%el(id), p, cd, cs)
166
167 ! Update the weighted sign, if the relative difference is small
168 if (abs(cd - distance) / distance .lt. tol) then
169 weighted_sign = weighted_sign + cs
170 else if (cd .lt. distance) then
171 weighted_sign = cs
172 end if
173
174 distance = min(cd, distance)
175 end do
176
177 distance = sign(min(distance, max_distance), weighted_sign)
178
179 end function tri_mesh_brute_force
180
194 function tri_mesh_aabb_tree(tree, object_list, p, max_distance) &
195 result(distance)
196 use aabb, only : aabb_t
198 use stack, only : stack_i4_t
199 implicit none
200
201 class(aabb_tree_t), intent(in) :: tree
202 class(tri_t), contiguous, dimension(:), intent(in) :: object_list
203 real(kind=dp), dimension(3), intent(in) :: p
204 real(kind=dp), intent(in) :: max_distance
205
206 real(kind=dp) :: distance
207 real(kind=dp) :: weighted_sign
208
209 real(kind=dp), parameter :: tol = 1.0e-6_dp
210
211 type(stack_i4_t) :: simple_stack
212 integer :: current_index
213
214 type(aabb_node_t) :: current_node
215 type(aabb_t) :: current_aabb
216 integer :: current_object_index
217 real(kind=dp) :: current_distance
218 real(kind=dp) :: current_sign
219
220 type(aabb_node_t) :: left_node
221 type(aabb_node_t) :: right_node
222
223 type(aabb_t) :: root_box
224 type(aabb_t) :: search_box
225
226 integer :: root_index, left_index, right_index
227 real(kind=dp) :: random_value
228
229 ! Initialize the stack and the search box
230 call simple_stack%init(size(object_list) * 2)
231 call search_box%init(p - max_distance, p + max_distance)
232
233 ! Check if the root node overlaps the search box, if it does, push it to
234 ! the stack and update the search box to a randomly selected object.
235 root_index = tree%get_root_index()
236 root_box = tree%get_aabb(root_index)
237
238 if (.not. root_box%overlaps(search_box)) then
239 distance = max_distance
240 weighted_sign = 1.0_dp
241 return
242 end if
243
244 ! Grab a random object and compute the distance to it
245 call random_number(random_value)
246 current_object_index = floor(random_value * size(object_list) + 1)
247 call element_distance(object_list(current_object_index), p, distance, &
248 weighted_sign)
249 distance = distance + object_list(current_object_index)%diameter()
250
251 ! Update the search box to the new distance and push the root node
252 call search_box%init(p - distance, p + distance)
253 call simple_stack%push(root_index)
254
255 ! Traverse the tree and compute the signed distance to the elements
256 do while (.not. simple_stack%is_empty())
257 current_index = simple_stack%pop()
258 if (current_index .eq. aabb_null_node) cycle
259
260 current_node = tree%get_node(current_index)
261 current_aabb = current_node%get_aabb()
262
263 if (current_node%is_leaf()) then
264 if (distance .lt. current_node%min_distance(p)) then
265 cycle
266 end if
267
268 current_object_index = current_node%get_object_index()
269 call element_distance(object_list(current_object_index), p, &
270 current_distance, current_sign)
271
272 ! Update the weighted sign, if the relative difference is small
273 if (abs(current_distance - distance) / distance .lt. tol) then
274 weighted_sign = weighted_sign + current_sign
275 else if (current_distance .lt. distance) then
276 weighted_sign = current_sign
277 end if
278
279 distance = min(distance, current_distance)
280
281 ! Update the search box to the new distance
282 if (distance .gt. current_aabb%get_diameter()) then
283 call search_box%init(p - distance, p + distance)
284 end if
285 else
286
287 left_index = tree%get_left_index(current_index)
288 if (left_index .ne. aabb_null_node) then
289 left_node = tree%get_left_node(current_index)
290 if (left_node%aabb%overlaps(search_box)) then
291 call simple_stack%push(left_index)
292 end if
293 end if
294
295 right_index = tree%get_right_index(current_index)
296 if (right_index .ne. aabb_null_node) then
297 right_node = tree%get_right_node(current_index)
298 if (right_node%aabb%overlaps(search_box)) then
299 call simple_stack%push(right_index)
300 end if
301 end if
302 end if
303 end do
304
305 if (distance .gt. max_distance) then
306 distance = max_distance
307 end if
308 distance = sign(distance, weighted_sign)
309
310 end function tri_mesh_aabb_tree
311
312 ! ========================================================================== !
313 ! Element distance functions
314 ! ========================================================================== !
315
317 subroutine element_distance(element, p, distance, weighted_sign)
318 class(*), intent(in) :: element
319 real(kind=dp), dimension(3), intent(in) :: p
320 real(kind=dp), intent(out) :: distance
321 real(kind=dp), intent(out), optional :: weighted_sign
322
323 select type (element)
324 type is (tri_t)
325 call element_distance_triangle(element, p, distance, weighted_sign)
326
327 class default
328 print *, "Error: Element type not supported."
329 stop
330 end select
331 end subroutine element_distance
332
333 ! -------------------------------------------------------------------------- !
334 ! Specific element distance functions
335
353 subroutine element_distance_triangle(triangle, p, distance, weighted_sign)
354 type(tri_t), intent(in) :: triangle
355 real(kind=dp), dimension(3), intent(in) :: p
356
357 real(kind=dp), intent(out) :: distance
358 real(kind=dp), intent(out), optional :: weighted_sign
359
360 real(kind=dp), dimension(3) :: v1, v2, v3
361 real(kind=dp), dimension(3) :: normal
362 real(kind=dp) :: normal_length
363 real(kind=dp) :: b1, b2, b3
364
365 real(kind=dp), dimension(3) :: projection
366 real(kind=dp), dimension(3) :: edge
367 real(kind=dp) :: tol = 1e-10_dp
368
369 real(kind=dp) :: face_distance
370 real(kind=dp) :: t
371
372 ! Get vertices and the normal vector
373 v1 = triangle%pts(1)%p%x
374 v2 = triangle%pts(2)%p%x
375 v3 = triangle%pts(3)%p%x
376
377 normal = cross(v2 - v1, v3 - v1)
378 normal_length = norm2(normal)
379
380 if (normal_length .lt. tol) then
381 distance = huge(1.0_dp)
382 weighted_sign = 0.0_dp
383 return
384 end if
385 normal = normal / normal_length
386
387 ! Compute Barycentric coordinates to determine if the point is inside the
388 ! triangular prism, of along an edge or by a face.
389 face_distance = dot_product(p - v1, normal)
390
391 projection = p - normal * face_distance
392 b1 = dot_product(normal, cross(v2 - v1, projection - v1)) / normal_length
393 b2 = dot_product(normal, cross(v3 - v2, projection - v2)) / normal_length
394 b3 = dot_product(normal, cross(v1 - v3, projection - v3)) / normal_length
395
396 if (b1 .le. tol) then
397 edge = v2 - v1
398 t = dot_product(p - v1, edge) / norm2(edge)
399 t = max(0.0_dp, min(1.0_dp, t))
400
401 projection = v1 + t * edge
402 else if (b2 .le. tol) then
403 edge = v3 - v2
404 t = dot_product(p - v2, edge) / norm2(edge)
405 t = max(0.0_dp, min(1.0_dp, t))
406
407 projection = v2 + t * edge
408 else if (b3 .le. tol) then
409 edge = v1 - v3
410 t = dot_product(p - v3, edge) / norm2(edge)
411 t = max(0.0_dp, min(1.0_dp, t))
412
413 projection = v3 + t * edge
414 end if
415
416 distance = norm2(projection - p)
417 if (present(weighted_sign)) then
418 weighted_sign = face_distance / distance
419 else
420 distance = sign(distance, face_distance)
421 end if
422
423 end subroutine element_distance_triangle
424
429 pure function cross(a, b) result(c)
430 real(kind=dp), dimension(3), intent(in) :: a
431 real(kind=dp), dimension(3), intent(in) :: b
432 real(kind=dp), dimension(3) :: c
433
434 c(1) = a(2) * b(3) - a(3) * b(2)
435 c(2) = a(3) * b(1) - a(1) * b(3)
436 c(3) = a(1) * b(2) - a(2) * b(1)
437
438 end function cross
439
440end module signed_distance
double real
Copy data between host and device (or device and device)
Definition device.F90:71
Axis Aligned Bounding Box (aabb) Tree data structure.
Definition aabb_tree.f90:70
integer, parameter, public aabb_null_node
Definition aabb_tree.f90:79
Axis Aligned Bounding Box (aabb) implementation in Fortran.
Definition aabb.f90:71
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
Defines a field.
Definition field.f90:34
Defines a mesh.
Definition mesh.f90:34
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements a point.
Definition point.f90:35
Project x onto X, the space of old solutions and back again.
Module containing Signed Distance Functions.
subroutine signed_distance_field_tri_mesh(field_data, mesh, max_distance)
Signed distance field for a triangular mesh.
pure real(kind=dp) function, dimension(3) cross(a, b)
Compute cross product of two vectors.
subroutine element_distance_triangle(triangle, p, distance, weighted_sign)
Signed distance function for a triangle.
real(kind=dp) function tri_mesh_aabb_tree(tree, object_list, p, max_distance)
Signed distance function using an AABB tree.
subroutine, public signed_distance_field(field_data, object, max_distance)
Signed distance field.
subroutine element_distance(element, p, distance, weighted_sign)
Main interface for the signed distance function for an element.
real(kind=dp) function tri_mesh_brute_force(mesh, p, max_distance)
Signed distance function.
Implements a dynamic stack ADT.
Definition stack.f90:35
Defines a triangular surface mesh.
Definition tri_mesh.f90:35
Defines a triangular element.
Definition tri.f90:34
Utilities.
Definition utils.f90:35
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition utils.f90:346
Axis Aligned Bounding Box (aabb) data structure.
Definition aabb.f90:107
Node type for the Axis Aligned Bounding Box (aabb) Tree.
Definition aabb_tree.f90:86
Axis Aligned Bounding Box (aabb) Tree.
A point in with coordinates .
Definition point.f90:43
Integer based stack.
Definition stack.f90:63
Triangular element.
Definition tri.f90:61
#define max(a, b)
Definition tensor.cu:40