94 type(
field_t),
intent(inout) :: field_data
96 real(kind=
dp),
intent(in) :: max_distance
101 real(kind=
dp),
dimension(3) :: p
102 real(kind=
dp) :: distance
105 field_data%x = 0.0_dp
106 total_size = field_data%dof%size()
108 call search_tree%init(
mesh%nelv)
109 call search_tree%build(
mesh%el)
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.")
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)
123 field_data%x(id, 1, 1, 1) =
real(distance, kind=
rp)
128 & Device version not implemented.")
129 call device_memcpy(field_data%x, field_data%x_d, field_data%size(), &
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
206 real(kind=dp) :: distance
207 real(kind=dp) :: weighted_sign
209 real(kind=dp),
parameter :: tol = 1.0e-6_dp
212 integer :: current_index
215 type(
aabb_t) :: current_aabb
216 integer :: current_object_index
217 real(kind=dp) :: current_distance
218 real(kind=dp) :: current_sign
224 type(
aabb_t) :: search_box
226 integer :: root_index, left_index, right_index
227 real(kind=dp) :: random_value
230 call simple_stack%init(
size(object_list) * 2)
231 call search_box%init(p - max_distance, p + max_distance)
235 root_index = tree%get_root_index()
236 root_box = tree%get_aabb(root_index)
238 if (.not. root_box%overlaps(search_box))
then
239 distance = max_distance
240 weighted_sign = 1.0_dp
245 call random_number(random_value)
246 current_object_index = floor(random_value *
size(object_list) + 1)
249 distance = distance + object_list(current_object_index)%diameter()
252 call search_box%init(p - distance, p + distance)
253 call simple_stack%push(root_index)
256 do while (.not. simple_stack%is_empty())
257 current_index = simple_stack%pop()
260 current_node = tree%get_node(current_index)
261 current_aabb = current_node%get_aabb()
263 if (current_node%is_leaf())
then
264 if (distance .lt. current_node%min_distance(p))
then
268 current_object_index = current_node%get_object_index()
270 current_distance, current_sign)
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
279 distance = min(distance, current_distance)
282 if (distance .gt. current_aabb%get_diameter())
then
283 call search_box%init(p - distance, p + distance)
287 left_index = tree%get_left_index(current_index)
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)
295 right_index = tree%get_right_index(current_index)
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)
305 if (distance .gt. max_distance)
then
306 distance = max_distance
308 distance = sign(distance, weighted_sign)
354 type(tri_t),
intent(in) :: triangle
355 real(kind=dp),
dimension(3),
intent(in) :: p
357 real(kind=dp),
intent(out) :: distance
358 real(kind=dp),
intent(out),
optional :: weighted_sign
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
366 real(kind=dp),
dimension(3) :: edge
367 real(kind=dp) :: tol = 1e-10_dp
369 real(kind=dp) :: face_distance
373 v1 = triangle%pts(1)%p%x
374 v2 = triangle%pts(2)%p%x
375 v3 = triangle%pts(3)%p%x
377 normal =
cross(v2 - v1, v3 - v1)
378 normal_length = norm2(normal)
380 if (normal_length .lt. tol)
then
381 distance = huge(1.0_dp)
382 weighted_sign = 0.0_dp
385 normal = normal / normal_length
389 face_distance = dot_product(p - v1, normal)
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
396 if (b1 .le. tol)
then
398 t = dot_product(p - v1, edge) / norm2(edge)
399 t =
max(0.0_dp, min(1.0_dp, t))
402 else if (b2 .le. tol)
then
404 t = dot_product(p - v2, edge) / norm2(edge)
405 t =
max(0.0_dp, min(1.0_dp, t))
408 else if (b3 .le. tol)
then
410 t = dot_product(p - v3, edge) / norm2(edge)
411 t =
max(0.0_dp, min(1.0_dp, t))
417 if (
present(weighted_sign))
then
418 weighted_sign = face_distance / distance
420 distance = sign(distance, face_distance)