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 811 for branches/dev_001_SBC/NEMO/NST_SRC/agrif_opa_sponge.F90 – NEMO

Ignore:
Timestamp:
2008-02-07T17:00:12+01:00 (16 years ago)
Author:
ctlod
Message:

dev_001_SBC: merge with the trunk last changesets: #780, 782, 783, 784, 785, 788, 789, 793, 794

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/dev_001_SBC/NEMO/NST_SRC/agrif_opa_sponge.F90

    r699 r811  
    1010   USE dom_oce 
    1111   USE in_out_manager 
     12   USE agrif_oce 
    1213 
    1314   IMPLICIT NONE 
     
    1617   PUBLIC Agrif_Sponge_Tra, Agrif_Sponge_Dyn, interptn, interpsn, interpun, interpvn 
    1718 
    18    !! * Namelist (namagrif) 
    19    REAL(wp) :: visc_tra = rdt 
    20    REAL(wp) :: visc_dyn = rdt 
    2119 
    2220   CONTAINS 
    2321 
    24    SUBROUTINE Agrif_Sponge_Tra( kt ) 
     22   SUBROUTINE Agrif_Sponge_Tra 
    2523      !!--------------------------------------------- 
    2624      !!   *** ROUTINE Agrif_Sponge_Tra *** 
     
    2826#include "domzgr_substitute.h90" 
    2927 
    30       INTEGER, INTENT(in) :: kt 
    31  
    3228      INTEGER :: ji,jj,jk 
    33       REAL(wp), DIMENSION(jpi,jpj,jpk) :: umasktemp,vmasktemp 
    3429      INTEGER :: spongearea 
    3530      REAL(wp) :: timecoeff 
    3631      REAL(wp) :: zta, zsa, zabe1, zabe2, zbtr 
    3732      REAL(wp), DIMENSION(jpi,jpj) :: localviscsponge 
    38       REAL(wp), DIMENSION(jpi,jpj,jpk) :: ztab, tbdiff, sbdiff 
     33      REAL(wp), DIMENSION(jpi,jpj,jpk) :: tbdiff, sbdiff 
    3934      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ztu ,ztv, zsu ,zsv 
     35      REAL(wp), DIMENSION(jpi,jpj,jpk) ::  ztab 
    4036 
    4137#if defined SPONGE 
    42  
    43       IF( kt == nit000  )   CALL agrif_sponge_init 
    4438 
    4539      timecoeff = REAL(Agrif_NbStepint(),wp)/Agrif_rhot() 
     
    6155      sbdiff(:,:,:) = sb(:,:,:) - ztab(:,:,:) 
    6256 
     57 
    6358      spongearea = 2 + 2 * Agrif_irhox() 
    6459 
    6560      localviscsponge = 0. 
    66       umasktemp = 0. 
    67       vmasktemp = 0. 
     61       
     62      IF (.NOT. spongedoneT) THEN 
     63         spe1ur(:,:) = 0. 
     64         spe2vr(:,:) = 0. 
    6865 
    6966      IF ((nbondi == -1).OR.(nbondi == 2)) THEN 
     
    7168            localviscsponge(ji,:) = visc_tra * (spongearea-ji)/real(spongearea-2) 
    7269         ENDDO 
    73          DO jk = 1, jpkm1 
    74             umasktemp(2:spongearea-1,:,jk) = umask(2:spongearea-1,:,jk) & 
    75                * 0.5 * (localviscsponge(2:spongearea-1,:) + localviscsponge(3:spongearea,:)) 
    76          ENDDO 
    77          DO jk = 1, jpkm1 
    78             vmasktemp(2:spongearea,1:jpjm1,jk) = vmask(2:spongearea,1:jpjm1,jk) & 
    79                * 0.5 * (localviscsponge(2:spongearea,1:jpjm1) + localviscsponge(2:spongearea,2:jpj)) 
    80          ENDDO 
     70     
     71    spe1ur(2:spongearea-1,:)=0.5 * (localviscsponge(2:spongearea-1,:) + localviscsponge(3:spongearea,:)) & 
     72          * e2u(2:spongearea-1,:) / e1u(2:spongearea-1,:) 
     73 
     74         spe2vr(2:spongearea,1:jpjm1) = 0.5 * (localviscsponge(2:spongearea,1:jpjm1) + & 
     75             localviscsponge(2:spongearea,2:jpj)) & 
     76           * e1v(2:spongearea,1:jpjm1) / e2v(2:spongearea,1:jpjm1) 
    8177      ENDIF 
    8278 
     
    8581            localviscsponge(ji,:) = visc_tra * (ji - (nlci-spongearea+1))/real(spongearea-2) 
    8682         ENDDO 
    87          DO jk = 1, jpkm1 
    88             umasktemp(nlci-spongearea + 1:nlci-2,:,jk) = umask(nlci-spongearea + 1:nlci-2,:,jk) & 
    89                * 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-2,:) + localviscsponge(nlci-spongearea + 2:nlci-1,:)) 
    90          ENDDO 
    91          DO jk = 1, jpkm1 
    92             vmasktemp(nlci-spongearea + 1:nlci-1,1:jpjm1,jk) = vmask(nlci-spongearea + 1:nlci-1,1:jpjm1,jk) & 
    93                * 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-1,1:jpjm1) + localviscsponge(nlci-spongearea + 1:nlci-1,2:jpj)) 
    94          ENDDO 
     83     
     84    spe1ur(nlci-spongearea + 1:nlci-2,:)=0.5 * (localviscsponge(nlci-spongearea + 1:nlci-2,:) + & 
     85           localviscsponge(nlci-spongearea + 2:nlci-1,:)) & 
     86          * e2u(nlci-spongearea + 1:nlci-2,:) / e1u(nlci-spongearea + 1:nlci-2,:) 
     87 
     88         spe2vr(nlci-spongearea + 1:nlci-1,1:jpjm1) = 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-1,1:jpjm1) & 
     89              + localviscsponge(nlci-spongearea + 1:nlci-1,2:jpj)) & 
     90           * e1v(nlci-spongearea + 1:nlci-1,1:jpjm1) / e2v(nlci-spongearea + 1:nlci-1,1:jpjm1) 
    9591      ENDIF 
    9692 
     
    10096            localviscsponge(:,jj) = visc_tra * (spongearea-jj)/real(spongearea-2) 
    10197         ENDDO 
    102          DO jk = 1, jpkm1 
    103             vmasktemp(:,2:spongearea-1,jk) = vmask(:,2:spongearea-1,jk) & 
    104                * 0.5 * (localviscsponge(:,2:spongearea-1) + localviscsponge(:,3:spongearea)) 
    105          ENDDO 
    106          DO jk = 1, jpkm1 
    107             umasktemp(1:jpim1,2:spongearea,jk) = umask(1:jpim1,2:spongearea,jk) & 
    108                * 0.5 * (localviscsponge(1:jpim1,2:spongearea) + localviscsponge(2:jpi,2:spongearea)) 
    109          ENDDO 
     98     
     99    spe1ur(1:jpim1,2:spongearea)=0.5 * (localviscsponge(1:jpim1,2:spongearea) + & 
     100           localviscsponge(2:jpi,2:spongearea)) & 
     101          * e2u(1:jpim1,2:spongearea) / e1u(1:jpim1,2:spongearea) 
     102 
     103         spe2vr(:,2:spongearea-1) = 0.5 * (localviscsponge(:,2:spongearea-1) + & 
     104             localviscsponge(:,3:spongearea)) & 
     105           * e1v(:,2:spongearea-1) / e2v(:,2:spongearea-1) 
    110106      ENDIF 
    111107 
     
    114110            localviscsponge(:,jj) = visc_tra * (jj - (nlcj-spongearea+1))/real(spongearea-2) 
    115111         ENDDO 
    116          DO jk = 1, jpkm1 
    117             vmasktemp(:,nlcj-spongearea + 1:nlcj-2,jk) = vmask(:,nlcj-spongearea + 1:nlcj-2,jk) & 
    118                * 0.5 * (localviscsponge(:,nlcj-spongearea + 1:nlcj-2) + localviscsponge(:,nlcj-spongearea + 2:nlcj-1)) 
    119          ENDDO 
    120          DO jk = 1, jpkm1 
    121             umasktemp(1:jpim1,nlcj-spongearea + 1:nlcj-1,jk) = umask(1:jpim1,nlcj-spongearea + 1:nlcj-1,jk) & 
    122                * 0.5 * (localviscsponge(1:jpim1,nlcj-spongearea + 1:nlcj-1) + localviscsponge(2:jpi,nlcj-spongearea + 1:nlcj-1)) 
    123          ENDDO 
    124       ENDIF 
    125  
    126       IF (.NOT. spongedoneT) THEN 
    127          zspe1ur(:,:) = e2u(:,:) / e1u(:,:) 
    128          zspe2vr(:,:) = e1v(:,:) / e2v(:,:) 
    129          zspbtr2(:,:) = 1. / ( e1t(:,:) * e2t(:,:)) 
     112     
     113    spe1ur(1:jpim1,nlcj-spongearea + 1:nlcj-1)=0.5 * (localviscsponge(1:jpim1,nlcj-spongearea + 1:nlcj-1) + & 
     114            localviscsponge(2:jpi,nlcj-spongearea + 1:nlcj-1)) & 
     115          * e2u(1:jpim1,nlcj-spongearea + 1:nlcj-1) / e1u(1:jpim1,nlcj-spongearea + 1:nlcj-1) 
     116 
     117         spe2vr(:,nlcj-spongearea + 1:nlcj-2) = 0.5 * (localviscsponge(:,nlcj-spongearea + 1:nlcj-2) + & 
     118            localviscsponge(:,nlcj-spongearea + 2:nlcj-1)) & 
     119           * e1v(:,nlcj-spongearea + 1:nlcj-2) / e2v(:,nlcj-spongearea + 1:nlcj-2) 
     120      ENDIF 
     121       
     122         spbtr2(:,:) = 1. / ( e1t(:,:) * e2t(:,:)) 
    130123 
    131124         spongedoneT = .TRUE. 
     
    136129            DO ji = 1, jpim1 
    137130#if defined key_zco 
    138                zabe1 = umasktemp(ji,jj,jk) * zspe1ur(ji,jj) 
    139                zabe2 = vmasktemp(ji,jj,jk) * zspe2vr(ji,jj) 
    140 #else 
    141                zabe1 = umasktemp(ji,jj,jk) * zspe1ur(ji,jj) * fse3u(ji,jj,jk) 
    142                zabe2 = vmasktemp(ji,jj,jk) * zspe2vr(ji,jj) * fse3v(ji,jj,jk) 
     131               zabe1 = umask(ji,jj,jk) * spe1ur(ji,jj) 
     132               zabe2 = vmask(ji,jj,jk) * spe2vr(ji,jj) 
     133#else 
     134               zabe1 = umask(ji,jj,jk) * spe1ur(ji,jj) * fse3u(ji,jj,jk) 
     135               zabe2 = vmask(ji,jj,jk) * spe2vr(ji,jj) * fse3v(ji,jj,jk) 
    143136#endif 
    144137               ztu(ji,jj,jk) = zabe1 * ( tbdiff(ji+1,jj  ,jk) - tbdiff(ji,jj,jk) ) 
     
    152145            DO ji = 2,jpim1 
    153146#if defined key_zco 
    154                zbtr = zspbtr2(ji,jj) 
    155 #else 
    156                zbtr = zspbtr2(ji,jj) / fse3t(ji,jj,jk) 
     147               zbtr = spbtr2(ji,jj) 
     148#else 
     149               zbtr = spbtr2(ji,jj) / fse3t(ji,jj,jk) 
    157150#endif 
    158151               ! horizontal diffusive trends 
     
    173166   END SUBROUTINE Agrif_Sponge_Tra 
    174167 
    175    SUBROUTINE Agrif_Sponge_dyn( kt ) 
     168   SUBROUTINE Agrif_Sponge_dyn 
    176169      !!--------------------------------------------- 
    177170      !!   *** ROUTINE Agrif_Sponge_dyn *** 
    178171      !!--------------------------------------------- 
    179172#include "domzgr_substitute.h90" 
    180  
    181       INTEGER,INTENT(in) :: kt 
    182173 
    183174      INTEGER :: ji,jj,jk 
    184175      INTEGER :: spongearea 
    185176      REAL(wp) :: timecoeff 
    186       REAL(wp) :: ze2u, ze1v, zua, zva 
     177      REAL(wp) :: ze2u, ze1v, zua, zva, zbtr 
    187178      REAL(wp), DIMENSION(jpi,jpj) :: localviscsponge 
    188179      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ztab, ubdiff, vbdiff,rotdiff,hdivdiff 
    189       REAL(wp), DIMENSION(jpi,jpj,jpk) :: umasktemp,vmasktemp 
    190180 
    191181#if defined SPONGE 
     
    194184 
    195185      Agrif_SpecialValue=0. 
    196       Agrif_UseSpecialValue = .TRUE. 
     186      Agrif_UseSpecialValue = ln_spc_dyn 
    197187      ztab = 0.e0 
    198188      CALL Agrif_Bc_Variable(ztab, ua,calledweight=timecoeff,procname=interpun) 
    199189      Agrif_UseSpecialValue = .FALSE. 
    200190 
    201       ubdiff(:,:,:) = ub(:,:,:) - ztab(:,:,:) 
     191      ubdiff(:,:,:) = (ub(:,:,:) - ztab(:,:,:))*umask(:,:,:) 
    202192 
    203193      ztab = 0.e0 
    204194      Agrif_SpecialValue=0. 
    205       Agrif_UseSpecialValue = .TRUE. 
     195      Agrif_UseSpecialValue = ln_spc_dyn 
    206196      CALL Agrif_Bc_Variable(ztab, va,calledweight=timecoeff,procname=interpvn) 
    207197      Agrif_UseSpecialValue = .FALSE. 
    208198 
    209       vbdiff(:,:,:) = vb(:,:,:) - ztab(:,:,:) 
     199      vbdiff(:,:,:) = (vb(:,:,:) - ztab(:,:,:))*vmask(:,:,:) 
    210200 
    211201      spongearea = 2 + 2 * Agrif_irhox() 
    212202 
    213203      localviscsponge = 0. 
    214       umasktemp = 0. 
    215       vmasktemp = 0. 
     204 
     205      IF (.NOT. spongedoneU) THEN 
     206         spe1ur2(:,:) = 0. 
     207         spe2vr2(:,:) = 0. 
    216208 
    217209      IF ((nbondi == -1).OR.(nbondi == 2)) THEN 
     
    219211            localviscsponge(ji,:) = visc_dyn * (spongearea-ji)/real(spongearea-2) 
    220212         ENDDO 
    221          DO jk = 1, jpkm1 
    222             umasktemp(2:spongearea-1,:,jk) = umask(2:spongearea-1,:,jk) & 
    223                * 0.5 * (localviscsponge(2:spongearea-1,:) + localviscsponge(3:spongearea,:)) 
    224          ENDDO 
    225          DO jk = 1, jpkm1 
    226             vmasktemp(2:spongearea,1:jpjm1,jk) = vmask(2:spongearea,1:jpjm1,jk) & 
    227                * 0.5 * (localviscsponge(2:spongearea,1:jpjm1) + localviscsponge(2:spongearea,2:jpj)) 
    228          ENDDO 
     213     
     214    spe1ur2(2:spongearea-1,:)=0.5 * (localviscsponge(2:spongearea-1,:) + localviscsponge(3:spongearea,:)) 
     215 
     216         spe2vr2(2:spongearea,1:jpjm1) = 0.5 * (localviscsponge(2:spongearea,1:jpjm1) + & 
     217             localviscsponge(2:spongearea,2:jpj)) 
    229218      ENDIF 
    230219 
     
    233222            localviscsponge(ji,:) = visc_dyn * (ji - (nlci-spongearea+1))/real(spongearea-2) 
    234223         ENDDO 
    235          DO jk = 1, jpkm1 
    236             umasktemp(nlci-spongearea + 1:nlci-2,:,jk) = umask(nlci-spongearea + 1:nlci-2,:,jk) & 
    237                * 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-2,:) + localviscsponge(nlci-spongearea + 2:nlci-1,:)) 
    238          ENDDO 
    239          DO jk = 1, jpkm1 
    240             vmasktemp(nlci-spongearea + 1:nlci-1,1:jpjm1,jk) = vmask(nlci-spongearea + 1:nlci-1,1:jpjm1,jk) & 
    241                * 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-1,1:jpjm1) + localviscsponge(nlci-spongearea + 1:nlci-1,2:jpj)) 
    242          ENDDO 
    243       ENDIF 
     224     
     225    spe1ur2(nlci-spongearea + 1:nlci-2,:)=0.5 * (localviscsponge(nlci-spongearea + 1:nlci-2,:) + & 
     226           localviscsponge(nlci-spongearea + 2:nlci-1,:)) 
     227 
     228         spe2vr2(nlci-spongearea + 1:nlci-1,1:jpjm1) = 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-1,1:jpjm1) & 
     229              + localviscsponge(nlci-spongearea + 1:nlci-1,2:jpj)) 
     230      ENDIF 
     231 
    244232 
    245233      IF ((nbondj == -1).OR.(nbondj == 2)) THEN 
     
    247235            localviscsponge(:,jj) = visc_dyn * (spongearea-jj)/real(spongearea-2) 
    248236         ENDDO 
    249          DO jk = 1, jpkm1 
    250             vmasktemp(:,2:spongearea-1,jk) = vmask(:,2:spongearea-1,jk) & 
    251                * 0.5 * (localviscsponge(:,2:spongearea-1) + localviscsponge(:,3:spongearea)) 
    252          ENDDO 
    253          DO jk = 1, jpkm1 
    254             umasktemp(1:jpim1,2:spongearea,jk) = umask(1:jpim1,2:spongearea,jk) & 
    255                * 0.5 * (localviscsponge(1:jpim1,2:spongearea) + localviscsponge(2:jpi,2:spongearea)) 
    256          ENDDO 
     237     
     238    spe1ur2(1:jpim1,2:spongearea)=0.5 * (localviscsponge(1:jpim1,2:spongearea) + & 
     239           localviscsponge(2:jpi,2:spongearea)) 
     240 
     241         spe2vr2(:,2:spongearea-1) = 0.5 * (localviscsponge(:,2:spongearea-1) + & 
     242             localviscsponge(:,3:spongearea)) 
    257243      ENDIF 
    258244 
     
    261247            localviscsponge(:,jj) = visc_dyn * (jj - (nlcj-spongearea+1))/real(spongearea-2) 
    262248         ENDDO 
    263          DO jk = 1, jpkm1 
    264             vmasktemp(:,nlcj-spongearea + 1:nlcj-2,jk) = vmask(:,nlcj-spongearea + 1:nlcj-2,jk) & 
    265                * 0.5 * (localviscsponge(:,nlcj-spongearea + 1:nlcj-2) + localviscsponge(:,nlcj-spongearea + 2:nlcj-1)) 
    266          ENDDO 
    267          DO jk = 1, jpkm1 
    268             umasktemp(1:jpim1,nlcj-spongearea + 1:nlcj-1,jk) = umask(1:jpim1,nlcj-spongearea + 1:nlcj-1,jk) & 
    269                * 0.5 * (localviscsponge(1:jpim1,nlcj-spongearea + 1:nlcj-1) + localviscsponge(2:jpi,nlcj-spongearea + 1:nlcj-1)) 
    270          ENDDO 
    271       ENDIF 
    272  
    273       ubdiff = ubdiff * umasktemp 
    274       vbdiff = vbdiff * vmasktemp 
    275  
     249     
     250    spe1ur2(1:jpim1,nlcj-spongearea + 1:nlcj-1)=0.5 * (localviscsponge(1:jpim1,nlcj-spongearea + 1:nlcj-1) + & 
     251            localviscsponge(2:jpi,nlcj-spongearea + 1:nlcj-1)) 
     252 
     253         spe2vr2(:,nlcj-spongearea + 1:nlcj-2) = 0.5 * (localviscsponge(:,nlcj-spongearea + 1:nlcj-2) + & 
     254            localviscsponge(:,nlcj-spongearea + 2:nlcj-1)) 
     255      ENDIF 
     256 
     257         spongedoneU = .TRUE. 
     258     
     259    spbtr3(:,:) = 1./( e1f(:,:) * e2f(:,:)) 
     260      ENDIF 
     261       
     262      IF (.NOT. spongedoneT) THEN 
     263        spbtr2(:,:) = 1. / ( e1t(:,:) * e2t(:,:))       
     264      ENDIF 
     265       
     266      DO jk=1,jpkm1 
     267      ubdiff(:,:,jk) = ubdiff(:,:,jk) * spe1ur2(:,:) 
     268      vbdiff(:,:,jk) = vbdiff(:,:,jk) * spe2vr2(:,:) 
     269      ENDDO 
     270       
    276271      hdivdiff = 0. 
    277272      rotdiff = 0. 
     
    286281            DO ji = 2, jpim1   ! vector opt. 
    287282#if defined key_zco 
     283               zbtr = spbtr2(ji,jj) 
    288284               hdivdiff(ji,jj,jk) = (  e2u(ji,jj) * ubdiff(ji,jj,jk) & 
    289285                  - e2u(ji-1,jj  ) * ubdiff(ji-1,jj  ,jk)      & 
    290286                  &               + e1v(ji,jj) * vbdiff(ji,jj,jk) - & 
    291                   &              e1v(ji  ,jj-1) * vbdiff(ji  ,jj-1,jk)  )   & 
    292                   &            / ( e1t(ji,jj) * e2t(ji,jj) ) 
    293 #else 
     287                  &              e1v(ji  ,jj-1) * vbdiff(ji  ,jj-1,jk)  ) * zbtr 
     288#else 
     289               zbtr = spbtr2(ji,jj) / fse3t(ji,jj,jk) 
    294290               hdivdiff(ji,jj,jk) =   & 
    295291                  (  e2u(ji,jj)*fse3u(ji,jj,jk) * &  
     
    298294                  + e1v(ji,jj)*fse3v(ji,jj,jk) * & 
    299295                  vbdiff(ji,jj,jk) - e1v(ji  ,jj-1)* & 
    300                   fse3v(ji  ,jj-1,jk)  * vbdiff(ji  ,jj-1,jk)  )    & 
    301                   / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) 
     296                  fse3v(ji  ,jj-1,jk)  * vbdiff(ji  ,jj-1,jk)  ) * zbtr 
    302297#endif 
    303298            END DO 
     
    306301         DO jj = 1, jpjm1 
    307302            DO ji = 1, jpim1   ! vector opt. 
     303#if defined key_zco       
     304               zbtr = spbtr3(ji,jj) 
    308305               rotdiff(ji,jj,jk) = (  e2v(ji+1,jj  ) * vbdiff(ji+1,jj  ,jk) - e2v(ji,jj) * vbdiff(ji,jj,jk)    & 
    309306                  &              - e1u(ji  ,jj+1) * ubdiff(ji  ,jj+1,jk) + e1u(ji,jj) * ubdiff(ji,jj,jk)  ) & 
    310                   &           * fmask(ji,jj,jk) / ( e1f(ji,jj) * e2f(ji,jj) ) 
     307                  &           * fmask(ji,jj,jk) * zbtr 
     308#else 
     309               zbtr = spbtr3(ji,jj) * fse3f(ji,jj,jk) 
     310               rotdiff(ji,jj,jk) = (  e2v(ji+1,jj  ) * vbdiff(ji+1,jj  ,jk) - e2v(ji,jj) * vbdiff(ji,jj,jk)    & 
     311                  &              - e1u(ji  ,jj+1) * ubdiff(ji  ,jj+1,jk) + e1u(ji,jj) * ubdiff(ji,jj,jk)  ) & 
     312                  &           * fmask(ji,jj,jk) * zbtr 
     313#endif         
    311314            END DO 
    312315         END DO 
     
    333336                  ze1v                  ) / e2v(ji,jj) 
    334337#else 
    335                ze2u = rotdiff (ji,jj,jk)*fse3f(ji,jj,jk) 
     338               ze2u = rotdiff (ji,jj,jk) 
    336339               ze1v = hdivdiff(ji,jj,jk) 
    337340               ! horizontal diffusive trends 
    338                zua = - ( ze2u - rotdiff (ji,jj-1,jk)* & 
    339                   fse3f(ji,jj-1,jk) ) / ( e2u(ji,jj) * fse3u(ji,jj,jk) )   & 
     341               zua = - ( ze2u - rotdiff (ji,jj-1,jk)) / ( e2u(ji,jj) * fse3u(ji,jj,jk) )   & 
    340342                  + ( hdivdiff(ji+1,jj,jk) - ze1v      & 
    341343                  ) / e1u(ji,jj) 
    342344 
    343                zva = + ( ze2u - rotdiff (ji-1,jj,jk)* & 
    344                   fse3f(ji-1,jj,jk) ) / ( e1v(ji,jj) * fse3v(ji,jj,jk) )   & 
     345               zva = + ( ze2u - rotdiff (ji-1,jj,jk)) / ( e1v(ji,jj) * fse3v(ji,jj,jk) )   & 
    345346                  + ( hdivdiff(ji,jj+1,jk) - ze1v    & 
    346347                  ) / e2v(ji,jj) 
     
    360361   END SUBROUTINE Agrif_Sponge_dyn 
    361362 
    362    SUBROUTINE agrif_sponge_init 
    363       !!--------------------------------------------- 
    364       !!   *** ROUTINE agrif_sponge_init *** 
    365       !!--------------------------------------------- 
    366       NAMELIST/namagrif/ visc_tra, visc_dyn 
    367       REWIND ( numnam ) 
    368       READ   ( numnam, namagrif ) 
    369  
    370       IF(lwp) THEN 
    371          WRITE(numout,*) 
    372          WRITE(numout,*) 'agrif_sponge_init : agrif sponge parameters' 
    373          WRITE(numout,*) '~~~~~~~~~~~~' 
    374          WRITE(numout,*) '          Namelist namagrif : set sponge parameters' 
    375          WRITE(numout,*) '             sponge coefficient for tracers =  ', visc_tra 
    376          WRITE(numout,*) '             sponge coefficient for dynamics = ', visc_dyn 
    377       ENDIF 
    378  
    379    END SUBROUTINE agrif_sponge_init 
    380  
    381363   SUBROUTINE interptn(tabres,i1,i2,j1,j2,k1,k2) 
    382364      !!--------------------------------------------- 
     
    405387   END SUBROUTINE interpsn 
    406388 
    407  
    408389   SUBROUTINE interpun(tabres,i1,i2,j1,j2,k1,k2) 
    409390      !!--------------------------------------------- 
Note: See TracChangeset for help on using the changeset viewer.