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          select type (sb=>gr%box)
 
  183            if (bounds(1,2) > sb%radius) 
then 
  184              bounds(1,2) = sb%radius
 
  186            message(1) = 
"Info: using spherical absorbing boundaries." 
  189            do imdim = 1, space%dim
 
  190              if (bounds(imdim,2) > sb%half_length(imdim)) 
then 
  191                bounds(imdim,2) = sb%half_length(imdim)
 
  194            message(1) = 
"Info: using cubic absorbing boundaries." 
  197            if (bounds(1,2) > sb%half_length) 
then 
  198              bounds(1,2) = sb%half_length
 
  200            if (bounds(1,2) > sb%radius) 
then 
  201              bounds(1,2) = sb%radius
 
  203            message(1) = 
"Info: using cylindrical absorbing boundaries." 
  208          this%ab_user_def = .
true.
 
  209          safe_allocate(this%ab_ufn(1:gr%np))
 
  218            this%ab_ufn(ip) = ufn_re
 
  220          message(1) = 
"Input: using user-defined function from expression:" 
  221          write(
message(2),
'(a,a)') 
'   F(x,y,z) = ', trim(user_def_expr)
 
  224          message(1) = 
"Input: ABShape block must have at least 2 columns." 
  239      abwidth_def = bounds(1,2) - bounds(1,1)
 
  241      bounds(:, 1) = bounds(:, 2) - abwidth
 
  243      select type (sb=>gr%box)
 
  249        if (this%ab_user_def) 
then 
  255      write(
message(1),
'(a,2a,9f12.6)') &
 
  258      write(
message(2),
'(a,2a,9f12.6)') &
 
  264      safe_allocate(mf(1:gr%np))
 
  268      safe_allocate(this%mf(1:gr%np))
 
  270      select case (this%abtype)
 
  274        this%mf(:) = abheight * mf(:)
 
  281      safe_deallocate_a(mf)
 
  293    if (this%abtype /= not_absorbing) 
then 
  294      safe_deallocate_a(this%mf)
 
  303    class(
mesh_t),                
intent(in) :: mesh
 
  305    class(
space_t),               
intent(in) :: space
 
  333    class(
space_t),               
intent(in)    :: space
 
  334    type(
grid_t),                 
intent(in)    :: gr
 
  335    real(real64),                 
intent(in)    :: bounds(1:space%dim, 1:2)
 
  336    real(real64),                 
intent(inout) :: mf(:)
 
  339    real(real64)   :: width(space%dim)
 
  340    real(real64)   :: xx(space%dim), rr, dd, ddv(space%dim), tmp(space%dim)
 
  347    width = bounds(:, 2) - bounds(:,1)
 
  352      if (this%ab_user_def) 
then 
  353        dd = this%ab_ufn(ip) - bounds(1,1)
 
  355          if (this%ab_ufn(ip) < bounds(1,2)) 
then 
  364        select type (sb => gr%box)
 
  366          rr = norm2(xx - sb%center)
 
  367          dd = rr -  bounds(1,1)
 
  369            if (dd  <  width(1)) 
then 
  380          ddv = abs(xx - sb%center) -  bounds(:, 1)
 
  381          do dir = 1, space%dim
 
  382            if (ddv(dir) > 
m_zero) 
then 
  383              if (ddv(dir)  <  width(dir)) 
then 
  389            mf(ip) = mf(ip) * tmp(dir)
 
  391          mf(ip) = 
m_one - mf(ip)
 
  394          rr = norm2(xx(2:space%dim) - sb%center(2:space%dim))
 
  397          ddv(1) = abs(xx(1) - sb%center(1)) - bounds(1,1)
 
  398          ddv(2) = rr         - bounds(2,1)
 
  400            if (ddv(dir) > 
m_zero) 
then 
  401              if (ddv(dir)  <  width(dir)) 
then 
  407            mf(ip) = mf(ip) * tmp(dir)
 
  409          mf(ip) = 
m_one - mf(ip)
 
  418    if (this%ab_user_def) 
then 
  419      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)
 
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
 
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
 
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)
 
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.