MODULE ad_rho_eos_mod ! !git $Id$ !svn $Id: ad_rho_eos.F 1188 2023-08-03 19:26:47Z 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 computes "in situ" density and other associated ! ! quantitites as a function of potential temperature, salinity, ! ! and pressure from a polynomial expression (Jackett and McDougall, ! ! 1992). The polynomial expression was found from fitting to 248 ! ! values in the oceanographic ranges of salinity, potential ! ! temperature, and pressure. It assumes no pressure variation ! ! along geopotential surfaces, that is, depth (meters; negative) ! ! and pressure (dbar; assumed negative here) are interchangeable. ! ! ! ! Check Values: (T=3 C, S=35.5 PSU, Z=-5000 m) ! ! ! ! alpha = 2.1014611551470d-04 (1/Celsius) ! ! beta = 7.2575037309946d-04 (1/PSU) ! ! gamma = 3.9684764511766d-06 (1/Pa) ! ! den = 1050.3639165364 (kg/m3) ! ! den1 = 1028.2845117925 (kg/m3) ! ! sound = 1548.8815240223 (m/s) ! ! bulk = 23786.056026320 (Pa) ! ! ! ! Reference: ! ! ! ! Jackett, D. R. and T. J. McDougall, 1995, Minimal Adjustment of ! ! Hydrostatic Profiles to Achieve Static Stability, J. of Atmos. ! ! and Oceanic Techn., vol. 12, pp. 381-389. ! ! ! !======================================================================= ! implicit none ! PRIVATE PUBLIC :: ad_rho_eos ! CONTAINS ! !*********************************************************************** SUBROUTINE ad_rho_eos (ng, tile, model) !*********************************************************************** ! USE mod_param USE mod_parallel USE mod_coupling USE mod_grid USE mod_mixing 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_rho_eos.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, model, 14, 72, MyFile) CALL ad_rho_eos_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nrhs(ng), & & GRID(ng) % rmask, & & GRID(ng) % z_r, & & GRID(ng) % ad_z_r, & & OCEAN(ng) % t, & & OCEAN(ng) % ad_t, & & OCEAN(ng) % rho, & & OCEAN(ng) % ad_rho) CALL wclock_off (ng, model, 14, 114, MyFile) ! RETURN END SUBROUTINE ad_rho_eos ! !*********************************************************************** SUBROUTINE ad_rho_eos_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nrhs, & & rmask, & & z_r, ad_z_r, & & t, ad_t, & & rho, ad_rho) !*********************************************************************** ! USE mod_param USE mod_eoscoef USE mod_scalars ! USE ad_exchange_2d_mod USE ad_exchange_3d_mod USE mp_exchange_mod, ONLY : ad_mp_exchange2d, ad_mp_exchange3d ! ! 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) :: nrhs ! real(r8), intent(in) :: rmask(LBi:,LBj:) real(r8), intent(in) :: z_r(LBi:,LBj:,:) real(r8), intent(in) :: t(LBi:,LBj:,:,:,:) real(r8), intent(in) :: rho(LBi:,LBj:,:) real(r8), intent(inout) :: ad_z_r(LBi:,LBj:,:) real(r8), intent(inout) :: ad_t(LBi:,LBj:,:,:,:) real(r8), intent(inout) :: ad_rho(LBi:,LBj:,:) ! ! Local variable declarations. ! integer :: i, ised, itrc, j, k real(r8) :: SedDen, Tp, Tpr10, Ts, Tt, sqrtTs real(r8) :: ad_SedDen, ad_Tp, ad_Tpr10, ad_Ts, ad_Tt, ad_sqrtTs real(r8) :: cff, cff1, cff2, cff3 real(r8) :: ad_cff, ad_cff1, ad_cff2, ad_cff3 real(r8) :: adfac, adfac1, adfac2, adfac3 real(r8), dimension(0:9) :: C real(r8), dimension(0:9) :: ad_C real(r8), dimension(0:9) :: dCdT(0:9) real(r8), dimension(0:9) :: ad_dCdT(0:9) real(r8), dimension(0:9) :: d2Cd2T(0:9) real(r8), dimension(IminS:ImaxS,N(ng)) :: DbulkDS real(r8), dimension(IminS:ImaxS,N(ng)) :: DbulkDT real(r8), dimension(IminS:ImaxS,N(ng)) :: Dden1DS real(r8), dimension(IminS:ImaxS,N(ng)) :: Dden1DT real(r8), dimension(IminS:ImaxS,N(ng)) :: Scof real(r8), dimension(IminS:ImaxS,N(ng)) :: Tcof real(r8), dimension(IminS:ImaxS,N(ng)) :: wrk real(r8), dimension(IminS:ImaxS,N(ng)) :: ad_DbulkDS real(r8), dimension(IminS:ImaxS,N(ng)) :: ad_DbulkDT real(r8), dimension(IminS:ImaxS,N(ng)) :: ad_Dden1DS real(r8), dimension(IminS:ImaxS,N(ng)) :: ad_Dden1DT real(r8), dimension(IminS:ImaxS,N(ng)) :: ad_Scof real(r8), dimension(IminS:ImaxS,N(ng)) :: ad_Tcof real(r8), dimension(IminS:ImaxS,N(ng)) :: ad_wrk real(r8), dimension(IminS:ImaxS,N(ng)) :: bulk real(r8), dimension(IminS:ImaxS,N(ng)) :: bulk0 real(r8), dimension(IminS:ImaxS,N(ng)) :: bulk1 real(r8), dimension(IminS:ImaxS,N(ng)) :: bulk2 real(r8), dimension(IminS:ImaxS,N(ng)) :: den real(r8), dimension(IminS:ImaxS,N(ng)) :: den1 real(r8), dimension(IminS:ImaxS,N(ng)) :: ad_bulk real(r8), dimension(IminS:ImaxS,N(ng)) :: ad_bulk0 real(r8), dimension(IminS:ImaxS,N(ng)) :: ad_bulk1 real(r8), dimension(IminS:ImaxS,N(ng)) :: ad_bulk2 real(r8), dimension(IminS:ImaxS,N(ng)) :: ad_den real(r8), dimension(IminS:ImaxS,N(ng)) :: ad_den1 ! !----------------------------------------------------------------------- ! 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_Tt=0.0_r8 ad_Ts=0.0_r8 ad_Tp=0.0_r8 ad_Tpr10=0.0_r8 ad_sqrtTs=0.0_r8 ad_cff=0.0_r8 ad_cff1=0.0_r8 ad_cff2=0.0_r8 ad_cff3=0.0_r8 ad_C=0.0_r8 ad_dCdT=0.0_r8 DO k=1,N(ng) DO i=IminS,ImaxS ad_DbulkDS(i,k)=0.0_r8 ad_DbulkDT(i,k)=0.0_r8 ad_Dden1DS(i,k)=0.0_r8 ad_Dden1DT(i,k)=0.0_r8 ad_Scof(i,k)=0.0_r8 ad_Tcof(i,k)=0.0_r8 ad_wrk(i,k)=0.0_r8 ad_bulk(i,k)=0.0_r8 ad_bulk0(i,k)=0.0_r8 ad_bulk1(i,k)=0.0_r8 ad_bulk2(i,k)=0.0_r8 ad_den(i,k)=0.0_r8 ad_den1(i,k)=0.0_r8 END DO END DO ! !======================================================================= ! Adjoint nonlinear equation of state. Notice that this equation ! of state is only valid for potential temperature range of -2C to 40C ! and a salinity range of 0 PSU to 42 PSU. !======================================================================= ! !----------------------------------------------------------------------- ! Exchange boundary data. !----------------------------------------------------------------------- ! !^ CALL mp_exchange3d (ng, tile, model, 1, & !^ & LBi, UBi, LBj, UBj, 1, N(ng), & !^ & NghostPoints, & !^ & EWperiodic(ng), NSperiodic(ng), & !^ & tl_rho) !^ CALL ad_mp_exchange3d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, 1, N(ng), & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_rho) ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN !^ CALL exchange_r3d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, 1, N(ng), & !^ & tl_rho) !^ CALL ad_exchange_r3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ad_rho) END IF ! !----------------------------------------------------------------------- ! Compute BASIC STATE related variables. !----------------------------------------------------------------------- ! DO j=JstrT,JendT DO k=1,N(ng) DO i=IstrT,IendT Tt=MAX(-2.0_r8,t(i,j,k,nrhs,itemp)) Ts=MAX(0.0_r8,t(i,j,k,nrhs,isalt)) sqrtTs=SQRT(Ts) Tp=z_r(i,j,k) Tpr10=0.1_r8*Tp ! ! Compute local nonlinear equation of state coefficients and their ! derivatives when appropriate. ! C(0)=Q00+Tt*(Q01+Tt*(Q02+Tt*(Q03+Tt*(Q04+Tt*Q05)))) C(1)=U00+Tt*(U01+Tt*(U02+Tt*(U03+Tt*U04))) C(2)=V00+Tt*(V01+Tt*V02) C(3)=A00+Tt*(A01+Tt*(A02+Tt*(A03+Tt*A04))) C(4)=B00+Tt*(B01+Tt*(B02+Tt*B03)) C(5)=D00+Tt*(D01+Tt*D02) C(6)=E00+Tt*(E01+Tt*(E02+Tt*E03)) C(7)=F00+Tt*(F01+Tt*F02) C(8)=G01+Tt*(G02+Tt*G03) C(9)=H00+Tt*(H01+Tt*H02) ! dCdT(0)=Q01+Tt*(2.0_r8*Q02+Tt*(3.0_r8*Q03+Tt*(4.0_r8*Q04+ & & Tt*5.0_r8*Q05))) dCdT(1)=U01+Tt*(2.0_r8*U02+Tt*(3.0_r8*U03+Tt*4.0_r8*U04)) dCdT(2)=V01+Tt*2.0_r8*V02 dCdT(3)=A01+Tt*(2.0_r8*A02+Tt*(3.0_r8*A03+Tt*4.0_r8*A04)) dCdT(4)=B01+Tt*(2.0_r8*B02+Tt*3.0_r8*B03) dCdT(5)=D01+Tt*2.0_r8*D02 dCdT(6)=E01+Tt*(2.0_r8*E02+Tt*3.0_r8*E03) dCdT(7)=F01+Tt*2.0_r8*F02 dCdT(8)=G02+Tt*2.0_r8*G03 dCdT(9)=H01+Tt*2.0_r8*H02 ! d2Cd2T(0)=2.0_r8*Q02+Tt*(6.0_r8*Q03+Tt*(12.0_r8*Q04+ & & Tt*20.0_r8*Q05)) d2Cd2T(1)=2.0_r8*U02+Tt*(6.0_r8*U03+Tt*12.0_r8*U04) d2Cd2T(2)=2.0_r8*V02 d2Cd2T(3)=2.0_r8*A02+Tt*(6.0_r8*A03+Tt*12.0_r8*A04) d2Cd2T(4)=2.0_r8*B02+Tt*6.0_r8*B03 d2Cd2T(5)=2.0_r8*D02 d2Cd2T(6)=2.0_r8*E02+Tt*6.0_r8*E03 d2Cd2T(7)=2.0_r8*F02 d2Cd2T(8)=2.0_r8*G03 d2Cd2T(9)=2.0_r8*H02 ! ! Compute BASIC STATE density (kg/m3) at standard one atmosphere ! pressure. ! den1(i,k)=C(0)+Ts*(C(1)+sqrtTs*C(2)+Ts*W00) ! ! Compute BASIC STATE d(den1)/d(S) and d(den1)/d(T) derivatives used ! in the computation of thermal expansion and saline contraction ! coefficients. ! Dden1DS(i,k)=C(1)+1.5_r8*C(2)*sqrtTs+2.0_r8*W00*Ts Dden1DT(i,k)=dCdT(0)+Ts*(dCdT(1)+sqrtTs*dCdT(2)) ! ! Compute BASIC STATE secant bulk modulus. ! bulk0(i,k)=C(3)+Ts*(C(4)+sqrtTs*C(5)) bulk1(i,k)=C(6)+Ts*(C(7)+sqrtTs*G00) bulk2(i,k)=C(8)+Ts*C(9) bulk (i,k)=bulk0(i,k)-Tp*(bulk1(i,k)-Tp*bulk2(i,k)) ! ! Compute local "in situ" density anomaly (kg/m3 - 1000). The (i,k) ! DO-loop is closed here because of the adjoint to facilitate vertical ! integrals of the BASIC STATE. ! cff=1.0_r8/(bulk(i,k)+Tpr10) den(i,k)=den1(i,k)*bulk(i,k)*cff den(i,k)=den(i,k)-1000.0_r8 den(i,k)=den(i,k)*rmask(i,j) END DO END DO ! !----------------------------------------------------------------------- ! Load adjoint "in situ" density anomaly (kg/m3 - 1000) and adjoint ! potential density anomaly (kg/m3 - 1000) referenced to the surface ! into global arrays. !----------------------------------------------------------------------- ! DO k=1,N(ng) DO i=IstrT,IendT !^ tl_rho(i,j,k)=tl_den(i,k) !^ ad_den(i,k)=ad_den(i,k)+ad_rho(i,j,k) ad_rho(i,j,k)=0.0_r8 END DO END DO ! !----------------------------------------------------------------------- ! Adjoint nonlinear equation of state. !----------------------------------------------------------------------- ! DO k=1,N(ng) DO i=IstrT,IendT ! ! Check temperature and salinity minimum valid values. Assign depth ! to the pressure. ! Tt=MAX(-2.0_r8,t(i,j,k,nrhs,itemp)) Ts=MAX(0.0_r8,t(i,j,k,nrhs,isalt)) sqrtTs=SQRT(Ts) Tp=z_r(i,j,k) Tpr10=0.1_r8*Tp ! ! Compute local nonlinear equation of state coefficients and their ! derivatives when appropriate. These coefficients can be stored ! in slab (i,k) arrays to avoid recompute them twice. However, the ! equivalent of 50 slabs arrays are required. ! C(0)=Q00+Tt*(Q01+Tt*(Q02+Tt*(Q03+Tt*(Q04+Tt*Q05)))) C(1)=U00+Tt*(U01+Tt*(U02+Tt*(U03+Tt*U04))) C(2)=V00+Tt*(V01+Tt*V02) C(3)=A00+Tt*(A01+Tt*(A02+Tt*(A03+Tt*A04))) C(4)=B00+Tt*(B01+Tt*(B02+Tt*B03)) C(5)=D00+Tt*(D01+Tt*D02) C(6)=E00+Tt*(E01+Tt*(E02+Tt*E03)) C(7)=F00+Tt*(F01+Tt*F02) C(8)=G01+Tt*(G02+Tt*G03) C(9)=H00+Tt*(H01+Tt*H02) ! dCdT(0)=Q01+Tt*(2.0_r8*Q02+Tt*(3.0_r8*Q03+Tt*(4.0_r8*Q04+ & & Tt*5.0_r8*Q05))) dCdT(1)=U01+Tt*(2.0_r8*U02+Tt*(3.0_r8*U03+Tt*4.0_r8*U04)) dCdT(2)=V01+Tt*2.0_r8*V02 dCdT(3)=A01+Tt*(2.0_r8*A02+Tt*(3.0_r8*A03+Tt*4.0_r8*A04)) dCdT(4)=B01+Tt*(2.0_r8*B02+Tt*3.0_r8*B03) dCdT(5)=D01+Tt*2.0_r8*D02 dCdT(6)=E01+Tt*(2.0_r8*E02+Tt*3.0_r8*E03) dCdT(7)=F01+Tt*2.0_r8*F02 dCdT(8)=G02+Tt*2.0_r8*G03 dCdT(9)=H01+Tt*2.0_r8*H02 ! d2Cd2T(0)=2.0_r8*Q02+Tt*(6.0_r8*Q03+Tt*(12.0_r8*Q04+ & & Tt*20.0_r8*Q05)) d2Cd2T(1)=2.0_r8*U02+Tt*(6.0_r8*U03+Tt*12.0_r8*U04) d2Cd2T(2)=2.0_r8*V02 d2Cd2T(3)=2.0_r8*A02+Tt*(6.0_r8*A03+Tt*12.0_r8*A04) d2Cd2T(4)=2.0_r8*B02+Tt*6.0_r8*B03 d2Cd2T(5)=2.0_r8*D02 d2Cd2T(6)=2.0_r8*E02+Tt*6.0_r8*E03 d2Cd2T(7)=2.0_r8*F02 d2Cd2T(8)=2.0_r8*G03 d2Cd2T(9)=2.0_r8*H02 ! !----------------------------------------------------------------------- ! Compute local adjoint "in situ" density anomaly (kg/m3 - 1000). !----------------------------------------------------------------------- ! cff=1.0_r8/(bulk(i,k)+Tpr10) !^ tl_den(i,k)=tl_den(i,k)*rmask(i,j) !^ ad_den(i,k)=ad_den(i,k)*rmask(i,j) !^ tl_den(i,k)=tl_den1(i,k)*bulk(i,k)*cff+ & !^ & den1(i,k)*(tl_bulk(i,k)*cff+ & !^ & bulk(i,k)*tl_cff) !^ adfac1=den1(i,k)*ad_den(i,k) ad_den1(i,k)=ad_den1(i,k)+bulk(i,k)*cff*ad_den(i,k) ad_bulk(i,k)=ad_bulk(i,k)+cff*adfac1 ad_cff=ad_cff+bulk(i,k)*adfac1 ad_den(i,k)=0.0_r8 !^ tl_cff=-cff*cff*(tl_bulk(i,k)+tl_Tpr10) !^ adfac=-cff*cff*ad_cff ad_bulk(i,k)=ad_bulk(i,k)+adfac ad_Tpr10=ad_Tpr10+adfac ad_cff=0.0_r8 ! ! Compute adjoint secant bulk modulus. ! !^ tl_bulk (i,k)=tl_bulk0(i,k)- & !^ & tl_Tp*(bulk1(i,k)-Tp*bulk2(i,k))- & !^ & Tp*(tl_bulk1(i,k)- & !^ & tl_Tp*bulk2(i,k)- & !^ & Tp*tl_bulk2(i,k)) !^ adfac=Tp*ad_bulk(i,k) ad_bulk0(i,k)=ad_bulk0(i,k)+ad_bulk(i,k) ad_bulk1(i,k)=ad_bulk1(i,k)-adfac ad_bulk2(i,k)=ad_bulk2(i,k)+adfac*Tp ad_Tp=ad_Tp- & & ad_bulk(i,k)*(bulk1(i,k)-Tp*bulk2(i,k))+ & & adfac*bulk2(i,k) ad_bulk(i,k)=0.0_r8 !^ tl_bulk2(i,k)=tl_C(8)+tl_Ts*C(9)+Ts*tl_C(9) !^ ad_C(8)=ad_C(8)+ad_bulk2(i,k) ad_C(9)=ad_C(9)+Ts*ad_bulk2(i,k) ad_Ts=ad_Ts+ad_bulk2(i,k)*C(9) ad_bulk2(i,k)=0.0_r8 !^ tl_bulk1(i,k)=tl_C(6)+ & !^ & tl_Ts*(C(7)+sqrtTs*G00)+ & !^ & Ts*(tl_C(7)+tl_sqrtTs*G00) !^ adfac=Ts*ad_bulk1(i,k) ad_C(6)=ad_C(6)+ad_bulk1(i,k) ad_C(7)=ad_C(7)+adfac ad_Ts=ad_Ts+ad_bulk1(i,k)*(C(7)+sqrtTs*G00) ad_sqrtTs=ad_sqrtTs+adfac*G00 ad_bulk1(i,k)=0.0_r8 !^ tl_bulk0(i,k)=tl_C(3)+ & !^ & tl_Ts*(C(4)+sqrtTs*C(5))+ & !^ & Ts*(tl_C(4)+tl_sqrtTs*C(5)+ & !^ & sqrtTs*tl_C(5)) !^ adfac=Ts*ad_bulk0(i,k) ad_C(3)=ad_C(3)+ad_bulk0(i,k) ad_C(4)=ad_C(4)+adfac ad_C(5)=ad_C(5)+sqrtTs*adfac ad_Ts=ad_Ts+ad_bulk0(i,k)*(C(4)+sqrtTs*C(5)) ad_sqrtTs=ad_sqrtTs+C(5)*adfac ad_bulk0(i,k)=0.0_r8 !^ tl_C(9)=tl_Tt*dCdT(9) !^ tl_C(8)=tl_Tt*dCdT(8) !^ tl_C(7)=tl_Tt*dCdT(7) !^ tl_C(6)=tl_Tt*dCdT(6) !^ tl_C(5)=tl_Tt*dCdT(5) !^ tl_C(4)=tl_Tt*dCdT(4) !^ tl_C(3)=tl_Tt*dCdT(3) !^ ad_Tt=ad_Tt+ad_C(9)*dCdT(9)+ & & ad_C(8)*dCdT(8)+ & & ad_C(7)*dCdT(7)+ & & ad_C(6)*dCdT(6)+ & & ad_C(5)*dCdT(5)+ & & ad_C(4)*dCdT(4)+ & & ad_C(3)*dCdT(3) ad_C(9)=0.0_r8 ad_C(8)=0.0_r8 ad_C(7)=0.0_r8 ad_C(6)=0.0_r8 ad_C(5)=0.0_r8 ad_C(4)=0.0_r8 ad_C(3)=0.0_r8 ! ! Compute d(den1)/d(S) and d(den1)/d(T) derivatives used in the ! computation of thermal expansion and saline contraction ! coefficients. ! !^ tl_Dden1DT(i,k)=tl_dCdT(0)+ & !^ & tl_Ts*(dCdT(1)+sqrtTs*dCdT(2))+ & !^ & Ts*(tl_dCdT(1)+tl_sqrtTs*dCdT(2)+ & !^ & sqrtTs*tl_dCdT(2)) !^ adfac1=Ts*ad_Dden1DT(i,k) ad_dCdT(0)=ad_dCdT(0)+ad_Dden1DT(i,k) ad_dCdT(1)=ad_dCdT(1)+adfac1 ad_dCdT(2)=ad_dCdT(2)+sqrtTs*adfac1 ad_Ts=ad_Ts+ & & (dCdT(1)+sqrtTs*dCdT(2))*ad_Dden1DT(i,k) ad_sqrtTs=ad_sqrtTs+dCdT(2)*adfac1 ad_Dden1DT(i,k)=0.0_r8 !^ tl_Dden1DS(i,k)=tl_C(1)+ & !^ & 1.5_r8*(tl_C(2)*sqrtTs+ & !^ & C(2)*tl_sqrtTs)+ & !^ & 2.0_r8*W00*tl_Ts !^ adfac1=1.5_r8*ad_Dden1DS(i,k) ad_C(1)=ad_C(1)+ad_Dden1DS(i,k) ad_C(2)=ad_C(2)+sqrtTs*adfac1 ad_Ts=ad_Ts+2.0_r8*W00*ad_Dden1DS(i,k) ad_sqrtTs=ad_sqrtTs+C(2)*adfac1 ad_Dden1DS(i,k)=0.0_r8 !^ tl_dCdT(2)=tl_Tt*d2Cd2T(2) !^ tl_dCdT(1)=tl_Tt*d2Cd2T(1) !^ tl_dCdT(0)=tl_Tt*d2Cd2T(0) !^ ad_Tt=ad_Tt+d2Cd2T(2)*ad_dCdT(2)+ & & d2Cd2T(1)*ad_dCdT(1)+ & & d2Cd2T(0)*ad_dCdT(0) ad_dCdT(2)=0.0_r8 ad_dCdT(1)=0.0_r8 ad_dCdT(0)=0.0_r8 ! ! Compute basic state and tangent linear density (kg/m3) at standard ! one atmosphere pressure. ! !^ tl_den1(i,k)=tl_C(0)+ & !^ & tl_Ts*(C(1)+sqrtTs*C(2)+Ts*W00)+ & !^ & Ts*(tl_C(1)+tl_sqrtTs*C(2)+ & !^ & sqrtTs*tl_C(2)+tl_Ts*W00) !^ adfac=Ts*ad_den1(i,k) ad_C(0)=ad_C(0)+ad_den1(i,k) ad_C(1)=ad_C(1)+adfac ad_C(2)=ad_C(2)+adfac*sqrtTs ad_Ts=ad_Ts+ & & ad_den1(i,k)*(C(1)+sqrtTs*C(2)+Ts*W00)+ & & adfac*W00 ad_sqrtTs=ad_sqrtTs+adfac*C(2) ad_den1(i,k)=0.0_r8 !^ tl_C(2)=tl_Tt*dCdT(2) !^ tl_C(1)=tl_Tt*dCdT(1) !^ tl_C(0)=tl_Tt*dCdT(0) !^ ad_Tt=ad_Tt+ad_C(2)*dCdT(2)+ & & ad_C(1)*dCdT(1)+ & & ad_C(0)*dCdT(0) ad_C(2)=0.0_r8 ad_C(1)=0.0_r8 ad_C(0)=0.0_r8 ! ! Check temperature and salinity minimum valid values. Assign depth ! to the pressure. ! !^ tl_Tpr10=0.1_r8*tl_Tp !^ ad_Tp=ad_Tp+0.1_r8*ad_Tpr10 ad_Tpr10=0.0_r8 !^ tl_Tp=tl_z_r(i,j,k) !^ ad_z_r(i,j,k)=ad_z_r(i,j,k)+ad_Tp ad_Tp=0.0_r8 IF (Ts.ne.0.0_r8) THEN !^ tl_sqrtTs=0.5_r8*tl_Ts/SQRT(Ts) !^ ad_Ts=ad_Ts+0.5_r8*ad_sqrtTs/SQRT(Ts) ad_sqrtTs=0.0_r8 ELSE !^ tl_sqrtTs=0.0_r8 !^ ad_sqrtTs=0.0_r8 END IF !^ tl_Ts=(0.5_r8-SIGN(0.5_r8,-t(i,j,k,nrhs,isalt)))* !^ & tl_t(i,j,k,nrhs,isalt) !^ ad_t(i,j,k,nrhs,isalt)=ad_t(i,j,k,nrhs,isalt)+ & & (0.5_r8-SIGN(0.5_r8, & & -t(i,j,k,nrhs,isalt)))* & & ad_Ts ad_Ts=0.0_r8 !^ tl_Tt=(0.5_r8-SIGN(0.5_r8,-2.0_r8-t(i,j,k,nrhs,itemp)))* !^ & tl_t(i,j,k,nrhs,itemp) !^ ad_t(i,j,k,nrhs,itemp)=ad_t(i,j,k,nrhs,itemp)+ & & (0.5_r8-SIGN(0.5_r8,-2.0_r8- & & t(i,j,k,nrhs,itemp)))* & & ad_Tt ad_Tt=0.0_r8 END DO END DO END DO ! RETURN END SUBROUTINE ad_rho_eos_tile END MODULE ad_rho_eos_mod