Octopus
box_minimum.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
23 use box_oct_m
25 use debug_oct_m
26 use global_oct_m
27 use lookup_oct_m
28 use math_oct_m
32 use unit_oct_m
34
35 implicit none
36
37 private
38 public :: box_minimum_t
39
43 type, extends(box_shape_t) :: box_minimum_t
44 private
45 real(real64), public :: radius = m_zero
46 integer :: n_sites
47 real(real64), allocatable :: site_position(:,:)
48 type(lookup_t) :: site_lookup
49 contains
50 procedure :: shape_contains_points => box_minimum_shape_contains_points
51 procedure :: write_info => box_minimum_write_info
52 procedure :: short_info => box_minimum_short_info
54 end type box_minimum_t
55
56 interface box_minimum_t
57 procedure box_minimum_constructor
58 end interface box_minimum_t
59
60contains
61
62 !--------------------------------------------------------------
63 function box_minimum_constructor(dim, radius, n_sites, site_position, namespace) result(box)
64 integer, intent(in) :: dim
65 real(real64), intent(in) :: radius
66 integer, intent(in) :: n_sites
67 real(real64), intent(in) :: site_position(1:dim,1:n_sites)
68 type(namespace_t), intent(in) :: namespace
69 class(box_minimum_t), pointer :: box
70
71 real(real64) :: center(dim)
72
74
75 ! Allocate memory
76 safe_allocate(box)
77
78 ! Initialize box
79 center = m_zero
80 call box_shape_init(box, namespace, dim, center, bounding_box_min=minval(abs(site_position), dim=2) - radius, &
81 bounding_box_max=maxval(abs(site_position), dim=2) + radius)
82 box%radius = radius
83 box%n_sites = n_sites
84 safe_allocate_source(box%site_position, site_position)
85
86 box%bounding_box_l = maxval(abs(site_position), dim=2) + box%radius
87
88 call lookup_init(box%site_lookup, box%dim, box%n_sites, box%site_position)
89
91 end function box_minimum_constructor
92
93 !--------------------------------------------------------------
94 subroutine box_minimum_finalize(this)
95 type(box_minimum_t), intent(inout) :: this
96
97 push_sub(box_minimum_finalize)
98
99 call lookup_end(this%site_lookup)
100
101 safe_deallocate_a(this%site_position)
102
103 call box_shape_end(this)
104
105 pop_sub(box_minimum_finalize)
106 end subroutine box_minimum_finalize
107
108 !--------------------------------------------------------------
109 function box_minimum_shape_contains_points(this, nn, xx) result(contained)
110 class(box_minimum_t), intent(in) :: this
111 integer, intent(in) :: nn
112 real(real64), contiguous, intent(in) :: xx(:,:)
113 logical :: contained(1:nn)
114
115 integer :: ip
116 integer, allocatable :: nlist(:)
117
118 safe_allocate(nlist(1:nn))
119
120 call lookup_get_list(this%site_lookup, nn, xx, this%radius + box_boundary_delta, nlist)
121
122 do ip = 1, nn
123 contained(ip) = nlist(ip) /= 0 .neqv. this%is_inside_out()
124 end do
125
126 safe_deallocate_a(nlist)
127
129
130 !--------------------------------------------------------------
131 subroutine box_minimum_write_info(this, iunit, namespace)
132 class(box_minimum_t), intent(in) :: this
133 integer, optional, intent(in) :: iunit
134 type(namespace_t), optional, intent(in) :: namespace
135
137
138 write(message(1), '(2x,a)') 'Type = minimum'
139 call messages_info(1, iunit, namespace=namespace)
140 write(message(1),'(2x,3a,f7.3)') 'Radius [', trim(units_abbrev(units_out%length)), '] = ', &
141 units_from_atomic(units_out%length, this%radius)
142 call messages_info(1, iunit, namespace=namespace)
147 !--------------------------------------------------------------
148 character(len=BOX_INFO_LEN) function box_minimum_short_info(this, unit_length) result(info)
149 class(box_minimum_t), intent(in) :: this
150 type(unit_t), intent(in) :: unit_length
151
152 push_sub(box_minimum_short_info)
153
154 write(info,'(a,f11.6,a,a)') 'BoxShape = minimum; Radius =', units_from_atomic(unit_length, this%radius), ' ', &
155 trim(units_abbrev(unit_length))
158 end function box_minimum_short_info
159
160end module box_minimum_oct_m
161
162!! Local Variables:
163!! mode: f90
164!! coding: utf-8
165!! End:
subroutine info()
Definition: em_resp.F90:1096
subroutine box_minimum_finalize(this)
class(box_minimum_t) function, pointer box_minimum_constructor(dim, radius, n_sites, site_position, namespace)
subroutine box_minimum_write_info(this, iunit, namespace)
character(len=box_info_len) function box_minimum_short_info(this, unit_length)
logical function, dimension(1:nn) box_minimum_shape_contains_points(this, nn, xx)
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:175
subroutine, public box_shape_end(this)
Definition: box_shape.F90:271
real(real64), parameter, public m_zero
Definition: global.F90:187
subroutine, public lookup_end(this)
Definition: lookup.F90:158
subroutine, public lookup_init(this, dim, nobjs, pos)
Definition: lookup.F90:140
subroutine, public lookup_get_list(this, npoint, points, radius, nlist, list)
Definition: lookup.F90:184
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
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 box that is a union of spheres. We do this in a specific class instead of using ...
Base class for more specialized boxes that are defined by a shape and have a center and basis vectors...
Definition: box_shape.F90:134