MODULE ad_ini_fields_mod ! !git $Id$ !svn $Id: ad_ini_fields.F 1180 2023-07-13 02:42:10Z 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 routine initializes other time levels for 2D adjoint fields. ! ! It also couples 3D and 2D momentum equations: it initializes 2D ! ! momentum (ubar,vbar) to the vertical integral of initial 3D ! ! momentum (u,v). ! ! ! !======================================================================= ! implicit none ! PRIVATE PUBLIC :: ad_ini_fields PUBLIC :: ad_ini_zeta PUBLIC :: ad_out_fields PUBLIC :: ad_out_zeta ! CONTAINS ! !*********************************************************************** SUBROUTINE ad_ini_fields (ng, tile, model) !*********************************************************************** ! USE mod_param USE mod_grid USE mod_coupling USE mod_ocean USE mod_stepping ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model ! ! Local variable declarations. ! character (len=*), parameter :: MyFile = & & "ROMS/Adjoint/ad_ini_fields.F" ! 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 wclock_on (ng, iADM, 2, 57, MyFile) CALL ad_ini_fields_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & kstp(ng), krhs(ng), knew(ng), & & nstp(ng), nnew(ng), & & GRID(ng) % rmask, & & GRID(ng) % umask, & & GRID(ng) % vmask, & & GRID(ng) % Hz, & & GRID(ng) % ad_Hz, & & OCEAN(ng) % ad_t, & & OCEAN(ng) % u, & & OCEAN(ng) % ad_u, & & OCEAN(ng) % v, & & OCEAN(ng) % ad_v, & & OCEAN(ng) % ubar, & & OCEAN(ng) % ad_ubar_sol, & & OCEAN(ng) % ad_ubar, & & OCEAN(ng) % vbar, & & OCEAN(ng) % ad_vbar_sol, & & OCEAN(ng) % ad_vbar, & & OCEAN(ng) % zeta, & & OCEAN(ng) % ad_zeta) CALL wclock_off (ng, iADM, 2, 89, MyFile) ! RETURN END SUBROUTINE ad_ini_fields ! !*********************************************************************** SUBROUTINE ad_ini_fields_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & kstp, krhs, knew, & & nstp, nnew, & & rmask, umask, vmask, & & Hz, ad_Hz, & & ad_t, & & u, ad_u, & & v, ad_v, & & ubar, ad_ubar_sol, ad_ubar, & & vbar, ad_vbar_sol, ad_vbar, & & zeta, ad_zeta) !*********************************************************************** ! USE mod_param USE mod_ncparam USE mod_scalars ! USE ad_exchange_2d_mod USE ad_exchange_3d_mod USE mp_exchange_mod, ONLY : ad_mp_exchange2d USE mp_exchange_mod, ONLY : ad_mp_exchange3d, ad_mp_exchange4d USE ad_set_depth_mod, ONLY : ad_set_depth_tile USE ad_t3dbc_mod, ONLY : ad_t3dbc_tile USE ad_u3dbc_mod, ONLY : ad_u3dbc_tile USE ad_v3dbc_mod, ONLY : ad_v3dbc_tile USE ad_u2dbc_mod, ONLY : ad_u2dbc_tile USE ad_v2dbc_mod, ONLY : ad_v2dbc_tile ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model integer, intent(in) :: LBi, UBi, LBj, UBj integer, intent(in) :: IminS, ImaxS, JminS, JmaxS integer, intent(in) :: kstp, krhs, knew integer, intent(in) :: nstp, nnew ! real(r8), intent(in) :: rmask(LBi:,LBj:) real(r8), intent(in) :: umask(LBi:,LBj:) real(r8), intent(in) :: vmask(LBi:,LBj:) real(r8), intent(in) :: Hz(LBi:,LBj:,:) real(r8), intent(in) :: u(LBi:,LBj:,:,:) real(r8), intent(in) :: v(LBi:,LBj:,:,:) real(r8), intent(in) :: ubar(LBi:,LBj:,:) real(r8), intent(in) :: vbar(LBi:,LBj:,:) real(r8), intent(in) :: zeta(LBi:,LBj:,:) real(r8), intent(inout) :: ad_Hz(LBi:,LBj:,:) real(r8), intent(inout) :: ad_t(LBi:,LBj:,:,:,:) real(r8), intent(inout) :: ad_u(LBi:,LBj:,:,:) real(r8), intent(inout) :: ad_v(LBi:,LBj:,:,:) real(r8), intent(inout) :: ad_ubar_sol(LBi:,LBj:) real(r8), intent(inout) :: ad_vbar_sol(LBi:,LBj:) real(r8), intent(inout) :: ad_ubar(LBi:,LBj:,:) real(r8), intent(inout) :: ad_vbar(LBi:,LBj:,:) real(r8), intent(inout) :: ad_zeta(LBi:,LBj:,:) ! ! Local variable declarations. ! integer :: i, ic, itrc, j, k, kbed real(r8) :: cff1 real(r8) :: adfac, ad_cff1, ad_cff2 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: CF real(r8), dimension(IminS:ImaxS,0:N(ng)) :: DC real(r8), dimension(IminS:ImaxS,0:N(ng)) :: ad_CF real(r8), dimension(IminS:ImaxS,0:N(ng)) :: ad_DC ! !----------------------------------------------------------------------- ! 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_cff1=0.0_r8 ad_cff2=0.0_r8 DO k=0,N(ng) DO i=IminS,ImaxS ad_CF(i,k)=0.0_r8 ad_DC(i,k)=0.0_r8 END DO END DO ! !----------------------------------------------------------------------- ! Adjoint of initialize other time levels for tracers. Only time level ! "nstp" is needed in the adjoint. !----------------------------------------------------------------------- ! ! Apply boundary conditions. ! !^ CALL mp_exchange4d (ng, tile, model, 1, & !^ & LBi, UBi, LBj, UBj, 1, N(ng), 1, NT(ng), & !^ & NghostPoints, & !^ & EWperiodic(ng), NSperiodic(ng), & !^ & tl_t(:,:,:,nstp,:), & !^ & tl_t(:,:,:,nnew,:)) !^ only nstp needed CALL ad_mp_exchange4d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, 1, N(ng), 1, NT(ng), & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_t(:,:,:,nstp,:)) ! ic=0 DO itrc=1,NT(ng) IF (LtracerCLM(itrc,ng).and.LnudgeTCLM(itrc,ng)) THEN ic=ic+1 ! OBC nudging coefficient index END IF IF (EWperiodic(ng).or.NSperiodic(ng)) THEN !^ CALL exchange_r3d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, 1, N(ng), & !^ & tl_t(:,:,:,nstp,itrc)) !^ only nstp needed CALL ad_exchange_r3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ad_t(:,:,:,nstp,itrc)) END IF ! !^ CALL tl_t3dbc_tile (ng, tile, itrc, ic, & !^ & LBi, UBi, LBj, UBj, N(ng), NT(ng), & !^ & IminS, ImaxS, JminS, JmaxS, & !^ & nstp, nstp, & !^ & tl_t) !^ only nstp needed CALL ad_t3dbc_tile (ng, tile, itrc, ic, & & LBi, UBi, LBj, UBj, N(ng), NT(ng), & & IminS, ImaxS, JminS, JmaxS, & & nstp, nstp, & & ad_t) ! ! Adjoint of tracers initialization. ! DO k=1,N(ng) DO j=JstrB,JendB DO i=IstrB,IendB !^ tl_t(i,j,k,nstp,itrc)=tl_cff1 !^ ad_cff1=ad_cff1+ad_t(i,j,k,nstp,itrc) ad_t(i,j,k,nstp,itrc)=0.0_r8 !^ tl_cff1=tl_cff1*rmask(i,j) !^ ad_cff1=ad_cff1*rmask(i,j) !^ tl_cff1=tl_t(i,j,k,nstp,itrc) !^ ad_t(i,j,k,nstp,itrc)=ad_t(i,j,k,nstp,itrc)+ad_cff1 ad_cff1=0.0_r8 END DO END DO END DO END DO ! !----------------------------------------------------------------------- ! Adjoint of compute vertically-integrated momentum (tl_ubar, tl_vbar) ! from initial 3D momentum (tl_u, tl_v). !----------------------------------------------------------------------- ! ! Apply boundary conditions. ! !^ CALL mp_exchange2d (ng, tile, model, 4, & !^ & LBi, UBi, LBj, UBj, & !^ & NghostPoints, & !^ & EWperiodic(ng), NSperiodic(ng), & !^ & tl_ubar(:,:,kstp), tl_vbar(:,:,kstp), & !^ & tl_ubar(:,:,knew), tl_vbar(:,:,knew)) !^ CALL ad_mp_exchange2d (ng, tile, model, 4, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_ubar(:,:,kstp), ad_vbar(:,:,kstp), & & ad_ubar(:,:,knew), ad_vbar(:,:,knew)) ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN !^ CALL exchange_v2d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & tl_vbar(:,:,knew)) !^ CALL ad_exchange_v2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_vbar(:,:,knew)) !^ CALL exchange_u2d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & tl_ubar(:,:,knew)) !^ CALL ad_exchange_u2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_ubar(:,:,knew)) !^ CALL exchange_v2d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & tl_vbar(:,:,kstp)) !^ CALL ad_exchange_v2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_vbar(:,:,kstp)) !^ CALL exchange_u2d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & tl_ubar(:,:,kstp)) !^ CALL ad_exchange_u2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_ubar(:,:,kstp)) END IF ! IF (.not.(ANY(ad_LBC(:,isUbar,ng)%radiation).or. & & ANY(ad_LBC(:,isVbar,ng)%radiation).or. & & ANY(ad_LBC(:,isUbar,ng)%Flather).or. & & ANY(ad_LBC(:,isVbar,ng)%Flather))) THEN !^ CALL tl_v2dbc_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & IminS, ImaxS, JminS, JmaxS, & !^ & krhs, kstp, knew, & !^ & ubar, vbar, zeta, & !^ & tl_ubar, tl_vbar, tl_zeta) !^ CALL ad_v2dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, knew, & & ubar, vbar, zeta, & & ad_ubar, ad_vbar, ad_zeta) !^ CALL tl_u2dbc_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & IminS, ImaxS, JminS, JmaxS, & !^ & krhs, kstp, knew, & !^ & ubar, vbar, zeta, & !^ & tl_ubar, tl_vbar, tl_zeta) !^ CALL ad_u2dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, knew, & & ubar, vbar, zeta, & & ad_ubar, ad_vbar, ad_zeta) !^ CALL tl_v2dbc_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & IminS, ImaxS, JminS, JmaxS, & !^ & krhs, kstp, kstp, & !^ & ubar, vbar, zeta, & !^ & tl_ubar, tl_vbar, tl_zeta) !^ CALL ad_v2dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, kstp, & & ubar, vbar, zeta, & & ad_ubar, ad_vbar, ad_zeta) !^ CALL tl_u2dbc_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & IminS, ImaxS, JminS, JmaxS, & !^ & krhs, kstp, kstp, & !^ & ubar, vbar, zeta, & !^ & tl_ubar, tl_vbar, tl_zeta) !^ CALL ad_u2dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, kstp, & & ubar, vbar, zeta, & & ad_ubar, ad_vbar, ad_zeta) END IF ! ! Load 2D adjoint momentum solution into IO arrays. ! DO j=JstrT,JendT DO i=IstrP,IendT ad_ubar_sol(i,j)=ad_ubar(i,j,kstp)+ad_ubar(i,j,knew) END DO IF (j.ge.JstrP) THEN DO i=IstrT,IendT ad_vbar_sol(i,j)=ad_vbar(i,j,kstp)+ad_vbar(i,j,knew) END DO END IF END DO ! ! Adjoint of compute vertically-integrated momentum (tl_ubar, tl_vbar) ! from initial 3D momentum (tl_u, tl_v). Here DC(i,1:N) are the grid ! cell thicknesses, DC(i,0) is total depth of the water column, and ! CF(i,0) is the vertical integral. ! DO j=JstrB,JendB IF (j.ge.JstrM) THEN DO i=IstrB,IendB DC(i,0)=0.0_r8 CF(i,0)=0.0_r8 END DO DO k=1,N(ng) DO i=IstrB,IendB DC(i,k)=0.5_r8*(Hz(i,j,k)+Hz(i,j-1,k)) DC(i,0)=DC(i,0)+DC(i,k) CF(i,0)=CF(i,0)+DC(i,k)*v(i,j,k,nstp) END DO END DO DO i=IstrB,IendB cff1=1.0_r8/DC(i,0) !^ tl_vbar(i,j,kstp)=tl_cff2 !^ tl_vbar(i,j,knew)=tl_cff2 !^ ad_cff2=ad_cff2+ad_vbar(i,j,knew)+ & & ad_vbar(i,j,kstp) ad_vbar(i,j,knew)=0.0_r8 ad_vbar(i,j,kstp)=0.0_r8 !^ tl_cff2=tl_cff2*vmask(i,j) !^ ad_cff2=ad_cff2*vmask(i,j) !^ tl_cff2=tl_CF(i,0)*cff1+CF(i,0)*tl_cff1 !^ ad_cff1=ad_cff1+CF(i,0)*ad_cff2 ad_CF(i,0)=ad_CF(i,0)+cff1*ad_cff2 ad_cff2=0.0_r8 !^ tl_cff1=-cff1*cff1*tl_DC(i,0) !^ ad_DC(i,0)=ad_DC(i,0)-cff1*cff1*ad_cff1 ad_cff1=0.0_r8 END DO DO k=1,N(ng) DO i=IstrB,IendB !^ tl_CF(i,0)=tl_CF(i,0)+tl_DC(i,k)*v(i,j,k,nstp)+ & !^ & DC(i,k)*tl_v(i,j,k,nstp) !^ ad_DC(i,k)=ad_DC(i,k)+v(i,j,k,nstp)*ad_CF(i,0) ad_v(i,j,k,nstp)=ad_v(i,j,k,nstp)+DC(i,k)*ad_CF(i,0) !^ tl_DC(i,0)=tl_DC(i,0)+tl_DC(i,k) !^ ad_DC(i,k)=ad_DC(i,k)+ad_DC(i,0) !^ tl_DC(i,k)=0.5_r8*(tl_Hz(i,j,k)+tl_Hz(i,j-1,k)) !^ adfac=0.5_r8*ad_DC(i,k) ad_Hz(i,j-1,k)=ad_Hz(i,j-1,k)+adfac ad_Hz(i,j ,k)=ad_Hz(i,j ,k)+adfac ad_DC(i,k)=0.0_r8 END DO END DO DO i=IstrB,IendB !^ tl_CF(i,0)=0.0_r8 !^ ad_CF(i,0)=0.0_r8 !^ tl_DC(i,0)=0.0_r8 !^ ad_DC(i,0)=0.0_r8 END DO END IF ! DO i=IstrM,IendB DC(i,0)=0.0_r8 CF(i,0)=0.0_r8 END DO DO k=1,N(ng) DO i=IstrM,IendB DC(i,k)=0.5_r8*(Hz(i,j,k)+Hz(i-1,j,k)) DC(i,0)=DC(i,0)+DC(i,k) CF(i,0)=CF(i,0)+DC(i,k)*u(i,j,k,nstp) END DO END DO DO i=IstrM,IendB cff1=1.0_r8/DC(i,0) !^ tl_ubar(i,j,kstp)=tl_cff2 !^ tl_ubar(i,j,knew)=tl_cff2 !^ ad_cff2=ad_cff2+ad_ubar(i,j,knew)+ & & ad_ubar(i,j,kstp) ad_ubar(i,j,knew)=0.0_r8 ad_ubar(i,j,kstp)=0.0_r8 !^ tl_cff2=tl_cff2*umask(i,j) !^ ad_cff2=ad_cff2*umask(i,j) !^ tl_cff2=tl_CF(i,0)*cff1+CF(i,0)*tl_cff1 !^ ad_cff1=ad_cff1+CF(i,0)*ad_cff2 ad_CF(i,0)=ad_CF(i,0)+cff1*ad_cff2 ad_cff2=0.0_r8 !^ tl_cff1=-cff1*cff1*tl_DC(i,0) !^ ad_DC(i,0)=ad_DC(i,0)-cff1*cff1*ad_cff1 ad_cff1=0.0_r8 END DO DO k=1,N(ng) DO i=IstrM,IendB !^ tl_CF(i,0)=tl_CF(i,0)+tl_DC(i,k)*u(i,j,k,nstp)+ & !^ & DC(i,k)*tl_u(i,j,k,nstp) !^ ad_DC(i,k)=ad_DC(i,k)+u(i,j,k,nstp)*ad_CF(i,0) ad_u(i,j,k,nstp)=ad_u(i,j,k,nstp)+DC(i,k)*ad_CF(i,0) !^ tl_DC(i,0)=tl_DC(i,0)+tl_DC(i,k) !^ ad_DC(i,k)=ad_DC(i,k)+ad_DC(i,0) !^ tl_DC(i,k)=0.5_r8*(tl_Hz(i,j,k)+tl_Hz(i-1,j,k)) !^ adfac=0.5_r8*ad_DC(i,k) ad_Hz(i-1,j,k)=ad_Hz(i-1,j,k)+adfac ad_Hz(i ,j,k)=ad_Hz(i ,j,k)+adfac ad_DC(i,k)=0.0_r8 END DO END DO DO i=IstrM,IendB !^ tl_CF(i,0)=0.0_r8 !^ ad_CF(i,0)=0.0_r8 !^ tl_DC(i,0)=0.0_r8 !^ ad_DC(i,0)=0.0_r8 END DO END DO ! !----------------------------------------------------------------------- ! Adjoint of initialize other time levels for 3D momentum. !----------------------------------------------------------------------- ! ! Apply boundary conditions. ! !^ CALL mp_exchange3d (ng, tile, model, 4, & !^ & LBi, UBi, LBj, UBj, 1, N(ng), & !^ & NghostPoints, & !^ & EWperiodic(ng), NSperiodic(ng), & !^ & tl_u(:,:,:,nstp), tl_v(:,:,:,nstp), & !^ & tl_u(:,:,:,nnew), tl_v(:,:,:,nnew)) !^ only nstp needed CALL ad_mp_exchange3d (ng, tile, model, 2, & & LBi, UBi, LBj, UBj, 1, N(ng), & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_u(:,:,:,nstp), ad_v(:,:,:,nstp)) ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN !^ CALL exchange_v3d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, 1, N(ng), & !^ & tl_v(:,:,:,nstp)) !^ only nstp needed CALL ad_exchange_v3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ad_v(:,:,:,nstp)) !^ CALL exchange_u3d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, 1, N(ng), & !^ & tl_u(:,:,:,nstp)) !^ only nstp needed CALL ad_exchange_u3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ad_u(:,:,:,nstp)) END IF ! !^ CALL tl_v3dbc_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, N(ng), & !^ & IminS, ImaxS, JminS, JmaxS, & !^ & nstp, nstp, & !^ & tl_v) !^ only nstp needed CALL ad_v3dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, N(ng), & & IminS, ImaxS, JminS, JmaxS, & & nstp, nstp, & & ad_v) !^ CALL tl_u3dbc_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, N(ng), & !^ & IminS, ImaxS, JminS, JmaxS, & !^ & nstp, nstp, & !^ & tl_u) !^ only nstp needed CALL ad_u3dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, N(ng), & & IminS, ImaxS, JminS, JmaxS, & & nstp, nstp, & & ad_u) ! ! Adjoint of initialize other time levels for momentum. ! DO j=JstrB,JendB IF (j.ge.JstrM) THEN DO k=1,N(ng) DO i=IstrB,IendB !^ tl_v(i,j,k,nstp)=tl_cff2 !^ ad_cff2=ad_cff2+ad_v(i,j,k,nstp) ad_v(i,j,k,nstp)=0.0_r8 !^ tl_cff2=tl_cff2*vmask(i,j) !^ ad_cff2=ad_cff2*vmask(i,j) !^ tl_cff2=tl_v(i,j,k,nstp) !^ ad_v(i,j,k,nstp)=ad_v(i,j,k,nstp)+ad_cff2 ad_cff2=0.0_r8 END DO END DO END IF DO k=1,N(ng) DO i=IstrM,IendB !^ tl_u(i,j,k,nstp)=tl_cff1 !^ ad_cff1=ad_cff1+ad_u(i,j,k,nstp) ad_u(i,j,k,nstp)=0.0_r8 !^ tl_cff1=tl_cff1*umask(i,j) !^ ad_cff1=ad_cff1*umask(i,j) !^ tl_cff1=tl_u(i,j,k,nstp) !^ ad_u(i,j,k,nstp)=ad_u(i,j,k,nstp)+ad_cff1 ad_cff1=0.0_r8 END DO END DO END DO ! RETURN END SUBROUTINE ad_ini_fields_tile ! !*********************************************************************** SUBROUTINE ad_ini_zeta (ng, tile, model) !*********************************************************************** ! USE mod_param USE mod_grid USE mod_coupling USE mod_ocean USE mod_stepping ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model ! ! Local variable declarations. ! character (len=*), parameter :: MyFile = & & "ROMS/Adjoint/ad_ini_fields.F"//", ad_ini_zeta" ! 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 wclock_on (ng, iADM, 2, 990, MyFile) CALL ad_ini_zeta_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & kstp(ng), krhs(ng), knew(ng), & & GRID(ng) % rmask, & & COUPLING(ng) % ad_Zt_avg1, & & OCEAN(ng) % zeta, & & OCEAN(ng) % ad_zeta_sol, & & OCEAN(ng) % ad_zeta) CALL wclock_off (ng, iADM, 2, 1014, MyFile) ! RETURN END SUBROUTINE ad_ini_zeta ! !*********************************************************************** SUBROUTINE ad_ini_zeta_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & kstp, krhs, knew, & & rmask, & & ad_Zt_avg1, & & zeta, ad_zeta_sol, ad_zeta) !*********************************************************************** ! USE mod_param USE mod_ncparam USE mod_scalars ! USE ad_exchange_2d_mod, ONLY : ad_exchange_r2d_tile USE mp_exchange_mod, ONLY : ad_mp_exchange2d USE ad_zetabc_mod, ONLY : ad_zetabc_tile ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model integer, intent(in) :: LBi, UBi, LBj, UBj integer, intent(in) :: IminS, ImaxS, JminS, JmaxS integer, intent(in) :: kstp, krhs, knew ! real(r8), intent(in) :: rmask(LBi:,LBj:) real(r8), intent(in) :: zeta(LBi:,LBj:,:) real(r8), intent(inout) :: ad_Zt_avg1(LBi:,LBj:) real(r8), intent(inout) :: ad_zeta_sol(LBi:,LBj:) real(r8), intent(inout) :: ad_zeta(LBi:,LBj:,:) ! ! Local variable declarations. ! integer :: Imin, Imax, Jmin, Jmax integer :: i, j, kbed real(r8) :: ad_cff1 ! !----------------------------------------------------------------------- ! 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_cff1=0.0_r8 ! !----------------------------------------------------------------------- ! Initialize fast-time averaged free-surface (Zt_avg1) with the inital ! free-surface !----------------------------------------------------------------------- ! !^ CALL mp_exchange2d (ng, tile, model, 1, & !^ & LBi, UBi, LBj, UBj, & !^ & NghostPoints, & !^ & EWperiodic(ng), NSperiodic(ng), & !^ & tl_Zt_avg1) !^ CALL ad_mp_exchange2d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_Zt_avg1) ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN !^ CALL exchange_r2d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & tl_Zt_avg1) !^ CALL ad_exchange_r2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_Zt_avg1) END IF ! DO j=JstrT,JendT DO i=IstrT,IendT !^ tl_Zt_avg1(i,j)=tl_zeta(i,j,kstp) !^ ad_zeta(i,j,kstp)=ad_zeta(i,j,kstp)+ad_Zt_avg1(i,j) ad_Zt_avg1(i,j)=0.0_r8 END DO END DO ! !----------------------------------------------------------------------- ! Adjoint of initialize other time levels for free-surface. !----------------------------------------------------------------------- ! ! Apply boundary conditions. ! !^ CALL mp_exchange2d (ng, tile, model, 2, & !^ & LBi, UBi, LBj, UBj, & !^ & NghostPoints, & !^ & EWperiodic(ng), NSperiodic(ng), & !^ & tl_zeta(:,:,kstp), & !^ & tl_zeta(:,:,knew)) !^ CALL ad_mp_exchange2d (ng, tile, model, 2, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_zeta(:,:,kstp), & & ad_zeta(:,:,knew)) ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN !^ CALL exchange_r2d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & tl_zeta(:,:,knew)) !^ CALL ad_exchange_r2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_zeta(:,:,knew)) !^ CALL exchange_r2d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & tl_zeta(:,:,kstp)) !^ CALL ad_exchange_r2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_zeta(:,:,kstp)) END IF ! IF (.not.(ANY(ad_LBC(:,isFsur,ng)%radiation).or. & & ANY(ad_LBC(:,isFsur,ng)%Chapman_explicit).or. & & ANY(ad_LBC(:,isFsur,ng)%Chapman_implicit))) THEN !^ CALL tl_zetabc_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & IminS, ImaxS, JminS, JmaxS, & !^ & krhs, kstp, knew, & !^ & zeta, & !^ & tl_zeta) !^ CALL ad_zetabc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, knew, & & zeta, & & ad_zeta) !^ CALL tl_zetabc_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & IminS, ImaxS, JminS, JmaxS, & !^ & krhs, kstp, kstp, & !^ & zeta, & !^ & tl_zeta) !^ CALL ad_zetabc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, kstp, & & zeta, & & ad_zeta) END IF ! ! Adjoint of free-surface initialization. ! IF (.not.(ANY(ad_LBC(:,isFsur,ng)%radiation).or. & & ANY(ad_LBC(:,isFsur,ng)%Chapman_explicit).or. & & ANY(ad_LBC(:,isFsur,ng)%Chapman_implicit))) THEN Imin=IstrB Imax=IendB Jmin=JstrB Jmax=JendB ELSE Imin=IstrT Imax=IendT Jmin=JstrT Jmax=JendT END IF DO j=Jmin,Jmax DO i=Imin,Imax !^ tl_zeta(i,j,kstp)=tl_cff1 !^ tl_zeta(i,j,knew)=tl_cff1 !^ ad_cff1=ad_cff1+ad_zeta(i,j,knew)+ad_zeta(i,j,kstp) ad_zeta(i,j,knew)=0.0_r8 ad_zeta(i,j,kstp)=0.0_r8 !^ tl_cff1=tl_cff1*rmask(i,j) !^ ad_cff1=ad_cff1*rmask(i,j) !^ tl_cff1=tl_zeta(i,j,kstp) !^ ad_zeta(i,j,kstp)=ad_zeta(i,j,kstp)+ad_cff1 ad_cff1=0.0_r8 END DO END DO ! !----------------------------------------------------------------------- ! Load free-surface adjoint solution into IO arrays. !----------------------------------------------------------------------- ! DO j=JstrT,JendT DO i=IstrT,IendT ad_zeta_sol(i,j)=ad_zeta(i,j,kstp) END DO END DO ! RETURN END SUBROUTINE ad_ini_zeta_tile ! !*********************************************************************** SUBROUTINE ad_out_fields (ng, tile, model) !*********************************************************************** ! USE mod_param USE mod_grid USE mod_coupling USE mod_ocean USE mod_stepping ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model ! ! Local variable declarations. ! character (len=*), parameter :: MyFile = & & "ROMS/Adjoint/ad_ini_fields.F"//", ad_out_fields" ! 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 wclock_on (ng, iADM, 2, 1453, MyFile) CALL ad_out_fields_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & kstp(ng), krhs(ng), knew(ng), & & nstp(ng), nnew(ng), & & GRID(ng) % rmask, & & GRID(ng) % umask, & & GRID(ng) % vmask, & & OCEAN(ng) % ad_u, & & OCEAN(ng) % ad_u_sol, & & OCEAN(ng) % ad_v, & & OCEAN(ng) % ad_v_sol, & & OCEAN(ng) % ad_t, & & OCEAN(ng) % ad_t_sol, & & OCEAN(ng) % zeta, & & OCEAN(ng) % ubar, & & OCEAN(ng) % vbar, & & OCEAN(ng) % ad_ubar_sol, & & OCEAN(ng) % ad_ubar, & & OCEAN(ng) % ad_vbar_sol, & & OCEAN(ng) % ad_vbar, & & OCEAN(ng) % ad_zeta) CALL wclock_off (ng, iADM, 2, 1482, MyFile) ! RETURN END SUBROUTINE ad_out_fields ! !*********************************************************************** SUBROUTINE ad_out_fields_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & kstp, krhs, knew, & & nstp, nnew, & & rmask, umask, vmask, & & ad_u, ad_u_sol, & & ad_v, ad_v_sol, & & ad_t, ad_t_sol, & & zeta, ubar, vbar, & & ad_ubar_sol, ad_ubar, & & ad_vbar_sol, ad_vbar, ad_zeta) !*********************************************************************** ! USE mod_param USE mod_ncparam USE mod_scalars ! USE ad_exchange_2d_mod USE ad_exchange_3d_mod USE mp_exchange_mod, ONLY : ad_mp_exchange2d USE mp_exchange_mod, ONLY : ad_mp_exchange3d, ad_mp_exchange4d USE ad_set_depth_mod, ONLY : ad_set_depth_tile USE ad_t3dbc_mod, ONLY : ad_t3dbc_tile USE ad_u3dbc_mod, ONLY : ad_u3dbc_tile USE ad_v3dbc_mod, ONLY : ad_v3dbc_tile USE ad_u2dbc_mod, ONLY : ad_u2dbc_tile USE ad_v2dbc_mod, ONLY : ad_v2dbc_tile ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model integer, intent(in) :: LBi, UBi, LBj, UBj integer, intent(in) :: IminS, ImaxS, JminS, JmaxS integer, intent(in) :: kstp, krhs, knew integer, intent(in) :: nstp, nnew ! real(r8), intent(in) :: rmask(LBi:,LBj:) real(r8), intent(in) :: umask(LBi:,LBj:) real(r8), intent(in) :: vmask(LBi:,LBj:) real(r8), intent(in) :: ubar(LBi:,LBj:,:) real(r8), intent(in) :: vbar(LBi:,LBj:,:) real(r8), intent(in) :: zeta(LBi:,LBj:,:) real(r8), intent(inout) :: ad_ubar_sol(LBi:,LBj:) real(r8), intent(inout) :: ad_vbar_sol(LBi:,LBj:) real(r8), intent(inout) :: ad_zeta(LBi:,LBj:,:) real(r8), intent(inout) :: ad_ubar(LBi:,LBj:,:) real(r8), intent(inout) :: ad_vbar(LBi:,LBj:,:) real(r8), intent(inout) :: ad_t(LBi:,LBj:,:,:,:) real(r8), intent(inout) :: ad_u(LBi:,LBj:,:,:) real(r8), intent(inout) :: ad_v(LBi:,LBj:,:,:) real(r8), intent(inout) :: ad_t_sol(LBi:,LBj:,:,:) real(r8), intent(inout) :: ad_u_sol(LBi:,LBj:,:) real(r8), intent(inout) :: ad_v_sol(LBi:,LBj:,:) ! ! Local variable declarations. ! integer :: i, ic, itrc, j, k, kbed, kout real(r8) :: cff1 real(r8) :: adfac, ad_cff1, ad_cff2 ! !----------------------------------------------------------------------- ! 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 ! kout=knew ! CALL ad_mp_exchange2d (ng, tile, model, 4, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_ubar(:,:,kstp), ad_vbar(:,:,kstp), & & ad_ubar(:,:,kout), ad_vbar(:,:,kout)) ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN CALL ad_exchange_v2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_vbar(:,:,kout)) CALL ad_exchange_u2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_ubar(:,:,kout)) CALL ad_exchange_v2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_vbar(:,:,kstp)) CALL ad_exchange_u2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_ubar(:,:,kstp)) END IF ! IF (.not.(ANY(ad_LBC(:,isUbar,ng)%radiation).or. & & ANY(ad_LBC(:,isVbar,ng)%radiation).or. & & ANY(ad_LBC(:,isUbar,ng)%Flather).or. & & ANY(ad_LBC(:,isVbar,ng)%Flather))) THEN CALL ad_v2dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, knew, & & ubar, vbar, zeta, & & ad_ubar, ad_vbar, ad_zeta) CALL ad_u2dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, knew, & & ubar, vbar, zeta, & & ad_ubar, ad_vbar, ad_zeta) CALL ad_v2dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, kstp, & & ubar, vbar, zeta, & & ad_ubar, ad_vbar, ad_zeta) CALL ad_u2dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, kstp, & & ubar, vbar, zeta, & & ad_ubar, ad_vbar, ad_zeta) END IF ! ! Load 2D adjoint momentum solution into IO arrays. ! DO j=JstrT,JendT DO i=IstrP,IendT ad_ubar_sol(i,j)=ad_ubar(i,j,kstp)+ad_ubar(i,j,kout) END DO IF (j.ge.JstrP) THEN DO i=IstrT,IendT ad_vbar_sol(i,j)=ad_vbar(i,j,kstp)+ad_vbar(i,j,kout) END DO END IF END DO ! CALL ad_mp_exchange3d (ng, tile, model, 4, & & LBi, UBi, LBj, UBj, 1, N(ng), & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_u(:,:,:,nstp), ad_v(:,:,:,nstp), & & ad_u(:,:,:,nnew), ad_v(:,:,:,nnew)) IF (EWperiodic(ng).or.NSperiodic(ng)) THEN CALL ad_exchange_v3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ad_v(:,:,:,nstp)) CALL ad_exchange_u3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ad_u(:,:,:,nstp)) CALL ad_exchange_v3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ad_v(:,:,:,nnew)) CALL ad_exchange_u3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ad_u(:,:,:,nnew)) END IF CALL ad_v3dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, N(ng), & & IminS, ImaxS, JminS, JmaxS, & & nstp, nnew, & & ad_v) CALL ad_v3dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, N(ng), & & IminS, ImaxS, JminS, JmaxS, & & nstp, nstp, & & ad_v) CALL ad_u3dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, N(ng), & & IminS, ImaxS, JminS, JmaxS, & & nstp, nnew, & & ad_u) CALL ad_u3dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, N(ng), & & IminS, ImaxS, JminS, JmaxS, & & nstp, nstp, & & ad_u) ! ! Add second piece of the 3D adjoint momentum solution into IO arrays. ! The first is loaded in "ad_step3d_uv". ! DO j=JstrB,JendB IF (j.ge.JstrM) THEN DO k=1,N(ng) DO i=IstrB,IendB ad_v_sol(i,j,k)=ad_v_sol(i,j,k)+ad_v(i,j,k,nnew) END DO END DO END IF DO k=1,N(ng) DO i=IstrM,IendB ad_u_sol(i,j,k)=ad_u_sol(i,j,k)+ad_u(i,j,k,nnew) END DO END DO END DO CALL ad_mp_exchange4d (ng, tile, model, 2, & & LBi, UBi, LBj, UBj, 1, N(ng), 1, NT(ng), & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_t(:,:,:,nstp,:), & & ad_t(:,:,:,nnew,:)) ! ic=0 DO itrc=1,NT(ng) IF (LtracerCLM(itrc,ng).and.LnudgeTCLM(itrc,ng)) THEN ic=ic+1 ! OBC nudging coefficient index END IF IF (EWperiodic(ng).or.NSperiodic(ng)) THEN CALL ad_exchange_r3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ad_t(:,:,:,nstp,itrc)) CALL ad_exchange_r3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ad_t(:,:,:,nnew,itrc)) END IF CALL ad_t3dbc_tile (ng, tile, itrc, ic, & & LBi, UBi, LBj, UBj, N(ng), NT(ng), & & IminS, ImaxS, JminS, JmaxS, & & nstp, nnew, & & ad_t) CALL ad_t3dbc_tile (ng, tile, itrc, ic, & & LBi, UBi, LBj, UBj, N(ng), NT(ng), & & IminS, ImaxS, JminS, JmaxS, & & nstp, nstp, & & ad_t) DO k=1,N(ng) DO j=JstrB,JendB DO i=IstrB,IendB ad_t_sol(i,j,k,itrc)=ad_t(i,j,k,nstp,itrc)+ & & ad_t(i,j,k,nnew,itrc) END DO END DO END DO END DO ! RETURN END SUBROUTINE ad_out_fields_tile ! !*********************************************************************** SUBROUTINE ad_out_zeta (ng, tile, model) !*********************************************************************** ! USE mod_param USE mod_grid USE mod_coupling USE mod_ocean USE mod_stepping ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model ! ! Local variable declarations. ! character (len=*), parameter :: MyFile = & & "ROMS/Adjoint/ad_ini_fields.F"//", ad_out_zeta" ! 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 wclock_on (ng, iADM, 2, 1942, MyFile) CALL ad_out_zeta_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & kstp(ng), krhs(ng), knew(ng), & & GRID(ng) % rmask, & & COUPLING(ng) % ad_Zt_avg1, & & OCEAN(ng) % zeta, & & OCEAN(ng) % ad_zeta_sol, & & OCEAN(ng) % ad_zeta) CALL wclock_off (ng, iADM, 2, 1966, MyFile) ! RETURN END SUBROUTINE ad_out_zeta ! !*********************************************************************** SUBROUTINE ad_out_zeta_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & kstp, krhs, knew, & & rmask, & & ad_Zt_avg1, & & zeta, ad_zeta_sol, ad_zeta) !*********************************************************************** ! USE mod_param USE mod_ncparam USE mod_scalars ! USE ad_exchange_2d_mod, ONLY : ad_exchange_r2d_tile USE mp_exchange_mod, ONLY : ad_mp_exchange2d USE ad_zetabc_mod, ONLY : ad_zetabc_tile ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model integer, intent(in) :: LBi, UBi, LBj, UBj integer, intent(in) :: IminS, ImaxS, JminS, JmaxS integer, intent(in) :: kstp, krhs, knew ! real(r8), intent(in) :: rmask(LBi:,LBj:) real(r8), intent(in) :: zeta(LBi:,LBj:,:) real(r8), intent(inout) :: ad_Zt_avg1(LBi:,LBj:) real(r8), intent(inout) :: ad_zeta_sol(LBi:,LBj:) real(r8), intent(inout) :: ad_zeta(LBi:,LBj:,:) ! ! Local variable declarations. ! integer :: Imin, Imax, Jmin, Jmax integer :: i, j, kbed, kout real(r8) :: ad_cff1 ! !----------------------------------------------------------------------- ! 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 ! kout=knew ! !----------------------------------------------------------------------- ! Load free-surface adjoint solution into IO arrays. !----------------------------------------------------------------------- ! CALL ad_mp_exchange2d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_Zt_avg1) ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN CALL ad_exchange_r2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_Zt_avg1) END IF ! CALL ad_mp_exchange2d (ng, tile, model, 2, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_zeta(:,:,kstp), & & ad_zeta(:,:,kout)) ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN CALL ad_exchange_r2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_zeta(:,:,kout)) CALL ad_exchange_r2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_zeta(:,:,kstp)) END IF ! IF (.not.(ANY(ad_LBC(:,isFsur,ng)%radiation).or. & & ANY(ad_LBC(:,isFsur,ng)%Chapman_explicit).or. & & ANY(ad_LBC(:,isFsur,ng)%Chapman_implicit))) THEN CALL ad_zetabc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, kout, & & zeta, & & ad_zeta) CALL ad_zetabc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, kstp, & & zeta, & & ad_zeta) END IF ! DO j=JstrT,JendT DO i=IstrT,IendT ad_zeta_sol(i,j)=ad_zeta(i,j,kout)+ad_zeta(i,j,kstp)+ & & ad_Zt_avg1(i,j) END DO END DO ! RETURN END SUBROUTINE ad_out_zeta_tile END MODULE ad_ini_fields_mod