Octopus
box_sphere.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
24 use debug_oct_m
25 use global_oct_m
26 use, intrinsic :: iso_fortran_env
30 use unit_oct_m
32
33 implicit none
34
35 private
36 public :: box_sphere_t
37
39 type, extends(box_shape_t) :: box_sphere_t
40 real(real64) :: radius
41 contains
42 procedure :: shape_contains_points => box_sphere_shape_contains_points
43 procedure :: get_surface_points => box_sphere_shape_get_surface_points
44 procedure :: get_surface_point_info => box_sphere_shape_get_surface_point_info
45 procedure :: write_info => box_sphere_write_info
46 procedure :: short_info => box_sphere_short_info
47 final :: box_sphere_finalize
48 end type box_sphere_t
49
50 interface box_sphere_t
51 procedure box_sphere_constructor
52 end interface box_sphere_t
53
54contains
55
56 !--------------------------------------------------------------
57 function box_sphere_constructor(dim, center, radius, namespace) result(box)
58 integer, intent(in) :: dim
59 real(real64), intent(in) :: center(dim)
60 real(real64), intent(in) :: radius
61 type(namespace_t), intent(in) :: namespace
62 class(box_sphere_t), pointer :: box
63
64 integer :: idir
65
67
68 ! Allocate memory
69 safe_allocate(box)
70
71 ! Initialize box
72 call box_shape_init(box, namespace, dim, center, bounding_box_min=[(-radius, idir=1, dim)], &
73 bounding_box_max=[(radius, idir=1, dim)])
74 box%radius = radius
75 box%bounding_box_l = abs(center) + radius
76
78 end function box_sphere_constructor
79
80 !--------------------------------------------------------------
81 subroutine box_sphere_finalize(this)
82 type(box_sphere_t), intent(inout) :: this
83
84 push_sub(box_sphere_finalize)
85
86 call box_shape_end(this)
87
88 pop_sub(box_sphere_finalize)
89 end subroutine box_sphere_finalize
90
91 !--------------------------------------------------------------
92 function box_sphere_shape_contains_points(this, nn, xx, tol) result(contained)
93 class(box_sphere_t), intent(in) :: this
94 integer, intent(in) :: nn
95 real(real64), contiguous, intent(in) :: xx(:,:)
96 real(real64), optional, intent(in) :: tol
97 logical :: contained(1:nn)
98
99 integer :: ip
100 real(real64) :: tol_
101
103
104 do ip = 1, nn
105 contained(ip) = sum((xx(ip, 1:this%dim) - this%center(1:this%dim))**2) <= (this%radius + tol_)**2 &
106 .neqv. this%is_inside_out()
107 end do
108
110
111 !--------------------------------------------------------------
115 function box_sphere_shape_get_surface_points(this, namespace, mesh_spacing, nn, xx, number_of_layers) result(surface_points)
116 class(box_sphere_t), intent(in) :: this
117 type(namespace_t), intent(in) :: namespace
118 real(real64), intent(in) :: mesh_spacing(:)
119 integer, intent(in) :: nn
120 real(real64), intent(in) :: xx(:,:)
121 integer, optional, intent(in) :: number_of_layers
122
123 logical :: surface_points(1:nn)
124 real(real64) :: shrink
125 real(real64), parameter :: compression = 0.94_real64
126 class(box_sphere_t), pointer :: shrinked_box
127 ! Check that the spacing is the same in each direction
128 assert(abs(mesh_spacing(1) - mesh_spacing(2)) < m_epsilon .and. abs(mesh_spacing(1) - mesh_spacing(3)) < m_epsilon)
129
130 ! First of all determine how much the shrink should be based on the spacing
131 shrink = 1 - compression * (mesh_spacing(1) / this%radius)
133 ! Then, we have to initialize the new sphere to be slightly smaller that the original one
134 shrinked_box => box_sphere_t(this%dim, this%center, this%radius * shrink, namespace)
136 ! Then we look for points of the old sphere that are not containted in the smaller sphere
137 ! box_sphere_shape_contains_points will return the point contained in the shrinked box, but we are interested in the ones not
138 ! contained, as they will the surface points of the original sphere
139 surface_points = .not. box_sphere_shape_contains_points(shrinked_box, nn, xx)
141 call box_shape_end(shrinked_box)
142 safe_deallocate_p(shrinked_box)
143
145
146 !--------------------------------------------------------------
147 subroutine box_sphere_shape_get_surface_point_info(this, point_coordinates, mesh_spacing, normal_vector, surface_element)
148 class(box_sphere_t), intent(in) :: this
149 real(real64), intent(in) :: point_coordinates(:)
150 real(real64), intent(in) :: mesh_spacing(:)
151 real(real64), intent(out) :: normal_vector(:)
152 real(real64), intent(out) :: surface_element
153
154 real(real64) :: dtheta, dphi, rr
155
157 dtheta = m_zero
158 dphi = m_zero
159
160 ! Compute the normal vector to the surface point
161
162 ! For a sphere, the normal vector to a point is simply the vector representing its radius, but normalized.
163 ! Therefore it is enough to normalize the vector containing the coordinates of the vector
164
165 ! So first of all we have to determine in which face we are. To do that let us divide the vector containing the coordinates
166 ! of the point by the maximum length of the parallelepiped
167 normal_vector(:) = point_coordinates(:) - this%center(:)
168
169 rr = sqrt(normal_vector(1)**2 + normal_vector(2)**2)
170 if (rr > m_zero) then
171 dtheta = normal_vector(1) * normal_vector(3) / (rr * this%radius**2) * mesh_spacing(1) + &
172 normal_vector(2) * normal_vector(3) / (rr * this%radius**2) * mesh_spacing(2) - rr / this%radius**2 * mesh_spacing(3)
173
174 dphi = (- normal_vector(2) * mesh_spacing(1) + normal_vector(1) * mesh_spacing(2)) / rr**2
175 end if
176
177 normal_vector(:) = normal_vector(:) / norm2(normal_vector)
178
179 ! Compute the surface element, which for a sphere is r^2 * sin(theta) dxdy, where theta = acos(z / r)
180 surface_element = abs(dtheta * dphi) * this%radius**2 * sin(acos(normal_vector(3)))
181
184
185 !--------------------------------------------------------------
186 subroutine box_sphere_write_info(this, iunit, namespace)
187 class(box_sphere_t), intent(in) :: this
188 integer, optional, intent(in) :: iunit
189 type(namespace_t), optional, intent(in) :: namespace
190
191 push_sub(box_sphere_write_info)
192
193 write(message(1),'(2x,a)') 'Type = sphere'
194 write(message(2),'(2x,3a,f7.3)') 'Radius [', trim(units_abbrev(units_out%length)), '] = ', &
195 units_from_atomic(units_out%length, this%radius)
196 call messages_info(2, iunit=iunit, namespace=namespace)
197
198 pop_sub(box_sphere_write_info)
199 end subroutine box_sphere_write_info
200
201 !--------------------------------------------------------------
202 character(len=BOX_INFO_LEN) function box_sphere_short_info(this, unit_length) result(info)
203 class(box_sphere_t), intent(in) :: this
204 type(unit_t), intent(in) :: unit_length
205
206 push_sub(box_sphere_short_info)
207
208 write(info,'(a,f11.6,a,a)') 'BoxShape = sphere; Radius =', units_from_atomic(unit_length, this%radius), ' ', &
209 trim(units_abbrev(unit_length))
210
211 pop_sub(box_sphere_short_info)
212 end function box_sphere_short_info
213
214end module box_sphere_oct_m
215
216!! Local Variables:
217!! mode: f90
218!! coding: utf-8
219!! End:
subroutine info()
Definition: em_resp.F90:1096
double acos(double __x) __attribute__((__nothrow__
double sin(double __x) __attribute__((__nothrow__
real(real64), parameter, public box_boundary_delta
Definition: box.F90:132
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
subroutine box_sphere_write_info(this, iunit, namespace)
Definition: box_sphere.F90:280
character(len=box_info_len) function box_sphere_short_info(this, unit_length)
Definition: box_sphere.F90:296
subroutine box_sphere_finalize(this)
Definition: box_sphere.F90:175
logical function, dimension(1:nn) box_sphere_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.
Definition: box_sphere.F90:209
logical function, dimension(1:nn) box_sphere_shape_contains_points(this, nn, xx, tol)
Definition: box_sphere.F90:186
subroutine box_sphere_shape_get_surface_point_info(this, point_coordinates, mesh_spacing, normal_vector, surface_element)
Definition: box_sphere.F90:241
class(box_sphere_t) function, pointer box_sphere_constructor(dim, center, radius, namespace)
Definition: box_sphere.F90:151
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_epsilon
Definition: global.F90:204
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
Base class for more specialized boxes that are defined by a shape and have a center and basis vectors...
Definition: box_shape.F90:134
Class implementing a spherical box.
Definition: box_sphere.F90:132