Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
aabb.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!
33! ============================================================================ !
34! Original C++ Implementation from:
35! https://github.com/JamesRandall/SimpleVoxelEngine/blob/master/voxelEngine/include/AABB.h
36!
37! Translated to Fortran by:
38! @author Tim Felle Olsen
39! @date 9 Feb 2024
40!
41! C++ Code License:
42! The MIT License (MIT)
43!
44! Copyright (c) 2017 James Randall
45!
46! Permission is hereby granted, free of charge, to any person obtaining a copy of
47! this software and associated documentation files (the "Software"), to deal in
48! the Software without restriction, including without limitation the rights to
49! use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of
50! the Software, and to permit persons to whom the Software is furnished to do so,
51! subject to the following conditions:
52!
53! The above copyright notice and this permission notice shall be included in all
54! copies or substantial portions of the Software.
55!
56! THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
57! IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS
58! FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR
59! COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER
60! IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN
61! CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
62! ============================================================================ !
63
71module aabb
72 use num_types, only: dp
73 use element, only: element_t
74 use point, only: point_t
75 use tri, only: tri_t
76 use quad, only: quad_t
77 use tet, only: tet_t
78 use hex, only: hex_t
79 use mesh, only: mesh_t
80 use tri_mesh, only: tri_mesh_t
81 use tet_mesh, only: tet_mesh_t
82 use utils, only: neko_error
83
84 implicit none
85 private
87
88 ! ========================================================================== !
89 ! Public interface for free functions
90 ! ========================================================================== !
91
93 interface merge
94 module procedure merge_aabb
95 end interface merge
96
98 interface intersection
99 module procedure intersection_aabb
100 end interface intersection
101
107 type :: aabb_t
108 private
109
110 logical :: initialized = .false.
111 real(kind=dp) :: box_min(3) = huge(0.0_dp)
112 real(kind=dp) :: box_max(3) = -huge(0.0_dp)
113 real(kind=dp) :: center(3) = 0.0_dp
114 real(kind=dp) :: diameter = huge(0.0_dp)
115 real(kind=dp) :: surface_area = 0.0_dp
116
117 contains
118
119 ! Initializers
120 procedure, pass(this), public :: init => aabb_init
121
122 ! Getters
123 procedure, pass(this), public :: get_min => aabb_get_min
124 procedure, pass(this), public :: get_max => aabb_get_max
125 procedure, pass(this), public :: get_width => aabb_get_width
126 procedure, pass(this), public :: get_height => aabb_get_height
127 procedure, pass(this), public :: get_depth => aabb_get_depth
128 procedure, pass(this), public :: get_diameter => aabb_get_diameter
129 procedure, pass(this), public :: get_surface_area => aabb_get_surface_area
130 procedure, pass(this), public :: get_center => aabb_get_center
131 procedure, pass(this), public :: get_diagonal => aabb_get_diagonal
132
133 procedure, pass(this), public :: add_padding
134
135 ! Comparison operators
136 generic :: operator(.lt.) => less
137 generic :: operator(.gt.) => greater
138
140 procedure, pass(this), public :: overlaps => aabb_overlaps
142 procedure, pass(this), public :: contains => aabb_contains_other
144 procedure, pass(this), public :: contains_point => aabb_contains_point
145
146 ! Private comparison operators
147 procedure, pass(this) :: less => aabb_less
148 procedure, pass(this) :: greater => aabb_greater
149
150 ! Private operations
151 procedure, pass(this), private :: calculate_surface_area
152
153 end type aabb_t
154
155contains
156
157 ! ========================================================================== !
158 ! Constructors
159 ! ========================================================================== !
160
180 function get_aabb(object, padding) result(box)
181 use utils, only: neko_error
182 implicit none
183
184 class(*), intent(in) :: object
185 real(kind=dp), intent(in), optional :: padding
186 type(aabb_t) :: box
187
188 select type(object)
189
190 type is (aabb_t)
191 call box%init(object%get_min(), object%get_max())
192 type is (tri_t)
193 box = get_aabb_element(object)
194 type is (hex_t)
195 box = get_aabb_element(object)
196 type is (tet_t)
197 box = get_aabb_element(object)
198 type is (quad_t)
199 box = get_aabb_element(object)
200
201 type is (mesh_t)
202 box = get_aabb_mesh(object)
203 type is (tri_mesh_t)
204 box = get_aabb_tri_mesh(object)
205 type is (tet_mesh_t)
206 box = get_aabb_tet_mesh(object)
207
208 class default
209 call neko_error("get_aabb: Unsupported object type")
210 end select
211
212 if (present(padding)) then
213 call box%add_padding(padding)
214 end if
215
216 end function get_aabb
217
223 subroutine add_padding(this, padding)
224 class(aabb_t), intent(inout) :: this
225 real(kind=dp), intent(in) :: padding
226 real(kind=dp) :: box_min(3), box_max(3)
227
228 box_min = this%box_min - padding * (this%diameter)
229 box_max = this%box_max + padding * (this%diameter)
230
231 call this%init(box_min, box_max)
232 end subroutine add_padding
233
243 function get_aabb_element(object) result(box)
244 class(element_t), intent(in) :: object
245 type(aabb_t) :: box
246
247 integer :: i
248 type(point_t), pointer :: pi
249 real(kind=dp) :: box_min(3), box_max(3)
250
251 box_min = huge(0.0_dp); box_max = -huge(0.0_dp)
252
253 do i = 1, object%n_points()
254 pi => object%p(i)
255 box_min = min(box_min, pi%x)
256 box_max = max(box_max, pi%x)
257 end do
258
259 call box%init(box_min, box_max)
260 end function get_aabb_element
261
270 function get_aabb_mesh(object) result(box)
271 type(mesh_t), intent(in) :: object
272 type(aabb_t) :: box
273
274 integer :: i
275 type(aabb_t) :: temp_box
276
277 do i = 1, object%nelv
278 temp_box = get_aabb(object%elements(i))
279 box = merge(box, temp_box)
280 end do
281
282 end function get_aabb_mesh
283
292 function get_aabb_tri_mesh(object) result(box)
293 type(tri_mesh_t), intent(in) :: object
294 type(aabb_t) :: box
295
296 integer :: i
297 type(aabb_t) :: temp_box
298
299 do i = 1, object%nelv
300 temp_box = get_aabb(object%el(i))
301 box = merge(box, temp_box)
302 end do
303
304 end function get_aabb_tri_mesh
305
314 function get_aabb_tet_mesh(object) result(box)
315 type(tet_mesh_t), intent(in) :: object
316 type(aabb_t) :: box
317
318 integer :: i
319 type(aabb_t) :: temp_box
320
321 do i = 1, object%nelv
322 temp_box = get_aabb(object%el(i))
323 box = merge(box, temp_box)
324 end do
325
326 end function get_aabb_tet_mesh
327
328
329 ! ========================================================================== !
330 ! Initializers
331 ! ========================================================================== !
332
336 subroutine aabb_init(this, lower_left_front, upper_right_back)
337 class(aabb_t), intent(inout) :: this
338 real(kind=dp), dimension(3), intent(in) :: lower_left_front
339 real(kind=dp), dimension(3), intent(in) :: upper_right_back
340
341 this%initialized = .true.
342 this%box_min = lower_left_front
343 this%box_max = upper_right_back
344 this%center = (this%box_min + this%box_max) / 2.0_dp
345 this%diameter = norm2(this%box_max - this%box_min)
346
347 this%surface_area = this%calculate_surface_area()
348
349 end subroutine aabb_init
350
351 ! ========================================================================== !
352 ! Getters
353 ! ========================================================================== !
354
356 pure function aabb_get_min(this) result(min)
357 class(aabb_t), intent(in) :: this
358 real(kind=dp), dimension(3) :: min
359
360 min = this%box_min
361 end function aabb_get_min
362
364 pure function aabb_get_max(this) result(max)
365 class(aabb_t), intent(in) :: this
366 real(kind=dp), dimension(3) :: max
367
368 max = this%box_max
369 end function aabb_get_max
370
372 pure function aabb_get_width(this) result(width)
373 class(aabb_t), intent(in) :: this
374 real(kind=dp) :: width
375
376 width = this%box_max(1) - this%box_min(1)
377 end function aabb_get_width
378
380 pure function aabb_get_depth(this) result(depth)
381 class(aabb_t), intent(in) :: this
382 real(kind=dp) :: depth
383
384 depth = this%box_max(2) - this%box_min(2)
385 end function aabb_get_depth
386
388 pure function aabb_get_height(this) result(height)
389 class(aabb_t), intent(in) :: this
390 real(kind=dp) :: height
391
392 height = this%box_max(3) - this%box_min(3)
393 end function aabb_get_height
394
396 pure function aabb_get_diameter(this) result(diameter)
397 class(aabb_t), intent(in) :: this
398 real(kind=dp) :: diameter
399
400 diameter = this%diameter
401 end function aabb_get_diameter
402
404 pure function aabb_get_surface_area(this) result(surface_area)
405 class(aabb_t), intent(in) :: this
406 real(kind=dp) :: surface_area
407
408 surface_area = this%surface_area
409 end function aabb_get_surface_area
410
412 pure function aabb_get_center(this) result(center)
413 class(aabb_t), intent(in) :: this
414 real(kind=dp), dimension(3) :: center
415
416 center = this%center
417 end function aabb_get_center
418
420 pure function aabb_get_diagonal(this) result(diagonal)
421 class(aabb_t), intent(in) :: this
422 real(kind=dp), dimension(3) :: diagonal
423
424 diagonal = this%box_max - this%box_min
425 end function aabb_get_diagonal
426
427 ! ========================================================================== !
428 ! Operations
429 ! ========================================================================== !
430
432 function aabb_overlaps(this, other) result(is_overlapping)
433 class(aabb_t), intent(in) :: this
434 class(aabb_t), intent(in) :: other
435 logical :: is_overlapping
436
437 if (.not. this%initialized .or. .not. other%initialized) then
438 ! call neko_error("aabb_overlaps: One or both aabbs are not initialized")
439 is_overlapping = .false.
440 else
441
442 is_overlapping = all(this%box_min .le. other%box_max) .and. &
443 all(this%box_max .ge. other%box_min)
444 end if
445
446 end function aabb_overlaps
447
449 function aabb_contains_other(this, other) result(is_contained)
450 class(aabb_t), intent(in) :: this
451 class(aabb_t), intent(in) :: other
452 logical :: is_contained
453
454 ! if (.not. this%initialized .or. .not. other%initialized) then
455 ! call neko_error("aabb_contains: One or both aabbs are not initialized")
456 ! end if
457
458 is_contained = all(this%box_min .le. other%box_min) .and. &
459 all(this%box_max .ge. other%box_max)
460
461 end function aabb_contains_other
462
464 function aabb_contains_point(this, p) result(is_contained)
465 class(aabb_t), intent(in) :: this
466 real(kind=dp), dimension(3), intent(in) :: p
467 logical :: is_contained
468
469 ! if (.not. this%initialized) then
470 ! call neko_error("aabb_contains_point: One or both aabbs are not initialized")
471 ! end if
472
473 is_contained = all(p .ge. this%box_min) .and. all(p .le. this%box_max)
474 end function aabb_contains_point
475
476 ! ========================================================================== !
477 ! Binary operations
478 ! ========================================================================== !
479
481 function merge_aabb(box1, box2) result(merged)
482 class(aabb_t), intent(in) :: box1
483 class(aabb_t), intent(in) :: box2
484 type(aabb_t) :: merged
485
486 real(kind=dp), dimension(3) :: box_min, box_max
487
488 box_min = min(box1%box_min, box2%box_min)
489 box_max = max(box1%box_max, box2%box_max)
490
491 call merged%init(box_min, box_max)
492 end function merge_aabb
493
495 function intersection_aabb(box1, box2) result(intersected)
496 class(aabb_t), intent(in) :: box1
497 class(aabb_t), intent(in) :: box2
498 type(aabb_t) :: intersected
499
500 real(kind=dp), dimension(3) :: box_min, box_max
501
502 box_min = max(box1%box_min, box2%box_min)
503 box_max = min(box1%box_max, box2%box_max)
504
505 call intersected%init(box_min, box_max)
506 end function intersection_aabb
507
508 ! ========================================================================== !
509 ! Private operations
510 ! ========================================================================== !
511
513 pure function calculate_surface_area(this) result(surface_area)
514 class(aabb_t), intent(in) :: this
515 real(kind=dp) :: surface_area
516
517 surface_area = 2.0 * (&
518 & this%get_width() * this%get_height() &
519 & + this%get_width() * this%get_depth() &
520 & + this%get_height() * this%get_depth() &
521 &)
522 end function calculate_surface_area
523
524 ! ========================================================================== !
525 ! Comparison operators
526 ! ========================================================================== !
527
529 pure function aabb_less(this, other)
530 class(aabb_t), intent(in) :: this
531 class(aabb_t), intent(in) :: other
532 logical :: aabb_less
533 logical :: equal
534
535 if (.not. this%initialized .or. .not. other%initialized) then
536 aabb_less = .false.
537 return
538 end if
539
540 aabb_less = this%box_min(1) .lt. other%box_min(1)
541 equal = this%box_min(1) .le. other%box_min(1)
542
543 if (.not. aabb_less .and. equal) then
544 aabb_less = this%box_min(2) .lt. other%box_min(2)
545 equal = this%box_min(2) .le. other%box_min(2)
546 end if
547 if (.not. aabb_less .and. equal) then
548 aabb_less = this%box_min(3) .lt. other%box_min(3)
549 equal = this%box_min(3) .le. other%box_min(3)
550 end if
551 if (.not. aabb_less .and. equal) then
552 aabb_less = this%box_max(1) .lt. other%box_max(1)
553 equal = this%box_max(1) .le. other%box_max(1)
554 end if
555 if (.not. aabb_less .and. equal) then
556 aabb_less = this%box_max(2) .lt. other%box_max(2)
557 equal = this%box_max(2) .le. other%box_max(2)
558 end if
559 if (.not. aabb_less .and. equal) then
560 aabb_less = this%box_max(3) .lt. other%box_max(3)
561 end if
562
563 end function aabb_less
564
566 pure function aabb_greater(this, other)
567 class(aabb_t), intent(in) :: this
568 class(aabb_t), intent(in) :: other
569 logical :: aabb_greater
570 logical :: equal
571
572 if (.not. this%initialized .or. .not. other%initialized) then
573 aabb_greater = .false.
574 return
575 end if
576
577 aabb_greater = this%box_min(1) .gt. other%box_min(1)
578 equal = this%box_min(1) .ge. other%box_min(1)
579
580 if (.not. aabb_greater .and. equal) then
581 aabb_greater = this%box_min(2) .gt. other%box_min(2)
582 equal = this%box_min(2) .ge. other%box_min(2)
583 end if
584 if (.not. aabb_greater .and. equal) then
585 aabb_greater = this%box_min(3) .gt. other%box_min(3)
586 equal = this%box_min(3) .ge. other%box_min(3)
587 end if
588 if (.not. aabb_greater .and. equal) then
589 aabb_greater = this%box_max(1) .gt. other%box_max(1)
590 equal = this%box_max(1) .ge. other%box_max(1)
591 end if
592 if (.not. aabb_greater .and. equal) then
593 aabb_greater = this%box_max(2) .gt. other%box_max(2)
594 equal = this%box_max(2) .ge. other%box_max(2)
595 end if
596 if (.not. aabb_greater .and. equal) then
597 aabb_greater = this%box_max(3) .gt. other%box_max(3)
598 end if
599
600 end function aabb_greater
601
602end module aabb
Intersect two aabbs.
Definition aabb.f90:98
Merge two aabbs.
Definition aabb.f90:93
Axis Aligned Bounding Box (aabb) implementation in Fortran.
Definition aabb.f90:71
pure real(kind=dp) function, dimension(3) aabb_get_max(this)
Get the maximum point of the aabb.
Definition aabb.f90:365
pure real(kind=dp) function aabb_get_surface_area(this)
Get the surface area of the aabb.
Definition aabb.f90:405
pure real(kind=dp) function aabb_get_diameter(this)
Get the diameter length of the aabb.
Definition aabb.f90:397
subroutine add_padding(this, padding)
Add padding to the aabb.
Definition aabb.f90:224
logical function aabb_contains_point(this, p)
Check if this aabb contains a point.
Definition aabb.f90:465
subroutine aabb_init(this, lower_left_front, upper_right_back)
Initialize the aabb.
Definition aabb.f90:337
type(aabb_t) function get_aabb_tri_mesh(object)
Get the aabb of a triangular mesh.
Definition aabb.f90:293
logical function aabb_overlaps(this, other)
Check if two aabbs are overlapping.
Definition aabb.f90:433
type(aabb_t) function get_aabb_mesh(object)
Get the aabb of a mesh.
Definition aabb.f90:271
type(aabb_t) function intersection_aabb(box1, box2)
Get the intersection of two aabbs.
Definition aabb.f90:496
pure real(kind=dp) function aabb_get_height(this)
Get the height of the aabb. Also known as the z-axis length.
Definition aabb.f90:389
pure logical function aabb_greater(this, other)
Greater than comparison operator.
Definition aabb.f90:567
type(aabb_t) function, public get_aabb(object, padding)
Construct the aabb of a predefined object.
Definition aabb.f90:181
logical function aabb_contains_other(this, other)
Check if this aabb contains another aabb.
Definition aabb.f90:450
pure real(kind=dp) function calculate_surface_area(this)
Calculate the surface area of the aabb.
Definition aabb.f90:514
pure logical function aabb_less(this, other)
Less than comparison operator.
Definition aabb.f90:530
type(aabb_t) function get_aabb_element(object)
Get the aabb of an arbitrary element.
Definition aabb.f90:244
type(aabb_t) function get_aabb_tet_mesh(object)
Get the aabb of a tetrahedral mesh.
Definition aabb.f90:315
pure real(kind=dp) function, dimension(3) aabb_get_center(this)
Get the center of the aabb.
Definition aabb.f90:413
pure real(kind=dp) function aabb_get_width(this)
Get the width of the aabb. Also known as the x-axis length.
Definition aabb.f90:373
pure real(kind=dp) function, dimension(3) aabb_get_diagonal(this)
Get the diagonal of the aabb.
Definition aabb.f90:421
type(aabb_t) function merge_aabb(box1, box2)
Merge two aabbs.
Definition aabb.f90:482
pure real(kind=dp) function aabb_get_depth(this)
Get the depth of the aabb. Also known as the y-axis length.
Definition aabb.f90:381
pure real(kind=dp) function, dimension(3) aabb_get_min(this)
Get the minimum point of the aabb.
Definition aabb.f90:357
Defines a hexahedron element.
Definition hex.f90:34
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public dp
Definition num_types.f90:9
Implements a point.
Definition point.f90:35
Defines a quadrilateral element.
Definition quad.f90:34
Defines a tetrahedral mesh.
Definition tet_mesh.f90:35
Defines a tetrahedral element.
Definition tet.f90:34
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
Base type for an element.
Definition element.f90:44
Hexahedron element.
Definition hex.f90:63
A point in with coordinates .
Definition point.f90:43
Quadrilateral element.
Definition quad.f90:58
Tetrahedral element.
Definition tet.f90:63
Triangular element.
Definition tri.f90:60
#define max(a, b)
Definition tensor.cu:40