Octopus
absorbing_boundaries.F90
Go to the documentation of this file.
1!! Copyright (C) 2015 U. De Giovannini
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#include "global.h"
19
21 use box_oct_m
26 use debug_oct_m
27 use global_oct_m
28 use grid_oct_m
30 use, intrinsic :: iso_fortran_env
31 use mesh_oct_m
34 use parser_oct_m
36 use space_oct_m
37 use unit_oct_m
40
41 implicit none
42
43 private
44 public :: &
48
50 private
51 integer, public :: abtype
52 real(real64), allocatable, public :: mf(:)
53 type(cube_function_t) :: cf
54 logical :: ab_user_def
55 real(real64), allocatable :: ab_ufn(:)
57
58 integer, public, parameter :: &
59 NOT_ABSORBING = 0, &
60 mask_absorbing = 1, &
62 exterior = 3
63
64contains
65
66 ! ---------------------------------------------------------
67 subroutine absorbing_boundaries_init(this, namespace, space, gr)
68 type(absorbing_boundaries_t), intent(out) :: this
69 type(namespace_t), intent(in) :: namespace
70 class(space_t), intent(in) :: space
71 type(grid_t), intent(in) :: gr
72
73 integer :: ip
74 real(real64) :: bounds(space%dim, 2)
75 integer :: cols_abshape_block, imdim, maxdim
76
77 real(real64) :: xx(space%dim), rr
78 real(real64) :: ufn_re, ufn_im
79 character(len=1024) :: user_def_expr
80
81 real(real64), allocatable :: mf(:)
82 real(real64) :: abheight, abwidth, abwidth_def
83 type(block_t) :: blk
84
86
87 bounds = m_zero
88
89 this%ab_user_def = .false.
90
91 !%Variable AbsorbingBoundaries
92 !%Type flag
93 !%Default not_absorbing
94 !%Section Time-Dependent::Absorbing Boundaries
95 !%Description
96 !% To improve the quality of the spectra by avoiding the formation of
97 !% standing density waves, one can make the boundaries of the simulation
98 !% box absorbing and use exterior complex scaling.
99 !%Option not_absorbing 0
100 !% Reflecting boundaries.
101 !%Option mask 1
102 !% Absorbing boundaries with a mask function.
103 !%Option cap 2
104 !% Absorbing boundaries with a complex absorbing potential.
105 !%Option exterior 3
106 !% Exterior complex scaling (not yet implemented).
107 !%End
108 call parse_variable(namespace, 'AbsorbingBoundaries', not_absorbing, this%abtype)
109 if (.not. varinfo_valid_option('AbsorbingBoundaries', this%abtype, is_flag = .true.)) then
110 call messages_input_error(namespace, 'AbsorbingBoundaries')
111 end if
112
113 if (this%abtype == exterior) then
114 call messages_not_implemented('Exterior complex scaling', namespace=namespace)
115 end if
116
117 if (this%abtype /= not_absorbing) then
118 call messages_print_with_emphasis(msg='Absorbing Boundaries', namespace=namespace)
119
120 !%Variable ABCapHeight
121 !%Type float
122 !%Default -0.2 a.u.
123 !%Section Time-Dependent::Absorbing Boundaries
124 !%Description
125 !% When <tt>AbsorbingBoundaries = cap</tt>, this is the height of the imaginary potential.
126 !%End
127 if (this%abtype == imaginary_absorbing) then
128 call parse_variable(namespace, 'ABCapHeight', -0.2_real64, abheight, units_inp%energy)
129 end if
130
131 !%Variable ABShape
132 !%Type block
133 !%Section Time-Dependent::Absorbing Boundaries
134 !%Description
135 !% Set the shape of the absorbing boundaries. Here you can set the inner
136 !% and outer bounds by setting the block as follows:
137 !%
138 !% <tt>%ABShape
139 !% <br>&nbsp;&nbsp; inner | outer | "user-defined"
140 !% <br>%</tt>
141 !%
142 !% The optional 3rd column is a user-defined expression for the absorbing
143 !% boundaries. For example, <math>r</math> creates a spherical absorbing zone for
144 !% coordinates with <math>{\tt inner} < r < {\tt outer}</math>, and <math>z</math> creates an
145 !% absorbing plane.
146 !% Note, values <tt>outer</tt> larger than the box size may lead in these cases to
147 !% unexpected reflection behaviours.
148 !% If no expression is given, the absorbing zone follows the edges of the
149 !% box (not valid for user-defined box).
150 !%End
152 cols_abshape_block = 0
153 if (parse_block(namespace, 'ABShape', blk) < 0) then
154 message(1) = "Input: ABShape not specified. Using default values for absorbing boundaries."
155 call messages_info(1, namespace=namespace)
156
157 select type (sb=>gr%box)
158 type is (box_sphere_t)
159 bounds(1,1) = sb%radius/m_two
160 bounds(1,2) = sb%radius
161 type is (box_parallelepiped_t)
162 bounds(:,1) = sb%half_length/m_two
163 bounds(:,2) = sb%half_length
164 type is (box_cylinder_t)
165 bounds(1,2) = sb%half_length
166 bounds(2,2) = sb%radius
167 bounds(1,1) = bounds(1,1)/m_two
168 bounds(2,1) = bounds(2,2)/m_two
169 class default
170 bounds(:,1) = m_zero
171 bounds(:,2) = 0.4_real64
172 end select
173 else
174 cols_abshape_block = parse_block_cols(blk, 0)
175
176 select case (cols_abshape_block)
177 case (2)
178 call parse_block_float(blk, 0, 0, bounds(1,1), units_inp%length)
179 call parse_block_float(blk, 0, 1, bounds(1,2), units_inp%length)
180 ! Apply single inner/outer pair to all dimensions.
181 bounds(:,1) = bounds(1,1)
182 bounds(:,2) = bounds(1,2)
183
184 select type (sb=>gr%box)
185 type is (box_sphere_t)
186 if (bounds(1,2) > sb%radius) then
187 bounds(1,2) = sb%radius
188 end if
189 message(1) = "Info: using spherical absorbing boundaries."
190
191 type is (box_parallelepiped_t)
192 do imdim = 1, space%dim
193 if (bounds(imdim,2) > sb%half_length(imdim)) then
194 bounds(imdim,2) = sb%half_length(imdim)
195 end if
196 end do
197 message(1) = "Info: using cubic absorbing boundaries."
198
199 type is (box_cylinder_t)
200 if (bounds(1,2) > sb%half_length) then
201 bounds(1,2) = sb%half_length
202 end if
203 if (bounds(2,2) > sb%radius) then
204 bounds(2,2) = sb%radius
205 end if
206 message(1) = "Info: using cylindrical absorbing boundaries."
207 end select
208 call messages_info(1, namespace=namespace)
209
210 case (3)
211 this%ab_user_def = .true.
212 safe_allocate(this%ab_ufn(1:gr%np))
213 this%ab_ufn = m_zero
214 call parse_block_float( blk, 0, 0, bounds(1,1), units_inp%length)
215 call parse_block_float( blk, 0, 1, bounds(1,2), units_inp%length)
216 call parse_block_string(blk, 0, 2, user_def_expr)
217 do ip = 1, gr%np
218 xx = units_from_atomic(units_inp%length, gr%x(ip, :))
219 rr = norm2(xx)
220 call parse_expression(ufn_re, ufn_im, space%dim, xx, rr, m_zero, user_def_expr)
221 this%ab_ufn(ip) = ufn_re
222 end do
223 message(1) = "Input: using user-defined function from expression:"
224 write(message(2),'(a,a)') ' F(x,y,z) = ', trim(user_def_expr)
225 call messages_info(2, namespace=namespace)
226 case default
227 message(1) = "Input: ABShape block must have at least 2 columns."
228 call messages_fatal(1, namespace=namespace)
229 end select
230
231 call parse_block_end(blk)
232 end if
233
234 !%Variable ABWidth
235 !%Type float
236 !%Section Time-Dependent::Absorbing Boundaries
237 !%Description
238 !% Specifies the boundary width. For a finer control over the absorbing boundary
239 !% shape use ABShape.
240 !%End
241 abwidth_def = bounds(1,2) - bounds(1,1)
242 call parse_variable(namespace, 'ABWidth', abwidth_def, abwidth, units_inp%length)
243 bounds(:, 1) = bounds(:, 2) - abwidth
244
245 select type (sb=>gr%box)
246 type is (box_sphere_t)
247 maxdim = 1
248 type is (box_cylinder_t)
249 maxdim = 2
250 class default
251 if (this%ab_user_def) then
252 maxdim = 1
253 else
254 maxdim = space%dim
255 end if
256 end select
257 write(message(1),'(a,2a,9f12.6)') &
258 " Absorbing boundary spans from [", trim(units_abbrev(units_inp%length)), '] =', &
259 (units_from_atomic(units_inp%length, bounds(imdim,1)), imdim=1,maxdim)
260 write(message(2),'(a,2a,9f12.6)') &
261 " to [", trim(units_abbrev(units_inp%length)), '] =', &
262 (units_from_atomic(units_inp%length, bounds(imdim,2)), imdim=1,maxdim)
263 call messages_info(2, namespace=namespace)
264
265 ! generate boundary function
266 safe_allocate(mf(1:gr%np))
267 call absorbing_boundaries_generate_mf(this, space, gr, bounds, mf)
268
269 ! mask or cap
270 safe_allocate(this%mf(1:gr%np))
271
272 select case (this%abtype)
273 case (mask_absorbing)
274 this%mf = m_one - mf
276 this%mf(:) = abheight * mf(:)
277 end select
278
279 if (debug%info) call absorbing_boundaries_write_info(this, gr, namespace, space)
280
281 call messages_print_with_emphasis(namespace=namespace)
282
283 safe_deallocate_a(mf)
284 end if
285
287 end subroutine absorbing_boundaries_init
288
289 ! ---------------------------------------------------------
290 subroutine absorbing_boundaries_end(this)
291 type(absorbing_boundaries_t), intent(inout) :: this
292
294
295 if (this%abtype /= not_absorbing) then
296 safe_deallocate_a(this%mf)
297 end if
298
300 end subroutine absorbing_boundaries_end
301
302 ! ---------------------------------------------------------
303 subroutine absorbing_boundaries_write_info(this, mesh, namespace, space)
304 type(absorbing_boundaries_t), intent(in) :: this
305 class(mesh_t), intent(in) :: mesh
306 type(namespace_t), intent(in) :: namespace
307 class(space_t), intent(in) :: space
308
309 integer :: err
310
311 if (space%dim < 3) return
312
314
315 if (this%abtype == mask_absorbing) then
316 call dio_function_output(io_function_fill_how("VTK"), "./td.general", "mask", namespace, space, mesh, &
317 this%mf(1:mesh%np), unit_one, err)
318 call dio_function_output(io_function_fill_how("PlaneZ"), "./td.general", "mask", namespace, space, mesh, &
319 this%mf(1:mesh%np), unit_one, err)
320 else if (this%abtype == imaginary_absorbing) then
321 call dio_function_output(io_function_fill_how("VTK"), "./td.general", "cap", namespace, space, mesh, &
322 this%mf(1:mesh%np), unit_one, err)
323 call dio_function_output(io_function_fill_how("PlaneZ"), "./td.general", "cap", namespace, space, mesh, &
324 this%mf(1:mesh%np), unit_one, err)
325 call dio_function_output(io_function_fill_how("PlaneX"), "./td.general", "cap", namespace, space, mesh, &
326 this%mf(1:mesh%np), unit_one, err)
327 call dio_function_output(io_function_fill_how("PlaneY"), "./td.general", "cap", namespace, space, mesh, &
328 this%mf(1:mesh%np), unit_one, err)
329 end if
330
333
334 ! ---------------------------------------------------------
335 subroutine absorbing_boundaries_generate_mf(this, space, gr, bounds, mf)
336 type(absorbing_boundaries_t), intent(inout) :: this
337 class(space_t), intent(in) :: space
338 type(grid_t), intent(in) :: gr
339 real(real64), intent(in) :: bounds(1:space%dim, 1:2)
340 real(real64), intent(inout) :: mf(:)
341
342 integer :: ip, dir
343 real(real64) :: width(space%dim)
344 real(real64) :: xx(space%dim), rr, dd, ddv(space%dim), tmp(space%dim)
345
347
348 mf = m_zero
349
350 ! generate the boundaries on the mesh
351 width = bounds(:, 2) - bounds(:,1)
352
353 do ip = 1, gr%np
354 xx = gr%x(ip, :)
355
356 if (this%ab_user_def) then
357 dd = this%ab_ufn(ip) - bounds(1,1)
358 if (dd > m_zero) then
359 if (this%ab_ufn(ip) < bounds(1,2)) then
360 mf(ip) = sin(dd * m_pi / (m_two * (width(1))))**2
361 else
362 mf(ip) = m_one
363 end if
364 end if
365
366 else ! this%ab_user_def == .false.
367
368 select type (sb => gr%box)
369 type is (box_sphere_t)
370 rr = norm2(xx - sb%center)
371 dd = rr - bounds(1,1)
372 if (dd > m_zero) then
373 if (dd < width(1)) then
374 mf(ip) = sin(dd * m_pi / (m_two * (width(1))))**2
375 else
376 mf(ip) = m_one
377 end if
378 end if
379
380 type is (box_parallelepiped_t)
381 ! We are filling from the center opposite to the spherical case
382 tmp = m_one
383 mf(ip) = m_one
384 ddv = abs(xx - sb%center) - bounds(:, 1)
385 do dir = 1, space%dim
386 if (ddv(dir) > m_zero) then
387 if (ddv(dir) < width(dir)) then
388 tmp(dir) = m_one - sin(ddv(dir) * m_pi / (m_two * (width(dir))))**2
389 else
390 tmp(dir) = m_zero
391 end if
392 end if
393 mf(ip) = mf(ip) * tmp(dir)
394 end do
395 mf(ip) = m_one - mf(ip)
397 type is (box_cylinder_t)
398 rr = norm2(xx(2:space%dim) - sb%center(2:space%dim))
399 tmp = m_one
400 mf(ip) = m_one
401 ddv(1) = abs(xx(1) - sb%center(1)) - bounds(1,1)
402 ddv(2) = rr - bounds(2,1)
403 do dir = 1, 2
404 if (ddv(dir) > m_zero) then
405 if (ddv(dir) < width(dir)) then
406 tmp(dir) = m_one - sin(ddv(dir) * m_pi / (m_two * (width(dir))))**2
407 else
408 tmp(dir) = m_zero
409 end if
410 end if
411 mf(ip) = mf(ip) * tmp(dir)
412 end do
413 mf(ip) = m_one - mf(ip)
414
415 class default
416 ! Other box shapes are not implemented
417 assert(.false.)
418 end select
419 end if
420 end do
421
422 if (this%ab_user_def) then
423 safe_deallocate_a(this%ab_ufn)
424 end if
425
430
431!! Local Variables:
432!! mode: f90
433!! coding: utf-8
434!! End:
double sin(double __x) __attribute__((__nothrow__
integer, parameter, public imaginary_absorbing
subroutine absorbing_boundaries_write_info(this, mesh, namespace, space)
integer, parameter, public mask_absorbing
subroutine absorbing_boundaries_generate_mf(this, space, gr, bounds, mf)
subroutine, public absorbing_boundaries_end(this)
integer, parameter, public exterior
subroutine, public absorbing_boundaries_init(this, namespace, space, gr)
type(debug_t), save, public debug
Definition: debug.F90:156
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_pi
some mathematical constants
Definition: global.F90:186
real(real64), parameter, public m_one
Definition: global.F90:189
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public dio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
integer(int64) function, public io_function_fill_how(where)
Use this function to quickly plot functions for debugging purposes: call dio_function_output(io_funct...
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:921
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1114
character(len=512), private msg
Definition: messages.F90:166
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:714
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
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_inp
the units systems for reading and writing
type(unit_t), public unit_one
some special units required for particular quantities
Class implementing a cylinder box. The cylinder axis is always along the first direction defined by t...
Class implementing a parallelepiped box. Currently this is restricted to a rectangular cuboid (all th...
Class implementing a spherical box.
Definition: box_sphere.F90:132
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
Describes mesh distribution to nodes.
Definition: mesh.F90:186
int true(void)