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.

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

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_pl, 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().

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 1564 of file scale_comm_icoA.F90.

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

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().

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 2043 of file scale_comm_icoA.F90.

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

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().

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 2521 of file scale_comm_icoA.F90.

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

References scale_atmos_grid_icoa_index::adm_kall, copy_info_r2r, copy_list_r2r, scale_atmos_grid_icoa_index::imax, 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.

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 2702 of file scale_comm_icoA.F90.

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

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.

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 2910 of file scale_comm_icoA.F90.

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

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.

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 3115 of file scale_comm_icoA.F90.

3115  implicit none
3116 
3117  integer :: suffix
3118  integer :: i, j
3119  !---------------------------------------------------------------------------
3120 
3121  suffix = adm_gall_1d * (j-1) + i
3122 

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_prc_icoa::i_n, scale_prc_icoa::i_npl, scale_prc_icoa::i_s, scale_prc_icoa::i_spl, 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, scale_prc_icoa::prc_rgn_vert_pl, and scale_precision::rp.

Referenced by comm_setup().

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 3483 of file scale_comm_icoA.F90.

3483  use scale_prc, only: &
3485  prc_nprocs
3486  implicit none
3487 
3488  real(SP), intent(in) :: localsum
3489  real(SP), intent(out) :: globalsum
3490 
3491  real(SP) :: sendbuf(1)
3492  real(SP) :: recvbuf(PRC_nprocs)
3493 
3494  integer :: ierr
3495  !---------------------------------------------------------------------------
3496 
3497  if ( comm_pl ) then
3498  sendbuf(1) = localsum
3499 
3500  call mpi_allgather( sendbuf, &
3501  1, &
3502  mpi_real, &
3503  recvbuf, &
3504  1, &
3505  mpi_real, &
3507  ierr )
3508 
3509  globalsum = sum( recvbuf(:) )
3510  else
3511  globalsum = localsum
3512  endif
3513 
3514  return

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

◆ 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 3519 of file scale_comm_icoA.F90.

3519  use scale_prc, only: &
3521  prc_nprocs
3522  implicit none
3523 
3524  real(DP), intent(in) :: localsum
3525  real(DP), intent(out) :: globalsum
3526 
3527  real(DP) :: sendbuf(1)
3528  real(DP) :: recvbuf(PRC_nprocs)
3529 
3530  integer :: ierr
3531  !---------------------------------------------------------------------------
3532 
3533  if ( comm_pl ) then
3534  sendbuf(1) = localsum
3535 
3536  call mpi_allgather( sendbuf, &
3537  1, &
3538  mpi_double_precision, &
3539  recvbuf, &
3540  1, &
3541  mpi_double_precision, &
3543  ierr )
3544 
3545  globalsum = sum( recvbuf(:) )
3546  else
3547  globalsum = localsum
3548  endif
3549 
3550  return

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

◆ 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 3555 of file scale_comm_icoA.F90.

3555  use scale_prc, only: &
3557  prc_nprocs
3558  implicit none
3559 
3560  integer, intent(in) :: kall
3561  real(SP), intent(in) :: localsum (kall)
3562  real(SP), intent(out) :: globalsum(kall)
3563 
3564  real(SP) :: sendbuf(kall)
3565  integer :: displs (PRC_nprocs)
3566  integer :: counts (PRC_nprocs)
3567  real(SP) :: recvbuf(kall,PRC_nprocs)
3568 
3569  integer :: ierr
3570  integer :: k, p
3571  !---------------------------------------------------------------------------
3572 
3573  do p = 1, prc_nprocs
3574  displs(p) = (p-1) * kall
3575  counts(p) = kall
3576  enddo
3577 
3578  if ( comm_pl ) then
3579  sendbuf(:) = localsum(:)
3580 
3581  call mpi_allgatherv( sendbuf, &
3582  kall, &
3583  mpi_real, &
3584  recvbuf, &
3585  counts, &
3586  displs, &
3587  mpi_real, &
3589  ierr )
3590 
3591  do k = 1, kall
3592  globalsum(k) = sum( recvbuf(k,:) )
3593  enddo
3594  else
3595  do k = 1, kall
3596  globalsum(k) = localsum(k)
3597  enddo
3598  endif
3599 
3600  return

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

◆ 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 3605 of file scale_comm_icoA.F90.

3605  use scale_prc, only: &
3607  prc_nprocs
3608  implicit none
3609 
3610  integer, intent(in) :: kall
3611  real(DP), intent(in) :: localsum (kall)
3612  real(DP), intent(out) :: globalsum(kall)
3613 
3614  real(DP) :: sendbuf(kall)
3615  integer :: displs (PRC_nprocs)
3616  integer :: counts (PRC_nprocs)
3617  real(DP) :: recvbuf(kall,PRC_nprocs)
3618 
3619  integer :: ierr
3620  integer :: k, p
3621  !---------------------------------------------------------------------------
3622 
3623  do p = 1, prc_nprocs
3624  displs(p) = (p-1) * kall
3625  counts(p) = kall
3626  enddo
3627 
3628  if ( comm_pl ) then
3629  sendbuf(:) = localsum(:)
3630 
3631  call mpi_allgatherv( sendbuf, &
3632  kall, &
3633  mpi_double_precision, &
3634  recvbuf, &
3635  counts, &
3636  displs, &
3637  mpi_double_precision, &
3639  ierr )
3640 
3641  do k = 1, kall
3642  globalsum(k) = sum( recvbuf(k,:) )
3643  enddo
3644  else
3645  do k = 1, kall
3646  globalsum(k) = localsum(k)
3647  enddo
3648  endif
3649 
3650  return

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

◆ 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 3655 of file scale_comm_icoA.F90.

3655  use scale_prc, only: &
3657  prc_nprocs
3658  implicit none
3659 
3660  real(SP), intent(in) :: localavg
3661  real(SP), intent(out) :: globalavg
3662 
3663  real(SP) :: sendbuf(1)
3664  real(SP) :: recvbuf(PRC_nprocs)
3665 
3666  integer :: ierr
3667  !---------------------------------------------------------------------------
3668 
3669  if ( comm_pl ) then
3670  sendbuf(1) = localavg
3671 
3672  call mpi_allgather( sendbuf, &
3673  1, &
3674  mpi_real, &
3675  recvbuf, &
3676  1, &
3677  mpi_real, &
3679  ierr )
3680 
3681  globalavg = sum( recvbuf(:) ) / real(prc_nprocs,kind=sp)
3682  else
3683  globalavg = localavg
3684  endif
3685 
3686  return

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

◆ 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 3691 of file scale_comm_icoA.F90.

3691  use scale_prc, only: &
3693  prc_nprocs
3694  implicit none
3695 
3696  real(DP), intent(in) :: localavg
3697  real(DP), intent(out) :: globalavg
3698 
3699  real(DP) :: sendbuf(1)
3700  real(DP) :: recvbuf(PRC_nprocs)
3701 
3702  integer :: ierr
3703  !---------------------------------------------------------------------------
3704 
3705  if ( comm_pl ) then
3706  sendbuf(1) = localavg
3707 
3708  call mpi_allgather( sendbuf, &
3709  1, &
3710  mpi_double_precision, &
3711  recvbuf, &
3712  1, &
3713  mpi_double_precision, &
3715  ierr )
3716 
3717  globalavg = sum( recvbuf(:) ) / real(prc_nprocs,kind=dp)
3718  else
3719  globalavg = localavg
3720  endif
3721 
3722  return

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

◆ 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 3727 of file scale_comm_icoA.F90.

3727  use scale_prc, only: &
3729  prc_nprocs
3730  implicit none
3731 
3732  real(SP), intent(in) :: localmax
3733  real(SP), intent(out) :: globalmax
3734 
3735  real(SP) :: sendbuf(1)
3736  real(SP) :: recvbuf(PRC_nprocs)
3737 
3738  integer :: ierr
3739  !---------------------------------------------------------------------------
3740 
3741  sendbuf(1) = localmax
3742 
3743  call mpi_allgather( sendbuf, &
3744  1, &
3745  mpi_real, &
3746  recvbuf, &
3747  1, &
3748  mpi_real, &
3750  ierr )
3751 
3752  globalmax = maxval( recvbuf(:) )
3753 
3754  return

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

◆ 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 3759 of file scale_comm_icoA.F90.

3759  use scale_prc, only: &
3761  prc_nprocs
3762  implicit none
3763 
3764  real(DP), intent(in) :: localmax
3765  real(DP), intent(out) :: globalmax
3766 
3767  real(DP) :: sendbuf(1)
3768  real(DP) :: recvbuf(PRC_nprocs)
3769 
3770  integer :: ierr
3771  !---------------------------------------------------------------------------
3772 
3773  sendbuf(1) = localmax
3774 
3775  call mpi_allgather( sendbuf, &
3776  1, &
3777  mpi_double_precision, &
3778  recvbuf, &
3779  1, &
3780  mpi_double_precision, &
3782  ierr )
3783 
3784  globalmax = maxval( recvbuf(:) )
3785 
3786  return

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

◆ 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 3791 of file scale_comm_icoA.F90.

3791  use scale_prc, only: &
3793  prc_nprocs
3794  implicit none
3795 
3796  real(SP), intent(in) :: localmin
3797  real(SP), intent(out) :: globalmin
3798 
3799  real(SP) :: sendbuf(1)
3800  real(SP) :: recvbuf(PRC_nprocs)
3801 
3802  integer :: ierr
3803  !---------------------------------------------------------------------------
3804 
3805  sendbuf(1) = localmin
3806 
3807  call mpi_allgather( sendbuf, &
3808  1, &
3809  mpi_real, &
3810  recvbuf, &
3811  1, &
3812  mpi_real, &
3814  ierr )
3815 
3816  globalmin = minval( recvbuf(:) )
3817 
3818  return

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

◆ 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 3823 of file scale_comm_icoA.F90.

3823  use scale_prc, only: &
3825  prc_nprocs
3826  implicit none
3827 
3828  real(DP), intent(in) :: localmin
3829  real(DP), intent(out) :: globalmin
3830 
3831  real(DP) :: sendbuf(1)
3832  real(DP) :: recvbuf(PRC_nprocs)
3833 
3834  integer :: ierr
3835  !---------------------------------------------------------------------------
3836 
3837  sendbuf(1) = localmin
3838 
3839  call mpi_allgather( sendbuf, &
3840  1, &
3841  mpi_double_precision, &
3842  recvbuf, &
3843  1, &
3844  mpi_double_precision, &
3846  ierr )
3847 
3848  globalmin = minval( recvbuf(:) )
3849 
3850  return

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

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.

79  integer, public :: COMM_datatype

Referenced by comm_setup().

◆ copy_info_r2r

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

Definition at line 82 of file scale_comm_icoA.F90.

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

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

◆ 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.

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

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

◆ send_info_r2r

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

Definition at line 84 of file scale_comm_icoA.F90.

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

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

◆ copy_info_p2r

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

Definition at line 86 of file scale_comm_icoA.F90.

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

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

◆ recv_info_p2r

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

Definition at line 87 of file scale_comm_icoA.F90.

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

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

◆ send_info_p2r

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

Definition at line 88 of file scale_comm_icoA.F90.

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

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

◆ copy_info_r2p

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

Definition at line 90 of file scale_comm_icoA.F90.

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

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

◆ recv_info_r2p

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

Definition at line 91 of file scale_comm_icoA.F90.

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

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

◆ send_info_r2p

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

Definition at line 92 of file scale_comm_icoA.F90.

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

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

◆ singular_info

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

Definition at line 94 of file scale_comm_icoA.F90.

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

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

◆ copy_list_r2r

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

Definition at line 97 of file scale_comm_icoA.F90.

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

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

◆ 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.

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

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

◆ send_list_r2r

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

Definition at line 99 of file scale_comm_icoA.F90.

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

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

◆ copy_list_p2r

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

Definition at line 101 of file scale_comm_icoA.F90.

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

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

◆ recv_list_p2r

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

Definition at line 102 of file scale_comm_icoA.F90.

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

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

◆ send_list_p2r

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

Definition at line 103 of file scale_comm_icoA.F90.

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

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

◆ copy_list_r2p

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

Definition at line 105 of file scale_comm_icoA.F90.

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

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

◆ recv_list_r2p

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

Definition at line 106 of file scale_comm_icoA.F90.

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

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

◆ send_list_r2p

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

Definition at line 107 of file scale_comm_icoA.F90.

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

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

◆ singular_list

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

Definition at line 109 of file scale_comm_icoA.F90.

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

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

◆ req_list

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

Definition at line 112 of file scale_comm_icoA.F90.

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

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

◆ sendbuf_r2r_sp

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

Definition at line 114 of file scale_comm_icoA.F90.

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

Referenced by comm_data_transfer_sp(), and comm_setup().

◆ recvbuf_r2r_sp

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

Definition at line 115 of file scale_comm_icoA.F90.

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

Referenced by comm_data_transfer_sp(), and comm_setup().

◆ sendbuf_p2r_sp

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

Definition at line 116 of file scale_comm_icoA.F90.

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

Referenced by comm_data_transfer_sp(), and comm_setup().

◆ recvbuf_p2r_sp

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

Definition at line 117 of file scale_comm_icoA.F90.

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

Referenced by comm_data_transfer_sp(), and comm_setup().

◆ sendbuf_r2p_sp

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

Definition at line 118 of file scale_comm_icoA.F90.

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

Referenced by comm_data_transfer_sp(), and comm_setup().

◆ recvbuf_r2p_sp

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

Definition at line 119 of file scale_comm_icoA.F90.

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

Referenced by comm_data_transfer_sp(), and comm_setup().

◆ sendbuf_r2r_dp

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

Definition at line 121 of file scale_comm_icoA.F90.

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

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

◆ recvbuf_r2r_dp

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

Definition at line 122 of file scale_comm_icoA.F90.

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

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

◆ sendbuf_p2r_dp

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

Definition at line 123 of file scale_comm_icoA.F90.

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

Referenced by comm_data_transfer_dp(), and comm_setup().

◆ recvbuf_p2r_dp

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

Definition at line 124 of file scale_comm_icoA.F90.

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

Referenced by comm_data_transfer_dp(), and comm_setup().

◆ sendbuf_r2p_dp

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

Definition at line 125 of file scale_comm_icoA.F90.

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

Referenced by comm_data_transfer_dp(), and comm_setup().

◆ recvbuf_r2p_dp

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

Definition at line 126 of file scale_comm_icoA.F90.

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

Referenced by comm_data_transfer_dp(), and comm_setup().

scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_prc_icoa::prc_rgn_total_pl
integer, parameter, public prc_rgn_total_pl
number of pole region
Definition: scale_prc_icoA.F90:73
scale_prc_icoa::prc_rgn_r2p_pl
integer, dimension(prc_rgn_total_pl), public prc_rgn_r2p_pl
process ID which have the pole regions
Definition: scale_prc_icoA.F90:88
scale_prc::prc_local_comm_world
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:88
scale_prc_icoa::i_npl
integer, parameter, public i_npl
north pole
Definition: scale_prc_icoA.F90:62
scale_prc_icoa
module process / icoA
Definition: scale_prc_icoA.F90:11
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_prc_icoa::prc_rgn_lp2r
integer, dimension(:,:), allocatable, public prc_rgn_lp2r
l,prc => rgn
Definition: scale_prc_icoA.F90:84
scale_prc_icoa::i_spl
integer, parameter, public i_spl
south pole
Definition: scale_prc_icoA.F90:63
scale_prc_icoa::prc_rgn_rgn4pl
integer, dimension(prc_rgn_total_pl), public prc_rgn_rgn4pl
region, having pole data in the halo
Definition: scale_prc_icoA.F90:89