Neko 1.99.3
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 num_types, only : dp
147
148 implicit none
149
150 type(tri_mesh_t), intent(in) :: mesh
151 real(kind=dp), intent(in) :: p(3)
152 real(kind=dp), intent(in) :: max_distance
153
154 integer :: id
155 real(kind=dp) :: distance, weighted_sign
156 real(kind=dp) :: cd, cs
157 real(kind=dp) :: tol = 1e-6_dp
158
159 distance = 1e10_dp
160 weighted_sign = 0.0_dp
161
162 do id = 1, mesh%nelv
163 call element_distance(mesh%el(id), p, cd, cs)
164
165 ! Update the weighted sign, if the relative difference is small
166 if (abs(cd - distance) / distance .lt. tol) then
167 weighted_sign = weighted_sign + cs
168 else if (cd .lt. distance) then
169 weighted_sign = cs
170 end if
171
172 distance = min(cd, distance)
173 end do
174
175 distance = sign(min(distance, max_distance), weighted_sign)
176
177 end function tri_mesh_brute_force
178
192 function tri_mesh_aabb_tree(tree, object_list, p, max_distance) &
193 result(distance)
194 use aabb, only : aabb_t
196 use stack, only : stack_i4_t
197 implicit none
198
199 class(aabb_tree_t), intent(in) :: tree
200 class(tri_t), contiguous, dimension(:), intent(in) :: object_list
201 real(kind=dp), dimension(3), intent(in) :: p
202 real(kind=dp), intent(in) :: max_distance
203
204 real(kind=dp) :: distance
205 real(kind=dp) :: weighted_sign
206
207 real(kind=dp), parameter :: tol = 1.0e-6_dp
208
209 type(stack_i4_t) :: simple_stack
210 integer :: current_index
211
212 type(aabb_node_t) :: current_node
213 type(aabb_t) :: current_aabb
214 integer :: current_object_index
215 real(kind=dp) :: current_distance
216 real(kind=dp) :: current_sign
217
218 type(aabb_node_t) :: left_node
219 type(aabb_node_t) :: right_node
220
221 type(aabb_t) :: root_box
222 type(aabb_t) :: search_box
223
224 integer :: root_index, left_index, right_index
225 real(kind=dp) :: random_value
226
227 ! Initialize the stack and the search box
228 call simple_stack%init(size(object_list) * 2)
229 call search_box%init(p - max_distance, p + max_distance)
230
231 ! Check if the root node overlaps the search box, if it does, push it to
232 ! the stack and update the search box to a randomly selected object.
233 root_index = tree%get_root_index()
234 root_box = tree%get_aabb(root_index)
235
236 if (.not. root_box%overlaps(search_box)) then
237 distance = max_distance
238 weighted_sign = 1.0_dp
239 return
240 end if
241
242 ! Grab a random object and compute the distance to it
243 call random_number(random_value)
244 current_object_index = floor(random_value * size(object_list) + 1)
245 call element_distance(object_list(current_object_index), p, distance, &
246 weighted_sign)
247 distance = distance + object_list(current_object_index)%diameter()
248
249 ! Update the search box to the new distance and push the root node
250 call search_box%init(p - distance, p + distance)
251 call simple_stack%push(root_index)
252
253 ! Traverse the tree and compute the signed distance to the elements
254 do while (.not. simple_stack%is_empty())
255 current_index = simple_stack%pop()
256 if (current_index .eq. aabb_null_node) cycle
257
258 current_node = tree%get_node(current_index)
259 current_aabb = current_node%get_aabb()
260
261 if (current_node%is_leaf()) then
262 if (distance .lt. current_node%min_distance(p)) then
263 cycle
264 end if
265
266 current_object_index = current_node%get_object_index()
267 call element_distance(object_list(current_object_index), p, &
268 current_distance, current_sign)
269
270 ! Update the weighted sign, if the relative difference is small
271 if (abs(current_distance - distance) / distance .lt. tol) then
272 weighted_sign = weighted_sign + current_sign
273 else if (current_distance .lt. distance) then
274 weighted_sign = current_sign
275 end if
276
277 distance = min(distance, current_distance)
278
279 ! Update the search box to the new distance
280 if (distance .gt. current_aabb%get_diameter()) then
281 call search_box%init(p - distance, p + distance)
282 end if
283 else
284
285 left_index = tree%get_left_index(current_index)
286 if (left_index .ne. aabb_null_node) then
287 left_node = tree%get_left_node(current_index)
288 if (left_node%aabb%overlaps(search_box)) then
289 call simple_stack%push(left_index)
290 end if
291 end if
292
293 right_index = tree%get_right_index(current_index)
294 if (right_index .ne. aabb_null_node) then
295 right_node = tree%get_right_node(current_index)
296 if (right_node%aabb%overlaps(search_box)) then
297 call simple_stack%push(right_index)
298 end if
299 end if
300 end if
301 end do
302
303 if (distance .gt. max_distance) then
304 distance = max_distance
305 end if
306 distance = sign(distance, weighted_sign)
307
308 end function tri_mesh_aabb_tree
309
310 ! ========================================================================== !
311 ! Element distance functions
312 ! ========================================================================== !
313
315 subroutine element_distance(element, p, distance, weighted_sign)
316 class(*), intent(in) :: element
317 real(kind=dp), dimension(3), intent(in) :: p
318 real(kind=dp), intent(out) :: distance
319 real(kind=dp), intent(out), optional :: weighted_sign
320
321 select type (element)
322 type is (tri_t)
323 call element_distance_triangle(element, p, distance, weighted_sign)
324
325 class default
326 print *, "Error: Element type not supported."
327 stop
328 end select
329 end subroutine element_distance
330
331 ! -------------------------------------------------------------------------- !
332 ! Specific element distance functions
333
351 subroutine element_distance_triangle(triangle, p, distance, weighted_sign)
352 type(tri_t), intent(in) :: triangle
353 real(kind=dp), dimension(3), intent(in) :: p
354
355 real(kind=dp), intent(out) :: distance
356 real(kind=dp), intent(out), optional :: weighted_sign
357
358 real(kind=dp), dimension(3) :: v1, v2, v3
359 real(kind=dp), dimension(3) :: normal
360 real(kind=dp) :: normal_length
361 real(kind=dp) :: b1, b2, b3
362
363 real(kind=dp), dimension(3) :: projection
364 real(kind=dp), dimension(3) :: edge
365 real(kind=dp) :: tol = 1e-10_dp
366
367 real(kind=dp) :: face_distance
368 real(kind=dp) :: t
369
370 ! Get vertices and the normal vector
371 v1 = triangle%pts(1)%p%x
372 v2 = triangle%pts(2)%p%x
373 v3 = triangle%pts(3)%p%x
374
375 normal = cross(v2 - v1, v3 - v1)
376 normal_length = norm2(normal)
377
378 if (normal_length .lt. tol) then
379 distance = huge(1.0_dp)
380 weighted_sign = 0.0_dp
381 return
382 end if
383 normal = normal / normal_length
384
385 ! Compute Barycentric coordinates to determine if the point is inside the
386 ! triangular prism, of along an edge or by a face.
387 face_distance = dot_product(p - v1, normal)
388
389 projection = p - normal * face_distance
390 b1 = dot_product(normal, cross(v2 - v1, projection - v1)) / normal_length
391 b2 = dot_product(normal, cross(v3 - v2, projection - v2)) / normal_length
392 b3 = dot_product(normal, cross(v1 - v3, projection - v3)) / normal_length
393
394 if (b1 .le. tol) then
395 edge = v2 - v1
396 t = dot_product(p - v1, edge) / norm2(edge)
397 t = max(0.0_dp, min(1.0_dp, t))
398
399 projection = v1 + t * edge
400 else if (b2 .le. tol) then
401 edge = v3 - v2
402 t = dot_product(p - v2, edge) / norm2(edge)
403 t = max(0.0_dp, min(1.0_dp, t))
404
405 projection = v2 + t * edge
406 else if (b3 .le. tol) then
407 edge = v1 - v3
408 t = dot_product(p - v3, edge) / norm2(edge)
409 t = max(0.0_dp, min(1.0_dp, t))
410
411 projection = v3 + t * edge
412 end if
413
414 distance = norm2(projection - p)
415 if (present(weighted_sign)) then
416 weighted_sign = face_distance / distance
417 else
418 distance = sign(distance, face_distance)
419 end if
420
421 end subroutine element_distance_triangle
422
427 pure function cross(a, b) result(c)
428 real(kind=dp), dimension(3), intent(in) :: a
429 real(kind=dp), dimension(3), intent(in) :: b
430 real(kind=dp), dimension(3) :: c
431
432 c(1) = a(2) * b(3) - a(3) * b(2)
433 c(2) = a(3) * b(1) - a(1) * b(3)
434 c(3) = a(1) * b(2) - a(2) * b(1)
435
436 end function cross
437
438end 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
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.
Integer based stack.
Definition stack.f90:63
Triangular element.
Definition tri.f90:61
#define max(a, b)
Definition tensor.cu:40