Octopus
multibox.F90
Go to the documentation of this file.
1!! Copyright (C) 2021 M. Oliveira, K. Lively, A. Obzhirov, I. Albar
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 multibox_oct_m
23 use box_oct_m
24 use debug_oct_m
25 use global_oct_m
29
30 implicit none
31
32 private
33 public :: &
34 multibox_t, &
36
38 type, abstract, extends(box_t) :: multibox_t
39 private
40 type(box_list_t), public :: list
41 contains
42 procedure :: bounds => multibox_bounds
43 procedure :: add_box => multibox_add_box
44 end type multibox_t
45
46contains
47
48 !--------------------------------------------------------------
49 subroutine multibox_end(this)
50 class(multibox_t), intent(inout) :: this
51
52 type(box_iterator_t) :: iter
53 class(box_t), pointer :: box
54
55 push_sub(multibox_end)
56
57 safe_deallocate_a(this%bounding_box_l)
58
59 call iter%start(this%list)
60 do while (iter%has_next())
61 box => iter%get_next()
62 safe_deallocate_p(box)
63 end do
64
65 pop_sub(multibox_end)
66 end subroutine multibox_end
67
68 !--------------------------------------------------------------
69 ! Although formally this will return the bounds of a union of boxes, we will
70 ! also use it for an intersection of boxes. This should be a reasonable (and
71 ! much simpler to implement) approximation to the bounds of an intersection of
72 ! boxes, as the intersection should be contained in the union (unless the box
73 ! is unbound, but in that case it does not matter how we handle this).
74 function multibox_bounds(this, axes) result(bounds)
75 class(multibox_t), intent(in) :: this
76 class(basis_vectors_t), optional, intent(in) :: axes
77 real(real64) :: bounds(2, this%dim)
78
79 real(real64) :: box_bounds(2, this%dim)
80 type(box_iterator_t) :: iter
81 class(box_t), pointer :: box
82
83 push_sub(multibox_bounds)
84
85 bounds(1, :) = m_huge
86 bounds(2, :) = -m_huge
87 call iter%start(this%list)
88 do while (iter%has_next())
89 box => iter%get_next()
90 if (this%is_inside_out()) then
91 call box%turn_inside_out()
92 end if
93 box_bounds = box%bounds(axes)
94
95 where (box_bounds(1,:) < bounds(1,:))
96 bounds(1, :) = box_bounds(1,:)
97 elsewhere (box_bounds(2,:) > bounds(2,:))
98 bounds(2, :) = box_bounds(2,:)
99 end where
100 if (this%is_inside_out()) then
101 call box%turn_inside_out()
102 end if
103 end do
104
105 pop_sub(multibox_bounds)
106 end function multibox_bounds
107
108 !--------------------------------------------------------------
109 subroutine multibox_add_box(this, new_box)
110 class(multibox_t), intent(inout) :: this
111 class(box_t), intent(in) :: new_box
112
113 integer :: idir
115 push_sub(multibox_add_box)
116
117 assert(this%dim == new_box%dim)
118
119 do idir = 1, this%dim
120 if (new_box%bounding_box_l(idir) > this%bounding_box_l(idir)) then
121 this%bounding_box_l(idir) = new_box%bounding_box_l(idir)
122 end if
123 end do
124
125 call this%list%add(new_box)
126
127 pop_sub(multibox_add_box)
128 end subroutine multibox_add_box
129
130end module multibox_oct_m
132!! Local Variables:
133!! mode: f90
134!! coding: utf-8
135!! End:
Box bounds along some axes.
Definition: box.F90:177
real(real64), parameter, public m_huge
Definition: global.F90:205
This module implements fully polymorphic linked lists, and some specializations thereof.
real(real64) function, dimension(2, this%dim) multibox_bounds(this, axes)
Definition: multibox.F90:168
subroutine multibox_add_box(this, new_box)
Definition: multibox.F90:203
subroutine, public multibox_end(this)
Definition: multibox.F90:143
class to tell whether a point is inside or outside
Definition: box.F90:141
Abstract class for boxes that are made up of a list of boxes.
Definition: multibox.F90:131