Octopus
box_intersection.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
22 use box_oct_m
23 use debug_oct_m
26 use global_oct_m
31 use unit_oct_m
32
33 implicit none
34
35 private
36 public :: &
38
46 type, extends(multibox_t) :: box_intersection_t
47 private
48 contains
49 procedure :: contains_points => box_intersection_contains_points
50 procedure :: write_info => box_intersection_write_info
51 procedure :: short_info => box_intersection_short_info
53 end type box_intersection_t
54
55 interface box_intersection_t
56 procedure box_intersection_constructor
57 end interface box_intersection_t
58
59contains
60
61 !--------------------------------------------------------------
62 function box_intersection_constructor(dim) result(box)
63 integer, intent(in) :: dim
64 class(box_intersection_t), pointer :: box
65
67
68 ! Allocate memory
69 safe_allocate(box)
70 safe_allocate(box%bounding_box_l(1:dim))
71
72 ! Initialize box
73 box%dim = dim
74 box%bounding_box_l = m_zero
75
78
79 !--------------------------------------------------------------
80 subroutine box_intersection_finalize(this)
81 type(box_intersection_t), intent(inout) :: this
82
84
85 call multibox_end(this)
86
88 end subroutine box_intersection_finalize
89
90 !--------------------------------------------------------------
91 recursive function box_intersection_contains_points(this, nn, xx) result(contained)
92 class(box_intersection_t), intent(in) :: this
93 integer, intent(in) :: nn
94 real(real64), contiguous, intent(in) :: xx(:,:)
95 logical :: contained(nn)
96
97 integer :: ip
98 real(real64) :: point(1:this%dim)
99 type(box_iterator_t) :: iter
100 class(box_t), pointer :: box
101
102 ! A point must be inside all boxes to be considered inside an intersection of boxes
103 do ip = 1, nn
104 point(1:this%dim) = xx(ip, 1:this%dim)
105 contained(ip) = .true.
106
107 call iter%start(this%list)
108 do while (iter%has_next())
109 box => iter%get_next()
110 contained(ip) = box%contains_point(point)
111 if (.not. contained(ip)) exit
112 end do
113
114 contained(ip) = contained(ip) .neqv. this%is_inside_out()
115 end do
116
118
119 !--------------------------------------------------------------
120 subroutine box_intersection_write_info(this, iunit, namespace)
121 class(box_intersection_t), intent(in) :: this
122 integer, optional, intent(in) :: iunit
123 type(namespace_t), optional, intent(in) :: namespace
124
126
127 ! Todo: need to decide how best to display the information of the boxes that make the intersection
128
130 end subroutine box_intersection_write_info
131
132 !--------------------------------------------------------------
133 character(len=BOX_INFO_LEN) function box_intersection_short_info(this, unit_length) result(info)
134 class(box_intersection_t), intent(in) :: this
135 type(unit_t), intent(in) :: unit_length
136
138
139 ! Todo: need to decide how best to display the information of the boxes that make the intersection
140 info = ''
141
146
147!! Local Variables:
148!! mode: f90
149!! coding: utf-8
150!! End:
subroutine info()
Definition: em_resp.F90:1096
subroutine box_intersection_finalize(this)
class(box_intersection_t) function, pointer box_intersection_constructor(dim)
recursive logical function, dimension(nn) box_intersection_contains_points(this, nn, xx)
character(len=box_info_len) function box_intersection_short_info(this, unit_length)
subroutine box_intersection_write_info(this, iunit, namespace)
real(real64), parameter, public m_zero
Definition: global.F90:187
This module implements fully polymorphic linked lists, and some specializations thereof.
subroutine, public multibox_end(this)
Definition: multibox.F90:143
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
Class implementing a box that is an intersection other boxes. Note that the we do not override the bo...
Abstract class for boxes that are made up of a list of boxes.
Definition: multibox.F90:131
int true(void)