SCALE-RM
Data Types | Functions/Subroutines | Variables
scale_comm_icoa Module Reference

module COMMUNICATION More...

Functions/Subroutines

subroutine, public comm_setup
 Setup. More...
 
subroutine comm_data_transfer_sp (var, var_pl)
 Data transfer kernel. More...
 
subroutine comm_data_transfer_dp (var, var_pl)
 Data transfer kernel. More...
 
subroutine, public comm_data_transfer_nopl (var)
 Data transfer kernel. More...
 
subroutine comm_var_sp (var, var_pl, kmax, vmax)
 Data transfer with region halo => pole center. More...
 
subroutine comm_var_dp (var, var_pl, kmax, vmax)
 Data transfer with region halo => pole center. More...
 
integer function suf (i, j)
 suffix calculation More...
 
subroutine comm_stat_sum_sp (localsum, globalsum)
 
subroutine comm_stat_sum_dp (localsum, globalsum)
 
subroutine comm_stat_sum_eachlayer_sp (kall, localsum, globalsum)
 
subroutine comm_stat_sum_eachlayer_dp (kall, localsum, globalsum)
 
subroutine comm_stat_avg_sp (localavg, globalavg)
 
subroutine comm_stat_avg_dp (localavg, globalavg)
 
subroutine comm_stat_max_sp (localmax, globalmax)
 
subroutine comm_stat_max_dp (localmax, globalmax)
 
subroutine comm_stat_min_sp (localmin, globalmin)
 
subroutine comm_stat_min_dp (localmin, globalmin)
 

Variables

logical, public comm_pl = .true.
 
integer, public comm_datatype
 
integer, dimension(:), allocatable, public copy_info_r2r
 
integer, allocatable, public node
 
integer, dimension(singular point), allocatable, public information
 
integer, dimension(:,:), allocatable, public recv_info_r2r
 
integer, dimension(:,:), allocatable, public send_info_r2r
 
integer, dimension(:), allocatable, public copy_info_p2r
 
integer, dimension(:,:), allocatable, public recv_info_p2r
 
integer, dimension(:,:), allocatable, public send_info_p2r
 
integer, dimension(:), allocatable, public copy_info_r2p
 
integer, dimension(:,:), allocatable, public recv_info_r2p
 
integer, dimension(:,:), allocatable, public send_info_r2p
 
integer, dimension(:), allocatable, public singular_info
 
integer, dimension(:,:), allocatable, public copy_list_r2r
 
integer, allocatable, public grid
 
integer, allocatable, public local
 
integer, allocatable, public region
 
integer, dimension(:,:,:), allocatable, public recv_list_r2r
 
integer, dimension(:,:,:), allocatable, public send_list_r2r
 
integer, dimension(:,:), allocatable, public copy_list_p2r
 
integer, dimension(:,:,:), allocatable, public recv_list_p2r
 
integer, dimension(:,:,:), allocatable, public send_list_p2r
 
integer, dimension(:,:), allocatable, public copy_list_r2p
 
integer, dimension(:,:,:), allocatable, public recv_list_r2p
 
integer, dimension(:,:,:), allocatable, public send_list_r2p
 
integer, dimension(:,:), allocatable, public singular_list
 
integer, dimension(:), allocatable, public req_list
 
real(sp), dimension(:,:), allocatable, public sendbuf_r2r_sp
 
real(sp), dimension(:,:), allocatable, public recvbuf_r2r_sp
 
real(sp), dimension(:,:), allocatable, public sendbuf_p2r_sp
 
real(sp), dimension(:,:), allocatable, public recvbuf_p2r_sp
 
real(sp), dimension(:,:), allocatable, public sendbuf_r2p_sp
 
real(sp), dimension(:,:), allocatable, public recvbuf_r2p_sp
 
real(dp), dimension(:,:), allocatable, public sendbuf_r2r_dp
 
real(dp), dimension(:,:), allocatable, public recvbuf_r2r_dp
 
real(dp), dimension(:,:), allocatable, public sendbuf_p2r_dp
 
real(dp), dimension(:,:), allocatable, public recvbuf_p2r_dp
 
real(dp), dimension(:,:), allocatable, public sendbuf_r2p_dp
 
real(dp), dimension(:,:), allocatable, public recvbuf_r2p_dp
 

Detailed Description

module COMMUNICATION

Description
MPI Communication module for Icosahedral A-grid
Author
NICAM developers, Team SCALE
NAMELIST
  • PARAM_COMM_ICOA
    nametypedefault valuecomment
    COMM_APPLY_BARRIER logical .false. barrier option
    COMM_VARMAX integer 50 maximum number of 3D variable for the communication buffer
    DEBUG logical .false. debug option
    TESTONLY logical .false. test option

History Output
No history output

Function/Subroutine Documentation

◆ comm_setup()

subroutine, public scale_comm_icoa::comm_setup ( )

Setup.

Definition at line 198 of file scale_comm_icoA.F90.

References scale_atmos_grid_icoa_index::adm_gall, scale_atmos_grid_icoa_index::adm_gall_1d, scale_atmos_grid_icoa_index::adm_gmax, scale_atmos_grid_icoa_index::adm_gmax_pl, scale_atmos_grid_icoa_index::adm_gmin, scale_atmos_grid_icoa_index::adm_gmin_pl, scale_atmos_grid_icoa_index::adm_gslf_pl, scale_atmos_grid_icoa_index::adm_kall, scale_atmos_grid_icoa_index::adm_lall, scale_atmos_grid_icoa_index::adm_vlink, comm_datatype, comm_pl, copy_info_p2r, copy_info_r2p, copy_info_r2r, copy_list_p2r, copy_list_r2p, copy_list_r2r, scale_precision::dp, scale_prc_icoa::i_dir, scale_prc_icoa::i_e, scale_prc_icoa::i_l, scale_prc_icoa::i_n, scale_prc_icoa::i_ne, scale_prc_icoa::i_npl, scale_prc_icoa::i_nw, scale_prc_icoa::i_prc, scale_prc_icoa::i_rgnid, scale_prc_icoa::i_s, scale_prc_icoa::i_se, scale_prc_icoa::i_spl, scale_prc_icoa::i_sw, scale_prc_icoa::i_w, scale_io::io_fid_conf, scale_io::io_fid_log, scale_io::io_l, scale_io::io_nml, scale_prc::prc_abort(), scale_prc_icoa::prc_have_pl, scale_prc::prc_local_comm_world, scale_prc::prc_myrank, scale_prc::prc_nprocs, scale_prc_icoa::prc_rgn_edge_tab, scale_prc_icoa::prc_rgn_l2r, scale_prc_icoa::prc_rgn_local, scale_prc_icoa::prc_rgn_lp2r, scale_prc_icoa::prc_rgn_r2lp, scale_prc_icoa::prc_rgn_r2p_pl, scale_prc_icoa::prc_rgn_vert_num, scale_prc_icoa::prc_rgn_vert_tab, scale_prc_icoa::prc_rgn_vert_tab_pl, recv_info_p2r, recv_info_r2p, recv_info_r2r, recv_list_p2r, recv_list_r2p, recv_list_r2r, recvbuf_p2r_dp, recvbuf_p2r_sp, recvbuf_r2p_dp, recvbuf_r2p_sp, recvbuf_r2r_dp, recvbuf_r2r_sp, req_list, scale_precision::rp, send_info_p2r, send_info_r2p, send_info_r2r, send_list_p2r, send_list_r2p, send_list_r2r, sendbuf_p2r_dp, sendbuf_p2r_sp, sendbuf_r2p_dp, sendbuf_r2p_sp, sendbuf_r2r_dp, sendbuf_r2r_sp, singular_info, singular_list, scale_precision::sp, and suf().

198  use scale_prc, only: &
199  prc_abort
200  use scale_prc_icoa, only: &
201  prc_rgn_r2p_pl, &
202  i_npl, &
203  i_spl
204  implicit none
205 
206  namelist / param_comm_icoa / &
207  comm_apply_barrier, &
208  comm_varmax, &
209  debug, &
210  testonly
211 
212  integer :: ierr
213  !---------------------------------------------------------------------------
214 
215  !--- read parameters
216  if( io_l ) write(io_fid_log,*)
217  if( io_l ) write(io_fid_log,*) '+++ Module[comm]/Category[common share]'
218  rewind(io_fid_conf)
219  read(io_fid_conf,nml=param_comm_icoa,iostat=ierr)
220  if ( ierr < 0 ) then
221  if( io_l ) write(io_fid_log,*) '*** PARAM_COMM_ICOA is not specified. use default.'
222  elseif( ierr > 0 ) then
223  write(*,*) 'xxx Not appropriate names in namelist PARAM_COMM_ICOA. STOP.'
224  call prc_abort
225  endif
226  if( io_nml ) write(io_fid_log,nml=param_comm_icoa)
227 
228  if ( rp == dp ) then
229  comm_datatype = mpi_double_precision
230  elseif( rp == sp ) then
231  comm_datatype = mpi_real
232  else
233  write(*,*) 'xxx precision is not supportd'
234  call prc_abort
235  endif
236 
237  if ( prc_rgn_r2p_pl(i_npl) < 0 &
238  .AND. prc_rgn_r2p_pl(i_spl) < 0 ) then
239  comm_pl = .false.
240  endif
241 
242  if( io_l ) write(io_fid_log,*)
243  if( io_l ) write(io_fid_log,*) '====== communication information ======'
244 
245  call comm_list_generate
246 
247  call comm_sortdest
248  call comm_sortdest_pl
249  call comm_sortdest_singular
250 
251  allocate( req_list( recv_nmax_r2r + send_nmax_r2r &
252  + recv_nmax_p2r + send_nmax_p2r &
253  + recv_nmax_r2p + send_nmax_r2p ) )
254 
255  if( testonly ) call comm_debugtest
256 
257  return
integer, public comm_datatype
datatype of variable
integer, parameter, public i_spl
south pole
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
logical, public io_l
output log or not? (this process)
Definition: scale_io.F90:61
module PROCESS
Definition: scale_prc.F90:11
integer, dimension(prc_rgn_total_pl), public prc_rgn_r2p_pl
process ID which have the pole regions
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
integer, parameter, public i_npl
north pole
integer, parameter, public sp
logical, public io_nml
output log or not? (for namelist, this process)
Definition: scale_io.F90:62
integer, public io_fid_log
Log file ID.
Definition: scale_io.F90:56
integer, parameter, public rp
module process / icoA
integer, parameter, public dp
Here is the call graph for this function:

◆ comm_data_transfer_sp()

subroutine scale_comm_icoa::comm_data_transfer_sp ( real(sp), dimension (:,:,:,:), intent(inout)  var,
real(sp), dimension(:,:,:,:), intent(inout)  var_pl 
)

Data transfer kernel.

Definition at line 1558 of file scale_comm_icoA.F90.

References scale_atmos_grid_icoa_index::adm_kall, copy_info_p2r, copy_info_r2p, copy_info_r2r, copy_list_p2r, copy_list_r2p, copy_list_r2r, scale_prc::prc_abort(), scale_prc::prc_local_comm_world, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), recv_info_p2r, recv_info_r2p, recv_info_r2r, recv_list_p2r, recv_list_r2p, recv_list_r2r, recvbuf_p2r_sp, recvbuf_r2p_sp, recvbuf_r2r_sp, req_list, send_info_p2r, send_info_r2p, send_info_r2r, send_list_p2r, send_list_r2p, send_list_r2r, sendbuf_p2r_sp, sendbuf_r2p_sp, sendbuf_r2r_sp, singular_info, and singular_list.

Referenced by comm_var_sp().

1558  use scale_prc, only: &
1560  prc_abort
1561  implicit none
1562 
1563  real(SP), intent(inout) :: var (:,:,:,:)
1564  real(SP), intent(inout) :: var_pl(:,:,:,:)
1565 
1566  integer :: shp(4), kmax, vmax
1567  integer :: totalsize, rank, tag
1568  integer :: irank, ipos, imax
1569  integer :: ij_from, l_from, ij_to, l_to
1570 
1571  integer :: k, v, ikv
1572  integer :: ierr
1573  !---------------------------------------------------------------------------
1574  !$acc data &
1575  !$acc pcopy(var) &
1576  !$acc pcopy(var_pl)
1577 
1578  if ( comm_apply_barrier ) then
1579  call prof_rapstart('COMM_barrier',2)
1580  call mpi_barrier(prc_local_comm_world,ierr)
1581  call prof_rapend ('COMM_barrier',2)
1582  endif
1583 
1584  !$acc wait
1585 
1586  call prof_rapstart('COMM_data_transfer',2)
1587 
1588  shp = shape(var)
1589  kmax = shp(2)
1590  vmax = shp(4)
1591 
1592  if ( kmax * vmax > adm_kall * comm_varmax ) then
1593  write(*,*) 'xxx [COMM_data_transfer] kmax * vmax exceeds ADM_kall * COMM_varmax, stop!'
1594  write(*,*) 'xxx kmax * vmax = ', kmax * vmax
1595  write(*,*) 'xxx ADM_kall * COMM_varmax = ', adm_kall * comm_varmax
1596  call prc_abort
1597  endif
1598 
1599  !---< start communication >---
1600  ! Theres no p2r & r2p communication without calling COMM_sortdest_pl.
1601  ! receive pole => region
1602  ! receive region => pole
1603  ! receive region => region
1604  ! pack and send pole => region
1605  ! pack and send region => pole
1606  ! pack and send region => region
1607  ! copy pole => region
1608  ! copy region => pole
1609  ! copy region => region
1610  ! wait all
1611  ! unpack pole => region
1612  ! unpack region => pole
1613  ! unpack region => region
1614  ! copy region halo => region halo (singular point)
1615 
1616  req_count = 0
1617 
1618  !--- receive r2r
1619  do irank = 1, recv_nmax_r2r
1620  req_count = req_count + 1
1621  totalsize = recv_info_r2r(i_size ,irank) * kmax * vmax
1622  rank = recv_info_r2r(i_prc_from,irank)
1623  tag = recv_info_r2r(i_prc_from,irank)
1624 
1625  call mpi_irecv( recvbuf_r2r_sp(1,irank), & ! [OUT]
1626  totalsize, & ! [IN]
1627  mpi_real, & ! [IN]
1628  rank, & ! [IN]
1629  tag, & ! [IN]
1630  prc_local_comm_world, & ! [IN]
1631  req_list(req_count), & ! [OUT]
1632  ierr ) ! [OUT]
1633  enddo
1634 
1635  !--- receive p2r
1636  do irank = 1, recv_nmax_p2r
1637  req_count = req_count + 1
1638  totalsize = recv_info_p2r(i_size ,irank) * kmax * vmax
1639  rank = recv_info_p2r(i_prc_from,irank)
1640  tag = recv_info_p2r(i_prc_from,irank) + 1000000
1641 
1642  call mpi_irecv( recvbuf_p2r_sp(1,irank), & ! [OUT]
1643  totalsize, & ! [IN]
1644  mpi_real, & ! [IN]
1645  rank, & ! [IN]
1646  tag, & ! [IN]
1647  prc_local_comm_world, & ! [IN]
1648  req_list(req_count), & ! [OUT]
1649  ierr ) ! [OUT]
1650  enddo
1651 
1652  !--- receive r2p
1653  do irank = 1, recv_nmax_r2p
1654  req_count = req_count + 1
1655  totalsize = recv_info_r2p(i_size ,irank) * kmax * vmax
1656  rank = recv_info_r2p(i_prc_from,irank)
1657  tag = recv_info_r2p(i_prc_from,irank) + 2000000
1658 
1659  call mpi_irecv( recvbuf_r2p_sp(1,irank), & ! [OUT]
1660  totalsize, & ! [IN]
1661  mpi_real, & ! [IN]
1662  rank, & ! [IN]
1663  tag, & ! [IN]
1664  prc_local_comm_world, & ! [IN]
1665  req_list(req_count), & ! [OUT]
1666  ierr ) ! [OUT]
1667  enddo
1668 
1669  !$acc data &
1670  !$acc pcopyin(Send_list_r2r,Send_list_p2r,Send_list_r2p)
1671  !$acc wait
1672 
1673  !--- pack and send r2r
1674  do irank = 1, send_nmax_r2r
1675  imax = send_info_r2r(i_size,irank)
1676 
1677  !$acc kernels pcopy(sendbuf_r2r_SP)
1678  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
1679  !$omp shared(irank,imax,kmax,vmax,Send_list_r2r,sendbuf_r2r_SP,var) &
1680  !$omp collapse(2)
1681  !$acc loop independent
1682  do v = 1, vmax
1683  !$acc loop independent
1684  do k = 1, kmax
1685  !$acc loop independent
1686  do ipos = 1, imax
1687  ij_from = send_list_r2r(i_grid_from,ipos,irank)
1688  l_from = send_list_r2r(i_l_from ,ipos,irank)
1689 
1690  ikv = (v-1) * imax * kmax &
1691  + (k-1) * imax &
1692  + ipos
1693 
1694  sendbuf_r2r_sp(ikv,irank) = var(ij_from,k,l_from,v)
1695  enddo
1696  enddo
1697  enddo
1698  !$omp end parallel do
1699  !$acc end kernels
1700 
1701  req_count = req_count + 1
1702  totalsize = imax * kmax * vmax
1703  rank = send_info_r2r(i_prc_to ,irank)
1704  tag = send_info_r2r(i_prc_from,irank)
1705 
1706  !$acc wait
1707  call mpi_isend( sendbuf_r2r_sp(1,irank), & ! [IN]
1708  totalsize, & ! [IN]
1709  mpi_real, & ! [IN]
1710  rank, & ! [IN]
1711  tag, & ! [IN]
1712  prc_local_comm_world, & ! [IN]
1713  req_list(req_count), & ! [OUT]
1714  ierr ) ! [OUT]
1715  enddo
1716 
1717  !--- pack and send p2r
1718  do irank = 1, send_nmax_p2r
1719  imax = send_info_p2r(i_size,irank)
1720 
1721  !$acc kernels pcopy(sendbuf_p2r_SP)
1722  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
1723  !$omp shared(irank,imax,kmax,vmax,Send_list_p2r,sendbuf_p2r_SP,var_pl) &
1724  !$omp collapse(2)
1725  !$acc loop independent
1726  do v = 1, vmax
1727  !$acc loop independent
1728  do k = 1, kmax
1729  !$acc loop independent
1730  do ipos = 1, imax
1731  ij_from = send_list_p2r(i_grid_from,ipos,irank)
1732  l_from = send_list_p2r(i_l_from ,ipos,irank)
1733 
1734  ikv = (v-1) * imax * kmax &
1735  + (k-1) * imax &
1736  + ipos
1737 
1738  sendbuf_p2r_sp(ikv,irank) = var_pl(ij_from,k,l_from,v)
1739  enddo
1740  enddo
1741  enddo
1742  !$omp end parallel do
1743  !$acc end kernels
1744 
1745  req_count = req_count + 1
1746  totalsize = imax * kmax * vmax
1747  rank = send_info_p2r(i_prc_to ,irank)
1748  tag = send_info_p2r(i_prc_from,irank) + 1000000
1749 
1750  !$acc wait
1751  call mpi_isend( sendbuf_p2r_sp(1,irank), & ! [IN]
1752  totalsize, & ! [IN]
1753  mpi_real, & ! [IN]
1754  rank, & ! [IN]
1755  tag, & ! [IN]
1756  prc_local_comm_world, & ! [IN]
1757  req_list(req_count), & ! [OUT]
1758  ierr ) ! [OUT]
1759  enddo
1760 
1761  !--- pack and send r2p
1762  do irank = 1, send_nmax_r2p
1763  imax = send_info_r2p(i_size,irank)
1764 
1765  !$acc kernels pcopy(sendbuf_r2p_SP)
1766  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
1767  !$omp shared(irank,imax,kmax,vmax,Send_list_r2p,sendbuf_r2p_SP,var) &
1768  !$omp collapse(2)
1769  !$acc loop independent
1770  do v = 1, vmax
1771  !$acc loop independent
1772  do k = 1, kmax
1773  !$acc loop independent
1774  do ipos = 1, imax
1775  ij_from = send_list_r2p(i_grid_from,ipos,irank)
1776  l_from = send_list_r2p(i_l_from ,ipos,irank)
1777 
1778  ikv = (v-1) * imax * kmax &
1779  + (k-1) * imax &
1780  + ipos
1781 
1782  sendbuf_r2p_sp(ikv,irank) = var(ij_from,k,l_from,v)
1783  enddo
1784  enddo
1785  enddo
1786  !$omp end parallel do
1787  !$acc end kernels
1788 
1789  req_count = req_count + 1
1790  totalsize = imax * kmax * vmax
1791  rank = send_info_r2p(i_prc_to ,irank)
1792  tag = send_info_r2p(i_prc_from,irank) + 2000000
1793 
1794  !$acc wait
1795  call mpi_isend( sendbuf_r2p_sp(1,irank), & ! [IN]
1796  totalsize, & ! [IN]
1797  mpi_real, & ! [IN]
1798  rank, & ! [IN]
1799  tag, & ! [IN]
1800  prc_local_comm_world, & ! [IN]
1801  req_list(req_count), & ! [OUT]
1802  ierr ) ! [OUT]
1803  enddo
1804 
1805  !$acc end data
1806  !$acc data &
1807  !$acc pcopyin(Copy_list_r2r,Copy_list_p2r,Copy_list_r2p)
1808  !$acc wait
1809 
1810  !$omp parallel default(none),private(ipos,k,v,imax,irank,ij_from,l_from,ij_to,l_to) &
1811  !$omp shared(kmax,vmax,var,var_pl, &
1812  !$omp Copy_nmax_p2r,Copy_info_p2r,Copy_list_p2r, &
1813  !$omp Copy_nmax_r2p,Copy_info_r2p,Copy_list_r2p, &
1814  !$omp Copy_nmax_r2r,Copy_info_r2r,Copy_list_r2r )
1815 
1816  !--- copy r2r
1817  do irank = 1, copy_nmax_r2r
1818  imax = copy_info_r2r(i_size)
1819 
1820  !$acc kernels
1821  !$omp do collapse(2)
1822  !$acc loop independent
1823  do v = 1, vmax
1824  !$acc loop independent
1825  do k = 1, kmax
1826  !$acc loop independent
1827  do ipos = 1, imax
1828  ij_from = copy_list_r2r(i_grid_from,ipos)
1829  l_from = copy_list_r2r(i_l_from ,ipos)
1830  ij_to = copy_list_r2r(i_grid_to ,ipos)
1831  l_to = copy_list_r2r(i_l_to ,ipos)
1832 
1833  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
1834  enddo
1835  enddo
1836  enddo
1837  !$omp end do nowait
1838  !$acc end kernels
1839  enddo
1840 
1841  !--- copy p2r
1842  do irank = 1, copy_nmax_p2r
1843  imax = copy_info_p2r(i_size)
1844 
1845  !$acc kernels
1846  !$omp do collapse(2)
1847  !$acc loop independent
1848  do v = 1, vmax
1849  !$acc loop independent
1850  do k = 1, kmax
1851  !$acc loop independent
1852  do ipos = 1, imax
1853  ij_from = copy_list_p2r(i_grid_from,ipos)
1854  l_from = copy_list_p2r(i_l_from ,ipos)
1855  ij_to = copy_list_p2r(i_grid_to ,ipos)
1856  l_to = copy_list_p2r(i_l_to ,ipos)
1857 
1858  var(ij_to,k,l_to,v) = var_pl(ij_from,k,l_from,v)
1859  enddo
1860  enddo
1861  enddo
1862  !$omp end do nowait
1863  !$acc end kernels
1864  enddo
1865 
1866  !--- copy r2p
1867  do irank = 1, copy_nmax_r2p
1868  imax = copy_info_r2p(i_size)
1869 
1870  !$acc kernels
1871  !$omp do collapse(2)
1872  !$acc loop independent
1873  do v = 1, vmax
1874  !$acc loop independent
1875  do k = 1, kmax
1876  !$acc loop independent
1877  do ipos = 1, imax
1878  ij_from = copy_list_r2p(i_grid_from,ipos)
1879  l_from = copy_list_r2p(i_l_from ,ipos)
1880  ij_to = copy_list_r2p(i_grid_to ,ipos)
1881  l_to = copy_list_r2p(i_l_to ,ipos)
1882 
1883  var_pl(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
1884  enddo
1885  enddo
1886  enddo
1887  !$omp end do nowait
1888  !$acc end kernels
1889  enddo
1890 
1891  !$omp end parallel
1892 
1893  !--- wait all
1894  if ( req_count > 0 ) then
1895  call mpi_waitall( req_count, & ! [IN]
1896  req_list(:), & ! [IN]
1897  mpi_statuses_ignore, & ! [OUT]
1898  ierr ) ! [OUT]
1899  endif
1900 
1901  !$acc end data
1902  !$acc data &
1903  !$acc pcopyin(Recv_list_r2r,Recv_list_p2r,Recv_list_r2p)
1904  !$acc wait
1905 
1906  !$omp parallel default(none),private(ipos,k,v,imax,irank,ikv,ij_from,l_from,ij_to,l_to) &
1907  !$omp shared(kmax,vmax,var,var_pl, &
1908  !$omp Recv_nmax_p2r,Recv_info_p2r,Recv_list_p2r,recvbuf_p2r_SP, &
1909  !$omp Recv_nmax_r2p,Recv_info_r2p,Recv_list_r2p,recvbuf_r2p_SP, &
1910  !$omp Recv_nmax_r2r,Recv_info_r2r,Recv_list_r2r,recvbuf_r2r_SP, &
1911  !$omp Singular_nmax,Singular_info,Singular_list )
1912 
1913  !--- unpack r2r
1914  do irank = 1, recv_nmax_r2r
1915  imax = recv_info_r2r(i_size,irank)
1916 
1917  !$acc kernels pcopyin(recvbuf_r2r_SP)
1918  !$omp do collapse(2)
1919  !$acc loop independent
1920  do v = 1, vmax
1921  !$acc loop independent
1922  do k = 1, kmax
1923  !$acc loop independent
1924  do ipos = 1, imax
1925  ij_to = recv_list_r2r(i_grid_to,ipos,irank)
1926  l_to = recv_list_r2r(i_l_to ,ipos,irank)
1927 
1928  ikv = (v-1) * imax * kmax &
1929  + (k-1) * imax &
1930  + ipos
1931 
1932  var(ij_to,k,l_to,v) = recvbuf_r2r_sp(ikv,irank)
1933  enddo
1934  enddo
1935  enddo
1936  !$omp end do nowait
1937  !$acc end kernels
1938  enddo
1939 
1940  !--- unpack p2r
1941  do irank = 1, recv_nmax_p2r
1942  imax = recv_info_p2r(i_size,irank)
1943 
1944  !$acc kernels pcopyin(recvbuf_p2r_SP)
1945  !$omp do collapse(2)
1946  !$acc loop independent
1947  do v = 1, vmax
1948  !$acc loop independent
1949  do k = 1, kmax
1950  !$acc loop independent
1951  do ipos = 1, imax
1952  ij_to = recv_list_p2r(i_grid_to,ipos,irank)
1953  l_to = recv_list_p2r(i_l_to ,ipos,irank)
1954 
1955  ikv = (v-1) * imax * kmax &
1956  + (k-1) * imax &
1957  + ipos
1958 
1959  var(ij_to,k,l_to,v) = recvbuf_p2r_sp(ikv,irank)
1960  enddo
1961  enddo
1962  enddo
1963  !$omp end do nowait
1964  !$acc end kernels
1965  enddo
1966 
1967  !--- unpack r2p
1968  do irank = 1, recv_nmax_r2p
1969  imax = recv_info_r2p(i_size,irank)
1970 
1971  !$acc kernels pcopyin(recvbuf_r2p_SP)
1972  !$omp do collapse(2)
1973  !$acc loop independent
1974  do v = 1, vmax
1975  !$acc loop independent
1976  do k = 1, kmax
1977  !$acc loop independent
1978  do ipos = 1, imax
1979  ij_to = recv_list_r2p(i_grid_to,ipos,irank)
1980  l_to = recv_list_r2p(i_l_to ,ipos,irank)
1981 
1982  ikv = (v-1) * imax * kmax &
1983  + (k-1) * imax &
1984  + ipos
1985 
1986  var_pl(ij_to,k,l_to,v) = recvbuf_r2p_sp(ikv,irank)
1987  enddo
1988  enddo
1989  enddo
1990  !$omp end do nowait
1991  !$acc end kernels
1992  enddo
1993 
1994  !$acc end data
1995  !$acc wait
1996 
1997  !--- singular point (halo to halo)
1998  do irank = 1, singular_nmax
1999  imax = singular_info(i_size)
2000 
2001  !$acc kernels pcopyin(Singular_list)
2002  !$omp do collapse(2)
2003  !$acc loop independent
2004  do v = 1, vmax
2005  !$acc loop independent
2006  do k = 1, kmax
2007  !$acc loop independent
2008  do ipos = 1, imax
2009  ij_from = singular_list(i_grid_from,ipos)
2010  l_from = singular_list(i_l_from ,ipos)
2011  ij_to = singular_list(i_grid_to ,ipos)
2012  l_to = singular_list(i_l_to ,ipos)
2013 
2014  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2015  enddo
2016  enddo
2017  enddo
2018  !$omp end do nowait
2019  !$acc end kernels
2020  enddo
2021 
2022  !$omp end parallel
2023 
2024  call prof_rapend('COMM_data_transfer',2)
2025 
2026  !$acc end data
2027  !$acc wait
2028 
2029  return
integer, public imax
of computational cells: x, local
module PROCESS
Definition: scale_prc.F90:11
integer, public kmax
of computational cells: z, local
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87
Here is the call graph for this function:
Here is the caller graph for this function:

◆ comm_data_transfer_dp()

subroutine scale_comm_icoa::comm_data_transfer_dp ( real(dp), dimension (:,:,:,:), intent(inout)  var,
real(dp), dimension(:,:,:,:), intent(inout)  var_pl 
)

Data transfer kernel.

Definition at line 2037 of file scale_comm_icoA.F90.

References scale_atmos_grid_icoa_index::adm_kall, copy_info_p2r, copy_info_r2p, copy_info_r2r, copy_list_p2r, copy_list_r2p, copy_list_r2r, scale_prc::prc_abort(), scale_prc::prc_local_comm_world, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), recv_info_p2r, recv_info_r2p, recv_info_r2r, recv_list_p2r, recv_list_r2p, recv_list_r2r, recvbuf_p2r_dp, recvbuf_r2p_dp, recvbuf_r2r_dp, req_list, send_info_p2r, send_info_r2p, send_info_r2r, send_list_p2r, send_list_r2p, send_list_r2r, sendbuf_p2r_dp, sendbuf_r2p_dp, sendbuf_r2r_dp, singular_info, and singular_list.

Referenced by comm_var_dp().

2037  use scale_prc, only: &
2039  prc_abort
2040  implicit none
2041 
2042  real(DP), intent(inout) :: var (:,:,:,:)
2043  real(DP), intent(inout) :: var_pl(:,:,:,:)
2044 
2045  integer :: shp(4), kmax, vmax
2046  integer :: totalsize, rank, tag
2047  integer :: irank, ipos, imax
2048  integer :: ij_from, l_from, ij_to, l_to
2049 
2050  integer :: k, v, ikv
2051  integer :: ierr
2052  !---------------------------------------------------------------------------
2053  !$acc data &
2054  !$acc pcopy(var) &
2055  !$acc pcopy(var_pl)
2056 
2057  if ( comm_apply_barrier ) then
2058  call prof_rapstart('COMM_barrier',2)
2059  call mpi_barrier(prc_local_comm_world,ierr)
2060  call prof_rapend ('COMM_barrier',2)
2061  endif
2062 
2063  !$acc wait
2064 
2065  call prof_rapstart('COMM_data_transfer',2)
2066 
2067  shp = shape(var)
2068  kmax = shp(2)
2069  vmax = shp(4)
2070 
2071  if ( kmax * vmax > adm_kall * comm_varmax ) then
2072  write(*,*) 'xxx [COMM_data_transfer] kmax * vmax exceeds ADM_kall * COMM_varmax, stop!'
2073  write(*,*) 'xxx kmax * vmax = ', kmax * vmax
2074  write(*,*) 'xxx ADM_kall * COMM_varmax = ', adm_kall * comm_varmax
2075  call prc_abort
2076  endif
2077 
2078  !---< start communication >---
2079  ! Theres no p2r & r2p communication without calling COMM_sortdest_pl.
2080  ! receive pole => region
2081  ! receive region => pole
2082  ! receive region => region
2083  ! pack and send pole => region
2084  ! pack and send region => pole
2085  ! pack and send region => region
2086  ! copy pole => region
2087  ! copy region => pole
2088  ! copy region => region
2089  ! wait all
2090  ! unpack pole => region
2091  ! unpack region => pole
2092  ! unpack region => region
2093  ! copy region halo => region halo (singular point)
2094 
2095  req_count = 0
2096 
2097  !--- receive r2r
2098  do irank = 1, recv_nmax_r2r
2099  req_count = req_count + 1
2100  totalsize = recv_info_r2r(i_size ,irank) * kmax * vmax
2101  rank = recv_info_r2r(i_prc_from,irank)
2102  tag = recv_info_r2r(i_prc_from,irank)
2103 
2104  call mpi_irecv( recvbuf_r2r_dp(1,irank), & ! [OUT]
2105  totalsize, & ! [IN]
2106  mpi_double_precision, & ! [IN]
2107  rank, & ! [IN]
2108  tag, & ! [IN]
2109  prc_local_comm_world, & ! [IN]
2110  req_list(req_count), & ! [OUT]
2111  ierr ) ! [OUT]
2112  enddo
2113 
2114  !--- receive p2r
2115  do irank = 1, recv_nmax_p2r
2116  req_count = req_count + 1
2117  totalsize = recv_info_p2r(i_size ,irank) * kmax * vmax
2118  rank = recv_info_p2r(i_prc_from,irank)
2119  tag = recv_info_p2r(i_prc_from,irank) + 1000000
2120 
2121  call mpi_irecv( recvbuf_p2r_dp(1,irank), & ! [OUT]
2122  totalsize, & ! [IN]
2123  mpi_double_precision, & ! [IN]
2124  rank, & ! [IN]
2125  tag, & ! [IN]
2126  prc_local_comm_world, & ! [IN]
2127  req_list(req_count), & ! [OUT]
2128  ierr ) ! [OUT]
2129  enddo
2130 
2131  !--- receive r2p
2132  do irank = 1, recv_nmax_r2p
2133  req_count = req_count + 1
2134  totalsize = recv_info_r2p(i_size ,irank) * kmax * vmax
2135  rank = recv_info_r2p(i_prc_from,irank)
2136  tag = recv_info_r2p(i_prc_from,irank) + 2000000
2137 
2138  call mpi_irecv( recvbuf_r2p_dp(1,irank), & ! [OUT]
2139  totalsize, & ! [IN]
2140  mpi_double_precision, & ! [IN]
2141  rank, & ! [IN]
2142  tag, & ! [IN]
2143  prc_local_comm_world, & ! [IN]
2144  req_list(req_count), & ! [OUT]
2145  ierr ) ! [OUT]
2146  enddo
2147 
2148  !$acc data &
2149  !$acc pcopyin(Send_list_r2r,Send_list_p2r,Send_list_r2p)
2150  !$acc wait
2151 
2152  !--- pack and send r2r
2153  do irank = 1, send_nmax_r2r
2154  imax = send_info_r2r(i_size,irank)
2155 
2156  !$acc kernels pcopy(sendbuf_r2r_DP)
2157  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
2158  !$omp shared(irank,imax,kmax,vmax,Send_list_r2r,sendbuf_r2r_DP,var) &
2159  !$omp collapse(2)
2160  !$acc loop independent
2161  do v = 1, vmax
2162  !$acc loop independent
2163  do k = 1, kmax
2164  !$acc loop independent
2165  do ipos = 1, imax
2166  ij_from = send_list_r2r(i_grid_from,ipos,irank)
2167  l_from = send_list_r2r(i_l_from ,ipos,irank)
2168 
2169  ikv = (v-1) * imax * kmax &
2170  + (k-1) * imax &
2171  + ipos
2172 
2173  sendbuf_r2r_dp(ikv,irank) = var(ij_from,k,l_from,v)
2174  enddo
2175  enddo
2176  enddo
2177  !$omp end parallel do
2178  !$acc end kernels
2179 
2180  req_count = req_count + 1
2181  totalsize = imax * kmax * vmax
2182  rank = send_info_r2r(i_prc_to ,irank)
2183  tag = send_info_r2r(i_prc_from,irank)
2184 
2185  !$acc wait
2186  call mpi_isend( sendbuf_r2r_dp(1,irank), & ! [IN]
2187  totalsize, & ! [IN]
2188  mpi_double_precision, & ! [IN]
2189  rank, & ! [IN]
2190  tag, & ! [IN]
2191  prc_local_comm_world, & ! [IN]
2192  req_list(req_count), & ! [OUT]
2193  ierr ) ! [OUT]
2194  enddo
2195 
2196  !--- pack and send p2r
2197  do irank = 1, send_nmax_p2r
2198  imax = send_info_p2r(i_size,irank)
2199 
2200  !$acc kernels pcopy(sendbuf_p2r_DP)
2201  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
2202  !$omp shared(irank,imax,kmax,vmax,Send_list_p2r,sendbuf_p2r_DP,var_pl) &
2203  !$omp collapse(2)
2204  !$acc loop independent
2205  do v = 1, vmax
2206  !$acc loop independent
2207  do k = 1, kmax
2208  !$acc loop independent
2209  do ipos = 1, imax
2210  ij_from = send_list_p2r(i_grid_from,ipos,irank)
2211  l_from = send_list_p2r(i_l_from ,ipos,irank)
2212 
2213  ikv = (v-1) * imax * kmax &
2214  + (k-1) * imax &
2215  + ipos
2216 
2217  sendbuf_p2r_dp(ikv,irank) = var_pl(ij_from,k,l_from,v)
2218  enddo
2219  enddo
2220  enddo
2221  !$omp end parallel do
2222  !$acc end kernels
2223 
2224  req_count = req_count + 1
2225  totalsize = imax * kmax * vmax
2226  rank = send_info_p2r(i_prc_to ,irank)
2227  tag = send_info_p2r(i_prc_from,irank) + 1000000
2228 
2229  !$acc wait
2230  call mpi_isend( sendbuf_p2r_dp(1,irank), & ! [IN]
2231  totalsize, & ! [IN]
2232  mpi_double_precision, & ! [IN]
2233  rank, & ! [IN]
2234  tag, & ! [IN]
2235  prc_local_comm_world, & ! [IN]
2236  req_list(req_count), & ! [OUT]
2237  ierr ) ! [OUT]
2238  enddo
2239 
2240  !--- pack and send r2p
2241  do irank = 1, send_nmax_r2p
2242  imax = send_info_r2p(i_size,irank)
2243 
2244  !$acc kernels pcopy(sendbuf_r2p_DP)
2245  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
2246  !$omp shared(irank,imax,kmax,vmax,Send_list_r2p,sendbuf_r2p_DP,var) &
2247  !$omp collapse(2)
2248  !$acc loop independent
2249  do v = 1, vmax
2250  !$acc loop independent
2251  do k = 1, kmax
2252  !$acc loop independent
2253  do ipos = 1, imax
2254  ij_from = send_list_r2p(i_grid_from,ipos,irank)
2255  l_from = send_list_r2p(i_l_from ,ipos,irank)
2256 
2257  ikv = (v-1) * imax * kmax &
2258  + (k-1) * imax &
2259  + ipos
2260 
2261  sendbuf_r2p_dp(ikv,irank) = var(ij_from,k,l_from,v)
2262  enddo
2263  enddo
2264  enddo
2265  !$omp end parallel do
2266  !$acc end kernels
2267 
2268  req_count = req_count + 1
2269  totalsize = imax * kmax * vmax
2270  rank = send_info_r2p(i_prc_to ,irank)
2271  tag = send_info_r2p(i_prc_from,irank) + 2000000
2272 
2273  !$acc wait
2274  call mpi_isend( sendbuf_r2p_dp(1,irank), & ! [IN]
2275  totalsize, & ! [IN]
2276  mpi_double_precision, & ! [IN]
2277  rank, & ! [IN]
2278  tag, & ! [IN]
2279  prc_local_comm_world, & ! [IN]
2280  req_list(req_count), & ! [OUT]
2281  ierr ) ! [OUT]
2282  enddo
2283 
2284  !$acc end data
2285  !$acc data &
2286  !$acc pcopyin(Copy_list_r2r,Copy_list_p2r,Copy_list_r2p)
2287  !$acc wait
2288 
2289  !$omp parallel default(none),private(ipos,k,v,imax,irank,ij_from,l_from,ij_to,l_to) &
2290  !$omp shared(kmax,vmax,var,var_pl, &
2291  !$omp Copy_nmax_p2r,Copy_info_p2r,Copy_list_p2r, &
2292  !$omp Copy_nmax_r2p,Copy_info_r2p,Copy_list_r2p, &
2293  !$omp Copy_nmax_r2r,Copy_info_r2r,Copy_list_r2r )
2294 
2295  !--- copy r2r
2296  do irank = 1, copy_nmax_r2r
2297  imax = copy_info_r2r(i_size)
2298 
2299  !$acc kernels
2300  !$omp do collapse(2)
2301  !$acc loop independent
2302  do v = 1, vmax
2303  !$acc loop independent
2304  do k = 1, kmax
2305  !$acc loop independent
2306  do ipos = 1, imax
2307  ij_from = copy_list_r2r(i_grid_from,ipos)
2308  l_from = copy_list_r2r(i_l_from ,ipos)
2309  ij_to = copy_list_r2r(i_grid_to ,ipos)
2310  l_to = copy_list_r2r(i_l_to ,ipos)
2311 
2312  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2313  enddo
2314  enddo
2315  enddo
2316  !$omp end do nowait
2317  !$acc end kernels
2318  enddo
2319 
2320  !--- copy p2r
2321  do irank = 1, copy_nmax_p2r
2322  imax = copy_info_p2r(i_size)
2323 
2324  !$acc kernels
2325  !$omp do collapse(2)
2326  !$acc loop independent
2327  do v = 1, vmax
2328  !$acc loop independent
2329  do k = 1, kmax
2330  !$acc loop independent
2331  do ipos = 1, imax
2332  ij_from = copy_list_p2r(i_grid_from,ipos)
2333  l_from = copy_list_p2r(i_l_from ,ipos)
2334  ij_to = copy_list_p2r(i_grid_to ,ipos)
2335  l_to = copy_list_p2r(i_l_to ,ipos)
2336 
2337  var(ij_to,k,l_to,v) = var_pl(ij_from,k,l_from,v)
2338  enddo
2339  enddo
2340  enddo
2341  !$omp end do nowait
2342  !$acc end kernels
2343  enddo
2344 
2345  !--- copy r2p
2346  do irank = 1, copy_nmax_r2p
2347  imax = copy_info_r2p(i_size)
2348 
2349  !$acc kernels
2350  !$omp do collapse(2)
2351  !$acc loop independent
2352  do v = 1, vmax
2353  !$acc loop independent
2354  do k = 1, kmax
2355  !$acc loop independent
2356  do ipos = 1, imax
2357  ij_from = copy_list_r2p(i_grid_from,ipos)
2358  l_from = copy_list_r2p(i_l_from ,ipos)
2359  ij_to = copy_list_r2p(i_grid_to ,ipos)
2360  l_to = copy_list_r2p(i_l_to ,ipos)
2361 
2362  var_pl(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2363  enddo
2364  enddo
2365  enddo
2366  !$omp end do nowait
2367  !$acc end kernels
2368  enddo
2369 
2370  !$omp end parallel
2371 
2372  !--- wait all
2373  if ( req_count > 0 ) then
2374  call mpi_waitall( req_count, & ! [IN]
2375  req_list(:), & ! [IN]
2376  mpi_statuses_ignore, & ! [OUT]
2377  ierr ) ! [OUT]
2378  endif
2379 
2380  !$acc end data
2381  !$acc data &
2382  !$acc pcopyin(Recv_list_r2r,Recv_list_p2r,Recv_list_r2p)
2383  !$acc wait
2384 
2385  !$omp parallel default(none),private(ipos,k,v,imax,irank,ikv,ij_from,l_from,ij_to,l_to) &
2386  !$omp shared(kmax,vmax,var,var_pl, &
2387  !$omp Recv_nmax_p2r,Recv_info_p2r,Recv_list_p2r,recvbuf_p2r_DP, &
2388  !$omp Recv_nmax_r2p,Recv_info_r2p,Recv_list_r2p,recvbuf_r2p_DP, &
2389  !$omp Recv_nmax_r2r,Recv_info_r2r,Recv_list_r2r,recvbuf_r2r_DP, &
2390  !$omp Singular_nmax,Singular_info,Singular_list )
2391 
2392  !--- unpack r2r
2393  do irank = 1, recv_nmax_r2r
2394  imax = recv_info_r2r(i_size,irank)
2395 
2396  !$acc kernels pcopyin(recvbuf_r2r_DP)
2397  !$omp do collapse(2)
2398  !$acc loop independent
2399  do v = 1, vmax
2400  !$acc loop independent
2401  do k = 1, kmax
2402  !$acc loop independent
2403  do ipos = 1, imax
2404  ij_to = recv_list_r2r(i_grid_to,ipos,irank)
2405  l_to = recv_list_r2r(i_l_to ,ipos,irank)
2406 
2407  ikv = (v-1) * imax * kmax &
2408  + (k-1) * imax &
2409  + ipos
2410 
2411  var(ij_to,k,l_to,v) = recvbuf_r2r_dp(ikv,irank)
2412  enddo
2413  enddo
2414  enddo
2415  !$omp end do nowait
2416  !$acc end kernels
2417  enddo
2418 
2419  !--- unpack p2r
2420  do irank = 1, recv_nmax_p2r
2421  imax = recv_info_p2r(i_size,irank)
2422 
2423  !$acc kernels pcopyin(recvbuf_p2r_DP)
2424  !$omp do collapse(2)
2425  !$acc loop independent
2426  do v = 1, vmax
2427  !$acc loop independent
2428  do k = 1, kmax
2429  !$acc loop independent
2430  do ipos = 1, imax
2431  ij_to = recv_list_p2r(i_grid_to,ipos,irank)
2432  l_to = recv_list_p2r(i_l_to ,ipos,irank)
2433 
2434  ikv = (v-1) * imax * kmax &
2435  + (k-1) * imax &
2436  + ipos
2437 
2438  var(ij_to,k,l_to,v) = recvbuf_p2r_dp(ikv,irank)
2439  enddo
2440  enddo
2441  enddo
2442  !$omp end do nowait
2443  !$acc end kernels
2444  enddo
2445 
2446  !--- unpack r2p
2447  do irank = 1, recv_nmax_r2p
2448  imax = recv_info_r2p(i_size,irank)
2449 
2450  !$acc kernels pcopyin(recvbuf_r2p_DP)
2451  !$omp do collapse(2)
2452  !$acc loop independent
2453  do v = 1, vmax
2454  !$acc loop independent
2455  do k = 1, kmax
2456  !$acc loop independent
2457  do ipos = 1, imax
2458  ij_to = recv_list_r2p(i_grid_to,ipos,irank)
2459  l_to = recv_list_r2p(i_l_to ,ipos,irank)
2460 
2461  ikv = (v-1) * imax * kmax &
2462  + (k-1) * imax &
2463  + ipos
2464 
2465  var_pl(ij_to,k,l_to,v) = recvbuf_r2p_dp(ikv,irank)
2466  enddo
2467  enddo
2468  enddo
2469  !$omp end do nowait
2470  !$acc end kernels
2471  enddo
2472 
2473  !$acc end data
2474  !$acc wait
2475 
2476  !--- singular point (halo to halo)
2477  do irank = 1, singular_nmax
2478  imax = singular_info(i_size)
2479 
2480  !$acc kernels pcopyin(Singular_list)
2481  !$omp do collapse(2)
2482  !$acc loop independent
2483  do v = 1, vmax
2484  !$acc loop independent
2485  do k = 1, kmax
2486  !$acc loop independent
2487  do ipos = 1, imax
2488  ij_from = singular_list(i_grid_from,ipos)
2489  l_from = singular_list(i_l_from ,ipos)
2490  ij_to = singular_list(i_grid_to ,ipos)
2491  l_to = singular_list(i_l_to ,ipos)
2492 
2493  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2494  enddo
2495  enddo
2496  enddo
2497  !$omp end do nowait
2498  !$acc end kernels
2499  enddo
2500 
2501  !$omp end parallel
2502 
2503  call prof_rapend('COMM_data_transfer',2)
2504 
2505  !$acc end data
2506  !$acc wait
2507 
2508  return
integer, public imax
of computational cells: x, local
module PROCESS
Definition: scale_prc.F90:11
integer, public kmax
of computational cells: z, local
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87
Here is the call graph for this function:
Here is the caller graph for this function:

◆ comm_data_transfer_nopl()

subroutine, public scale_comm_icoa::comm_data_transfer_nopl ( real(dp), dimension (:,:,:,:), intent(inout)  var)

Data transfer kernel.

Definition at line 2515 of file scale_comm_icoA.F90.

References scale_atmos_grid_icoa_index::adm_kall, copy_info_r2r, copy_list_r2r, scale_prc::prc_abort(), scale_prc::prc_local_comm_world, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), recv_info_r2r, recv_list_r2r, recvbuf_r2r_dp, req_list, send_info_r2r, send_list_r2r, sendbuf_r2r_dp, singular_info, and singular_list.

2515  use scale_prc, only: &
2517  prc_abort
2518  implicit none
2519 
2520  real(DP), intent(inout) :: var (:,:,:,:)
2521 
2522  integer :: shp(4), kmax, vmax
2523  integer :: totalsize, rank, tag
2524  integer :: irank, ipos, imax
2525  integer :: ij_from, l_from, ij_to, l_to
2526 
2527  integer :: k, v, ikv
2528  integer :: ierr
2529  !---------------------------------------------------------------------------
2530 
2531  if ( comm_apply_barrier ) then
2532  call prof_rapstart('COMM_barrier',2)
2533  call mpi_barrier(prc_local_comm_world,ierr)
2534  call prof_rapend ('COMM_barrier',2)
2535  endif
2536 
2537  call prof_rapstart('COMM_data_transfer',2)
2538 
2539  shp = shape(var)
2540  kmax = shp(2)
2541  vmax = shp(4)
2542 
2543  if ( kmax * vmax > adm_kall * comm_varmax ) then
2544  write(*,*) 'xxx [COMM_data_transfer_nopl] kmax * vmax exceeds ADM_kall * COMM_varmax, stop!'
2545  write(*,*) 'xxx kmax * vmax = ', kmax * vmax
2546  write(*,*) 'xxx ADM_kall * COMM_varmax = ', adm_kall * comm_varmax
2547  call prc_abort
2548  endif
2549 
2550  !---< start communication >---
2551  ! receive region => region
2552  ! pack and send region => region
2553  ! copy region => region
2554  ! wait all
2555  ! unpack region => region
2556  ! copy region halo => region halo (singular point)
2557 
2558  req_count = 0
2559 
2560  !--- receive r2r
2561  do irank = 1, recv_nmax_r2r
2562  req_count = req_count + 1
2563  totalsize = recv_info_r2r(i_size ,irank) * kmax * vmax
2564  rank = recv_info_r2r(i_prc_from,irank)
2565  tag = recv_info_r2r(i_prc_from,irank)
2566 
2567  call mpi_irecv( recvbuf_r2r_dp(1,irank), & ! [OUT]
2568  totalsize, & ! [IN]
2569  mpi_double_precision, & ! [IN]
2570  rank, & ! [IN]
2571  tag, & ! [IN]
2572  prc_local_comm_world, & ! [IN]
2573  req_list(req_count), & ! [OUT]
2574  ierr ) ! [OUT]
2575  enddo
2576 
2577  !--- pack and send r2r
2578  do irank = 1, send_nmax_r2r
2579  imax = send_info_r2r(i_size,irank)
2580 
2581  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
2582  !$omp shared(irank,imax,kmax,vmax,Send_list_r2r,sendbuf_r2r_DP,var)
2583  do v = 1, vmax
2584  do k = 1, kmax
2585  do ipos = 1, imax
2586  ij_from = send_list_r2r(i_grid_from,ipos,irank)
2587  l_from = send_list_r2r(i_l_from ,ipos,irank)
2588 
2589  ikv = (v-1) * imax * kmax &
2590  + (k-1) * imax &
2591  + ipos
2592 
2593  sendbuf_r2r_dp(ikv,irank) = var(ij_from,k,l_from,v)
2594  enddo
2595  enddo
2596  enddo
2597  !$omp end parallel do
2598 
2599  req_count = req_count + 1
2600  totalsize = imax * kmax * vmax
2601  rank = send_info_r2r(i_prc_to ,irank)
2602  tag = send_info_r2r(i_prc_from,irank)
2603 
2604  call mpi_isend( sendbuf_r2r_dp(1,irank), & ! [IN]
2605  totalsize, & ! [IN]
2606  mpi_double_precision, & ! [IN]
2607  rank, & ! [IN]
2608  tag, & ! [IN]
2609  prc_local_comm_world, & ! [IN]
2610  req_list(req_count), & ! [OUT]
2611  ierr ) ! [OUT]
2612  enddo
2613 
2614  !--- copy r2r
2615  !$omp parallel do default(none),private(ipos,k,v,imax,irank,ij_from,l_from,ij_to,l_to) &
2616  !$omp shared(kmax,vmax,Copy_nmax_r2r,Copy_info_r2r,Copy_list_r2r,var)
2617  do irank = 1, copy_nmax_r2r
2618  imax = copy_info_r2r(i_size)
2619 
2620  do v = 1, vmax
2621  do k = 1, kmax
2622  do ipos = 1, imax
2623  ij_from = copy_list_r2r(i_grid_from,ipos)
2624  l_from = copy_list_r2r(i_l_from ,ipos)
2625  ij_to = copy_list_r2r(i_grid_to ,ipos)
2626  l_to = copy_list_r2r(i_l_to ,ipos)
2627 
2628  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2629  enddo
2630  enddo
2631  enddo
2632  enddo
2633  !$omp end parallel do
2634 
2635  !--- wait all
2636  if ( req_count > 0 ) then
2637  call mpi_waitall( req_count, & ! [IN]
2638  req_list(:), & ! [IN]
2639  mpi_statuses_ignore, & ! [OUT]
2640  ierr ) ! [OUT]
2641  endif
2642 
2643  !--- unpack r2r
2644  !$omp parallel do default(none),private(ipos,k,v,imax,irank,ij_to,l_to,ikv) &
2645  !$omp shared(kmax,vmax,Recv_nmax_r2r,Recv_info_r2r,Recv_list_r2r,recvbuf_r2r_DP,var)
2646  do irank = 1, recv_nmax_r2r
2647  imax = recv_info_r2r(i_size,irank)
2648 
2649  do v = 1, vmax
2650  do k = 1, kmax
2651  do ipos = 1, imax
2652  ij_to = recv_list_r2r(i_grid_to,ipos,irank)
2653  l_to = recv_list_r2r(i_l_to ,ipos,irank)
2654 
2655  ikv = (v-1) * imax * kmax &
2656  + (k-1) * imax &
2657  + ipos
2658 
2659  var(ij_to,k,l_to,v) = recvbuf_r2r_dp(ikv,irank)
2660  enddo
2661  enddo
2662  enddo
2663  enddo
2664  !$omp end parallel do
2665 
2666  !--- singular point (halo to halo)
2667  do irank = 1, singular_nmax
2668  imax = singular_info(i_size)
2669 
2670  do v = 1, vmax
2671  do k = 1, kmax
2672  do ipos = 1, imax
2673  ij_from = singular_list(i_grid_from,ipos)
2674  l_from = singular_list(i_l_from ,ipos)
2675  ij_to = singular_list(i_grid_to ,ipos)
2676  l_to = singular_list(i_l_to ,ipos)
2677 
2678  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2679  enddo
2680  enddo
2681  enddo
2682  enddo
2683 
2684  call prof_rapend('COMM_data_transfer',2)
2685 
2686  return
integer, public imax
of computational cells: x, local
module PROCESS
Definition: scale_prc.F90:11
integer, public kmax
of computational cells: z, local
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87
Here is the call graph for this function:

◆ comm_var_sp()

subroutine scale_comm_icoa::comm_var_sp ( real(sp), dimension (adm_gall, kmax,adm_lall, vmax), intent(inout)  var,
real(sp), dimension(adm_gall_pl,kmax,adm_lall_pl,vmax), intent(inout)  var_pl,
integer, intent(in)  kmax,
integer, intent(in)  vmax 
)

Data transfer with region halo => pole center.

Definition at line 2696 of file scale_comm_icoA.F90.

References comm_data_transfer_sp(), comm_pl, copy_info_p2r, copy_list_p2r, scale_prc_icoa::i_npl, scale_prc_icoa::i_spl, scale_prc::prc_local_comm_world, scale_prc_icoa::prc_rgn_lp2r, scale_prc_icoa::prc_rgn_rgn4pl, scale_prc_icoa::prc_rgn_total_pl, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), recv_info_p2r, recv_list_p2r, req_list, send_info_p2r, and send_list_p2r.

2696  use scale_prc, only: &
2698  use scale_prc_icoa, only: &
2699  prc_rgn_total_pl, &
2700  prc_rgn_rgn4pl, &
2701  prc_rgn_lp2r, &
2702  i_npl, &
2703  i_spl
2704  implicit none
2705 
2706  integer, intent(in) :: kmax
2707  integer, intent(in) :: vmax
2708  real(SP), intent(inout) :: var (adm_gall, kmax,adm_lall, vmax)
2709  real(SP), intent(inout) :: var_pl(adm_gall_pl,kmax,adm_lall_pl,vmax)
2710 
2711  real(SP) :: sendbuf_h2p_sp(kmax*vmax,prc_rgn_total_pl)
2712  real(SP) :: recvbuf_h2p_sp(kmax*vmax,prc_rgn_total_pl)
2713 
2714  integer :: totalsize, rank, tag
2715  integer :: irank, ipos
2716  integer :: ij_from, l_from, ij_to, l_to
2717  integer :: r_from, r_to
2718 
2719  integer :: k, v, kk
2720  integer :: ierr
2721  !---------------------------------------------------------------------------
2722 
2723  if ( comm_apply_barrier ) then
2724  call prof_rapstart('COMM_barrier',2)
2725  call mpi_barrier(prc_local_comm_world,ierr)
2726  call prof_rapend ('COMM_barrier',2)
2727  endif
2728 
2729  call prof_rapstart('COMM_var',2)
2730 
2731  if ( comm_pl ) then
2732 
2733  req_count = 0
2734 
2735  !--- receive p2r-reverse
2736  do irank = 1, send_nmax_p2r
2737  do ipos = 1, send_info_p2r(i_size,irank)
2738  l_from = send_list_p2r(i_l_to ,ipos,irank)
2739  r_from = prc_rgn_lp2r(l_from,send_info_p2r(i_prc_to,irank))
2740 
2741  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
2742  req_count = req_count + 1
2743  totalsize = kmax * vmax
2744  rank = send_info_p2r(i_prc_to ,irank)
2745  tag = send_info_p2r(i_prc_from,irank) + 1000000
2746 
2747  call mpi_irecv( recvbuf_h2p_sp(1,i_npl), & ! [OUT]
2748  totalsize, & ! [IN]
2749  mpi_real, & ! [IN]
2750  rank, & ! [IN]
2751  tag, & ! [IN]
2752  prc_local_comm_world, & ! [IN]
2753  req_list(req_count), & ! [OUT]
2754  ierr ) ! [OUT]
2755  endif
2756 
2757  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
2758  req_count = req_count + 1
2759  totalsize = kmax * vmax
2760  rank = send_info_p2r(i_prc_to ,irank)
2761  tag = send_info_p2r(i_prc_from,irank) + 2000000
2762 
2763  call mpi_irecv( recvbuf_h2p_sp(1,i_spl), & ! [OUT]
2764  totalsize, & ! [IN]
2765  mpi_real, & ! [IN]
2766  rank, & ! [IN]
2767  tag, & ! [IN]
2768  prc_local_comm_world, & ! [IN]
2769  req_list(req_count), & ! [OUT]
2770  ierr ) ! [OUT]
2771  endif
2772  enddo
2773  enddo
2774 
2775  !--- pack and send p2r-reverse
2776  do irank = 1, recv_nmax_p2r
2777  do ipos = 1, recv_info_p2r(i_size,irank)
2778  ij_from = recv_list_p2r(i_grid_to,ipos,irank)
2779  l_from = recv_list_p2r(i_l_to ,ipos,irank)
2780  r_from = prc_rgn_lp2r(l_from,recv_info_p2r(i_prc_to,irank))
2781 
2782  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
2783  do v = 1, vmax
2784  do k = 1, kmax
2785  kk = (v-1) * kmax + k
2786  sendbuf_h2p_sp(kk,i_npl) = var(ij_from,k,l_from,v)
2787  enddo
2788  enddo
2789 
2790  req_count = req_count + 1
2791  totalsize = kmax * vmax
2792  rank = recv_info_p2r(i_prc_from,irank)
2793  tag = recv_info_p2r(i_prc_from,irank) + 1000000
2794 
2795  call mpi_isend( sendbuf_h2p_sp(1,i_npl), & ! [IN]
2796  totalsize, & ! [IN]
2797  mpi_real, & ! [IN]
2798  rank, & ! [IN]
2799  tag, & ! [IN]
2800  prc_local_comm_world, & ! [IN]
2801  req_list(req_count), & ! [OUT]
2802  ierr ) ! [OUT]
2803  endif
2804 
2805  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
2806  do v = 1, vmax
2807  do k = 1, kmax
2808  kk = (v-1) * kmax + k
2809  sendbuf_h2p_sp(kk,i_spl) = var(ij_from,k,l_from,v)
2810  enddo
2811  enddo
2812 
2813  req_count = req_count + 1
2814  totalsize = kmax * vmax
2815  rank = recv_info_p2r(i_prc_from,irank)
2816  tag = recv_info_p2r(i_prc_from,irank) + 2000000
2817 
2818  call mpi_isend( sendbuf_h2p_sp(1,i_spl), & ! [IN]
2819  totalsize, & ! [IN]
2820  mpi_real, & ! [IN]
2821  rank, & ! [IN]
2822  tag, & ! [IN]
2823  prc_local_comm_world, & ! [IN]
2824  req_list(req_count), & ! [OUT]
2825  ierr ) ! [OUT]
2826  endif
2827  enddo
2828  enddo
2829 
2830  !--- copy p2r-reverse
2831  do irank = 1, copy_nmax_p2r
2832  do ipos = 1, copy_info_p2r(i_size)
2833  ij_from = copy_list_p2r(i_grid_to ,ipos)
2834  l_from = copy_list_p2r(i_l_to ,ipos)
2835  r_from = prc_rgn_lp2r(l_from,copy_info_p2r(i_prc_to))
2836  ij_to = copy_list_p2r(i_grid_from,ipos)
2837  l_to = copy_list_p2r(i_l_from ,ipos)
2838  r_to = prc_rgn_lp2r(l_to,copy_info_p2r(i_prc_from))
2839 
2840  if ( r_from == prc_rgn_rgn4pl(i_npl) &
2841  .OR. r_from == prc_rgn_rgn4pl(i_spl) ) then
2842  do v = 1, vmax
2843  do k = 1, kmax
2844  var_pl(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2845  enddo
2846  enddo
2847  endif
2848  enddo
2849  enddo
2850 
2851  !--- wait all
2852  if ( req_count > 0 ) then
2853  call mpi_waitall( req_count, & ! [IN]
2854  req_list(:), & ! [IN]
2855  mpi_statuses_ignore, & ! [OUT]
2856  ierr ) ! [OUT]
2857  endif
2858 
2859  !--- unpack p2r-reverse
2860  do irank = 1, send_nmax_p2r
2861  do ipos = 1, send_info_p2r(i_size,irank)
2862  l_from = send_list_p2r(i_l_to ,ipos,irank)
2863  r_from = prc_rgn_lp2r(l_from,send_info_p2r(i_prc_to,irank))
2864 
2865  ij_to = send_list_p2r(i_grid_from,ipos,irank)
2866  l_to = send_list_p2r(i_l_from ,ipos,irank)
2867 
2868  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
2869  do v = 1, vmax
2870  do k = 1, kmax
2871  kk = (v-1) * kmax + k
2872  var_pl(ij_to,k,l_to,v) = recvbuf_h2p_sp(kk,i_npl)
2873  enddo
2874  enddo
2875  endif
2876 
2877  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
2878  do v = 1, vmax
2879  do k = 1, kmax
2880  kk = (v-1) * kmax + k
2881  var_pl(ij_to,k,l_to,v) = recvbuf_h2p_sp(kk,i_spl)
2882  enddo
2883  enddo
2884  endif
2885  enddo
2886  enddo
2887 
2888  endif
2889 
2890  call comm_data_transfer_sp(var,var_pl)
2891 
2892  call prof_rapend('COMM_var',2)
2893 
2894  return
integer, parameter, public i_spl
south pole
module PROCESS
Definition: scale_prc.F90:11
integer, dimension(prc_rgn_total_pl), public prc_rgn_rgn4pl
region, having pole data in the halo
integer, public kmax
of computational cells: z, local
integer, parameter, public i_npl
north pole
integer, parameter, public prc_rgn_total_pl
number of pole region
integer, dimension(:,:), allocatable, public prc_rgn_lp2r
l,prc => rgn
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87
module process / icoA
Here is the call graph for this function:

◆ comm_var_dp()

subroutine scale_comm_icoa::comm_var_dp ( real(dp), dimension (adm_gall, kmax,adm_lall, vmax), intent(inout)  var,
real(dp), dimension(adm_gall_pl,kmax,adm_lall_pl,vmax), intent(inout)  var_pl,
integer, intent(in)  kmax,
integer, intent(in)  vmax 
)

Data transfer with region halo => pole center.

Definition at line 2904 of file scale_comm_icoA.F90.

References comm_data_transfer_dp(), comm_pl, copy_info_p2r, copy_list_p2r, scale_prc_icoa::i_npl, scale_prc_icoa::i_spl, scale_prc::prc_local_comm_world, scale_prc_icoa::prc_rgn_lp2r, scale_prc_icoa::prc_rgn_rgn4pl, scale_prc_icoa::prc_rgn_total_pl, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), recv_info_p2r, recv_list_p2r, req_list, send_info_p2r, and send_list_p2r.

2904  use scale_prc, only: &
2906  use scale_prc_icoa, only: &
2907  prc_rgn_total_pl, &
2908  prc_rgn_rgn4pl, &
2909  prc_rgn_lp2r, &
2910  i_npl, &
2911  i_spl
2912  implicit none
2913 
2914  integer, intent(in) :: kmax
2915  integer, intent(in) :: vmax
2916  real(DP), intent(inout) :: var (adm_gall, kmax,adm_lall, vmax)
2917  real(DP), intent(inout) :: var_pl(adm_gall_pl,kmax,adm_lall_pl,vmax)
2918 
2919  real(DP) :: sendbuf_h2p_dp(kmax*vmax,prc_rgn_total_pl)
2920  real(DP) :: recvbuf_h2p_dp(kmax*vmax,prc_rgn_total_pl)
2921 
2922  integer :: totalsize, rank, tag
2923  integer :: irank, ipos
2924  integer :: ij_from, l_from, ij_to, l_to
2925  integer :: r_from, r_to
2926 
2927  integer :: k, v, kk
2928  integer :: ierr
2929  !---------------------------------------------------------------------------
2930 
2931  if ( comm_apply_barrier ) then
2932  call prof_rapstart('COMM_barrier',2)
2933  call mpi_barrier(prc_local_comm_world,ierr)
2934  call prof_rapend ('COMM_barrier',2)
2935  endif
2936 
2937  call prof_rapstart('COMM_var',2)
2938 
2939  if ( comm_pl ) then
2940 
2941  req_count = 0
2942 
2943  !--- receive p2r-reverse
2944  do irank = 1, send_nmax_p2r
2945  do ipos = 1, send_info_p2r(i_size,irank)
2946  l_from = send_list_p2r(i_l_to ,ipos,irank)
2947  r_from = prc_rgn_lp2r(l_from,send_info_p2r(i_prc_to,irank))
2948 
2949  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
2950  req_count = req_count + 1
2951  totalsize = kmax * vmax
2952  rank = send_info_p2r(i_prc_to ,irank)
2953  tag = send_info_p2r(i_prc_from,irank) + 1000000
2954 
2955  call mpi_irecv( recvbuf_h2p_dp(1,i_npl), & ! [OUT]
2956  totalsize, & ! [IN]
2957  mpi_double_precision, & ! [IN]
2958  rank, & ! [IN]
2959  tag, & ! [IN]
2960  prc_local_comm_world, & ! [IN]
2961  req_list(req_count), & ! [OUT]
2962  ierr ) ! [OUT]
2963  endif
2964 
2965  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
2966  req_count = req_count + 1
2967  totalsize = kmax * vmax
2968  rank = send_info_p2r(i_prc_to ,irank)
2969  tag = send_info_p2r(i_prc_from,irank) + 2000000
2970 
2971  call mpi_irecv( recvbuf_h2p_dp(1,i_spl), & ! [OUT]
2972  totalsize, & ! [IN]
2973  mpi_double_precision, & ! [IN]
2974  rank, & ! [IN]
2975  tag, & ! [IN]
2976  prc_local_comm_world, & ! [IN]
2977  req_list(req_count), & ! [OUT]
2978  ierr ) ! [OUT]
2979  endif
2980  enddo
2981  enddo
2982 
2983  !--- pack and send p2r-reverse
2984  do irank = 1, recv_nmax_p2r
2985  do ipos = 1, recv_info_p2r(i_size,irank)
2986  ij_from = recv_list_p2r(i_grid_to,ipos,irank)
2987  l_from = recv_list_p2r(i_l_to ,ipos,irank)
2988  r_from = prc_rgn_lp2r(l_from,recv_info_p2r(i_prc_to,irank))
2989 
2990  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
2991  do v = 1, vmax
2992  do k = 1, kmax
2993  kk = (v-1) * kmax + k
2994  sendbuf_h2p_dp(kk,i_npl) = var(ij_from,k,l_from,v)
2995  enddo
2996  enddo
2997 
2998  req_count = req_count + 1
2999  totalsize = kmax * vmax
3000  rank = recv_info_p2r(i_prc_from,irank)
3001  tag = recv_info_p2r(i_prc_from,irank) + 1000000
3002 
3003  call mpi_isend( sendbuf_h2p_dp(1,i_npl), & ! [IN]
3004  totalsize, & ! [IN]
3005  mpi_double_precision, & ! [IN]
3006  rank, & ! [IN]
3007  tag, & ! [IN]
3008  prc_local_comm_world, & ! [IN]
3009  req_list(req_count), & ! [OUT]
3010  ierr ) ! [OUT]
3011  endif
3012 
3013  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
3014  do v = 1, vmax
3015  do k = 1, kmax
3016  kk = (v-1) * kmax + k
3017  sendbuf_h2p_dp(kk,i_spl) = var(ij_from,k,l_from,v)
3018  enddo
3019  enddo
3020 
3021  req_count = req_count + 1
3022  totalsize = kmax * vmax
3023  rank = recv_info_p2r(i_prc_from,irank)
3024  tag = recv_info_p2r(i_prc_from,irank) + 2000000
3025 
3026  call mpi_isend( sendbuf_h2p_dp(1,i_spl), & ! [IN]
3027  totalsize, & ! [IN]
3028  mpi_double_precision, & ! [IN]
3029  rank, & ! [IN]
3030  tag, & ! [IN]
3031  prc_local_comm_world, & ! [IN]
3032  req_list(req_count), & ! [OUT]
3033  ierr ) ! [OUT]
3034  endif
3035  enddo
3036  enddo
3037 
3038  !--- copy p2r-reverse
3039  do irank = 1, copy_nmax_p2r
3040  do ipos = 1, copy_info_p2r(i_size)
3041  ij_from = copy_list_p2r(i_grid_to ,ipos)
3042  l_from = copy_list_p2r(i_l_to ,ipos)
3043  r_from = prc_rgn_lp2r(l_from,copy_info_p2r(i_prc_to))
3044  ij_to = copy_list_p2r(i_grid_from,ipos)
3045  l_to = copy_list_p2r(i_l_from ,ipos)
3046  r_to = prc_rgn_lp2r(l_to,copy_info_p2r(i_prc_from))
3047 
3048  if ( r_from == prc_rgn_rgn4pl(i_npl) &
3049  .OR. r_from == prc_rgn_rgn4pl(i_spl) ) then
3050  do v = 1, vmax
3051  do k = 1, kmax
3052  var_pl(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
3053  enddo
3054  enddo
3055  endif
3056  enddo
3057  enddo
3058 
3059  !--- wait all
3060  if ( req_count > 0 ) then
3061  call mpi_waitall( req_count, & ! [IN]
3062  req_list(:), & ! [IN]
3063  mpi_statuses_ignore, & ! [OUT]
3064  ierr ) ! [OUT]
3065  endif
3066 
3067  !--- unpack p2r-reverse
3068  do irank = 1, send_nmax_p2r
3069  do ipos = 1, send_info_p2r(i_size,irank)
3070  l_from = send_list_p2r(i_l_to ,ipos,irank)
3071  r_from = prc_rgn_lp2r(l_from,send_info_p2r(i_prc_to,irank))
3072 
3073  ij_to = send_list_p2r(i_grid_from,ipos,irank)
3074  l_to = send_list_p2r(i_l_from ,ipos,irank)
3075 
3076  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
3077  do v = 1, vmax
3078  do k = 1, kmax
3079  kk = (v-1) * kmax + k
3080  var_pl(ij_to,k,l_to,v) = recvbuf_h2p_dp(kk,i_npl)
3081  enddo
3082  enddo
3083  endif
3084 
3085  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
3086  do v = 1, vmax
3087  do k = 1, kmax
3088  kk = (v-1) * kmax + k
3089  var_pl(ij_to,k,l_to,v) = recvbuf_h2p_dp(kk,i_spl)
3090  enddo
3091  enddo
3092  endif
3093  enddo
3094  enddo
3095 
3096  endif
3097 
3098  call comm_data_transfer_dp(var,var_pl)
3099 
3100  call prof_rapend('COMM_var',2)
3101 
3102  return
integer, parameter, public i_spl
south pole
module PROCESS
Definition: scale_prc.F90:11
integer, dimension(prc_rgn_total_pl), public prc_rgn_rgn4pl
region, having pole data in the halo
integer, public kmax
of computational cells: z, local
integer, parameter, public i_npl
north pole
integer, parameter, public prc_rgn_total_pl
number of pole region
integer, dimension(:,:), allocatable, public prc_rgn_lp2r
l,prc => rgn
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87
module process / icoA
Here is the call graph for this function:

◆ suf()

integer function scale_comm_icoa::suf ( integer  i,
integer  j 
)

suffix calculation

Returns
suf

Definition at line 3109 of file scale_comm_icoA.F90.

References scale_atmos_grid_icoa_index::adm_gall_1d, scale_atmos_grid_icoa_index::adm_gall_pl, scale_atmos_grid_icoa_index::adm_gmax, scale_atmos_grid_icoa_index::adm_gmin, scale_atmos_grid_icoa_index::adm_kall, scale_atmos_grid_icoa_index::adm_kmax, scale_atmos_grid_icoa_index::adm_kmin, scale_atmos_grid_icoa_index::adm_lall, scale_atmos_grid_icoa_index::adm_lall_pl, scale_atmos_grid_icoa_index::adm_vlink, scale_prc_icoa::i_n, scale_prc_icoa::i_s, scale_io::io_fid_log, scale_io::io_l, scale_prc::prc_mpifinish(), scale_prc::prc_myrank, scale_prc_icoa::prc_rgn_have_sgp, scale_prc_icoa::prc_rgn_l2r, and scale_prc_icoa::prc_rgn_vert_num.

Referenced by comm_setup().

3109  implicit none
3110 
3111  integer :: suffix
3112  integer :: i, j
3113  !---------------------------------------------------------------------------
3114 
3115  suffix = adm_gall_1d * (j-1) + i
3116 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ comm_stat_sum_sp()

subroutine scale_comm_icoa::comm_stat_sum_sp ( real(sp), intent(in)  localsum,
real(sp), intent(out)  globalsum 
)

Definition at line 3475 of file scale_comm_icoA.F90.

References comm_pl, scale_prc::prc_local_comm_world, and scale_prc::prc_nprocs.

3475  use scale_prc, only: &
3477  prc_nprocs
3478  implicit none
3479 
3480  real(SP), intent(in) :: localsum
3481  real(SP), intent(out) :: globalsum
3482 
3483  real(SP) :: sendbuf(1)
3484  real(SP) :: recvbuf(prc_nprocs)
3485 
3486  integer :: ierr
3487  !---------------------------------------------------------------------------
3488 
3489  if ( comm_pl ) then
3490  sendbuf(1) = localsum
3491 
3492  call mpi_allgather( sendbuf, &
3493  1, &
3494  mpi_real, &
3495  recvbuf, &
3496  1, &
3497  mpi_real, &
3499  ierr )
3500 
3501  globalsum = sum( recvbuf(:) )
3502  else
3503  globalsum = localsum
3504  endif
3505 
3506  return
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:88
module PROCESS
Definition: scale_prc.F90:11
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87

◆ comm_stat_sum_dp()

subroutine scale_comm_icoa::comm_stat_sum_dp ( real(dp), intent(in)  localsum,
real(dp), intent(out)  globalsum 
)

Definition at line 3511 of file scale_comm_icoA.F90.

References comm_pl, scale_prc::prc_local_comm_world, and scale_prc::prc_nprocs.

3511  use scale_prc, only: &
3513  prc_nprocs
3514  implicit none
3515 
3516  real(DP), intent(in) :: localsum
3517  real(DP), intent(out) :: globalsum
3518 
3519  real(DP) :: sendbuf(1)
3520  real(DP) :: recvbuf(prc_nprocs)
3521 
3522  integer :: ierr
3523  !---------------------------------------------------------------------------
3524 
3525  if ( comm_pl ) then
3526  sendbuf(1) = localsum
3527 
3528  call mpi_allgather( sendbuf, &
3529  1, &
3530  mpi_double_precision, &
3531  recvbuf, &
3532  1, &
3533  mpi_double_precision, &
3535  ierr )
3536 
3537  globalsum = sum( recvbuf(:) )
3538  else
3539  globalsum = localsum
3540  endif
3541 
3542  return
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:88
module PROCESS
Definition: scale_prc.F90:11
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87

◆ comm_stat_sum_eachlayer_sp()

subroutine scale_comm_icoa::comm_stat_sum_eachlayer_sp ( integer, intent(in)  kall,
real(sp), dimension (kall), intent(in)  localsum,
real(sp), dimension(kall), intent(out)  globalsum 
)

Definition at line 3547 of file scale_comm_icoA.F90.

References comm_pl, scale_prc::prc_local_comm_world, and scale_prc::prc_nprocs.

3547  use scale_prc, only: &
3549  prc_nprocs
3550  implicit none
3551 
3552  integer, intent(in) :: kall
3553  real(SP), intent(in) :: localsum (kall)
3554  real(SP), intent(out) :: globalsum(kall)
3555 
3556  real(SP) :: sendbuf(kall)
3557  integer :: displs (prc_nprocs)
3558  integer :: counts (prc_nprocs)
3559  real(SP) :: recvbuf(kall,prc_nprocs)
3560 
3561  integer :: ierr
3562  integer :: k, p
3563  !---------------------------------------------------------------------------
3564 
3565  do p = 1, prc_nprocs
3566  displs(p) = (p-1) * kall
3567  counts(p) = kall
3568  enddo
3569 
3570  if ( comm_pl ) then
3571  sendbuf(:) = localsum(:)
3572 
3573  call mpi_allgatherv( sendbuf, &
3574  kall, &
3575  mpi_real, &
3576  recvbuf, &
3577  counts, &
3578  displs, &
3579  mpi_real, &
3581  ierr )
3582 
3583  do k = 1, kall
3584  globalsum(k) = sum( recvbuf(k,:) )
3585  enddo
3586  else
3587  do k = 1, kall
3588  globalsum(k) = localsum(k)
3589  enddo
3590  endif
3591 
3592  return
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:88
module PROCESS
Definition: scale_prc.F90:11
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87

◆ comm_stat_sum_eachlayer_dp()

subroutine scale_comm_icoa::comm_stat_sum_eachlayer_dp ( integer, intent(in)  kall,
real(dp), dimension (kall), intent(in)  localsum,
real(dp), dimension(kall), intent(out)  globalsum 
)

Definition at line 3597 of file scale_comm_icoA.F90.

References comm_pl, scale_prc::prc_local_comm_world, and scale_prc::prc_nprocs.

3597  use scale_prc, only: &
3599  prc_nprocs
3600  implicit none
3601 
3602  integer, intent(in) :: kall
3603  real(DP), intent(in) :: localsum (kall)
3604  real(DP), intent(out) :: globalsum(kall)
3605 
3606  real(DP) :: sendbuf(kall)
3607  integer :: displs (prc_nprocs)
3608  integer :: counts (prc_nprocs)
3609  real(DP) :: recvbuf(kall,prc_nprocs)
3610 
3611  integer :: ierr
3612  integer :: k, p
3613  !---------------------------------------------------------------------------
3614 
3615  do p = 1, prc_nprocs
3616  displs(p) = (p-1) * kall
3617  counts(p) = kall
3618  enddo
3619 
3620  if ( comm_pl ) then
3621  sendbuf(:) = localsum(:)
3622 
3623  call mpi_allgatherv( sendbuf, &
3624  kall, &
3625  mpi_double_precision, &
3626  recvbuf, &
3627  counts, &
3628  displs, &
3629  mpi_double_precision, &
3631  ierr )
3632 
3633  do k = 1, kall
3634  globalsum(k) = sum( recvbuf(k,:) )
3635  enddo
3636  else
3637  do k = 1, kall
3638  globalsum(k) = localsum(k)
3639  enddo
3640  endif
3641 
3642  return
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:88
module PROCESS
Definition: scale_prc.F90:11
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87

◆ comm_stat_avg_sp()

subroutine scale_comm_icoa::comm_stat_avg_sp ( real(sp), intent(in)  localavg,
real(sp), intent(out)  globalavg 
)

Definition at line 3647 of file scale_comm_icoA.F90.

References comm_pl, scale_prc::prc_local_comm_world, and scale_prc::prc_nprocs.

3647  use scale_prc, only: &
3649  prc_nprocs
3650  implicit none
3651 
3652  real(SP), intent(in) :: localavg
3653  real(SP), intent(out) :: globalavg
3654 
3655  real(SP) :: sendbuf(1)
3656  real(SP) :: recvbuf(prc_nprocs)
3657 
3658  integer :: ierr
3659  !---------------------------------------------------------------------------
3660 
3661  if ( comm_pl ) then
3662  sendbuf(1) = localavg
3663 
3664  call mpi_allgather( sendbuf, &
3665  1, &
3666  mpi_real, &
3667  recvbuf, &
3668  1, &
3669  mpi_real, &
3671  ierr )
3672 
3673  globalavg = sum( recvbuf(:) ) / real(prc_nprocs,kind=sp)
3674  else
3675  globalavg = localavg
3676  endif
3677 
3678  return
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:88
module PROCESS
Definition: scale_prc.F90:11
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87
integer, parameter, public sp

◆ comm_stat_avg_dp()

subroutine scale_comm_icoa::comm_stat_avg_dp ( real(dp), intent(in)  localavg,
real(dp), intent(out)  globalavg 
)

Definition at line 3683 of file scale_comm_icoA.F90.

References comm_pl, scale_prc::prc_local_comm_world, and scale_prc::prc_nprocs.

3683  use scale_prc, only: &
3685  prc_nprocs
3686  implicit none
3687 
3688  real(DP), intent(in) :: localavg
3689  real(DP), intent(out) :: globalavg
3690 
3691  real(DP) :: sendbuf(1)
3692  real(DP) :: recvbuf(prc_nprocs)
3693 
3694  integer :: ierr
3695  !---------------------------------------------------------------------------
3696 
3697  if ( comm_pl ) then
3698  sendbuf(1) = localavg
3699 
3700  call mpi_allgather( sendbuf, &
3701  1, &
3702  mpi_double_precision, &
3703  recvbuf, &
3704  1, &
3705  mpi_double_precision, &
3707  ierr )
3708 
3709  globalavg = sum( recvbuf(:) ) / real(prc_nprocs,kind=dp)
3710  else
3711  globalavg = localavg
3712  endif
3713 
3714  return
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:88
module PROCESS
Definition: scale_prc.F90:11
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87
integer, parameter, public dp

◆ comm_stat_max_sp()

subroutine scale_comm_icoa::comm_stat_max_sp ( real(sp), intent(in)  localmax,
real(sp), intent(out)  globalmax 
)

Definition at line 3719 of file scale_comm_icoA.F90.

References scale_prc::prc_local_comm_world, and scale_prc::prc_nprocs.

3719  use scale_prc, only: &
3721  prc_nprocs
3722  implicit none
3723 
3724  real(SP), intent(in) :: localmax
3725  real(SP), intent(out) :: globalmax
3726 
3727  real(SP) :: sendbuf(1)
3728  real(SP) :: recvbuf(prc_nprocs)
3729 
3730  integer :: ierr
3731  !---------------------------------------------------------------------------
3732 
3733  sendbuf(1) = localmax
3734 
3735  call mpi_allgather( sendbuf, &
3736  1, &
3737  mpi_real, &
3738  recvbuf, &
3739  1, &
3740  mpi_real, &
3742  ierr )
3743 
3744  globalmax = maxval( recvbuf(:) )
3745 
3746  return
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:88
module PROCESS
Definition: scale_prc.F90:11
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87

◆ comm_stat_max_dp()

subroutine scale_comm_icoa::comm_stat_max_dp ( real(dp), intent(in)  localmax,
real(dp), intent(out)  globalmax 
)

Definition at line 3751 of file scale_comm_icoA.F90.

References scale_prc::prc_local_comm_world, and scale_prc::prc_nprocs.

3751  use scale_prc, only: &
3753  prc_nprocs
3754  implicit none
3755 
3756  real(DP), intent(in) :: localmax
3757  real(DP), intent(out) :: globalmax
3758 
3759  real(DP) :: sendbuf(1)
3760  real(DP) :: recvbuf(prc_nprocs)
3761 
3762  integer :: ierr
3763  !---------------------------------------------------------------------------
3764 
3765  sendbuf(1) = localmax
3766 
3767  call mpi_allgather( sendbuf, &
3768  1, &
3769  mpi_double_precision, &
3770  recvbuf, &
3771  1, &
3772  mpi_double_precision, &
3774  ierr )
3775 
3776  globalmax = maxval( recvbuf(:) )
3777 
3778  return
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:88
module PROCESS
Definition: scale_prc.F90:11
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87

◆ comm_stat_min_sp()

subroutine scale_comm_icoa::comm_stat_min_sp ( real(sp), intent(in)  localmin,
real(sp), intent(out)  globalmin 
)

Definition at line 3783 of file scale_comm_icoA.F90.

References scale_prc::prc_local_comm_world, and scale_prc::prc_nprocs.

3783  use scale_prc, only: &
3785  prc_nprocs
3786  implicit none
3787 
3788  real(SP), intent(in) :: localmin
3789  real(SP), intent(out) :: globalmin
3790 
3791  real(SP) :: sendbuf(1)
3792  real(SP) :: recvbuf(prc_nprocs)
3793 
3794  integer :: ierr
3795  !---------------------------------------------------------------------------
3796 
3797  sendbuf(1) = localmin
3798 
3799  call mpi_allgather( sendbuf, &
3800  1, &
3801  mpi_real, &
3802  recvbuf, &
3803  1, &
3804  mpi_real, &
3806  ierr )
3807 
3808  globalmin = minval( recvbuf(:) )
3809 
3810  return
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:88
module PROCESS
Definition: scale_prc.F90:11
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87

◆ comm_stat_min_dp()

subroutine scale_comm_icoa::comm_stat_min_dp ( real(dp), intent(in)  localmin,
real(dp), intent(out)  globalmin 
)

Definition at line 3815 of file scale_comm_icoA.F90.

References scale_prc::prc_local_comm_world, and scale_prc::prc_nprocs.

3815  use scale_prc, only: &
3817  prc_nprocs
3818  implicit none
3819 
3820  real(DP), intent(in) :: localmin
3821  real(DP), intent(out) :: globalmin
3822 
3823  real(DP) :: sendbuf(1)
3824  real(DP) :: recvbuf(prc_nprocs)
3825 
3826  integer :: ierr
3827  !---------------------------------------------------------------------------
3828 
3829  sendbuf(1) = localmin
3830 
3831  call mpi_allgather( sendbuf, &
3832  1, &
3833  mpi_double_precision, &
3834  recvbuf, &
3835  1, &
3836  mpi_double_precision, &
3838  ierr )
3839 
3840  globalmin = minval( recvbuf(:) )
3841 
3842  return
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:88
module PROCESS
Definition: scale_prc.F90:11
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87

Variable Documentation

◆ comm_pl

logical, public scale_comm_icoa::comm_pl = .true.

◆ comm_datatype

integer, public scale_comm_icoa::comm_datatype

Definition at line 79 of file scale_comm_icoA.F90.

Referenced by comm_setup().

79  integer, public :: comm_datatype
integer, public comm_datatype
datatype of variable

◆ copy_info_r2r

integer, dimension(:), allocatable, public scale_comm_icoa::copy_info_r2r

Definition at line 82 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_nopl(), comm_data_transfer_sp(), and comm_setup().

82  integer, public, allocatable :: copy_info_r2r(:)

◆ node

integer allocatable public scale_comm_icoa::node

Definition at line 82 of file scale_comm_icoA.F90.

◆ information

integer, dimension (singular point), allocatable, public scale_comm_icoa::information

Definition at line 82 of file scale_comm_icoA.F90.

◆ recv_info_r2r

integer, dimension(:,:), allocatable, public scale_comm_icoa::recv_info_r2r

Definition at line 83 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_nopl(), comm_data_transfer_sp(), and comm_setup().

83  integer, public, allocatable :: recv_info_r2r(:,:)

◆ send_info_r2r

integer, dimension(:,:), allocatable, public scale_comm_icoa::send_info_r2r

Definition at line 84 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_nopl(), comm_data_transfer_sp(), and comm_setup().

84  integer, public, allocatable :: send_info_r2r(:,:)

◆ copy_info_p2r

integer, dimension(:), allocatable, public scale_comm_icoa::copy_info_p2r

Definition at line 86 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_sp(), comm_setup(), comm_var_dp(), and comm_var_sp().

86  integer, public, allocatable :: copy_info_p2r(:)

◆ recv_info_p2r

integer, dimension(:,:), allocatable, public scale_comm_icoa::recv_info_p2r

Definition at line 87 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_sp(), comm_setup(), comm_var_dp(), and comm_var_sp().

87  integer, public, allocatable :: recv_info_p2r(:,:)

◆ send_info_p2r

integer, dimension(:,:), allocatable, public scale_comm_icoa::send_info_p2r

Definition at line 88 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_sp(), comm_setup(), comm_var_dp(), and comm_var_sp().

88  integer, public, allocatable :: send_info_p2r(:,:)

◆ copy_info_r2p

integer, dimension(:), allocatable, public scale_comm_icoa::copy_info_r2p

Definition at line 90 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_sp(), and comm_setup().

90  integer, public, allocatable :: copy_info_r2p(:)

◆ recv_info_r2p

integer, dimension(:,:), allocatable, public scale_comm_icoa::recv_info_r2p

Definition at line 91 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_sp(), and comm_setup().

91  integer, public, allocatable :: recv_info_r2p(:,:)

◆ send_info_r2p

integer, dimension(:,:), allocatable, public scale_comm_icoa::send_info_r2p

Definition at line 92 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_sp(), and comm_setup().

92  integer, public, allocatable :: send_info_r2p(:,:)

◆ singular_info

integer, dimension(:), allocatable, public scale_comm_icoa::singular_info

Definition at line 94 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_nopl(), comm_data_transfer_sp(), and comm_setup().

94  integer, public, allocatable :: singular_info(:)

◆ copy_list_r2r

integer, dimension(:,:), allocatable, public scale_comm_icoa::copy_list_r2r

Definition at line 97 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_nopl(), comm_data_transfer_sp(), and comm_setup().

97  integer, public, allocatable :: copy_list_r2r(:,:)

◆ grid

integer allocatable public scale_comm_icoa::grid

Definition at line 97 of file scale_comm_icoA.F90.

◆ local

integer allocatable public scale_comm_icoa::local

Definition at line 97 of file scale_comm_icoA.F90.

◆ region

integer allocatable public scale_comm_icoa::region

Definition at line 97 of file scale_comm_icoA.F90.

◆ recv_list_r2r

integer, dimension(:,:,:), allocatable, public scale_comm_icoa::recv_list_r2r

Definition at line 98 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_nopl(), comm_data_transfer_sp(), and comm_setup().

98  integer, public, allocatable :: recv_list_r2r(:,:,:)

◆ send_list_r2r

integer, dimension(:,:,:), allocatable, public scale_comm_icoa::send_list_r2r

Definition at line 99 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_nopl(), comm_data_transfer_sp(), and comm_setup().

99  integer, public, allocatable :: send_list_r2r(:,:,:)

◆ copy_list_p2r

integer, dimension(:,:), allocatable, public scale_comm_icoa::copy_list_p2r

Definition at line 101 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_sp(), comm_setup(), comm_var_dp(), and comm_var_sp().

101  integer, public, allocatable :: copy_list_p2r(:,:)

◆ recv_list_p2r

integer, dimension(:,:,:), allocatable, public scale_comm_icoa::recv_list_p2r

Definition at line 102 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_sp(), comm_setup(), comm_var_dp(), and comm_var_sp().

102  integer, public, allocatable :: recv_list_p2r(:,:,:)

◆ send_list_p2r

integer, dimension(:,:,:), allocatable, public scale_comm_icoa::send_list_p2r

Definition at line 103 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_sp(), comm_setup(), comm_var_dp(), and comm_var_sp().

103  integer, public, allocatable :: send_list_p2r(:,:,:)

◆ copy_list_r2p

integer, dimension(:,:), allocatable, public scale_comm_icoa::copy_list_r2p

Definition at line 105 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_sp(), and comm_setup().

105  integer, public, allocatable :: copy_list_r2p(:,:)

◆ recv_list_r2p

integer, dimension(:,:,:), allocatable, public scale_comm_icoa::recv_list_r2p

Definition at line 106 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_sp(), and comm_setup().

106  integer, public, allocatable :: recv_list_r2p(:,:,:)

◆ send_list_r2p

integer, dimension(:,:,:), allocatable, public scale_comm_icoa::send_list_r2p

Definition at line 107 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_sp(), and comm_setup().

107  integer, public, allocatable :: send_list_r2p(:,:,:)

◆ singular_list

integer, dimension(:,:), allocatable, public scale_comm_icoa::singular_list

Definition at line 109 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_nopl(), comm_data_transfer_sp(), and comm_setup().

109  integer, public, allocatable :: singular_list(:,:)

◆ req_list

integer, dimension(:), allocatable, public scale_comm_icoa::req_list

Definition at line 112 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_nopl(), comm_data_transfer_sp(), comm_setup(), comm_var_dp(), and comm_var_sp().

112  integer, public, allocatable :: req_list(:)

◆ sendbuf_r2r_sp

real(sp), dimension(:,:), allocatable, public scale_comm_icoa::sendbuf_r2r_sp

Definition at line 114 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_sp(), and comm_setup().

114  real(SP), public, allocatable :: sendbuf_r2r_sp(:,:)

◆ recvbuf_r2r_sp

real(sp), dimension(:,:), allocatable, public scale_comm_icoa::recvbuf_r2r_sp

Definition at line 115 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_sp(), and comm_setup().

115  real(SP), public, allocatable :: recvbuf_r2r_sp(:,:)

◆ sendbuf_p2r_sp

real(sp), dimension(:,:), allocatable, public scale_comm_icoa::sendbuf_p2r_sp

Definition at line 116 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_sp(), and comm_setup().

116  real(SP), public, allocatable :: sendbuf_p2r_sp(:,:)

◆ recvbuf_p2r_sp

real(sp), dimension(:,:), allocatable, public scale_comm_icoa::recvbuf_p2r_sp

Definition at line 117 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_sp(), and comm_setup().

117  real(SP), public, allocatable :: recvbuf_p2r_sp(:,:)

◆ sendbuf_r2p_sp

real(sp), dimension(:,:), allocatable, public scale_comm_icoa::sendbuf_r2p_sp

Definition at line 118 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_sp(), and comm_setup().

118  real(SP), public, allocatable :: sendbuf_r2p_sp(:,:)

◆ recvbuf_r2p_sp

real(sp), dimension(:,:), allocatable, public scale_comm_icoa::recvbuf_r2p_sp

Definition at line 119 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_sp(), and comm_setup().

119  real(SP), public, allocatable :: recvbuf_r2p_sp(:,:)

◆ sendbuf_r2r_dp

real(dp), dimension(:,:), allocatable, public scale_comm_icoa::sendbuf_r2r_dp

Definition at line 121 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_nopl(), and comm_setup().

121  real(DP), public, allocatable :: sendbuf_r2r_dp(:,:)

◆ recvbuf_r2r_dp

real(dp), dimension(:,:), allocatable, public scale_comm_icoa::recvbuf_r2r_dp

Definition at line 122 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), comm_data_transfer_nopl(), and comm_setup().

122  real(DP), public, allocatable :: recvbuf_r2r_dp(:,:)

◆ sendbuf_p2r_dp

real(dp), dimension(:,:), allocatable, public scale_comm_icoa::sendbuf_p2r_dp

Definition at line 123 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), and comm_setup().

123  real(DP), public, allocatable :: sendbuf_p2r_dp(:,:)

◆ recvbuf_p2r_dp

real(dp), dimension(:,:), allocatable, public scale_comm_icoa::recvbuf_p2r_dp

Definition at line 124 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), and comm_setup().

124  real(DP), public, allocatable :: recvbuf_p2r_dp(:,:)

◆ sendbuf_r2p_dp

real(dp), dimension(:,:), allocatable, public scale_comm_icoa::sendbuf_r2p_dp

Definition at line 125 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), and comm_setup().

125  real(DP), public, allocatable :: sendbuf_r2p_dp(:,:)

◆ recvbuf_r2p_dp

real(dp), dimension(:,:), allocatable, public scale_comm_icoa::recvbuf_r2p_dp

Definition at line 126 of file scale_comm_icoA.F90.

Referenced by comm_data_transfer_dp(), and comm_setup().

126  real(DP), public, allocatable :: recvbuf_r2p_dp(:,:)