Populates the redistribute types describing the communication pattern required to deal with linked boundary conditions for the mixed wfb. This is used when passing_wfb = trapped_wfb = false. This is the default.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in), | dimension(:, :) | :: | l_links | ||
integer, | intent(in), | dimension(:, :) | :: | r_links | ||
integer, | intent(in) | :: | n_links_max | |||
type(redist_type), | intent(out) | :: | wfb_p | |||
type(redist_type), | intent(out) | :: | wfb_h |
subroutine setup_connected_bc_redistribute_mixed_wfb(l_links, r_links, n_links_max, &
wfb_p, wfb_h)
use gs2_layouts, only: g_lo, il_idx, ik_idx, it_idx, proc_id, ie_idx, is_idx
use le_grids, only: il_is_wfb
use redistribute, only: index_list_type, init_fill, delete_list, redist_type
use theta_grid, only: ntgrid
use mp, only: nproc, iproc
implicit none
integer, dimension(:, :), intent(in) :: l_links, r_links
integer, intent(in) :: n_links_max
type(redist_type), intent(out) ::wfb_p, wfb_h
type(index_list_type), dimension(0:nproc-1) :: from, to_p, to_h
integer, dimension (0:nproc-1) :: nn_from, nn_to
integer, dimension (3) :: to_low, from_low, to_high, from_high
integer :: il, ik, it, ncell, ip, iglo_star, j, n, iglo
integer :: iglo_right, ipright, iglo_left, ipleft, ie, is
nn_to = 0
nn_from = 0
!NOTE: No special restriction/counting for wfb_h unlike links_h
! Whilst we loop over the entire domain here, we should not have any work
! to do for ik, ie, is and il indices which do not belong to this processor's
! local g_lo range. We may be better off writing this as an explicit loop over
! the xyles dimensions with yles range being set my the min/max that this
! processor sees. For now we can just get the indices and cycle.
!$OMP PARALLEL DO DEFAULT(none) &
!$OMP PRIVATE(iglo, il, ik, ie, is, it, ncell, ip, iglo_right, iglo_left, &
!$OMP ipleft, j, ipright, iglo_star) &
!$OMP SHARED(g_lo, iproc, l_links, r_links) &
!$OMP REDUCTION(+: nn_to, nn_from) &
!$OMP SCHEDULE(static)
do iglo = g_lo%llim_world, g_lo%ulim_world
!CMR, 20/10/2013:
! Communicate pattern sends wfb endpoint to ALL linked cells
! (ntgrid,1,iglo) => ALL connected cells (j,1,iglo_conn)
! where j index in receive buffer = r_links(iglo)+1
! (-ntgrid,2,iglo) => ALL connected cells (j,2,iglo_conn)
! where j index in receive buffer = l_links(iglo)+1
!
il = il_idx(g_lo,iglo)
if (.not. il_is_wfb(il)) cycle
if (il > g_lo%il_max .or. il < g_lo%il_min) cycle
ik = ik_idx(g_lo,iglo)
if (ik > g_lo%ik_max .or. ik < g_lo%ik_min) cycle
ie = ie_idx(g_lo, iglo)
if (ie > g_lo%ie_max .or. ie < g_lo%ie_min) cycle
is = is_idx(g_lo, iglo)
if (is > g_lo%is_max .or. is < g_lo%is_min) cycle
it = it_idx(g_lo,iglo)
ncell = r_links(it, ik) + l_links(it, ik) + 1
if (ncell == 1) cycle
ip = proc_id(g_lo,iglo)
iglo_right = iglo ; iglo_left = iglo ; ipright = ip ; ipleft = ip
! v_par > 0:
!CMR: introduced iglo_star to make notation below less confusing
!
call find_leftmost_link (iglo, iglo_star, ipright)
do j = 1, ncell
! sender
if (ip == iproc) nn_from(ipright) = nn_from(ipright) + 1
! receiver
if (ipright == iproc) nn_to(ip) = nn_to(ip) + 1
call get_right_connection (iglo_star, iglo_right, ipright)
iglo_star=iglo_right
end do
! v_par < 0:
call find_rightmost_link (iglo, iglo_star, ipleft)
do j = 1, ncell
! sender
if (ip == iproc) nn_from(ipleft) = nn_from(ipleft) + 1
! receiver
if (ipleft == iproc) nn_to(ip) = nn_to(ip) + 1
call get_left_connection (iglo_star, iglo_left, ipleft)
iglo_star=iglo_left
end do
end do
!$OMP END PARALLEL DO
do ip = 0, nproc-1
if (nn_from(ip) > 0) then
allocate (from(ip)%first(nn_from(ip)))
allocate (from(ip)%second(nn_from(ip)))
allocate (from(ip)%third(nn_from(ip)))
endif
if (nn_to(ip) > 0) then
allocate (to_p(ip)%first(nn_to(ip)))
allocate (to_p(ip)%second(nn_to(ip)))
allocate (to_p(ip)%third(nn_to(ip)))
allocate (to_h(ip)%first(nn_to(ip)))
allocate (to_h(ip)%second(nn_to(ip)))
allocate (to_h(ip)%third(nn_to(ip)))
endif
end do
nn_from = 0
nn_to = 0
do iglo = g_lo%llim_world, g_lo%ulim_world
il = il_idx(g_lo,iglo)
if (.not. il_is_wfb(il)) cycle
if (il > g_lo%il_max .or. il < g_lo%il_min) cycle
ik = ik_idx(g_lo,iglo)
if (ik > g_lo%ik_max .or. ik < g_lo%ik_min) cycle
ie = ie_idx(g_lo, iglo)
if (ie > g_lo%ie_max .or. ie < g_lo%ie_min) cycle
is = is_idx(g_lo, iglo)
if (is > g_lo%is_max .or. is < g_lo%is_min) cycle
it = it_idx(g_lo,iglo)
ncell = r_links(it, ik) + l_links(it, ik) + 1
if (ncell == 1) cycle
ip = proc_id(g_lo,iglo)
iglo_right = iglo ; iglo_left = iglo ; ipright = ip ; ipleft = ip
! v_par > 0:
call find_leftmost_link (iglo, iglo_star, ipright)
do j = 1, ncell
! sender
if (ip == iproc) then
n = nn_from(ipright) + 1
nn_from(ipright) = n
from(ipright)%first(n) = ntgrid
from(ipright)%second(n) = 1
from(ipright)%third(n) = iglo
end if
! receiver
if (ipright == iproc) then
n = nn_to(ip) + 1
nn_to(ip) = n
to_p(ip)%first(n) = r_links(it, ik) + 1
to_p(ip)%second(n) = 1
to_p(ip)%third(n) = iglo_star
to_h(ip)%first(n) = 2*ncell-r_links(it, ik)
to_h(ip)%second(n) = 1
to_h(ip)%third(n) = iglo_star
end if
call get_right_connection (iglo_star, iglo_right, ipright)
iglo_star=iglo_right
end do
! v_par < 0:
call find_rightmost_link (iglo, iglo_star, ipleft)
do j = 1, ncell
! sender
if (ip == iproc) then
n = nn_from(ipleft) + 1
nn_from(ipleft) = n
from(ipleft)%first(n) = -ntgrid
from(ipleft)%second(n) = 2
from(ipleft)%third(n) = iglo
end if
! receiver
if (ipleft == iproc) then
n = nn_to(ip) + 1
nn_to(ip) = n
to_p(ip)%first(n) = l_links(it, ik) + 1
to_p(ip)%second(n) = 2
to_p(ip)%third(n) = iglo_star
to_h(ip)%first(n) = 2*ncell-l_links(it, ik)
to_h(ip)%second(n) = 2
to_h(ip)%third(n) = iglo_star
end if
call get_left_connection (iglo_star, iglo_left, ipleft)
iglo_star=iglo_left
end do
end do
from_low = [-ntgrid, 1, g_lo%llim_proc]
from_high = [ntgrid, 2, g_lo%ulim_alloc]
to_low = [1, 1, g_lo%llim_proc]
to_high = [n_links_max, 2, g_lo%ulim_alloc]
call init_fill (wfb_p, 'c', to_low, to_high, to_p, from_low, from_high, from)
call init_fill (wfb_h, 'c', to_low, to_high, to_h, from_low, from_high, from)
call delete_list (from)
call delete_list (to_p)
call delete_list (to_h)
end subroutine setup_connected_bc_redistribute_mixed_wfb