New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 5619 for branches/NERC/dev_r5589_is_oce_cpl/NEMOGCM/NEMO/OPA_SRC/LBC/lib_mpp.F90 – NEMO

Ignore:
Timestamp:
2015-07-20T19:43:15+02:00 (9 years ago)
Author:
mathiot
Message:

ocean/ice sheet coupling: initial commit

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/NERC/dev_r5589_is_oce_cpl/NEMOGCM/NEMO/OPA_SRC/LBC/lib_mpp.F90

    r5429 r5619  
    7272   PUBLIC   mpp_lnk_3d, mpp_lnk_3d_gather, mpp_lnk_2d, mpp_lnk_2d_e 
    7373   PUBLIC   mpp_lnk_2d_9  
     74   PUBLIC   mpp_sum_3d, mpp_sum_2d 
    7475   PUBLIC   mppscatter, mppgather 
    7576   PUBLIC   mpp_ini_ice, mpp_ini_znl 
     
    13941395   END SUBROUTINE mpp_lnk_2d_e 
    13951396 
     1397   SUBROUTINE mpp_sum_3d( ptab, cd_type, psgn, cd_mpp, pval ) 
     1398      !!---------------------------------------------------------------------- 
     1399      !!                  ***  routine mpp_sum_3d  *** 
     1400      !! 
     1401      !! ** Purpose :   Message passing manadgement 
     1402      !! 
     1403      !! ** Method  :   Use mppsend and mpprecv function for passing mask 
     1404      !!      between processors following neighboring subdomains. 
     1405      !!            domain parameters 
     1406      !!                    nlci   : first dimension of the local subdomain 
     1407      !!                    nlcj   : second dimension of the local subdomain 
     1408      !!                    nbondi : mark for "east-west local boundary" 
     1409      !!                    nbondj : mark for "north-south local boundary" 
     1410      !!                    noea   : number for local neighboring processors 
     1411      !!                    nowe   : number for local neighboring processors 
     1412      !!                    noso   : number for local neighboring processors 
     1413      !!                    nono   : number for local neighboring processors 
     1414      !! 
     1415      !! ** Action  :   ptab with update value at its periphery 
     1416      !! 
     1417      !!---------------------------------------------------------------------- 
     1418      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   ptab     ! 3D array on which the boundary condition is applied 
     1419      CHARACTER(len=1)                , INTENT(in   ) ::   cd_type  ! define the nature of ptab array grid-points 
     1420      !                                                             ! = T , U , V , F , W points 
     1421      REAL(wp)                        , INTENT(in   ) ::   psgn     ! =-1 the sign change across the north fold boundary 
     1422      !                                                             ! =  1. , the sign is kept 
     1423      CHARACTER(len=3), OPTIONAL      , INTENT(in   ) ::   cd_mpp   ! fill the overlap area only 
     1424      REAL(wp)        , OPTIONAL      , INTENT(in   ) ::   pval     ! background value (used at closed boundaries) 
     1425      !! 
     1426      INTEGER  ::   ji, jj, jk, jl             ! dummy loop indices 
     1427      INTEGER  ::   imigr, iihom, ijhom        ! temporary integers 
     1428      INTEGER  ::   ml_req1, ml_req2, ml_err   ! for key_mpi_isend 
     1429      REAL(wp) ::   zland 
     1430      INTEGER, DIMENSION(MPI_STATUS_SIZE) ::   ml_stat   ! for key_mpi_isend 
     1431      ! 
     1432      REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE ::   zt3ns, zt3sn   ! 3d for north-south & south-north 
     1433      REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE ::   zt3ew, zt3we   ! 3d for east-west & west-east 
     1434 
     1435      !!---------------------------------------------------------------------- 
     1436       
     1437      ALLOCATE( zt3ns(jpi,jprecj,jpk,2), zt3sn(jpi,jprecj,jpk,2),   & 
     1438         &      zt3ew(jpj,jpreci,jpk,2), zt3we(jpj,jpreci,jpk,2)  ) 
     1439 
     1440      ! 
     1441      IF( PRESENT( pval ) ) THEN   ;   zland = pval      ! set land value 
     1442      ELSE                         ;   zland = 0.e0      ! zero by default 
     1443      ENDIF 
     1444 
     1445      ! 1. standard boundary treatment 
     1446      ! ------------------------------ 
     1447      IF( PRESENT( cd_mpp ) ) THEN      ! only fill added line/raw with existing values 
     1448         ! 
     1449         ! WARNING ptab is defined only between nld and nle 
     1450!         DO jk = 1, jpk 
     1451!            DO jj = nlcj+1, jpj                 ! added line(s)   (inner only) 
     1452!               ptab(nldi  :nlei  , jj          ,jk) = ptab(nldi:nlei,     nlej,jk) 
     1453!               ptab(1     :nldi-1, jj          ,jk) = ptab(nldi     ,     nlej,jk) 
     1454!               ptab(nlei+1:nlci  , jj          ,jk) = ptab(     nlei,     nlej,jk) 
     1455!            END DO 
     1456!            DO ji = nlci+1, jpi                 ! added column(s) (full) 
     1457!               ptab(ji           ,nldj  :nlej  ,jk) = ptab(     nlei,nldj:nlej,jk) 
     1458!               ptab(ji           ,1     :nldj-1,jk) = ptab(     nlei,nldj     ,jk) 
     1459!               ptab(ji           ,nlej+1:jpj   ,jk) = ptab(     nlei,     nlej,jk) 
     1460!            END DO 
     1461!         END DO 
     1462         ! 
     1463      ELSE                              ! standard close or cyclic treatment 
     1464         ! 
     1465         !                                   ! East-West boundaries 
     1466         !                                        !* Cyclic east-west 
     1467         IF( nbondi == 2 .AND. (nperio == 1 .OR. nperio == 4 .OR. nperio == 6) ) THEN 
     1468!            ptab( 1 ,:,:) = ptab(jpim1,:,:) 
     1469!            ptab(jpi,:,:) = ptab(  2  ,:,:) 
     1470         ELSE                                     !* closed 
     1471!            IF( .NOT. cd_type == 'F' )   ptab(     1       :jpreci,:,:) = zland    ! south except F-point 
     1472!                                         ptab(nlci-jpreci+1:jpi   ,:,:) = zland    ! north 
     1473         ENDIF 
     1474         !                                   ! North-South boundaries (always closed) 
     1475!         IF( .NOT. cd_type == 'F' )   ptab(:,     1       :jprecj,:) = zland       ! south except F-point 
     1476!                                      ptab(:,nlcj-jprecj+1:jpj   ,:) = zland       ! north 
     1477         ! 
     1478      ENDIF 
     1479 
     1480      ! 2. East and west directions exchange 
     1481      ! ------------------------------------ 
     1482      ! we play with the neigbours AND the row number because of the periodicity 
     1483      ! 
     1484      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions 
     1485      CASE ( -1, 0, 1 )                ! all exept 2 (i.e. close case) 
     1486      iihom = nlci-jpreci 
     1487         DO jl = 1, jpreci 
     1488            zt3ew(:,jl,:,1) = ptab(jl      ,:,:) ; ptab(jl      ,:,:) = 0.0_wp 
     1489            zt3we(:,jl,:,1) = ptab(iihom+jl,:,:) ; ptab(iihom+jl,:,:) = 0.0_wp  
     1490         END DO 
     1491      END SELECT 
     1492      ! 
     1493      !                           ! Migrations 
     1494      imigr = jpreci * jpj * jpk 
     1495      ! 
     1496      SELECT CASE ( nbondi ) 
     1497      CASE ( -1 ) 
     1498         CALL mppsend( 2, zt3we(1,1,1,1), imigr, noea, ml_req1 ) 
     1499         CALL mpprecv( 1, zt3ew(1,1,1,2), imigr, noea ) 
     1500         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err) 
     1501      CASE ( 0 ) 
     1502         CALL mppsend( 1, zt3ew(1,1,1,1), imigr, nowe, ml_req1 ) 
     1503         CALL mppsend( 2, zt3we(1,1,1,1), imigr, noea, ml_req2 ) 
     1504         CALL mpprecv( 1, zt3ew(1,1,1,2), imigr, noea ) 
     1505         CALL mpprecv( 2, zt3we(1,1,1,2), imigr, nowe ) 
     1506         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err) 
     1507         IF(l_isend) CALL mpi_wait(ml_req2, ml_stat, ml_err) 
     1508      CASE ( 1 ) 
     1509         CALL mppsend( 1, zt3ew(1,1,1,1), imigr, nowe, ml_req1 ) 
     1510         CALL mpprecv( 2, zt3we(1,1,1,2), imigr, nowe ) 
     1511         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err) 
     1512      END SELECT 
     1513      ! 
     1514      !                           ! Write Dirichlet lateral conditions 
     1515      iihom = nlci-nreci 
     1516      ! 
     1517      SELECT CASE ( nbondi ) 
     1518      CASE ( -1 ) 
     1519         DO jl = 1, jpreci 
     1520            ptab(iihom+jl,:,:) = ptab(iihom+jl,:,:) + zt3ew(:,jl,:,2) 
     1521         END DO 
     1522      CASE ( 0 ) 
     1523         DO jl = 1, jpreci 
     1524            ptab(jpreci+jl,:,:) = ptab(jpreci+jl,:,:) + zt3we(:,jl,:,2) 
     1525            ptab(iihom +jl,:,:) = ptab(iihom +jl,:,:) + zt3ew(:,jl,:,2) 
     1526         END DO 
     1527      CASE ( 1 ) 
     1528         DO jl = 1, jpreci 
     1529            ptab(jpreci+jl,:,:) = ptab(jpreci+jl,:,:) + zt3we(:,jl,:,2) 
     1530         END DO 
     1531      END SELECT 
     1532 
     1533 
     1534      ! 3. North and south directions 
     1535      ! ----------------------------- 
     1536      ! always closed : we play only with the neigbours 
     1537      ! 
     1538      IF( nbondj /= 2 ) THEN      ! Read Dirichlet lateral conditions 
     1539         ijhom = nlcj-jprecj 
     1540         DO jl = 1, jprecj 
     1541            zt3sn(:,jl,:,1) = ptab(:,ijhom+jl,:) ; ptab(:,ijhom+jl,:) = 0.0_wp 
     1542            zt3ns(:,jl,:,1) = ptab(:,jl      ,:) ; ptab(:,jl      ,:) = 0.0_wp 
     1543         END DO 
     1544      ENDIF 
     1545      ! 
     1546      !                           ! Migrations 
     1547      imigr = jprecj * jpi * jpk 
     1548      ! 
     1549      SELECT CASE ( nbondj ) 
     1550      CASE ( -1 ) 
     1551         CALL mppsend( 4, zt3sn(1,1,1,1), imigr, nono, ml_req1 ) 
     1552         CALL mpprecv( 3, zt3ns(1,1,1,2), imigr, nono ) 
     1553         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err) 
     1554      CASE ( 0 ) 
     1555         CALL mppsend( 3, zt3ns(1,1,1,1), imigr, noso, ml_req1 ) 
     1556         CALL mppsend( 4, zt3sn(1,1,1,1), imigr, nono, ml_req2 ) 
     1557         CALL mpprecv( 3, zt3ns(1,1,1,2), imigr, nono ) 
     1558         CALL mpprecv( 4, zt3sn(1,1,1,2), imigr, noso ) 
     1559         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err) 
     1560         IF(l_isend) CALL mpi_wait(ml_req2, ml_stat, ml_err) 
     1561      CASE ( 1 ) 
     1562         CALL mppsend( 3, zt3ns(1,1,1,1), imigr, noso, ml_req1 ) 
     1563         CALL mpprecv( 4, zt3sn(1,1,1,2), imigr, noso ) 
     1564         IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err) 
     1565      END SELECT 
     1566      ! 
     1567      !                           ! Write Dirichlet lateral conditions 
     1568      ijhom = nlcj-nrecj 
     1569      ! 
     1570      SELECT CASE ( nbondj ) 
     1571      CASE ( -1 ) 
     1572         DO jl = 1, jprecj 
     1573            ptab(:,ijhom+jl,:) = ptab(:,ijhom+jl,:) + zt3ns(:,jl,:,2) 
     1574         END DO 
     1575      CASE ( 0 ) 
     1576         DO jl = 1, jprecj 
     1577            ptab(:,jprecj+jl,:) = ptab(:,jprecj+jl,:) + zt3sn(:,jl,:,2) 
     1578            ptab(:,ijhom +jl,:) = ptab(:,ijhom +jl,:) + zt3ns(:,jl,:,2) 
     1579         END DO 
     1580      CASE ( 1 ) 
     1581         DO jl = 1, jprecj 
     1582            ptab(:,jprecj+jl,:) = ptab(:,jprecj+jl,:) + zt3sn(:,jl   ,:,2) 
     1583         END DO 
     1584      END SELECT 
     1585 
     1586 
     1587      ! 4. north fold treatment 
     1588      ! ----------------------- 
     1589      ! 
     1590      IF( npolj /= 0 .AND. .NOT. PRESENT(cd_mpp) ) THEN 
     1591         ! 
     1592         SELECT CASE ( jpni ) 
     1593         CASE ( 1 )     ;   CALL lbc_nfd      ( ptab, cd_type, psgn )   ! only 1 northern proc, no mpp 
     1594         CASE DEFAULT   ;   CALL mpp_lbc_north( ptab, cd_type, psgn )   ! for all northern procs. 
     1595         END SELECT 
     1596         ! 
     1597      ENDIF 
     1598      ! 
     1599      DEALLOCATE( zt3ns, zt3sn, zt3ew, zt3we ) 
     1600      ! 
     1601   END SUBROUTINE mpp_sum_3d 
     1602 
     1603   SUBROUTINE mpp_sum_2d( pt2d, cd_type, psgn, cd_mpp, pval ) 
     1604      !!---------------------------------------------------------------------- 
     1605      !!                  ***  routine mpp_sum_2d  *** 
     1606      !! 
     1607      !! ** Purpose :   Message passing manadgement for 2d array 
     1608      !! 
     1609      !! ** Method  :   Use mppsend and mpprecv function for passing mask 
     1610      !!      between processors following neighboring subdomains. 
     1611      !!            domain parameters 
     1612      !!                    nlci   : first dimension of the local subdomain 
     1613      !!                    nlcj   : second dimension of the local subdomain 
     1614      !!                    nbondi : mark for "east-west local boundary" 
     1615      !!                    nbondj : mark for "north-south local boundary" 
     1616      !!                    noea   : number for local neighboring processors 
     1617      !!                    nowe   : number for local neighboring processors 
     1618      !!                    noso   : number for local neighboring processors 
     1619      !!                    nono   : number for local neighboring processors 
     1620      !! 
     1621      !!---------------------------------------------------------------------- 
     1622      REAL(wp), DIMENSION(jpi,jpj), INTENT(inout) ::   pt2d     ! 2D array on which the boundary condition is applied 
     1623      CHARACTER(len=1)            , INTENT(in   ) ::   cd_type  ! define the nature of ptab array grid-points 
     1624      !                                                         ! = T , U , V , F , W and I points 
     1625      REAL(wp)                    , INTENT(in   ) ::   psgn     ! =-1 the sign change across the north fold boundary 
     1626      !                                                         ! =  1. , the sign is kept 
     1627      CHARACTER(len=3), OPTIONAL  , INTENT(in   ) ::   cd_mpp   ! fill the overlap area only 
     1628      REAL(wp)        , OPTIONAL  , INTENT(in   ) ::   pval     ! background value (used at closed boundaries) 
     1629      !! 
     1630      INTEGER  ::   ji, jj, jl   ! dummy loop indices 
     1631      INTEGER  ::   imigr, iihom, ijhom        ! temporary integers 
     1632      INTEGER  ::   ml_req1, ml_req2, ml_err   ! for key_mpi_isend 
     1633      REAL(wp) ::   zland 
     1634      INTEGER, DIMENSION(MPI_STATUS_SIZE) ::   ml_stat   ! for key_mpi_isend 
     1635      ! 
     1636      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  zt2ns, zt2sn   ! 2d for north-south & south-north 
     1637      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  zt2ew, zt2we   ! 2d for east-west & west-east 
     1638 
     1639      !!---------------------------------------------------------------------- 
     1640 
     1641      ALLOCATE( zt2ns(jpi,jprecj,2), zt2sn(jpi,jprecj,2),  & 
     1642         &      zt2ew(jpj,jpreci,2), zt2we(jpj,jpreci,2)   ) 
     1643 
     1644      ! 
     1645      IF( PRESENT( pval ) ) THEN   ;   zland = pval      ! set land value 
     1646      ELSE                         ;   zland = 0.e0      ! zero by default 
     1647      ENDIF 
     1648 
     1649      ! 1. standard boundary treatment 
     1650      ! ------------------------------ 
     1651      ! 
     1652!      IF( PRESENT( cd_mpp ) ) THEN      ! only fill added line/raw with existing values 
     1653!         ! 
     1654!         ! WARNING pt2d is defined only between nld and nle 
     1655!         DO jj = nlcj+1, jpj                 ! added line(s)   (inner only) 
     1656!            pt2d(nldi  :nlei  , jj          ) = pt2d(nldi:nlei,     nlej) 
     1657!            pt2d(1     :nldi-1, jj          ) = pt2d(nldi     ,     nlej) 
     1658!            pt2d(nlei+1:nlci  , jj          ) = pt2d(     nlei,     nlej) 
     1659!         END DO 
     1660!         DO ji = nlci+1, jpi                 ! added column(s) (full) 
     1661!            pt2d(ji           ,nldj  :nlej  ) = pt2d(     nlei,nldj:nlej) 
     1662!            pt2d(ji           ,1     :nldj-1) = pt2d(     nlei,nldj     ) 
     1663!            pt2d(ji           ,nlej+1:jpj   ) = pt2d(     nlei,     nlej) 
     1664!         END DO 
     1665!         ! 
     1666!      ELSE                              ! standard close or cyclic treatment 
     1667!         ! 
     1668!         !                                   ! East-West boundaries 
     1669!         IF( nbondi == 2 .AND.   &                ! Cyclic east-west 
     1670!            &    (nperio == 1 .OR. nperio == 4 .OR. nperio == 6) ) THEN 
     1671!            pt2d( 1 ,:) = pt2d(jpim1,:)                                    ! west 
     1672!            pt2d(jpi,:) = pt2d(  2  ,:)                                    ! east 
     1673!         ELSE                                     ! closed 
     1674!            IF( .NOT. cd_type == 'F' )   pt2d(     1       :jpreci,:) = zland    ! south except F-point 
     1675!                                         pt2d(nlci-jpreci+1:jpi   ,:) = zland    ! north 
     1676!         ENDIF 
     1677!         !                                   ! North-South boundaries (always closed) 
     1678!            IF( .NOT. cd_type == 'F' )   pt2d(:,     1       :jprecj) = zland    !south except F-point 
     1679!                                         pt2d(:,nlcj-jprecj+1:jpj   ) = zland    ! north 
     1680!         ! 
     1681!      ENDIF 
     1682 
     1683      ! 2. East and west directions exchange 
     1684      ! ------------------------------------ 
     1685      ! we play with the neigbours AND the row number because of the periodicity 
     1686      ! 
     1687      SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions 
     1688      CASE ( -1, 0, 1 )                ! all exept 2 (i.e. close case) 
     1689         iihom = nlci - jpreci 
     1690         DO jl = 1, jpreci 
     1691            zt2ew(:,jl,1) = pt2d(jl       ,:) ; pt2d(jl       ,:) = 0.0_wp 
     1692            zt2we(:,jl,1) = pt2d(iihom +jl,:) ; pt2d(iihom +jl,:) = 0.0_wp 
     1693         END DO 
     1694      END SELECT 
     1695      ! 
     1696      !                           ! Migrations 
     1697      imigr = jpreci * jpj 
     1698      ! 
     1699      SELECT CASE ( nbondi ) 
     1700      CASE ( -1 ) 
     1701         CALL mppsend( 2, zt2we(1,1,1), imigr, noea, ml_req1 ) 
     1702         CALL mpprecv( 1, zt2ew(1,1,2), imigr, noea ) 
     1703         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 
     1704      CASE ( 0 ) 
     1705         CALL mppsend( 1, zt2ew(1,1,1), imigr, nowe, ml_req1 ) 
     1706         CALL mppsend( 2, zt2we(1,1,1), imigr, noea, ml_req2 ) 
     1707         CALL mpprecv( 1, zt2ew(1,1,2), imigr, noea ) 
     1708         CALL mpprecv( 2, zt2we(1,1,2), imigr, nowe ) 
     1709         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 
     1710         IF(l_isend) CALL mpi_wait(ml_req2,ml_stat,ml_err) 
     1711      CASE ( 1 ) 
     1712         CALL mppsend( 1, zt2ew(1,1,1), imigr, nowe, ml_req1 ) 
     1713         CALL mpprecv( 2, zt2we(1,1,2), imigr, nowe ) 
     1714         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 
     1715      END SELECT 
     1716      ! 
     1717      !                           ! Write Dirichlet lateral conditions 
     1718      iihom = nlci-nreci 
     1719      ! 
     1720      SELECT CASE ( nbondi ) 
     1721      CASE ( -1 ) 
     1722         DO jl = 1, jpreci 
     1723            pt2d(iihom+jl,:) = pt2d(iihom+jl,:) + zt2ew(:,jl,2) 
     1724         END DO 
     1725      CASE ( 0 ) 
     1726         DO jl = 1, jpreci 
     1727            pt2d(jpreci+jl,:) = pt2d(jpreci+jl,:) + zt2we(:,jl,2) 
     1728            pt2d(iihom +jl,:) = pt2d(iihom +jl,:) + zt2ew(:,jl,2) 
     1729         END DO 
     1730      CASE ( 1 ) 
     1731         DO jl = 1, jpreci 
     1732            pt2d(jpreci+jl,:) = pt2d(jpreci+jl,:) + zt2we(:,jl,2) 
     1733         END DO 
     1734      END SELECT 
     1735 
     1736 
     1737      ! 3. North and south directions 
     1738      ! ----------------------------- 
     1739      ! always closed : we play only with the neigbours 
     1740      ! 
     1741      IF( nbondj /= 2 ) THEN      ! Read Dirichlet lateral conditions 
     1742         ijhom = nlcj - jprecj 
     1743         DO jl = 1, jprecj 
     1744            zt2sn(:,jl,1) = pt2d(:,ijhom +jl) ; pt2d(:,ijhom +jl) = 0.0_wp 
     1745            zt2ns(:,jl,1) = pt2d(:,jl       ) ; pt2d(:,jl       ) = 0.0_wp 
     1746         END DO 
     1747      ENDIF 
     1748      ! 
     1749      !                           ! Migrations 
     1750      imigr = jprecj * jpi 
     1751      ! 
     1752      SELECT CASE ( nbondj ) 
     1753      CASE ( -1 ) 
     1754         CALL mppsend( 4, zt2sn(1,1,1), imigr, nono, ml_req1 ) 
     1755         CALL mpprecv( 3, zt2ns(1,1,2), imigr, nono ) 
     1756         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 
     1757      CASE ( 0 ) 
     1758         CALL mppsend( 3, zt2ns(1,1,1), imigr, noso, ml_req1 ) 
     1759         CALL mppsend( 4, zt2sn(1,1,1), imigr, nono, ml_req2 ) 
     1760         CALL mpprecv( 3, zt2ns(1,1,2), imigr, nono ) 
     1761         CALL mpprecv( 4, zt2sn(1,1,2), imigr, noso ) 
     1762         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 
     1763         IF(l_isend) CALL mpi_wait(ml_req2,ml_stat,ml_err) 
     1764      CASE ( 1 ) 
     1765         CALL mppsend( 3, zt2ns(1,1,1), imigr, noso, ml_req1 ) 
     1766         CALL mpprecv( 4, zt2sn(1,1,2), imigr, noso ) 
     1767         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 
     1768      END SELECT 
     1769      ! 
     1770      !                           ! Write Dirichlet lateral conditions 
     1771      ijhom = nlcj-nrecj 
     1772      ! 
     1773      SELECT CASE ( nbondj ) 
     1774      CASE ( -1 ) 
     1775         DO jl = 1, jprecj 
     1776            pt2d(:,ijhom+jl) = pt2d(:,ijhom+jl) + zt2ns(:,jl,2) 
     1777         END DO 
     1778      CASE ( 0 ) 
     1779         DO jl = 1, jprecj 
     1780            pt2d(:,jprecj+jl) = pt2d(:,jprecj+jl) + zt2sn(:,jl,2) 
     1781            pt2d(:,ijhom +jl) = pt2d(:,ijhom +jl) + zt2ns(:,jl,2) 
     1782         END DO 
     1783      CASE ( 1 ) 
     1784         DO jl = 1, jprecj 
     1785            pt2d(:,jprecj+jl) = pt2d(:,jprecj+jl) + zt2sn(:,jl,2) 
     1786         END DO 
     1787      END SELECT 
     1788 
     1789 
     1790      ! 4. north fold treatment 
     1791      ! ----------------------- 
     1792      ! 
     1793      IF( npolj /= 0 .AND. .NOT. PRESENT(cd_mpp) ) THEN 
     1794         ! 
     1795         SELECT CASE ( jpni ) 
     1796         CASE ( 1 )     ;   CALL lbc_nfd      ( pt2d, cd_type, psgn )   ! only 1 northern proc, no mpp 
     1797         CASE DEFAULT   ;   CALL mpp_lbc_north( pt2d, cd_type, psgn )   ! for all northern procs. 
     1798         END SELECT 
     1799         ! 
     1800      ENDIF 
     1801      ! 
     1802      DEALLOCATE( zt2ns, zt2sn, zt2ew, zt2we ) 
     1803      ! 
     1804   END SUBROUTINE mpp_sum_2d 
    13961805 
    13971806   SUBROUTINE mppsend( ktyp, pmess, kbytes, kdest, md_req ) 
Note: See TracChangeset for help on using the changeset viewer.