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 1218 for trunk/NEMO/OPA_SRC/SBC/cpl_oasis3.F90 – NEMO

Ignore:
Timestamp:
2008-10-28T10:12:16+01:00 (16 years ago)
Author:
smasson
Message:

first implementation of the new coupling interface in the trunk, see ticket:155

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMO/OPA_SRC/SBC/cpl_oasis3.F90

    r1146 r1218  
    2121   !!   cpl_prism_init     : initialization of coupled mode communication 
    2222   !!   cpl_prism_define   : definition of grid and fields 
    23    !!   cpl_prism_send     : send out fields in coupled mode 
    24    !!   cpl_prism_recv     : receive fields in coupled mode 
     23   !!   cpl_prism_snd     : snd out fields in coupled mode 
     24   !!   cpl_prism_rcv     : receive fields in coupled mode 
    2525   !!   cpl_prism_finalize : finalize the coupled mode communication 
    2626   !!---------------------------------------------------------------------- 
    27    !! * Modules used 
    28 !##################### WARNING coupled mode ############################### 
    29 !##################### WARNING coupled mode ############################### 
    30 !   Following lines must be enabled if coupling with OASIS 
     27   USE mod_prism_proto              ! OASIS3 prism module 
     28   USE mod_prism_def_partition_proto! OASIS3 prism module for partitioning 
     29   USE mod_prism_grids_writing      ! OASIS3 prism module for writing grid files 
     30   USE mod_prism_put_proto          ! OASIS3 prism module for snding 
     31   USE mod_prism_get_proto          ! OASIS3 prism module for receiving 
     32   USE mod_prism_grids_writing      ! OASIS3 prism module for writing grids 
     33   USE par_oce                      ! ocean parameters 
     34   USE dom_oce                      ! ocean space and time domain 
     35   USE in_out_manager               ! I/O manager 
     36   USE lib_mpp 
     37   USE lbclnk          ! ocean lateral boundary conditions (or mpp link) 
     38   IMPLICIT NONE 
     39   PRIVATE 
    3140! 
    32 !   USE mod_prism_proto              ! OASIS3 prism module 
    33 !   USE mod_prism_def_partition_proto! OASIS3 prism module for partitioning 
    34 !   USE mod_prism_grids_writing      ! OASIS3 prism module for writing grid files 
    35 !   USE mod_prism_put_proto          ! OASIS3 prism module for sending 
    36 !   USE mod_prism_get_proto          ! OASIS3 prism module for receiving 
    37 !   USE mod_prism_grids_writing      ! OASIS3 prism module for writing grids 
    38 !##################### WARNING coupled mode ############################### 
    39 !##################### WARNING coupled mode ############################### 
    40 #if defined key_mpp_mpi 
    41    USE lib_mpp, only : mppsize, mpprank ! message passing 
    42    USE lib_mpp, only : mppsend          ! message passing 
    43    USE lib_mpp, only : mpprecv          ! message passing 
    44 #endif 
    45    USE daymod                       ! date and time info 
    46    USE dom_oce                      ! ocean space and time domain 
    47    USE sbc_ice                      ! surface boundary condition: ice 
    48    USE in_out_manager               ! I/O manager 
    49    USE par_oce                      ! 
    50    USE phycst, only : rt0           ! freezing point of sea water 
    51  
    52    USE oce, only: tn, un, vn 
    53 #if defined key_lim2 
    54    USE ice_2, only: frld, hicif, hsnif 
    55 #endif 
    56  
    57    IMPLICIT NONE 
    58 ! 
    59 ! Exchange parameters for coupling ORCA-LIM with ECHAM5 
    60 ! 
    61 #if defined key_cpl_ocevel 
    62    INTEGER, PARAMETER         :: nsend =  6 
    63 #else 
    64    INTEGER, PARAMETER         :: nsend =  4 
    65 #endif 
    66  
    67 #if defined key_cpl_discharge 
    68    INTEGER, PARAMETER         :: nrecv = 20 
    69 #else 
    70    INTEGER, PARAMETER         :: nrecv = 17 
    71 #endif 
    72  
    73    INTEGER, DIMENSION(nsend)  :: send_id 
    74    INTEGER, DIMENSION(nrecv)  :: recv_id 
    75  
    76    CHARACTER(len=32)          :: cpl_send (nsend) 
    77    CHARACTER(len=32)          :: cpl_recv (nrecv) 
    78  
    79    PRIVATE 
    80  
    81    INTEGER                    :: localRank      ! local MPI rank 
    82    INTEGER                    :: comp_id        ! id returned by prism_init_comp 
    83  
    84    INTEGER                    :: range(5) 
    85  
    86    INTEGER, PARAMETER         :: localRoot  = 0 
    87    INTEGER                    :: localSize      ! local MPI size 
    88    INTEGER                    :: localComm      ! local MPI size 
    89    LOGICAL                    :: commRank       ! true for ranks doing OASIS communication 
    90  
    91    LOGICAL, SAVE              :: prism_was_initialized 
    92    LOGICAL, SAVE              :: prism_was_terminated 
    93    INTEGER, SAVE              :: write_grid 
    94  
    95    INTEGER                    :: ierror         ! return error code 
     41   INTEGER, PUBLIC            :: nlocalComm 
     42   LOGICAL, PUBLIC, PARAMETER :: lk_cpl = .TRUE.   ! coupled flag 
     43   INTEGER                    :: ncomp_id          ! id returned by prism_init_comp 
     44   INTEGER                    :: nerror            ! return error code 
     45 
     46   INTEGER, PUBLIC :: nrcv, nsnd    ! Number of received and sent coupling fields 
     47 
     48   INTEGER, PARAMETER :: nmaxfld=40    ! Maximum number of coupling fields 
     49    
     50   TYPE, PUBLIC ::   FLD_CPL            ! Type for coupling field information 
     51      LOGICAL            ::   laction   ! To be coupled or not 
     52      CHARACTER(len = 8) ::   clname    ! Name of the coupling field    
     53      CHARACTER(len = 1) ::   clgrid    ! Grid type   
     54      REAL(wp)           ::   nsgn      ! Control of the sign change 
     55      INTEGER            ::   nid       ! Id of the field 
     56   END TYPE FLD_CPL 
     57 
     58   TYPE(FLD_CPL), DIMENSION(nmaxfld), PUBLIC :: srcv, ssnd   ! Coupling fields 
    9659 
    9760   REAL(wp), DIMENSION(:,:), ALLOCATABLE :: exfld  ! Temporary buffer for receiving 
    98  
    99 #ifdef key_cpl_rootexchg 
    100    LOGICAL                               :: rootexchg =.true.     ! logical switch  
    101 #else 
    102    LOGICAL                               :: rootexchg =.false.    ! logical switch 
    103 #endif 
    104  
    105    REAL(wp), DIMENSION(:),   ALLOCATABLE :: buffer ! Temporary buffer for exchange 
    106    INTEGER, DIMENSION(:,:),  ALLOCATABLE :: ranges ! Temporary buffer for exchange 
    10761 
    10862   !! Routine accessibility 
    10963   PUBLIC cpl_prism_init 
    11064   PUBLIC cpl_prism_define 
    111    PUBLIC cpl_prism_send 
    112    PUBLIC cpl_prism_recv 
     65   PUBLIC cpl_prism_snd 
     66   PUBLIC cpl_prism_rcv 
    11367   PUBLIC cpl_prism_finalize 
    11468 
    115    PUBLIC send_id, recv_id 
    116  
    11769   !!---------------------------------------------------------------------- 
    11870   !!   OPA 9.0 , LOCEAN-IPSL (2006)  
    119    !! $Id$ 
     71   !! $Header$  
    12072   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt) 
    12173   !!---------------------------------------------------------------------- 
     
    12375CONTAINS 
    12476 
    125    SUBROUTINE cpl_prism_init( localCommunicator ) 
    126  
    127       IMPLICIT NONE 
     77   SUBROUTINE cpl_prism_init 
    12878 
    12979      !!------------------------------------------------------------------- 
     
    13585      !! ** Method  :   OASIS3 MPI communication  
    13686      !!-------------------------------------------------------------------- 
    137       !! * Arguments 
    138       !! 
    139       INTEGER, INTENT(OUT)       :: localCommunicator 
    140       !! 
    141       !! * Local declarations 
    142       !! 
    143       CHARACTER(len=4)           :: comp_name      ! name of this PRISM component 
    144       !! 
    145       !!-------------------------------------------------------------------- 
    146       !! 
    147       IF(lwp) WRITE(numout,*) 
     87      !! 
     88 
    14889      IF(lwp) WRITE(numout,*) 'cpl_prism_init : initialization in coupled ocean/atmosphere case' 
    14990      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~~~' 
    15091      IF(lwp) WRITE(numout,*) 
    151       
    152 #if defined key_flx_bulk_monthly || defined key_flx_bulk_daily || defined key_flx_forced_daily 
    153       IF(lwp)WRITE(numout,cform_err) 
    154       IF(lwp)WRITE(numout,*) ' key_coupled and key_flx_bulk_* key_flx_forced_daily are incompatible' 
    155       nstop = nstop + 1 
    156 #endif 
    157  
    158       comp_name = 'opa9' 
    159  
    16092      !------------------------------------------------------------------ 
    16193      ! 1st Initialize the PRISM system for the application 
    16294      !------------------------------------------------------------------ 
    163  
    164       CALL prism_init_comp_proto ( comp_id, comp_name, ierror ) 
    165       IF ( ierror /= PRISM_Ok ) & 
    166          CALL prism_abort_proto (comp_id, 'cpl_prism_init', 'Failure in prism_init_comp_proto') 
    167       prism_was_initialized = .true. 
     95      CALL prism_init_comp_proto ( ncomp_id, 'oceanx', nerror ) 
     96      IF ( nerror /= PRISM_Ok ) & 
     97         CALL prism_abort_proto (ncomp_id, 'cpl_prism_init', 'Failure in prism_init_comp_proto') 
    16898 
    16999      !------------------------------------------------------------------ 
     
    171101      !------------------------------------------------------------------ 
    172102 
    173       CALL prism_get_localcomm_proto ( localComm, ierror ) 
    174       IF ( ierror /= PRISM_Ok ) & 
    175          CALL prism_abort_proto (comp_id, 'cpl_prism_init','Failure in prism_get_localcomm_proto' ) 
    176  
    177       localCommunicator = localComm 
     103      CALL prism_get_localcomm_proto ( nlocalComm, nerror ) 
     104      IF ( nerror /= PRISM_Ok ) & 
     105         CALL prism_abort_proto (ncomp_id, 'cpl_prism_init','Failure in prism_get_localcomm_proto' ) 
    178106 
    179107   END SUBROUTINE cpl_prism_init 
     
    181109 
    182110   SUBROUTINE cpl_prism_define () 
    183  
    184       IMPLICIT NONE 
    185111 
    186112      !!------------------------------------------------------------------- 
     
    196122      !! * Local declarations 
    197123      !! 
    198       INTEGER                    :: grid_id(2)     ! id returned by prism_def_grid 
    199       INTEGER                    :: part_id 
    200  
     124      INTEGER                    :: id_part 
    201125      INTEGER                    :: paral(5)       ! OASIS3 box partition 
    202  
    203       INTEGER                    :: shape(2,3)     ! shape of arrays passed to PSMILe 
    204       INTEGER                    :: nodim(2) 
    205       INTEGER                    :: data_type      ! data type of transients 
    206  
    207       INTEGER                    :: ji, jj         ! local loop indicees 
    208       INTEGER                    :: nx, ny, nc     ! local variables 
    209       INTEGER                    :: im1, ip1 
    210       INTEGER                    :: jm1, jp1 
    211       INTEGER                    :: i_grid         ! loop index 
    212       INTEGER                    :: info 
    213       INTEGER                    :: maxlen 
    214       INTEGER                    :: mask(jpi,jpj) 
    215       REAL(kind=wp)              :: area(jpi,jpj) 
    216  
    217       CHARACTER(len=4)           :: point_name     ! name of the grid points 
    218  
    219       REAL(kind=wp)              :: rclam(jpi,jpj,4) 
    220       REAL(kind=wp)              :: rcphi(jpi,jpj,4) 
    221  
    222       REAL(kind=wp)              :: glam_b(jpi,jpj) ! buffer for orca2 grid correction 
    223       REAL(kind=wp)              :: gphi_b(jpi,jpj) ! buffer for orca2 grid correction 
    224       !! 
    225       !!-------------------------------------------------------------------- 
    226       
     126      INTEGER                    :: ishape(2,2)    ! shape of arrays passed to PSMILe 
     127      INTEGER                    :: ji             ! local loop indicees 
     128      !! 
     129      !!-------------------------------------------------------------------- 
     130 
    227131      IF(lwp) WRITE(numout,*) 
    228132      IF(lwp) WRITE(numout,*) 'cpl_prism_define : initialization in coupled ocean/atmosphere case' 
    229133      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~~~~' 
    230134      IF(lwp) WRITE(numout,*) 
    231       
    232 #if defined key_flx_bulk_monthly || defined key_flx_bulk_daily || defined key_flx_forced_daily 
    233       IF(lwp)WRITE(numout,cform_err) 
    234       IF(lwp)WRITE(numout,*) ' key_coupled and key_flx_bulk_... are incompatible' 
    235       nstop = nstop + 1 
    236 #endif 
    237  
    238       ! ----------------------------------------------------------------- 
    239       ! ... Some initialisation 
    240       ! ----------------------------------------------------------------- 
    241  
    242       send_id = 0 
    243       recv_id = 0 
    244  
    245 #if defined key_mpp_mpi 
    246  
    247       ! ----------------------------------------------------------------- 
    248       ! ... Some MPI stuff relevant for optional exchange via root only 
    249       ! ----------------------------------------------------------------- 
    250  
    251       commRank = .false. 
    252  
    253       localRank = mpprank ! from lib_mpp 
    254       localSize = mppsize ! from lib_mpp 
    255  
    256       IF ( rootexchg ) THEN 
    257          IF ( localRank == localRoot ) commRank = .true. 
    258       ELSE 
    259          commRank = .true. 
     135 
     136      ! 
     137      ! ... Define the shape for the area that excludes the halo 
     138      !     For serial configuration (key_mpp_mpi not being active) 
     139      !     nl* is set to the global values 1 and jp*glo. 
     140      ! 
     141      ishape(:,1) = (/ 1, nlei-nldi+1 /) 
     142      ishape(:,2) = (/ 1, nlej-nldj+1 /) 
     143      ! 
     144      ! ... Allocate memory for data exchange 
     145      ! 
     146      ALLOCATE(exfld(nlei-nldi+1, nlej-nldj+1), stat = nerror) 
     147      IF (nerror > 0) THEN 
     148         CALL prism_abort_proto ( ncomp_id, 'cpl_prism_define', 'Failure in allocating exfld') 
     149         RETURN 
    260150      ENDIF 
    261  
    262       IF ( rootexchg .and. localRank == localRoot ) THEN 
    263          ALLOCATE(ranges(5,0:localSize-1), stat = ierror) 
    264          IF (ierror > 0) THEN 
    265             CALL prism_abort_proto ( comp_id, 'cpl_prism_define', 'Failure in allocating Integer') 
    266             RETURN 
    267          ENDIF 
    268       ENDIF 
    269  
    270 #else 
    271       ! 
    272       ! For non-parallel configurations the one and only process ("localRoot") 
    273       ! takes part in the communication 
    274       !  
    275       localRank = localRoot 
    276       commRank = .true. 
    277  
    278 #endif 
    279  
    280       ! ----------------------------------------------------------------- 
    281       ! ... If necessary the root process writes the global grid info 
    282       ! ----------------------------------------------------------------- 
    283  
    284       IF ( localRank == localRoot ) THEN 
    285  
    286          WRITE(numout,*)'Opening file SSTOCEAN, unit= 199' 
    287  
    288          OPEN (199,STATUS='NEW',FILE="sstocean",FORM='UNFORMATTED',err=310) 
    289  
    290          ! In case the sstocean of OASIS3 from a previous run exists 
    291          ! the programs jumps to the end of the if-block 
    292 !     
    293 !*    2.0    Write exchange fields to OASIS data file. 
    294 !            ----------------------------------------- 
    295  
    296          WHERE (tmask(:,:,1) > 0.5 ) 
    297             mask(:,:) = 0 
    298          ELSE WHERE 
    299             mask(:,:) = 1 
    300          END WHERE 
    301  
    302 ! Initialise ice mask at the very first start only 
    303          frld = 1. 
    304  
    305          WRITE(199) 'SSTOCEAN' 
    306          WRITE(199) (tn(:,:,1)*mask(:,:))+rt0 
    307  
    308          WRITE(199) 'SICOCEAN' 
    309          WRITE(199) (1.-frld(:,:))*mask(:,:) 
    310  
    311 #if defined key_cpl_albedo 
    312 # if defined key_lim3 
    313          Must be adapted for LIM3 
    314 # endif 
    315          tn_ice  = 271.285 
    316     alb_ice =   0.75 
    317  
    318          WRITE(199) 'STIOCEAN' 
    319          WRITE(199) tn_ice(:,:) 
    320  
    321          WRITE(199) 'SAIOCEAN' 
    322          WRITE(199) alb_ice(:,:) 
    323 #else 
    324          hicit = 0. 
    325          hsnit = 0. 
    326          WRITE(199) 'SITOCEAN' 
    327          WRITE(199) hicif(:,:)*mask(:,:) 
    328  
    329          WRITE(199) 'SNTOCEAN' 
    330          WRITE(199) hsnif(:,:)*mask(:,:) 
    331 #endif 
    332  
    333 #if defined key_cpl_ocevel 
    334          un(:,:,1) = 0. 
    335          vn(:,:,1) = 0. 
    336  
    337          WHERE (umask(:,:,1) > 0.5 ) 
    338             mask(:,:) = 0 
    339          ELSE WHERE 
    340             mask(:,:) = 1 
    341          END WHERE 
    342  
    343          WRITE(199) 'SUNOCEAN' 
    344          WRITE(199) un(:,:,1)*mask(:,:) 
    345  
    346          WHERE (vmask(:,:,1) > 0.5 ) 
    347             mask(:,:) = 0 
    348          ELSE WHERE 
    349             mask(:,:) = 1 
    350          END WHERE 
    351  
    352          WRITE(199) 'SVNOCEAN' 
    353          WRITE(199) vn(:,:,1)*mask(:,:) 
    354 #endif 
    355  
    356          WRITE(numout,*) 
    357          WRITE(numout,*)' sstocean written' 
    358          WRITE(numout,*)' ***************' 
    359  
    360          CLOSE(199) 
    361  
    362  310     CONTINUE 
    363  
    364          CALL prism_start_grids_writing ( write_grid ) 
    365  
    366       ENDIF  ! localRank == localRoot 
    367  
    368       IF ( localRank == localRoot .and. write_grid == 1 ) THEN 
    369  
    370          !------------------------------------------------------------------ 
    371          ! 1st write global grid information (ORCA tripolar) characteristics 
    372          !     for surface coupling into a OASIS3 specific grid file. For 
    373          !     surface coupling it is sufficient to specify only one vertical 
    374          !     z-level. 
    375          !------------------------------------------------------------------ 
    376          ! 
    377          ! ... Treat corners in the horizontal plane 
    378          ! 
    379          nx = jpi 
    380          ny = jpj 
    381          nc = 4 
    382  
    383          DO i_grid = 1, 3 
    384  
    385             IF ( i_grid == 1 ) THEN 
    386  
    387                ! -------------------------------------------------------- 
    388                ! ... Write the grid info for T points 
    389                ! -------------------------------------------------------- 
    390  
    391                point_name = 'opat' 
    392  
    393                glam_b = glamt 
    394                gphi_b = gphit 
    395  
    396                DO ji = 1, jpi 
    397                   DO jj = 1, jpj 
    398  
    399                      im1 = ji-1 
    400                      jm1 = jj-1 
    401                      IF (ji == 1) im1 = jpi-2 
    402                      IF (jj == 1) jm1 = jj 
    403  
    404                      rclam(ji,jj,1) = glamf(ji,jj) 
    405                      rclam(ji,jj,2) = glamf(im1,jj) 
    406                      rclam(ji,jj,3) = glamf(im1,jm1) 
    407                      rclam(ji,jj,4) = glamf(ji,jm1) 
    408  
    409                      rcphi(ji,jj,1) = gphif(ji,jj) 
    410                      rcphi(ji,jj,2) = gphif(im1,jj) 
    411                      rcphi(ji,jj,3) = gphif(im1,jm1) 
    412                      rcphi(ji,jj,4) = gphif(ji,jm1) 
    413  
    414                   END DO 
    415                END DO 
    416  
    417                ! Correction of one (land) grid cell of the orca2 grid. 
    418                ! It was causing problems with the SCRIP interpolation. 
    419  
    420                IF (jpiglo == 182 .AND. jpjglo == 149) THEN 
    421                   rclam(145,106,2) = -1.0 
    422                   rcphi(145,106,2) = 41.0 
    423                ENDIF 
    424  
    425                WHERE (tmask(:,:,1) > 0.5 ) 
    426                   mask(:,:) = 0 
    427                ELSE WHERE 
    428                   mask(:,:) = 1 
    429                END WHERE 
    430  
    431                area = e1t * e2t 
    432  
    433             ELSE IF ( i_grid == 2 ) THEN 
    434  
    435                ! -------------------------------------------------------- 
    436                ! ... Write the grid info for u points 
    437                ! -------------------------------------------------------- 
    438  
    439                point_name = 'opau' 
    440  
    441                glam_b = glamu 
    442                gphi_b = gphiu 
    443  
    444                DO ji = 1, jpi 
    445                   DO jj = 1, jpj 
    446  
    447                      ip1 = ji+1 
    448                      jm1 = jj-1 
    449  
    450                      IF (ji == jpiglo) ip1 = 3 
    451                      IF (jj == 1) jm1 = jj 
    452  
    453                      rclam(ji,jj,1) = glamv(ip1,jj) 
    454                      rclam(ji,jj,2) = glamv(ji,jj) 
    455                      rclam(ji,jj,3) = glamv(ji,jm1) 
    456                      rclam(ji,jj,4) = glamv(ip1,jm1) 
    457  
    458                      rcphi(ji,jj,1) = gphiv(ip1,jj) 
    459                      rcphi(ji,jj,2) = gphiv(ji,jj) 
    460                      rcphi(ji,jj,3) = gphiv(ji,jm1) 
    461                      rcphi(ji,jj,4) = gphiv(ip1,jm1) 
    462  
    463                   END DO 
    464                END DO 
    465  
    466                ! Correction of three (land) grid cell of the orca2 grid. 
    467                ! It was causing problems with the SCRIP interpolation. 
    468  
    469                IF (jpiglo == 182 .AND. jpjglo == 149) THEN 
    470                   glam_b(144,106)   = -1.0 
    471                   gphi_b(144,106)   = 40.5 
    472                   rclam (144,106,2) = -1.5   
    473                   rcphi (144,106,2) = 41.0 
    474  
    475                   glam_b(144,107)   = -1.0 
    476                   gphi_b(144,107)   = 41.5 
    477                   rclam (144,107,2) = -1.5   
    478                   rcphi (144,107,2) = 42.0 
    479                   rclam (144,107,3) = -1.5   
    480                   rcphi (144,107,3) = 41.0 
    481  
    482                   glam_b(144,108)   = -1.0 
    483                   gphi_b(144,108)   = 42.5 
    484                   rclam (144,108,2) = -1.5   
    485                   rcphi (144,108,2) = 43.0 
    486                   rclam (144,108,3) = -1.5   
    487                   rcphi (144,108,3) = 42.0 
    488                ENDIF 
    489  
    490                WHERE (umask(:,:,1) > 0.5 ) 
    491                   mask(:,:) = 0 
    492                ELSE WHERE 
    493                   mask(:,:) = 1 
    494                END WHERE 
    495  
    496                area = e1u * e2u 
    497  
    498             ELSE IF ( i_grid == 3 ) THEN 
    499  
    500                ! -------------------------------------------------------- 
    501                ! ... Write the grid info for v points 
    502                ! -------------------------------------------------------- 
    503  
    504                point_name = 'opav' 
    505  
    506                glam_b = glamv 
    507                gphi_b = gphiv 
    508  
    509                DO ji = 1, jpi 
    510                   DO jj = 1, jpj 
    511  
    512                      im1 = ji-1 
    513                      jp1 = jj+1 
    514                      IF (ji == 1) im1 = jpiglo-2 
    515                      IF (jj == jpjglo) jp1 = jj 
    516  
    517                      rclam(ji,jj,1) = glamu(ji,jp1) 
    518                      rclam(ji,jj,2) = glamu(im1,jp1) 
    519                      rclam(ji,jj,3) = glamu(im1,jj) 
    520                      rclam(ji,jj,4) = glamu(ji,jj) 
    521  
    522                      rcphi(ji,jj,1) = gphiu(ji,jp1) 
    523                      rcphi(ji,jj,2) = gphiu(im1,jp1) 
    524                      rcphi(ji,jj,3) = gphiu(im1,jj) 
    525                      rcphi(ji,jj,4) = gphiu(ji,jj) 
    526  
    527                   END DO 
    528                END DO 
    529  
    530                ! Correction of one (land) grid cell of the orca2 grid. 
    531                ! It was causing problems with the SCRIP interpolation. 
    532  
    533                IF (jpiglo == 182 .AND. jpjglo == 149) THEN 
    534                   rclam(145,105,2) = -1.0   
    535                   rcphi(145,105,2) = 40.5 
    536                ENDIF 
    537  
    538                WHERE (vmask(:,:,1) > 0.5 ) 
    539                   mask(:,:) = 0 
    540                ELSE WHERE 
    541                   mask(:,:) = 1 
    542                END WHERE 
    543  
    544                area = e1v * e2v 
    545  
    546             ENDIF ! i_grid 
    547  
    548             WHERE (glam_b(:,:) < 0.) 
    549                glam_b(:,:) = glam_b(:,:) + 360. 
    550             END WHERE 
    551             WHERE (glam_b(:,:) > 360.) 
    552                glam_b(:,:) = glam_b(:,:) - 360. 
    553             END WHERE 
    554  
    555             WHERE (rclam(:,:,:) < 0.) 
    556                rclam(:,:,:) = rclam(:,:,:) + 360. 
    557             END WHERE 
    558             WHERE (rclam(:,:,:) > 360.) 
    559                rclam(:,:,:) = rclam(:,:,:) - 360. 
    560             END WHERE 
    561  
    562             mask(:,jpjglo)=1 
    563  
    564             CALL prism_write_grid   ( point_name, nx, ny, glam_b, gphi_b )  
    565             CALL prism_write_corner ( point_name, nx, ny, nc, rclam, rcphi ) 
    566             CALL prism_write_mask   ( point_name, nx, ny, mask ) 
    567             CALL prism_write_area   ( point_name, nx, ny, area ) 
    568  
    569          END DO ! i_grid 
    570  
    571          CALL prism_terminate_grids_writing () 
    572  
    573       ENDIF ! localRank == localRoot .and. write_grid == 1 
    574  
     151      ! 
    575152      ! ----------------------------------------------------------------- 
    576153      ! ... Define the partition  
    577154      ! ----------------------------------------------------------------- 
    578  
    579       IF ( rootexchg ) THEN 
    580  
    581          paral(1) = 2              ! box partitioning 
    582          paral(2) = 0              ! NEMO lower left corner global offset     
    583          paral(3) = jpiglo         ! local extent in i  
    584          paral(4) = jpjglo         ! local extent in j 
    585          paral(5) = jpiglo         ! global extent in x 
    586  
    587          range(1) = nimpp-1+nldi   ! global start in i 
    588          range(2) = nlei-nldi+1    ! local size in i of valid region 
    589          range(3) = njmpp-1+nldj   ! global start in j 
    590          range(4) = nlej-nldj+1    ! local size in j of valid region 
    591          range(5) = range(2) & 
    592                   * range(4)       ! local horizontal size 
    593  
    594          IF(ln_ctl) THEN 
    595          write(numout,*) ' rootexchg: range(1:5)', range 
     155       
     156      paral(1) = 2                                              ! box partitioning 
     157      paral(2) = jpiglo * (nldj-1+njmpp-1) + (nldi-1+nimpp-1)   ! NEMO lower left corner global offset     
     158      paral(3) = nlei-nldi+1                                    ! local extent in i  
     159      paral(4) = nlej-nldj+1                                    ! local extent in j 
     160      paral(5) = jpiglo                                         ! global extent in x 
     161       
     162      IF( ln_ctl ) THEN 
     163         WRITE(numout,*) ' multiexchg: paral (1:5)', paral 
     164         WRITE(numout,*) ' multiexchg: jpi, jpj =', jpi, jpj 
     165         WRITE(numout,*) ' multiexchg: nldi, nlei, nimpp =', nldi, nlei, nimpp 
     166         WRITE(numout,*) ' multiexchg: nldj, nlej, njmpp =', nldj, nlej, njmpp 
     167      ENDIF 
     168       
     169      CALL prism_def_partition_proto ( id_part, paral, nerror ) 
     170      ! 
     171      ! ... Announce send variables.  
     172      ! 
     173      DO ji = 1, nsnd 
     174         IF ( ssnd(ji)%laction ) THEN  
     175            CALL prism_def_var_proto (ssnd(ji)%nid, ssnd(ji)%clname, id_part, (/ 2, 0/),  & 
     176               &                      PRISM_Out   , ishape   , PRISM_REAL, nerror) 
     177            IF ( nerror /= PRISM_Ok ) THEN 
     178               WRITE(numout,*) 'Failed to define transient ', ji, TRIM(ssnd(ji)%clname) 
     179               CALL prism_abort_proto ( ssnd(ji)%nid, 'cpl_prism_define', 'Failure in prism_def_var') 
     180            ENDIF 
    596181         ENDIF 
    597  
    598          ! 
    599          ! Collect ranges from all NEMO procs on the local root process 
    600          ! 
    601          CALL mpi_gather(range,  5, MPI_INTEGER, & 
    602                          ranges, 5, MPI_INTEGER, localRoot, localComm, ierror) 
    603  
    604          IF ( localRank == localRoot ) THEN 
    605  
    606             maxlen = maxval(ranges(5,:)) 
    607              
    608             ALLOCATE(buffer(1:maxlen), stat = ierror) 
    609             IF (ierror > 0) THEN 
    610                CALL prism_abort_proto ( comp_id, 'cpl_prism_define', 'Failure in allocating buffer') 
    611                RETURN 
     182      END DO 
     183      ! 
     184      ! ... Announce received variables.  
     185      ! 
     186      DO ji = 1, nrcv 
     187         IF ( srcv(ji)%laction ) THEN  
     188            CALL prism_def_var_proto ( srcv(ji)%nid, srcv(ji)%clname, id_part, (/ 2, 0/),   & 
     189               &                      PRISM_In    , ishape   , PRISM_REAL, nerror) 
     190            IF ( nerror /= PRISM_Ok ) THEN 
     191               WRITE(numout,*) 'Failed to define transient ', ji, TRIM(srcv(ji)%clname) 
     192               CALL prism_abort_proto ( srcv(ji)%nid, 'cpl_prism_define', 'Failure in prism_def_var') 
    612193            ENDIF 
    613  
    614           ENDIF 
    615  
    616       ELSE 
    617  
    618          paral(1) = 2                  ! box partitioning 
    619 !2dtest         paral(2) = jpiglo           & 
    620 !2dtest                  * (nldj-1+njmpp-1) & 
    621 !2dtest                  + (nldi-1+nimpp-1)   ! NEMO lower left corner global offset     
    622          paral(2) = jpiglo & 
    623                   * (nldj-1+njmpp-1)   ! NEMO lower left corner global offset     
    624          paral(3) = nlei-nldi+1        ! local extent in i  
    625          paral(4) = nlej-nldj+1        ! local extent in j 
    626          paral(5) = jpiglo             ! global extent in x 
    627  
    628          IF(ln_ctl) THEN 
    629             print*, ' multiexchg: paral (1:5)', paral 
    630             print*, ' multiexchg: jpi, jpj =', jpi, jpj 
    631             print*, ' multiexchg: nldi, nlei, nimpp =', nldi, nlei, nimpp 
    632             print*, ' multiexchg: nldj, nlej, njmpp =', nldj, nlej, njmpp 
    633194         ENDIF 
    634  
    635          IF ( paral(3) /= nlei-nldi+1 ) THEN 
    636               print*, 'WARNING!!! in cpl_oasis3 - cpl_prism_define' 
    637               print*, 'cpl_prism_define: local extend in i is ', paral(3), ' should equal ', nlei-nldi+1 
    638          ENDIF 
    639          IF ( paral(4) /= nlej-nldj+1 ) THEN 
    640               print*, 'WARNING!!! in cpl_oasis3 - cpl_prism_define' 
    641               print*, 'cpl_prism_define: local extend in j is ', paral(4), ' should equal ', nlej-nldj+1 
    642          ENDIF 
    643  
    644       ENDIF 
    645  
    646       IF ( commRank ) & 
    647       CALL prism_def_partition_proto ( part_id, paral, ierror ) 
    648  
    649       grid_id(1)= part_id 
    650  
    651       !------------------------------------------------------------------ 
    652       ! 3rd Declare the transient variables 
    653       !------------------------------------------------------------------ 
    654       ! 
    655       ! ... Define symbolic names for the transient fields send by the ocean 
    656       !     These must be identical to the names specified in the SMIOC file. 
    657       ! 
    658       cpl_send( 1)='SSTOCEAN' ! sea surface temperature              -> sst_io 
    659       cpl_send( 2)='SICOCEAN' ! sea ice area fraction                -> 1.-frld 
    660 #if defined key_cpl_albedo 
    661       cpl_send( 3)='STIOCEAN' ! surface temperature over sea ice     -> tn_ice 
    662       cpl_send( 4)='SAIOCEAN' ! albedo over sea ice                  -> alb_ice 
    663 #else 
    664       cpl_send( 3)='SITOCEAN' ! sea ice thickness                    -> hicif (only 1 layer available!) 
    665       cpl_send( 4)='SNTOCEAN' ! surface snow thickness over sea ice  -> hsnif 
    666 #endif 
    667 #if defined key_cpl_ocevel 
    668       cpl_send( 5)='SUNOCEAN' ! U-velocity                           -> un 
    669       cpl_send( 6)='SVNOCEAN' ! V-velocity                           -> vn 
    670 #endif 
    671       ! 
    672       ! ...  Define symbolic names for transient fields received by the ocean. 
    673       !      These must be identical to the names specified in the SMIOC file. 
    674       ! 
    675       ! ...  a) U-Grid fields 
    676       ! 
    677       cpl_recv( 1)='TXWOCEWU' ! weighted surface downward eastward stress 
    678       cpl_recv( 2)='TYWOCEWU' ! weighted surface downward northward stress 
    679       cpl_recv( 3)='TXIOCEWU' ! weighted surface downward eastward stress over ice 
    680       cpl_recv( 4)='TYIOCEWU' ! weighted surface downward northward stress over ice 
    681       ! 
    682       ! ...  a) V-Grid fields 
    683       ! 
    684       cpl_recv( 5)='TXWOCEWV' ! weighted surface downward eastward stress 
    685       cpl_recv( 6)='TYWOCEWV' ! weighted surface downward northward stress 
    686       cpl_recv( 7)='TXIOCEWV' ! weighted surface downward eastward stress over ice 
    687       cpl_recv( 8)='TYIOCEWV' ! weighted surface downward northward stress over ice 
    688       ! 
    689       ! ...  a) T-Grid fields 
    690       ! 
    691       cpl_recv( 9)='FRWOCEPE' ! P-E over water                               -> zpew 
    692       cpl_recv(10)='FRIOCEPE' ! P-E over ice                                 -> zpei 
    693       cpl_recv(11)='FRROCESN' ! surface downward snow fall                   -> zpsol 
    694       cpl_recv(12)='FRIOCEEV' ! surface upward snow flux where sea ice       -> zevice 
    695  
    696       cpl_recv(13)='QSWOCESR' ! surface net downward shortwave flux          -> qsr_oce 
    697       cpl_recv(14)='QSWOCENS' ! surface downward non-solar heat flux in air  -> qnsr_oce 
    698       cpl_recv(15)='QSIOCESR' ! solar heat flux on sea ice                   -> qsr_ice 
    699       cpl_recv(16)='QSIOCENS' ! non-solar heat flux on sea ice               -> qnsr_ice 
    700       cpl_recv(17)='QSIOCEDQ' ! non-solar heat flux derivative               -> dqns_ice 
    701  
    702 #ifdef key_cpl_discharge 
    703       cpl_recv(18)='FRWOCEIS' ! ice discharge into ocean                     -> calving 
    704       cpl_recv(19)='FRWOCERD' ! river discharge into ocean                   -> zrunriv 
    705       cpl_recv(20)='FRWOCECD' ! continental discharge into ocean             -> zruncot 
    706 #endif 
    707       ! 
    708       ! data_type has to be PRISM_REAL as PRISM_DOUBLE is not supported. 
    709       ! For exchange of double precision fields the OASIS3 has to be compiled 
    710       ! with use_realtype_single. (see OASIS3 User Guide prism_2-4, 5th Ed., 
    711       ! p. 13 and p. 53 for further explanation.) 
    712       ! 
    713       data_type = PRISM_REAL 
    714  
    715       nodim(1) = 3 ! check 
    716       nodim(2) = 0 
    717  
    718       ! 
    719       ! ... Define the shape for the area that excludes the halo 
    720       !     For serial configuration (key_mpp_mpi not being active) 
    721       !     nl* is set to the global values 1 and jp*glo. 
    722       ! 
    723       IF ( rootexchg ) THEN 
    724          shape(1,1) = 1 
    725          shape(2,1) = jpiglo 
    726          shape(1,2) = 1 
    727          shape(2,2) = jpjglo 
    728          shape(1,3) = 1 
    729          shape(2,3) = 1 
    730        ELSE 
    731          shape(1,1) = 1 
    732          shape(2,1) = nlei-nldi+1 ! jpi 
    733          shape(1,2) = 1 
    734          shape(2,2) = nlej-nldj+1 ! jpj 
    735          shape(1,3) = 1 
    736          shape(2,3) = 1 
    737       ENDIF 
    738       ! 
    739       ! ----------------------------------------------------------------- 
    740       ! ... Allocate memory for data exchange 
    741       ! ----------------------------------------------------------------- 
    742       ! 
    743       ALLOCATE(exfld(shape(1,1):shape(2,1),shape(1,2):shape(2,2)), stat = ierror) 
    744       IF (ierror > 0) THEN 
    745          CALL prism_abort_proto ( comp_id, 'cpl_prism_define', 'Failure in allocating exfld') 
    746          RETURN 
    747       ENDIF 
    748       ! 
    749       ! ... Announce send variables, all on T points.  
    750       ! 
    751       info = PRISM_Out 
    752       ! 
    753  
    754       IF ( commRank ) THEN 
    755  
    756          DO ji = 1, nsend 
    757             !        if ( ji == 2 ) ; then ; nodim(2) = 2 ; else ; nodim(2) = 0 ; endif 
    758             CALL prism_def_var_proto (send_id(ji), cpl_send(ji), grid_id(1), & 
    759                  nodim, info, shape, data_type, ierror) 
    760             IF ( ierror /= PRISM_Ok ) THEN 
    761                PRINT *, 'Failed to define transient ', ji, TRIM(cpl_send(ji)) 
    762                CALL prism_abort_proto ( comp_id, 'cpl_prism_define', 'Failure in prism_def_var') 
    763             ENDIF 
    764          ENDDO 
    765          ! 
    766          nodim(1) = 3 ! check 
    767          nodim(2) = 0 
    768          ! 
    769          ! ... Announce recv variables.  
    770          ! 
    771          info = PRISM_In 
    772          ! 
    773          ! ... a) on U points 
    774          ! 
    775          DO ji = 1, 4 
    776             CALL prism_def_var_proto (recv_id(ji), cpl_recv(ji), grid_id(1), & 
    777                  nodim, info, shape, data_type, ierror) 
    778             IF ( ierror /= PRISM_Ok ) THEN 
    779                PRINT *, 'Failed to define transient ', ji, TRIM(cpl_recv(ji)) 
    780                CALL prism_abort_proto ( comp_id, 'cpl_prism_define', 'Failure in prism_def_var') 
    781             ENDIF 
    782          ENDDO 
    783          ! 
    784          ! ... b) on V points 
    785          ! 
    786          DO ji = 5, 8 
    787             CALL prism_def_var_proto (recv_id(ji), cpl_recv(ji), grid_id(1), & 
    788                  nodim, info, shape, data_type, ierror) 
    789             IF ( ierror /= PRISM_Ok ) THEN 
    790                PRINT *, 'Failed to define transient ', TRIM(cpl_recv(ji)) 
    791                CALL prism_abort_proto ( comp_id, 'cpl_prism_define', 'Failure in prism_def_var') 
    792             ENDIF 
    793          ENDDO 
    794          ! 
    795          ! ... c) on T points 
    796          ! 
    797          DO ji = 9, nrecv 
    798             CALL prism_def_var_proto (recv_id(ji), cpl_recv(ji), grid_id(1), & 
    799                  nodim, info, shape, data_type, ierror) 
    800             IF ( ierror /= PRISM_Ok ) THEN 
    801                PRINT *, 'Failed to define transient ', TRIM(cpl_recv(ji)) 
    802                CALL prism_abort_proto ( comp_id, 'cpl_prism_define', 'Failure in prism_def_var') 
    803             ENDIF 
    804          ENDDO 
    805  
    806       ENDIF ! commRank 
    807  
    808       !------------------------------------------------------------------ 
    809       ! 4th End of definition phase 
    810       !------------------------------------------------------------------ 
    811  
    812       IF ( commRank ) THEN 
    813          CALL prism_enddef_proto(ierror) 
    814          IF ( ierror /= PRISM_Ok ) & 
    815               CALL prism_abort_proto ( comp_id, 'cpl_prism_define', 'Failure in prism_enddef') 
    816       ENDIF 
    817  
     195      END DO 
     196       
     197      !------------------------------------------------------------------ 
     198      ! End of definition phase 
     199      !------------------------------------------------------------------ 
     200       
     201      CALL prism_enddef_proto(nerror) 
     202      IF ( nerror /= PRISM_Ok )   CALL prism_abort_proto ( ncomp_id, 'cpl_prism_define', 'Failure in prism_enddef') 
     203       
    818204   END SUBROUTINE cpl_prism_define 
    819  
    820  
    821  
    822    SUBROUTINE cpl_prism_send( var_id, date, data_array, info ) 
    823  
    824       IMPLICIT NONE 
     205    
     206    
     207   SUBROUTINE cpl_prism_snd( kid, kstep, pdata, kinfo ) 
    825208 
    826209      !!--------------------------------------------------------------------- 
    827       !!              ***  ROUTINE cpl_prism_send  *** 
     210      !!              ***  ROUTINE cpl_prism_snd  *** 
    828211      !! 
    829212      !! ** Purpose : - At each coupling time-step,this routine sends fields 
     
    832215      !! * Arguments 
    833216      !! 
    834       INTEGER, INTENT( IN )  :: var_id    ! variable Id 
    835       INTEGER, INTENT( OUT ) :: info      ! OASIS3 info argument 
    836       INTEGER, INTENT( IN )  :: date      ! ocean time-step in seconds 
    837       REAL(wp)               :: data_array(:,:) 
    838       !! 
    839       !! * Local declarations 
    840       !! 
    841 #if defined key_mpp_mpi 
    842       REAL(wp)               :: global_array(jpiglo,jpjglo) 
    843       ! 
    844 !mpi  INTEGER                :: status(MPI_STATUS_SIZE) 
    845 !mpi  INTEGER                :: type       ! MPI data type 
    846       INTEGER                :: request    ! MPI isend request 
    847       INTEGER                :: ji, jj, jn ! local loop indicees 
    848 #else 
    849       INTEGER                :: ji 
    850 #endif 
    851       !! 
    852       !!-------------------------------------------------------------------- 
    853       !! 
    854  
    855 #if defined key_mpp_mpi 
    856  
    857       request = 0 
    858  
    859       IF ( rootexchg ) THEN 
    860          ! 
    861 !mpi     IF ( wp == 4 ) type = MPI_REAL 
    862 !mpi     IF ( wp == 8 ) type = MPI_DOUBLE_PRECISION 
    863          ! 
    864          ! collect data on the local root process 
    865          ! 
    866  
    867          if ( var_id == 1 .and. localRank == localRoot .and. ln_ctl )  then 
    868              do ji = 0, localSize-1 
    869                 WRITE(numout,*) ' rootexchg: ranges for rank ', ji, ' are ', ranges(:,ji)  
    870              enddo 
    871          endif 
    872  
    873          IF ( localRank /= localRoot ) THEN 
    874  
    875             DO jj = nldj, nlej 
    876                DO ji = nldi, nlei 
    877                   exfld(ji-nldi+1,jj-nldj+1)=data_array(ji,jj) 
    878                ENDDO 
    879             ENDDO 
    880  
    881 !mpi        CALL mpi_send(exfld, range(5), type, localRoot, localRank, localComm, ierror) 
    882             CALL mppsend (localRank, exfld, range(5), localRoot, request)   
    883  
    884             if ( var_id == 1 .and. ln_ctl )  then 
    885                WRITE(numout,*) ' rootexchg: This is process       ', localRank 
    886                WRITE(numout,*) ' rootexchg: We have a range of    ', range  
    887 !               WRITE(numout,*) ' rootexchg: We got SST to process ', data_array  
    888             endif 
    889  
    890          ENDIF 
    891  
    892          IF ( localRank == localRoot ) THEN 
    893  
    894             DO jj = ranges(3,localRoot), ranges(3,localRoot)+ranges(4,localRoot)-1 
    895                DO ji = ranges(1,localRoot), ranges(1,localRoot)+ranges(2,localRoot)-1 
    896                   global_array(ji,jj) = data_array(ji,jj) ! workaround 
    897                ENDDO 
    898             ENDDO 
    899  
    900             DO jn = 1, localSize-1 
    901  
    902 !mpi           CALL mpi_recv(buffer, ranges(5,jn), type, localRoot, jn, localComm, status, ierror) 
    903                CALL mpprecv(jn, buffer, ranges(5,jn)) 
    904  
    905                if ( var_id == 1 .and. ln_ctl )  then 
    906                    WRITE(numout,*) ' rootexchg: Handling data from process ', jn 
    907 !                   WRITE(numout,*) ' rootexchg: We got SST to process      ', buffer 
    908                endif 
    909  
    910  
    911                DO jj = ranges(3,jn), ranges(3,jn)+ranges(4,jn)-1 
    912                   DO ji = ranges(1,jn), ranges(1,jn)+ranges(2,jn)-1 
    913                      global_array(ji,jj) = buffer((jj-ranges(3,jn))*ranges(2,jn) + ji-ranges(1,jn)+1) 
    914                   ENDDO 
    915                ENDDO 
    916  
    917             ENDDO 
    918  
    919             CALL prism_put_proto ( var_id, date, global_array, info ) 
    920  
    921          ENDIF 
    922  
    923       ELSE 
    924  
    925          DO jj = nldj, nlej 
    926             DO ji = nldi, nlei 
    927                exfld(ji-nldi+1,jj-nldj+1)=data_array(ji,jj) 
    928             ENDDO 
    929          ENDDO 
    930  
    931          CALL prism_put_proto ( var_id, date, exfld, info ) 
    932  
     217      INTEGER,                      INTENT( IN    )   :: kid       ! variable intex in the array 
     218      INTEGER,                      INTENT(   OUT )   :: kinfo     ! OASIS3 info argument 
     219      INTEGER,                      INTENT( IN    )   :: kstep     ! ocean time-step in seconds 
     220      REAL(wp), DIMENSION(jpi,jpj), INTENT( IN    )   :: pdata 
     221      !! 
     222      !! 
     223      !!-------------------------------------------------------------------- 
     224      ! 
     225      ! snd data to OASIS3 
     226      ! 
     227      IF( lk_mpp ) THEN   ;   CALL prism_put_proto ( ssnd(kid)%nid, kstep, pdata(nldi:nlei, nldj:nlej), kinfo ) 
     228      ELSE                ;   CALL prism_put_proto ( ssnd(kid)%nid, kstep, pdata                      , kinfo ) 
    933229      ENDIF 
    934  
    935 #else 
    936  
    937       ! 
    938       ! send local data from every process to OASIS3 
    939       ! 
    940       IF ( commRank ) & 
    941       CALL prism_put_proto ( var_id, date, data_array, info ) 
    942  
    943 #endif 
    944  
    945       IF ( commRank ) THEN 
    946  
    947          IF (ln_ctl .and. lwp) THEN         
    948  
    949             IF ( info == PRISM_Sent     .OR. & 
    950                  info == PRISM_ToRest   .OR. & 
    951                  info == PRISM_SentOut  .OR. & 
    952                  info == PRISM_ToRestOut       ) THEN 
    953                WRITE(numout,*) '****************' 
    954                DO ji = 1, nsend 
    955                   IF (var_id == send_id(ji) ) THEN 
    956                      WRITE(numout,*) 'prism_put_proto: Outgoing ', cpl_send(ji) 
    957                      EXIT 
    958                   ENDIF 
    959                ENDDO 
    960                WRITE(numout,*) 'prism_put_proto: var_id ', var_id 
    961                WRITE(numout,*) 'prism_put_proto:   date ', date 
    962                WRITE(numout,*) 'prism_put_proto:   info ', info 
    963                WRITE(numout,*) '     - Minimum value is ', MINVAL(data_array) 
    964                WRITE(numout,*) '     - Maximum value is ', MAXVAL(data_array) 
    965                WRITE(numout,*) '     -     Sum value is ', SUM(data_array) 
    966                WRITE(numout,*) '****************' 
    967             ENDIF 
    968  
    969          ENDIF 
    970  
    971       ENDIF 
    972  
    973    END SUBROUTINE cpl_prism_send 
    974  
    975  
    976  
    977    SUBROUTINE cpl_prism_recv( var_id, date, data_array, info ) 
    978  
    979       IMPLICIT NONE 
     230       
     231      IF ( ln_ctl ) THEN         
     232         IF ( kinfo == PRISM_Sent     .OR. kinfo == PRISM_ToRest .OR.   & 
     233            & kinfo == PRISM_SentOut  .OR. kinfo == PRISM_ToRestOut ) THEN 
     234            WRITE(numout,*) '****************' 
     235            WRITE(numout,*) 'prism_put_proto: Outgoing ', ssnd(kid)%clname 
     236            WRITE(numout,*) 'prism_put_proto: ivarid ', ssnd(kid)%nid 
     237            WRITE(numout,*) 'prism_put_proto:  kstep ', kstep 
     238            WRITE(numout,*) 'prism_put_proto:   info ', kinfo 
     239            WRITE(numout,*) '     - Minimum value is ', MINVAL(pdata) 
     240            WRITE(numout,*) '     - Maximum value is ', MAXVAL(pdata) 
     241            WRITE(numout,*) '     -     Sum value is ', SUM(pdata) 
     242            WRITE(numout,*) '****************' 
     243        ENDIF 
     244     ENDIF 
     245    END SUBROUTINE cpl_prism_snd 
     246 
     247 
     248   SUBROUTINE cpl_prism_rcv( kid, kstep, pdata, kinfo ) 
    980249 
    981250      !!--------------------------------------------------------------------- 
    982       !!              ***  ROUTINE cpl_prism_recv  *** 
     251      !!              ***  ROUTINE cpl_prism_rcv  *** 
    983252      !! 
    984253      !! ** Purpose : - At each coupling time-step,this routine receives fields 
    985254      !!      like stresses and fluxes from the coupler or remote application. 
    986255      !!---------------------------------------------------------------------- 
    987       !! * Arguments 
    988       !! 
    989       INTEGER, INTENT( IN )  :: var_id    ! variable Id 
    990       INTEGER, INTENT( OUT ) :: info      ! variable Id 
    991       INTEGER, INTENT( IN )  :: date      ! ocean time-step in seconds 
    992       REAL(wp),INTENT( OUT ) :: data_array(:,:) 
    993       !! 
    994       !! * Local declarations 
    995       !! 
    996 #if defined key_mpp_mpi 
    997       REAL(wp)               :: global_array(jpiglo,jpjglo) 
    998       ! 
    999 !      LOGICAL                :: action = .false. 
    1000       LOGICAL                :: action 
    1001 !mpi  INTEGER                :: status(MPI_STATUS_SIZE) 
    1002 !mpi  INTEGER                :: type       ! MPI data type 
    1003       INTEGER                :: request    ! MPI isend request 
    1004       INTEGER                :: ji, jj, jn ! local loop indices 
    1005 #else 
    1006       INTEGER                :: ji 
    1007 #endif 
    1008       !! 
    1009       !!-------------------------------------------------------------------- 
    1010       !! 
    1011 #ifdef key_mpp_mpi 
    1012       action = .false. 
    1013       request = 0 
    1014  
    1015       IF ( rootexchg ) THEN 
    1016          ! 
    1017          ! receive data from OASIS3 on local root 
    1018          ! 
    1019          IF ( commRank ) & 
    1020               CALL prism_get_proto ( var_id, date, global_array, info ) 
    1021  
    1022          CALL MPI_BCAST ( info, 1, MPI_INTEGER, localRoot, localComm, ierror ) 
    1023  
    1024       ELSE 
    1025          ! 
    1026          ! receive local data from OASIS3 on every process 
    1027          ! 
    1028          CALL prism_get_proto ( var_id, date, exfld, info ) 
    1029  
     256      INTEGER,                      INTENT( IN    )   :: kid       ! variable intex in the array 
     257      INTEGER,                      INTENT( IN    )   :: kstep     ! ocean time-step in seconds 
     258      REAL(wp), DIMENSION(jpi,jpj), INTENT( INOUT )   :: pdata     ! IN to keep the value if nothing is done 
     259      INTEGER,                      INTENT(   OUT )   :: kinfo     ! OASIS3 info argument 
     260      !! 
     261      LOGICAL                :: llaction 
     262      !!-------------------------------------------------------------------- 
     263      ! 
     264      ! receive local data from OASIS3 on every process 
     265      ! 
     266      CALL prism_get_proto ( srcv(kid)%nid, kstep, exfld, kinfo )          
     267 
     268      llaction = .false. 
     269      IF( kinfo == PRISM_Recvd   .OR. kinfo == PRISM_FromRest .OR.   & 
     270          kinfo == PRISM_RecvOut .OR. kinfo == PRISM_FromRestOut )   llaction = .TRUE. 
     271 
     272      IF ( ln_ctl )   WRITE(numout,*) "llaction, kinfo, kstep, ivarid: " , llaction, kinfo, kstep, srcv(kid)%nid 
     273 
     274      IF ( llaction ) THEN 
     275 
     276         IF( lk_mpp ) THEN   ;   pdata(nldi:nlei, nldj:nlej) = exfld(:,:) 
     277         ELSE                ;   pdata(    :    ,     :    ) = exfld(:,:) 
     278         ENDIF 
     279          
     280         !--- Fill the overlap areas and extra hallows (mpp) 
     281         !--- check periodicity conditions (all cases) 
     282         CALL lbc_lnk( pdata, srcv(kid)%clgrid, srcv(kid)%nsgn )    
     283          
     284         IF ( ln_ctl ) THEN         
     285            WRITE(numout,*) '****************' 
     286            WRITE(numout,*) 'prism_get_proto: Incoming ', srcv(kid)%clname 
     287            WRITE(numout,*) 'prism_get_proto: ivarid '  , srcv(kid)%nid 
     288            WRITE(numout,*) 'prism_get_proto:   kstep', kstep 
     289            WRITE(numout,*) 'prism_get_proto:   info ', kinfo 
     290            WRITE(numout,*) '     - Minimum value is ', MINVAL(pdata) 
     291            WRITE(numout,*) '     - Maximum value is ', MAXVAL(pdata) 
     292            WRITE(numout,*) '     -     Sum value is ', SUM(pdata) 
     293            WRITE(numout,*) '****************' 
     294            call flush(numout) 
     295         ENDIF 
     296       
    1030297      ENDIF 
    1031298 
    1032       IF ( info == PRISM_Recvd        .OR. & 
    1033            info == PRISM_FromRest     .OR. & 
    1034            info == PRISM_RecvOut      .OR. & 
    1035            info == PRISM_FromRestOut ) action = .true. 
    1036  
    1037       IF (ln_ctl .and. lwp) THEN         
    1038          WRITE(numout,*) "info", info, var_id 
    1039          WRITE(numout,*) "date", date, var_id 
    1040          WRITE(numout,*) "action", action, var_id 
    1041       ENDIF 
    1042  
    1043       IF ( rootexchg .and. action ) THEN 
    1044          ! 
    1045 !mpi     IF ( wp == 4 ) type = MPI_REAL 
    1046 !mpi     IF ( wp == 8 ) type = MPI_DOUBLE_PRECISION 
    1047          ! 
    1048          ! distribute data to processes 
    1049          ! 
    1050          IF ( localRank == localRoot ) THEN 
    1051  
    1052             DO jj = ranges(3,localRoot), ranges(3,localRoot)+ranges(4,localRoot)-1 
    1053                DO ji = ranges(1,localRoot), ranges(1,localRoot)+ranges(2,localRoot)-1 
    1054                   exfld(ji-ranges(1,localRoot)+1,jj-ranges(3,localRoot)+1) = global_array(ji,jj) 
    1055                ENDDO 
    1056             ENDDO 
    1057  
    1058             DO jn = 1, localSize-1 
    1059  
    1060                DO jj = ranges(3,jn), ranges(3,jn)+ranges(4,jn)-1 
    1061                   DO ji = ranges(1,jn), ranges(1,jn)+ranges(2,jn)-1 
    1062                      buffer((jj-ranges(3,jn))*ranges(2,jn) + ji-ranges(1,jn)+1) = global_array(ji,jj) 
    1063                   ENDDO 
    1064                ENDDO 
    1065  
    1066 !mpi           CALL mpi_send(buffer, ranges(5,jn), type, jn, jn, localComm, ierror) 
    1067                CALL mppsend (jn, buffer, ranges(5,jn), jn, request)   
    1068  
    1069             ENDDO 
    1070  
    1071          ENDIF 
    1072  
    1073          IF ( localRank /= localRoot ) THEN 
    1074 !mpi         CALL mpi_recv(exfld, range(5), type, localRoot, localRank, localComm, status, ierror) 
    1075              CALL mpprecv(localRank, exfld, range(5)) 
    1076          ENDIF 
    1077  
    1078       ENDIF 
    1079  
    1080       IF ( action ) THEN 
    1081  
    1082          data_array = 0.0 
    1083  
    1084          DO jj = nldj, nlej 
    1085             DO ji = nldi, nlei 
    1086                data_array(ji,jj)=exfld(ji-nldi+1,jj-nldj+1) 
    1087             ENDDO 
    1088          ENDDO 
    1089  
    1090          IF (ln_ctl .and. lwp) THEN         
    1091             WRITE(numout,*) '****************' 
    1092             DO ji = 1, nrecv 
    1093                IF (var_id == recv_id(ji) ) THEN 
    1094                   WRITE(numout,*) 'prism_get_proto: Incoming ', cpl_recv(ji) 
    1095                   EXIT 
    1096                ENDIF 
    1097             ENDDO 
    1098             WRITE(numout,*) 'prism_get_proto: var_id ', var_id 
    1099             WRITE(numout,*) 'prism_get_proto:   date ', date 
    1100             WRITE(numout,*) 'prism_get_proto:   info ', info 
    1101             WRITE(numout,*) '     - Minimum value is ', MINVAL(data_array) 
    1102             WRITE(numout,*) '     - Maximum value is ', MAXVAL(data_array) 
    1103             WRITE(numout,*) '     -     Sum value is ', SUM(data_array) 
    1104             WRITE(numout,*) '****************' 
    1105          ENDIF 
    1106  
    1107       ENDIF 
    1108 #else 
    1109       CALL prism_get_proto ( var_id, date, exfld, info) 
    1110        
    1111       IF (info == PRISM_Recvd        .OR. & 
    1112           info == PRISM_FromRest     .OR. & 
    1113           info == PRISM_RecvOut      .OR. & 
    1114           info == PRISM_FromRestOut )      THEN 
    1115              data_array = exfld 
    1116  
    1117          IF (ln_ctl .and. lwp ) THEN         
    1118             WRITE(numout,*) '****************' 
    1119             DO ji = 1, nrecv 
    1120                IF (var_id == recv_id(ji) ) THEN 
    1121                   WRITE(numout,*) 'prism_get_proto: Incoming ', cpl_recv(ji) 
    1122                   EXIT 
    1123                ENDIF 
    1124             ENDDO 
    1125             WRITE(numout,*) 'prism_get_proto: var_id ', var_id 
    1126             WRITE(numout,*) 'prism_get_proto:   date ', date 
    1127             WRITE(numout,*) 'prism_get_proto:   info ', info 
    1128             WRITE(numout,*) '     - Minimum value is ', MINVAL(data_array) 
    1129             WRITE(numout,*) '     - Maximum value is ', MAXVAL(data_array) 
    1130             WRITE(numout,*) '     -     Sum value is ', SUM(data_array) 
    1131             WRITE(numout,*) '****************' 
    1132          ENDIF 
    1133  
    1134        ENDIF 
    1135 #endif 
    1136  
    1137    END SUBROUTINE cpl_prism_recv 
    1138  
     299   END SUBROUTINE cpl_prism_rcv 
    1139300 
    1140301 
    1141302   SUBROUTINE cpl_prism_finalize 
    1142  
    1143       IMPLICIT NONE 
    1144303 
    1145304      !!--------------------------------------------------------------------- 
     
    1152311 
    1153312      DEALLOCATE(exfld) 
    1154  
    1155       if ( prism_was_initialized ) then 
    1156  
    1157          if ( prism_was_terminated ) then 
    1158             print *, 'prism has already been terminated.' 
    1159          else 
    1160             call prism_terminate_proto ( ierror ) 
    1161             prism_was_terminated = .true. 
    1162          endif 
    1163  
    1164       else 
    1165  
    1166          print *, 'Initialize prism before terminating it.' 
    1167  
    1168       endif 
    1169  
     313      CALL prism_terminate_proto ( nerror )          
    1170314 
    1171315   END SUBROUTINE cpl_prism_finalize 
    1172316 
     317#else 
     318 
     319   !!---------------------------------------------------------------------- 
     320   !!   Default case                                Forced Ocean/Atmosphere 
     321   !!---------------------------------------------------------------------- 
     322   !!   Empty module 
     323   !!---------------------------------------------------------------------- 
     324   USE in_out_manager               ! I/O manager 
     325   LOGICAL, PUBLIC, PARAMETER :: lk_cpl = .FALSE.   !: coupled flag 
     326   PUBLIC cpl_prism_init 
     327   PUBLIC cpl_prism_finalize 
     328 
     329CONTAINS 
     330 
     331   SUBROUTINE cpl_prism_init 
     332      WRITE(numout,*) 'cpl_prism_init: Error you sould not be there...' 
     333   END SUBROUTINE cpl_prism_init 
     334 
     335   SUBROUTINE cpl_prism_finalize 
     336      WRITE(numout,*) 'cpl_prism_finalize: Error you sould not be there...' 
     337   END SUBROUTINE cpl_prism_finalize 
     338 
    1173339#endif 
    1174340 
Note: See TracChangeset for help on using the changeset viewer.