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 5429 – NEMO

Changeset 5429


Ignore:
Timestamp:
2015-06-16T11:57:07+02:00 (9 years ago)
Author:
smasson
Message:

merge dev_r5302_CNRS18_HPC_scalability into the trunk, see #1523

Location:
trunk/NEMOGCM
Files:
8 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMOGCM/CONFIG/SHARED/namelist_ice_lim3_ref

    r5344 r5429  
    7676/ 
    7777!------------------------------------------------------------------------------ 
     78&namicehdf     !   Ice horizontal diffusion 
     79!------------------------------------------------------------------------------ 
     80   nn_convfrq     = 5              !  convergence check frequency of the Crant-Nicholson scheme (perf. optimization) 
     81/ 
     82!------------------------------------------------------------------------------ 
    7883&namicethd     !   Ice thermodynamics 
    7984!------------------------------------------------------------------------------ 
  • trunk/NEMOGCM/NEMO/LIM_SRC_3/limadv.F90

    r5123 r5429  
    207207 
    208208      !-- Lateral boundary conditions 
    209       CALL lbc_lnk( psm , 'T',  1. )   ;   CALL lbc_lnk( ps0 , 'T',  1. ) 
    210       CALL lbc_lnk( psx , 'T', -1. )   ;   CALL lbc_lnk( psy , 'T', -1. )      ! caution gradient ==> the sign changes 
    211       CALL lbc_lnk( psxx, 'T',  1. )   ;   CALL lbc_lnk( psyy, 'T',  1. ) 
    212       CALL lbc_lnk( psxy, 'T',  1. ) 
     209      CALL lbc_lnk_multi( psm , 'T',  1., ps0 , 'T',  1.   & 
     210         &              , psx , 'T', -1., psy , 'T', -1.   &   ! caution gradient ==> the sign changes 
     211         &              , psxx, 'T',  1., psyy, 'T',  1.   & 
     212         &              , psxy, 'T',  1. ) 
    213213 
    214214      IF(ln_ctl) THEN 
     
    393393 
    394394      !-- Lateral boundary conditions 
    395       CALL lbc_lnk( psm , 'T',  1. )   ;   CALL lbc_lnk( ps0 , 'T',  1. ) 
    396       CALL lbc_lnk( psx , 'T', -1. )   ;   CALL lbc_lnk( psy , 'T', -1. )      ! caution gradient ==> the sign changes 
    397       CALL lbc_lnk( psxx, 'T',  1. )   ;   CALL lbc_lnk( psyy, 'T',  1. ) 
    398       CALL lbc_lnk( psxy, 'T',  1. ) 
     395      CALL lbc_lnk_multi( psm , 'T',  1.,  ps0 , 'T',  1.   & 
     396         &              , psx , 'T', -1.,  psy , 'T', -1.   &   ! caution gradient ==> the sign changes 
     397         &              , psxx, 'T',  1.,  psyy, 'T',  1.   & 
     398         &              , psxy, 'T',  1. ) 
    399399 
    400400      IF(ln_ctl) THEN 
  • trunk/NEMOGCM/NEMO/LIM_SRC_3/limhdf.F90

    r5123 r5429  
    1313   !!---------------------------------------------------------------------- 
    1414   !!   lim_hdf       : diffusion trend on sea-ice variable 
     15   !!   lim_hdf_init  : initialisation of diffusion trend on sea-ice variable 
    1516   !!---------------------------------------------------------------------- 
    1617   USE dom_oce        ! ocean domain 
     
    2627   PRIVATE 
    2728 
    28    PUBLIC   lim_hdf     ! called by lim_trp 
     29   PUBLIC   lim_hdf         ! called by lim_trp 
     30   PUBLIC   lim_hdf_init    ! called by sbc_lim_init 
    2931 
    3032   LOGICAL  ::   linit = .TRUE.                             ! initialization flag (set to flase after the 1st call) 
     33   INTEGER  ::   nn_convfrq                                 !:  convergence check frequency of the Crant-Nicholson scheme 
    3134   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   efact   ! metric coefficient 
    3235 
     
    121124         CALL lbc_lnk( zrlx, 'T', 1. )                   ! lateral boundary condition 
    122125         ! 
    123          zconv = 0._wp                                   ! convergence test 
    124          DO jj = 2, jpjm1 
    125             DO ji = fs_2, fs_jpim1 
    126                zconv = MAX( zconv, ABS( zrlx(ji,jj) - ptab(ji,jj) )  ) 
    127             END DO 
    128          END DO 
    129          IF( lk_mpp )   CALL mpp_max( zconv )            ! max over the global domain 
     126         IF ( MOD( iter, nn_convfrq ) == 0 )  THEN    ! convergence test every nn_convfrq iterations (perf. optimization) 
     127            zconv = 0._wp 
     128            DO jj = 2, jpjm1 
     129               DO ji = fs_2, fs_jpim1 
     130                  zconv = MAX( zconv, ABS( zrlx(ji,jj) - ptab(ji,jj) )  ) 
     131               END DO 
     132            END DO 
     133            IF( lk_mpp )   CALL mpp_max( zconv )      ! max over the global domain 
     134         ENDIF 
    130135         ! 
    131136         ptab(:,:) = zrlx(:,:) 
     
    162167   END SUBROUTINE lim_hdf 
    163168 
     169    
     170   SUBROUTINE lim_hdf_init 
     171      !!------------------------------------------------------------------- 
     172      !!                  ***  ROUTINE lim_hdf_init  *** 
     173      !! 
     174      !! ** Purpose : Initialisation of horizontal diffusion of sea-ice  
     175      !! 
     176      !! ** Method  : Read the namicehdf namelist 
     177      !! 
     178      !! ** input   : Namelist namicehdf 
     179      !!------------------------------------------------------------------- 
     180      INTEGER  ::   ios                 ! Local integer output status for namelist read 
     181      NAMELIST/namicehdf/ nn_convfrq 
     182      !!------------------------------------------------------------------- 
     183      ! 
     184      IF(lwp) THEN 
     185         WRITE(numout,*) 
     186         WRITE(numout,*) 'lim_hdf : Ice horizontal diffusion' 
     187         WRITE(numout,*) '~~~~~~~' 
     188      ENDIF 
     189      ! 
     190      REWIND( numnam_ice_ref )              ! Namelist namicehdf in reference namelist : Ice horizontal diffusion 
     191      READ  ( numnam_ice_ref, namicehdf, IOSTAT = ios, ERR = 901) 
     192901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicehdf in reference namelist', lwp ) 
     193 
     194      REWIND( numnam_ice_cfg )              ! Namelist namicehdf in configuration namelist : Ice horizontal diffusion 
     195      READ  ( numnam_ice_cfg, namicehdf, IOSTAT = ios, ERR = 902 ) 
     196902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicehdf in configuration namelist', lwp ) 
     197      IF(lwm) WRITE ( numoni, namicehdf ) 
     198      ! 
     199      IF(lwp) THEN                          ! control print 
     200         WRITE(numout,*) 
     201         WRITE(numout,*)'   Namelist of ice parameters for ice horizontal diffusion computation ' 
     202         WRITE(numout,*)'      convergence check frequency of the Crant-Nicholson scheme   nn_convfrq   = ', nn_convfrq 
     203      ENDIF 
     204      ! 
     205   END SUBROUTINE lim_hdf_init 
    164206#else 
    165207   !!---------------------------------------------------------------------- 
    166208   !!   Default option          Dummy module           NO LIM sea-ice model 
    167209   !!---------------------------------------------------------------------- 
    168 CONTAINS 
    169    SUBROUTINE lim_hdf         ! Empty routine 
    170    END SUBROUTINE lim_hdf 
    171210#endif 
    172211 
  • trunk/NEMOGCM/NEMO/LIM_SRC_3/limrhg.F90

    r5123 r5429  
    377377            END DO 
    378378         END DO 
    379          CALL lbc_lnk( v_ice1  , 'U', -1. )   ;   CALL lbc_lnk( u_ice2  , 'V', -1. )      ! lateral boundary cond. 
    380  
     379 
     380         CALL lbc_lnk_multi( v_ice1, 'U', -1., u_ice2, 'V', -1. )      ! lateral boundary cond. 
     381          
    381382         DO jj = k_j1+1, k_jpj-1 
    382383            DO ji = fs_2, fs_jpim1 
     
    412413            END DO 
    413414         END DO 
    414          CALL lbc_lnk( zs1 , 'T', 1. )   ;   CALL lbc_lnk( zs2, 'T', 1. ) 
    415          CALL lbc_lnk( zs12, 'F', 1. ) 
    416  
     415 
     416         CALL lbc_lnk_multi( zs1 , 'T', 1., zs2, 'T', 1., zs12, 'F', 1. ) 
     417  
    417418         ! Ice internal stresses (Appendix C of Hunke and Dukowicz, 2002) 
    418419         DO jj = k_j1+1, k_jpj-1 
     
    570571      END DO 
    571572 
    572       CALL lbc_lnk( u_ice(:,:), 'U', -1. )  
    573       CALL lbc_lnk( v_ice(:,:), 'V', -1. )  
     573      CALL lbc_lnk_multi( u_ice(:,:), 'U', -1., v_ice(:,:), 'V', -1. ) 
     574 
    574575#if defined key_agrif && defined key_lim2 
    575576      CALL agrif_rhg_lim2( nn_nevp , nn_nevp, 'U' ) 
     
    595596      END DO 
    596597 
    597       CALL lbc_lnk( u_ice2(:,:), 'V', -1. )  
    598       CALL lbc_lnk( v_ice1(:,:), 'U', -1. ) 
     598      CALL lbc_lnk_multi( u_ice2(:,:), 'V', -1., v_ice1(:,:), 'U', -1. ) 
    599599 
    600600      ! Recompute delta, shear and div, inputs for mechanical redistribution  
     
    643643 
    644644      ! Lateral boundary condition 
    645       CALL lbc_lnk( divu_i (:,:), 'T', 1. ) 
    646       CALL lbc_lnk( delta_i(:,:), 'T', 1. ) 
    647       ! CALL lbc_lnk( shear_i(:,:), 'F', 1. ) 
    648       CALL lbc_lnk( shear_i(:,:), 'T', 1. ) 
     645      CALL lbc_lnk_multi( divu_i (:,:), 'T', 1., delta_i(:,:), 'T', 1.,  shear_i(:,:), 'T', 1. ) 
    649646 
    650647      ! * Store the stress tensor for the next time step 
  • trunk/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg_ts.F90

    r5217 r5429  
    557557               END DO 
    558558            END DO 
    559             CALL lbc_lnk( zwx, 'U', 1._wp )    ;   CALL lbc_lnk( zwy, 'V', 1._wp ) 
     559            CALL lbc_lnk_multi( zwx, 'U', 1._wp, zwy, 'V', 1._wp ) 
    560560            ! 
    561561            zhup2_e (:,:) = hu_0(:,:) + zwx(:,:)                ! Ocean depth at U- and V-points 
     
    635635               END DO 
    636636            END DO 
    637             CALL lbc_lnk( zsshu_a, 'U', 1._wp )   ;   CALL lbc_lnk( zsshv_a, 'V', 1._wp ) 
     637            CALL lbc_lnk_multi( zsshu_a, 'U', 1._wp, zsshv_a, 'V', 1._wp ) 
    638638         ENDIF    
    639639         !                                  
     
    803803         !                                                 !  ----------------------- 
    804804         ! 
    805          CALL lbc_lnk( ua_e  , 'U', -1._wp )               ! local domain boundaries  
    806          CALL lbc_lnk( va_e  , 'V', -1._wp ) 
     805         CALL lbc_lnk_multi( ua_e, 'U', -1._wp, va_e , 'V', -1._wp ) 
    807806 
    808807#if defined key_bdy   
     
    859858            END DO 
    860859         END DO 
    861          CALL lbc_lnk( zsshu_a, 'U', 1._wp )   ;   CALL lbc_lnk( zsshv_a, 'V', 1._wp ) ! Boundary conditions 
     860         CALL lbc_lnk_multi( zsshu_a, 'U', 1._wp, zsshv_a, 'V', 1._wp ) ! Boundary conditions 
    862861      ENDIF 
    863862      ! 
  • trunk/NEMOGCM/NEMO/OPA_SRC/LBC/lbclnk.F90

    r4990 r5429  
    2222   USE lib_mpp          ! distributed memory computing library 
    2323 
     24 
     25   INTERFACE lbc_lnk_multi 
     26      MODULE PROCEDURE mpp_lnk_2d_9 
     27   END INTERFACE 
     28 
    2429   INTERFACE lbc_lnk 
    2530      MODULE PROCEDURE mpp_lnk_3d_gather, mpp_lnk_3d, mpp_lnk_2d 
     
    3944 
    4045   PUBLIC lbc_lnk       ! ocean lateral boundary conditions 
     46   PUBLIC lbc_lnk_multi ! modified ocean lateral boundary conditions 
    4147   PUBLIC lbc_lnk_e 
    4248   PUBLIC lbc_bdy_lnk   ! ocean lateral BDY boundary conditions 
  • trunk/NEMOGCM/NEMO/OPA_SRC/LBC/lib_mpp.F90

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

    r5410 r5429  
    3737   USE limdyn          ! Ice dynamics 
    3838   USE limtrp          ! Ice transport 
     39   USE limhdf          ! Ice horizontal diffusion 
    3940   USE limthd          ! Ice thermodynamics 
    4041   USE limitd_me       ! Mechanics on ice thickness distribution 
     
    293294      CALL lim_itd_init                ! ice thickness distribution initialization 
    294295      ! 
     296      CALL lim_hdf_init                ! set ice horizontal diffusion computation parameters 
     297      ! 
    295298      CALL lim_thd_init                ! set ice thermodynics parameters 
    296299      ! 
Note: See TracChangeset for help on using the changeset viewer.