MODULE ad_zetabc_mod ! !git $Id$ !svn $Id: ad_zetabc.F 1151 2023-02-09 03:08:53Z arango $ !================================================== Hernan G. Arango === ! Copyright (c) 2002-2023 The ROMS/TOMS Group Andrew M. Moore ! ! Licensed under a MIT/X style license ! ! See License_ROMS.md ! !======================================================================= ! ! ! This subroutine sets adjoint lateral boundary conditions for ! ! free-surface. It updates the specified "kout" index. ! ! ! ! BASIC STATE variables needed: zeta ! ! ! !======================================================================= ! implicit none PRIVATE PUBLIC :: ad_zetabc, ad_zetabc_tile CONTAINS ! !*********************************************************************** SUBROUTINE ad_zetabc (ng, tile, kout) !*********************************************************************** ! USE mod_param USE mod_ocean USE mod_stepping ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, kout ! ! Local variable declarations. ! integer :: IminS, ImaxS, JminS, JmaxS integer :: LBi, UBi, LBj, UBj, LBij, UBij ! ! Set horizontal starting and ending indices for automatic private ! storage arrays. ! IminS=BOUNDS(ng)%Istr(tile)-3 ImaxS=BOUNDS(ng)%Iend(tile)+3 JminS=BOUNDS(ng)%Jstr(tile)-3 JmaxS=BOUNDS(ng)%Jend(tile)+3 ! ! Determine array lower and upper bounds in the I- and J-directions. ! LBi=BOUNDS(ng)%LBi(tile) UBi=BOUNDS(ng)%UBi(tile) LBj=BOUNDS(ng)%LBj(tile) UBj=BOUNDS(ng)%UBj(tile) ! ! Set array lower and upper bounds for MIN(I,J) directions and ! MAX(I,J) directions. ! LBij=BOUNDS(ng)%LBij UBij=BOUNDS(ng)%UBij ! CALL ad_zetabc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs(ng), kstp(ng), kout, & & OCEAN(ng) % zeta, & & OCEAN(ng) % ad_zeta) RETURN END SUBROUTINE ad_zetabc ! !*********************************************************************** SUBROUTINE ad_zetabc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, kout, & & zeta, ad_zeta) !*********************************************************************** ! USE mod_param USE mod_boundary USE mod_grid USE mod_ncparam USE mod_scalars ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile integer, intent(in) :: LBi, UBi, LBj, UBj integer, intent(in) :: IminS, ImaxS, JminS, JmaxS integer, intent(in) :: krhs, kstp, kout ! real(r8), intent(in) :: zeta(LBi:,LBj:,:) real(r8), intent(inout) :: ad_zeta(LBi:,LBj:,:) ! ! Local variable declarations. ! integer :: i, j, know real(r8) :: Ce, Cx real(r8) :: cff, cff1, cff2, dt2d, tau real(r8) :: ad_Ce, ad_Cx real(r8) :: ad_cff1, ad_cff2 real(r8) :: adfac real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_grad ! !----------------------------------------------------------------------- ! Set lower and upper tile bounds and staggered variables bounds for ! this horizontal domain partition. Notice that if tile=-1, it will ! set the values for the global grid. !----------------------------------------------------------------------- ! integer :: Istr, IstrB, IstrP, IstrR, IstrT, IstrM, IstrU integer :: Iend, IendB, IendP, IendR, IendT integer :: Jstr, JstrB, JstrP, JstrR, JstrT, JstrM, JstrV integer :: Jend, JendB, JendP, JendR, JendT integer :: Istrm3, Istrm2, Istrm1, IstrUm2, IstrUm1 integer :: Iendp1, Iendp2, Iendp2i, Iendp3 integer :: Jstrm3, Jstrm2, Jstrm1, JstrVm2, JstrVm1 integer :: Jendp1, Jendp2, Jendp2i, Jendp3 ! Istr =BOUNDS(ng) % Istr (tile) IstrB =BOUNDS(ng) % IstrB (tile) IstrM =BOUNDS(ng) % IstrM (tile) IstrP =BOUNDS(ng) % IstrP (tile) IstrR =BOUNDS(ng) % IstrR (tile) IstrT =BOUNDS(ng) % IstrT (tile) IstrU =BOUNDS(ng) % IstrU (tile) Iend =BOUNDS(ng) % Iend (tile) IendB =BOUNDS(ng) % IendB (tile) IendP =BOUNDS(ng) % IendP (tile) IendR =BOUNDS(ng) % IendR (tile) IendT =BOUNDS(ng) % IendT (tile) Jstr =BOUNDS(ng) % Jstr (tile) JstrB =BOUNDS(ng) % JstrB (tile) JstrM =BOUNDS(ng) % JstrM (tile) JstrP =BOUNDS(ng) % JstrP (tile) JstrR =BOUNDS(ng) % JstrR (tile) JstrT =BOUNDS(ng) % JstrT (tile) JstrV =BOUNDS(ng) % JstrV (tile) Jend =BOUNDS(ng) % Jend (tile) JendB =BOUNDS(ng) % JendB (tile) JendP =BOUNDS(ng) % JendP (tile) JendR =BOUNDS(ng) % JendR (tile) JendT =BOUNDS(ng) % JendT (tile) ! Istrm3 =BOUNDS(ng) % Istrm3 (tile) ! Istr-3 Istrm2 =BOUNDS(ng) % Istrm2 (tile) ! Istr-2 Istrm1 =BOUNDS(ng) % Istrm1 (tile) ! Istr-1 IstrUm2=BOUNDS(ng) % IstrUm2(tile) ! IstrU-2 IstrUm1=BOUNDS(ng) % IstrUm1(tile) ! IstrU-1 Iendp1 =BOUNDS(ng) % Iendp1 (tile) ! Iend+1 Iendp2 =BOUNDS(ng) % Iendp2 (tile) ! Iend+2 Iendp2i=BOUNDS(ng) % Iendp2i(tile) ! Iend+2 interior Iendp3 =BOUNDS(ng) % Iendp3 (tile) ! Iend+3 Jstrm3 =BOUNDS(ng) % Jstrm3 (tile) ! Jstr-3 Jstrm2 =BOUNDS(ng) % Jstrm2 (tile) ! Jstr-2 Jstrm1 =BOUNDS(ng) % Jstrm1 (tile) ! Jstr-1 JstrVm2=BOUNDS(ng) % JstrVm2(tile) ! JstrV-2 JstrVm1=BOUNDS(ng) % JstrVm1(tile) ! JstrV-1 Jendp1 =BOUNDS(ng) % Jendp1 (tile) ! Jend+1 Jendp2 =BOUNDS(ng) % Jendp2 (tile) ! Jend+2 Jendp2i=BOUNDS(ng) % Jendp2i(tile) ! Jend+2 interior Jendp3 =BOUNDS(ng) % Jendp3 (tile) ! Jend+3 ! !----------------------------------------------------------------------- ! Initialize adjoint private variables. !----------------------------------------------------------------------- ! ad_Ce=0.0_r8 ad_Cx=0.0_r8 ad_cff1=0.0_r8 ad_cff2=0.0_r8 ad_grad(LBi:UBi,LBj:UBj)=0.0_r8 ! !----------------------------------------------------------------------- ! Set time-indices !----------------------------------------------------------------------- ! IF (iif(ng).eq.1) THEN know=krhs dt2d=dtfast(ng) ELSE IF (PREDICTOR_2D_STEP(ng)) THEN know=krhs dt2d=2.0_r8*dtfast(ng) ELSE know=kstp dt2d=dtfast(ng) END IF ! !----------------------------------------------------------------------- ! Boundary corners. !----------------------------------------------------------------------- ! IF (.not.(EWperiodic(ng).or.NSperiodic(ng))) THEN IF (DOMAIN(ng)%NorthEast_Corner(tile)) THEN IF (LBC_apply(ng)%north(Iend+1).and. & & LBC_apply(ng)%east (Jend+1)) THEN !^ tl_zeta(Iend+1,Jend+1,kout)=0.5_r8* & !^ & (tl_zeta(Iend+1,Jend ,kout)+ & !^ & tl_zeta(Iend ,Jend+1,kout)) !^ adfac=0.5_r8*ad_zeta(Iend+1,Jend+1,kout) ad_zeta(Iend ,Jend+1,kout)=ad_zeta(Iend ,Jend+1,kout)+ & & adfac ad_zeta(Iend+1,Jend ,kout)=ad_zeta(Iend+1,Jend ,kout)+ & & adfac ad_zeta(Iend+1,Jend+1,kout)=0.0_r8 END IF END IF IF (DOMAIN(ng)%NorthWest_Corner(tile)) THEN IF (LBC_apply(ng)%north(Istr-1).and. & & LBC_apply(ng)%west (Jend+1)) THEN !^ tl_zeta(Istr-1,Jend+1,kout)=0.5_r8* & !^ & (tl_zeta(Istr-1,Jend ,kout)+ & !^ & tl_zeta(Istr ,Jend+1,kout)) !^ adfac=0.5_r8*ad_zeta(Istr-1,Jend+1,kout) ad_zeta(Istr-1,Jend ,kout)=ad_zeta(Istr-1,Jend ,kout)+ & & adfac ad_zeta(Istr ,Jend+1,kout)=ad_zeta(Istr ,Jend+1,kout)+ & & adfac ad_zeta(Istr-1,Jend+1,kout)=0.0_r8 END IF END IF IF (DOMAIN(ng)%SouthEast_Corner(tile)) THEN IF (LBC_apply(ng)%south(Iend+1).and. & & LBC_apply(ng)%east (Jstr-1)) THEN !^ tl_zeta(Iend+1,Jstr-1,kout)=0.5_r8* & !^ & (tl_zeta(Iend ,Jstr-1,kout)+ & !^ & tl_zeta(Iend+1,Jstr ,kout)) !^ adfac=0.5_r8*ad_zeta(Iend+1,Jstr-1,kout) ad_zeta(Iend ,Jstr-1,kout)=ad_zeta(Iend ,Jstr-1,kout)+ & & adfac ad_zeta(Iend+1,Jstr ,kout)=ad_zeta(Iend+1,Jstr ,kout)+ & & adfac ad_zeta(Iend+1,Jstr-1,kout)=0.0_r8 END IF END IF IF (DOMAIN(ng)%SouthWest_Corner(tile)) THEN IF (LBC_apply(ng)%south(Istr-1).and. & & LBC_apply(ng)%west (Jstr-1)) THEN !^ tl_zeta(Istr-1,Jstr-1,kout)=0.5_r8* & !^ & (tl_zeta(Istr ,Jstr-1,kout)+ & !^ & tl_zeta(Istr-1,Jstr ,kout)) !^ adfac=0.5_r8*ad_zeta(Istr-1,Jstr-1,kout) ad_zeta(Istr ,Jstr-1,kout)=ad_zeta(Istr ,Jstr-1,kout)+ & & adfac ad_zeta(Istr-1,Jstr ,kout)=ad_zeta(Istr-1,Jstr ,kout)+ & & adfac ad_zeta(Istr-1,Jstr-1,kout)=0.0_r8 END IF END IF END IF ! !----------------------------------------------------------------------- ! Lateral boundary conditions at the northern edge. !----------------------------------------------------------------------- ! IF (DOMAIN(ng)%Northern_Edge(tile)) THEN ! ! Northern edge, implicit upstream radiation condition. ! IF (ad_LBC(inorth,isFsur,ng)%radiation) THEN IF (iic(ng).ne.0) THEN DO i=Istr,Iend IF (LBC_apply(ng)%north(i)) THEN !^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend+1,kout)* & !^ & GRID(ng)%rmask(i,Jend+1) !^ ad_zeta(i,Jend+1,kout)=ad_zeta(i,Jend+1,kout)* & & GRID(ng)%rmask(i,Jend+1) IF (ad_LBC(inorth,isFsur,ng)%nudging) THEN !^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend+1,kout)- & !^ & tau*tl_zeta(i,Jend+1,know) !^ ad_zeta(i,Jend+1,know)=ad_zeta(i,Jend+1,know)- & & tau*ad_zeta(i,Jend+1,kout) END IF !^ tl_zeta(i,Jend+1,kout)=(cff*tl_zeta(i,Jend+1,know)+ & !^ & Ce *tl_zeta(i,Jend ,kout)- & !^ & MAX(Cx,0.0_r8)* & !^ & tl_grad(i ,Jend+1)- & !^ & MIN(Cx,0.0_r8)* & !^ & tl_grad(i+1,Jend+1))/ & !^ & (cff+Ce) !^ adfac=ad_zeta(i,Jend+1,kout)/(cff+Ce) ad_grad(i ,Jend+1)=ad_grad(i ,Jend+1)- & & MAX(Cx,0.0_r8)*adfac ad_grad(i+1,Jend+1)=ad_grad(i+1,Jend+1)- & & MIN(Cx,0.0_r8)*adfac ad_zeta(i,Jend ,kout)=ad_zeta(i,Jend ,kout)+Ce *adfac ad_zeta(i,Jend+1,know)=ad_zeta(i,Jend+1,know)+cff*adfac ad_zeta(i,Jend+1,kout)=0.0_r8 END IF END DO END IF ! ! Northern edge, explicit Chapman boundary condition. ! ELSE IF (ad_LBC(inorth,isFsur,ng)%Chapman_explicit) THEN DO i=Istr,Iend IF (LBC_apply(ng)%north(i)) THEN cff=dt2d*GRID(ng)%pn(i,Jend) cff1=SQRT(g*(GRID(ng)%h(i,Jend)+ & & zeta(i,Jend,know))) Ce=cff*cff1 !^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend+1,kout)* & !^ & GRID(ng)%rmask(i,Jend+1) !^ ad_zeta(i,Jend+1,kout)=ad_zeta(i,Jend+1,kout)* & & GRID(ng)%rmask(i,Jend+1) !^ tl_zeta(i,Jend+1,kout)=(1.0_r8-Ce)*tl_zeta(i,Jend+1,know)+& !^ & tl_Ce*(zeta(i,Jend+1,know)+ & !^ & zeta(i,Jend ,know))+ & !^ & Ce*tl_zeta(i,Jend,know) !^ ad_zeta(i,Jend+1,know)=ad_zeta(i,Jend+1,know)+ & & (1.0_r8-Ce)*ad_zeta(i,Jend+1,kout) ad_Ce=ad_Ce+(zeta(i,Jend+1,know)+ & & zeta(i,Jend ,know))*ad_zeta(i,Jend+1,kout) ad_zeta(i,Jend,know)=ad_zeta(i,Jend,know)+ & & Ce*ad_zeta(i,Jend+1,kout) ad_zeta(i,Jend+1,kout)=0.0_r8 !^ tl_Ce=cff*tl_cff1 !^ ad_cff1=ad_cff1+cff*ad_Ce ad_Ce=0.0_r8 !^ tl_cff1=0.5_r8*g*(GRID(ng)%tl_h(i,Jend)+ & !^ & tl_zeta(i,Jend,know))/cff1 !^ adfac=0.5_r8*g*ad_cff1/cff1 GRID(ng)%ad_h(i,Jend)=GRID(ng)%ad_h(i,Jend)+adfac ad_zeta(i,Jend,know)=ad_zeta(i,Jend,know)+adfac ad_cff1=0.0_r8 END IF END DO ! ! Northern edge, implicit Chapman boundary condition. ! ELSE IF (ad_LBC(inorth,isFsur,ng)%Chapman_implicit) THEN DO i=Istr,Iend IF (LBC_apply(ng)%north(i)) THEN cff=dt2d*GRID(ng)%pn(i,Jend) cff1=SQRT(g*(GRID(ng)%h(i,Jend)+ & & zeta(i,Jend,know))) Ce=cff*cff1 cff2=1.0_r8/(1.0_r8+Ce) !^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend+1,kout)* & !^ & GRID(ng)%rmask(i,Jend+1) !^ ad_zeta(i,Jend+1,kout)=ad_zeta(i,Jend+1,kout)* & & GRID(ng)%rmask(i,Jend+1) !^ tl_zeta(i,Jend+1,kout)=tl_cff2*(zeta(i,Jend+1,know)+ & !^ & Ce*zeta(i,Jend,kout))+ & !^ & cff2*(tl_zeta(i,Jend+1,know)+ & !^ & tl_Ce*zeta(i,Jend,kout)+ & !^ & Ce*tl_zeta(i,Jend,kout)) !^ adfac=cff2*ad_zeta(i,Jend+1,kout) ad_zeta(i,Jend ,kout)=ad_zeta(i,Jend ,kout)+Ce*adfac ad_zeta(i,Jend+1,know)=ad_zeta(i,Jend+1,know)+adfac ad_Ce=ad_Ce+zeta(i,Jend,kout)*adfac ad_cff2=ad_cff2+ & & (zeta(i,Jend+1,know)+ & & Ce*zeta(i,Jend,kout))*ad_zeta(i,Jend+1,kout) ad_zeta(i,Jend+1,kout)=0.0_r8 !^ tl_cff2=-cff2*cff2*tl_Ce !^ ad_Ce=ad_Ce-cff2*cff2*ad_cff2 ad_cff2=0.0_r8 !^ tl_Ce=cff*tl_cff1 !^ ad_cff1=ad_cff1+cff*ad_Ce ad_Ce=0.0_r8 !^ tl_cff1=0.5_r8*g*(GRID(ng)%tl_h(i,Jend)+ & !^ & tl_zeta(i,Jend,know))/cff1 !^ adfac=0.5_r8*g*ad_cff1/cff1 GRID(ng)%ad_h(i,Jend)=GRID(ng)%ad_h(i,Jend)+adfac ad_zeta(i,Jend,know)=ad_zeta(i,Jend,know)+adfac ad_cff1=0.0_r8 END IF END DO ! ! Northern edge, clamped boundary condition. ! ELSE IF (ad_LBC(inorth,isFsur,ng)%clamped) THEN DO i=Istr,Iend IF (LBC_apply(ng)%north(i)) THEN !^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend+1,kout)* & !^ & GRID(ng)%rmask(i,Jend+1) !^ ad_zeta(i,Jend+1,kout)=ad_zeta(i,Jend+1,kout)* & & GRID(ng)%rmask(i,Jend+1) !^ tl_zeta(i,Jend+1,kout)=0.0_r8 !^ ad_zeta(i,Jend+1,kout)=0.0_r8 END IF END DO ! ! Northern edge, gradient boundary condition. ! ELSE IF (ad_LBC(inorth,isFsur,ng)%gradient) THEN DO i=Istr,Iend IF (LBC_apply(ng)%north(i)) THEN !^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend+1,kout)* & !^ & GRID(ng)%rmask(i,Jend+1) !^ ad_zeta(i,Jend+1,kout)=ad_zeta(i,Jend+1,kout)* & & GRID(ng)%rmask(i,Jend+1) !^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend,kout) !^ ad_zeta(i,Jend ,kout)=ad_zeta(i,Jend ,kout)+ & & ad_zeta(i,Jend+1,kout) ad_zeta(i,Jend+1,kout)=0.0_r8 END IF END DO ! ! Northern edge, closed boundary condition. ! ELSE IF (ad_LBC(inorth,isFsur,ng)%closed) THEN DO i=Istr,Iend IF (LBC_apply(ng)%north(i)) THEN !^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend+1,kout)* & !^ & GRID(ng)%rmask(i,Jend+1) !^ ad_zeta(i,Jend+1,kout)=ad_zeta(i,Jend+1,kout)* & & GRID(ng)%rmask(i,Jend+1) !^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend,kout) !^ ad_zeta(i,Jend ,kout)=ad_zeta(i,Jend ,kout)+ & & ad_zeta(i,Jend+1,kout) ad_zeta(i,Jend+1,kout)=0.0_r8 END IF END DO END IF END IF ! !----------------------------------------------------------------------- ! Lateral boundary conditions at the southern edge. !----------------------------------------------------------------------- ! IF (DOMAIN(ng)%Southern_Edge(tile)) THEN ! ! Southern edge, implicit upstream radiation condition. ! IF (ad_LBC(isouth,isFsur,ng)%radiation) THEN IF (iic(ng).ne.0) THEN DO i=Istr,Iend IF (LBC_apply(ng)%south(i)) THEN !^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr-1,kout)* & !^ & GRID(ng)%rmask(i,Jstr-1) !^ ad_zeta(i,Jstr-1,kout)=ad_zeta(i,Jstr-1,kout)* & & GRID(ng)%rmask(i,Jstr-1) IF (ad_LBC(isouth,isFsur,ng)%nudging) THEN !^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr-1,kout)- & !^ & tau*tl_zeta(i,Jstr-1,know) !^ ad_zeta(i,Jstr-1,know)=ad_zeta(i,Jstr-1,know)- & & tau*ad_zeta(i,Jstr-1,kout) END IF !^ tl_zeta(i,Jstr-1,kout)=(cff*tl_zeta(i,Jstr-1,know)+ & !^ & Ce *tl_zeta(i,Jstr ,kout)- & !^ & MAX(Cx,0.0_r8)* & !^ & tl_grad(i ,Jstr-1)- & !^ & MIN(Cx,0.0_r8)* & !^ & tl_grad(i+1,Jstr-1))/ & !^ & (cff+Ce) !^ adfac=ad_zeta(i,Jstr-1,kout)/(cff+Ce) ad_grad(i ,Jstr-1)=ad_grad(i ,Jstr-1)- & & MAX(Cx,0.0_r8)*adfac ad_grad(i+1,Jstr-1)=ad_grad(i+1,Jstr-1)- & & MIN(Cx,0.0_r8)*adfac ad_zeta(i,Jstr-1,know)=ad_zeta(i,Jstr-1,know)+cff*adfac ad_zeta(i,Jstr ,kout)=ad_zeta(i,Jstr ,kout)+Ce *adfac ad_zeta(i,Jstr-1,kout)=0.0_r8 END IF END DO END IF ! ! Southern edge, explicit Chapman boundary condition. ! ELSE IF (ad_LBC(isouth,isFsur,ng)%Chapman_explicit) THEN DO i=Istr,Iend IF (LBC_apply(ng)%south(i)) THEN cff=dt2d*GRID(ng)%pn(i,Jstr) cff1=SQRT(g*(GRID(ng)%h(i,Jstr)+ & & zeta(i,Jstr,know))) Ce=cff*cff1 !^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr-1,kout)* & !^ & GRID(ng)%rmask(i,Jstr-1) !^ ad_zeta(i,Jstr-1,kout)=ad_zeta(i,Jstr-1,kout)* & & GRID(ng)%rmask(i,Jstr-1) !^ tl_zeta(i,Jstr-1,kout)=(1.0_r8-Ce)*tl_zeta(i,Jstr-1,know)+& !^ & tl_Ce*(zeta(i,Jstr-1,know)+ & !^ & zeta(i,Jstr ,know))+ & !^ & Ce*tl_zeta(i,Jstr,know) !^ ad_zeta(i,Jstr-1,know)=ad_zeta(i,Jstr-1,know)+ & & (1.0_r8-Ce)*ad_zeta(i,Jstr-1,kout) ad_Ce=ad_Ce+(zeta(i,Jstr-1,know)+ & & zeta(i,Jstr ,know))*ad_zeta(i,Jstr-1,kout) ad_zeta(i,Jstr,know)=ad_zeta(i,Jstr,know)+ & & Ce*ad_zeta(i,Jstr-1,kout) ad_zeta(i,Jstr-1,kout)=0.0_r8 !^ tl_Ce=cff*tl_cff1 !^ ad_cff1=ad_cff1+cff*ad_Ce ad_Ce=0.0_r8 !^ tl_cff1=0.5_r8*g*(GRID(ng)%tl_h(i,Jstr)+ & !^ & tl_zeta(i,Jstr,know))/cff1 !^ adfac=0.5_r8*g*ad_cff1/cff1 GRID(ng)%ad_h(i,Jstr)=GRID(ng)%ad_h(i,Jstr)+adfac ad_zeta(i,Jstr,know)=ad_zeta(i,Jstr,know)+adfac ad_cff1=0.0_r8 END IF END DO ! ! Southern edge, implicit Chapman boundary condition. ! ELSE IF (ad_LBC(isouth,isFsur,ng)%Chapman_implicit) THEN DO i=Istr,Iend IF (LBC_apply(ng)%south(i)) THEN cff=dt2d*GRID(ng)%pn(i,Jstr) cff1=SQRT(g*(GRID(ng)%h(i,Jstr)+ & & zeta(i,Jstr,know))) Ce=cff*cff1 cff2=1.0_r8/(1.0_r8+Ce) !^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr-1,kout)* & !^ & GRID(ng)%rmask(i,Jstr-1) !^ ad_zeta(i,Jstr-1,kout)=ad_zeta(i,Jstr-1,kout)* & & GRID(ng)%rmask(i,Jstr-1) !^ tl_zeta(i,Jstr-1,kout)=tl_cff2*(zeta(i,Jstr-1,know)+ & !^ & Ce*zeta(i,Jstr,kout))+ & !^ & cff2*(tl_zeta(i,Jstr-1,know)+ & !^ & tl_Ce*zeta(i,Jstr,kout)+ & !^ & Ce*tl_zeta(i,Jstr,kout)) !^ adfac=cff2*ad_zeta(i,Jstr-1,kout) ad_zeta(i,Jstr-1,know)=ad_zeta(i,Jstr-1,know)+adfac ad_zeta(i,Jstr ,kout)=ad_zeta(i,Jstr ,kout)+Ce*adfac ad_Ce=ad_Ce+zeta(i,Jstr,kout)*adfac ad_cff2=ad_cff2+ & & (zeta(i,Jstr-1,know)+ & & Ce*zeta(i,Jstr,kout))*ad_zeta(i,Jstr-1,kout) ad_zeta(i,Jstr-1,kout)=0.0_r8 !^ tl_cff2=-cff2*cff2*tl_Ce !^ ad_Ce=ad_Ce-cff2*cff2*ad_cff2 ad_cff2=0.0_r8 !^ tl_Ce=cff*tl_cff1 !^ ad_cff1=ad_cff1+cff*ad_Ce ad_Ce=0.0_r8 !^ tl_cff1=0.5_r8*g*(GRID(ng)%tl_h(i,Jstr)+ & !^ & tl_zeta(i,Jstr,know))/cff1 !^ adfac=0.5_r8*g*ad_cff1/cff1 GRID(ng)%ad_h(i,Jstr)=GRID(ng)%ad_h(i,Jstr)+adfac ad_zeta(i,Jstr,know)=ad_zeta(i,Jstr,know)+adfac ad_cff1=0.0_r8 END IF END DO ! ! Southern edge, clamped boundary condition. ! ELSE IF (ad_LBC(isouth,isFsur,ng)%clamped) THEN DO i=Istr,Iend IF (LBC_apply(ng)%south(i)) THEN !^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr-1,kout)* & !^ & GRID(ng)%rmask(i,Jstr-1) !^ ad_zeta(i,Jstr-1,kout)=ad_zeta(i,Jstr-1,kout)* & & GRID(ng)%rmask(i,Jstr-1) !^ tl_zeta(i,Jstr-1,kout)=0.0_r8 !^ ad_zeta(i,Jstr-1,kout)=0.0_r8 END IF END DO ! ! Southern edge, gradient boundary condition. ! ELSE IF (ad_LBC(isouth,isFsur,ng)%gradient) THEN DO i=Istr,Iend IF (LBC_apply(ng)%south(i)) THEN !^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr-1,kout)* & !^ & GRID(ng)%rmask(i,Jstr-1) !^ ad_zeta(i,Jstr-1,kout)=ad_zeta(i,Jstr-1,kout)* & & GRID(ng)%rmask(i,Jstr-1) !^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr,kout) !^ ad_zeta(i,Jstr ,kout)=ad_zeta(i,Jstr ,kout)+ & & ad_zeta(i,Jstr-1,kout) ad_zeta(i,Jstr-1,kout)=0.0_r8 END IF END DO ! ! Southern edge, closed boundary condition. ! ELSE IF (ad_LBC(isouth,isFsur,ng)%closed) THEN DO i=Istr,Iend IF (LBC_apply(ng)%south(i)) THEN !^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr-1,kout)* & !^ & GRID(ng)%rmask(i,Jstr-1) !^ ad_zeta(i,Jstr-1,kout)=ad_zeta(i,Jstr-1,kout)* & & GRID(ng)%rmask(i,Jstr-1) !^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr,kout) !^ ad_zeta(i,Jstr ,kout)=ad_zeta(i,Jstr ,kout)+ & & ad_zeta(i,Jstr-1,kout) ad_zeta(i,Jstr-1,kout)=0.0_r8 END IF END DO END IF END IF ! !----------------------------------------------------------------------- ! Lateral boundary conditions at the eastern edge. !----------------------------------------------------------------------- ! IF (DOMAIN(ng)%Eastern_Edge(tile)) THEN ! ! Eastern edge, implicit upstream radiation condition. ! IF (ad_LBC(ieast,isFsur,ng)%radiation) THEN IF (iic(ng).ne.0) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%east(j)) THEN !^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend+1,j,kout)* & !^ & GRID(ng)%rmask(Iend+1,j) !^ ad_zeta(Iend+1,j,kout)=ad_zeta(Iend+1,j,kout)* & & GRID(ng)%rmask(Iend+1,j) IF (ad_LBC(ieast,isFsur,ng)%nudging) THEN !^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend+1,j,kout)- & !^ & tau*tl_zeta(Iend+1,j,know) !^ ad_zeta(Iend+1,j,know)=ad_zeta(Iend+1,j,know)- & & tau*ad_zeta(Iend+1,j,kout) END IF !^ tl_zeta(Iend+1,j,kout)=(cff*tl_zeta(Iend+1,j,know)+ & !^ & Cx *tl_zeta(Iend ,j,kout)- & !^ & MAX(Ce,0.0_r8)* & !^ & tl_grad(Iend+1,j )- & !^ & MIN(Ce,0.0_r8)* & !^ & tl_grad(Iend+1,j+1))/ & !^ & (cff+Cx) !^ adfac=ad_zeta(Iend+1,j,kout)/(cff+Cx) ad_grad(Iend+1,j )=ad_grad(Iend+1,j )- & & MAX(Ce,0.0_r8)*adfac ad_grad(Iend+1,j+1)=ad_grad(Iend+1,j+1)- & & MIN(Ce,0.0_r8)*adfac ad_zeta(Iend ,j,kout)=ad_zeta(Iend ,j,kout)+Cx *adfac ad_zeta(Iend+1,j,know)=ad_zeta(Iend+1,j,know)+cff*adfac ad_zeta(Iend+1,j,kout)=0.0_r8 END IF END DO END IF ! ! Eastern edge, explicit Chapman boundary condition. ! ELSE IF (ad_LBC(ieast,isFsur,ng)%Chapman_explicit) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%east(j)) THEN cff=dt2d*GRID(ng)%pm(Iend,j) cff1=SQRT(g*(GRID(ng)%h(Iend,j)+ & & zeta(Iend,j,know))) Cx=cff*cff1 !^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend+1,j,kout)* & !^ & GRID(ng)%rmask(Iend+1,j) !^ ad_zeta(Iend+1,j,kout)=ad_zeta(Iend+1,j,kout)* & & GRID(ng)%rmask(Iend+1,j) !^ tl_zeta(Iend+1,j,kout)=(1.0_r8-Cx)*tl_zeta(Iend+1,j,know)+& !^ & tl_Cx*(zeta(Iend+1,j,know)+ & !^ & zeta(Iend ,j,know))+ & !^ & Cx*tl_zeta(Iend,j,know) !^ ad_zeta(Iend+1,j,know)=ad_zeta(Iend+1,j,know)+ & & (1.0_r8-Cx)*ad_zeta(Iend+1,j,kout) ad_Cx=ad_Cx+(zeta(Iend+1,j,know)+ & & zeta(Iend ,j,know))*ad_zeta(Iend+1,j,kout) ad_zeta(Iend,j,know)=ad_zeta(Iend,j,know)+ & & Cx*ad_zeta(Iend+1,j,kout) ad_zeta(Iend+1,j,kout)=0.0_r8 !^ tl_Cx=cff*tl_cff1 !^ ad_cff1=ad_cff1+cff*ad_Cx ad_Cx=0.0_r8 !^ tl_cff1=0.5_r8*g*(GRID(ng)%tl_h(Iend,j)+ & !^ & tl_zeta(Iend,j,know))/cff1 !^ adfac=0.5_r8*g*ad_cff1/cff1 GRID(ng)%ad_h(Iend,j)=GRID(ng)%ad_h(Iend,j)+adfac ad_zeta(Iend,j,know)=ad_zeta(Iend,j,know)+adfac ad_cff1=0.0_r8 END IF END DO ! ! Eastern edge, implicit Chapman boundary condition. ! ELSE IF (ad_LBC(ieast,isFsur,ng)%Chapman_implicit) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%east(j)) THEN cff=dt2d*GRID(ng)%pm(Iend,j) cff1=SQRT(g*(GRID(ng)%h(Iend,j)+ & & zeta(Iend,j,know))) Cx=cff*cff1 cff2=1.0_r8/(1.0_r8+Cx) !^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend+1,j,kout)* & !^ & GRID(ng)%rmask(Iend+1,j) !^ ad_zeta(Iend+1,j,kout)=ad_zeta(Iend+1,j,kout)* & & GRID(ng)%rmask(Iend+1,j) !^ tl_zeta(Iend+1,j,kout)=tl_cff2*(zeta(Iend+1,j,know)+ & !^ & Cx*zeta(Iend,j,kout))+ & !^ & cff2*(tl_zeta(Iend+1,j,know)+ & !^ & tl_Cx*zeta(Iend,j,kout)+ & !^ & Cx*tl_zeta(Iend,j,kout)) !^ adfac=cff2*ad_zeta(Iend+1,j,kout) ad_zeta(Iend ,j,kout)=ad_zeta(Iend ,j,kout)+Cx*adfac ad_zeta(Iend+1,j,know)=ad_zeta(Iend+1,j,know)+adfac ad_Cx=ad_Cx+zeta(Iend,j,kout)*adfac ad_cff2=ad_cff2+ & & (zeta(Iend+1,j,know)+ & & Cx*zeta(Iend,j,kout))*ad_zeta(Iend+1,j,kout) ad_zeta(Iend+1,j,kout)=0.0_r8 !^ tl_cff2=-cff2*cff2*tl_Cx !^ ad_Cx=ad_Cx-cff2*cff2*ad_cff2 ad_cff2=0.0_r8 !^ tl_Cx=cff*tl_cff1 !^ ad_cff1=ad_cff1+cff*ad_Cx ad_Cx=0.0_r8 !^ tl_cff1=0.5_r8*g*(GRID(ng)%tl_h(Iend,j)+ & !^ & tl_zeta(Iend,j,know))/cff1 !^ adfac=0.5_r8*g*ad_cff1/cff1 GRID(ng)%ad_h(Iend,j)=GRID(ng)%ad_h(Iend,j)+adfac ad_zeta(Iend,j,know)=ad_zeta(Iend,j,know)+adfac ad_cff1=0.0_r8 END IF END DO ! ! Eastern edge, clamped boundary condition. ! ELSE IF (ad_LBC(ieast,isFsur,ng)%clamped) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%east(j)) THEN !^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend+1,j,kout)* & !^ & GRID(ng)%rmask(Iend+1,j) !^ ad_zeta(Iend+1,j,kout)=ad_zeta(Iend+1,j,kout)* & & GRID(ng)%rmask(Iend+1,j) !^ tl_zeta(Iend+1,j,kout)=0.0_r8 !^ ad_zeta(Iend+1,j,kout)=0.0_r8 END IF END DO ! ! Eastern edge, gradient boundary condition. ! ELSE IF (ad_LBC(ieast,isFsur,ng)%gradient) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%east(j)) THEN !^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend+1,j,kout)* & !^ & GRID(ng)%rmask(Iend+1,j) !^ ad_zeta(Iend+1,j,kout)=ad_zeta(Iend+1,j,kout)* & & GRID(ng)%rmask(Iend+1,j) !^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend,j,kout) !^ ad_zeta(Iend ,j,kout)=ad_zeta(Iend ,j,kout)+ & & ad_zeta(Iend+1,j,kout) ad_zeta(Iend+1,j,kout)=0.0_r8 END IF END DO ! ! Eastern edge, closed boundary condition. ! ELSE IF (ad_LBC(ieast,isFsur,ng)%closed) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%east(j)) THEN !^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend+1,j,kout)* & !^ & GRID(ng)%rmask(Iend+1,j) !^ ad_zeta(Iend+1,j,kout)=ad_zeta(Iend+1,j,kout)* & & GRID(ng)%rmask(Iend+1,j) !^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend,j,kout) !^ ad_zeta(Iend ,j,kout)=ad_zeta(Iend ,j,kout)+ & & ad_zeta(Iend+1,j,kout) ad_zeta(Iend+1,j,kout)=0.0_r8 END IF END DO END IF END IF ! !----------------------------------------------------------------------- ! Lateral boundary conditions at the western edge. !----------------------------------------------------------------------- ! IF (DOMAIN(ng)%Western_Edge(tile)) THEN ! ! Western edge, implicit upstream radiation condition. ! IF (ad_LBC(iwest,isFsur,ng)%radiation) THEN IF (iic(ng).ne.0) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%west(j)) THEN !^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr-1,j,kout)* & !^ & GRID(ng)%rmask(Istr-1,j) !^ ad_zeta(Istr-1,j,kout)=ad_zeta(Istr-1,j,kout)* & & GRID(ng)%rmask(Istr-1,j) IF (ad_LBC(iwest,isFsur,ng)%nudging) THEN !^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr-1,j,kout)- & !^ & tau*tl_zeta(Istr-1,j,know) !^ ad_zeta(Istr-1,j,know)=ad_zeta(Istr-1,j,know)- & & tau*ad_zeta(Istr-1,j,kout) END IF !^ tl_zeta(Istr-1,j,kout)=(cff*tl_zeta(Istr-1,j,know)+ & !^ & Cx *tl_zeta(Istr ,j,kout)- & !^ & MAX(Ce,0.0_r8)* & !^ & tl_grad(Istr-1,j )- & !^ & MIN(Ce,0.0_r8)* & !^ & tl_grad(Istr-1,j+1))/ & !^ & (cff+Cx) !^ adfac=ad_zeta(Istr-1,j,kout)/(cff+Cx) ad_grad(Istr-1,j )=ad_grad(Istr-1,j )- & & MAX(Ce,0.0_r8)*adfac ad_grad(Istr-1,j+1)=ad_grad(Istr-1,j+1)- & & MIN(Ce,0.0_r8)*adfac ad_zeta(Istr-1,j,know)=ad_zeta(Istr-1,j,know)+cff*adfac ad_zeta(Istr ,j,kout)=ad_zeta(Istr ,j,kout)+Cx *adfac ad_zeta(Istr-1,j,kout)=0.0_r8 END IF END DO END IF ! ! Western edge, explicit Chapman boundary condition. ! ELSE IF (ad_LBC(iwest,isFsur,ng)%Chapman_explicit) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%west(j)) THEN cff=dt2d*GRID(ng)%pm(Istr,j) cff1=SQRT(g*(GRID(ng)%h(Istr,j)+ & & zeta(Istr,j,know))) Cx=cff*cff1 !^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr-1,j,kout)* & !^ & GRID(ng)%rmask(Istr-1,j) !^ ad_zeta(Istr-1,j,kout)=ad_zeta(Istr-1,j,kout)* & & GRID(ng)%rmask(Istr-1,j) !^ tl_zeta(Istr-1,j,kout)=(1.0_r8-Cx)*tl_zeta(Istr-1,j,know)+& !^ & tl_Cx*(zeta(Istr-1,j,know)+ & !^ & zeta(Istr ,j,know))+ & !^ & Cx*tl_zeta(Istr,j,know) !^ ad_zeta(Istr-1,j,know)=ad_zeta(Istr-1,j,know)+ & & (1.0_r8-Cx)*ad_zeta(Istr-1,j,kout) ad_Cx=ad_Cx+(zeta(Istr-1,j,know)+ & & zeta(Istr ,j,know))*ad_zeta(Istr-1,j,kout) ad_zeta(Istr,j,know)=ad_zeta(Istr,j,know)+ & & Cx*ad_zeta(Istr-1,j,kout) ad_zeta(Istr-1,j,kout)=0.0_r8 !^ tl_Cx=cff*tl_cff1 !^ ad_cff1=ad_cff1+cff*ad_Cx ad_Cx=0.0_r8 !^ tl_cff1=0.5_r8*g*(GRID(ng)%tl_h(Istr,j)+ & !^ & tl_zeta(Istr,j,know))/cff1 !^ adfac=0.5_r8*g*ad_cff1/cff1 GRID(ng)%ad_h(Istr,j)=GRID(ng)%ad_h(Istr,j)+adfac ad_zeta(Istr,j,know)=ad_zeta(Istr,j,know)+adfac ad_cff1=0.0_r8 END IF END DO ! ! Western edge, implicit Chapman boundary condition. ! ELSE IF (ad_LBC(iwest,isFsur,ng)%Chapman_implicit) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%west(j)) THEN cff=dt2d*GRID(ng)%pm(Istr,j) cff1=SQRT(g*(GRID(ng)%h(Istr,j)+ & & zeta(Istr,j,know))) Cx=cff*cff1 cff2=1.0_r8/(1.0_r8+Cx) !^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr-1,j,kout)* & !^ & GRID(ng)%rmask(Istr-1,j) !^ ad_zeta(Istr-1,j,kout)=ad_zeta(Istr-1,j,kout)* & & GRID(ng)%rmask(Istr-1,j) !^ tl_zeta(Istr-1,j,kout)=tl_cff2*(zeta(Istr-1,j,know)+ & !^ & Cx*zeta(Istr,j,kout))+ & !^ & cff2*(tl_zeta(Istr-1,j,know)+ & !^ & tl_Cx*zeta(Istr,j,kout)+ & !^ & Cx*tl_zeta(Istr,j,kout)) !^ adfac=cff2*ad_zeta(Istr-1,j,kout) ad_zeta(Istr-1,j,know)=ad_zeta(Istr-1,j,know)+adfac ad_zeta(Istr ,j,kout)=ad_zeta(Istr ,j,kout)+Cx*adfac ad_Cx=ad_Cx+zeta(Istr,j,kout)*adfac ad_cff2=ad_cff2+ & & (zeta(Istr-1,j,know)+ & & Cx*zeta(Istr,j,kout))*ad_zeta(Istr-1,j,kout) ad_zeta(Istr-1,j,kout)=0.0_r8 !^ tl_cff2=-cff2*cff2*tl_Cx !^ ad_Cx=ad_Cx-cff2*cff2*ad_cff2 ad_cff2=0.0_r8 !^ tl_Cx=cff*tl_cff1 !^ ad_cff1=ad_cff1+cff*ad_Cx ad_Cx=0.0_r8 !^ tl_cff1=0.5_r8*g*(GRID(ng)%tl_h(Istr,j)+ & !^ & tl_zeta(Istr,j,know))/cff1 !^ adfac=0.5_r8*g*ad_cff1/cff1 GRID(ng)%ad_h(Istr,j)=GRID(ng)%ad_h(Istr,j)+adfac ad_zeta(Istr,j,know)=ad_zeta(Istr,j,know)+adfac ad_cff1=0.0_r8 END IF END DO ! ! Western edge, clamped boundary condition. ! ELSE IF (ad_LBC(iwest,isFsur,ng)%gradient) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%west(j)) THEN !^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr-1,j,kout)* & !^ & GRID(ng)%rmask(Istr-1,j) !^ ad_zeta(Istr-1,j,kout)=ad_zeta(Istr-1,j,kout)* & & GRID(ng)%rmask(Istr-1,j) !^ tl_zeta(Istr-1,j,kout)=0.0_r8 !^ ad_zeta(Istr-1,j,kout)=0.0_r8 END IF END DO ! ! Western edge, gradient boundary condition. ! ELSE IF (ad_LBC(iwest,isFsur,ng)%gradient) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%west(j)) THEN !^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr-1,j,kout)* & !^ & GRID(ng)%rmask(Istr-1,j) !^ ad_zeta(Istr-1,j,kout)=ad_zeta(Istr-1,j,kout)* & & GRID(ng)%rmask(Istr-1,j) !^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr,j,kout) !^ ad_zeta(Istr ,j,kout)=ad_zeta(Istr ,j,kout)+ & & ad_zeta(Istr-1,j,kout) ad_zeta(Istr-1,j,kout)=0.0_r8 END IF END DO ! ! Western edge, closed boundary condition. ! ELSE IF (ad_LBC(iwest,isFsur,ng)%closed) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%west(j)) THEN !^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr-1,j,kout)* & !^ & GRID(ng)%rmask(Istr-1,j) !^ ad_zeta(Istr-1,j,kout)=ad_zeta(Istr-1,j,kout)* & & GRID(ng)%rmask(Istr-1,j) !^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr,j,kout) !^ ad_zeta(Istr ,j,kout)=ad_zeta(Istr ,j,kout)+ & & ad_zeta(Istr-1,j,kout) ad_zeta(Istr-1,j,kout)=0.0_r8 END IF END DO END IF END IF RETURN END SUBROUTINE ad_zetabc_tile END MODULE ad_zetabc_mod