MODULE ad_bc_2d_mod ! !git $Id$ !svn $Id: ad_bc_2d.F 1151 2023-02-09 03:08:53Z arango $ !================================================== Hernan G. Arango === ! Copyright (c) 2002-2023 The ROMS/TOMS Group ! ! Licensed under a MIT/X style license ! ! See License_ROMS.md ! !======================================================================= ! ! ! These routines apply close, gradient or periodic boundary ! ! conditions to generic 2D adjoint fields. ! ! ! ! On Input: ! ! ! ! ng Nested grid number. ! ! tile Domain partition. ! ! LBi I-dimension Lower bound. ! ! UBi I-dimension Upper bound. ! ! LBj J-dimension Lower bound. ! ! UBj J-dimension Upper bound. ! ! ad_A 2D adjoint field. ! ! ! ! On Output: ! ! ! ! ad_A Processed 2D adjoint field. ! ! ! ! Routines: ! ! ! ! ad_bc_r2d_tile Boundary conditions for field at RHO-points ! ! ad_bc_u2d_tile Boundary conditions for field at U-points ! ! ad_bc_v2d_tile Boundary conditions for field at V-points ! ! ! !======================================================================= ! implicit none ! CONTAINS ! !*********************************************************************** SUBROUTINE ad_bc_r2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_A) !*********************************************************************** ! USE mod_param USE mod_boundary USE mod_ncparam USE mod_scalars ! USE ad_exchange_2d_mod, ONLY : ad_exchange_r2d_tile ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile integer, intent(in) :: LBi, UBi, LBj, UBj ! real(r8), intent(inout) :: ad_A(LBi:,LBj:) ! ! Local variable declarations. ! integer :: i, j real(r8) :: adfac ! !----------------------------------------------------------------------- ! 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 ! !----------------------------------------------------------------------- ! Apply adjoint periodic boundary conditons. !----------------------------------------------------------------------- ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN CALL ad_exchange_r2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_A) 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_A(Iend+1,Jend+1)=0.5_r8*(tl_A(Iend+1,Jend )+ & !^ & tl_A(Iend ,Jend+1)) !^ adfac=0.5_r8*ad_A(Iend+1,Jend+1) ad_A(Iend+1,Jend )=ad_A(Iend+1,Jend )+adfac ad_A(Iend ,Jend+1)=ad_A(Iend ,Jend+1)+adfac ad_A(Iend+1,Jend+1)=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_A(Istr-1,Jend+1)=0.5_r8*(tl_A(Istr-1,Jend )+ & !^ & tl_A(Istr ,Jend+1)) !^ adfac=0.5_r8*ad_A(Istr-1,Jend+1) ad_A(Istr-1,Jend )=ad_A(Istr-1,Jend )+adfac ad_A(Istr ,Jend+1)=ad_A(Istr ,Jend+1)+adfac ad_A(Istr-1,Jend+1)=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_A(Iend+1,Jstr-1)=0.5_r8*(tl_A(Iend ,Jstr-1)+ & !^ & tl_A(Iend+1,Jstr )) !^ adfac=0.5_r8*ad_A(Iend+1,Jstr-1) ad_A(Iend ,Jstr-1)=ad_A(Iend ,Jstr-1)+adfac ad_A(Iend+1,Jstr )=ad_A(Iend+1,Jstr )+adfac ad_A(Iend+1,Jstr-1)=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_A(Istr-1,Jstr-1)=0.5_r8*(tl_A(Istr ,Jstr-1)+ & !^ & tl_A(Istr-1,Jstr )) !^ adfac=0.5_r8*ad_A(Istr-1,Jstr-1) ad_A(Istr ,Jstr-1)=ad_A(Istr ,Jstr-1)+adfac ad_A(Istr-1,Jstr )=ad_A(Istr-1,Jstr )+adfac ad_A(Istr-1,Jstr-1)=0.0_r8 END IF END IF END IF ! !----------------------------------------------------------------------- ! Adjoint North-South gradient boundary conditions. !----------------------------------------------------------------------- ! IF (.not.NSperiodic(ng)) THEN IF (DOMAIN(ng)%Southern_Edge(tile)) THEN DO i=Istr,Iend IF (LBC_apply(ng)%south(i)) THEN !^ tl_A(i,Jstr-1)=tl_A(i,Jstr) !^ ad_A(i,Jstr )=ad_A(i,Jstr)+ad_A(i,Jstr-1) ad_A(i,Jstr-1)=0.0_r8 END IF END DO END IF IF (DOMAIN(ng)%Northern_Edge(tile)) THEN DO i=Istr,Iend IF (LBC_apply(ng)%north(i)) THEN !^ tl_A(i,Jend+1)=tl_A(i,Jend) !^ ad_A(i,Jend )=ad_A(i,Jend)+ad_A(i,Jend+1) ad_A(i,Jend+1)=0.0_r8 END IF END DO END IF END IF ! !----------------------------------------------------------------------- ! Adjoint East-West gradient boundary conditions. !----------------------------------------------------------------------- ! IF (.not.EWperiodic(ng)) THEN IF (DOMAIN(ng)%Western_Edge(tile)) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%west(j)) THEN !^ tl_A(Istr-1,j)=tl_A(Istr,j) !^ ad_A(Istr ,j)=ad_A(Istr,j)+ad_A(Istr-1,j) ad_A(Istr-1,j)=0.0_r8 END IF END DO END IF IF (DOMAIN(ng)%Eastern_Edge(tile)) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%east(j)) THEN !^ tl_A(Iend+1,j)=tl_A(Iend,j) !^ ad_A(Iend ,j)=ad_A(Iend,j)+ad_A(Iend+1,j) ad_A(Iend+1,j)=0.0_r8 END IF END DO END IF END IF RETURN END SUBROUTINE ad_bc_r2d_tile ! !*********************************************************************** SUBROUTINE ad_bc_u2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_A) !*********************************************************************** ! USE mod_param USE mod_boundary USE mod_grid USE mod_ncparam USE mod_scalars ! USE ad_exchange_2d_mod, ONLY : ad_exchange_u2d_tile ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile integer, intent(in) :: LBi, UBi, LBj, UBj ! real(r8), intent(inout) :: ad_A(LBi:,LBj:) ! ! Local variable declarations. ! integer :: Imin, Imax integer :: i, j real(r8) :: adfac ! !----------------------------------------------------------------------- ! 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 ! !----------------------------------------------------------------------- ! Apply adjoint periodic boundary conditons. !----------------------------------------------------------------------- ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN CALL ad_exchange_u2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_A) 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_A(Iend+1,Jend+1)=0.5_r8*(tl_A(Iend+1,Jend )+ & !^ & tl_A(Iend ,Jend+1)) !^ adfac=0.5_r8*ad_A(Iend+1,Jend+1) ad_A(Iend+1,Jend )=ad_A(Iend+1,Jend )+adfac ad_A(Iend ,Jend+1)=ad_A(Iend ,Jend+1)+adfac ad_A(Iend+1,Jend+1)=0.0_r8 END IF END IF IF (DOMAIN(ng)%NorthWest_Corner(tile)) THEN IF (LBC_apply(ng)%north(Istr ).and. & & LBC_apply(ng)%west (Jend+1)) THEN !^ tl_A(Istr ,Jend+1)=0.5_r8*(tl_A(Istr ,Jend )+ & !^ & tl_A(Istr+1,Jend+1)) !^ adfac=0.5_r8*ad_A(Istr,Jend+1) ad_A(Istr ,Jend )=ad_A(Istr ,Jend )+adfac ad_A(Istr+1,Jend+1)=ad_A(Istr+1,Jend+1)+adfac ad_A(Istr ,Jend+1)=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_A(Iend+1,Jstr-1)=0.5_r8*(tl_A(Iend ,Jstr-1)+ & !^ & tl_A(Iend+1,Jstr )) !^ adfac=0.5_r8*ad_A(Iend+1,Jstr-1) ad_A(Iend ,Jstr-1)=ad_A(Iend ,Jstr-1)+adfac ad_A(Iend+1,Jstr )=ad_A(Iend+1,Jstr )+adfac ad_A(Iend+1,Jstr-1)=0.0_r8 END IF END IF IF (DOMAIN(ng)%SouthWest_Corner(tile)) THEN IF (LBC_apply(ng)%south(Istr ).and. & & LBC_apply(ng)%west (Jstr-1)) THEN !^ tl_A(Istr ,Jstr-1)=0.5_r8*(tl_A(Istr+1,Jstr-1)+ & !^ & tl_A(Istr ,Jstr )) !^ adfac=0.5_r8*ad_A(Istr,Jstr-1) ad_A(Istr+1,Jstr-1)=ad_A(Istr+1,Jstr-1)+adfac ad_A(Istr ,Jstr )=ad_A(Istr ,Jstr )+adfac ad_A(Istr ,Jstr-1)=0.0_r8 END IF END IF END IF ! !----------------------------------------------------------------------- ! Adjoint North-South boundary conditions: Closed (free-slip/no-slip) ! or gradient. !----------------------------------------------------------------------- ! IF (.not.NSperiodic(ng)) THEN IF (DOMAIN(ng)%Southern_Edge(tile)) THEN IF (ad_LBC(isouth,isBu2d,ng)%closed) THEN IF (EWperiodic(ng)) THEN Imin=IstrU Imax=Iend ELSE Imin=Istr Imax=IendR END IF DO i=Imin,Imax IF (LBC_apply(ng)%south(i)) THEN !^ tl_A(i,Jstr-1)=tl_A(i,Jstr-1)*GRID(ng)%umask(i,Jstr-1) !^ ad_A(i,Jstr-1)=ad_A(i,Jstr-1)*GRID(ng)%umask(i,Jstr-1) !^ tl_A(i,Jstr-1)=gamma2(ng)*tl_A(i,Jstr) !^ ad_A(i,Jstr)=ad_A(i,Jstr)+gamma2(ng)*ad_A(i,Jstr-1) ad_A(i,Jstr-1)=0.0_r8 END IF END DO ELSE DO i=IstrU,Iend IF (LBC_apply(ng)%south(i)) THEN !^ tl_A(i,Jstr-1)=tl_A(i,Jstr) !^ ad_A(i,Jstr )=ad_A(i,Jstr)+ad_A(i,Jstr-1) ad_A(i,Jstr-1)=0.0_r8 END IF END DO END IF END IF IF (DOMAIN(ng)%Northern_Edge(tile)) THEN IF (ad_LBC(inorth,isBu2d,ng)%closed) THEN IF (EWperiodic(ng)) THEN Imin=IstrU Imax=Iend ELSE Imin=Istr Imax=IendR END IF DO i=Imin,Imax IF (LBC_apply(ng)%north(i)) THEN !^ tl_A(i,Jend+1)=tl_A(i,Jend+1)*GRID(ng)%umask(i,Jend+1) !^ ad_A(i,Jend+1)=ad_A(i,Jend+1)*GRID(ng)%umask(i,Jend+1) !^ tl_A(i,Jend+1)=gamma2(ng)*tl_A(i,Jend) !^ ad_A(i,Jend)=ad_A(i,Jend)+gamma2(ng)*ad_A(i,Jend+1) ad_A(i,Jend+1)=0.0_r8 END IF END DO ELSE DO i=IstrU,Iend IF (LBC_apply(ng)%north(i)) THEN !^ tl_A(i,Jend+1)=tl_A(i,Jend) !^ ad_A(i,Jend )=ad_A(i,Jend)+ad_A(i,Jend+1) ad_A(i,Jend+1)=0.0_r8 END IF END DO END IF END IF END IF ! !----------------------------------------------------------------------- ! Adjoint East-West boundary conditions: Closed or gradient. !----------------------------------------------------------------------- ! IF (.not.EWperiodic(ng)) THEN IF (DOMAIN(ng)%Western_Edge(tile)) THEN IF (ad_LBC(iwest,isBu2d,ng)%closed) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%west(j)) THEN !^ tl_A(Istr,j)=0.0_r8 !^ ad_A(Istr,j)=0.0_r8 END IF END DO ELSE DO j=Jstr,Jend IF (LBC_apply(ng)%west(j)) THEN !^ tl_A(Istr,j)=tl_A(Istr+1,j) !^ ad_A(Istr+1,j)=ad_A(Istr+1,j)+ad_A(Istr,j) ad_A(Istr ,j)=0.0_r8 END IF END DO END IF END IF IF (DOMAIN(ng)%Eastern_Edge(tile)) THEN IF (ad_LBC(ieast,isBu2d,ng)%closed) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%east(j)) THEN !^ A(Iend+1,j)=0.0_r8 !^ ad_A(Iend+1,j)=0.0_r8 END IF END DO ELSE DO j=Jstr,Jend IF (LBC_apply(ng)%east(j)) THEN !^ tl_A(Iend+1,j)=tl_A(Iend,j) !^ ad_A(Iend ,j)=ad_A(Iend,j)+ad_A(Iend+1,j) ad_A(Iend+1,j)=0.0_r8 END IF END DO END IF END IF END IF RETURN END SUBROUTINE ad_bc_u2d_tile ! !*********************************************************************** SUBROUTINE ad_bc_v2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_A) !*********************************************************************** ! USE mod_param USE mod_boundary USE mod_grid USE mod_ncparam USE mod_scalars ! USE ad_exchange_2d_mod, ONLY : ad_exchange_v2d_tile ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile integer, intent(in) :: LBi, UBi, LBj, UBj ! real(r8), intent(inout) :: ad_A(LBi:,LBj:) ! ! Local variable declarations. ! integer :: Jmin, Jmax integer :: i, j real(r8) :: adfac ! !----------------------------------------------------------------------- ! 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 ! !----------------------------------------------------------------------- ! Apply adjoint periodic boundary conditons. !----------------------------------------------------------------------- ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN CALL ad_exchange_v2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_A) 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_A(Iend+1,Jend+1)=0.5_r8*(tl_A(Iend+1,Jend )+ & !^ & tl_A(Iend ,Jend+1)) !^ adfac=0.5_r8*ad_A(Iend+1,Jend+1) ad_A(Iend+1,Jend )=ad_A(Iend+1,Jend )+adfac ad_A(Iend ,Jend+1)=ad_A(Iend ,Jend+1)+adfac ad_A(Iend+1,Jend+1)=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_A(Istr-1,Jend+1)=0.5_r8*(tl_A(Istr-1,Jend )+ & !^ & tl_A(Istr ,Jend+1)) !^ adfac=0.5_r8*ad_A(Istr-1,Jend+1) ad_A(Istr-1,Jend )=ad_A(Istr-1,Jend )+adfac ad_A(Istr ,Jend+1)=ad_A(Istr ,Jend+1)+adfac ad_A(Istr-1,Jend+1)=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 )) THEN !^ tl_A(Iend+1,Jstr )=0.5_r8*(tl_A(Iend ,Jstr )+ & !^ & tl_A(Iend+1,Jstr+1)) !^ adfac=0.5_r8*ad_A(Iend+1,Jstr ) ad_A(Iend ,Jstr )=ad_A(Iend ,Jstr )+adfac ad_A(Iend+1,Jstr+1)=ad_A(Iend+1,Jstr+1)+adfac ad_A(Iend+1,Jstr )=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 )) THEN !^ tl_A(Istr-1,Jstr )=0.5_r8*(tl_A(Istr ,Jstr )+ & !^ & tl_A(Istr-1,Jstr+1)) !^ adfac=0.5_r8*ad_A(Istr-1,Jstr ) ad_A(Istr ,Jstr )=ad_A(Istr ,Jstr )+adfac ad_A(Istr-1,Jstr+1)=ad_A(Istr-1,Jstr+1)+adfac ad_A(Istr-1,Jstr )=0.0_r8 END IF END IF END IF ! !----------------------------------------------------------------------- ! Adjoint North-South boundary conditions: Closed or Gradient. !----------------------------------------------------------------------- ! IF (.not.NSperiodic(ng)) THEN IF (DOMAIN(ng)%Southern_Edge(tile)) THEN IF (ad_LBC(isouth,isBv2d,ng)%closed) THEN DO i=Istr,Iend IF (LBC_apply(ng)%south(i)) THEN !^ tl_A(i,Jstr)=0.0_r8 !^ ad_A(i,Jstr)=0.0_r8 END IF END DO ELSE DO i=Istr,Iend IF (LBC_apply(ng)%south(i)) THEN !^ tl_A(i,Jstr)=tl_A(i,Jstr+1) !^ ad_A(i,Jstr+1)=ad_A(i,Jstr+1)+ad_A(i,Jstr) ad_A(i,Jstr)=0.0_r8 END IF END DO END IF END IF IF (DOMAIN(ng)%Northern_Edge(tile)) THEN IF (ad_LBC(inorth,isBv2d,ng)%closed) THEN DO i=Istr,Iend IF (LBC_apply(ng)%north(i)) THEN !^ tl_A(i,Jend+1)=0.0_r8 !^ ad_A(i,Jend+1)=0.0_r8 END IF END DO ELSE DO i=Istr,Iend IF (LBC_apply(ng)%north(i)) THEN !^ tl_A(i,Jend+1)=tl_A(i,Jend) !^ ad_A(i,Jend)=ad_A(i,Jend)+ad_A(i,Jend+1) ad_A(i,Jend+1)=0.0_r8 END IF END DO END IF END IF END IF ! !----------------------------------------------------------------------- ! Adjoint East-West boundary conditions. !----------------------------------------------------------------------- ! IF (.not.EWperiodic(ng)) THEN IF (DOMAIN(ng)%Western_Edge(tile)) THEN IF (ad_LBC(iwest,isBv2d,ng)%closed) THEN IF (NSperiodic(ng)) THEN Jmin=JstrV Jmax=Jend ELSE Jmin=Jstr Jmax=JendR END IF DO j=Jmin,Jmax IF (LBC_apply(ng)%west(j)) THEN !^ tl_A(Istr-1,j)=tl_A(Istr-1,j)*GRID(ng)%vmask(Istr-1,j) !^ ad_A(Istr-1,j)=ad_A(Istr-1,j)*GRID(ng)%vmask(Istr-1,j) !^ tl_A(Istr-1,j)=gamma2(ng)*tl_A(Istr,j) !^ ad_A(Istr ,j)=ad_A(Istr,j)+gamma2(ng)*ad_A(Istr-1,j) ad_A(Istr-1,j)=0.0_r8 END IF END DO ELSE DO j=JstrV,Jend IF (LBC_apply(ng)%west(j)) THEN !^ tl_A(Istr-1,j)=tl_A(Istr,j) !^ ad_A(Istr ,j)=ad_A(Istr,j)+ad_A(Istr-1,j) ad_A(Istr-1,j)=0.0_r8 END IF END DO END IF END IF IF (DOMAIN(ng)%Eastern_Edge(tile)) THEN IF (ad_LBC(ieast,isBv2d,ng)%closed) THEN IF (NSperiodic(ng)) THEN Jmin=JstrV Jmax=Jend ELSE Jmin=Jstr Jmax=JendR END IF DO j=Jmin,Jmax IF (LBC_apply(ng)%east(j)) THEN !^ tl_A(Iend+1,j)=tl_A(Iend+1,j)*GRID(ng)%vmask(Iend+1,j) !^ ad_A(Iend+1,j)=ad_A(Iend+1,j)*GRID(ng)%vmask(Iend+1,j) !^ tl_A(Iend+1,j)=gamma2(ng)*tl_A(Iend,j) !^ ad_A(Iend ,j)=ad_A(Iend,j)+ad_A(Iend+1,j) ad_A(Iend+1,j)=0.0_r8 END IF END DO ELSE DO j=JstrV,Jend IF (LBC_apply(ng)%east(j)) THEN !^ tl_A(Iend+1,j)=tl_A(Iend,j) !^ ad_A(Iend ,j)=ad_A(Iend,j)+ad_A(Iend+1,j) ad_A(Iend+1,j)=0.0_r8 END IF END DO END IF END IF END IF RETURN END SUBROUTINE ad_bc_v2d_tile ! !*********************************************************************** SUBROUTINE ad_dabc_r2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_A) !*********************************************************************** ! USE mod_param USE mod_boundary USE mod_ncparam USE mod_scalars ! USE ad_exchange_2d_mod, ONLY : ad_exchange_r2d_tile ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile integer, intent(in) :: LBi, UBi, LBj, UBj ! real(r8), intent(inout) :: ad_A(LBi:,LBj:) ! ! Local variable declarations. ! integer :: i, j real(r8) :: adfac ! !----------------------------------------------------------------------- ! 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 ! !----------------------------------------------------------------------- ! Apply adjoint periodic boundary conditons. !----------------------------------------------------------------------- ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN CALL ad_exchange_r2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_A) 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_A(Iend+1,Jend+1)=0.5_r8*(tl_A(Iend+1,Jend )+ & !^ & tl_A(Iend ,Jend+1)) !^ adfac=0.5_r8*ad_A(Iend+1,Jend+1) ad_A(Iend+1,Jend )=ad_A(Iend+1,Jend )+adfac ad_A(Iend ,Jend+1)=ad_A(Iend ,Jend+1)+adfac ad_A(Iend+1,Jend+1)=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_A(Istr-1,Jend+1)=0.5_r8*(tl_A(Istr-1,Jend )+ & !^ & tl_A(Istr ,Jend+1)) !^ adfac=0.5_r8*ad_A(Istr-1,Jend+1) ad_A(Istr-1,Jend )=ad_A(Istr-1,Jend )+adfac ad_A(Istr ,Jend+1)=ad_A(Istr ,Jend+1)+adfac ad_A(Istr-1,Jend+1)=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_A(Iend+1,Jstr-1)=0.5_r8*(tl_A(Iend ,Jstr-1)+ & !^ & tl_A(Iend+1,Jstr )) !^ adfac=0.5_r8*ad_A(Iend+1,Jstr-1) ad_A(Iend ,Jstr-1)=ad_A(Iend ,Jstr-1)+adfac ad_A(Iend+1,Jstr )=ad_A(Iend+1,Jstr )+adfac ad_A(Iend+1,Jstr-1)=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_A(Istr-1,Jstr-1)=0.5_r8*(tl_A(Istr ,Jstr-1)+ & !^ & tl_A(Istr-1,Jstr )) !^ adfac=0.5_r8*ad_A(Istr-1,Jstr-1) ad_A(Istr ,Jstr-1)=ad_A(Istr ,Jstr-1)+adfac ad_A(Istr-1,Jstr )=ad_A(Istr-1,Jstr )+adfac ad_A(Istr-1,Jstr-1)=0.0_r8 END IF END IF END IF ! !----------------------------------------------------------------------- ! Adjoint North-South gradient boundary conditions. !----------------------------------------------------------------------- ! IF (.not.NSperiodic(ng)) THEN IF (DOMAIN(ng)%Southern_Edge(tile)) THEN DO i=Istr,Iend IF (LBC_apply(ng)%south(i)) THEN !^ tl_A(i,Jstr-1)=tl_A(i,Jstr) !^ ad_A(i,Jstr )=ad_A(i,Jstr)+ad_A(i,Jstr-1) ad_A(i,Jstr-1)=0.0_r8 END IF END DO END IF IF (DOMAIN(ng)%Northern_Edge(tile)) THEN DO i=Istr,Iend IF (LBC_apply(ng)%north(i)) THEN !^ tl_A(i,Jend+1)=tl_A(i,Jend) !^ ad_A(i,Jend )=ad_A(i,Jend)+ad_A(i,Jend+1) ad_A(i,Jend+1)=0.0_r8 END IF END DO END IF END IF ! !----------------------------------------------------------------------- ! Adjoint East-West gradient boundary conditions. !----------------------------------------------------------------------- ! IF (.not.EWperiodic(ng)) THEN IF (DOMAIN(ng)%Western_Edge(tile)) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%west(j)) THEN !^ tl_A(Istr-1,j)=tl_A(Istr,j) !^ ad_A(Istr ,j)=ad_A(Istr,j)+ad_A(Istr-1,j) ad_A(Istr-1,j)=0.0_r8 END IF END DO END IF IF (DOMAIN(ng)%Eastern_Edge(tile)) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%east(j)) THEN !^ tl_A(Iend+1,j)=tl_A(Iend,j) !^ ad_A(Iend ,j)=ad_A(Iend,j)+ad_A(Iend+1,j) ad_A(Iend+1,j)=0.0_r8 END IF END DO END IF END IF RETURN END SUBROUTINE ad_dabc_r2d_tile ! !*********************************************************************** SUBROUTINE ad_dabc_u2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_A) !*********************************************************************** ! USE mod_param USE mod_boundary USE mod_grid USE mod_ncparam USE mod_scalars ! USE ad_exchange_2d_mod, ONLY : ad_exchange_u2d_tile ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile integer, intent(in) :: LBi, UBi, LBj, UBj ! real(r8), intent(inout) :: ad_A(LBi:,LBj:) ! ! Local variable declarations. ! integer :: Imin, Imax integer :: i, j real(r8) :: adfac ! !----------------------------------------------------------------------- ! 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 ! !----------------------------------------------------------------------- ! Apply adjoint periodic boundary conditons. !----------------------------------------------------------------------- ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN CALL ad_exchange_u2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_A) 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_A(Iend+1,Jend+1)=0.5_r8*(tl_A(Iend+1,Jend )+ & !^ & tl_A(Iend ,Jend+1)) !^ adfac=0.5_r8*ad_A(Iend+1,Jend+1) ad_A(Iend+1,Jend )=ad_A(Iend+1,Jend )+adfac ad_A(Iend ,Jend+1)=ad_A(Iend ,Jend+1)+adfac ad_A(Iend+1,Jend+1)=0.0_r8 END IF END IF IF (DOMAIN(ng)%NorthWest_Corner(tile)) THEN IF (LBC_apply(ng)%north(Istr ).and. & & LBC_apply(ng)%west (Jend+1)) THEN !^ tl_A(Istr ,Jend+1)=0.5_r8*(tl_A(Istr ,Jend )+ & !^ & tl_A(Istr+1,Jend+1)) !^ adfac=0.5_r8*ad_A(Istr,Jend+1) ad_A(Istr ,Jend )=ad_A(Istr ,Jend )+adfac ad_A(Istr+1,Jend+1)=ad_A(Istr+1,Jend+1)+adfac ad_A(Istr ,Jend+1)=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_A(Iend+1,Jstr-1)=0.5_r8*(tl_A(Iend ,Jstr-1)+ & !^ & tl_A(Iend+1,Jstr )) !^ adfac=0.5_r8*ad_A(Iend+1,Jstr-1) ad_A(Iend ,Jstr-1)=ad_A(Iend ,Jstr-1)+adfac ad_A(Iend+1,Jstr )=ad_A(Iend+1,Jstr )+adfac ad_A(Iend+1,Jstr-1)=0.0_r8 END IF END IF IF (DOMAIN(ng)%SouthWest_Corner(tile)) THEN IF (LBC_apply(ng)%south(Istr ).and. & & LBC_apply(ng)%west (Jstr-1)) THEN !^ tl_A(Istr ,Jstr-1)=0.5_r8*(tl_A(Istr+1,Jstr-1)+ & !^ & tl_A(Istr ,Jstr )) !^ adfac=0.5_r8*ad_A(Istr,Jstr-1) ad_A(Istr+1,Jstr-1)=ad_A(Istr+1,Jstr-1)+adfac ad_A(Istr ,Jstr )=ad_A(Istr ,Jstr )+adfac ad_A(Istr ,Jstr-1)=0.0_r8 END IF END IF END IF ! !----------------------------------------------------------------------- ! Adjoint North-South gradient boundary conditions. !----------------------------------------------------------------------- ! IF (.not.NSperiodic(ng)) THEN IF (DOMAIN(ng)%Southern_Edge(tile)) THEN DO i=IstrU,Iend IF (LBC_apply(ng)%south(i)) THEN !^ tl_A(i,Jstr-1)=tl_A(i,Jstr) !^ ad_A(i,Jstr )=ad_A(i,Jstr)+ad_A(i,Jstr-1) ad_A(i,Jstr-1)=0.0_r8 END IF END DO END IF IF (DOMAIN(ng)%Northern_Edge(tile)) THEN DO i=IstrU,Iend IF (LBC_apply(ng)%north(i)) THEN !^ tl_A(i,Jend+1)=tl_A(i,Jend) !^ ad_A(i,Jend )=ad_A(i,Jend)+ad_A(i,Jend+1) ad_A(i,Jend+1)=0.0_r8 END IF END DO END IF END IF ! !----------------------------------------------------------------------- ! Adjoint East-West gradient boundary conditions. !----------------------------------------------------------------------- ! IF (.not.EWperiodic(ng)) THEN IF (DOMAIN(ng)%Western_Edge(tile)) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%west(j)) THEN !^ tl_A(Istr,j)=tl_A(Istr+1,j) !^ ad_A(Istr+1,j)=ad_A(Istr+1,j)+ad_A(Istr,j) ad_A(Istr ,j)=0.0_r8 END IF END DO END IF IF (DOMAIN(ng)%Eastern_Edge(tile)) THEN DO j=Jstr,Jend IF (LBC_apply(ng)%east(j)) THEN !^ tl_A(Iend+1,j)=tl_A(Iend,j) !^ ad_A(Iend ,j)=ad_A(Iend,j)+ad_A(Iend+1,j) ad_A(Iend+1,j)=0.0_r8 END IF END DO END IF END IF RETURN END SUBROUTINE ad_dabc_u2d_tile ! !*********************************************************************** SUBROUTINE ad_dabc_v2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_A) !*********************************************************************** ! USE mod_param USE mod_boundary USE mod_grid USE mod_ncparam USE mod_scalars ! USE ad_exchange_2d_mod, ONLY : ad_exchange_v2d_tile ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile integer, intent(in) :: LBi, UBi, LBj, UBj ! real(r8), intent(inout) :: ad_A(LBi:,LBj:) ! ! Local variable declarations. ! integer :: Jmin, Jmax integer :: i, j real(r8) :: adfac ! !----------------------------------------------------------------------- ! 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 ! !----------------------------------------------------------------------- ! Apply adjoint periodic boundary conditons. !----------------------------------------------------------------------- ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN CALL ad_exchange_v2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_A) 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_A(Iend+1,Jend+1)=0.5_r8*(tl_A(Iend+1,Jend )+ & !^ & tl_A(Iend ,Jend+1)) !^ adfac=0.5_r8*ad_A(Iend+1,Jend+1) ad_A(Iend+1,Jend )=ad_A(Iend+1,Jend )+adfac ad_A(Iend ,Jend+1)=ad_A(Iend ,Jend+1)+adfac ad_A(Iend+1,Jend+1)=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_A(Istr-1,Jend+1)=0.5_r8*(tl_A(Istr-1,Jend )+ & !^ & tl_A(Istr ,Jend+1)) !^ adfac=0.5_r8*ad_A(Istr-1,Jend+1) ad_A(Istr-1,Jend )=ad_A(Istr-1,Jend )+adfac ad_A(Istr ,Jend+1)=ad_A(Istr ,Jend+1)+adfac ad_A(Istr-1,Jend+1)=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 )) THEN !^ tl_A(Iend+1,Jstr )=0.5_r8*(tl_A(Iend ,Jstr )+ & !^ & tl_A(Iend+1,Jstr+1)) !^ adfac=0.5_r8*ad_A(Iend+1,Jstr ) ad_A(Iend ,Jstr )=ad_A(Iend ,Jstr )+adfac ad_A(Iend+1,Jstr+1)=ad_A(Iend+1,Jstr+1)+adfac ad_A(Iend+1,Jstr )=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 )) THEN !^ tl_A(Istr-1,Jstr )=0.5_r8*(tl_A(Istr ,Jstr )+ & !^ & tl_A(Istr-1,Jstr+1)) !^ adfac=0.5_r8*ad_A(Istr-1,Jstr ) ad_A(Istr ,Jstr )=ad_A(Istr ,Jstr )+adfac ad_A(Istr-1,Jstr+1)=ad_A(Istr-1,Jstr+1)+adfac ad_A(Istr-1,Jstr )=0.0_r8 END IF END IF END IF ! !----------------------------------------------------------------------- ! Adjoint North-South gradient boundary conditions. !----------------------------------------------------------------------- ! IF (.not.NSperiodic(ng)) THEN IF (DOMAIN(ng)%Southern_Edge(tile)) THEN DO i=Istr,Iend IF (LBC_apply(ng)%south(i)) THEN !^ tl_A(i,Jstr)=tl_A(i,Jstr+1) !^ ad_A(i,Jstr+1)=ad_A(i,Jstr+1)+ad_A(i,Jstr) ad_A(i,Jstr)=0.0_r8 END IF END DO END IF IF (DOMAIN(ng)%Northern_Edge(tile)) THEN DO i=Istr,Iend IF (LBC_apply(ng)%north(i)) THEN !^ tl_A(i,Jend+1)=tl_A(i,Jend) !^ ad_A(i,Jend)=ad_A(i,Jend)+ad_A(i,Jend+1) ad_A(i,Jend+1)=0.0_r8 END IF END DO END IF END IF ! !----------------------------------------------------------------------- ! Adjoint East-West gradient boundary conditions. !----------------------------------------------------------------------- ! IF (.not.EWperiodic(ng)) THEN IF (DOMAIN(ng)%Western_Edge(tile)) THEN DO j=JstrV,Jend IF (LBC_apply(ng)%west(j)) THEN !^ tl_A(Istr-1,j)=tl_A(Istr,j) !^ ad_A(Istr ,j)=ad_A(Istr,j)+ad_A(Istr-1,j) ad_A(Istr-1,j)=0.0_r8 END IF END DO END IF IF (DOMAIN(ng)%Eastern_Edge(tile)) THEN DO j=JstrV,Jend IF (LBC_apply(ng)%east(j)) THEN !^ tl_A(Iend+1,j)=tl_A(Iend,j) !^ ad_A(Iend ,j)=ad_A(Iend,j)+ad_A(Iend+1,j) ad_A(Iend+1,j)=0.0_r8 END IF END DO END IF END IF RETURN END SUBROUTINE ad_dabc_v2d_tile END MODULE ad_bc_2d_mod