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, only : aabb_t
43 use utils, only : neko_error, neko_warning
44 use math, only : rzero
45 use stack, only : stack_i4_t
46
47 implicit none
48 private
49 public :: signed_distance_field
50
51contains
52
63 subroutine signed_distance_field(field_data, object, max_distance)
64 type(field_t), intent(inout) :: field_data
65 class(*), intent(in) :: object
66 real(kind=dp), intent(in), optional :: max_distance
67
68 real(kind=dp) :: max_dist
69
70 if (present(max_distance)) then
71 max_dist = max_distance
72 else
73 max_dist = huge(0.0_dp)
74 end if
75
76 select type (object)
77 type is (tri_mesh_t)
78 call signed_distance_field_tri_mesh(field_data, object, max_dist)
79
80 class default
81 call neko_error("signed_distance_field: Object type not supported.")
82 end select
83
84 end subroutine signed_distance_field
85
96 subroutine signed_distance_field_tri_mesh(field_data, mesh, max_distance)
97 type(field_t), intent(inout) :: field_data
98 type(tri_mesh_t), intent(in) :: mesh
99 real(kind=dp), intent(in) :: max_distance
100
101 integer :: total_size
102 integer :: id
103 type(aabb_tree_t) :: search_tree
104 real(kind=dp), dimension(3) :: p
105 real(kind=dp) :: distance
106
107 ! Zero the field
108 call rzero(field_data%x, field_data%size())
109 total_size = field_data%dof%size()
110
111 call search_tree%init(mesh%nelv)
112 call search_tree%build(mesh%el)
113
114 if (search_tree%get_size() .ne. mesh%nelv) then
115 call neko_error("signed_distance_field_tri_mesh: " // &
116 "Error building the search tree.")
117 end if
118
119 !$omp parallel do private(id, p, distance)
120 do id = 1, total_size
121 p(1) = field_data%dof%x(id, 1, 1, 1)
122 p(2) = field_data%dof%y(id, 1, 1, 1)
123 p(3) = field_data%dof%z(id, 1, 1, 1)
124
125 distance = tri_mesh_aabb_tree(search_tree, mesh%el, p, max_distance)
126
127 field_data%x(id, 1, 1, 1) = real(distance, kind=rp)
128 end do
129 !$omp end parallel do
130
131 if (neko_bcknd_device .eq. 1) then
132 call neko_warning("signed_distance_field_tri_mesh:&
133 & Device version not implemented.")
134 call device_memcpy(field_data%x, field_data%x_d, field_data%size(), &
135 host_to_device, sync = .false.)
136 end if
137
138 end subroutine signed_distance_field_tri_mesh
139
152 function tri_mesh_brute_force(mesh, p, max_distance) result(distance)
153 type(tri_mesh_t), intent(in) :: mesh
154 real(kind=dp), intent(in) :: p(3)
155 real(kind=dp), intent(in) :: max_distance
156
157 integer :: id
158 real(kind=dp) :: distance, weighted_sign
159 real(kind=dp) :: cd, cs
160 real(kind=dp) :: tol = 1e-6_dp
161
162 distance = huge(0.0_dp)
163 weighted_sign = 0.0_dp
164
165 do id = 1, mesh%nelv
166 call element_distance(mesh%el(id), p, cd, cs)
167
168 ! Update the weighted sign, if the relative difference is small
169 if (abs(cd - distance) / distance .lt. tol) then
170 weighted_sign = weighted_sign + cs
171 else if (cd .lt. distance) then
172 weighted_sign = cs
173 end if
174
175 distance = min(cd, distance)
176 end do
177
178 distance = sign(min(distance, max_distance), weighted_sign)
179
180 end function tri_mesh_brute_force
181
195 function tri_mesh_aabb_tree(tree, object_list, p, max_distance) &
196 result(distance)
197 class(aabb_tree_t), intent(in) :: tree
198 class(tri_t), contiguous, dimension(:), intent(in) :: object_list
199 real(kind=dp), dimension(3), intent(in) :: p
200 real(kind=dp), intent(in) :: max_distance
201
202 real(kind=dp) :: distance
203 real(kind=dp) :: weighted_sign
204
205 real(kind=dp), parameter :: tol = 1.0e-6_dp
206
207 type(stack_i4_t) :: simple_stack
208 integer :: current_index
209
210 type(aabb_node_t) :: current_node
211 type(aabb_t) :: current_aabb
212 integer :: current_object_index
213 real(kind=dp) :: current_distance
214 real(kind=dp) :: current_sign
215
216 type(aabb_node_t) :: left_node
217 type(aabb_node_t) :: right_node
218
219 type(aabb_t) :: root_box
220 type(aabb_t) :: search_box
221
222 integer :: root_index, left_index, right_index
223 real(kind=dp) :: random_value
224
225 ! Initialize the stack and the search box
226 call simple_stack%init(size(object_list) * 2)
227 call search_box%init(p - max_distance, p + max_distance)
228
229 ! Check if the root node overlaps the search box, if it does, push it to
230 ! the stack and update the search box to a randomly selected object.
231 root_index = tree%get_root_index()
232 root_box = tree%get_aabb(root_index)
233
234 if (.not. root_box%overlaps(search_box)) then
235 distance = max_distance
236 weighted_sign = 1.0_dp
237 return
238 end if
239
240 ! Grab a random object and compute the distance to it
241 call random_number(random_value)
242 current_object_index = floor(random_value * size(object_list) + 1)
243 call element_distance(object_list(current_object_index), p, distance, &
244 weighted_sign)
245 distance = distance + object_list(current_object_index)%diameter()
246
247 ! Update the search box to the new distance and push the root node
248 call search_box%init(p - distance, p + distance)
249 call simple_stack%push(root_index)
250
251 ! Traverse the tree and compute the signed distance to the elements
252 do while (.not. simple_stack%is_empty())
253 current_index = simple_stack%pop()
254 if (current_index .eq. aabb_null_node) cycle
255
256 current_node = tree%get_node(current_index)
257 current_aabb = current_node%get_aabb()
258
259 if (current_node%is_leaf()) then
260 if (distance .lt. current_node%min_distance(p)) then
261 cycle
262 end if
263
264 current_object_index = current_node%get_object_index()
265 call element_distance(object_list(current_object_index), p, &
266 current_distance, current_sign)
267
268 ! Update the weighted sign, if the relative difference is small
269 if (abs(current_distance - distance) / distance .lt. tol) then
270 weighted_sign = weighted_sign + current_sign
271 else if (current_distance .lt. distance) then
272 weighted_sign = current_sign
273 end if
274
275 distance = min(distance, current_distance)
276
277 ! Update the search box to the new distance
278 if (distance .gt. current_aabb%get_diameter()) then
279 call search_box%init(p - distance, p + distance)
280 end if
281 else
282
283 left_index = tree%get_left_index(current_index)
284 if (left_index .ne. aabb_null_node) then
285 left_node = tree%get_left_node(current_index)
286 if (left_node%aabb%overlaps(search_box)) then
287 call simple_stack%push(left_index)
288 end if
289 end if
290
291 right_index = tree%get_right_index(current_index)
292 if (right_index .ne. aabb_null_node) then
293 right_node = tree%get_right_node(current_index)
294 if (right_node%aabb%overlaps(search_box)) then
295 call simple_stack%push(right_index)
296 end if
297 end if
298 end if
299 end do
300
301 if (distance .gt. max_distance) then
302 distance = max_distance
303 end if
304 distance = sign(distance, weighted_sign)
305
306 end function tri_mesh_aabb_tree
307
308 ! ========================================================================== !
309 ! Element distance functions
310 ! ========================================================================== !
311
313 subroutine element_distance(element, p, distance, weighted_sign)
314 class(*), intent(in) :: element
315 real(kind=dp), dimension(3), intent(in) :: p
316 real(kind=dp), intent(out) :: distance
317 real(kind=dp), intent(out), optional :: weighted_sign
318
319 select type (element)
320 type is (tri_t)
321 call element_distance_triangle(element, p, distance, weighted_sign)
322
323 class default
324 print *, "Error: Element type not supported."
325 stop
326 end select
327 end subroutine element_distance
328
329 ! -------------------------------------------------------------------------- !
330 ! Specific element distance functions
331
349 subroutine element_distance_triangle(triangle, p, distance, weighted_sign)
350 type(tri_t), intent(in) :: triangle
351 real(kind=dp), dimension(3), intent(in) :: p
352
353 real(kind=dp), intent(out) :: distance
354 real(kind=dp), intent(out), optional :: weighted_sign
355
356 real(kind=dp), dimension(3) :: v1, v2, v3
357 real(kind=dp), dimension(3) :: normal
358 real(kind=dp) :: normal_length
359 real(kind=dp) :: b1, b2, b3
360
361 real(kind=dp), dimension(3) :: projection
362 real(kind=dp), dimension(3) :: edge
363 real(kind=dp) :: tol = 1e-10_dp
364
365 real(kind=dp) :: face_distance
366 real(kind=dp) :: t
367
368 ! Get vertices and the normal vector
369 v1 = triangle%pts(1)%p%x
370 v2 = triangle%pts(2)%p%x
371 v3 = triangle%pts(3)%p%x
372
373 normal = cross(v2 - v1, v3 - v1)
374 normal_length = norm2(normal)
375
376 if (normal_length .lt. tol) then
377 distance = huge(1.0_dp)
378 weighted_sign = 0.0_dp
379 return
380 end if
381 normal = normal / normal_length
382
383 ! Compute Barycentric coordinates to determine if the point is inside the
384 ! triangular prism, of along an edge or by a face.
385 face_distance = dot_product(p - v1, normal)
386
387 projection = p - normal * face_distance
388 b1 = dot_product(normal, cross(v2 - v1, projection - v1)) / normal_length
389 b2 = dot_product(normal, cross(v3 - v2, projection - v2)) / normal_length
390 b3 = dot_product(normal, cross(v1 - v3, projection - v3)) / normal_length
391
392 if (b1 .le. tol) then
393 edge = v2 - v1
394 t = dot_product(p - v1, edge) / norm2(edge)
395 t = max(0.0_dp, min(1.0_dp, t))
396
397 projection = v1 + t * edge
398 else if (b2 .le. tol) then
399 edge = v3 - v2
400 t = dot_product(p - v2, edge) / norm2(edge)
401 t = max(0.0_dp, min(1.0_dp, t))
402
403 projection = v2 + t * edge
404 else if (b3 .le. tol) then
405 edge = v1 - v3
406 t = dot_product(p - v3, edge) / norm2(edge)
407 t = max(0.0_dp, min(1.0_dp, t))
408
409 projection = v3 + t * edge
410 end if
411
412 distance = norm2(projection - p)
413 if (present(weighted_sign)) then
414 weighted_sign = face_distance / distance
415 else
416 distance = sign(distance, face_distance)
417 end if
418
419 end subroutine element_distance_triangle
420
425 pure function cross(a, b) result(c)
426 real(kind=dp), dimension(3), intent(in) :: a
427 real(kind=dp), dimension(3), intent(in) :: b
428 real(kind=dp), dimension(3) :: c
429
430 c(1) = a(2) * b(3) - a(3) * b(2)
431 c(2) = a(3) * b(1) - a(1) * b(3)
432 c(3) = a(1) * b(2) - a(2) * b(1)
433
434 end function cross
435
436end module signed_distance
double real
Copy data between host and device (or device and device)
Definition device.F90:72
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:48
Defines a field.
Definition field.f90:34
Definition math.f90:60
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:233
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:49
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:398
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:77
Triangular element.
Definition tri.f90:61
#define max(a, b)
Definition tensor.cu:40