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 5372 for branches/2015 – NEMO

Changeset 5372 for branches/2015


Ignore:
Timestamp:
2015-06-05T21:14:36+02:00 (9 years ago)
Author:
mcastril
Message:

ticket #1523 Message Packing

Location:
branches/2015/dev_r5302_CNRS18_HPC_scalability/NEMOGCM/NEMO
Files:
5 edited

Legend:

Unmodified
Added
Removed
  • branches/2015/dev_r5302_CNRS18_HPC_scalability/NEMOGCM/NEMO/LIM_SRC_3/limadv.F90

    r5123 r5372  
    206206      END DO 
    207207 
     208 
     209#if defined key_multisend 
     210      CALL lbc_lnk_multi(  psm , 'T',  1. ,  ps0 , 'T',  1. ,  psx , 'T', -1. ,  psy , 'T', -1. ,  psxx, 'T',  1.  ,  psyy, 'T',  1. ,  psxy, 'T',  1. ) 
     211#else 
    208212      !-- Lateral boundary conditions 
    209213      CALL lbc_lnk( psm , 'T',  1. )   ;   CALL lbc_lnk( ps0 , 'T',  1. ) 
     
    211215      CALL lbc_lnk( psxx, 'T',  1. )   ;   CALL lbc_lnk( psyy, 'T',  1. ) 
    212216      CALL lbc_lnk( psxy, 'T',  1. ) 
     217#endif 
    213218 
    214219      IF(ln_ctl) THEN 
     
    392397      END DO 
    393398 
     399#if defined key_multisend 
     400      CALL lbc_lnk_multi(  psm , 'T',  1. ,  ps0 , 'T',  1. ,  psx , 'T', -1. ,  psy , 'T', -1. ,  psxx, 'T',  1.  ,  psyy, 'T',  1. ,  psxy, 'T',  1. ) 
     401#else 
    394402      !-- Lateral boundary conditions 
    395403      CALL lbc_lnk( psm , 'T',  1. )   ;   CALL lbc_lnk( ps0 , 'T',  1. ) 
     
    397405      CALL lbc_lnk( psxx, 'T',  1. )   ;   CALL lbc_lnk( psyy, 'T',  1. ) 
    398406      CALL lbc_lnk( psxy, 'T',  1. ) 
     407#endif 
     408 
    399409 
    400410      IF(ln_ctl) THEN 
  • branches/2015/dev_r5302_CNRS18_HPC_scalability/NEMOGCM/NEMO/LIM_SRC_3/limrhg.F90

    r5123 r5372  
    377377            END DO 
    378378         END DO 
    379          CALL lbc_lnk( v_ice1  , 'U', -1. )   ;   CALL lbc_lnk( u_ice2  , 'V', -1. )      ! lateral boundary cond. 
     379 
     380 
     381#if defined key_multisend 
     382         CALL lbc_lnk_multi( v_ice1, 'U', -1. , u_ice2, 'V', -1. )      ! lateral boundary cond. 
     383#else 
     384         CALL lbc_lnk( v_ice1, 'U', -1. )   ;   CALL lbc_lnk( u_ice2, 'V', -1. )      ! lateral boundary cond. 
     385#endif 
    380386 
    381387         DO jj = k_j1+1, k_jpj-1 
     
    412418            END DO 
    413419         END DO 
     420 
     421#if defined key_multisend 
     422         CALL lbc_lnk_multi( zs1 , 'T', 1. , zs2, 'T', 1. , zs12, 'F', 1. ) 
     423          
     424#else 
    414425         CALL lbc_lnk( zs1 , 'T', 1. )   ;   CALL lbc_lnk( zs2, 'T', 1. ) 
    415426         CALL lbc_lnk( zs12, 'F', 1. ) 
     427#endif 
     428 
    416429 
    417430         ! Ice internal stresses (Appendix C of Hunke and Dukowicz, 2002) 
     
    570583      END DO 
    571584 
    572       CALL lbc_lnk( u_ice(:,:), 'U', -1. )  
    573       CALL lbc_lnk( v_ice(:,:), 'V', -1. )  
     585#if defined key_multisend 
     586      CALL lbc_lnk_multi( u_ice(:,:), 'U', -1. , v_ice(:,:), 'V', -1. ) 
     587 
     588#else 
     589      CALL lbc_lnk( u_ice(:,:), 'U', -1. ) 
     590      CALL lbc_lnk( v_ice(:,:), 'V', -1. ) 
     591#endif 
     592 
    574593#if defined key_agrif && defined key_lim2 
    575594      CALL agrif_rhg_lim2( nn_nevp , nn_nevp, 'U' ) 
     
    595614      END DO 
    596615 
    597       CALL lbc_lnk( u_ice2(:,:), 'V', -1. )  
     616#if defined key_multisend 
     617      CALL lbc_lnk_multi( u_ice2(:,:), 'V', -1. , v_ice1(:,:), 'U', -1. ) 
     618 
     619#else 
     620      CALL lbc_lnk( u_ice2(:,:), 'V', -1. ) 
    598621      CALL lbc_lnk( v_ice1(:,:), 'U', -1. ) 
     622#endif 
    599623 
    600624      ! Recompute delta, shear and div, inputs for mechanical redistribution  
     
    642666      END DO 
    643667 
     668#if defined key_multisend 
     669      ! Lateral boundary condition 
     670      CALL lbc_lnk_multi( divu_i (:,:), 'T', 1. , delta_i(:,:), 'T', 1. ,  shear_i(:,:), 'T', 1. ) 
     671      ! CALL lbc_lnk( shear_i(:,:), 'F', 1. ) 
     672 
     673#else 
    644674      ! Lateral boundary condition 
    645675      CALL lbc_lnk( divu_i (:,:), 'T', 1. ) 
     
    647677      ! CALL lbc_lnk( shear_i(:,:), 'F', 1. ) 
    648678      CALL lbc_lnk( shear_i(:,:), 'T', 1. ) 
     679#endif 
    649680 
    650681      ! * Store the stress tensor for the next time step 
  • branches/2015/dev_r5302_CNRS18_HPC_scalability/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg_ts.F90

    r5217 r5372  
    557557               END DO 
    558558            END DO 
     559#if defined key_multisend 
     560            CALL lbc_lnk_multi( zwx, 'U', 1._wp , zwy, 'V', 1._wp ) 
     561#else 
    559562            CALL lbc_lnk( zwx, 'U', 1._wp )    ;   CALL lbc_lnk( zwy, 'V', 1._wp ) 
     563#endif 
    560564            ! 
    561565            zhup2_e (:,:) = hu_0(:,:) + zwx(:,:)                ! Ocean depth at U- and V-points 
     
    635639               END DO 
    636640            END DO 
     641#if defined key_multisend 
     642            CALL lbc_lnk_multi( zsshu_a, 'U', 1._wp , zsshv_a, 'V', 1._wp ) 
     643#else 
    637644            CALL lbc_lnk( zsshu_a, 'U', 1._wp )   ;   CALL lbc_lnk( zsshv_a, 'V', 1._wp ) 
     645#endif 
    638646         ENDIF    
    639647         !                                  
     
    803811         !                                                 !  ----------------------- 
    804812         ! 
     813#if defined key_multisend 
     814         CALL lbc_lnk_multi( ua_e, 'U', -1._wp , va_e , 'V', -1._wp ) 
     815#else 
    805816         CALL lbc_lnk( ua_e  , 'U', -1._wp )               ! local domain boundaries  
    806817         CALL lbc_lnk( va_e  , 'V', -1._wp ) 
     818#endif 
     819 
     820 
    807821 
    808822#if defined key_bdy   
     
    859873            END DO 
    860874         END DO 
     875#if defined key_multisend 
     876         CALL lbc_lnk_multi( zsshu_a, 'U', 1._wp , zsshv_a, 'V', 1._wp ) ! Boundary conditions 
     877#else 
    861878         CALL lbc_lnk( zsshu_a, 'U', 1._wp )   ;   CALL lbc_lnk( zsshv_a, 'V', 1._wp ) ! Boundary conditions 
     879#endif 
    862880      ENDIF 
    863881      ! 
  • branches/2015/dev_r5302_CNRS18_HPC_scalability/NEMOGCM/NEMO/OPA_SRC/LBC/lbclnk.F90

    r4990 r5372  
    2222   USE lib_mpp          ! distributed memory computing library 
    2323 
     24 
     25#if defined key_multisend 
     26   INTERFACE lbc_lnk_multi 
     27      MODULE PROCEDURE mpp_lnk_2d_9 
     28   END INTERFACE 
     29#endif 
     30 
    2431   INTERFACE lbc_lnk 
    2532      MODULE PROCEDURE mpp_lnk_3d_gather, mpp_lnk_3d, mpp_lnk_2d 
     
    3946 
    4047   PUBLIC lbc_lnk       ! ocean lateral boundary conditions 
     48#if defined key_multisend 
     49   PUBLIC lbc_lnk_multi  ! modified ocean lateral boundary conditions 
     50#endif 
    4151   PUBLIC lbc_lnk_e 
    4252   PUBLIC lbc_bdy_lnk   ! ocean lateral BDY boundary conditions 
  • branches/2015/dev_r5302_CNRS18_HPC_scalability/NEMOGCM/NEMO/OPA_SRC/LBC/lib_mpp.F90

    r4990 r5372  
    7171   PUBLIC   mpp_min, mpp_max, mpp_sum, mpp_minloc, mpp_maxloc 
    7272   PUBLIC   mpp_lnk_3d, mpp_lnk_3d_gather, mpp_lnk_2d, mpp_lnk_2d_e 
     73   !MODIFICATION 
     74   PUBLIC   mpp_lnk_2d_9  
    7375   PUBLIC   mppscatter, mppgather 
    7476   PUBLIC   mpp_ini_ice, mpp_ini_znl 
     
    7779   PUBLIC   mpp_lnk_bdy_2d, mpp_lnk_bdy_3d 
    7880   PUBLIC   mpp_lbc_north_icb, mpp_lnk_2d_icb 
     81 
     82   type arrayptr 
     83      real , dimension (:,:),  pointer :: pt2d 
     84   end type arrayptr 
     85 
    7986 
    8087   !! * Interfaces 
     
    511518   END SUBROUTINE mpp_lnk_3d 
    512519 
     520   SUBROUTINE mpp_lnk_2d_multiple( pt2d_array , type_array , psgn_array , num_fields , cd_mpp, pval ) 
     521      !!---------------------------------------------------------------------- 
     522      !!                  ***  routine mpp_lnk_2d_multiple  *** 
     523      !! 
     524      !! ** Purpose :   Message passing management for multiple 2d arrays 
     525      !! 
     526      !! ** Method  :   Use mppsend and mpprecv function for passing mask 
     527      !!      between processors following neighboring subdomains. 
     528      !!            domain parameters 
     529      !!                    nlci   : first dimension of the local subdomain 
     530      !!                    nlcj   : second dimension of the local subdomain 
     531      !!                    nbondi : mark for "east-west local boundary" 
     532      !!                    nbondj : mark for "north-south local boundary" 
     533      !!                    noea   : number for local neighboring processors 
     534      !!                    nowe   : number for local neighboring processors 
     535      !!                    noso   : number for local neighboring processors 
     536      !!                    nono   : number for local neighboring processors 
     537      !! 
     538      !!---------------------------------------------------------------------- 
     539 
     540      INTEGER :: num_fields 
     541      TYPE( arrayptr ) , DIMENSION(:) :: pt2d_array 
     542      !REAL(wp), DIMENSION(jpi,jpj), INTENT(inout) ::   pt2d_array    ! Array of 2D arrays on which the boundary condition is applied 
     543      CHARACTER(len=1) , DIMENSION( : ), INTENT(in   ) ::   type_array ! define the nature of ptab array grid-points 
     544      !                                                         ! = T , U , V , F , W and I points 
     545      REAL(wp)        , DIMENSION( : ),  INTENT(in   ) ::   psgn_array    ! =-1 the sign change across the north fold boundary 
     546      !                                                         ! =  1. , the sign is kept 
     547      CHARACTER(len=3), OPTIONAL  , INTENT(in   ) ::   cd_mpp   ! fill the overlap area only 
     548      REAL(wp)        , OPTIONAL  , INTENT(in   ) ::   pval     ! background value (used at closed boundaries) 
     549      !! 
     550      INTEGER  ::   ji, jj, jl   ! dummy loop indices 
     551      INTEGER  ::   ii    !!MULTI SEND DUMMY LOOP INDICES 
     552      INTEGER  ::   imigr, iihom, ijhom        ! temporary integers 
     553      INTEGER  ::   ml_req1, ml_req2, ml_err   ! for key_mpi_isend 
     554 
     555      REAL(wp) ::   zland 
     556      INTEGER, DIMENSION(MPI_STATUS_SIZE) ::   ml_stat   ! for key_mpi_isend 
     557      ! 
     558      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  zt2ns, zt2sn   ! 2d for north-south & south-north 
     559      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  zt2ew, zt2we   ! 2d for east-west & west-east 
     560 
     561      !!---------------------------------------------------------------------- 
     562 
     563      ALLOCATE( zt2ns(jpi,jprecj,2*num_fields), zt2sn(jpi,jprecj,2*num_fields),  & 
     564         &      zt2ew(jpj,jpreci,2*num_fields), zt2we(jpj,jpreci,2*num_fields)   ) 
     565 
     566      ! 
     567      IF( PRESENT( pval ) ) THEN   ;   zland = pval      ! set land value 
     568      ELSE                         ;   zland = 0.e0      ! zero by default 
     569      ENDIF 
     570 
     571      ! 1. standard boundary treatment 
     572      ! ------------------------------ 
     573      ! 
     574      !First Array 
     575      DO ii = 1 , num_fields 
     576         IF( PRESENT( cd_mpp ) ) THEN      ! only fill added line/raw with existing values 
     577            ! 
     578            ! WARNING pt2d is defined only between nld and nle 
     579            DO jj = nlcj+1, jpj                 ! added line(s)   (inner only) 
     580 
     581 
     582               pt2d_array(ii)%pt2d(nldi  :nlei  , jj ) = pt2d_array(ii)%pt2d(nldi:nlei,     nlej ) 
     583 
     584 
     585               pt2d_array(ii)%pt2d(1     :nldi-1, jj ) = pt2d_array(ii)%pt2d(nldi     ,     nlej ) 
     586               pt2d_array(ii)%pt2d(nlei+1:nlci  , jj ) = pt2d_array(ii)%pt2d(     nlei,     nlej )  
     587            END DO 
     588            DO ji = nlci+1, jpi                 ! added column(s) (full) 
     589               pt2d_array(ii)%pt2d( ji ,nldj  :nlej  ) = pt2d_array(ii)%pt2d( nlei , nldj:nlej) 
     590               pt2d_array(ii)%pt2d( ji ,1     :nldj-1 ) = pt2d_array(ii)%pt2d( nlei , nldj ) 
     591               pt2d_array(ii)%pt2d( ji ,nlej+1:jpj ) = pt2d_array(ii)%pt2d( nlei , nlej ) 
     592            END DO 
     593            ! 
     594         ELSE                              ! standard close or cyclic treatment 
     595            ! 
     596            !                                   ! East-West boundaries 
     597            IF( nbondi == 2 .AND.   &                ! Cyclic east-west 
     598               &    (nperio == 1 .OR. nperio == 4 .OR. nperio == 6) ) THEN 
     599               pt2d_array(ii)%pt2d( 1 , : ) = pt2d_array(ii)%pt2d( jpim1 , : )                                    ! west 
     600               pt2d_array(ii)%pt2d( jpi , : ) = pt2d_array(ii)%pt2d(  2  , : )                                    ! east 
     601            ELSE                                     ! closed 
     602               IF( .NOT. type_array(ii) == 'F' )   pt2d_array(ii)%pt2d( 1 : jpreci , : ) = zland    ! south except F-point 
     603                                            pt2d_array(ii)%pt2d( nlci-jpreci+1 : jpi , : ) = zland    ! north 
     604            ENDIF 
     605            !                                   ! North-South boundaries (always closed) 
     606               IF( .NOT. type_array(ii) == 'F' )   pt2d_array(ii)%pt2d( : , 1 : jprecj ) = zland    !south except F-point 
     607                                            pt2d_array(ii)%pt2d( : , nlcj-jprecj+1:jpj ) = zland    ! north 
     608            ! 
     609         ENDIF 
     610      END DO 
     611 
     612 
     613      ! 2. East and west directions exchange 
     614      ! ------------------------------------ 
     615      ! we play with the neigbours AND the row number because of the periodicity 
     616      ! 
     617      DO ii = 1 , num_fields 
     618         SELECT CASE ( nbondi )      ! Read Dirichlet lateral conditions 
     619         CASE ( -1, 0, 1 )                ! all exept 2 (i.e. close case) 
     620            iihom = nlci-nreci 
     621            DO jl = 1, jpreci 
     622               zt2ew( : , jl , ii ) = pt2d_array(ii)%pt2d( jpreci+jl , : ) 
     623               zt2we( : , jl , ii ) = pt2d_array(ii)%pt2d( iihom +jl , : ) 
     624            END DO 
     625         END SELECT 
     626      END DO 
     627      ! 
     628      !                           ! Migrations 
     629      imigr = jpreci * jpj 
     630      ! 
     631      SELECT CASE ( nbondi ) 
     632      CASE ( -1 ) 
     633         CALL mppsend( 2, zt2we(1,1,1), num_fields*imigr, noea, ml_req1 ) 
     634         CALL mpprecv( 1, zt2ew(1,1,num_fields+1), num_fields*imigr, noea ) 
     635         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 
     636      CASE ( 0 ) 
     637         CALL mppsend( 1, zt2ew(1,1,1), num_fields*imigr, nowe, ml_req1 ) 
     638         CALL mppsend( 2, zt2we(1,1,1), num_fields*imigr, noea, ml_req2 ) 
     639         CALL mpprecv( 1, zt2ew(1,1,num_fields+1), num_fields*imigr, noea ) 
     640         CALL mpprecv( 2, zt2we(1,1,num_fields+1), num_fields*imigr, nowe ) 
     641         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 
     642         IF(l_isend) CALL mpi_wait(ml_req2,ml_stat,ml_err) 
     643      CASE ( 1 ) 
     644         CALL mppsend( 1, zt2ew(1,1,1), num_fields*imigr, nowe, ml_req1 ) 
     645         CALL mpprecv( 2, zt2we(1,1,num_fields+1), num_fields*imigr, nowe ) 
     646         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 
     647      END SELECT 
     648      ! 
     649      !                           ! Write Dirichlet lateral conditions 
     650      iihom = nlci - jpreci 
     651      ! 
     652 
     653   DO ii = 1 , num_fields 
     654      SELECT CASE ( nbondi ) 
     655      CASE ( -1 ) 
     656         DO jl = 1, jpreci 
     657            pt2d_array(ii)%pt2d( iihom+jl , : ) = zt2ew(:,jl,num_fields+ii) 
     658         END DO 
     659      CASE ( 0 ) 
     660         DO jl = 1, jpreci 
     661            pt2d_array(ii)%pt2d( jl , : ) = zt2we(:,jl,num_fields+ii) 
     662            pt2d_array(ii)%pt2d( iihom+jl , : ) = zt2ew(:,jl,num_fields+ii) 
     663         END DO 
     664      CASE ( 1 ) 
     665         DO jl = 1, jpreci 
     666            pt2d_array(ii)%pt2d( jl , : )= zt2we(:,jl,num_fields+ii) 
     667         END DO 
     668      END SELECT 
     669   END DO 
     670 
     671      ! 3. North and south directions 
     672      ! ----------------------------- 
     673      ! always closed : we play only with the neigbours 
     674      ! 
     675   !First Array 
     676   DO ii = 1 , num_fields 
     677      IF( nbondj /= 2 ) THEN      ! Read Dirichlet lateral conditions 
     678         ijhom = nlcj-nrecj 
     679         DO jl = 1, jprecj 
     680            zt2sn(:,jl , ii) = pt2d_array(ii)%pt2d( : , ijhom +jl ) 
     681            zt2ns(:,jl , ii) = pt2d_array(ii)%pt2d( : , jprecj+jl ) 
     682         END DO 
     683      ENDIF 
     684   END DO 
     685      ! 
     686      !                           ! Migrations 
     687      imigr = jprecj * jpi 
     688      ! 
     689      SELECT CASE ( nbondj ) 
     690      CASE ( -1 ) 
     691         CALL mppsend( 4, zt2sn(1,1,1), num_fields*imigr, nono, ml_req1 ) 
     692         CALL mpprecv( 3, zt2ns(1,1,num_fields+1), num_fields*imigr, nono ) 
     693         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 
     694      CASE ( 0 ) 
     695         CALL mppsend( 3, zt2ns(1,1,1), num_fields*imigr, noso, ml_req1 ) 
     696         CALL mppsend( 4, zt2sn(1,1,1), num_fields*imigr, nono, ml_req2 ) 
     697         CALL mpprecv( 3, zt2ns(1,1,num_fields+1), num_fields*imigr, nono ) 
     698         CALL mpprecv( 4, zt2sn(1,1,num_fields+1), num_fields*imigr, noso ) 
     699         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 
     700         IF(l_isend) CALL mpi_wait(ml_req2,ml_stat,ml_err) 
     701      CASE ( 1 ) 
     702         CALL mppsend( 3, zt2ns(1,1,1), num_fields*imigr, noso, ml_req1 ) 
     703         CALL mpprecv( 4, zt2sn(1,1,num_fields+1), num_fields*imigr, noso ) 
     704         IF(l_isend) CALL mpi_wait(ml_req1,ml_stat,ml_err) 
     705      END SELECT 
     706      ! 
     707      !                           ! Write Dirichlet lateral conditions 
     708      ijhom = nlcj - jprecj 
     709      ! 
     710 
     711   DO ii = 1 , num_fields 
     712   !First Array 
     713      SELECT CASE ( nbondj ) 
     714      CASE ( -1 ) 
     715         DO jl = 1, jprecj 
     716            pt2d_array(ii)%pt2d( : , ijhom+jl ) = zt2ns( : , jl , num_fields+ii ) 
     717         END DO 
     718      CASE ( 0 ) 
     719         DO jl = 1, jprecj 
     720            pt2d_array(ii)%pt2d( : , jl ) = zt2sn( : , jl , num_fields + ii) 
     721            pt2d_array(ii)%pt2d( : , ijhom + jl ) = zt2ns( : , jl , num_fields + ii ) 
     722         END DO 
     723      CASE ( 1 ) 
     724         DO jl = 1, jprecj 
     725            pt2d_array(ii)%pt2d( : , jl ) = zt2sn( : , jl , num_fields + ii ) 
     726         END DO 
     727      END SELECT 
     728   END DO 
     729 
     730      ! 4. north fold treatment 
     731      ! ----------------------- 
     732      ! 
     733 
     734 
     735   DO ii = 1 , num_fields 
     736   !First Array 
     737      IF( npolj /= 0 .AND. .NOT. PRESENT(cd_mpp) ) THEN 
     738         ! 
     739         SELECT CASE ( jpni ) 
     740         CASE ( 1 )     ;   CALL lbc_nfd      ( pt2d_array(ii)%pt2d( : , : ), type_array(ii) , psgn_array(ii) )   ! only 1 northern proc, no mpp 
     741         CASE DEFAULT   ;   CALL mpp_lbc_north( pt2d_array(ii)%pt2d( : , : ), type_array(ii), psgn_array(ii) )   ! for all northern procs. 
     742         END SELECT 
     743         ! 
     744      ENDIF 
     745      ! 
     746   END DO 
     747 
     748      DEALLOCATE( zt2ns, zt2sn, zt2ew, zt2we ) 
     749      ! 
     750   END SUBROUTINE mpp_lnk_2d_multiple 
     751 
     752   SUBROUTINE load_array(pt2d,cd_type,psgn,pt2d_array, type_array, psgn_array,num_fields) 
     753      REAL(wp), DIMENSION(jpi,jpj), TARGET   ,   INTENT(inout) ::   pt2d    ! Second 2D array on which the boundary condition is applied 
     754      CHARACTER(len=1)            , INTENT(in   ) ::   cd_type ! define the nature of ptab array grid-points 
     755      REAL(wp)         ,               INTENT(in   ) ::   psgn    ! =-1 the sign change across the north fold boundary 
     756      TYPE(arrayptr) , DIMENSION(9)               ::   pt2d_array 
     757      CHARACTER(len=1) , DIMENSION( 9 ) ::   type_array ! define the nature of ptab array grid-points 
     758      REAL(wp)        , DIMENSION( 9 ) ::   psgn_array    ! =-1 the sign change across the north fold boundary 
     759      INTEGER , INTENT (inout):: num_fields  
     760      num_fields=num_fields+1 
     761      pt2d_array(num_fields)%pt2d=>pt2d 
     762      type_array(num_fields)=cd_type 
     763      psgn_array(num_fields)=psgn 
     764   END SUBROUTINE load_array 
     765 
     766   SUBROUTINE mpp_lnk_2d_9( pt2dA, cd_typeA, psgnA, pt2dB, cd_typeB , psgnB, pt2dC, cd_typeC, psgnC, pt2dD, cd_typeD, psgnD, pt2dE, cd_typeE, psgnE, pt2dF, cd_typeF, psgnF, pt2dG, cd_typeG,  psgnG, pt2dH, cd_typeH, psgnH, pt2dI, cd_typeI, psgnI, cd_mpp, pval) 
     767      REAL(wp), DIMENSION(jpi,jpj), TARGET   ,   INTENT(inout) ::   pt2dA    ! Second 2D array on which the boundary condition is applied 
     768      REAL(wp), DIMENSION(jpi,jpj), TARGET   , OPTIONAL  ,  INTENT(inout) ::   pt2dB , pt2dC , pt2dD , pt2dE ,pt2dF , pt2dG , pt2dH , pt2dI  
     769      CHARACTER(len=1)            , INTENT(in   ) ::   cd_typeA ! define the nature of ptab array grid-points 
     770      CHARACTER(len=1)            , OPTIONAL  , INTENT(in   ) ::   cd_typeB , cd_typeC , cd_typeD , cd_typeE , cd_typeF , cd_typeG , cd_typeH , cd_typeI  
     771      REAL(wp)         ,               INTENT(in   ) ::   psgnA    ! =-1 the sign change across the north fold boundary 
     772      REAL(wp)         ,   OPTIONAL  ,            INTENT(in   ) ::   psgnB , psgnC , psgnD , psgnE , psgnF , psgnG , psgnH , psgnI    
     773      CHARACTER(len=3), OPTIONAL  , INTENT(in   ) ::   cd_mpp   ! fill the overlap area only 
     774      REAL(wp)        , OPTIONAL  , INTENT(in   ) ::   pval     ! background value (used at closed boundaries) 
     775      TYPE(arrayptr) , DIMENSION(9)               ::   pt2d_array  
     776      CHARACTER(len=1) , DIMENSION( 9 ) ::   type_array ! define the nature of ptab array grid-points 
     777      !                                                         ! = T , U , V , F , W and I points 
     778      REAL(wp)        , DIMENSION( 9 ) ::   psgn_array    ! =-1 the sign change across the north fold boundary 
     779      INTEGER :: num_fields 
     780       
     781      num_fields = 0 
     782 
     783      !! Load the first array 
     784      call load_array(pt2dA,cd_typeA,psgnA,pt2d_array, type_array, psgn_array,num_fields) 
     785 
     786      !! Look if more arrays are added 
     787      if(present (psgnB) )call load_array(pt2dB,cd_typeB,psgnB,pt2d_array, type_array, psgn_array,num_fields) 
     788      if(present (psgnC) )call load_array(pt2dC,cd_typeC,psgnC,pt2d_array, type_array, psgn_array,num_fields) 
     789      if(present (psgnD) )call load_array(pt2dD,cd_typeD,psgnD,pt2d_array, type_array, psgn_array,num_fields) 
     790      if(present (psgnE) )call load_array(pt2dE,cd_typeE,psgnE,pt2d_array, type_array, psgn_array,num_fields) 
     791      if(present (psgnF) )call load_array(pt2dF,cd_typeF,psgnF,pt2d_array, type_array, psgn_array,num_fields) 
     792      if(present (psgnG) )call load_array(pt2dG,cd_typeG,psgnG,pt2d_array, type_array, psgn_array,num_fields) 
     793      if(present (psgnH) )call load_array(pt2dH,cd_typeH,psgnH,pt2d_array, type_array, psgn_array,num_fields) 
     794      if(present (psgnI) )call load_array(pt2dI,cd_typeI,psgnI,pt2d_array, type_array, psgn_array,num_fields) 
     795          
     796      CALL mpp_lnk_2d_multiple(pt2d_array,type_array,psgn_array,num_fields,cd_mpp,pval) 
     797   END SUBROUTINE mpp_lnk_2d_9 
     798 
    513799 
    514800   SUBROUTINE mpp_lnk_2d( pt2d, cd_type, psgn, cd_mpp, pval ) 
Note: See TracChangeset for help on using the changeset viewer.