Octopus
box_shape.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 box_shape_oct_m
22 use box_oct_m
24 use debug_oct_m
25 use global_oct_m
26 use math_oct_m
30
31 implicit none
32
33 private
34 public :: &
38
41 type, abstract, extends(box_t) :: box_shape_t
42 private
43 real(real64), allocatable, public :: center(:)
44 type(basis_vectors_t), public :: axes
45 real(real64), allocatable :: bounding_box(:,:)
48 contains
49 procedure :: contains_points => box_shape_contains_points
50 procedure :: bounds => box_shape_bounds
51 procedure(box_shape_shape_contains_points), deferred :: shape_contains_points
52 end type box_shape_t
53
54 abstract interface
55 function box_shape_shape_contains_points(this, nn, xx, tol) result(contained)
56 import :: box_shape_t
57 import :: real64
58 class(box_shape_t), intent(in) :: this
59 integer, intent(in) :: nn
60 real(real64), contiguous, intent(in) :: xx(:,:)
61 real(real64), optional, intent(in) :: tol
62 logical :: contained(1:nn)
64 end interface
65
66contains
67
68 !--------------------------------------------------------------
69 recursive function box_shape_contains_points(this, nn, xx, tol) result(contained)
70 class(box_shape_t), intent(in) :: this
71 integer, intent(in) :: nn
72 real(real64), contiguous, intent(in) :: xx(:,:)
73 real(real64), optional, intent(in) :: tol
74 logical :: contained(1:nn)
75
76 ! no push_sub because this function is called very frequently
77
78 contained = this%shape_contains_points(nn, this%axes%from_cartesian(nn, xx), tol)
79
80 end function box_shape_contains_points
81
82 !--------------------------------------------------------------
83 subroutine box_shape_init(this, namespace, dim, center, bounding_box_min, bounding_box_max, axes)
84 class(box_shape_t), intent(inout) :: this
85 type(namespace_t), intent(in) :: namespace
86 integer, intent(in) :: dim
87 real(real64), intent(in) :: center(dim)
88 real(real64), intent(in) :: bounding_box_min(dim)
89 real(real64), intent(in) :: bounding_box_max(dim)
90 real(real64), optional, intent(in) :: axes(dim, dim)
91
92 integer :: idir
93 real(real64) :: normalized_axes(dim, dim)
94
95 this%dim = dim
96
97 safe_allocate(this%center(1:dim))
98 safe_allocate(this%bounding_box_l(1:dim))
99 this%bounding_box_l = m_zero
100
101 this%center(1:dim) = center(1:dim)
102
103 if (present(axes)) then
104 do idir = 1, dim
105 normalized_axes(:, idir) = axes(:, idir)/norm2(axes(:, idir))
106 end do
107 else
108 normalized_axes = diagonal_matrix(dim, m_one)
109 end if
110 this%axes = basis_vectors_t(namespace, dim, normalized_axes)
111
112 if (allocated(this%bounding_box)) then
113 safe_deallocate_a(this%bounding_box)
114 end if
115 safe_allocate(this%bounding_box(1:2, 1:dim))
116 this%bounding_box(1,:) = bounding_box_min + this%center
117 this%bounding_box(2,:) = bounding_box_max + this%center
118
119 end subroutine box_shape_init
120
121 !--------------------------------------------------------------
129 function box_shape_bounds(this, axes) result(bounds)
130 class(box_shape_t), intent(in) :: this
131 class(basis_vectors_t), optional, intent(in) :: axes
132 real(real64) :: bounds(2, this%dim)
133
134 integer :: idir, ib
135 real(real64) :: vertex_red(this%dim), vertex_cart(this%dim)
139 if (this%is_inside_out()) then
140 ! The box is unbound.
141 bounds(1, :) = -m_huge
142 bounds(2, :) = m_huge
143 else
144 ! The bounds must be located at the axis-aligned bounding box vertices, so
145 ! we loop over the vertices, convert them to reduced coordinates, and find
146 ! the maximum and minimum coordinates. This algorithm might not be
147 ! optimal, but it will do for now.
148 bounds(1, :) = m_huge
149 bounds(2, :) = -m_huge
150 do idir = 1, this%dim
151 do ib = 1, 2
152 ! Vertex in box reduced coordinates
153 vertex_red = m_zero
154 vertex_red(idir) = this%bounding_box(ib, idir)
155
156 ! Transform vertex to Cartesian coordinates
157 vertex_cart = this%axes%to_cartesian(vertex_red)
158
159 ! Convert to reduced coordinates of the input axes
160 if (present(axes)) then
161 vertex_red = axes%from_cartesian(vertex_cart)
162 else
163 vertex_red = vertex_cart
164 end if
165
166 where (vertex_red < bounds(1,:))
167 bounds(1, :) = vertex_red
168 elsewhere (vertex_red > bounds(2,:))
169 bounds(2, :) = vertex_red
170 end where
171 end do
172 end do
173 end if
174
175 pop_sub(box_shape_bounds)
176 end function box_shape_bounds
177
178 !--------------------------------------------------------------
179 subroutine box_shape_end(this)
180 class(box_shape_t), intent(inout) :: this
181
182 safe_deallocate_a(this%center)
183 safe_deallocate_a(this%bounding_box)
184 safe_deallocate_a(this%bounding_box_l)
185
186 end subroutine box_shape_end
187
188end module box_shape_oct_m
189
190!! Local Variables:
191!! mode: f90
192!! coding: utf-8
193!! End:
recursive logical function, dimension(1:nn) box_shape_contains_points(this, nn, xx, tol)
Definition: box_shape.F90:163
subroutine, public box_shape_init(this, namespace, dim, center, bounding_box_min, bounding_box_max, axes)
Definition: box_shape.F90:177
subroutine, public box_shape_end(this)
Definition: box_shape.F90:273
real(real64) function, dimension(2, this%dim) box_shape_bounds(this, axes)
Returns the bounding box of the shape. This is a bounding box aligned with some given axes (if presen...
Definition: box_shape.F90:223
real(real64), parameter, public m_huge
Definition: global.F90:206
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_one
Definition: global.F90:189
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
Vectors defining a basis in a vector space. This class provides methods to convert vector coordinates...
class to tell whether a point is inside or outside
Definition: box.F90:141
Base class for more specialized boxes that are defined by a shape and have a center and basis vectors...
Definition: box_shape.F90:134