- Timestamp:
- 2015-12-04T17:56:07+01:00 (8 years ago)
- Location:
- branches/2015/dev_MetOffice_merge_2015/NEMOGCM/NEMO/OPA_SRC/LBC
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/2015/dev_MetOffice_merge_2015/NEMOGCM/NEMO/OPA_SRC/LBC/lbclnk.F90
r5429 r6006 17 17 !!---------------------------------------------------------------------- 18 18 !! lbc_lnk : generic interface for mpp_lnk_3d and mpp_lnk_2d routines defined in lib_mpp 19 !! lbc_sum : generic interface for mpp_lnk_sum_3d and mpp_lnk_sum_2d routines defined in lib_mpp 19 20 !! lbc_lnk_e : generic interface for mpp_lnk_2d_e routine defined in lib_mpp 20 21 !! lbc_bdy_lnk : generic interface for mpp_lnk_bdy_2d and mpp_lnk_bdy_3d routines defined in lib_mpp … … 31 32 END INTERFACE 32 33 34 !JMM interface not defined if not key_mpp_mpi : likely do not compile without this CPP key !!!! 35 INTERFACE lbc_sum 36 MODULE PROCEDURE mpp_lnk_sum_3d, mpp_lnk_sum_2d 37 END INTERFACE 38 33 39 INTERFACE lbc_bdy_lnk 34 40 MODULE PROCEDURE mpp_lnk_bdy_2d, mpp_lnk_bdy_3d … … 45 51 PUBLIC lbc_lnk ! ocean lateral boundary conditions 46 52 PUBLIC lbc_lnk_multi ! modified ocean lateral boundary conditions 53 PUBLIC lbc_sum 47 54 PUBLIC lbc_lnk_e 48 55 PUBLIC lbc_bdy_lnk ! ocean lateral BDY boundary conditions … … 59 66 !! Default option shared memory computing 60 67 !!---------------------------------------------------------------------- 61 !! lbc_lnk : generic interface for lbc_lnk_3d and lbc_lnk_2d 62 !! lbc_lnk_3d : set the lateral boundary condition on a 3D variable on ocean mesh 63 !! lbc_lnk_2d : set the lateral boundary condition on a 2D variable on ocean mesh 64 !! lbc_bdy_lnk : set the lateral BDY boundary condition 68 !! lbc_sum : generic interface for mpp_lnk_sum_3d and mpp_lnk_sum_2d 69 !! lbc_lnk_sum_3d: compute sum over the halos on a 3D variable on ocean mesh 70 !! lbc_lnk_sum_3d: compute sum over the halos on a 2D variable on ocean mesh 71 !! lbc_lnk : generic interface for lbc_lnk_3d and lbc_lnk_2d 72 !! lbc_lnk_3d : set the lateral boundary condition on a 3D variable on ocean mesh 73 !! lbc_lnk_2d : set the lateral boundary condition on a 2D variable on ocean mesh 74 !! lbc_bdy_lnk : set the lateral BDY boundary condition 65 75 !!---------------------------------------------------------------------- 66 76 USE oce ! ocean dynamics and tracers … … 74 84 INTERFACE lbc_lnk 75 85 MODULE PROCEDURE lbc_lnk_3d_gather, lbc_lnk_3d, lbc_lnk_2d 86 END INTERFACE 87 88 INTERFACE lbc_sum 89 MODULE PROCEDURE mpp_lnk_sum_3d, mpp_lnk_sum_2d 76 90 END INTERFACE 77 91 -
branches/2015/dev_MetOffice_merge_2015/NEMOGCM/NEMO/OPA_SRC/LBC/lib_mpp.F90
r5836 r6006 72 72 PUBLIC mpp_lnk_3d, mpp_lnk_3d_gather, mpp_lnk_2d, mpp_lnk_2d_e 73 73 PUBLIC mpp_lnk_2d_9 74 PUBLIC mpp_lnk_sum_3d, mpp_lnk_sum_2d 74 75 PUBLIC mppscatter, mppgather 75 76 PUBLIC mpp_ini_ice, mpp_ini_znl … … 1402 1403 END SUBROUTINE mpp_lnk_2d_e 1403 1404 1405 SUBROUTINE mpp_lnk_sum_3d( ptab, cd_type, psgn, cd_mpp, pval ) 1406 !!---------------------------------------------------------------------- 1407 !! *** routine mpp_lnk_sum_3d *** 1408 !! 1409 !! ** Purpose : Message passing manadgement (sum the overlap region) 1410 !! 1411 !! ** Method : Use mppsend and mpprecv function for passing mask 1412 !! between processors following neighboring subdomains. 1413 !! domain parameters 1414 !! nlci : first dimension of the local subdomain 1415 !! nlcj : second dimension of the local subdomain 1416 !! nbondi : mark for "east-west local boundary" 1417 !! nbondj : mark for "north-south local boundary" 1418 !! noea : number for local neighboring processors 1419 !! nowe : number for local neighboring processors 1420 !! noso : number for local neighboring processors 1421 !! nono : number for local neighboring processors 1422 !! 1423 !! ** Action : ptab with update value at its periphery 1424 !! 1425 !!---------------------------------------------------------------------- 1426 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) :: ptab ! 3D array on which the boundary condition is applied 1427 CHARACTER(len=1) , INTENT(in ) :: cd_type ! define the nature of ptab array grid-points 1428 ! ! = T , U , V , F , W points 1429 REAL(wp) , INTENT(in ) :: psgn ! =-1 the sign change across the north fold boundary 1430 ! ! = 1. , the sign is kept 1431 CHARACTER(len=3), OPTIONAL , INTENT(in ) :: cd_mpp ! fill the overlap area only 1432 REAL(wp) , OPTIONAL , INTENT(in ) :: pval ! background value (used at closed boundaries) 1433 !! 1434 INTEGER :: ji, jj, jk, jl ! dummy loop indices 1435 INTEGER :: imigr, iihom, ijhom ! temporary integers 1436 INTEGER :: ml_req1, ml_req2, ml_err ! for key_mpi_isend 1437 REAL(wp) :: zland 1438 INTEGER, DIMENSION(MPI_STATUS_SIZE) :: ml_stat ! for key_mpi_isend 1439 ! 1440 REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE :: zt3ns, zt3sn ! 3d for north-south & south-north 1441 REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE :: zt3ew, zt3we ! 3d for east-west & west-east 1442 1443 !!---------------------------------------------------------------------- 1444 1445 ALLOCATE( zt3ns(jpi,jprecj,jpk,2), zt3sn(jpi,jprecj,jpk,2), & 1446 & zt3ew(jpj,jpreci,jpk,2), zt3we(jpj,jpreci,jpk,2) ) 1447 1448 ! 1449 IF( PRESENT( pval ) ) THEN ; zland = pval ! set land value 1450 ELSE ; zland = 0.e0 ! zero by default 1451 ENDIF 1452 1453 ! 1. standard boundary treatment 1454 ! ------------------------------ 1455 ! 2. East and west directions exchange 1456 ! ------------------------------------ 1457 ! we play with the neigbours AND the row number because of the periodicity 1458 ! 1459 SELECT CASE ( nbondi ) ! Read lateral conditions 1460 CASE ( -1, 0, 1 ) ! all exept 2 (i.e. close case) 1461 iihom = nlci-jpreci 1462 DO jl = 1, jpreci 1463 zt3ew(:,jl,:,1) = ptab(jl ,:,:) ; ptab(jl ,:,:) = 0.0_wp 1464 zt3we(:,jl,:,1) = ptab(iihom+jl,:,:) ; ptab(iihom+jl,:,:) = 0.0_wp 1465 END DO 1466 END SELECT 1467 ! 1468 ! ! Migrations 1469 imigr = jpreci * jpj * jpk 1470 ! 1471 SELECT CASE ( nbondi ) 1472 CASE ( -1 ) 1473 CALL mppsend( 2, zt3we(1,1,1,1), imigr, noea, ml_req1 ) 1474 CALL mpprecv( 1, zt3ew(1,1,1,2), imigr, noea ) 1475 IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err) 1476 CASE ( 0 ) 1477 CALL mppsend( 1, zt3ew(1,1,1,1), imigr, nowe, ml_req1 ) 1478 CALL mppsend( 2, zt3we(1,1,1,1), imigr, noea, ml_req2 ) 1479 CALL mpprecv( 1, zt3ew(1,1,1,2), imigr, noea ) 1480 CALL mpprecv( 2, zt3we(1,1,1,2), imigr, nowe ) 1481 IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err) 1482 IF(l_isend) CALL mpi_wait(ml_req2, ml_stat, ml_err) 1483 CASE ( 1 ) 1484 CALL mppsend( 1, zt3ew(1,1,1,1), imigr, nowe, ml_req1 ) 1485 CALL mpprecv( 2, zt3we(1,1,1,2), imigr, nowe ) 1486 IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err) 1487 END SELECT 1488 ! 1489 ! ! Write lateral conditions 1490 iihom = nlci-nreci 1491 ! 1492 SELECT CASE ( nbondi ) 1493 CASE ( -1 ) 1494 DO jl = 1, jpreci 1495 ptab(iihom+jl,:,:) = ptab(iihom+jl,:,:) + zt3ew(:,jl,:,2) 1496 END DO 1497 CASE ( 0 ) 1498 DO jl = 1, jpreci 1499 ptab(jpreci+jl,:,:) = ptab(jpreci+jl,:,:) + zt3we(:,jl,:,2) 1500 ptab(iihom +jl,:,:) = ptab(iihom +jl,:,:) + zt3ew(:,jl,:,2) 1501 END DO 1502 CASE ( 1 ) 1503 DO jl = 1, jpreci 1504 ptab(jpreci+jl,:,:) = ptab(jpreci+jl,:,:) + zt3we(:,jl,:,2) 1505 END DO 1506 END SELECT 1507 1508 1509 ! 3. North and south directions 1510 ! ----------------------------- 1511 ! always closed : we play only with the neigbours 1512 ! 1513 IF( nbondj /= 2 ) THEN ! Read lateral conditions 1514 ijhom = nlcj-jprecj 1515 DO jl = 1, jprecj 1516 zt3sn(:,jl,:,1) = ptab(:,ijhom+jl,:) ; ptab(:,ijhom+jl,:) = 0.0_wp 1517 zt3ns(:,jl,:,1) = ptab(:,jl ,:) ; ptab(:,jl ,:) = 0.0_wp 1518 END DO 1519 ENDIF 1520 ! 1521 ! ! Migrations 1522 imigr = jprecj * jpi * jpk 1523 ! 1524 SELECT CASE ( nbondj ) 1525 CASE ( -1 ) 1526 CALL mppsend( 4, zt3sn(1,1,1,1), imigr, nono, ml_req1 ) 1527 CALL mpprecv( 3, zt3ns(1,1,1,2), imigr, nono ) 1528 IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err) 1529 CASE ( 0 ) 1530 CALL mppsend( 3, zt3ns(1,1,1,1), imigr, noso, ml_req1 ) 1531 CALL mppsend( 4, zt3sn(1,1,1,1), imigr, nono, ml_req2 ) 1532 CALL mpprecv( 3, zt3ns(1,1,1,2), imigr, nono ) 1533 CALL mpprecv( 4, zt3sn(1,1,1,2), imigr, noso ) 1534 IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err) 1535 IF(l_isend) CALL mpi_wait(ml_req2, ml_stat, ml_err) 1536 CASE ( 1 ) 1537 CALL mppsend( 3, zt3ns(1,1,1,1), imigr, noso, ml_req1 ) 1538 CALL mpprecv( 4, zt3sn(1,1,1,2), imigr, noso ) 1539 IF(l_isend) CALL mpi_wait(ml_req1, ml_stat, ml_err) 1540 END SELECT 1541 ! 1542 ! ! Write lateral conditions 1543 ijhom = nlcj-nrecj 1544 ! 1545 SELECT CASE ( nbondj ) 1546 CASE ( -1 ) 1547 DO jl = 1, jprecj 1548 ptab(:,ijhom+jl,:) = ptab(:,ijhom+jl,:) + zt3ns(:,jl,:,2) 1549 END DO 1550 CASE ( 0 ) 1551 DO jl = 1, jprecj 1552 ptab(:,jprecj+jl,:) = ptab(:,jprecj+jl,:) + zt3sn(:,jl,:,2) 1553 ptab(:,ijhom +jl,:) = ptab(:,ijhom +jl,:) + zt3ns(:,jl,:,2) 1554 END DO 1555 CASE ( 1 ) 1556 DO jl = 1, jprecj 1557 ptab(:,jprecj+jl,:) = ptab(:,jprecj+jl,:) + zt3sn(:,jl ,:,2) 1558 END DO 1559 END SELECT 1560 1561 1562 ! 4. north fold treatment 1563 ! ----------------------- 1564 ! 1565 IF( npolj /= 0 .AND. .NOT. PRESENT(cd_mpp) ) THEN 1566 ! 1567 SELECT CASE ( jpni ) 1568 CASE ( 1 ) ; CALL lbc_nfd ( ptab, cd_type, psgn ) ! only 1 northern proc, no mpp 1569 CASE DEFAULT ; CALL mpp_lbc_north( ptab, cd_type, psgn ) ! for all northern procs. 1570 END SELECT 1571 ! 1572 ENDIF 1573 ! 1574 DEALLOCATE( zt3ns, zt3sn, zt3ew, zt3we ) 1575 ! 1576 END SUBROUTINE mpp_lnk_sum_3d 1577 1578 SUBROUTINE mpp_lnk_sum_2d( pt2d, cd_type, psgn, cd_mpp, pval ) 1579 !!---------------------------------------------------------------------- 1580 !! *** routine mpp_lnk_sum_2d *** 1581 !! 1582 !! ** Purpose : Message passing manadgement for 2d array (sum the overlap region) 1583 !! 1584 !! ** Method : Use mppsend and mpprecv function for passing mask 1585 !! between processors following neighboring subdomains. 1586 !! domain parameters 1587 !! nlci : first dimension of the local subdomain 1588 !! nlcj : second dimension of the local subdomain 1589 !! nbondi : mark for "east-west local boundary" 1590 !! nbondj : mark for "north-south local boundary" 1591 !! noea : number for local neighboring processors 1592 !! nowe : number for local neighboring processors 1593 !! noso : number for local neighboring processors 1594 !! nono : number for local neighboring processors 1595 !! 1596 !!---------------------------------------------------------------------- 1597 REAL(wp), DIMENSION(jpi,jpj), INTENT(inout) :: pt2d ! 2D array on which the boundary condition is applied 1598 CHARACTER(len=1) , INTENT(in ) :: cd_type ! define the nature of ptab array grid-points 1599 ! ! = T , U , V , F , W and I points 1600 REAL(wp) , INTENT(in ) :: psgn ! =-1 the sign change across the north fold boundary 1601 ! ! = 1. , the sign is kept 1602 CHARACTER(len=3), OPTIONAL , INTENT(in ) :: cd_mpp ! fill the overlap area only 1603 REAL(wp) , OPTIONAL , INTENT(in ) :: pval ! background value (used at closed boundaries) 1604 !! 1605 INTEGER :: ji, jj, jl ! dummy loop indices 1606 INTEGER :: imigr, iihom, ijhom ! temporary integers 1607 INTEGER :: ml_req1, ml_req2, ml_err ! for key_mpi_isend 1608 REAL(wp) :: zland 1609 INTEGER, DIMENSION(MPI_STATUS_SIZE) :: ml_stat ! for key_mpi_isend 1610 ! 1611 REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: zt2ns, zt2sn ! 2d for north-south & south-north 1612 REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: zt2ew, zt2we ! 2d for east-west & west-east 1613 1614 !!---------------------------------------------------------------------- 1615 1616 ALLOCATE( zt2ns(jpi,jprecj,2), zt2sn(jpi,jprecj,2), & 1617 & zt2ew(jpj,jpreci,2), zt2we(jpj,jpreci,2) ) 1618 1619 ! 1620 IF( PRESENT( pval ) ) THEN ; zland = pval ! set land value 1621 ELSE ; zland = 0.e0 ! zero by default 1622 ENDIF 1623 1624 ! 1. standard boundary treatment 1625 ! ------------------------------ 1626 ! 2. East and west directions exchange 1627 ! ------------------------------------ 1628 ! we play with the neigbours AND the row number because of the periodicity 1629 ! 1630 SELECT CASE ( nbondi ) ! Read lateral conditions 1631 CASE ( -1, 0, 1 ) ! all exept 2 (i.e. close case) 1632 iihom = nlci - jpreci 1633 DO jl = 1, jpreci 1634 zt2ew(:,jl,1) = pt2d(jl ,:) ; pt2d(jl ,:) = 0.0_wp 1635 zt2we(:,jl,1) = pt2d(iihom +jl,:) ; pt2d(iihom +jl,:) = 0.0_wp 1636 END DO 1637 END SELECT 1638 ! 1639 ! ! Migrations 1640 imigr = jpreci * jpj 1641 ! 1642 SELECT CASE ( nbondi ) 1643 CASE ( -1 ) 1644 CALL mppsend( 2, zt2we(1,1,1), imigr, noea, ml_req1 ) 1645 CALL mpprecv( 1, zt2ew(1,1,2), imigr, noea ) 1646 IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 1647 CASE ( 0 ) 1648 CALL mppsend( 1, zt2ew(1,1,1), imigr, nowe, ml_req1 ) 1649 CALL mppsend( 2, zt2we(1,1,1), imigr, noea, ml_req2 ) 1650 CALL mpprecv( 1, zt2ew(1,1,2), imigr, noea ) 1651 CALL mpprecv( 2, zt2we(1,1,2), imigr, nowe ) 1652 IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 1653 IF(l_isend) CALL mpi_wait(ml_req2,ml_stat,ml_err) 1654 CASE ( 1 ) 1655 CALL mppsend( 1, zt2ew(1,1,1), imigr, nowe, ml_req1 ) 1656 CALL mpprecv( 2, zt2we(1,1,2), imigr, nowe ) 1657 IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 1658 END SELECT 1659 ! 1660 ! ! Write lateral conditions 1661 iihom = nlci-nreci 1662 ! 1663 SELECT CASE ( nbondi ) 1664 CASE ( -1 ) 1665 DO jl = 1, jpreci 1666 pt2d(iihom+jl,:) = pt2d(iihom+jl,:) + zt2ew(:,jl,2) 1667 END DO 1668 CASE ( 0 ) 1669 DO jl = 1, jpreci 1670 pt2d(jpreci+jl,:) = pt2d(jpreci+jl,:) + zt2we(:,jl,2) 1671 pt2d(iihom +jl,:) = pt2d(iihom +jl,:) + zt2ew(:,jl,2) 1672 END DO 1673 CASE ( 1 ) 1674 DO jl = 1, jpreci 1675 pt2d(jpreci+jl,:) = pt2d(jpreci+jl,:) + zt2we(:,jl,2) 1676 END DO 1677 END SELECT 1678 1679 1680 ! 3. North and south directions 1681 ! ----------------------------- 1682 ! always closed : we play only with the neigbours 1683 ! 1684 IF( nbondj /= 2 ) THEN ! Read lateral conditions 1685 ijhom = nlcj - jprecj 1686 DO jl = 1, jprecj 1687 zt2sn(:,jl,1) = pt2d(:,ijhom +jl) ; pt2d(:,ijhom +jl) = 0.0_wp 1688 zt2ns(:,jl,1) = pt2d(:,jl ) ; pt2d(:,jl ) = 0.0_wp 1689 END DO 1690 ENDIF 1691 ! 1692 ! ! Migrations 1693 imigr = jprecj * jpi 1694 ! 1695 SELECT CASE ( nbondj ) 1696 CASE ( -1 ) 1697 CALL mppsend( 4, zt2sn(1,1,1), imigr, nono, ml_req1 ) 1698 CALL mpprecv( 3, zt2ns(1,1,2), imigr, nono ) 1699 IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 1700 CASE ( 0 ) 1701 CALL mppsend( 3, zt2ns(1,1,1), imigr, noso, ml_req1 ) 1702 CALL mppsend( 4, zt2sn(1,1,1), imigr, nono, ml_req2 ) 1703 CALL mpprecv( 3, zt2ns(1,1,2), imigr, nono ) 1704 CALL mpprecv( 4, zt2sn(1,1,2), imigr, noso ) 1705 IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 1706 IF(l_isend) CALL mpi_wait(ml_req2,ml_stat,ml_err) 1707 CASE ( 1 ) 1708 CALL mppsend( 3, zt2ns(1,1,1), imigr, noso, ml_req1 ) 1709 CALL mpprecv( 4, zt2sn(1,1,2), imigr, noso ) 1710 IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 1711 END SELECT 1712 ! 1713 ! ! Write lateral conditions 1714 ijhom = nlcj-nrecj 1715 ! 1716 SELECT CASE ( nbondj ) 1717 CASE ( -1 ) 1718 DO jl = 1, jprecj 1719 pt2d(:,ijhom+jl) = pt2d(:,ijhom+jl) + zt2ns(:,jl,2) 1720 END DO 1721 CASE ( 0 ) 1722 DO jl = 1, jprecj 1723 pt2d(:,jprecj+jl) = pt2d(:,jprecj+jl) + zt2sn(:,jl,2) 1724 pt2d(:,ijhom +jl) = pt2d(:,ijhom +jl) + zt2ns(:,jl,2) 1725 END DO 1726 CASE ( 1 ) 1727 DO jl = 1, jprecj 1728 pt2d(:,jprecj+jl) = pt2d(:,jprecj+jl) + zt2sn(:,jl,2) 1729 END DO 1730 END SELECT 1731 1732 1733 ! 4. north fold treatment 1734 ! ----------------------- 1735 ! 1736 IF( npolj /= 0 .AND. .NOT. PRESENT(cd_mpp) ) THEN 1737 ! 1738 SELECT CASE ( jpni ) 1739 CASE ( 1 ) ; CALL lbc_nfd ( pt2d, cd_type, psgn ) ! only 1 northern proc, no mpp 1740 CASE DEFAULT ; CALL mpp_lbc_north( pt2d, cd_type, psgn ) ! for all northern procs. 1741 END SELECT 1742 ! 1743 ENDIF 1744 ! 1745 DEALLOCATE( zt2ns, zt2sn, zt2ew, zt2we ) 1746 ! 1747 END SUBROUTINE mpp_lnk_sum_2d 1404 1748 1405 1749 SUBROUTINE mppsend( ktyp, pmess, kbytes, kdest, md_req ) -
branches/2015/dev_MetOffice_merge_2015/NEMOGCM/NEMO/OPA_SRC/LBC/mppini_2.h90
r5130 r6006 136 136 137 137 imask(:,:)=1 138 WHERE ( zdta(:,:) - zdtaisf(:,:) <= 0.) imask = 0138 WHERE ( zdta(:,:) - zdtaisf(:,:) <= rn_isfhmin ) imask = 0 139 139 140 140 ! 1. Dimension arrays for subdomains
Note: See TracChangeset
for help on using the changeset viewer.