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(), &
151 real(kind=
dp),
intent(in) :: p(3)
152 real(kind=
dp),
intent(in) :: max_distance
155 real(kind=
dp) :: distance, weighted_sign
156 real(kind=
dp) :: cd, cs
157 real(kind=
dp) :: tol = 1e-6_dp
160 weighted_sign = 0.0_dp
166 if (abs(cd - distance) / distance .lt. tol)
then
167 weighted_sign = weighted_sign + cs
168 else if (cd .lt. distance)
then
172 distance = min(cd, distance)
175 distance = sign(min(distance, max_distance), weighted_sign)
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
204 real(kind=dp) :: distance
205 real(kind=dp) :: weighted_sign
207 real(kind=dp),
parameter :: tol = 1.0e-6_dp
210 integer :: current_index
213 type(
aabb_t) :: current_aabb
214 integer :: current_object_index
215 real(kind=dp) :: current_distance
216 real(kind=dp) :: current_sign
222 type(
aabb_t) :: search_box
224 integer :: root_index, left_index, right_index
225 real(kind=dp) :: random_value
228 call simple_stack%init(
size(object_list) * 2)
229 call search_box%init(p - max_distance, p + max_distance)
233 root_index = tree%get_root_index()
234 root_box = tree%get_aabb(root_index)
236 if (.not. root_box%overlaps(search_box))
then
237 distance = max_distance
238 weighted_sign = 1.0_dp
243 call random_number(random_value)
244 current_object_index = floor(random_value *
size(object_list) + 1)
247 distance = distance + object_list(current_object_index)%diameter()
250 call search_box%init(p - distance, p + distance)
251 call simple_stack%push(root_index)
254 do while (.not. simple_stack%is_empty())
255 current_index = simple_stack%pop()
258 current_node = tree%get_node(current_index)
259 current_aabb = current_node%get_aabb()
261 if (current_node%is_leaf())
then
262 if (distance .lt. current_node%min_distance(p))
then
266 current_object_index = current_node%get_object_index()
268 current_distance, current_sign)
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
277 distance = min(distance, current_distance)
280 if (distance .gt. current_aabb%get_diameter())
then
281 call search_box%init(p - distance, p + distance)
285 left_index = tree%get_left_index(current_index)
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)
293 right_index = tree%get_right_index(current_index)
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)
303 if (distance .gt. max_distance)
then
304 distance = max_distance
306 distance = sign(distance, weighted_sign)
352 type(
tri_t),
intent(in) :: triangle
353 real(kind=dp),
dimension(3),
intent(in) :: p
355 real(kind=dp),
intent(out) :: distance
356 real(kind=dp),
intent(out),
optional :: weighted_sign
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
364 real(kind=dp),
dimension(3) :: edge
365 real(kind=dp) :: tol = 1e-10_dp
367 real(kind=dp) :: face_distance
371 v1 = triangle%pts(1)%p%x
372 v2 = triangle%pts(2)%p%x
373 v3 = triangle%pts(3)%p%x
375 normal =
cross(v2 - v1, v3 - v1)
376 normal_length = norm2(normal)
378 if (normal_length .lt. tol)
then
379 distance = huge(1.0_dp)
380 weighted_sign = 0.0_dp
383 normal = normal / normal_length
387 face_distance = dot_product(p - v1, normal)
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
394 if (b1 .le. tol)
then
396 t = dot_product(p - v1, edge) / norm2(edge)
397 t =
max(0.0_dp, min(1.0_dp, t))
400 else if (b2 .le. tol)
then
402 t = dot_product(p - v2, edge) / norm2(edge)
403 t =
max(0.0_dp, min(1.0_dp, t))
406 else if (b3 .le. tol)
then
408 t = dot_product(p - v3, edge) / norm2(edge)
409 t =
max(0.0_dp, min(1.0_dp, t))
415 if (
present(weighted_sign))
then
416 weighted_sign = face_distance / distance
418 distance = sign(distance, face_distance)