Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
intersection_detector.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
36 use aabb_tree, only : aabb_tree_t
37 use mesh, only : mesh_t
38 use hex, only : hex_t
39 use point, only : point_t
40 use utils, only : neko_error
41 use stack, only : stack_i4_t
42 implicit none
43 private
44
45 type, public :: intersect_detector_t
46 type(mesh_t), private, pointer :: msh => null()
47 type(aabb_tree_t), private :: search_tree
48 contains
49 procedure, pass(this) :: init => intersect_detector_init
50 procedure, pass(this) :: free => intersect_detector_free
51 procedure, pass(this) :: overlap => intersect_detector_overlap
53
54contains
55
57 subroutine intersect_detector_init(this, msh, padding)
58 class(intersect_detector_t), intent(inout) :: this
59 type(mesh_t), target, intent(in) :: msh
60 real(kind=dp), intent(in), optional :: padding
61 type(hex_t), allocatable :: elements(:)
62 integer :: i
63
64 call this%free()
65 this%msh => msh
66
67 call this%search_tree%init(msh%nelv)
68
70 allocate(elements(msh%nelv))
71 do i = 1, msh%nelv
72 select type(el => msh%elements(i)%e)
73 type is (hex_t)
74 elements(i) = el
75 class default
76 call neko_error('Unsupported element type')
77 end select
78 end do
79
80 if (present(padding)) then
81 call this%search_tree%build(elements, padding)
82 else
83 call this%search_tree%build(elements)
84 end if
85
86 deallocate(elements)
87
88 if (this%search_tree%get_size() .ne. (msh%nelv)) then
89 call neko_error("Error building the search tree.")
90 end if
91
92 end subroutine intersect_detector_init
93
95 subroutine intersect_detector_free(this)
96 class(intersect_detector_t), intent(inout) :: this
97
98 if (associated(this%msh)) then
99 nullify(this%msh)
100 end if
101
103
104 end subroutine intersect_detector_free
105
107 subroutine intersect_detector_overlap(this, p, overlaps)
108 class(intersect_detector_t), intent(inout) :: this
109 type(point_t), intent(in) :: p
110 type(stack_i4_t), intent(inout) :: overlaps
111
112 call this%search_tree%query_overlaps(p, -1, overlaps)
113
114 end subroutine intersect_detector_overlap
115
116end module intersection_detector
Axis Aligned Bounding Box (aabb) Tree data structure.
Definition aabb_tree.f90:70
Defines a hexahedron element.
Definition hex.f90:34
Implements a mesh intersection detector.
subroutine intersect_detector_free(this)
Destroy an intersector detector.
subroutine intersect_detector_overlap(this, p, overlaps)
Computes the overlap between elements and a given point p.
subroutine intersect_detector_init(this, msh, padding)
Initialise an intersector detector for a given mesh msh.
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public dp
Definition num_types.f90:9
Implements a point.
Definition point.f90:35
Implements a dynamic stack ADT.
Definition stack.f90:49
Utilities.
Definition utils.f90:35
Axis Aligned Bounding Box (aabb) Tree.
Hexahedron element.
Definition hex.f90:63
A point in with coordinates .
Definition point.f90:43
Integer based stack.
Definition stack.f90:77