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

Ignore:
Timestamp:
2007-03-07T14:28:16+01:00 (17 years ago)
Author:
opalod
Message:

nemo_v2_update_008:RB: clean agrif routines and add sponge layer coefficient in namelist

File:
1 edited

Legend:

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

    r469 r636  
    11#define SPONGE 
    22 
    3       Module agrif_opa_sponge 
     3Module agrif_opa_sponge 
    44#if defined key_agrif 
    5       USE par_oce 
    6       USE oce 
    7       USE dom_oce 
    8        
    9       Contains 
    10  
    11  
    12       Subroutine Agrif_Sponge_Tra( kt ) 
    13  
    14       implicit none 
    15  
    16       INTEGER :: kt 
    17       REAL(wp), DIMENSION(jpi,jpj,jpk) :: tabtemp, tbdiff, sbdiff 
     5   USE par_oce 
     6   USE oce 
     7   USE dom_oce 
     8   USE in_out_manager 
     9 
     10   IMPLICIT NONE 
     11   PRIVATE 
     12 
     13   PUBLIC Agrif_Sponge_Tra, Agrif_Sponge_Dyn, interptn, interpsn, interpun, interpvn 
     14 
     15   !! * Namelist (namagrif) 
     16   REAL(wp) :: visc_tra = rdt 
     17   REAL(wp) :: visc_dyn = rdt 
     18 
     19   CONTAINS 
     20 
     21   SUBROUTINE Agrif_Sponge_Tra( kt ) 
     22      !!--------------------------------------------- 
     23      !!   *** ROUTINE Agrif_Sponge_Tra *** 
     24      !!--------------------------------------------- 
     25#include "domzgr_substitute.h90" 
     26 
     27      INTEGER, INTENT(in) :: kt 
     28 
    1829      INTEGER :: ji,jj,jk 
    19       REAL(wp) :: viscsponge 
    2030      REAL(wp), DIMENSION(jpi,jpj,jpk) :: umasktemp,vmasktemp 
    2131      INTEGER :: spongearea 
    22       integer ipt,jpt 
    23       real,dimension(:,:),pointer :: e1tparent,e2tparent 
     32      REAL(wp) :: timecoeff 
     33      REAL(wp) :: zta, zsa, zabe1, zabe2, zbtr 
    2434      REAL(wp), DIMENSION(jpi,jpj) :: localviscsponge 
    25       real(wp) :: timecoeff 
     35      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ztab, tbdiff, sbdiff 
    2636      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ztu ,ztv, zsu ,zsv 
    27       REAL(wp) :: zta, zsa, zabe1, zabe2, zbtr 
    28  
    29 #include "domzgr_substitute.h90" 
    30  
    31          
     37 
    3238#if defined SPONGE 
    3339 
    34       timecoeff = real(Agrif_NbStepint())/Agrif_rhot() 
     40      IF( kt == nit000  )   CALL agrif_sponge_init 
     41 
     42      timecoeff = REAL(Agrif_NbStepint(),wp)/Agrif_rhot() 
    3543 
    3644      Agrif_SpecialValue=0. 
    3745      Agrif_UseSpecialValue = .TRUE. 
    38       tabtemp = 0. 
    39       Call Agrif_Bc_Variable(tabtemp, ta,calledweight=timecoeff,procname=interptn) 
     46      ztab = 0.e0 
     47      CALL Agrif_Bc_Variable(ztab, ta,calledweight=timecoeff,procname=interptn) 
    4048      Agrif_UseSpecialValue = .FALSE. 
    4149 
    42       tbdiff(:,:,:) = tb(:,:,:) - tabtemp(:,:,:) 
    43  
    44       tabtemp = 0. 
     50      tbdiff(:,:,:) = tb(:,:,:) - ztab(:,:,:) 
     51 
     52      ztab = 0.e0 
    4553      Agrif_SpecialValue=0. 
    4654      Agrif_UseSpecialValue = .TRUE. 
    47       Call Agrif_Bc_Variable(tabtemp, sa,calledweight=timecoeff,procname=interpsn) 
     55      CALL Agrif_Bc_Variable(ztab, sa,calledweight=timecoeff,procname=interpsn) 
    4856      Agrif_UseSpecialValue = .FALSE. 
    4957 
    50       sbdiff(:,:,:) = sb(:,:,:) - tabtemp(:,:,:) 
    51         
    52       viscsponge = rdt 
     58      sbdiff(:,:,:) = sb(:,:,:) - ztab(:,:,:) 
    5359 
    5460      spongearea = 2 + 2 * Agrif_irhox() 
     
    5965 
    6066      IF ((nbondi == -1).OR.(nbondi == 2)) THEN 
    61  
    62         DO ji = 2, spongearea 
    63           localviscsponge(ji,:) = viscsponge * (spongearea-ji)/real(spongearea-2) 
    64         ENDDO 
    65  
    66         DO jk = 1, jpkm1 
    67           umasktemp(2:spongearea-1,:,jk) = umask(2:spongearea-1,:,jk) & 
    68            * 0.5 * (localviscsponge(2:spongearea-1,:) + localviscsponge(3:spongearea,:)) 
    69         ENDDO 
    70  
    71        DO jk = 1, jpkm1 
    72          vmasktemp(2:spongearea,1:jpjm1,jk) = vmask(2:spongearea,1:jpjm1,jk) & 
    73           * 0.5 * (localviscsponge(2:spongearea,1:jpjm1) + localviscsponge(2:spongearea,2:jpj)) 
    74        ENDDO 
    75  
     67         DO ji = 2, spongearea 
     68            localviscsponge(ji,:) = visc_tra * (spongearea-ji)/real(spongearea-2) 
     69         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 
    7678      ENDIF 
    7779 
    7880      IF ((nbondi == 1).OR.(nbondi == 2)) THEN 
    79  
    80         DO ji = nlci-spongearea + 1,nlci-1 
    81           localviscsponge(ji,:) = viscsponge * (ji - (nlci-spongearea+1))/real(spongearea-2) 
    82         ENDDO 
    83  
    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  
    89        DO jk = 1, jpkm1 
    90         vmasktemp(nlci-spongearea + 1:nlci-1,1:jpjm1,jk) = vmask(nlci-spongearea + 1:nlci-1,1:jpjm1,jk) & 
    91           * 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-1,1:jpjm1) + localviscsponge(nlci-spongearea + 1:nlci-1,2:jpj)) 
    92        ENDDO 
    93  
    94       ENDIF 
    95  
     81         DO ji = nlci-spongearea + 1,nlci-1 
     82            localviscsponge(ji,:) = visc_tra * (ji - (nlci-spongearea+1))/real(spongearea-2) 
     83         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 
     92      ENDIF 
    9693 
    9794 
    9895      IF ((nbondj == -1).OR.(nbondj == 2)) THEN 
    99  
    100         DO jj = 2, spongearea 
    101         localviscsponge(:,jj) = viscsponge * (spongearea-jj)/real(spongearea-2) 
    102         ENDDO 
    103  
    104       DO jk = 1, jpkm1 
    105       vmasktemp(:,2:spongearea-1,jk) = vmask(:,2:spongearea-1,jk) & 
    106          * 0.5 * (localviscsponge(:,2:spongearea-1) + localviscsponge(:,3:spongearea)) 
    107         ENDDO 
    108  
    109       DO jk = 1, jpkm1 
    110        umasktemp(1:jpim1,2:spongearea,jk) = umask(1:jpim1,2:spongearea,jk) & 
    111          * 0.5 * (localviscsponge(1:jpim1,2:spongearea) + localviscsponge(2:jpi,2:spongearea)) 
    112       ENDDO 
    113  
     96         DO jj = 2, spongearea 
     97            localviscsponge(:,jj) = visc_tra * (spongearea-jj)/real(spongearea-2) 
     98         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 
    114107      ENDIF 
    115108 
    116109      IF ((nbondj == 1).OR.(nbondj == 2)) THEN 
    117  
    118         DO jj = nlcj-spongearea + 1,nlcj-1 
    119        localviscsponge(:,jj) = viscsponge * (jj - (nlcj-spongearea+1))/real(spongearea-2) 
    120        ENDDO 
    121  
    122       DO jk = 1, jpkm1 
    123        vmasktemp(:,nlcj-spongearea + 1:nlcj-2,jk) = vmask(:,nlcj-spongearea + 1:nlcj-2,jk) & 
    124           * 0.5 * (localviscsponge(:,nlcj-spongearea + 1:nlcj-2) + localviscsponge(:,nlcj-spongearea + 2:nlcj-1)) 
    125        ENDDO 
    126  
    127       DO jk = 1, jpkm1 
    128         umasktemp(1:jpim1,nlcj-spongearea + 1:nlcj-1,jk) = umask(1:jpim1,nlcj-spongearea + 1:nlcj-1,jk) & 
    129          * 0.5 * (localviscsponge(1:jpim1,nlcj-spongearea + 1:nlcj-1) + localviscsponge(2:jpi,nlcj-spongearea + 1:nlcj-1)) 
    130       ENDDO 
    131  
    132       ENDIF 
    133  
    134       IF (.Not. spongedoneT) THEN 
     110         DO jj = nlcj-spongearea + 1,nlcj-1 
     111            localviscsponge(:,jj) = visc_tra * (jj - (nlcj-spongearea+1))/real(spongearea-2) 
     112         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 
    135124         zspe1ur(:,:) = e2u(:,:) / e1u(:,:) 
    136125         zspe2vr(:,:) = e1v(:,:) / e2v(:,:) 
    137126         zspbtr2(:,:) = 1. / ( e1t(:,:) * e2t(:,:)) 
    138           
     127 
    139128         spongedoneT = .TRUE. 
    140129      ENDIF 
    141130 
    142         DO jk = 1, jpkm1 
    143           DO jj = 1, jpjm1 
     131      DO jk = 1, jpkm1 
     132         DO jj = 1, jpjm1 
    144133            DO ji = 1, jpim1 
    145134#if defined key_zco 
     
    155144               zsv(ji,jj,jk) = zabe2 * ( sbdiff(ji  ,jj+1,jk) - sbdiff(ji,jj,jk) ) 
    156145            ENDDO 
    157           ENDDO 
     146         ENDDO 
    158147 
    159148         DO jj = 2,jpjm1 
     
    175164         END DO 
    176165 
    177         ENDDO 
    178  
    179 #endif 
    180  
    181       Return 
    182       End Subroutine Agrif_Sponge_Tra 
    183        
    184       Subroutine Agrif_Sponge_dyn( kt ) 
    185  
    186       implicit none 
    187  
    188       INTEGER :: kt 
    189       REAL(wp), DIMENSION(jpi,jpj,jpk) :: tabtemp, ubdiff, vbdiff,rotdiff,hdivdiff 
     166      ENDDO 
     167 
     168#endif 
     169 
     170   END SUBROUTINE Agrif_Sponge_Tra 
     171 
     172   SUBROUTINE Agrif_Sponge_dyn( kt ) 
     173      !!--------------------------------------------- 
     174      !!   *** ROUTINE Agrif_Sponge_dyn *** 
     175      !!--------------------------------------------- 
     176#include "domzgr_substitute.h90" 
     177 
     178      INTEGER,INTENT(in) :: kt 
     179 
    190180      INTEGER :: ji,jj,jk 
    191       REAL(wp) :: viscsponge 
     181      INTEGER :: spongearea 
     182      REAL(wp) :: timecoeff 
     183      REAL(wp) :: ze2u, ze1v, zua, zva 
     184      REAL(wp), DIMENSION(jpi,jpj) :: localviscsponge 
     185      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ztab, ubdiff, vbdiff,rotdiff,hdivdiff 
    192186      REAL(wp), DIMENSION(jpi,jpj,jpk) :: umasktemp,vmasktemp 
    193       INTEGER :: spongearea 
    194       integer ipt,jpt 
    195       real,dimension(:,:),pointer :: e1tparent,e2tparent 
    196       REAL(wp), DIMENSION(jpi,jpj) :: localviscsponge 
    197       real(wp) :: timecoeff 
    198       REAL(wp):: ze2u, ze1v, zua, zva 
    199  
    200 #include "domzgr_substitute.h90" 
    201187 
    202188#if defined SPONGE 
    203189 
    204       timecoeff = real(Agrif_NbStepint())/Agrif_rhot() 
     190      timecoeff = REAL(Agrif_NbStepint(),wp)/Agrif_rhot() 
    205191 
    206192      Agrif_SpecialValue=0. 
    207193      Agrif_UseSpecialValue = .TRUE. 
    208       tabtemp = 0. 
    209       Call Agrif_Bc_Variable(tabtemp, ua,calledweight=timecoeff,procname=interpun) 
     194      ztab = 0.e0 
     195      CALL Agrif_Bc_Variable(ztab, ua,calledweight=timecoeff,procname=interpun) 
    210196      Agrif_UseSpecialValue = .FALSE. 
    211197 
    212       ubdiff(:,:,:) = ub(:,:,:) - tabtemp(:,:,:) 
    213  
    214       tabtemp = 0. 
     198      ubdiff(:,:,:) = ub(:,:,:) - ztab(:,:,:) 
     199 
     200      ztab = 0.e0 
    215201      Agrif_SpecialValue=0. 
    216202      Agrif_UseSpecialValue = .TRUE. 
    217       Call Agrif_Bc_Variable(tabtemp, va,calledweight=timecoeff,procname=interpvn) 
     203      CALL Agrif_Bc_Variable(ztab, va,calledweight=timecoeff,procname=interpvn) 
    218204      Agrif_UseSpecialValue = .FALSE. 
    219205 
    220       vbdiff(:,:,:) = vb(:,:,:) - tabtemp(:,:,:) 
    221         
    222       viscsponge = rdt 
     206      vbdiff(:,:,:) = vb(:,:,:) - ztab(:,:,:) 
    223207 
    224208      spongearea = 2 + 2 * Agrif_irhox() 
     
    229213 
    230214      IF ((nbondi == -1).OR.(nbondi == 2)) THEN 
    231  
    232         DO ji = 2, spongearea 
    233           localviscsponge(ji,:) = viscsponge * (spongearea-ji)/real(spongearea-2) 
    234         ENDDO 
    235  
    236         DO jk = 1, jpkm1 
    237           umasktemp(2:spongearea-1,:,jk) = umask(2:spongearea-1,:,jk) & 
    238            * 0.5 * (localviscsponge(2:spongearea-1,:) + localviscsponge(3:spongearea,:)) 
    239         ENDDO 
    240  
    241        DO jk = 1, jpkm1 
    242          vmasktemp(2:spongearea,1:jpjm1,jk) = vmask(2:spongearea,1:jpjm1,jk) & 
    243           * 0.5 * (localviscsponge(2:spongearea,1:jpjm1) + localviscsponge(2:spongearea,2:jpj)) 
    244        ENDDO 
    245  
     215         DO ji = 2, spongearea 
     216            localviscsponge(ji,:) = visc_dyn * (spongearea-ji)/real(spongearea-2) 
     217         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 
    246226      ENDIF 
    247227 
    248228      IF ((nbondi == 1).OR.(nbondi == 2)) THEN 
    249  
    250         DO ji = nlci-spongearea + 1,nlci-1 
    251           localviscsponge(ji,:) = viscsponge * (ji - (nlci-spongearea+1))/real(spongearea-2) 
    252         ENDDO 
    253  
    254        DO jk = 1, jpkm1 
    255         umasktemp(nlci-spongearea + 1:nlci-2,:,jk) = umask(nlci-spongearea + 1:nlci-2,:,jk) & 
    256          * 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-2,:) + localviscsponge(nlci-spongearea + 2:nlci-1,:)) 
    257        ENDDO 
    258  
    259        DO jk = 1, jpkm1 
    260         vmasktemp(nlci-spongearea + 1:nlci-1,1:jpjm1,jk) = vmask(nlci-spongearea + 1:nlci-1,1:jpjm1,jk) & 
    261           * 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-1,1:jpjm1) + localviscsponge(nlci-spongearea + 1:nlci-1,2:jpj)) 
    262        ENDDO 
    263  
    264       ENDIF 
    265  
    266  
     229         DO ji = nlci-spongearea + 1,nlci-1 
     230            localviscsponge(ji,:) = visc_dyn * (ji - (nlci-spongearea+1))/real(spongearea-2) 
     231         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 
    267241 
    268242      IF ((nbondj == -1).OR.(nbondj == 2)) THEN 
    269  
    270         DO jj = 2, spongearea 
    271         localviscsponge(:,jj) = viscsponge * (spongearea-jj)/real(spongearea-2) 
    272         ENDDO 
    273  
    274       DO jk = 1, jpkm1 
    275       vmasktemp(:,2:spongearea-1,jk) = vmask(:,2:spongearea-1,jk) & 
    276          * 0.5 * (localviscsponge(:,2:spongearea-1) + localviscsponge(:,3:spongearea)) 
    277         ENDDO 
    278  
    279       DO jk = 1, jpkm1 
    280        umasktemp(1:jpim1,2:spongearea,jk) = umask(1:jpim1,2:spongearea,jk) & 
    281          * 0.5 * (localviscsponge(1:jpim1,2:spongearea) + localviscsponge(2:jpi,2:spongearea)) 
    282       ENDDO 
    283  
     243         DO jj = 2, spongearea 
     244            localviscsponge(:,jj) = visc_dyn * (spongearea-jj)/real(spongearea-2) 
     245         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 
    284254      ENDIF 
    285255 
    286256      IF ((nbondj == 1).OR.(nbondj == 2)) THEN 
    287  
    288         DO jj = nlcj-spongearea + 1,nlcj-1 
    289        localviscsponge(:,jj) = viscsponge * (jj - (nlcj-spongearea+1))/real(spongearea-2) 
    290        ENDDO 
    291  
    292       DO jk = 1, jpkm1 
    293        vmasktemp(:,nlcj-spongearea + 1:nlcj-2,jk) = vmask(:,nlcj-spongearea + 1:nlcj-2,jk) & 
    294           * 0.5 * (localviscsponge(:,nlcj-spongearea + 1:nlcj-2) + localviscsponge(:,nlcj-spongearea + 2:nlcj-1)) 
    295        ENDDO 
    296  
    297       DO jk = 1, jpkm1 
    298         umasktemp(1:jpim1,nlcj-spongearea + 1:nlcj-1,jk) = umask(1:jpim1,nlcj-spongearea + 1:nlcj-1,jk) & 
    299          * 0.5 * (localviscsponge(1:jpim1,nlcj-spongearea + 1:nlcj-1) + localviscsponge(2:jpi,nlcj-spongearea + 1:nlcj-1)) 
    300       ENDDO 
    301  
    302       ENDIF 
    303        
     257         DO jj = nlcj-spongearea + 1,nlcj-1 
     258            localviscsponge(:,jj) = visc_dyn * (jj - (nlcj-spongearea+1))/real(spongearea-2) 
     259         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 
    304270      ubdiff = ubdiff * umasktemp 
    305271      vbdiff = vbdiff * vmasktemp 
    306        
     272 
    307273      hdivdiff = 0. 
    308274      rotdiff = 0. 
     
    318284#if defined key_zco 
    319285               hdivdiff(ji,jj,jk) = (  e2u(ji,jj) * ubdiff(ji,jj,jk) & 
    320                - e2u(ji-1,jj  ) * ubdiff(ji-1,jj  ,jk)      & 
    321      &               + e1v(ji,jj) * vbdiff(ji,jj,jk) - & 
    322      &              e1v(ji  ,jj-1) * vbdiff(ji  ,jj-1,jk)  )   & 
    323      &            / ( e1t(ji,jj) * e2t(ji,jj) ) 
     286                  - e2u(ji-1,jj  ) * ubdiff(ji-1,jj  ,jk)      & 
     287                  &               + e1v(ji,jj) * vbdiff(ji,jj,jk) - & 
     288                  &              e1v(ji  ,jj-1) * vbdiff(ji  ,jj-1,jk)  )   & 
     289                  &            / ( e1t(ji,jj) * e2t(ji,jj) ) 
    324290#else 
    325291               hdivdiff(ji,jj,jk) =   & 
     
    327293                  ubdiff(ji,jj,jk) - e2u(ji-1,jj  )* & 
    328294                  fse3u(ji-1,jj  ,jk)  * ubdiff(ji-1,jj  ,jk)       & 
    329                    + e1v(ji,jj)*fse3v(ji,jj,jk) * & 
     295                  + e1v(ji,jj)*fse3v(ji,jj,jk) * & 
    330296                  vbdiff(ji,jj,jk) - e1v(ji  ,jj-1)* & 
    331297                  fse3v(ji  ,jj-1,jk)  * vbdiff(ji  ,jj-1,jk)  )    & 
     
    334300            END DO 
    335301         END DO 
    336           
     302 
    337303         DO jj = 1, jpjm1 
    338304            DO ji = 1, jpim1   ! vector opt. 
     
    342308            END DO 
    343309         END DO 
    344                    
    345          ENDDO 
    346           
     310 
     311      ENDDO 
     312 
    347313      !                                                ! =============== 
    348314      DO jk = 1, jpkm1                                 ! Horizontal slab 
     
    355321               ze1v = hdivdiff(ji,jj,jk) 
    356322               zua = - (                ze2u                  - & 
    357                rotdiff (ji,jj-1,jk) ) / e2u(ji,jj)   & 
    358                      + ( hdivdiff(ji+1,jj,jk) -     & 
    359                 ze1v                  ) / e1u(ji,jj) 
     323                  rotdiff (ji,jj-1,jk) ) / e2u(ji,jj)   & 
     324                  + ( hdivdiff(ji+1,jj,jk) -     & 
     325                  ze1v                  ) / e1u(ji,jj) 
    360326 
    361327               zva = + (                ze2u                  - & 
    362                rotdiff (ji-1,jj,jk) ) / e1v(ji,jj)   & 
    363                      + ( hdivdiff(ji,jj+1,jk) -       & 
    364                 ze1v                  ) / e2v(ji,jj) 
     328                  rotdiff (ji-1,jj,jk) ) / e1v(ji,jj)   & 
     329                  + ( hdivdiff(ji,jj+1,jk) -       & 
     330                  ze1v                  ) / e2v(ji,jj) 
    365331#else 
    366332               ze2u = rotdiff (ji,jj,jk)*fse3f(ji,jj,jk) 
     
    368334               ! horizontal diffusive trends 
    369335               zua = - ( ze2u - rotdiff (ji,jj-1,jk)* & 
    370                fse3f(ji,jj-1,jk) ) / ( e2u(ji,jj) * fse3u(ji,jj,jk) )   & 
    371                      + ( hdivdiff(ji+1,jj,jk) - ze1v      & 
    372                ) / e1u(ji,jj) 
     336                  fse3f(ji,jj-1,jk) ) / ( e2u(ji,jj) * fse3u(ji,jj,jk) )   & 
     337                  + ( hdivdiff(ji+1,jj,jk) - ze1v      & 
     338                  ) / e1u(ji,jj) 
    373339 
    374340               zva = + ( ze2u - rotdiff (ji-1,jj,jk)* & 
    375                fse3f(ji-1,jj,jk) ) / ( e1v(ji,jj) * fse3v(ji,jj,jk) )   & 
    376                      + ( hdivdiff(ji,jj+1,jk) - ze1v    & 
    377                     ) / e2v(ji,jj) 
     341                  fse3f(ji-1,jj,jk) ) / ( e1v(ji,jj) * fse3v(ji,jj,jk) )   & 
     342                  + ( hdivdiff(ji,jj+1,jk) - ze1v    & 
     343                  ) / e2v(ji,jj) 
    378344#endif 
    379345 
    380346               ! add it to the general momentum trends 
    381  
    382347               ua(ji,jj,jk) = ua(ji,jj,jk) + zua 
    383348               va(ji,jj,jk) = va(ji,jj,jk) + zva 
     
    390355#endif 
    391356 
    392       Return 
    393       End Subroutine Agrif_Sponge_dyn 
    394  
    395        subroutine interptn(tabres,i1,i2,j1,j2,k1,k2) 
    396        Implicit none 
     357   END SUBROUTINE Agrif_Sponge_dyn 
     358 
     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 
     378   SUBROUTINE interptn(tabres,i1,i2,j1,j2,k1,k2) 
     379      !!--------------------------------------------- 
     380      !!   *** ROUTINE interptn *** 
     381      !!--------------------------------------------- 
    397382#  include "domzgr_substitute.h90"        
    398        integer i1,i2,j1,j2,k1,k2 
    399        real,dimension(i1:i2,j1:j2,k1:k2) :: tabres 
    400  
    401        tabres(i1:i2,j1:j2,k1:k2) = tn(i1:i2,j1:j2,k1:k2) 
    402  
    403        end subroutine interptn     
    404         
    405        subroutine interpsn(tabres,i1,i2,j1,j2,k1,k2) 
    406        Implicit none 
     383       
     384      INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2 
     385      REAL(wp), DIMENSION(i1:i2,j1:j2,k1:k2), INTENT(inout) :: tabres 
     386 
     387      tabres(i1:i2,j1:j2,k1:k2) = tn(i1:i2,j1:j2,k1:k2) 
     388 
     389   END SUBROUTINE interptn 
     390 
     391   SUBROUTINE interpsn(tabres,i1,i2,j1,j2,k1,k2) 
     392      !!--------------------------------------------- 
     393      !!   *** ROUTINE interpsn *** 
     394      !!--------------------------------------------- 
    407395#  include "domzgr_substitute.h90"        
    408        integer i1,i2,j1,j2,k1,k2 
    409        real,dimension(i1:i2,j1:j2,k1:k2) :: tabres 
    410  
    411        tabres(i1:i2,j1:j2,k1:k2) = sn(i1:i2,j1:j2,k1:k2) 
    412  
    413        end subroutine interpsn   
    414                   
    415   
    416        subroutine interpun(tabres,i1,i2,j1,j2,k1,k2) 
    417        Implicit none 
     396       
     397      INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2 
     398      REAL(wp), DIMENSION(i1:i2,j1:j2,k1:k2), INTENT(inout) :: tabres 
     399 
     400      tabres(i1:i2,j1:j2,k1:k2) = sn(i1:i2,j1:j2,k1:k2) 
     401 
     402   END SUBROUTINE interpsn 
     403 
     404 
     405   SUBROUTINE interpun(tabres,i1,i2,j1,j2,k1,k2) 
     406      !!--------------------------------------------- 
     407      !!   *** ROUTINE interpun *** 
     408      !!--------------------------------------------- 
    418409#  include "domzgr_substitute.h90"        
    419        integer i1,i2,j1,j2,k1,k2 
    420        real,dimension(i1:i2,j1:j2,k1:k2) :: tabres 
    421  
    422        tabres(i1:i2,j1:j2,k1:k2) = un(i1:i2,j1:j2,k1:k2) 
    423  
    424        end subroutine interpun  
    425         
    426        subroutine interpvn(tabres,i1,i2,j1,j2,k1,k2) 
    427        Implicit none 
     410       
     411      INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2 
     412      REAL(wp), DIMENSION(i1:i2,j1:j2,k1:k2), INTENT(inout) :: tabres 
     413 
     414      tabres(i1:i2,j1:j2,k1:k2) = un(i1:i2,j1:j2,k1:k2) 
     415 
     416   END SUBROUTINE interpun 
     417 
     418   SUBROUTINE interpvn(tabres,i1,i2,j1,j2,k1,k2) 
     419      !!--------------------------------------------- 
     420      !!   *** ROUTINE interpvn *** 
     421      !!--------------------------------------------- 
    428422#  include "domzgr_substitute.h90"        
    429        integer i1,i2,j1,j2,k1,k2 
    430        real,dimension(i1:i2,j1:j2,k1:k2) :: tabres 
    431  
    432        tabres(i1:i2,j1:j2,k1:k2) = vn(i1:i2,j1:j2,k1:k2) 
    433  
    434        end subroutine interpvn  
     423       
     424      INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2 
     425      REAL(wp), DIMENSION(i1:i2,j1:j2,k1:k2), INTENT(inout) :: tabres 
     426 
     427      tabres(i1:i2,j1:j2,k1:k2) = vn(i1:i2,j1:j2,k1:k2) 
     428 
     429   END SUBROUTINE interpvn 
    435430 
    436431#else 
    437        CONTAINS 
    438        subroutine agrif_opa_sponge_empty 
    439        end subroutine agrif_opa_sponge_empty 
    440 #endif 
    441                      
    442        End Module agrif_opa_sponge 
     432CONTAINS 
     433 
     434   SUBROUTINE agrif_opa_sponge_empty 
     435      !!--------------------------------------------- 
     436      !!   *** ROUTINE agrif_OPA_sponge_empty *** 
     437      !!--------------------------------------------- 
     438      WRITE(*,*)  'agrif_opa_sponge : You should not have seen this print! error?' 
     439   END SUBROUTINE agrif_opa_sponge_empty 
     440#endif 
     441 
     442END MODULE agrif_opa_sponge 
Note: See TracChangeset for help on using the changeset viewer.