MODULE tl_rho_eos_mod ! !git $Id$ !svn $Id: tl_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 :: tl_rho_eos ! CONTAINS ! !*********************************************************************** SUBROUTINE tl_rho_eos (ng, tile, model) !*********************************************************************** ! USE mod_param 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/Tangent/tl_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, 70, MyFile) CALL tl_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) % tl_z_r, & & GRID(ng) % z_w, & & GRID(ng) % tl_z_w, & & OCEAN(ng) % t, & & OCEAN(ng) % tl_t, & & MIXING(ng) % alpha, & & MIXING(ng) % tl_alpha, & & MIXING(ng) % beta, & & MIXING(ng) % tl_beta, & & OCEAN(ng) % rho, & & OCEAN(ng) % tl_rho) CALL wclock_off (ng, model, 14, 114, MyFile) ! RETURN END SUBROUTINE tl_rho_eos ! !*********************************************************************** SUBROUTINE tl_rho_eos_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nrhs, & & rmask, & & z_r, tl_z_r, & & z_w, tl_z_w, & & t, tl_t, & & alpha, tl_alpha, & & beta, tl_beta, & & rho, tl_rho) !*********************************************************************** ! USE mod_param USE mod_eoscoef USE mod_scalars ! USE exchange_2d_mod USE exchange_3d_mod USE mp_exchange_mod, ONLY : mp_exchange2d, 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) :: z_w(LBi:,LBj:,0:) real(r8), intent(in) :: t(LBi:,LBj:,:,:,:) real(r8), intent(in) :: tl_z_r(LBi:,LBj:,:) real(r8), intent(in) :: tl_z_w(LBi:,LBj:,0:) real(r8), intent(in) :: tl_t(LBi:,LBj:,:,:,:) real(r8), intent(inout) :: alpha(LBi:,LBj:) real(r8), intent(inout) :: beta(LBi:,LBj:) real(r8), intent(out) :: tl_alpha(LBi:,LBj:) real(r8), intent(out) :: tl_beta(LBi:,LBj:) real(r8), intent(out) :: rho(LBi:,LBj:,:) real(r8), intent(out) :: tl_rho(LBi:,LBj:,:) ! ! Local variable declarations. ! integer :: i, ised, itrc, j, k real(r8) :: SedDen, Tp, Tpr10, Ts, Tt, sqrtTs real(r8) :: tl_SedDen, tl_Tp, tl_Tpr10, tl_Ts, tl_Tt, tl_sqrtTs real(r8) :: cff, cff1, cff2, cff3 real(r8) :: tl_cff, tl_cff1, tl_cff2, tl_cff3 real(r8), dimension(0:9) :: C real(r8), dimension(0:9) :: tl_C real(r8), dimension(0:9) :: dCdT(0:9) real(r8), dimension(0:9) :: tl_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)) :: tl_DbulkDS real(r8), dimension(IminS:ImaxS,N(ng)) :: tl_DbulkDT real(r8), dimension(IminS:ImaxS,N(ng)) :: tl_Dden1DS real(r8), dimension(IminS:ImaxS,N(ng)) :: tl_Dden1DT real(r8), dimension(IminS:ImaxS,N(ng)) :: tl_Scof real(r8), dimension(IminS:ImaxS,N(ng)) :: tl_Tcof real(r8), dimension(IminS:ImaxS,N(ng)) :: tl_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)) :: tl_bulk real(r8), dimension(IminS:ImaxS,N(ng)) :: tl_bulk0 real(r8), dimension(IminS:ImaxS,N(ng)) :: tl_bulk1 real(r8), dimension(IminS:ImaxS,N(ng)) :: tl_bulk2 real(r8), dimension(IminS:ImaxS,N(ng)) :: tl_den real(r8), dimension(IminS:ImaxS,N(ng)) :: tl_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 ! !======================================================================= ! 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. !======================================================================= ! DO j=JstrT,JendT DO k=1,N(ng) DO i=IstrT,IendT ! ! Check temperature and salinity lower values. Assign depth to the ! pressure. ! Tt=MAX(-2.0_r8,t(i,j,k,nrhs,itemp)) 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) Ts=MAX(0.0_r8,t(i,j,k,nrhs,isalt)) tl_Ts=(0.5_r8-SIGN(0.5_r8,-t(i,j,k,nrhs,isalt)))* & & tl_t(i,j,k,nrhs,isalt) sqrtTs=SQRT(Ts) IF (Ts.ne.0.0_r8) THEN tl_sqrtTs=0.5_r8*tl_Ts/SQRT(Ts) ELSE tl_sqrtTs=0.0_r8 END IF Tp=z_r(i,j,k) tl_Tp=tl_z_r(i,j,k) Tpr10=0.1_r8*Tp tl_Tpr10=0.1_r8*tl_Tp ! !----------------------------------------------------------------------- ! Compute BASIC STATE and tangent linear density (kg/m3) at standard ! one atmosphere pressure. !----------------------------------------------------------------------- ! 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) ! 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 tl_C(0)=tl_Tt*dCdT(0) tl_C(1)=tl_Tt*dCdT(1) tl_C(2)=tl_Tt*dCdT(2) ! den1(i,k)=C(0)+Ts*(C(1)+sqrtTs*C(2)+Ts*W00) 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) ! ! Compute d(den1)/d(S) and d(den1)/d(T) derivatives used in the ! computation of thermal expansion and saline contraction ! coefficients. ! 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 ! tl_dCdT(0)=tl_Tt*d2Cd2T(0) tl_dCdT(1)=tl_Tt*d2Cd2T(1) tl_dCdT(2)=tl_Tt*d2Cd2T(2) ! 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)) ! tl_Dden1DS(i,k)=tl_C(1)+ & & 1.5_r8*(tl_C(2)*sqrtTs+ & & C(2)*tl_sqrtTs)+ & & 2.0_r8*W00*tl_Ts 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)) ! !----------------------------------------------------------------------- ! Compute BASIC STATE and tangent linear secant bulk modulus. !----------------------------------------------------------------------- ! 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(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 ! tl_C(3)=tl_Tt*dCdT(3) tl_C(4)=tl_Tt*dCdT(4) tl_C(5)=tl_Tt*dCdT(5) tl_C(6)=tl_Tt*dCdT(6) tl_C(7)=tl_Tt*dCdT(7) tl_C(8)=tl_Tt*dCdT(8) tl_C(9)=tl_Tt*dCdT(9) ! 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)) ! 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)) tl_bulk1(i,k)=tl_C(6)+ & & tl_Ts*(C(7)+sqrtTs*G00)+ & & Ts*(tl_C(7)+tl_sqrtTs*G00) tl_bulk2(i,k)=tl_C(8)+tl_Ts*C(9)+Ts*tl_C(9) 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)) ! ! Compute d(bulk)/d(S) and d(bulk)/d(T) derivatives used ! in the computation of thermal expansion and saline contraction ! coefficients. ! 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 ! tl_dCdT(3)=tl_Tt*d2Cd2T(3) tl_dCdT(4)=tl_Tt*d2Cd2T(4) tl_dCdT(5)=tl_Tt*d2Cd2T(5) tl_dCdT(6)=tl_Tt*d2Cd2T(6) tl_dCdT(7)=tl_Tt*d2Cd2T(7) tl_dCdT(8)=tl_Tt*d2Cd2T(8) tl_dCdT(9)=tl_Tt*d2Cd2T(9) ! DbulkDS(i,k)=C(4)+sqrtTs*1.5_r8*C(5)- & & Tp*(C(7)+sqrtTs*1.5_r8*G00-Tp*C(9)) DbulkDT(i,k)=dCdT(3)+Ts*(dCdT(4)+sqrtTs*dCdT(5))- & & Tp*(dCdT(6)+Ts*dCdT(7)- & & Tp*(dCdT(8)+Ts*dCdT(9))) ! tl_DbulkDS(i,k)=tl_C(4)+ & & 1.5_r8*(tl_sqrtTs*C(5)+ & & sqrtTs*tl_C(5))- & & tl_Tp*(C(7)+sqrtTs*1.5_r8*G00- & & Tp*C(9))- & & Tp*(tl_C(7)+tl_sqrtTs*1.5_r8*G00- & & tl_Tp*C(9)-Tp*tl_C(9)) tl_DbulkDT(i,k)=tl_dCdT(3)+ & & tl_Ts*(dCdT(4)+sqrtTs*dCdT(5))+ & & Ts*(tl_dCdT(4)+tl_sqrtTs*dCdT(5)+ & & sqrtTs*tl_dCdT(5))- & & tl_Tp*(dCdT(6)+Ts*dCdT(7)- & & Tp*(dCdT(8)+Ts*dCdT(9)))- & & Tp*(tl_dCdT(6)+tl_Ts*dCdT(7)+Ts*tl_dCdT(7)- & & tl_Tp*(dCdT(8)+Ts*dCdT(9))- & & Tp*(tl_dCdT(8)+tl_Ts*dCdT(9)+ & & Ts*tl_dCdT(9))) ! !----------------------------------------------------------------------- ! Compute local "in situ" density anomaly (kg/m3 - 1000). !----------------------------------------------------------------------- ! cff=1.0_r8/(bulk(i,k)+Tpr10) tl_cff=-cff*cff*(tl_bulk(i,k)+tl_Tpr10) den(i,k)=den1(i,k)*bulk(i,k)*cff 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) den(i,k)=den(i,k)-1000.0_r8 den(i,k)=den(i,k)*rmask(i,j) tl_den(i,k)=tl_den(i,k)*rmask(i,j) END DO END DO ! !----------------------------------------------------------------------- ! Compute thermal expansion (1/Celsius) and saline contraction ! (1/PSU) coefficients. !----------------------------------------------------------------------- ! DO k=N(ng),N(ng) DO i=IstrT,IendT Tpr10=0.1_r8*z_r(i,j,k) tl_Tpr10=0.1_r8*tl_z_r(i,j,k) ! ! Compute thermal expansion and saline contraction coefficients. ! cff=bulk(i,k)+Tpr10 tl_cff=tl_bulk(i,k)+tl_Tpr10 cff1=Tpr10*den1(i,k) tl_cff1=tl_Tpr10*den1(i,k)+Tpr10*tl_den1(i,k) cff2=bulk(i,k)*cff tl_cff2=tl_bulk(i,k)*cff+bulk(i,k)*tl_cff wrk(i,k)=(den(i,k)+1000.0_r8)*cff*cff tl_wrk(i,k)=cff*(cff*tl_den(i,k)+ & & 2.0_r8*tl_cff*(den(i,k)+1000.0_r8)) Tcof(i,k)=-(DbulkDT(i,k)*cff1+ & & Dden1DT(i,k)*cff2) tl_Tcof(i,k)=-(tl_DbulkDT(i,k)*cff1+ & & DbulkDT(i,k)*tl_cff1+ & & tl_Dden1DT(i,k)*cff2+ & & Dden1DT(i,k)*tl_cff2) Scof(i,k)= (DbulkDS(i,k)*cff1+ & & Dden1DS(i,k)*cff2) tl_Scof(i,k)= (tl_DbulkDS(i,k)*cff1+ & & DbulkDS(i,k)*tl_cff1+ & & tl_Dden1DS(i,k)*cff2+ & & Dden1DS(i,k)*tl_cff2) END DO IF (k.eq.N(ng)) THEN DO i=IstrT,IendT cff=1.0_r8/wrk(i,N(ng)) tl_cff=-cff*cff*tl_wrk(i,N(ng)) alpha(i,j)=cff*Tcof(i,N(ng)) tl_alpha(i,j)=tl_cff*Tcof(i,N(ng))+cff*tl_Tcof(i,N(ng)) beta (i,j)=cff*Scof(i,N(ng)) tl_beta (i,j)=tl_cff*Scof(i,N(ng))+cff*tl_Scof(i,N(ng)) END DO END IF END DO ! !----------------------------------------------------------------------- ! Load "in situ" density anomaly (kg/m3 - 1000) and potential ! density anomaly (kg/m3 - 1000) referenced to the surface into global ! arrays. Notice that this is done in a separate (i,k) DO-loops to ! facilitate the adjoint. !----------------------------------------------------------------------- ! DO k=1,N(ng) DO i=IstrT,IendT rho(i,j,k)=den(i,k) tl_rho(i,j,k)=tl_den(i,k) END DO END DO END DO ! !----------------------------------------------------------------------- ! Exchange boundary data. !----------------------------------------------------------------------- ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN CALL exchange_r3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & rho) CALL exchange_r3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & tl_rho) END IF ! CALL mp_exchange3d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, 1, N(ng), & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & rho) CALL mp_exchange3d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, 1, N(ng), & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & tl_rho) CALL mp_exchange2d (ng, tile, model, 2, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & alpha, beta) CALL mp_exchange2d (ng, tile, model, 2, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & tl_alpha, tl_beta) ! RETURN END SUBROUTINE tl_rho_eos_tile END MODULE tl_rho_eos_mod