30  use, 
intrinsic :: iso_fortran_env
 
   52    type(stencil_t),            
intent(inout) :: this
 
   53    integer,                    
intent(in)    :: dim
 
   54    class(coordinate_system_t), 
intent(in)    :: coord_system
 
   56    real(real64)   :: theta(dim), basis_vectors(dim, dim)
 
   59    real(real64)   :: norm, min_norm(dim)
 
   61    real(real64), 
parameter :: tol = 1e-14_real64
 
   67    this%stargeneral%narms = 0
 
   68    this%stargeneral%arms = 0
 
   76    select type (coord_system)
 
   78      basis_vectors = coord_system%basis%vectors
 
   80      message(1) = 
"Stencil stargeneral can currently only be used with affine coordinate systems." 
   86      theta(1) = 
acos(dot_product(basis_vectors(:, 1), basis_vectors(:, 2)))
 
   89        this%stargeneral%narms = this%stargeneral%narms + 1
 
   91        this%stargeneral%arms(this%stargeneral%narms, 1:dim) = arm
 
   93        this%stargeneral%narms = this%stargeneral%narms + 1
 
   95        this%stargeneral%arms(this%stargeneral%narms, 1:dim) = arm
 
  108    theta(1) = 
acos(dot_product(basis_vectors(:, 1), basis_vectors(:, 2)))
 
  109    if (abs(theta(1) - 
m_pi*
m_half) > 1e-8_real64) this%stargeneral%narms = this%stargeneral%narms + 1
 
  111    theta(2) = 
acos(dot_product(basis_vectors(:, 2), basis_vectors(:, 3)))
 
  112    if (abs(theta(2)-
m_pi*
m_half) > 1e-8_real64) this%stargeneral%narms = this%stargeneral%narms + 1
 
  114    theta(3) = 
acos(dot_product(basis_vectors(:, 3), basis_vectors(:, 1)))
 
  115    if (abs(theta(3) - 
m_pi*
m_half) > 1e-8_real64) this%stargeneral%narms = this%stargeneral%narms + 1
 
  119    if (this%stargeneral%narms == 0) 
then 
  131          if (ii == 0 .and. jj == 0) cycle
 
  132          if (ii == 0 .and. kk == 0) cycle
 
  133          if (kk == 0 .and. jj == 0) cycle
 
  136          if (abs(theta(1) - 
m_pi*
m_half) <= 1e-8_real64 .and. kk == 0) cycle
 
  137          if (abs(theta(2) - 
m_pi*
m_half) <= 1e-8_real64 .and. ii == 0) cycle
 
  138          if (abs(theta(3) - 
m_pi*
m_half) <= 1e-8_real64 .and. jj == 0) cycle
 
  140          norm = sum((ii*basis_vectors(:, 1) + jj*basis_vectors(:, 2) + kk*basis_vectors(:, 3))**2)
 
  142          if (min_norm(1) - norm > tol) 
then 
  143            if (min_norm(2) - min_norm(1) > tol .and. this%stargeneral%narms >= 2) 
then 
  144              if (min_norm(3) - min_norm(2) > tol .and. this%stargeneral%narms == 3) 
then 
  145                this%stargeneral%arms(3, 1:3) = this%stargeneral%arms(2, 1:3)
 
  146                min_norm(3) =  min_norm(2)
 
  148              this%stargeneral%arms(2, 1:3) = this%stargeneral%arms(1, 1:3)
 
  149              min_norm(2) = min_norm(1)
 
  152            this%stargeneral%arms(1, 1:3) = (/ii, jj, kk/)
 
  157          if (this%stargeneral%arms(1, 1) == -ii &
 
  158            .and. this%stargeneral%arms(1, 2) == -jj &
 
  159            .and. this%stargeneral%arms(1, 3) == -kk) cycle
 
  161          if (this%stargeneral%narms == 1) cycle
 
  163          if (min_norm(2) - norm > tol) 
then 
  164            if (min_norm(3) - min_norm(2) > tol .and. this%stargeneral%narms == 3) 
then 
  165              this%stargeneral%arms(3, 1:3) = this%stargeneral%arms(2, 1:3)
 
  166              min_norm(3) =  min_norm(2)
 
  169            this%stargeneral%arms(2, 1:3) = (/ii, jj, kk/)
 
  173          if (this%stargeneral%narms == 2) cycle
 
  176          if (this%stargeneral%arms(2, 1) == -ii &
 
  177            .and. this%stargeneral%arms(2, 2) == -jj &
 
  178            .and. this%stargeneral%arms(2, 3) == -kk) cycle
 
  180          if (min_norm(3) - norm > tol) 
then 
  182            this%stargeneral%arms(3, 1:3) = (/ii, jj, kk/)
 
  195    integer, 
intent(in) :: dim
 
  196    integer, 
intent(in) :: order
 
  204    n = n + 2 * order * this%stargeneral%narms
 
  215    integer, 
intent(in) :: dir
 
  216    integer, 
intent(in) :: order
 
  223    if (dir >= 1 .or. dir <= 3) 
then 
  240    integer,         
intent(in)    :: dim
 
  241    integer,         
intent(in)    :: order
 
  244    logical :: got_center
 
  258          this%points(i, n) = j
 
  267          this%points(i, n) = j
 
  273        do i = 1, this%stargeneral%narms
 
  275          this%points(1:2, n) = this%stargeneral%arms(i, 1:2)*j
 
  296          this%points(i, n) = j
 
  302        do i = 1, this%stargeneral%narms
 
  304          this%points(1:3, n) = this%stargeneral%arms(i, 1:3)*j
 
  321    integer, 
intent(in)          :: dim
 
  322    integer, 
intent(in)          :: order
 
  323    integer, 
intent(out)         :: pol(:,:)
 
  325    integer :: i, j, n, j1
 
  326    logical :: zeros_treated(dim)
 
  352        do i = 1, this%stargeneral%narms
 
  354          if (sum(this%stargeneral%arms(i,1:dim)) == 0) 
then 
  355            pol(1:2, n) = (/j,1/)
 
  357            pol(1:2, n) = (/1,j/)
 
  374      zeros_treated = .false.
 
  375      do i = 1, this%stargeneral%narms
 
  376        if (this%stargeneral%arms(i, 1) == 0) 
then 
  379            pol(1:3, n) = (/0, 1, j1/)
 
  381          zeros_treated(1) = .
true.
 
  382        else if (this%stargeneral%arms(i, 2) == 0) 
then 
  385            pol(1:3, n) = (/j1, 0, 1/)
 
  387          zeros_treated(2) = .
true.
 
  388        else if (this%stargeneral%arms(i, 3) == 0) 
then 
  391            pol(1:3, n) = (/1, j1, 0/)
 
  393          zeros_treated(3) = .
true.
 
  397      do i = 1, this%stargeneral%narms - count(zeros_treated)
 
  398        if (.not. zeros_treated(1)) 
then 
  401            pol(1:3, n) = (/0, 1, j1/)
 
  403          zeros_treated(1) = .
true.
 
  404        else if (.not. zeros_treated(2)) 
then 
  407            pol(1:3, n) = (/j1, 0, 1/)
 
  409          zeros_treated(2) = .
true.
 
  410        else if (.not. zeros_treated(3)) 
then 
  413            pol(1:3, n) = (/1, j1, 0/)
 
  415          zeros_treated(3) = .
true.
 
double acos(double __x) __attribute__((__nothrow__
 
real(real64), parameter, public m_huge
 
real(real64), parameter, public m_pi
some mathematical constants
 
real(real64), parameter, public m_half
 
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
 
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
 
This module defines stencils used in Octopus.
 
subroutine, public stencil_allocate(this, dim, size)
 
subroutine, public stencil_init_center(this)
 
This module defines routines, generating operators for a generalized star stencil.
 
subroutine, public stencil_stargeneral_get_arms(this, dim, coord_system)
Finds the direction of the arms of the star-general stencil as described in Natan et al....
 
subroutine, public stencil_stargeneral_pol_lapl(this, dim, order, pol)
 
integer function, public stencil_stargeneral_extent(dir, order)
Returns maximum extension of the stencil in spatial direction dir = 1, 2, 3 for a given discretizatio...
 
subroutine, public stencil_stargeneral_get_lapl(this, dim, order)
 
integer function, public stencil_stargeneral_size_lapl(this, dim, order)
 
The class representing the stencil, which is used for non-local mesh operations.