MODULE tl_variability_mod ! !git $Id$ !svn $Id: tl_variability.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 converts tangent error correlations, C, to background ! ! error covariances, B, by multiplying by the standard deviations, ! ! S. ! ! ! ! ! ! The background error covariance is defined as: ! ! ! ! B = S C S ! ! ! ! C = C^(1/2) C^(T/2) ! ! ! ! C^(1/2) = G L^(1/2) W^(-1/2) TLM ! ! C^(T/2) = W^(-1/2) L^(T/2) G ADM ! ! ! ! where ! ! ! ! B : background-error covariance matrix ! ! S : diagonal matrix of background-error standard deviations ! ! C : symmetric matrix of background-error correlations ! ! G : normalization coefficients matrix (convert to correlations) ! ! L : tangent linear and adjoint diffusion operators ! ! W : diagonal matrix of local area or volume metrics ! ! ! ! Here, T/2 denote the transpose of a squared-root factor. ! ! ! ! Reference: ! ! ! ! Weaver, A. and P. Courtier, 2001: Correlation modeling on the ! ! sphere using a generalized diffusion equation, Q.J.R. Meteo. ! ! Soc, 127, 1815-1846. ! ! ! !======================================================================! ! USE mod_kinds implicit none PRIVATE PUBLIC :: tl_variability CONTAINS ! !*********************************************************************** SUBROUTINE tl_variability (ng, tile, Linp, Lweak) !*********************************************************************** ! USE mod_param USE mod_ocean ! ! Imported variable declarations. ! logical, intent(in) :: Lweak integer, intent(in) :: ng, tile, Linp ! ! Local variable declarations. ! 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 tl_variability_tile (ng, tile, & & LBi, UBi, LBj, UBj, LBij, UBij, & & IminS, ImaxS, JminS, JmaxS, & & Linp, Lweak, & & OCEAN(ng) % e_t, & & OCEAN(ng) % e_u, & & OCEAN(ng) % e_v, & & OCEAN(ng) % e_zeta, & & OCEAN(ng) % tl_t, & & OCEAN(ng) % tl_u, & & OCEAN(ng) % tl_v, & & OCEAN(ng) % tl_zeta) RETURN END SUBROUTINE tl_variability ! !*********************************************************************** SUBROUTINE tl_variability_tile (ng, tile, & & LBi, UBi, LBj, UBj, LBij, UBij, & & IminS, ImaxS, JminS, JmaxS, & & Linp, Lweak, & & t_std, u_std, v_std, & & zeta_std, & & tl_t, tl_u, tl_v, & & tl_zeta) !*********************************************************************** ! USE mod_param USE mod_ncparam USE mod_scalars ! USE mp_exchange_mod ! ! Imported variable declarations. ! logical, intent(in) :: Lweak integer, intent(in) :: ng, tile integer, intent(in) :: LBi, UBi, LBj, UBj, LBij, UBij integer, intent(in) :: IminS, ImaxS, JminS, JmaxS integer, intent(in) :: Linp ! real(r8), intent(in) :: t_std(LBi:,LBj:,:,:,:) real(r8), intent(in) :: u_std(LBi:,LBj:,:,:) real(r8), intent(in) :: v_std(LBi:,LBj:,:,:) real(r8), intent(in) :: zeta_std(LBi:,LBj:,:) real(r8), intent(inout) :: tl_t(LBi:,LBj:,:,:,:) real(r8), intent(inout) :: tl_u(LBi:,LBj:,:,:) real(r8), intent(inout) :: tl_v(LBi:,LBj:,:,:) real(r8), intent(inout) :: tl_zeta(LBi:,LBj:,:) ! ! Local variable declarations. ! integer :: i, ib, ir, j, rec integer :: it, k ! !----------------------------------------------------------------------- ! 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 ! !----------------------------------------------------------------------- ! Determine error covariance standard deviation factors to use. !----------------------------------------------------------------------- ! IF (Lweak) THEN rec=2 ! weak constraint ELSE rec=1 ! strong constraint END IF ! !----------------------------------------------------------------------- ! Initial conditions and model error covariance: Multiply tangent ! linear state by its corresponding standard deviation. !----------------------------------------------------------------------- ! ! Free-surface. ! DO j=JstrT,JendT DO i=IstrT,IendT tl_zeta(i,j,Linp)=tl_zeta(i,j,Linp)*zeta_std(i,j,rec) END DO END DO CALL mp_exchange2d (ng, tile, iTLM, 1, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & tl_zeta(:,:,Linp)) ! ! 3D U-momentum. ! DO k=1,N(ng) DO j=JstrT,JendT DO i=IstrP,IendT tl_u(i,j,k,Linp)=tl_u(i,j,k,Linp)*u_std(i,j,k,rec) END DO END DO DO j=JstrP,JendT DO i=IstrT,IendT tl_v(i,j,k,Linp)=tl_v(i,j,k,Linp)*v_std(i,j,k,rec) END DO END DO END DO CALL mp_exchange3d (ng, tile, iTLM, 2, & & LBi, UBi, LBj, UBj, 1, N(ng), & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & tl_u(:,:,:,Linp), & & tl_v(:,:,:,Linp)) ! ! Tangent linear tracers. ! DO it=1,NT(ng) DO k=1,N(ng) DO j=JstrT,JendT DO i=IstrT,IendT tl_t(i,j,k,Linp,it)=tl_t(i,j,k,Linp,it)* & & t_std(i,j,k,rec,it) END DO END DO END DO END DO CALL mp_exchange4d (ng, tile, iTLM, 1, & & LBi, UBi, LBj, UBj, 1, N(ng), 1, NT(ng), & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & tl_t(:,:,:,Linp,:)) RETURN END SUBROUTINE tl_variability_tile END MODULE tl_variability_mod