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_ndiamond, 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 1567 of file scale_comm_icoA.F90.

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

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

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

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

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

3118  implicit none
3119 
3120  integer :: suffix
3121  integer :: i, j
3122  !---------------------------------------------------------------------------
3123 
3124  suffix = adm_gall_1d * (j-1) + i
3125 

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

3486  use scale_prc, only: &
3488  prc_nprocs
3489  implicit none
3490 
3491  real(SP), intent(in) :: localsum
3492  real(SP), intent(out) :: globalsum
3493 
3494  real(SP) :: sendbuf(1)
3495  real(SP) :: recvbuf(PRC_nprocs)
3496 
3497  integer :: ierr
3498  !---------------------------------------------------------------------------
3499 
3500  if ( comm_pl ) then
3501  sendbuf(1) = localsum
3502 
3503  call mpi_allgather( sendbuf, &
3504  1, &
3505  mpi_real, &
3506  recvbuf, &
3507  1, &
3508  mpi_real, &
3510  ierr )
3511 
3512  globalsum = sum( recvbuf(:) )
3513  else
3514  globalsum = localsum
3515  endif
3516 
3517  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 3522 of file scale_comm_icoA.F90.

3522  use scale_prc, only: &
3524  prc_nprocs
3525  implicit none
3526 
3527  real(DP), intent(in) :: localsum
3528  real(DP), intent(out) :: globalsum
3529 
3530  real(DP) :: sendbuf(1)
3531  real(DP) :: recvbuf(PRC_nprocs)
3532 
3533  integer :: ierr
3534  !---------------------------------------------------------------------------
3535 
3536  if ( comm_pl ) then
3537  sendbuf(1) = localsum
3538 
3539  call mpi_allgather( sendbuf, &
3540  1, &
3541  mpi_double_precision, &
3542  recvbuf, &
3543  1, &
3544  mpi_double_precision, &
3546  ierr )
3547 
3548  globalsum = sum( recvbuf(:) )
3549  else
3550  globalsum = localsum
3551  endif
3552 
3553  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 3558 of file scale_comm_icoA.F90.

3558  use scale_prc, only: &
3560  prc_nprocs
3561  implicit none
3562 
3563  integer, intent(in) :: kall
3564  real(SP), intent(in) :: localsum (kall)
3565  real(SP), intent(out) :: globalsum(kall)
3566 
3567  real(SP) :: sendbuf(kall)
3568  integer :: displs (PRC_nprocs)
3569  integer :: counts (PRC_nprocs)
3570  real(SP) :: recvbuf(kall,PRC_nprocs)
3571 
3572  integer :: ierr
3573  integer :: k, p
3574  !---------------------------------------------------------------------------
3575 
3576  do p = 1, prc_nprocs
3577  displs(p) = (p-1) * kall
3578  counts(p) = kall
3579  enddo
3580 
3581  if ( comm_pl ) then
3582  sendbuf(:) = localsum(:)
3583 
3584  call mpi_allgatherv( sendbuf, &
3585  kall, &
3586  mpi_real, &
3587  recvbuf, &
3588  counts, &
3589  displs, &
3590  mpi_real, &
3592  ierr )
3593 
3594  do k = 1, kall
3595  globalsum(k) = sum( recvbuf(k,:) )
3596  enddo
3597  else
3598  do k = 1, kall
3599  globalsum(k) = localsum(k)
3600  enddo
3601  endif
3602 
3603  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 3608 of file scale_comm_icoA.F90.

3608  use scale_prc, only: &
3610  prc_nprocs
3611  implicit none
3612 
3613  integer, intent(in) :: kall
3614  real(DP), intent(in) :: localsum (kall)
3615  real(DP), intent(out) :: globalsum(kall)
3616 
3617  real(DP) :: sendbuf(kall)
3618  integer :: displs (PRC_nprocs)
3619  integer :: counts (PRC_nprocs)
3620  real(DP) :: recvbuf(kall,PRC_nprocs)
3621 
3622  integer :: ierr
3623  integer :: k, p
3624  !---------------------------------------------------------------------------
3625 
3626  do p = 1, prc_nprocs
3627  displs(p) = (p-1) * kall
3628  counts(p) = kall
3629  enddo
3630 
3631  if ( comm_pl ) then
3632  sendbuf(:) = localsum(:)
3633 
3634  call mpi_allgatherv( sendbuf, &
3635  kall, &
3636  mpi_double_precision, &
3637  recvbuf, &
3638  counts, &
3639  displs, &
3640  mpi_double_precision, &
3642  ierr )
3643 
3644  do k = 1, kall
3645  globalsum(k) = sum( recvbuf(k,:) )
3646  enddo
3647  else
3648  do k = 1, kall
3649  globalsum(k) = localsum(k)
3650  enddo
3651  endif
3652 
3653  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 3658 of file scale_comm_icoA.F90.

3658  use scale_prc, only: &
3660  prc_nprocs
3661  implicit none
3662 
3663  real(SP), intent(in) :: localavg
3664  real(SP), intent(out) :: globalavg
3665 
3666  real(SP) :: sendbuf(1)
3667  real(SP) :: recvbuf(PRC_nprocs)
3668 
3669  integer :: ierr
3670  !---------------------------------------------------------------------------
3671 
3672  if ( comm_pl ) then
3673  sendbuf(1) = localavg
3674 
3675  call mpi_allgather( sendbuf, &
3676  1, &
3677  mpi_real, &
3678  recvbuf, &
3679  1, &
3680  mpi_real, &
3682  ierr )
3683 
3684  globalavg = sum( recvbuf(:) ) / real(prc_nprocs,kind=sp)
3685  else
3686  globalavg = localavg
3687  endif
3688 
3689  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 3694 of file scale_comm_icoA.F90.

3694  use scale_prc, only: &
3696  prc_nprocs
3697  implicit none
3698 
3699  real(DP), intent(in) :: localavg
3700  real(DP), intent(out) :: globalavg
3701 
3702  real(DP) :: sendbuf(1)
3703  real(DP) :: recvbuf(PRC_nprocs)
3704 
3705  integer :: ierr
3706  !---------------------------------------------------------------------------
3707 
3708  if ( comm_pl ) then
3709  sendbuf(1) = localavg
3710 
3711  call mpi_allgather( sendbuf, &
3712  1, &
3713  mpi_double_precision, &
3714  recvbuf, &
3715  1, &
3716  mpi_double_precision, &
3718  ierr )
3719 
3720  globalavg = sum( recvbuf(:) ) / real(prc_nprocs,kind=dp)
3721  else
3722  globalavg = localavg
3723  endif
3724 
3725  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 3730 of file scale_comm_icoA.F90.

3730  use scale_prc, only: &
3732  prc_nprocs
3733  implicit none
3734 
3735  real(SP), intent(in) :: localmax
3736  real(SP), intent(out) :: globalmax
3737 
3738  real(SP) :: sendbuf(1)
3739  real(SP) :: recvbuf(PRC_nprocs)
3740 
3741  integer :: ierr
3742  !---------------------------------------------------------------------------
3743 
3744  sendbuf(1) = localmax
3745 
3746  call mpi_allgather( sendbuf, &
3747  1, &
3748  mpi_real, &
3749  recvbuf, &
3750  1, &
3751  mpi_real, &
3753  ierr )
3754 
3755  globalmax = maxval( recvbuf(:) )
3756 
3757  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 3762 of file scale_comm_icoA.F90.

3762  use scale_prc, only: &
3764  prc_nprocs
3765  implicit none
3766 
3767  real(DP), intent(in) :: localmax
3768  real(DP), intent(out) :: globalmax
3769 
3770  real(DP) :: sendbuf(1)
3771  real(DP) :: recvbuf(PRC_nprocs)
3772 
3773  integer :: ierr
3774  !---------------------------------------------------------------------------
3775 
3776  sendbuf(1) = localmax
3777 
3778  call mpi_allgather( sendbuf, &
3779  1, &
3780  mpi_double_precision, &
3781  recvbuf, &
3782  1, &
3783  mpi_double_precision, &
3785  ierr )
3786 
3787  globalmax = maxval( recvbuf(:) )
3788 
3789  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 3794 of file scale_comm_icoA.F90.

3794  use scale_prc, only: &
3796  prc_nprocs
3797  implicit none
3798 
3799  real(SP), intent(in) :: localmin
3800  real(SP), intent(out) :: globalmin
3801 
3802  real(SP) :: sendbuf(1)
3803  real(SP) :: recvbuf(PRC_nprocs)
3804 
3805  integer :: ierr
3806  !---------------------------------------------------------------------------
3807 
3808  sendbuf(1) = localmin
3809 
3810  call mpi_allgather( sendbuf, &
3811  1, &
3812  mpi_real, &
3813  recvbuf, &
3814  1, &
3815  mpi_real, &
3817  ierr )
3818 
3819  globalmin = minval( recvbuf(:) )
3820 
3821  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 3826 of file scale_comm_icoA.F90.

3826  use scale_prc, only: &
3828  prc_nprocs
3829  implicit none
3830 
3831  real(DP), intent(in) :: localmin
3832  real(DP), intent(out) :: globalmin
3833 
3834  real(DP) :: sendbuf(1)
3835  real(DP) :: recvbuf(PRC_nprocs)
3836 
3837  integer :: ierr
3838  !---------------------------------------------------------------------------
3839 
3840  sendbuf(1) = localmin
3841 
3842  call mpi_allgather( sendbuf, &
3843  1, &
3844  mpi_double_precision, &
3845  recvbuf, &
3846  1, &
3847  mpi_double_precision, &
3849  ierr )
3850 
3851  globalmin = minval( recvbuf(:) )
3852 
3853  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_precision::sp
integer, parameter, public sp
Definition: scale_precision.F90:31
scale_comm_cartesc::comm_datatype
integer, public comm_datatype
datatype of variable
Definition: scale_comm_cartesC.F90:105
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
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:89
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_atmos_grid_cartesc_index::imax
integer, public imax
Definition: scale_atmos_grid_cartesC_index.F90:37
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_atmos_grid_cartesc_index::kmax
integer, public kmax
Definition: scale_atmos_grid_cartesC_index.F90:36
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
scale_precision::dp
integer, parameter, public dp
Definition: scale_precision.F90:32
scale_prc::prc_nprocs
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:90