#define SPONGE Module agrif_opa_sponge #if defined key_agrif USE par_oce USE oce USE dom_oce USE in_out_manager IMPLICIT NONE PRIVATE PUBLIC Agrif_Sponge_Tra, Agrif_Sponge_Dyn, interptn, interpsn, interpun, interpvn !! * Namelist (namagrif) REAL(wp) :: visc_tra = rdt REAL(wp) :: visc_dyn = rdt CONTAINS SUBROUTINE Agrif_Sponge_Tra( kt ) !!--------------------------------------------- !! *** ROUTINE Agrif_Sponge_Tra *** !!--------------------------------------------- #include "domzgr_substitute.h90" INTEGER, INTENT(in) :: kt INTEGER :: ji,jj,jk REAL(wp), DIMENSION(jpi,jpj,jpk) :: umasktemp,vmasktemp INTEGER :: spongearea REAL(wp) :: timecoeff REAL(wp) :: zta, zsa, zabe1, zabe2, zbtr REAL(wp), DIMENSION(jpi,jpj) :: localviscsponge REAL(wp), DIMENSION(jpi,jpj,jpk) :: ztab, tbdiff, sbdiff REAL(wp), DIMENSION(jpi,jpj,jpk) :: ztu ,ztv, zsu ,zsv #if defined SPONGE IF( kt == nit000 ) CALL agrif_sponge_init timecoeff = REAL(Agrif_NbStepint(),wp)/Agrif_rhot() Agrif_SpecialValue=0. Agrif_UseSpecialValue = .TRUE. ztab = 0.e0 CALL Agrif_Bc_Variable(ztab, ta,calledweight=timecoeff,procname=interptn) Agrif_UseSpecialValue = .FALSE. tbdiff(:,:,:) = tb(:,:,:) - ztab(:,:,:) ztab = 0.e0 Agrif_SpecialValue=0. Agrif_UseSpecialValue = .TRUE. CALL Agrif_Bc_Variable(ztab, sa,calledweight=timecoeff,procname=interpsn) Agrif_UseSpecialValue = .FALSE. sbdiff(:,:,:) = sb(:,:,:) - ztab(:,:,:) spongearea = 2 + 2 * Agrif_irhox() localviscsponge = 0. umasktemp = 0. vmasktemp = 0. IF ((nbondi == -1).OR.(nbondi == 2)) THEN DO ji = 2, spongearea localviscsponge(ji,:) = visc_tra * (spongearea-ji)/real(spongearea-2) ENDDO DO jk = 1, jpkm1 umasktemp(2:spongearea-1,:,jk) = umask(2:spongearea-1,:,jk) & * 0.5 * (localviscsponge(2:spongearea-1,:) + localviscsponge(3:spongearea,:)) ENDDO DO jk = 1, jpkm1 vmasktemp(2:spongearea,1:jpjm1,jk) = vmask(2:spongearea,1:jpjm1,jk) & * 0.5 * (localviscsponge(2:spongearea,1:jpjm1) + localviscsponge(2:spongearea,2:jpj)) ENDDO ENDIF IF ((nbondi == 1).OR.(nbondi == 2)) THEN DO ji = nlci-spongearea + 1,nlci-1 localviscsponge(ji,:) = visc_tra * (ji - (nlci-spongearea+1))/real(spongearea-2) ENDDO DO jk = 1, jpkm1 umasktemp(nlci-spongearea + 1:nlci-2,:,jk) = umask(nlci-spongearea + 1:nlci-2,:,jk) & * 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-2,:) + localviscsponge(nlci-spongearea + 2:nlci-1,:)) ENDDO DO jk = 1, jpkm1 vmasktemp(nlci-spongearea + 1:nlci-1,1:jpjm1,jk) = vmask(nlci-spongearea + 1:nlci-1,1:jpjm1,jk) & * 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-1,1:jpjm1) + localviscsponge(nlci-spongearea + 1:nlci-1,2:jpj)) ENDDO ENDIF IF ((nbondj == -1).OR.(nbondj == 2)) THEN DO jj = 2, spongearea localviscsponge(:,jj) = visc_tra * (spongearea-jj)/real(spongearea-2) ENDDO DO jk = 1, jpkm1 vmasktemp(:,2:spongearea-1,jk) = vmask(:,2:spongearea-1,jk) & * 0.5 * (localviscsponge(:,2:spongearea-1) + localviscsponge(:,3:spongearea)) ENDDO DO jk = 1, jpkm1 umasktemp(1:jpim1,2:spongearea,jk) = umask(1:jpim1,2:spongearea,jk) & * 0.5 * (localviscsponge(1:jpim1,2:spongearea) + localviscsponge(2:jpi,2:spongearea)) ENDDO ENDIF IF ((nbondj == 1).OR.(nbondj == 2)) THEN DO jj = nlcj-spongearea + 1,nlcj-1 localviscsponge(:,jj) = visc_tra * (jj - (nlcj-spongearea+1))/real(spongearea-2) ENDDO DO jk = 1, jpkm1 vmasktemp(:,nlcj-spongearea + 1:nlcj-2,jk) = vmask(:,nlcj-spongearea + 1:nlcj-2,jk) & * 0.5 * (localviscsponge(:,nlcj-spongearea + 1:nlcj-2) + localviscsponge(:,nlcj-spongearea + 2:nlcj-1)) ENDDO DO jk = 1, jpkm1 umasktemp(1:jpim1,nlcj-spongearea + 1:nlcj-1,jk) = umask(1:jpim1,nlcj-spongearea + 1:nlcj-1,jk) & * 0.5 * (localviscsponge(1:jpim1,nlcj-spongearea + 1:nlcj-1) + localviscsponge(2:jpi,nlcj-spongearea + 1:nlcj-1)) ENDDO ENDIF IF (.NOT. spongedoneT) THEN zspe1ur(:,:) = e2u(:,:) / e1u(:,:) zspe2vr(:,:) = e1v(:,:) / e2v(:,:) zspbtr2(:,:) = 1. / ( e1t(:,:) * e2t(:,:)) spongedoneT = .TRUE. ENDIF DO jk = 1, jpkm1 DO jj = 1, jpjm1 DO ji = 1, jpim1 #if defined key_zco zabe1 = umasktemp(ji,jj,jk) * zspe1ur(ji,jj) zabe2 = vmasktemp(ji,jj,jk) * zspe2vr(ji,jj) #else zabe1 = umasktemp(ji,jj,jk) * zspe1ur(ji,jj) * fse3u(ji,jj,jk) zabe2 = vmasktemp(ji,jj,jk) * zspe2vr(ji,jj) * fse3v(ji,jj,jk) #endif ztu(ji,jj,jk) = zabe1 * ( tbdiff(ji+1,jj ,jk) - tbdiff(ji,jj,jk) ) zsu(ji,jj,jk) = zabe1 * ( sbdiff(ji+1,jj ,jk) - sbdiff(ji,jj,jk) ) ztv(ji,jj,jk) = zabe2 * ( tbdiff(ji ,jj+1,jk) - tbdiff(ji,jj,jk) ) zsv(ji,jj,jk) = zabe2 * ( sbdiff(ji ,jj+1,jk) - sbdiff(ji,jj,jk) ) ENDDO ENDDO DO jj = 2,jpjm1 DO ji = 2,jpim1 #if defined key_zco zbtr = zspbtr2(ji,jj) #else zbtr = zspbtr2(ji,jj) / fse3t(ji,jj,jk) #endif ! horizontal diffusive trends zta = zbtr * ( ztu(ji,jj,jk) - ztu(ji-1,jj,jk) & & + ztv(ji,jj,jk) - ztv(ji,jj-1,jk) ) zsa = zbtr * ( zsu(ji,jj,jk) - zsu(ji-1,jj,jk) & & + zsv(ji,jj,jk) - zsv(ji,jj-1,jk) ) ! add it to the general tracer trends ta(ji,jj,jk) = (ta(ji,jj,jk) + zta) sa(ji,jj,jk) = (sa(ji,jj,jk) + zsa) END DO END DO ENDDO #endif END SUBROUTINE Agrif_Sponge_Tra SUBROUTINE Agrif_Sponge_dyn( kt ) !!--------------------------------------------- !! *** ROUTINE Agrif_Sponge_dyn *** !!--------------------------------------------- #include "domzgr_substitute.h90" INTEGER,INTENT(in) :: kt INTEGER :: ji,jj,jk INTEGER :: spongearea REAL(wp) :: timecoeff REAL(wp) :: ze2u, ze1v, zua, zva REAL(wp), DIMENSION(jpi,jpj) :: localviscsponge REAL(wp), DIMENSION(jpi,jpj,jpk) :: ztab, ubdiff, vbdiff,rotdiff,hdivdiff REAL(wp), DIMENSION(jpi,jpj,jpk) :: umasktemp,vmasktemp #if defined SPONGE timecoeff = REAL(Agrif_NbStepint(),wp)/Agrif_rhot() Agrif_SpecialValue=0. Agrif_UseSpecialValue = .TRUE. ztab = 0.e0 CALL Agrif_Bc_Variable(ztab, ua,calledweight=timecoeff,procname=interpun) Agrif_UseSpecialValue = .FALSE. ubdiff(:,:,:) = ub(:,:,:) - ztab(:,:,:) ztab = 0.e0 Agrif_SpecialValue=0. Agrif_UseSpecialValue = .TRUE. CALL Agrif_Bc_Variable(ztab, va,calledweight=timecoeff,procname=interpvn) Agrif_UseSpecialValue = .FALSE. vbdiff(:,:,:) = vb(:,:,:) - ztab(:,:,:) spongearea = 2 + 2 * Agrif_irhox() localviscsponge = 0. umasktemp = 0. vmasktemp = 0. IF ((nbondi == -1).OR.(nbondi == 2)) THEN DO ji = 2, spongearea localviscsponge(ji,:) = visc_dyn * (spongearea-ji)/real(spongearea-2) ENDDO DO jk = 1, jpkm1 umasktemp(2:spongearea-1,:,jk) = umask(2:spongearea-1,:,jk) & * 0.5 * (localviscsponge(2:spongearea-1,:) + localviscsponge(3:spongearea,:)) ENDDO DO jk = 1, jpkm1 vmasktemp(2:spongearea,1:jpjm1,jk) = vmask(2:spongearea,1:jpjm1,jk) & * 0.5 * (localviscsponge(2:spongearea,1:jpjm1) + localviscsponge(2:spongearea,2:jpj)) ENDDO ENDIF IF ((nbondi == 1).OR.(nbondi == 2)) THEN DO ji = nlci-spongearea + 1,nlci-1 localviscsponge(ji,:) = visc_dyn * (ji - (nlci-spongearea+1))/real(spongearea-2) ENDDO DO jk = 1, jpkm1 umasktemp(nlci-spongearea + 1:nlci-2,:,jk) = umask(nlci-spongearea + 1:nlci-2,:,jk) & * 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-2,:) + localviscsponge(nlci-spongearea + 2:nlci-1,:)) ENDDO DO jk = 1, jpkm1 vmasktemp(nlci-spongearea + 1:nlci-1,1:jpjm1,jk) = vmask(nlci-spongearea + 1:nlci-1,1:jpjm1,jk) & * 0.5 * (localviscsponge(nlci-spongearea + 1:nlci-1,1:jpjm1) + localviscsponge(nlci-spongearea + 1:nlci-1,2:jpj)) ENDDO ENDIF IF ((nbondj == -1).OR.(nbondj == 2)) THEN DO jj = 2, spongearea localviscsponge(:,jj) = visc_dyn * (spongearea-jj)/real(spongearea-2) ENDDO DO jk = 1, jpkm1 vmasktemp(:,2:spongearea-1,jk) = vmask(:,2:spongearea-1,jk) & * 0.5 * (localviscsponge(:,2:spongearea-1) + localviscsponge(:,3:spongearea)) ENDDO DO jk = 1, jpkm1 umasktemp(1:jpim1,2:spongearea,jk) = umask(1:jpim1,2:spongearea,jk) & * 0.5 * (localviscsponge(1:jpim1,2:spongearea) + localviscsponge(2:jpi,2:spongearea)) ENDDO ENDIF IF ((nbondj == 1).OR.(nbondj == 2)) THEN DO jj = nlcj-spongearea + 1,nlcj-1 localviscsponge(:,jj) = visc_dyn * (jj - (nlcj-spongearea+1))/real(spongearea-2) ENDDO DO jk = 1, jpkm1 vmasktemp(:,nlcj-spongearea + 1:nlcj-2,jk) = vmask(:,nlcj-spongearea + 1:nlcj-2,jk) & * 0.5 * (localviscsponge(:,nlcj-spongearea + 1:nlcj-2) + localviscsponge(:,nlcj-spongearea + 2:nlcj-1)) ENDDO DO jk = 1, jpkm1 umasktemp(1:jpim1,nlcj-spongearea + 1:nlcj-1,jk) = umask(1:jpim1,nlcj-spongearea + 1:nlcj-1,jk) & * 0.5 * (localviscsponge(1:jpim1,nlcj-spongearea + 1:nlcj-1) + localviscsponge(2:jpi,nlcj-spongearea + 1:nlcj-1)) ENDDO ENDIF ubdiff = ubdiff * umasktemp vbdiff = vbdiff * vmasktemp hdivdiff = 0. rotdiff = 0. DO jk = 1, jpkm1 ! Horizontal slab ! ! =============== ! ! -------- ! Horizontal divergence ! div ! ! -------- DO jj = 2, jpjm1 DO ji = 2, jpim1 ! vector opt. #if defined key_zco hdivdiff(ji,jj,jk) = ( e2u(ji,jj) * ubdiff(ji,jj,jk) & - e2u(ji-1,jj ) * ubdiff(ji-1,jj ,jk) & & + e1v(ji,jj) * vbdiff(ji,jj,jk) - & & e1v(ji ,jj-1) * vbdiff(ji ,jj-1,jk) ) & & / ( e1t(ji,jj) * e2t(ji,jj) ) #else hdivdiff(ji,jj,jk) = & ( e2u(ji,jj)*fse3u(ji,jj,jk) * & ubdiff(ji,jj,jk) - e2u(ji-1,jj )* & fse3u(ji-1,jj ,jk) * ubdiff(ji-1,jj ,jk) & + e1v(ji,jj)*fse3v(ji,jj,jk) * & vbdiff(ji,jj,jk) - e1v(ji ,jj-1)* & fse3v(ji ,jj-1,jk) * vbdiff(ji ,jj-1,jk) ) & / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) #endif END DO END DO DO jj = 1, jpjm1 DO ji = 1, jpim1 ! vector opt. rotdiff(ji,jj,jk) = ( e2v(ji+1,jj ) * vbdiff(ji+1,jj ,jk) - e2v(ji,jj) * vbdiff(ji,jj,jk) & & - e1u(ji ,jj+1) * ubdiff(ji ,jj+1,jk) + e1u(ji,jj) * ubdiff(ji,jj,jk) ) & & * fmask(ji,jj,jk) / ( e1f(ji,jj) * e2f(ji,jj) ) END DO END DO ENDDO ! ! =============== DO jk = 1, jpkm1 ! Horizontal slab ! ! =============== DO jj = 2, jpjm1 DO ji = 2, jpim1 ! vector opt. #if defined key_zco ! horizontal diffusive trends ze2u = rotdiff (ji,jj,jk) ze1v = hdivdiff(ji,jj,jk) zua = - ( ze2u - & rotdiff (ji,jj-1,jk) ) / e2u(ji,jj) & + ( hdivdiff(ji+1,jj,jk) - & ze1v ) / e1u(ji,jj) zva = + ( ze2u - & rotdiff (ji-1,jj,jk) ) / e1v(ji,jj) & + ( hdivdiff(ji,jj+1,jk) - & ze1v ) / e2v(ji,jj) #else ze2u = rotdiff (ji,jj,jk)*fse3f(ji,jj,jk) ze1v = hdivdiff(ji,jj,jk) ! horizontal diffusive trends zua = - ( ze2u - rotdiff (ji,jj-1,jk)* & fse3f(ji,jj-1,jk) ) / ( e2u(ji,jj) * fse3u(ji,jj,jk) ) & + ( hdivdiff(ji+1,jj,jk) - ze1v & ) / e1u(ji,jj) zva = + ( ze2u - rotdiff (ji-1,jj,jk)* & fse3f(ji-1,jj,jk) ) / ( e1v(ji,jj) * fse3v(ji,jj,jk) ) & + ( hdivdiff(ji,jj+1,jk) - ze1v & ) / e2v(ji,jj) #endif ! add it to the general momentum trends ua(ji,jj,jk) = ua(ji,jj,jk) + zua va(ji,jj,jk) = va(ji,jj,jk) + zva END DO END DO ! ! =============== END DO ! End of slab ! ! =============== #endif END SUBROUTINE Agrif_Sponge_dyn SUBROUTINE agrif_sponge_init !!--------------------------------------------- !! *** ROUTINE agrif_sponge_init *** !!--------------------------------------------- NAMELIST/namagrif/ visc_tra, visc_dyn REWIND ( numnam ) READ ( numnam, namagrif ) IF(lwp) THEN WRITE(numout,*) WRITE(numout,*) 'agrif_sponge_init : agrif sponge parameters' WRITE(numout,*) '~~~~~~~~~~~~' WRITE(numout,*) ' Namelist namagrif : set sponge parameters' WRITE(numout,*) ' sponge coefficient for tracers = ', visc_tra WRITE(numout,*) ' sponge coefficient for dynamics = ', visc_dyn ENDIF END SUBROUTINE agrif_sponge_init SUBROUTINE interptn(tabres,i1,i2,j1,j2,k1,k2) !!--------------------------------------------- !! *** ROUTINE interptn *** !!--------------------------------------------- # include "domzgr_substitute.h90" INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2 REAL(wp), DIMENSION(i1:i2,j1:j2,k1:k2), INTENT(inout) :: tabres tabres(i1:i2,j1:j2,k1:k2) = tn(i1:i2,j1:j2,k1:k2) END SUBROUTINE interptn SUBROUTINE interpsn(tabres,i1,i2,j1,j2,k1,k2) !!--------------------------------------------- !! *** ROUTINE interpsn *** !!--------------------------------------------- # include "domzgr_substitute.h90" INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2 REAL(wp), DIMENSION(i1:i2,j1:j2,k1:k2), INTENT(inout) :: tabres tabres(i1:i2,j1:j2,k1:k2) = sn(i1:i2,j1:j2,k1:k2) END SUBROUTINE interpsn SUBROUTINE interpun(tabres,i1,i2,j1,j2,k1,k2) !!--------------------------------------------- !! *** ROUTINE interpun *** !!--------------------------------------------- # include "domzgr_substitute.h90" INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2 REAL(wp), DIMENSION(i1:i2,j1:j2,k1:k2), INTENT(inout) :: tabres tabres(i1:i2,j1:j2,k1:k2) = un(i1:i2,j1:j2,k1:k2) END SUBROUTINE interpun SUBROUTINE interpvn(tabres,i1,i2,j1,j2,k1,k2) !!--------------------------------------------- !! *** ROUTINE interpvn *** !!--------------------------------------------- # include "domzgr_substitute.h90" INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2 REAL(wp), DIMENSION(i1:i2,j1:j2,k1:k2), INTENT(inout) :: tabres tabres(i1:i2,j1:j2,k1:k2) = vn(i1:i2,j1:j2,k1:k2) END SUBROUTINE interpvn #else CONTAINS SUBROUTINE agrif_opa_sponge_empty !!--------------------------------------------- !! *** ROUTINE agrif_OPA_sponge_empty *** !!--------------------------------------------- WRITE(*,*) 'agrif_opa_sponge : You should not have seen this print! error?' END SUBROUTINE agrif_opa_sponge_empty #endif END MODULE agrif_opa_sponge