30 use,
intrinsic :: iso_fortran_env
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(:)
58 integer,
public,
parameter :: &
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
74 real(real64) :: bounds(space%dim, 2)
75 integer :: cols_abshape_block, imdim, maxdim
77 real(real64) :: xx(space%dim), rr
78 real(real64) :: ufn_re, ufn_im
79 character(len=1024) :: user_def_expr
81 real(real64),
allocatable :: mf(:)
82 real(real64) :: abheight, abwidth, abwidth_def
89 this%ab_user_def = .false.
108 call parse_variable(namespace,
'AbsorbingBoundaries', not_absorbing, this%abtype)
117 if (this%abtype /= not_absorbing)
then
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."
157 select type (sb=>gr%box)
159 bounds(1,1) = sb%radius/
m_two
160 bounds(1,2) = sb%radius
162 bounds(:,1) = sb%half_length/
m_two
163 bounds(:,2) = sb%half_length
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
171 bounds(:,2) = 0.4_real64
176 select case (cols_abshape_block)
181 bounds(:,1) = bounds(1,1)
182 bounds(:,2) = bounds(1,2)
184 select type (sb=>gr%box)
186 if (bounds(1,2) > sb%radius)
then
187 bounds(1,2) = sb%radius
189 message(1) =
"Info: using spherical absorbing boundaries."
192 do imdim = 1, space%dim
193 if (bounds(imdim,2) > sb%half_length(imdim))
then
194 bounds(imdim,2) = sb%half_length(imdim)
197 message(1) =
"Info: using cubic absorbing boundaries."
200 if (bounds(1,2) > sb%half_length)
then
201 bounds(1,2) = sb%half_length
203 if (bounds(2,2) > sb%radius)
then
204 bounds(2,2) = sb%radius
206 message(1) =
"Info: using cylindrical absorbing boundaries."
211 this%ab_user_def = .
true.
212 safe_allocate(this%ab_ufn(1:gr%np))
221 this%ab_ufn(ip) = ufn_re
223 message(1) =
"Input: using user-defined function from expression:"
224 write(
message(2),
'(a,a)')
' F(x,y,z) = ', trim(user_def_expr)
227 message(1) =
"Input: ABShape block must have at least 2 columns."
241 abwidth_def = bounds(1,2) - bounds(1,1)
243 bounds(:, 1) = bounds(:, 2) - abwidth
245 select type (sb=>gr%box)
251 if (this%ab_user_def)
then
257 write(
message(1),
'(a,2a,9f12.6)') &
260 write(
message(2),
'(a,2a,9f12.6)') &
266 safe_allocate(mf(1:gr%np))
270 safe_allocate(this%mf(1:gr%np))
272 select case (this%abtype)
276 this%mf(:) = abheight * mf(:)
283 safe_deallocate_a(mf)
295 if (this%abtype /= not_absorbing)
then
296 safe_deallocate_a(this%mf)
305 class(
mesh_t),
intent(in) :: mesh
307 class(
space_t),
intent(in) :: space
311 if (space%dim < 3)
return
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(:)
343 real(real64) :: width(space%dim)
344 real(real64) :: xx(space%dim), rr, dd, ddv(space%dim), tmp(space%dim)
351 width = bounds(:, 2) - bounds(:,1)
356 if (this%ab_user_def)
then
357 dd = this%ab_ufn(ip) - bounds(1,1)
359 if (this%ab_ufn(ip) < bounds(1,2))
then
368 select type (sb => gr%box)
370 rr = norm2(xx - sb%center)
371 dd = rr - bounds(1,1)
373 if (dd < width(1))
then
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
393 mf(ip) = mf(ip) * tmp(dir)
395 mf(ip) =
m_one - mf(ip)
398 rr = norm2(xx(2:space%dim) - sb%center(2:space%dim))
401 ddv(1) = abs(xx(1) - sb%center(1)) - bounds(1,1)
402 ddv(2) = rr - bounds(2,1)
404 if (ddv(dir) >
m_zero)
then
405 if (ddv(dir) < width(dir))
then
411 mf(ip) = mf(ip) * tmp(dir)
413 mf(ip) =
m_one - mf(ip)
422 if (this%ab_user_def)
then
423 safe_deallocate_a(this%ab_ufn)
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
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_pi
some mathematical constants
real(real64), parameter, public m_one
This module implements the underlying real-space grid.
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.
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
subroutine, public messages_not_implemented(feature, namespace)
character(len=512), private msg
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_input_error(namespace, var, details, row, column)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
integer function, public parse_block(namespace, name, blk, check_varinfo_)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
character(len=20) pure function, public units_abbrev(this)
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.
Description of the grid, containing information on derivatives, stencil, and symmetries.
Describes mesh distribution to nodes.