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 8898 for branches/2017/dev_r8624_AGRIF3_VVL/NEMOGCM/NEMO/NST_SRC/agrif_opa_update.F90 – NEMO

Ignore:
Timestamp:
2017-12-05T14:06:56+01:00 (6 years ago)
Author:
jchanut
Message:

AGRIF+VVL: Merge with free surface changes - #1965

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2017/dev_r8624_AGRIF3_VVL/NEMOGCM/NEMO/NST_SRC/agrif_opa_update.F90

    r8803 r8898  
    103103      ! 
    104104      nbcline = nbcline + 1 
     105      ! 
     106# if ! defined DECAL_FEEDBACK 
     107      ! Account for updated thicknesses at boundary edges 
     108      IF (.NOT.ln_linssh) THEN 
     109         CALL Agrif_Update_Variable(un_update_id,locupdate1=(/0,0/),locupdate2=(/0,0/),procname = correct_u_bdy) 
     110         CALL Agrif_Update_Variable(vn_update_id,locupdate1=(/0,0/),locupdate2=(/0,0/),procname = correct_v_bdy) 
     111      ENDIF 
     112# endif 
    105113      !  
    106114      IF ( ln_dynspg_ts.AND.ln_bt_fw ) THEN 
     
    379387                  ! 
    380388                  IF (.NOT.(lk_agrif_fstep.AND.(neuler==0))) THEN ! Add asselin part 
    381                      zub  = ub(ji,jj,jk) * e3u_b(ji,jj,jk) 
     389                     zub  = ub(ji,jj,jk) * e3u_b(ji,jj,jk)  ! fse3t_b prior update should be used 
    382390                     zuno = un(ji,jj,jk) * e3u_a(ji,jj,jk) 
    383391                     zunu = tabres(ji,jj,jk) 
     
    399407   END SUBROUTINE updateu 
    400408 
     409   SUBROUTINE correct_u_bdy( tabres, i1, i2, j1, j2, k1, k2, before, nb, ndir ) 
     410      !!--------------------------------------------- 
     411      !!           *** ROUTINE correct_u_bdy *** 
     412      !!--------------------------------------------- 
     413      INTEGER                               , INTENT(in   ) :: i1, i2, j1, j2, k1, k2 
     414      REAL(wp), DIMENSION(i1:i2,j1:j2,k1:k2), INTENT(inout) :: tabres 
     415      LOGICAL                               , INTENT(in   ) :: before 
     416      INTEGER                               , INTENT(in)    :: nb, ndir 
     417      !! 
     418      LOGICAL :: western_side, eastern_side  
     419      ! 
     420      INTEGER  :: jj, jk 
     421      REAL(wp) :: zcor 
     422      !!--------------------------------------------- 
     423      !  
     424      IF( .NOT.before ) THEN 
     425         ! 
     426         western_side  = (nb == 1).AND.(ndir == 1) 
     427         eastern_side  = (nb == 1).AND.(ndir == 2) 
     428         ! 
     429         IF (western_side) THEN 
     430            DO jj=j1,j2 
     431               zcor = un_b(i1-1,jj) * hu_a(i1-1,jj) * r1_hu_n(i1-1,jj) - un_b(i1-1,jj) 
     432               un_b(i1-1,jj) = un_b(i1-1,jj) + zcor 
     433               DO jk=1,jpkm1 
     434                  un(i1-1,jj,jk) = un(i1-1,jj,jk) + zcor * umask(i1-1,jj,jk) 
     435               END DO  
     436            END DO 
     437         ENDIF 
     438         ! 
     439         IF (eastern_side) THEN 
     440            DO jj=j1,j2 
     441               zcor = un_b(i2+1,jj) * hu_a(i2+1,jj) * r1_hu_n(i2+1,jj) - un_b(i2+1,jj) 
     442               un_b(i2+1,jj) = un_b(i2+1,jj) + zcor 
     443               DO jk=1,jpkm1 
     444                  un(i2+1,jj,jk) = un(i2+1,jj,jk) + zcor * umask(i2+1,jj,jk) 
     445               END DO  
     446            END DO 
     447         ENDIF 
     448         ! 
     449      ENDIF 
     450      !  
     451   END SUBROUTINE correct_u_bdy 
     452 
    401453 
    402454   SUBROUTINE updatev( tabres, i1, i2, j1, j2, k1, k2, before) 
     
    428480                  ! 
    429481                  IF (.NOT.(lk_agrif_fstep.AND.(neuler==0))) THEN ! Add asselin part 
    430                      zvb  = vb(ji,jj,jk) * e3v_b(ji,jj,jk) 
     482                     zvb  = vb(ji,jj,jk) * e3v_b(ji,jj,jk) ! fse3t_b prior update should be used 
    431483                     zvno = vn(ji,jj,jk) * e3v_a(ji,jj,jk) 
    432484                     zvnu = tabres(ji,jj,jk) 
     
    447499      !  
    448500   END SUBROUTINE updatev 
     501 
     502   SUBROUTINE correct_v_bdy( tabres, i1, i2, j1, j2, k1, k2, before, nb, ndir ) 
     503      !!--------------------------------------------- 
     504      !!           *** ROUTINE correct_u_bdy *** 
     505      !!--------------------------------------------- 
     506      INTEGER                               , INTENT(in   ) :: i1, i2, j1, j2, k1, k2 
     507      REAL(wp), DIMENSION(i1:i2,j1:j2,k1:k2), INTENT(inout) :: tabres 
     508      LOGICAL                               , INTENT(in   ) :: before 
     509      INTEGER                               , INTENT(in)    :: nb, ndir 
     510      !! 
     511      LOGICAL :: southern_side, northern_side  
     512      ! 
     513      INTEGER  :: ji, jk 
     514      REAL(wp) :: zcor 
     515      !!--------------------------------------------- 
     516      !  
     517      IF( .NOT.before ) THEN 
     518         ! 
     519         southern_side = (nb == 2).AND.(ndir == 1) 
     520         northern_side = (nb == 2).AND.(ndir == 2) 
     521         ! 
     522         IF (southern_side) THEN 
     523            DO ji=i1,i2 
     524               zcor = vn_b(ji,j1-1) * hv_a(ji,j1-1) * r1_hv_n(ji,j1-1) - vn_b(ji,j1-1) 
     525               vn_b(ji,j1-1) = vn_b(ji,j1-1) + zcor 
     526               DO jk=1,jpkm1 
     527                  vn(ji,j1-1,jk) = vn(ji,j1-1,jk) + zcor * vmask(ji,j1-1,jk) 
     528               END DO  
     529            END DO 
     530         ENDIF 
     531         ! 
     532         IF (northern_side) THEN 
     533            DO ji=i1,i2 
     534               zcor = vn_b(ji,j2+1) * hv_a(ji,j2+1) * r1_hv_n(ji,j2+1) - vn_b(ji,j2+1) 
     535               vn_b(ji,j2+1) = vn_b(ji,j2+1) + zcor 
     536               DO jk=1,jpkm1 
     537                  vn(ji,j2+1,jk) = vn(ji,j2+1,jk) + zcor * vmask(ji,j2+1,jk) 
     538               END DO  
     539            END DO 
     540         ENDIF 
     541         ! 
     542      ENDIF 
     543      !  
     544   END SUBROUTINE correct_v_bdy 
    449545 
    450546 
     
    602698         END DO 
    603699      ELSE 
    604          IF (     (.NOT.(lk_agrif_fstep.AND.(neuler==0)).AND.(ln_dynspg_exp)) & 
    605             & .OR.(.NOT.(lk_agrif_fstep.AND.(neuler==0)).AND.(ln_dynspg_ts    & 
    606             & .AND.(.NOT.ln_bt_fw)))) THEN 
    607 ! tsplit_new 
    608 !         IF (.NOT.(lk_agrif_fstep.AND.(neuler==0))) THEN 
     700         IF (.NOT.(lk_agrif_fstep.AND.(neuler==0))) THEN 
    609701            DO jj=j1,j2 
    610702               DO ji=i1,i2 
     
    662754               ub2_i_b(ji,jj) = ub2_i_b(ji,jj) + za1 * zcor  
    663755               ! Update corrective fluxes: 
    664 ! tsplit_new 
    665 !               un_bf(ji,jj)  = un_bf(ji,jj) + zcor 
     756               un_bf(ji,jj)  = un_bf(ji,jj) + zcor 
    666757               ! Update half step back fluxes: 
    667758               ub2_b(ji,jj) = tabres(ji,jj) 
     
    720811   END SUBROUTINE reflux_sshu 
    721812 
    722  
    723813   SUBROUTINE updatevb2b( tabres, i1, i2, j1, j2, before ) 
    724814      !!--------------------------------------------- 
     
    752842               vb2_i_b(ji,jj) = vb2_i_b(ji,jj) + za1 * zcor  
    753843               ! Update corrective fluxes: 
    754 ! tsplit_new 
    755 !               vn_bf(ji,jj)  = vn_bf(ji,jj) + zcor 
     844               vn_bf(ji,jj)  = vn_bf(ji,jj) + zcor 
    756845               ! Update half step back fluxes: 
    757846               vb2_b(ji,jj) = tabres(ji,jj) 
     
    9521041!         hdivn(i1:i2,j1:j2,1:jpkm1)   = e3t_b(i1:i2,j1:j2,1:jpkm1) 
    9531042 
    954          IF (     (.NOT.(lk_agrif_fstep.AND.(neuler==0)).AND.(ln_dynspg_exp)) & 
    955             & .OR.(.NOT.(lk_agrif_fstep.AND.(neuler==0)).AND.(ln_dynspg_ts    & 
    956             & .AND.(.NOT.ln_bt_fw)))) THEN 
    957 ! tsplit_new 
    958 !         IF (.NOT.(lk_agrif_fstep.AND.(neuler==0) )) THEN 
    959  
     1043         IF (.NOT.(lk_agrif_fstep.AND.(neuler==0) )) THEN 
    9601044            DO jk = 1, jpkm1 
    9611045               DO jj=j1,j2 
Note: See TracChangeset for help on using the changeset viewer.