1!! Copyright (C) 2014 M. Oliveira, J. Jornet-Somoza 2!! 3!! This program is free software; you can redistribute it and/or modify 4!! it under the terms of the GNU General Public License as published by 5!! the Free Software Foundation; either version 2, or (at your option) 6!! any later version. 7!! 8!! This program is distributed in the hope that it will be useful, 9!! but WITHOUT ANY WARRANTY; without even the implied warranty of 10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 11!! GNU General Public License for more details. 12!! 13!! You should have received a copy of the GNU General Public License 14!! along with this program; if not, write to the Free Software 15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 16!! 02110-1301, USA. 17!! 18 19#include "global.h" 20 21module box_union_oct_m 22 use global_oct_m 23 use messages_oct_m 24 use profiling_oct_m 25 use box_oct_m 26 27 implicit none 28 29 private 30 public :: & 31 box_union_t, & 32 box_union_init, & 33 box_union_end, & 34 box_union_inside_vec, & 35 box_union_inside, & 36 box_union_get_nboxes, & 37 box_union_get_center 38 39 type box_union_t 40 private 41 42 !> TODO: make this a linked list, so that boxes can be added and removed efficiently on-the-fly 43 integer :: n_boxes 44 type(box_t), pointer :: boxes(:) 45 end type box_union_t 46 47contains 48 49 !-------------------------------------------------------------- 50 subroutine box_union_init(union, n_boxes, boxes) 51 type(box_union_t), intent(out) :: union 52 integer, intent(in) :: n_boxes 53 type(box_t), intent(in) :: boxes(:) 54 55 integer :: ibox 56 57 PUSH_SUB(box_union_init) 58 59 union%n_boxes = n_boxes 60 SAFE_ALLOCATE(union%boxes(1:n_boxes)) 61 62 do ibox = 1,n_boxes 63 call box_copy(union%boxes(ibox), boxes(ibox)) 64 end do 65 66 POP_SUB(box_union_init) 67 end subroutine box_union_init 68 69 !-------------------------------------------------------------- 70 subroutine box_union_end(union) 71 type(box_union_t), intent(inout) :: union 72 73 integer :: ibox 74 75 PUSH_SUB(box_union_end) 76 77 do ibox = 1, union%n_boxes 78 call box_end(union%boxes(ibox)) 79 end do 80 SAFE_DEALLOCATE_P(union%boxes) 81 82 union%n_boxes = 0 83 84 POP_SUB(box_union_end) 85 end subroutine box_union_end 86 87 !-------------------------------------------------------------- 88 !> Checks if a vector of points are inside the box. 89 subroutine box_union_inside_vec(union, npoints, points, inside) 90 type(box_union_t), intent(in) :: union 91 integer, intent(in) :: npoints 92 FLOAT, intent(in) :: points(:, :) 93 logical, intent(out) :: inside(:) 94 95 integer :: ibox 96 logical, allocatable :: inside2(:) 97 98 ! no push_sub because this function is called very frequently 99 SAFE_ALLOCATE(inside2(1:npoints)) 100 101 inside = .false. 102 do ibox = 1, union%n_boxes 103 call box_inside_vec(union%boxes(ibox), npoints, points, inside2) 104 inside = inside .or. inside2 105 end do 106 107 SAFE_DEALLOCATE_A(inside2) 108 109 end subroutine box_union_inside_vec 110 111 !-------------------------------------------------------------- 112 !> Checks if a point are inside the union box. 113 logical function box_union_inside(union, point) result(inside) 114 type(box_union_t), intent(in) :: union 115 FLOAT, intent(in) :: point(:) 116 117 integer :: ibox 118 119 ! no push_sub because this function is called very frequently 120 121 inside = .false. 122 do ibox = 1, union%n_boxes 123 if(box_inside(union%boxes(ibox), point)) inside = .true. 124 end do 125 126 end function box_union_inside 127 128 !-------------------------------------------------------------- 129 !> Returns number of boxes inside domain 130 pure integer function box_union_get_nboxes(union) result(nbox) 131 type(box_union_t), intent(in) :: union 132 133 ! no push_sub because this function is called very frequently 134 135 nbox = union%n_boxes 136 137 end function box_union_get_nboxes 138 139 !-------------------------------------------------------------- 140 !> Returns number of boxes inside domain 141 pure function box_union_get_center(union, ibox) result(x) 142 type(box_union_t), intent(in) :: union 143 integer, intent(in) :: ibox 144 FLOAT, dimension(MAX_DIM) :: x 145 146 ! no push_sub because this function is called very frequently 147 148 x = box_get_center(union%boxes(ibox)) 149 150 end function box_union_get_center 151end module box_union_oct_m 152 153 154!! Local Variables: 155!! mode: f90 156!! coding: utf-8 157!! End: 158