Neko  0.8.1
A portable framework for high-order spectral element flow simulations
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
40 
41  implicit none
42 
43 contains
44 
55  subroutine signed_distance_field(field_data, object, max_distance)
56  use utils, only: neko_error
57  implicit none
58 
59  type(field_t), intent(inout) :: field_data
60  class(*), intent(in) :: object
61  real(kind=dp), intent(in), optional :: max_distance
62 
63  real(kind=dp) :: max_dist
64 
65  if (present(max_distance)) then
66  max_dist = max_distance
67  else
68  max_dist = huge(0.0_dp)
69  end if
70 
71  select type(object)
72  type is (tri_mesh_t)
73  call signed_distance_field_tri_mesh(field_data, object, max_dist)
74 
75  class default
76  call neko_error("signed_distance_field: Object type not supported.")
77  end select
78 
79  end subroutine signed_distance_field
80 
91  subroutine signed_distance_field_tri_mesh(field_data, mesh, max_distance)
92  use utils, only: neko_error
93  implicit none
94 
95  type(field_t), intent(inout) :: field_data
96  type(tri_mesh_t), intent(in) :: mesh
97  real(kind=dp), intent(in) :: max_distance
98 
99  integer :: total_size
100  integer :: id
101  type(aabb_tree_t) :: search_tree
102  real(kind=dp), dimension(3) :: p
103  real(kind=dp) :: distance
104 
105  ! Zero the field
106  field_data%x = 0.0_dp
107  total_size = field_data%dof%size()
108 
109  call search_tree%init(mesh%nelv)
110  call search_tree%build(mesh%el)
111 
112  if (search_tree%get_size() .ne. mesh%nelv) then
113  call neko_error("signed_distance_field_tri_mesh: &
114  & Error building the search tree.")
115  end if
116 
117  do id = 1, total_size
118  p(1) = field_data%dof%x(id, 1, 1, 1)
119  p(2) = field_data%dof%y(id, 1, 1, 1)
120  p(3) = field_data%dof%z(id, 1, 1, 1)
121 
122  distance = tri_mesh_aabb_tree(search_tree, mesh%el, p, max_distance)
123 
124  field_data%x(id, 1, 1, 1) = real(distance, kind=rp)
125  end do
126 
127  end subroutine signed_distance_field_tri_mesh
128 
139  function tri_mesh_brute_force(mesh, p, max_distance) result(distance)
140  use tri, only: tri_t
141  use point, only: point_t
142  use num_types, only: dp
143 
144  implicit none
145 
146  type(tri_mesh_t), intent(in) :: mesh
147  real(kind=dp), intent(in) :: p(3)
148  real(kind=dp), intent(in) :: max_distance
149 
150  integer :: id
151  real(kind=dp) :: distance, weighted_sign
152  real(kind=dp) :: cd, cs
153  real(kind=dp) :: tol = 1e-6_dp
154 
155  distance = 1e10_dp
156  weighted_sign = 0.0_dp
157 
158  do id = 1, mesh%nelv
159  call element_distance(mesh%el(id), p, cd, cs)
160 
161  ! Update the weighted sign, if the relative difference is small
162  if (abs(cd - distance) / distance .lt. tol) then
163  weighted_sign = weighted_sign + cs
164  else if (cd .lt. distance) then
165  weighted_sign = cs
166  end if
167 
168  distance = min(cd, distance)
169  end do
170 
171  distance = sign(min(distance, max_distance), weighted_sign)
172 
173  end function tri_mesh_brute_force
174 
188  function tri_mesh_aabb_tree(tree, object_list, p, max_distance) result(distance)
189  use aabb, only: aabb_t
191  use stack, only: stack_i4_t
192  implicit none
193 
194  class(aabb_tree_t), intent(in) :: tree
195  class(tri_t), dimension(:), intent(in) :: object_list
196  real(kind=dp), dimension(3), intent(in) :: p
197  real(kind=dp), intent(in) :: max_distance
198 
199  real(kind=dp) :: distance
200  real(kind=dp) :: weighted_sign
201 
202  real(kind=dp), parameter :: tol = 1.0e-6_dp
203 
204  type(stack_i4_t) :: simple_stack
205  integer :: current_index
206 
207  type(aabb_node_t) :: current_node
208  type(aabb_t) :: current_aabb
209  integer :: current_object_index
210  real(kind=dp) :: current_distance
211  real(kind=dp) :: current_sign
212 
213  type(aabb_node_t) :: left_node
214  type(aabb_node_t) :: right_node
215 
216  type(aabb_t) :: root_box
217  type(aabb_t) :: search_box
218 
219  integer :: root_index, left_index, right_index
220  real(kind=dp) :: random_value
221 
222  ! Initialize the stack and the search box
223  call simple_stack%init(size(object_list) * 2)
224  call search_box%init(p - max_distance, p + max_distance)
225 
226  ! Check if the root node overlaps the search box, if it does, push it to
227  ! the stack and update the search box to a randomly selected object.
228  root_index = tree%get_root_index()
229  root_box = tree%get_aabb(root_index)
230 
231  if (.not. root_box%overlaps(search_box)) then
232  distance = max_distance
233  weighted_sign = 1.0_dp
234  return
235  end if
236 
237  ! Grab a random object and compute the distance to it
238  call random_number(random_value)
239  current_object_index = floor(random_value * size(object_list) + 1)
240  call element_distance(object_list(current_object_index), p, distance, weighted_sign)
241  distance = distance + object_list(current_object_index)%diameter()
242 
243  ! Update the search box to the new distance and push the root node
244  call search_box%init(p - distance, p + distance)
245  call simple_stack%push(root_index)
246 
247  ! Traverse the tree and compute the signed distance to the elements
248  do while (.not. simple_stack%is_empty())
249  current_index = simple_stack%pop()
250  if (current_index .eq. aabb_null_node) cycle
251 
252  current_node = tree%get_node(current_index)
253  current_aabb = current_node%get_aabb()
254 
255  if (current_node%is_leaf()) then
256  if (distance .lt. current_node%min_distance(p)) then
257  cycle
258  end if
259 
260  current_object_index = current_node%get_object_index()
261  call element_distance(object_list(current_object_index), p, current_distance, current_sign)
262 
263  ! Update the weighted sign, if the relative difference is small
264  if (abs(current_distance - distance) / distance .lt. tol) then
265  weighted_sign = weighted_sign + current_sign
266  else if (current_distance .lt. distance) then
267  weighted_sign = current_sign
268  end if
269 
270  distance = min(distance, current_distance)
271 
272  ! Update the search box to the new distance
273  if (distance .gt. current_aabb%get_diameter()) then
274  call search_box%init(p - distance, p + distance)
275  end if
276  else
277 
278  left_index = tree%get_left_index(current_index)
279  if (left_index .ne. aabb_null_node) then
280  left_node = tree%get_left_node(current_index)
281  if (left_node%aabb%overlaps(search_box)) then
282  call simple_stack%push(left_index)
283  end if
284  end if
285 
286  right_index = tree%get_right_index(current_index)
287  if (right_index .ne. aabb_null_node) then
288  right_node = tree%get_right_node(current_index)
289  if (right_node%aabb%overlaps(search_box)) then
290  call simple_stack%push(right_index)
291  end if
292  end if
293  end if
294  end do
295 
296  if (distance .gt. max_distance) then
297  distance = max_distance
298  end if
299  distance = sign(distance, weighted_sign)
300 
301  end function tri_mesh_aabb_tree
302 
303  ! ========================================================================== !
304  ! Element distance functions
305  ! ========================================================================== !
306 
308  subroutine element_distance(element, p, distance, weighted_sign)
309  class(*), intent(in) :: element
310  real(kind=dp), dimension(3), intent(in) :: p
311  real(kind=dp), intent(out) :: distance
312  real(kind=dp), intent(out), optional :: weighted_sign
313 
314  select type(element)
315  type is (tri_t)
316  call element_distance_triangle(element, p, distance, weighted_sign)
317 
318  class default
319  print *, "Error: Element type not supported."
320  stop
321  end select
322  end subroutine element_distance
323 
324  ! -------------------------------------------------------------------------- !
325  ! Specific element distance functions
326 
344  subroutine element_distance_triangle(triangle, p, distance, weighted_sign)
345  type(tri_t), intent(in) :: triangle
346  real(kind=dp), dimension(3), intent(in) :: p
347 
348  real(kind=dp), intent(out) :: distance
349  real(kind=dp), intent(out), optional :: weighted_sign
350 
351  real(kind=dp), dimension(3) :: v1, v2, v3
352  real(kind=dp), dimension(3) :: normal
353  real(kind=dp) :: normal_length
354  real(kind=dp) :: b1, b2, b3
355 
356  real(kind=dp), dimension(3) :: projection
357  real(kind=dp), dimension(3) :: edge
358  real(kind=dp) :: tol = 1e-10_dp
359 
360  real(kind=dp) :: face_distance
361  real(kind=dp) :: t
362 
363  ! Get vertices and the normal vector
364  v1 = triangle%pts(1)%p%x
365  v2 = triangle%pts(2)%p%x
366  v3 = triangle%pts(3)%p%x
367 
368  normal = cross(v2 - v1, v3 - v1)
369  normal_length = norm2(normal)
370 
371  if (normal_length .lt. tol) then
372  distance = huge(1.0_dp)
373  weighted_sign = 0.0_dp
374  return
375  end if
376  normal = normal / normal_length
377 
378  ! Compute Barycentric coordinates to determine if the point is inside the
379  ! triangular prism, of along an edge or by a face.
380  face_distance = dot_product(p - v1, normal)
381 
382  projection = p - normal * face_distance
383  b1 = dot_product(normal, cross(v2 - v1, projection - v1)) / normal_length
384  b2 = dot_product(normal, cross(v3 - v2, projection - v2)) / normal_length
385  b3 = dot_product(normal, cross(v1 - v3, projection - v3)) / normal_length
386 
387  if (b1 .le. tol) then
388  edge = v2 - v1
389  t = dot_product(p - v1, edge) / norm2(edge)
390  t = max(0.0_dp, min(1.0_dp, t))
391 
392  projection = v1 + t * edge
393  else if (b2 .le. tol) then
394  edge = v3 - v2
395  t = dot_product(p - v2, edge) / norm2(edge)
396  t = max(0.0_dp, min(1.0_dp, t))
397 
398  projection = v2 + t * edge
399  else if (b3 .le. tol) then
400  edge = v1 - v3
401  t = dot_product(p - v3, edge) / norm2(edge)
402  t = max(0.0_dp, min(1.0_dp, t))
403 
404  projection = v3 + t * edge
405  end if
406 
407  distance = norm2(projection - p)
408  if (present(weighted_sign)) then
409  weighted_sign = face_distance / distance
410  else
411  distance = sign(distance, face_distance)
412  end if
413 
414  end subroutine element_distance_triangle
415 
420  pure function cross(a, b) result(c)
421  real(kind=dp), dimension(3), intent(in) :: a
422  real(kind=dp), dimension(3), intent(in) :: b
423  real(kind=dp), dimension(3) :: c
424 
425  c(1) = a(2) * b(3) - a(3) * b(2)
426  c(2) = a(3) * b(1) - a(1) * b(3)
427  c(3) = a(1) * b(2) - a(2) * b(1)
428 
429  end function cross
430 
431 end module signed_distance
double real
Definition: device_config.h:12
Axis Aligned Bounding Box (aabb) Tree data structure.
Definition: aabb_tree.f90:70
integer, parameter, public aabb_null_node
Definition: aabb_tree.f90:78
Axis Aligned Bounding Box (aabb) implementation in Fortran.
Definition: aabb.f90:71
Defines a field.
Definition: field.f90:34
Defines a mesh.
Definition: mesh.f90:34
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.
Definition: projection.f90:63
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.
subroutine signed_distance_field(field_data, object, max_distance)
Signed distance field.
real(kind=dp) function tri_mesh_aabb_tree(tree, object_list, p, max_distance)
Signed distance function using an AABB tree.
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
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:85
Axis Aligned Bounding Box (aabb) Tree.
Definition: aabb_tree.f90:125
A point in with coordinates .
Definition: point.f90:43
Integer based stack.
Definition: stack.f90:63
Triangular element.
Definition: tri.f90:60
#define max(a, b)
Definition: tensor.cu:40