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
181 select type (sb=>gr%box)
182 type is (box_sphere_t)
183 if (bounds(1,2) > sb%radius) then
184 bounds(1,2) = sb%radius
185 end if
186 message(1) = "Info: using spherical absorbing boundaries."
187
188 type is (box_parallelepiped_t)
189 do imdim = 1, space%dim
190 if (bounds(imdim,2) > sb%half_length(imdim)) then
191 bounds(imdim,2) = sb%half_length(imdim)
192 end if
193 end do
194 message(1) = "Info: using cubic absorbing boundaries."
195
196 type is (box_cylinder_t)
197 if (bounds(1,2) > sb%half_length) then
198 bounds(1,2) = sb%half_length
199 end if
200 if (bounds(1,2) > sb%radius) then
201 bounds(1,2) = sb%radius
202 end if
203 message(1) = "Info: using cylindrical absorbing boundaries."
204 end select
205 call messages_info(1, namespace=namespace)
206
207 case (3)
208 this%ab_user_def = .true.
209 safe_allocate(this%ab_ufn(1:gr%np))
210 this%ab_ufn = m_zero
211 call parse_block_float( blk, 0, 0, bounds(1,1), units_inp%length)
212 call parse_block_float( blk, 0, 1, bounds(1,2), units_inp%length)
213 call parse_block_string(blk, 0, 2, user_def_expr)
214 do ip = 1, gr%np
215 xx = units_from_atomic(units_inp%length, gr%x(ip, :))
216 rr = norm2(xx)
217 call parse_expression(ufn_re, ufn_im, space%dim, xx, rr, m_zero, user_def_expr)
218 this%ab_ufn(ip) = ufn_re
219 end do
220 message(1) = "Input: using user-defined function from expression:"
221 write(message(2),'(a,a)') ' F(x,y,z) = ', trim(user_def_expr)
222 call messages_info(2, namespace=namespace)
223 case default
224 message(1) = "Input: ABShape block must have at least 2 columns."
225 call messages_fatal(1, namespace=namespace)
226 end select
227
228 call parse_block_end(blk)
229 end if
230
231 !%Variable ABWidth
232 !%Type float
233 !%Section Time-Dependent::Absorbing Boundaries
234 !%Description
235 !% Specifies the boundary width. For a finer control over the absorbing boundary
236 !% shape use ABShape.
237 !%End
238! call messages_obsolete_variable('ABWidth', 'ABShape', namespace=namespace)
239 abwidth_def = bounds(1,2) - bounds(1,1)
240 call parse_variable(namespace, 'ABWidth', abwidth_def, abwidth, units_inp%length)
241 bounds(:, 1) = bounds(:, 2) - abwidth
242
243 select type (sb=>gr%box)
244 type is (box_sphere_t)
245 maxdim = 1
246 type is (box_cylinder_t)
247 maxdim = 2
248 class default
249 if (this%ab_user_def) then
250 maxdim = 1
251 else
252 maxdim = space%dim
253 end if
254 end select
255 write(message(1),'(a,2a,9f12.6)') &
256 " Absorbing boundary spans from [", trim(units_abbrev(units_inp%length)), '] =', &
257 (units_from_atomic(units_inp%length, bounds(imdim,1)), imdim=1,maxdim)
258 write(message(2),'(a,2a,9f12.6)') &
259 " to [", trim(units_abbrev(units_inp%length)), '] =', &
260 (units_from_atomic(units_inp%length, bounds(imdim,2)), imdim=1,maxdim)
261 call messages_info(2, namespace=namespace)
262
263 ! generate boundary function
264 safe_allocate(mf(1:gr%np))
265 call absorbing_boundaries_generate_mf(this, space, gr, bounds, mf)
266
267 ! mask or cap
268 safe_allocate(this%mf(1:gr%np))
269
270 select case (this%abtype)
271 case (mask_absorbing)
272 this%mf = m_one - mf
274 this%mf(:) = abheight * mf(:)
275 end select
276
277 if (debug%info) call absorbing_boundaries_write_info(this, gr, namespace, space)
278
279 call messages_print_with_emphasis(namespace=namespace)
280
281 safe_deallocate_a(mf)
282 end if
283
285 end subroutine absorbing_boundaries_init
286
287 ! ---------------------------------------------------------
288 subroutine absorbing_boundaries_end(this)
289 type(absorbing_boundaries_t), intent(inout) :: this
290
292
293 if (this%abtype /= not_absorbing) then
294 safe_deallocate_a(this%mf)
295 end if
296
298 end subroutine absorbing_boundaries_end
299
300 ! ---------------------------------------------------------
301 subroutine absorbing_boundaries_write_info(this, mesh, namespace, space)
302 type(absorbing_boundaries_t), intent(in) :: this
303 class(mesh_t), intent(in) :: mesh
304 type(namespace_t), intent(in) :: namespace
305 class(space_t), intent(in) :: space
306
307 integer :: err
308
309 if (space%dim < 3) return
310
312
313 if (this%abtype == mask_absorbing) then
314 call dio_function_output(io_function_fill_how("VTK"), "./td.general", "mask", namespace, space, mesh, &
315 this%mf(1:mesh%np), unit_one, err)
316 call dio_function_output(io_function_fill_how("PlaneZ"), "./td.general", "mask", namespace, space, mesh, &
317 this%mf(1:mesh%np), unit_one, err)
318 else if (this%abtype == imaginary_absorbing) then
319 call dio_function_output(io_function_fill_how("VTK"), "./td.general", "cap", namespace, space, mesh, &
320 this%mf(1:mesh%np), unit_one, err)
321 call dio_function_output(io_function_fill_how("PlaneZ"), "./td.general", "cap", namespace, space, mesh, &
322 this%mf(1:mesh%np), unit_one, err)
323 call dio_function_output(io_function_fill_how("PlaneX"), "./td.general", "cap", namespace, space, mesh, &
324 this%mf(1:mesh%np), unit_one, err)
325 call dio_function_output(io_function_fill_how("PlaneY"), "./td.general", "cap", namespace, space, mesh, &
326 this%mf(1:mesh%np), unit_one, err)
327 end if
328
331
332 ! ---------------------------------------------------------
333 subroutine absorbing_boundaries_generate_mf(this, space, gr, bounds, mf)
334 type(absorbing_boundaries_t), intent(inout) :: this
335 class(space_t), intent(in) :: space
336 type(grid_t), intent(in) :: gr
337 real(real64), intent(in) :: bounds(1:space%dim, 1:2)
338 real(real64), intent(inout) :: mf(:)
339
340 integer :: ip, dir
341 real(real64) :: width(space%dim)
342 real(real64) :: xx(space%dim), rr, dd, ddv(space%dim), tmp(space%dim)
343
345
346 mf = m_zero
347
348 ! generate the boundaries on the mesh
349 width = bounds(:, 2) - bounds(:,1)
350
351 do ip = 1, gr%np
352 xx = gr%x(ip, :)
353
354 if (this%ab_user_def) then
355 dd = this%ab_ufn(ip) - bounds(1,1)
356 if (dd > m_zero) then
357 if (this%ab_ufn(ip) < bounds(1,2)) then
358 mf(ip) = sin(dd * m_pi / (m_two * (width(1))))**2
359 else
360 mf(ip) = m_one
361 end if
362 end if
363
364 else ! this%ab_user_def == .false.
365
366 select type (sb => gr%box)
367 type is (box_sphere_t)
368 rr = norm2(xx - sb%center)
369 dd = rr - bounds(1,1)
370 if (dd > m_zero) then
371 if (dd < width(1)) then
372 mf(ip) = sin(dd * m_pi / (m_two * (width(1))))**2
373 else
374 mf(ip) = m_one
375 end if
376 end if
377
378 type is (box_parallelepiped_t)
379 ! We are filling from the center opposite to the spherical case
380 tmp = m_one
381 mf(ip) = m_one
382 ddv = abs(xx - sb%center) - bounds(:, 1)
383 do dir = 1, space%dim
384 if (ddv(dir) > m_zero) then
385 if (ddv(dir) < width(dir)) then
386 tmp(dir) = m_one - sin(ddv(dir) * m_pi / (m_two * (width(dir))))**2
387 else
388 tmp(dir) = m_zero
389 end if
390 end if
391 mf(ip) = mf(ip) * tmp(dir)
392 end do
393 mf(ip) = m_one - mf(ip)
395 type is (box_cylinder_t)
396 rr = norm2(xx(2:space%dim) - sb%center(2:space%dim))
397 tmp = m_one
398 mf(ip) = m_one
399 ddv(1) = abs(xx(1) - sb%center(1)) - bounds(1,1)
400 ddv(2) = rr - bounds(2,1)
401 do dir = 1, 2
402 if (ddv(dir) > m_zero) then
403 if (ddv(dir) < width(dir)) then
404 tmp(dir) = m_one - sin(ddv(dir) * m_pi / (m_two * (width(dir))))**2
405 else
406 tmp(dir) = m_zero
407 end if
408 end if
409 mf(ip) = mf(ip) * tmp(dir)
410 end do
411 mf(ip) = m_one - mf(ip)
412
413 class default
414 ! Other box shapes are not implemented
415 assert(.false.)
416 end select
417 end if
418 end do
419
420 if (this%ab_user_def) then
421 safe_deallocate_a(this%ab_ufn)
422 end if
423
428
429!! Local Variables:
430!! mode: f90
431!! coding: utf-8
432!! 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:903
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1096
character(len=512), private msg
Definition: messages.F90:165
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:414
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:696
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:599
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:169
Describes mesh distribution to nodes.
Definition: mesh.F90:186
int true(void)