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 782 for trunk/NEMO/NST_SRC/agrif_opa_sponge.F90 – NEMO

Ignore:
Timestamp:
2008-01-07T10:07:26+01:00 (16 years ago)
Author:
rblod
Message:

Improvment of AGRIF-NEMO routines, see ticket #42

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMO/NST_SRC/agrif_opa_sponge.F90

    r719 r782  
    77   USE dom_oce 
    88   USE in_out_manager 
     9   USE agrif_oce 
    910 
    1011   IMPLICIT NONE 
     
    1314   PUBLIC Agrif_Sponge_Tra, Agrif_Sponge_Dyn, interptn, interpsn, interpun, interpvn 
    1415 
    15    !! * Namelist (namagrif) 
    16    REAL(wp) :: visc_tra = rdt 
    17    REAL(wp) :: visc_dyn = rdt 
    1816 
    1917   CONTAINS 
    2018 
    21    SUBROUTINE Agrif_Sponge_Tra( kt ) 
     19   SUBROUTINE Agrif_Sponge_Tra 
    2220      !!--------------------------------------------- 
    2321      !!   *** ROUTINE Agrif_Sponge_Tra *** 
     
    2523#include "domzgr_substitute.h90" 
    2624 
    27       INTEGER, INTENT(in) :: kt 
    28  
    2925      INTEGER :: ji,jj,jk 
    30       REAL(wp), DIMENSION(jpi,jpj,jpk) :: umasktemp,vmasktemp 
    3126      INTEGER :: spongearea 
    3227      REAL(wp) :: timecoeff 
    3328      REAL(wp) :: zta, zsa, zabe1, zabe2, zbtr 
    3429      REAL(wp), DIMENSION(jpi,jpj) :: localviscsponge 
    35       REAL(wp), DIMENSION(jpi,jpj,jpk) :: ztab, tbdiff, sbdiff 
     30      REAL(wp), DIMENSION(jpi,jpj,jpk) :: tbdiff, sbdiff 
    3631      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ztu ,ztv, zsu ,zsv 
     32      REAL(wp), DIMENSION(jpi,jpj,jpk) ::  ztab 
    3733 
    3834#if defined SPONGE 
    39  
    40       IF( kt == nit000  )   CALL agrif_sponge_init 
    4135 
    4236      timecoeff = REAL(Agrif_NbStepint(),wp)/Agrif_rhot() 
     
    5852      sbdiff(:,:,:) = sb(:,:,:) - ztab(:,:,:) 
    5953 
     54 
    6055      spongearea = 2 + 2 * Agrif_irhox() 
    6156 
    6257      localviscsponge = 0. 
    63       umasktemp = 0. 
    64       vmasktemp = 0. 
     58       
     59      IF (.NOT. spongedoneT) THEN 
     60         spe1ur(:,:) = 0. 
     61         spe2vr(:,:) = 0. 
    6562 
    6663      IF ((nbondi == -1).OR.(nbondi == 2)) THEN 
     
    6865            localviscsponge(ji,:) = visc_tra * (spongearea-ji)/real(spongearea-2) 
    6966         ENDDO 
    70          DO jk = 1, jpkm1 
    71             umasktemp(2:spongearea-1,:,jk) = umask(2:spongearea-1,:,jk) & 
    72                * 0.5 * (localviscsponge(2:spongearea-1,:) + localviscsponge(3:spongearea,:)) 
    73          ENDDO 
    74          DO jk = 1, jpkm1 
    75             vmasktemp(2:spongearea,1:jpjm1,jk) = vmask(2:spongearea,1:jpjm1,jk) & 
    76                * 0.5 * (localviscsponge(2:spongearea,1:jpjm1) + localviscsponge(2:spongearea,2:jpj)) 
    77          ENDDO 
     67     
     68    spe1ur(2:spongearea-1,:)=0.5 * (localviscsponge(2:spongearea-1,:) + localviscsponge(3:spongearea,:)) & 
     69          * e2u(2:spongearea-1,:) / e1u(2:spongearea-1,:) 
     70 
     71         spe2vr(2:spongearea,1:jpjm1) = 0.5 * (localviscsponge(2:spongearea,1:jpjm1) + & 
     72             localviscsponge(2:spongearea,2:jpj)) & 
     73           * e1v(2:spongearea,1:jpjm1) / e2v(2:spongearea,1:jpjm1) 
    7874      ENDIF 
    7975 
     
    8278            localviscsponge(ji,:) = visc_tra * (ji - (nlci-spongearea+1))/real(spongearea-2) 
    8379         ENDDO 
    84          DO jk = 1, jpkm1 
    85             umasktemp(nlci-spongearea + 1:nlci-2,:,jk) = umask(nlci-spongearea + 1:nlci-2,:,jk) & 
    86                * 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-2,:) + localviscsponge(nlci-spongearea + 2:nlci-1,:)) 
    87          ENDDO 
    88          DO jk = 1, jpkm1 
    89             vmasktemp(nlci-spongearea + 1:nlci-1,1:jpjm1,jk) = vmask(nlci-spongearea + 1:nlci-1,1:jpjm1,jk) & 
    90                * 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-1,1:jpjm1) + localviscsponge(nlci-spongearea + 1:nlci-1,2:jpj)) 
    91          ENDDO 
     80     
     81    spe1ur(nlci-spongearea + 1:nlci-2,:)=0.5 * (localviscsponge(nlci-spongearea + 1:nlci-2,:) + & 
     82           localviscsponge(nlci-spongearea + 2:nlci-1,:)) & 
     83          * e2u(nlci-spongearea + 1:nlci-2,:) / e1u(nlci-spongearea + 1:nlci-2,:) 
     84 
     85         spe2vr(nlci-spongearea + 1:nlci-1,1:jpjm1) = 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-1,1:jpjm1) & 
     86              + localviscsponge(nlci-spongearea + 1:nlci-1,2:jpj)) & 
     87           * e1v(nlci-spongearea + 1:nlci-1,1:jpjm1) / e2v(nlci-spongearea + 1:nlci-1,1:jpjm1) 
    9288      ENDIF 
    9389 
     
    9793            localviscsponge(:,jj) = visc_tra * (spongearea-jj)/real(spongearea-2) 
    9894         ENDDO 
    99          DO jk = 1, jpkm1 
    100             vmasktemp(:,2:spongearea-1,jk) = vmask(:,2:spongearea-1,jk) & 
    101                * 0.5 * (localviscsponge(:,2:spongearea-1) + localviscsponge(:,3:spongearea)) 
    102          ENDDO 
    103          DO jk = 1, jpkm1 
    104             umasktemp(1:jpim1,2:spongearea,jk) = umask(1:jpim1,2:spongearea,jk) & 
    105                * 0.5 * (localviscsponge(1:jpim1,2:spongearea) + localviscsponge(2:jpi,2:spongearea)) 
    106          ENDDO 
     95     
     96    spe1ur(1:jpim1,2:spongearea)=0.5 * (localviscsponge(1:jpim1,2:spongearea) + & 
     97           localviscsponge(2:jpi,2:spongearea)) & 
     98          * e2u(1:jpim1,2:spongearea) / e1u(1:jpim1,2:spongearea) 
     99 
     100         spe2vr(:,2:spongearea-1) = 0.5 * (localviscsponge(:,2:spongearea-1) + & 
     101             localviscsponge(:,3:spongearea)) & 
     102           * e1v(:,2:spongearea-1) / e2v(:,2:spongearea-1) 
    107103      ENDIF 
    108104 
     
    111107            localviscsponge(:,jj) = visc_tra * (jj - (nlcj-spongearea+1))/real(spongearea-2) 
    112108         ENDDO 
    113          DO jk = 1, jpkm1 
    114             vmasktemp(:,nlcj-spongearea + 1:nlcj-2,jk) = vmask(:,nlcj-spongearea + 1:nlcj-2,jk) & 
    115                * 0.5 * (localviscsponge(:,nlcj-spongearea + 1:nlcj-2) + localviscsponge(:,nlcj-spongearea + 2:nlcj-1)) 
    116          ENDDO 
    117          DO jk = 1, jpkm1 
    118             umasktemp(1:jpim1,nlcj-spongearea + 1:nlcj-1,jk) = umask(1:jpim1,nlcj-spongearea + 1:nlcj-1,jk) & 
    119                * 0.5 * (localviscsponge(1:jpim1,nlcj-spongearea + 1:nlcj-1) + localviscsponge(2:jpi,nlcj-spongearea + 1:nlcj-1)) 
    120          ENDDO 
    121       ENDIF 
    122  
    123       IF (.NOT. spongedoneT) THEN 
    124          zspe1ur(:,:) = e2u(:,:) / e1u(:,:) 
    125          zspe2vr(:,:) = e1v(:,:) / e2v(:,:) 
    126          zspbtr2(:,:) = 1. / ( e1t(:,:) * e2t(:,:)) 
     109     
     110    spe1ur(1:jpim1,nlcj-spongearea + 1:nlcj-1)=0.5 * (localviscsponge(1:jpim1,nlcj-spongearea + 1:nlcj-1) + & 
     111            localviscsponge(2:jpi,nlcj-spongearea + 1:nlcj-1)) & 
     112          * e2u(1:jpim1,nlcj-spongearea + 1:nlcj-1) / e1u(1:jpim1,nlcj-spongearea + 1:nlcj-1) 
     113 
     114         spe2vr(:,nlcj-spongearea + 1:nlcj-2) = 0.5 * (localviscsponge(:,nlcj-spongearea + 1:nlcj-2) + & 
     115            localviscsponge(:,nlcj-spongearea + 2:nlcj-1)) & 
     116           * e1v(:,nlcj-spongearea + 1:nlcj-2) / e2v(:,nlcj-spongearea + 1:nlcj-2) 
     117      ENDIF 
     118       
     119         spbtr2(:,:) = 1. / ( e1t(:,:) * e2t(:,:)) 
    127120 
    128121         spongedoneT = .TRUE. 
     
    133126            DO ji = 1, jpim1 
    134127#if defined key_zco 
    135                zabe1 = umasktemp(ji,jj,jk) * zspe1ur(ji,jj) 
    136                zabe2 = vmasktemp(ji,jj,jk) * zspe2vr(ji,jj) 
    137 #else 
    138                zabe1 = umasktemp(ji,jj,jk) * zspe1ur(ji,jj) * fse3u(ji,jj,jk) 
    139                zabe2 = vmasktemp(ji,jj,jk) * zspe2vr(ji,jj) * fse3v(ji,jj,jk) 
     128               zabe1 = umask(ji,jj,jk) * spe1ur(ji,jj) 
     129               zabe2 = vmask(ji,jj,jk) * spe2vr(ji,jj) 
     130#else 
     131               zabe1 = umask(ji,jj,jk) * spe1ur(ji,jj) * fse3u(ji,jj,jk) 
     132               zabe2 = vmask(ji,jj,jk) * spe2vr(ji,jj) * fse3v(ji,jj,jk) 
    140133#endif 
    141134               ztu(ji,jj,jk) = zabe1 * ( tbdiff(ji+1,jj  ,jk) - tbdiff(ji,jj,jk) ) 
     
    149142            DO ji = 2,jpim1 
    150143#if defined key_zco 
    151                zbtr = zspbtr2(ji,jj) 
    152 #else 
    153                zbtr = zspbtr2(ji,jj) / fse3t(ji,jj,jk) 
     144               zbtr = spbtr2(ji,jj) 
     145#else 
     146               zbtr = spbtr2(ji,jj) / fse3t(ji,jj,jk) 
    154147#endif 
    155148               ! horizontal diffusive trends 
     
    170163   END SUBROUTINE Agrif_Sponge_Tra 
    171164 
    172    SUBROUTINE Agrif_Sponge_dyn( kt ) 
     165   SUBROUTINE Agrif_Sponge_dyn 
    173166      !!--------------------------------------------- 
    174167      !!   *** ROUTINE Agrif_Sponge_dyn *** 
    175168      !!--------------------------------------------- 
    176169#include "domzgr_substitute.h90" 
    177  
    178       INTEGER,INTENT(in) :: kt 
    179170 
    180171      INTEGER :: ji,jj,jk 
    181172      INTEGER :: spongearea 
    182173      REAL(wp) :: timecoeff 
    183       REAL(wp) :: ze2u, ze1v, zua, zva 
     174      REAL(wp) :: ze2u, ze1v, zua, zva, zbtr 
    184175      REAL(wp), DIMENSION(jpi,jpj) :: localviscsponge 
    185176      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ztab, ubdiff, vbdiff,rotdiff,hdivdiff 
    186       REAL(wp), DIMENSION(jpi,jpj,jpk) :: umasktemp,vmasktemp 
    187177 
    188178#if defined SPONGE 
     
    191181 
    192182      Agrif_SpecialValue=0. 
    193       Agrif_UseSpecialValue = .TRUE. 
     183      Agrif_UseSpecialValue = ln_spc_dyn 
    194184      ztab = 0.e0 
    195185      CALL Agrif_Bc_Variable(ztab, ua,calledweight=timecoeff,procname=interpun) 
    196186      Agrif_UseSpecialValue = .FALSE. 
    197187 
    198       ubdiff(:,:,:) = ub(:,:,:) - ztab(:,:,:) 
     188      ubdiff(:,:,:) = (ub(:,:,:) - ztab(:,:,:))*umask(:,:,:) 
    199189 
    200190      ztab = 0.e0 
    201191      Agrif_SpecialValue=0. 
    202       Agrif_UseSpecialValue = .TRUE. 
     192      Agrif_UseSpecialValue = ln_spc_dyn 
    203193      CALL Agrif_Bc_Variable(ztab, va,calledweight=timecoeff,procname=interpvn) 
    204194      Agrif_UseSpecialValue = .FALSE. 
    205195 
    206       vbdiff(:,:,:) = vb(:,:,:) - ztab(:,:,:) 
     196      vbdiff(:,:,:) = (vb(:,:,:) - ztab(:,:,:))*vmask(:,:,:) 
    207197 
    208198      spongearea = 2 + 2 * Agrif_irhox() 
    209199 
    210200      localviscsponge = 0. 
    211       umasktemp = 0. 
    212       vmasktemp = 0. 
     201 
     202      IF (.NOT. spongedoneU) THEN 
     203         spe1ur2(:,:) = 0. 
     204         spe2vr2(:,:) = 0. 
    213205 
    214206      IF ((nbondi == -1).OR.(nbondi == 2)) THEN 
     
    216208            localviscsponge(ji,:) = visc_dyn * (spongearea-ji)/real(spongearea-2) 
    217209         ENDDO 
    218          DO jk = 1, jpkm1 
    219             umasktemp(2:spongearea-1,:,jk) = umask(2:spongearea-1,:,jk) & 
    220                * 0.5 * (localviscsponge(2:spongearea-1,:) + localviscsponge(3:spongearea,:)) 
    221          ENDDO 
    222          DO jk = 1, jpkm1 
    223             vmasktemp(2:spongearea,1:jpjm1,jk) = vmask(2:spongearea,1:jpjm1,jk) & 
    224                * 0.5 * (localviscsponge(2:spongearea,1:jpjm1) + localviscsponge(2:spongearea,2:jpj)) 
    225          ENDDO 
     210     
     211    spe1ur2(2:spongearea-1,:)=0.5 * (localviscsponge(2:spongearea-1,:) + localviscsponge(3:spongearea,:)) 
     212 
     213         spe2vr2(2:spongearea,1:jpjm1) = 0.5 * (localviscsponge(2:spongearea,1:jpjm1) + & 
     214             localviscsponge(2:spongearea,2:jpj)) 
    226215      ENDIF 
    227216 
     
    230219            localviscsponge(ji,:) = visc_dyn * (ji - (nlci-spongearea+1))/real(spongearea-2) 
    231220         ENDDO 
    232          DO jk = 1, jpkm1 
    233             umasktemp(nlci-spongearea + 1:nlci-2,:,jk) = umask(nlci-spongearea + 1:nlci-2,:,jk) & 
    234                * 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-2,:) + localviscsponge(nlci-spongearea + 2:nlci-1,:)) 
    235          ENDDO 
    236          DO jk = 1, jpkm1 
    237             vmasktemp(nlci-spongearea + 1:nlci-1,1:jpjm1,jk) = vmask(nlci-spongearea + 1:nlci-1,1:jpjm1,jk) & 
    238                * 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-1,1:jpjm1) + localviscsponge(nlci-spongearea + 1:nlci-1,2:jpj)) 
    239          ENDDO 
    240       ENDIF 
     221     
     222    spe1ur2(nlci-spongearea + 1:nlci-2,:)=0.5 * (localviscsponge(nlci-spongearea + 1:nlci-2,:) + & 
     223           localviscsponge(nlci-spongearea + 2:nlci-1,:)) 
     224 
     225         spe2vr2(nlci-spongearea + 1:nlci-1,1:jpjm1) = 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-1,1:jpjm1) & 
     226              + localviscsponge(nlci-spongearea + 1:nlci-1,2:jpj)) 
     227      ENDIF 
     228 
    241229 
    242230      IF ((nbondj == -1).OR.(nbondj == 2)) THEN 
     
    244232            localviscsponge(:,jj) = visc_dyn * (spongearea-jj)/real(spongearea-2) 
    245233         ENDDO 
    246          DO jk = 1, jpkm1 
    247             vmasktemp(:,2:spongearea-1,jk) = vmask(:,2:spongearea-1,jk) & 
    248                * 0.5 * (localviscsponge(:,2:spongearea-1) + localviscsponge(:,3:spongearea)) 
    249          ENDDO 
    250          DO jk = 1, jpkm1 
    251             umasktemp(1:jpim1,2:spongearea,jk) = umask(1:jpim1,2:spongearea,jk) & 
    252                * 0.5 * (localviscsponge(1:jpim1,2:spongearea) + localviscsponge(2:jpi,2:spongearea)) 
    253          ENDDO 
     234     
     235    spe1ur2(1:jpim1,2:spongearea)=0.5 * (localviscsponge(1:jpim1,2:spongearea) + & 
     236           localviscsponge(2:jpi,2:spongearea)) 
     237 
     238         spe2vr2(:,2:spongearea-1) = 0.5 * (localviscsponge(:,2:spongearea-1) + & 
     239             localviscsponge(:,3:spongearea)) 
    254240      ENDIF 
    255241 
     
    258244            localviscsponge(:,jj) = visc_dyn * (jj - (nlcj-spongearea+1))/real(spongearea-2) 
    259245         ENDDO 
    260          DO jk = 1, jpkm1 
    261             vmasktemp(:,nlcj-spongearea + 1:nlcj-2,jk) = vmask(:,nlcj-spongearea + 1:nlcj-2,jk) & 
    262                * 0.5 * (localviscsponge(:,nlcj-spongearea + 1:nlcj-2) + localviscsponge(:,nlcj-spongearea + 2:nlcj-1)) 
    263          ENDDO 
    264          DO jk = 1, jpkm1 
    265             umasktemp(1:jpim1,nlcj-spongearea + 1:nlcj-1,jk) = umask(1:jpim1,nlcj-spongearea + 1:nlcj-1,jk) & 
    266                * 0.5 * (localviscsponge(1:jpim1,nlcj-spongearea + 1:nlcj-1) + localviscsponge(2:jpi,nlcj-spongearea + 1:nlcj-1)) 
    267          ENDDO 
    268       ENDIF 
    269  
    270       ubdiff = ubdiff * umasktemp 
    271       vbdiff = vbdiff * vmasktemp 
    272  
     246     
     247    spe1ur2(1:jpim1,nlcj-spongearea + 1:nlcj-1)=0.5 * (localviscsponge(1:jpim1,nlcj-spongearea + 1:nlcj-1) + & 
     248            localviscsponge(2:jpi,nlcj-spongearea + 1:nlcj-1)) 
     249 
     250         spe2vr2(:,nlcj-spongearea + 1:nlcj-2) = 0.5 * (localviscsponge(:,nlcj-spongearea + 1:nlcj-2) + & 
     251            localviscsponge(:,nlcj-spongearea + 2:nlcj-1)) 
     252      ENDIF 
     253 
     254         spongedoneU = .TRUE. 
     255     
     256    spbtr3(:,:) = 1./( e1f(:,:) * e2f(:,:)) 
     257      ENDIF 
     258       
     259      IF (.NOT. spongedoneT) THEN 
     260        spbtr2(:,:) = 1. / ( e1t(:,:) * e2t(:,:))       
     261      ENDIF 
     262       
     263      DO jk=1,jpkm1 
     264      ubdiff(:,:,jk) = ubdiff(:,:,jk) * spe1ur2(:,:) 
     265      vbdiff(:,:,jk) = vbdiff(:,:,jk) * spe2vr2(:,:) 
     266      ENDDO 
     267       
    273268      hdivdiff = 0. 
    274269      rotdiff = 0. 
     
    283278            DO ji = 2, jpim1   ! vector opt. 
    284279#if defined key_zco 
     280               zbtr = spbtr2(ji,jj) 
    285281               hdivdiff(ji,jj,jk) = (  e2u(ji,jj) * ubdiff(ji,jj,jk) & 
    286282                  - e2u(ji-1,jj  ) * ubdiff(ji-1,jj  ,jk)      & 
    287283                  &               + e1v(ji,jj) * vbdiff(ji,jj,jk) - & 
    288                   &              e1v(ji  ,jj-1) * vbdiff(ji  ,jj-1,jk)  )   & 
    289                   &            / ( e1t(ji,jj) * e2t(ji,jj) ) 
    290 #else 
     284                  &              e1v(ji  ,jj-1) * vbdiff(ji  ,jj-1,jk)  ) * zbtr 
     285#else 
     286               zbtr = spbtr2(ji,jj) / fse3t(ji,jj,jk) 
    291287               hdivdiff(ji,jj,jk) =   & 
    292288                  (  e2u(ji,jj)*fse3u(ji,jj,jk) * &  
     
    295291                  + e1v(ji,jj)*fse3v(ji,jj,jk) * & 
    296292                  vbdiff(ji,jj,jk) - e1v(ji  ,jj-1)* & 
    297                   fse3v(ji  ,jj-1,jk)  * vbdiff(ji  ,jj-1,jk)  )    & 
    298                   / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) 
     293                  fse3v(ji  ,jj-1,jk)  * vbdiff(ji  ,jj-1,jk)  ) * zbtr 
    299294#endif 
    300295            END DO 
     
    303298         DO jj = 1, jpjm1 
    304299            DO ji = 1, jpim1   ! vector opt. 
     300#if defined key_zco       
     301               zbtr = spbtr3(ji,jj) 
    305302               rotdiff(ji,jj,jk) = (  e2v(ji+1,jj  ) * vbdiff(ji+1,jj  ,jk) - e2v(ji,jj) * vbdiff(ji,jj,jk)    & 
    306303                  &              - e1u(ji  ,jj+1) * ubdiff(ji  ,jj+1,jk) + e1u(ji,jj) * ubdiff(ji,jj,jk)  ) & 
    307                   &           * fmask(ji,jj,jk) / ( e1f(ji,jj) * e2f(ji,jj) ) 
     304                  &           * fmask(ji,jj,jk) * zbtr 
     305#else 
     306               zbtr = spbtr3(ji,jj) * fse3f(ji,jj,jk) 
     307               rotdiff(ji,jj,jk) = (  e2v(ji+1,jj  ) * vbdiff(ji+1,jj  ,jk) - e2v(ji,jj) * vbdiff(ji,jj,jk)    & 
     308                  &              - e1u(ji  ,jj+1) * ubdiff(ji  ,jj+1,jk) + e1u(ji,jj) * ubdiff(ji,jj,jk)  ) & 
     309                  &           * fmask(ji,jj,jk) * zbtr 
     310#endif         
    308311            END DO 
    309312         END DO 
     
    330333                  ze1v                  ) / e2v(ji,jj) 
    331334#else 
    332                ze2u = rotdiff (ji,jj,jk)*fse3f(ji,jj,jk) 
     335               ze2u = rotdiff (ji,jj,jk) 
    333336               ze1v = hdivdiff(ji,jj,jk) 
    334337               ! horizontal diffusive trends 
    335                zua = - ( ze2u - rotdiff (ji,jj-1,jk)* & 
    336                   fse3f(ji,jj-1,jk) ) / ( e2u(ji,jj) * fse3u(ji,jj,jk) )   & 
     338               zua = - ( ze2u - rotdiff (ji,jj-1,jk)) / ( e2u(ji,jj) * fse3u(ji,jj,jk) )   & 
    337339                  + ( hdivdiff(ji+1,jj,jk) - ze1v      & 
    338340                  ) / e1u(ji,jj) 
    339341 
    340                zva = + ( ze2u - rotdiff (ji-1,jj,jk)* & 
    341                   fse3f(ji-1,jj,jk) ) / ( e1v(ji,jj) * fse3v(ji,jj,jk) )   & 
     342               zva = + ( ze2u - rotdiff (ji-1,jj,jk)) / ( e1v(ji,jj) * fse3v(ji,jj,jk) )   & 
    342343                  + ( hdivdiff(ji,jj+1,jk) - ze1v    & 
    343344                  ) / e2v(ji,jj) 
     
    357358   END SUBROUTINE Agrif_Sponge_dyn 
    358359 
    359    SUBROUTINE agrif_sponge_init 
    360       !!--------------------------------------------- 
    361       !!   *** ROUTINE agrif_sponge_init *** 
    362       !!--------------------------------------------- 
    363       NAMELIST/namagrif/ visc_tra, visc_dyn 
    364       REWIND ( numnam ) 
    365       READ   ( numnam, namagrif ) 
    366  
    367       IF(lwp) THEN 
    368          WRITE(numout,*) 
    369          WRITE(numout,*) 'agrif_sponge_init : agrif sponge parameters' 
    370          WRITE(numout,*) '~~~~~~~~~~~~' 
    371          WRITE(numout,*) '          Namelist namagrif : set sponge parameters' 
    372          WRITE(numout,*) '             sponge coefficient for tracers =  ', visc_tra 
    373          WRITE(numout,*) '             sponge coefficient for dynamics = ', visc_dyn 
    374       ENDIF 
    375  
    376    END SUBROUTINE agrif_sponge_init 
    377  
    378360   SUBROUTINE interptn(tabres,i1,i2,j1,j2,k1,k2) 
    379361      !!--------------------------------------------- 
     
    402384   END SUBROUTINE interpsn 
    403385 
    404  
    405386   SUBROUTINE interpun(tabres,i1,i2,j1,j2,k1,k2) 
    406387      !!--------------------------------------------- 
Note: See TracChangeset for help on using the changeset viewer.