MODULE ad_set_vbc_mod ! !git $Id$ !svn $Id: ad_set_vbc.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 module sets vertical boundary conditons for adjoint momentum ! ! and tracers. ! ! ! ! BASIC STATE variables needed: btflx, stflx, dqdt, t, sss, ! ! z_r, v, u ! ! (btflx and stflx are over written) ! ! ! !======================================================================= ! implicit none ! PRIVATE PUBLIC :: ad_set_vbc ! CONTAINS ! !*********************************************************************** SUBROUTINE ad_set_vbc (ng, tile) !*********************************************************************** ! USE mod_param USE mod_grid USE mod_forces USE mod_ocean USE mod_stepping ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile ! ! Local variable declarations. ! character (len=*), parameter :: MyFile = & & "ROMS/Adjoint/ad_set_vbc.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, 6, 53, MyFile) CALL ad_set_vbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & GRID(ng) % rdrag2, & & OCEAN(ng) % t, & & OCEAN(ng) % ad_t, & & OCEAN(ng) % u, & & OCEAN(ng) % ad_u, & & OCEAN(ng) % v, & & OCEAN(ng) % ad_v, & & FORCES(ng) % ad_bustr, & & FORCES(ng) % ad_bustr_sol, & & FORCES(ng) % ad_bvstr, & & FORCES(ng) % ad_bvstr_sol, & & FORCES(ng) % stflux, & & FORCES(ng) % btflux, & & FORCES(ng) % stflx, & & FORCES(ng) % ad_stflx, & & FORCES(ng) % btflx, & & FORCES(ng) % ad_btflx, & & nrhs(ng)) CALL wclock_off (ng, iADM, 6, 120, MyFile) ! RETURN END SUBROUTINE ad_set_vbc ! !*********************************************************************** SUBROUTINE ad_set_vbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & rdrag2, & & t, ad_t, & & u, ad_u, & & v, ad_v, & & ad_bustr, ad_bustr_sol, & & ad_bvstr, ad_bvstr_sol, & & stflux, btflux, & & stflx, ad_stflx, & & btflx, ad_btflx, & & nrhs) !*********************************************************************** ! USE mod_param USE mod_scalars ! USE ad_bc_2d_mod USE mp_exchange_mod, ONLY : ad_mp_exchange2d ! ! 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) :: nrhs ! real(r8), intent(in) :: rdrag2(LBi:,LBj:) real(r8), intent(in) :: t(LBi:,LBj:,:,:,:) real(r8), intent(in) :: u(LBi:,LBj:,:,:) real(r8), intent(in) :: v(LBi:,LBj:,:,:) real(r8), intent(in) :: stflx(LBi:,LBj:,:) real(r8), intent(in) :: btflx(LBi:,LBj:,:) real(r8), intent(in) :: stflux(LBi:,LBj:,:) real(r8), intent(in) :: btflux(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_bustr(LBi:,LBj:) real(r8), intent(inout) :: ad_bvstr(LBi:,LBj:) real(r8), intent(inout) :: ad_stflx(LBi:,LBj:,:) real(r8), intent(inout) :: ad_btflx(LBi:,LBj:,:) real(r8), intent(out) :: ad_bustr_sol(LBi:,LBj:) real(r8), intent(out) :: ad_bvstr_sol(LBi:,LBj:) ! ! Local variable declarations. ! integer :: i, j, itrc ! real(r8) :: EmP, ad_EmP real(r8) :: cff1, cff2, cff3 real(r8) :: ad_cff1, ad_cff2, ad_cff3, adfac, adfac1, adfac2 ! ! !----------------------------------------------------------------------- ! 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_EmP=0.0_r8 ad_cff1=0.0_r8 ad_cff2=0.0_r8 ad_cff3=0.0_r8 ! !----------------------------------------------------------------------- ! Set kinematic bottom momentum flux (m2/s2). !----------------------------------------------------------------------- ! ! Apply boundary conditions. ! !^ CALL mp_exchange2d (ng, tile, iTLM, 2, & !^ & LBi, UBi, LBj, UBj, & !^ & NghostPoints, & !^ & EWperiodic(ng), NSperiodic(ng), & !^ & tl_bustr, tl_bvstr) !^ CALL ad_mp_exchange2d (ng, tile, iADM, 2, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_bustr, ad_bvstr) !^ CALL bc_u2d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & tl_bustr) !^ CALL ad_bc_u2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_bustr) !^ CALL bc_v2d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & tl_bvstr) !^ CALL ad_bc_v2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_bvstr) ! ! Save adjoint bottom momentum flux for output purposes. ! DO j=JstrR,JendR DO i=Istr,IendR ad_bustr_sol(i,j)=ad_bustr(i,j) END DO IF (j.ge.Jstr) THEN DO i=Istr,IendR ad_bvstr_sol(i,j)=ad_bvstr(i,j) END DO END IF END DO ! !----------------------------------------------------------------------- ! Set kinematic bottom momentum flux (m2/s2). !----------------------------------------------------------------------- ! ! Set quadratic bottom stress. ! DO j=JstrV,Jend DO i=Istr,Iend cff1=0.25_r8*(u(i ,j ,1,nrhs)+ & & u(i+1,j ,1,nrhs)+ & & u(i ,j-1,1,nrhs)+ & & u(i+1,j-1,1,nrhs)) cff2=SQRT(cff1*cff1+v(i,j,1,nrhs)*v(i,j,1,nrhs)) !^ tl_bvstr(i,j)=0.5_r8*(rdrag2(i,j-1)+rdrag2(i,j))* & !^ & (tl_v(i,j,1,nrhs)*cff2+ & !^ & v(i,j,1,nrhs)*tl_cff2) !^ adfac=0.5_r8*(rdrag2(i,j-1)+rdrag2(i,j))* & & ad_bvstr(i,j) ad_v(i,j,1,nrhs)=ad_v(i,j,1,nrhs)+adfac*cff2 ad_cff2=ad_cff2+adfac*v(i,j,1,nrhs) ad_bvstr(i,j)=0.0_r8 IF (cff2.ne.0.0_r8) THEN !^ tl_cff2=(cff1*tl_cff1+v(i,j,1,nrhs)*tl_v(i,j,1,nrhs))/cff2 !^ adfac=ad_cff2/cff2 ad_v(i,j,1,nrhs)=ad_v(i,j,1,nrhs)+adfac*v(i,j,1,nrhs) ad_cff1=ad_cff1+adfac*cff1 ad_cff2=0.0_r8 ELSE !^ tl_cff2=0.0_r8 !^ ad_cff2=0.0_r8 END IF !^ tl_cff1=0.25_r8*(tl_u(i ,j ,1,nrhs)+ & !^ & tl_u(i+1,j ,1,nrhs)+ & !^ & tl_u(i ,j-1,1,nrhs)+ & !^ & tl_u(i+1,j-1,1,nrhs)) !^ adfac=0.25_r8*ad_cff1 ad_u(i ,j-1,1,nrhs)=ad_u(i ,j-1,1,nrhs)+adfac ad_u(i+1,j-1,1,nrhs)=ad_u(i+1,j-1,1,nrhs)+adfac ad_u(i ,j ,1,nrhs)=ad_u(i ,j ,1,nrhs)+adfac ad_u(i+1,j ,1,nrhs)=ad_u(i+1,j ,1,nrhs)+adfac ad_cff1=0.0_r8 END DO END DO DO j=Jstr,Jend DO i=IstrU,Iend cff1=0.25_r8*(v(i ,j ,1,nrhs)+ & & v(i ,j+1,1,nrhs)+ & & v(i-1,j ,1,nrhs)+ & & v(i-1,j+1,1,nrhs)) cff2=SQRT(u(i,j,1,nrhs)*u(i,j,1,nrhs)+cff1*cff1) !^ tl_bustr(i,j)=0.5_r8*(rdrag2(i-1,j)+rdrag2(i,j))* & !^ & (tl_u(i,j,1,nrhs)*cff2+ & !^ & u(i,j,1,nrhs)*tl_cff2) !^ adfac=0.5_r8*(rdrag2(i-1,j)+rdrag2(i,j))* & & ad_bustr(i,j) ad_u(i,j,1,nrhs)=ad_u(i,j,1,nrhs)+adfac*cff2 ad_cff2=ad_cff2+adfac*u(i,j,1,nrhs) ad_bustr(i,j)=0.0_r8 IF (cff2.ne.0.0_r8) THEN !^ tl_cff2=(u(i,j,1,nrhs)*tl_u(i,j,1,nrhs)+cff1*tl_cff1)/cff2 !^ adfac=ad_cff2/cff2 ad_u(i,j,1,nrhs)=ad_u(i,j,1,nrhs)+adfac*u(i,j,1,nrhs) ad_cff1=ad_cff1+adfac*cff1 ad_cff2=0.0_r8 ELSE !^ tl_cff2=0.0_r8 !^ ad_cff2=0.0_r8 END IF !^ tl_cff1=0.25_r8*(tl_v(i ,j ,1,nrhs)+ & !^ & tl_v(i ,j+1,1,nrhs)+ & !^ & tl_v(i-1,j ,1,nrhs)+ & !^ & tl_v(i-1,j+1,1,nrhs)) !^ adfac=0.25_r8*ad_cff1 ad_v(i-1,j ,1,nrhs)=ad_v(i-1,j ,1,nrhs)+adfac ad_v(i ,j ,1,nrhs)=ad_v(i ,j ,1,nrhs)+adfac ad_v(i-1,j+1,1,nrhs)=ad_v(i-1,j+1,1,nrhs)+adfac ad_v(i ,j+1,1,nrhs)=ad_v(i ,j+1,1,nrhs)+adfac ad_cff1=0.0_r8 END DO END DO ! !----------------------------------------------------------------------- ! Add salt flux correction or multiply flux by salinity. !----------------------------------------------------------------------- ! DO j=JstrR,JendR DO i=IstrR,IendR EmP=stflux(i,j,isalt) !^ tl_btflx(i,j,isalt)=btflx(i,j,isalt)* & !^ & tl_t(i,j,1,nrhs,isalt) !^ ad_t(i,j,1,nrhs,isalt)=ad_t(i,j,1,nrhs,isalt)+ & & btflx(i,j,isalt)*ad_btflx(i,j,isalt) ad_btflx(i,j,isalt)=0.0_r8 !^ tl_stflx(i,j,isalt)=EmP*tl_t(i,j,N(ng),nrhs,isalt)+ & !^ & tl_EmP*t(i,j,N(ng),nrhs,isalt) !^ ad_t(i,j,N(ng),nrhs,isalt)=ad_t(i,j,N(ng),nrhs,isalt)+ & & EmP*ad_stflx(i,j,isalt) ad_EmP=ad_EmP+ & t(i,j,N(ng),nrhs,isalt)*ad_stflx(i,j,isalt) ad_stflx(i,j,isalt)=0.0_r8 !^ tl_EmP=0.0_r8 !^ ad_EmP=0.0_r8 END DO END DO ! !----------------------------------------------------------------------- ! Adjoint of loading surface and bottom net heat flux (degC m/s) into ! state variables "stflx" and "btflx". ! ! Notice that the forcing net heat flux stflux(:,:,itemp) is processed ! elsewhere from data, coupling, bulk flux parameterization, ! or analytical formulas. ! ! During model coupling, we need to make sure that this forcing is ! unaltered during the coupling interval when ROMS timestep size is ! smaller. Notice that further manipulations to state variable ! stflx(:,:,itemp) are allowed below. !----------------------------------------------------------------------- ! DO j=JstrR,JendR DO i=IstrR,IendR !^ tl_btflx(i,j,itemp)=0.0_r8 ! not needed in TLM !^ !^ tl_stflx(i,j,itemp)=0.0_r8 ! not needed in TLM !^ END DO END DO ! RETURN END SUBROUTINE ad_set_vbc_tile END MODULE ad_set_vbc_mod