Octopus
box_parallelepiped.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2021 M. Oliveira, K. Lively, A. Obzhirov, I. Albar
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
24 use box_oct_m
26 use debug_oct_m
27 use global_oct_m
31 use unit_oct_m
33
34 implicit none
35
36 private
37 public :: box_parallelepiped_t
38
42 type, extends(box_shape_t) :: box_parallelepiped_t
43 private
44 real(real64), allocatable, public :: half_length(:)
45
46 integer, public :: n_periodic_boundaries = 0
47 contains
48 procedure :: shape_contains_points => box_parallelepiped_shape_contains_points
49 procedure :: get_surface_points => box_parallelepiped_shape_get_surface_points
50 procedure :: get_surface_point_info => box_parallelepiped_shape_get_surface_point_info
51 procedure :: write_info => box_parallelepiped_write_info
52 procedure :: short_info => box_parallelepiped_short_info
53 procedure :: regenerate => box_parallelepiped_regenerate
56
57 interface box_parallelepiped_t
58 procedure box_parallelepiped_constructor
59 end interface box_parallelepiped_t
60
61contains
62
63 !--------------------------------------------------------------
64 function box_parallelepiped_constructor(dim, center, axes, length, namespace, n_periodic_boundaries) result(box)
65 integer, intent(in) :: dim
66 real(real64), intent(in) :: center(dim)
67 real(real64), intent(in) :: axes(dim, dim)
68 real(real64), intent(in) :: length(dim)
69 type(namespace_t), intent(in) :: namespace
70 integer, optional, intent(in) :: n_periodic_boundaries
71 class(box_parallelepiped_t), pointer :: box
72
74
75 ! Allocate memory
76 safe_allocate(box)
77
78 ! Initialize box
79 safe_allocate(box%half_length(1:dim))
80 box%half_length = m_half * length
81
82 if (present(n_periodic_boundaries)) then
83 box%n_periodic_boundaries = n_periodic_boundaries
84 end if
85 call box_shape_init(box, namespace, dim, center, bounding_box_min=-box%half_length, bounding_box_max=box%half_length, axes=axes)
86
87 box%bounding_box_l = box%half_length + abs(center)
88
91
92 !--------------------------------------------------------------
93 subroutine box_parallelepiped_regenerate(box, dim, axes, length, namespace)
94 class(box_parallelepiped_t), intent(inout) :: box
95 integer, intent(in) :: dim
96 real(real64), intent(in) :: axes(dim, dim)
97 real(real64), intent(in) :: length(dim)
98 type(namespace_t), intent(in) :: namespace
99
100 real(real64) :: center(dim)
101
103
104 center = box%center
105 safe_deallocate_a(box%center)
106 safe_deallocate_a(box%bounding_box_l)
107
108 ! Initialize box
109 box%half_length = m_half * length
110 call box_shape_init(box, namespace, dim, center, bounding_box_min=-box%half_length, bounding_box_max=box%half_length, axes=axes)
111
112 box%bounding_box_l = box%half_length + abs(center)
113
116
117
118 !--------------------------------------------------------------
119 subroutine box_parallelepiped_finalize(this)
120 type(box_parallelepiped_t), intent(inout) :: this
121
123
124 call box_shape_end(this)
125 safe_deallocate_a(this%half_length)
126
128 end subroutine box_parallelepiped_finalize
129
130 !--------------------------------------------------------------
131 function box_parallelepiped_shape_contains_points(this, nn, xx, tol) result(contained)
132 class(box_parallelepiped_t), intent(in) :: this
133 integer, intent(in) :: nn
134 real(real64), contiguous, intent(in) :: xx(:,:)
135 real(real64), optional, intent(in) :: tol
136 logical :: contained(1:nn)
138 integer :: ip, idir
139 real(real64) :: ulimit(this%dim), llimit(this%dim), tol_
140
143 llimit = this%center - this%half_length - tol_
144 do idir = 1, this%dim
145 if (idir <= this%n_periodic_boundaries) then
146 ! When periodic, we exclude one of the faces from the box.
147 ulimit(idir) = this%center(idir) + this%half_length(idir) - tol_
148 else
149 ulimit(idir) = this%center(idir) + this%half_length(idir) + tol_
150 end if
151 end do
152
153 do ip = 1, nn
154 contained(ip) = .true.
155 do idir = 1, this%dim
156 contained(ip) = contained(ip) .and. xx(ip, idir) >= llimit(idir) .and. xx(ip, idir) <= ulimit(idir)
157 if (idir > this%n_periodic_boundaries) then
158 ! We only consider the box to be inside out along the non-periodic directions.
159 contained(ip) = contained(ip) .neqv. this%is_inside_out()
160 end if
161 end do
162 end do
163
165
166 !--------------------------------------------------------------
170 function box_parallelepiped_shape_get_surface_points(this,namespace,mesh_spacing,nn,xx, number_of_layers) result(surface_points)
171 class(box_parallelepiped_t), intent(in) :: this
172 type(namespace_t), intent(in) :: namespace
173 real(real64), intent(in) :: mesh_spacing(:)
174 integer, intent(in) :: nn
175 real(real64), intent(in) :: xx(:,:)
176 integer, optional, intent(in) :: number_of_layers
177
178 logical :: surface_points(1:nn)
179 integer :: idir, number_of_layers_
180 real(real64) :: shrink(this%dim), test_axis(this%dim)
181 class(box_parallelepiped_t), pointer :: shrinked_box
182
183 ! Check that each axis is along the Cartersian axis
184 do idir = 1, this%dim
185 test_axis = m_zero
186 test_axis(idir) = m_one
187 assert(all(abs(this%axes%vectors(:, idir) - test_axis(:)) < m_epsilon))
188 end do
189
190 number_of_layers_ = 1
191 if (present(number_of_layers)) number_of_layers_ = number_of_layers
192
193 ! Determine how much the shrink should be based on the spacing in each direction
194 shrink(:) = 1 - number_of_layers_ * (1 - box_boundary_delta) * (mesh_spacing(:) / this%half_length(:))
195
196 ! Then we create a shrunk parallelepiped box
197 shrinked_box => box_parallelepiped_t(this%dim, this%center, this%axes%vectors, &
198 m_two*this%half_length(:)*shrink(:), namespace)
199
200 ! Then we look for points of the old box that are not containted in the smaller box
201 ! box_parallelepiped_shape_contains_points will return the point contained in the shrinked box,
202 ! but we are interested in the ones not contained, as they will the surface points of the original parallelepiped
203 if (SIZE(xx, 1) == 3) then
204 surface_points = .not. box_parallelepiped_shape_contains_points(shrinked_box, nn, transpose(xx))
205 else
206 surface_points = .not. box_parallelepiped_shape_contains_points(shrinked_box, nn, xx)
207 end if
208
209 safe_deallocate_p(shrinked_box)
211
212 !--------------------------------------------------------------
213 subroutine box_parallelepiped_shape_get_surface_point_info(this, point_coordinates,mesh_spacing, normal_vector, surface_element)
214 class(box_parallelepiped_t), intent(in) :: this
215 real(real64), intent(in) :: point_coordinates(:)
216 real(real64), intent(in) :: mesh_spacing(:)
217 real(real64), intent(out) :: normal_vector(:)
218 real(real64), intent(out) :: surface_element
219
220 real(real64) :: vector_norm
221 integer :: idir, first_index, second_index, third_index, face_counter
222
225 ! Compute the normal vector to the surface point
226
227 ! For a parallelepiped, the normal vector to a point is simply a vector that has plus or minus one in the direction at
228 ! which the face is pointing, and zero in the others. For instance, if the top face is facing the z direction, then
229 ! the normal vector for all points in that face should be (0,0,1)
230
231 ! So first of all we have to determine in which face we are. To do that let us divide the vector containing the coordinates
232 ! of the point by the maximum length of the parallelepiped
233 normal_vector(:) = (point_coordinates(:) - this%center(:)) / this%half_length(:)
234 face_counter = m_zero
235 surface_element = m_zero
236
237 ! Loop over the three directions
238 do idir = 1, this%dim
239 ! Determine the indeces of the faces
240 first_index = this%dim - mod(idir, this%dim)
241 second_index = this%dim - mod(idir + 1, this%dim)
242 third_index = this%dim - mod(idir + 2, this%dim)
243 if (abs(normal_vector(first_index)) >= abs(normal_vector(second_index)) .and. &
244 abs(normal_vector(first_index)) >= abs(normal_vector(third_index))) then
245 ! Define the surface element for the first_index
246 surface_element = mesh_spacing(second_index) * mesh_spacing(third_index)
247 face_counter = face_counter + 1
248 ! Now we check whether we are on the face or on a vertex or corner
249 if (abs(normal_vector(second_index)) / abs(normal_vector(first_index)) < &
250 m_one - m_half * mesh_spacing(second_index) / this%half_length(first_index)) then
251 normal_vector(second_index) = m_zero
252 end if
253 if (abs(normal_vector(third_index)) / abs(normal_vector(first_index)) < &
254 m_one - m_half * mesh_spacing(third_index) / this%half_length(first_index)) then
255 normal_vector(third_index) = m_zero
256 end if
257 end if
258 end do
259
260 if (face_counter == 3) surface_element = surface_element * (m_three / m_four)
261
262 ! Now if the point lies in a face, then the normalized_coordinates has one component which is one, and all the others are
263 ! zero. If the point is on a vertex, more than one coordinate is non zero. the final step is to normalize the vector
264 vector_norm = norm2(normal_vector(:))
265 ! Normal vector should never be zero
266 assert(vector_norm > m_epsilon)
267 normal_vector(:) = normal_vector(:) / vector_norm
268
271
272 !--------------------------------------------------------------
273 subroutine box_parallelepiped_write_info(this, iunit, namespace)
274 class(box_parallelepiped_t), intent(in) :: this
275 integer, optional, intent(in) :: iunit
276 type(namespace_t), optional, intent(in) :: namespace
277
278 integer :: idir
279
281
282 write(message(1),'(2x,a)') 'Type = parallelepiped'
283 write(message(2),'(2x,3a, 99(f8.3,a))') 'Lengths [', trim(units_abbrev(units_out%length)), '] = (', &
284 (units_from_atomic(units_out%length, m_two*this%half_length(idir)), ',', idir = 1, this%dim - 1), &
285 units_from_atomic(units_out%length, m_two*this%half_length(this%dim)), ')'
286 call messages_info(2, iunit=iunit, namespace=namespace)
287
289 end subroutine box_parallelepiped_write_info
290
291 !--------------------------------------------------------------
292 character(len=BOX_INFO_LEN) function box_parallelepiped_short_info(this, unit_length) result(info)
293 class(box_parallelepiped_t), intent(in) :: this
294 type(unit_t), intent(in) :: unit_length
295
296 integer :: idir
297
299
300 write(info, '(a,a,a,99(f11.6,a))') 'BoxShape = parallelepiped; Lengths [', trim(units_abbrev(unit_length)),'] = [', &
301 (units_from_atomic(unit_length, m_two*this%half_length(idir)), ',', idir = 1, this%dim - 1), &
302 units_from_atomic(unit_length, m_two*this%half_length(this%dim)), ']'
303
308
309!! Local Variables:
310!! mode: f90
311!! coding: utf-8
312!! End:
subroutine info()
Definition: em_resp.F90:1096
real(real64), parameter, public box_boundary_delta
Definition: box.F90:132
class(box_parallelepiped_t) function, pointer box_parallelepiped_constructor(dim, center, axes, length, namespace, n_periodic_boundaries)
subroutine box_parallelepiped_shape_get_surface_point_info(this, point_coordinates, mesh_spacing, normal_vector, surface_element)
subroutine box_parallelepiped_regenerate(box, dim, axes, length, namespace)
subroutine box_parallelepiped_finalize(this)
subroutine box_parallelepiped_write_info(this, iunit, namespace)
logical function, dimension(1:nn) box_parallelepiped_shape_get_surface_points(this, namespace, mesh_spacing, nn, xx, number_of_layers)
Get a mask for the grid points telling which of them are surface points.
logical function, dimension(1:nn) box_parallelepiped_shape_contains_points(this, nn, xx, tol)
character(len=box_info_len) function box_parallelepiped_short_info(this, unit_length)
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), parameter, public m_two
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_four
Definition: global.F90:192
real(real64), parameter, public m_epsilon
Definition: global.F90:204
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
real(real64), parameter, public m_three
Definition: global.F90:191
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
Class implementing a parallelepiped box. Currently this is restricted to a rectangular cuboid (all th...
Base class for more specialized boxes that are defined by a shape and have a center and basis vectors...
Definition: box_shape.F90:134
int true(void)