MODULE conv_3d_mod ! !git $Id$ !svn $Id: conv_3d.F 1151 2023-02-09 03:08:53Z 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 ! !======================================================================= ! ! ! These routines applies the background error covariance to data ! ! assimilation fields via the space convolution of the diffusion ! ! equation (filter) for 3D state variables. The diffusion filter ! ! is solved using an implicit or explicit vertical algorithm. ! ! ! ! For Gaussian (bell-shaped) correlations, the space convolution ! ! of the diffusion operator is an efficient way to estimate the ! ! finite domain error covariances. ! ! ! ! Notice that "z_r" and "Hz" are assumed to be time invariant in ! ! the spatial convolution. That it, they are not adjointable. ! ! ! ! On Input: ! ! ! ! ng Nested grid number. ! ! model Calling model identifier. ! ! Istr Starting tile index in the I-direction. ! ! Iend Ending tile index in the I-direction. ! ! Jstr Starting tile index in the J-direction. ! ! Jend Ending tile index in the J-direction. ! ! LBi I-dimension lower bound. ! ! UBi I-dimension upper bound. ! ! LBj J-dimension lower bound. ! ! UBj J-dimension upper bound. ! ! LBk K-dimension lower bound. ! ! UBk K-dimension upper bound. ! ! Nghost Number of ghost points. ! ! NHsteps Number of horizontal diffusion integration steps. ! ! NVsteps Number of vertical diffusion integration steps. ! ! DTsizeH Horizontal diffusion pseudo time-step size. ! ! DTsizeV Vertical diffusion pseudo time-step size. ! ! Kh Horizontal diffusion coefficients. ! ! Kv Vertical diffusion coefficients. ! ! A 3D state variable to diffuse. ! ! ! ! On Output: ! ! ! ! A Convolved 3D state variable. ! ! ! ! Routines: ! ! ! ! conv_r3d_tile Nonlinear 3D convolution at RHO-points ! ! conv_u3d_tile Nonlinear 3D convolution at U-points ! ! conv_v3d_tile Nonlinear 3D convolution at V-points ! ! ! !======================================================================= ! implicit none ! PUBLIC ! CONTAINS ! !*********************************************************************** SUBROUTINE conv_r3d_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, LBk, UBk, & & IminS, ImaxS, JminS, JmaxS, & & Nghost, NHsteps, NVsteps, & & DTsizeH, DTsizeV, & & Kh, Kv, & & pm, pn, & & pmon_u, pnom_v, & & rmask, umask, vmask, & & Hz, z_r, & & A) !*********************************************************************** ! USE mod_param USE mod_scalars ! USE bc_3d_mod, ONLY: dabc_r3d_tile USE mp_exchange_mod, ONLY : mp_exchange3d ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model integer, intent(in) :: LBi, UBi, LBj, UBj, LBk, UBk integer, intent(in) :: IminS, ImaxS, JminS, JmaxS integer, intent(in) :: Nghost, NHsteps, NVsteps real(r8), intent(in) :: DTsizeH, DTsizeV ! real(r8), intent(in) :: pm(LBi:,LBj:) real(r8), intent(in) :: pn(LBi:,LBj:) real(r8), intent(in) :: pmon_u(LBi:,LBj:) real(r8), intent(in) :: pnom_v(LBi:,LBj:) 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) :: z_r(LBi:,LBj:,:) real(r8), intent(in) :: Kh(LBi:,LBj:) real(r8), intent(in) :: Kv(LBi:,LBj:,0:) real(r8), intent(inout) :: A(LBi:,LBj:,LBk:) ! ! Local variable declarations. ! integer :: Nnew, Nold, Nsav, i, j, k, k1, k2, step real(r8) :: cff, cff1, cff2, cff3, cff4 real(r8), dimension(LBi:UBi,LBj:UBj,LBk:UBk,2) :: Awrk real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FE real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FX real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Hfac real(r8), dimension(IminS:ImaxS,JminS:JmaxS,0:N(ng)) :: FC real(r8), dimension(IminS:ImaxS,0:N(ng)) :: BC real(r8), dimension(IminS:ImaxS,0:N(ng)) :: CF real(r8), dimension(IminS:ImaxS,0:N(ng)) :: 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 ! !----------------------------------------------------------------------- ! Space convolution of the diffusion equation for a 3D state variable ! at RHO-points. !----------------------------------------------------------------------- ! ! Compute metrics factors. Notice that "z_r" and "Hz" are assumed to ! be time invariant in the vertical convolution. Scratch array are ! used for efficiency. ! DO j=Jstr-1,Jend+1 DO i=Istr-1,Iend+1 Hfac(i,j)=DTsizeH*pm(i,j)*pn(i,j) FC(i,j,N(ng))=0.0_r8 DO k=1,N(ng)-1 FC(i,j,k)=-DTsizeV*Kv(i,j,k)/(z_r(i,j,k+1)-z_r(i,j,k)) END DO FC(i,j,0)=0.0_r8 END DO END DO ! ! Set integration indices and initial conditions. ! Nold=1 Nnew=2 CALL dabc_r3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, LBk, UBk, & & A) CALL mp_exchange3d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, LBk, UBk, & & Nghost, & & EWperiodic(ng), NSperiodic(ng), & & A) DO k=1,N(ng) DO j=Jstr-1,Jend+1 DO i=Istr-1,Iend+1 Awrk(i,j,k,Nold)=A(i,j,k) END DO END DO END DO ! !----------------------------------------------------------------------- ! Integrate horizontal diffusion equation. !----------------------------------------------------------------------- ! DO step=1,NHsteps ! ! Diffusion along S-coordinates: compute XI- and ETA-components of ! diffusive flux. ! DO k=1,N(ng) DO j=Jstr,Jend DO i=Istr,Iend+1 FX(i,j)=pmon_u(i,j)*0.5_r8*(Kh(i-1,j)+Kh(i,j))* & & (Awrk(i,j,k,Nold)-Awrk(i-1,j,k,Nold)) FX(i,j)=FX(i,j)*umask(i,j) END DO END DO DO j=Jstr,Jend+1 DO i=Istr,Iend FE(i,j)=pnom_v(i,j)*0.5_r8*(Kh(i,j-1)+Kh(i,j))* & & (Awrk(i,j,k,Nold)-Awrk(i,j-1,k,Nold)) FE(i,j)=FE(i,j)*vmask(i,j) END DO END DO ! ! Time-step horizontal diffusion equation. ! DO j=Jstr,Jend DO i=Istr,Iend Awrk(i,j,k,Nnew)=Awrk(i,j,k,Nold)+ & & Hfac(i,j)* & & (FX(i+1,j)-FX(i,j)+ & & FE(i,j+1)-FE(i,j)) END DO END DO END DO ! ! Apply boundary conditions. ! CALL dabc_r3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, LBk, UBk, & & Awrk(:,:,:,Nnew)) CALL mp_exchange3d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, LBk, UBk, & & Nghost, & & EWperiodic(ng), NSperiodic(ng), & & Awrk(:,:,:,Nnew)) ! ! Update integration indices. ! Nsav=Nold Nold=Nnew Nnew=Nsav END DO ! !----------------------------------------------------------------------- ! Integrate vertical diffusion equation implicitly. !----------------------------------------------------------------------- ! DO step=1,NVsteps ! ! Compute diagonal matrix coefficients BC and load right-hand-side ! terms for the diffusion equation into DC. ! DO j=Jstr,Jend DO k=1,N(ng) DO i=Istr,Iend BC(i,k)=Hz(i,j,k)-FC(i,j,k)-FC(i,j,k-1) DC(i,k)=Awrk(i,j,k,Nold)*Hz(i,j,k) END DO END DO ! ! Solve the tridiagonal system. ! DO i=Istr,Iend cff=1.0_r8/BC(i,1) CF(i,1)=cff*FC(i,j,1) DC(i,1)=cff*DC(i,1) END DO DO k=2,N(ng)-1 DO i=Istr,Iend cff=1.0_r8/(BC(i,k)-FC(i,j,k-1)*CF(i,k-1)) CF(i,k)=cff*FC(i,j,k) DC(i,k)=cff*(DC(i,k)-FC(i,j,k-1)*DC(i,k-1)) END DO END DO ! ! Compute new solution by back substitution. ! DO i=Istr,Iend DC(i,N(ng))=(DC(i,N(ng))- & & FC(i,j,N(ng)-1)*DC(i,N(ng)-1))/ & & (BC(i,N(ng))-FC(i,j,N(ng)-1)*CF(i,N(ng)-1)) Awrk(i,j,N(ng),Nnew)=DC(i,N(ng)) Awrk(i,j,N(ng),Nnew)=Awrk(i,j,N(ng),Nnew)*rmask(i,j) END DO DO k=N(ng)-1,1,-1 DO i=Istr,Iend DC(i,k)=DC(i,k)-CF(i,k)*DC(i,k+1) Awrk(i,j,k,Nnew)=DC(i,k) Awrk(i,j,k,Nnew)=Awrk(i,j,k,Nnew)*rmask(i,j) END DO END DO END DO ! ! Update integration indices. ! Nsav=Nold Nold=Nnew Nnew=Nsav END DO ! !----------------------------------------------------------------------- ! Load convolved solution. !----------------------------------------------------------------------- ! DO k=1,N(ng) DO j=Jstr,Jend DO i=Istr,Iend A(i,j,k)=Awrk(i,j,k,Nold) END DO END DO END DO CALL dabc_r3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, LBk, UBk, & & A) CALL mp_exchange3d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, LBk, UBk, & & Nghost, & & EWperiodic(ng), NSperiodic(ng), & & A) RETURN END SUBROUTINE conv_r3d_tile ! !*********************************************************************** SUBROUTINE conv_u3d_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, LBk, UBk, & & IminS, ImaxS, JminS, JmaxS, & & Nghost, NHsteps, NVsteps, & & DTsizeH, DTsizeV, & & Kh, Kv, & & pm, pn, & & pmon_r, pnom_p, & & umask, pmask, & & Hz, z_r, & & A) !*********************************************************************** ! USE mod_param USE mod_scalars ! USE bc_3d_mod, ONLY: dabc_u3d_tile USE mp_exchange_mod, ONLY : mp_exchange3d ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model integer, intent(in) :: LBi, UBi, LBj, UBj, LBk, UBk integer, intent(in) :: IminS, ImaxS, JminS, JmaxS integer, intent(in) :: Nghost, NHsteps, NVsteps real(r8), intent(in) :: DTsizeH, DTsizeV ! real(r8), intent(in) :: pm(LBi:,LBj:) real(r8), intent(in) :: pn(LBi:,LBj:) real(r8), intent(in) :: pmon_r(LBi:,LBj:) real(r8), intent(in) :: pnom_p(LBi:,LBj:) real(r8), intent(in) :: umask(LBi:,LBj:) real(r8), intent(in) :: pmask(LBi:,LBj:) real(r8), intent(in) :: Hz(LBi:,LBj:,:) real(r8), intent(in) :: z_r(LBi:,LBj:,:) real(r8), intent(in) :: Kh(LBi:,LBj:) real(r8), intent(in) :: Kv(LBi:,LBj:,0:) real(r8), intent(inout) :: A(LBi:,LBj:,LBk:) ! ! Local variable declarations. ! integer :: Nnew, Nold, Nsav, i, j, k, k1, k2, step real(r8) :: cff, cff1, cff2, cff3, cff4 real(r8), dimension(LBi:UBi,LBj:UBj,LBk:UBk,2) :: Awrk real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FE real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FX real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Hfac real(r8), dimension(IminS:ImaxS,JminS:JmaxS,0:N(ng)) :: FC real(r8), dimension(IminS:ImaxS,0:N(ng)) :: BC real(r8), dimension(IminS:ImaxS,0:N(ng)) :: CF real(r8), dimension(IminS:ImaxS,0:N(ng)) :: 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 ! !----------------------------------------------------------------------- ! Space convolution of the diffusion equation for a 3D state variable ! at U-points. !----------------------------------------------------------------------- ! ! Compute metrics factors. Notice that "z_r" and "Hz" are assumed to ! be time invariant in the vertical convolution. Scratch array are ! used for efficiency. ! cff=DTsizeH*0.25_r8 DO j=Jstr-1,Jend+1 DO i=IstrU-1,Iend+1 Hfac(i,j)=cff*(pm(i-1,j)+pm(i,j))*(pn(i-1,j)+pn(i,j)) FC(i,j,N(ng))=0.0_r8 DO k=1,N(ng)-1 FC(i,j,k)=-DTsizeV*(Kv(i-1,j,k)+Kv(i,j,k))/ & & (z_r(i-1,j,k+1)+z_r(i,j,k+1)- & & z_r(i-1,j,k )-z_r(i,j,k )) END DO FC(i,j,0)=0.0_r8 END DO END DO ! ! Set integration indices and initial conditions. ! Nold=1 Nnew=2 CALL dabc_u3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, LBk, UBk, & & A) CALL mp_exchange3d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, LBk, UBk, & & Nghost, & & EWperiodic(ng), NSperiodic(ng), & & A) DO k=1,N(ng) DO j=Jstr-1,Jend+1 DO i=IstrU-1,Iend+1 Awrk(i,j,k,Nold)=A(i,j,k) END DO END DO END DO ! !----------------------------------------------------------------------- ! Integrate horizontal diffusion equation. !----------------------------------------------------------------------- ! DO step=1,NHsteps ! ! Diffusion along S-coordinates: compute XI- and ETA-components of ! diffusive flux. ! DO k=1,N(ng) DO j=Jstr,Jend DO i=IstrU-1,Iend FX(i,j)=pmon_r(i,j)*Kh(i,j)* & & (Awrk(i+1,j,k,Nold)-Awrk(i,j,k,Nold)) END DO END DO DO j=Jstr,Jend+1 DO i=IstrU,Iend FE(i,j)=pnom_p(i,j)*0.25_r8*(Kh(i-1,j )+Kh(i,j )+ & & Kh(i-1,j-1)+Kh(i,j-1))* & & (Awrk(i,j,k,Nold)-Awrk(i,j-1,k,Nold)) FE(i,j)=FE(i,j)*pmask(i,j) END DO END DO ! ! Time-step horizontal diffusion equation. ! DO j=Jstr,Jend DO i=IstrU,Iend Awrk(i,j,k,Nnew)=Awrk(i,j,k,Nold)+ & & Hfac(i,j)* & & (FX(i,j)-FX(i-1,j)+ & & FE(i,j+1)-FE(i,j)) END DO END DO END DO ! ! Apply boundary conditions. ! CALL dabc_u3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, LBk, UBk, & & Awrk(:,:,:,Nnew)) CALL mp_exchange3d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, LBk, UBk, & & Nghost, & & EWperiodic(ng), NSperiodic(ng), & & Awrk(:,:,:,Nnew)) ! ! Update integration indices. ! Nsav=Nold Nold=Nnew Nnew=Nsav END DO ! !----------------------------------------------------------------------- ! Integrate vertical diffusion equation implicitly. !----------------------------------------------------------------------- ! DO step=1,NVsteps ! ! Compute diagonal matrix coefficients BC and load right-hand-side ! terms for the diffusion equation into DC. ! DO j=Jstr,Jend DO k=1,N(ng) DO i=IstrU,Iend cff=0.5_r8*(Hz(i-1,j,k)+Hz(i,j,k)) BC(i,k)=cff-FC(i,j,k)-FC(i,j,k-1) DC(i,k)=Awrk(i,j,k,Nold)*cff END DO END DO ! ! Solve the tridiagonal system. ! DO i=IstrU,Iend cff=1.0_r8/BC(i,1) CF(i,1)=cff*FC(i,j,1) DC(i,1)=cff*DC(i,1) END DO DO k=2,N(ng)-1 DO i=IstrU,Iend cff=1.0_r8/(BC(i,k)-FC(i,j,k-1)*CF(i,k-1)) CF(i,k)=cff*FC(i,j,k) DC(i,k)=cff*(DC(i,k)-FC(i,j,k-1)*DC(i,k-1)) END DO END DO ! ! Compute new solution by back substitution. ! DO i=IstrU,Iend DC(i,N(ng))=(DC(i,N(ng))- & & FC(i,j,N(ng)-1)*DC(i,N(ng)-1))/ & & (BC(i,N(ng))-FC(i,j,N(ng)-1)*CF(i,N(ng)-1)) Awrk(i,j,N(ng),Nnew)=DC(i,N(ng)) Awrk(i,j,N(ng),Nnew)=Awrk(i,j,N(ng),Nnew)*umask(i,j) END DO DO k=N(ng)-1,1,-1 DO i=IstrU,Iend DC(i,k)=DC(i,k)-CF(i,k)*DC(i,k+1) Awrk(i,j,k,Nnew)=DC(i,k) Awrk(i,j,k,Nnew)=Awrk(i,j,k,Nnew)*umask(i,j) END DO END DO END DO ! ! Update integration indices. ! Nsav=Nold Nold=Nnew Nnew=Nsav END DO ! !----------------------------------------------------------------------- ! Load convolved solution. !----------------------------------------------------------------------- ! DO k=1,N(ng) DO j=Jstr,Jend DO i=IstrU,Iend A(i,j,k)=Awrk(i,j,k,Nold) END DO END DO END DO CALL dabc_u3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, LBk, UBk, & & A) CALL mp_exchange3d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, LBk, UBk, & & Nghost, & & EWperiodic(ng), NSperiodic(ng), & & A) RETURN END SUBROUTINE conv_u3d_tile ! !*********************************************************************** SUBROUTINE conv_v3d_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, LBk, UBk, & & IminS, ImaxS, JminS, JmaxS, & & Nghost, NHsteps, NVsteps, & & DTsizeH, DTsizeV, & & Kh, Kv, & & pm, pn, & & pmon_p, pnom_r, & & vmask, pmask, & & Hz, z_r, & & A) !*********************************************************************** ! USE mod_param USE mod_scalars ! USE bc_3d_mod, ONLY: dabc_v3d_tile USE mp_exchange_mod, ONLY : mp_exchange3d ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model integer, intent(in) :: LBi, UBi, LBj, UBj, LBk, UBk integer, intent(in) :: IminS, ImaxS, JminS, JmaxS integer, intent(in) :: Nghost, NHsteps, NVsteps real(r8), intent(in) :: DTsizeH, DTsizeV ! real(r8), intent(in) :: pm(LBi:,LBj:) real(r8), intent(in) :: pn(LBi:,LBj:) real(r8), intent(in) :: pmon_p(LBi:,LBj:) real(r8), intent(in) :: pnom_r(LBi:,LBj:) real(r8), intent(in) :: vmask(LBi:,LBj:) real(r8), intent(in) :: pmask(LBi:,LBj:) real(r8), intent(in) :: Hz(LBi:,LBj:,:) real(r8), intent(in) :: z_r(LBi:,LBj:,:) real(r8), intent(in) :: Kh(LBi:,LBj:) real(r8), intent(in) :: Kv(LBi:,LBj:,0:) real(r8), intent(inout) :: A(LBi:,LBj:,LBk:) ! ! Local variable declarations. ! integer :: Nnew, Nold, Nsav, i, j, k, k1, k2, step real(r8) :: cff, cff1, cff2, cff3, cff4 real(r8), dimension(LBi:UBi,LBj:UBj,LBk:UBk,2) :: Awrk real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FE real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FX real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Hfac real(r8), dimension(IminS:ImaxS,JminS:JmaxS,0:N(ng)) :: FC real(r8), dimension(IminS:ImaxS,0:N(ng)) :: BC real(r8), dimension(IminS:ImaxS,0:N(ng)) :: CF real(r8), dimension(IminS:ImaxS,0:N(ng)) :: 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 ! !----------------------------------------------------------------------- ! Space convolution of the diffusion equation for a 3D state variable ! at V-points. !----------------------------------------------------------------------- ! ! Compute metrics factors. Notice that "z_r" and "Hz" are assumed to ! be time invariant in the vertical convolution. Scratch array are ! used for efficiency. ! cff=DTsizeH*0.25_r8 DO j=JstrV-1,Jend+1 DO i=Istr-1,Iend+1 Hfac(i,j)=cff*(pm(i,j-1)+pm(i,j))*(pn(i,j-1)+pn(i,j)) FC(i,j,N(ng))=0.0_r8 DO k=1,N(ng)-1 FC(i,j,k)=-DTsizeV*(Kv(i,j-1,k)+Kv(i,j,k))/ & & (z_r(i,j-1,k+1)+z_r(i,j,k+1)- & & z_r(i,j-1,k )-z_r(i,j,k )) END DO FC(i,j,0)=0.0_r8 END DO END DO ! ! Set integration indices and initial conditions. ! Nold=1 Nnew=2 CALL dabc_v3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, LBk, UBk, & & A) CALL mp_exchange3d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, LBk, UBk, & & Nghost, & & EWperiodic(ng), NSperiodic(ng), & & A) DO k=1,N(ng) DO j=JstrV-1,Jend+1 DO i=Istr-1,Iend+1 Awrk(i,j,k,Nold)=A(i,j,k) END DO END DO END DO ! !----------------------------------------------------------------------- ! Integrate horizontal diffusion equation. !----------------------------------------------------------------------- ! DO step=1,NHsteps ! ! Compute XI- and ETA-components of diffusive flux. ! DO k=1,N(ng) DO j=JstrV,Jend DO i=Istr,Iend+1 FX(i,j)=pmon_p(i,j)*0.25_r8*(Kh(i-1,j )+Kh(i,j )+ & & Kh(i-1,j-1)+Kh(i,j-1))* & & (Awrk(i,j,k,Nold)-Awrk(i-1,j,k,Nold)) FX(i,j)=FX(i,j)*pmask(i,j) END DO END DO DO j=JstrV-1,Jend DO i=Istr,Iend FE(i,j)=pnom_r(i,j)*Kh(i,j)* & & (Awrk(i,j+1,k,Nold)-Awrk(i,j,k,Nold)) END DO END DO ! ! Time-step horizontal diffusion equation. ! DO j=JstrV,Jend DO i=Istr,Iend Awrk(i,j,k,Nnew)=Awrk(i,j,k,Nold)+ & & Hfac(i,j)* & & (FX(i+1,j)-FX(i,j)+ & & FE(i,j)-FE(i,j-1)) END DO END DO END DO ! ! Apply boundary conditions. ! CALL dabc_v3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, LBk, UBk, & & Awrk(:,:,:,Nnew)) CALL mp_exchange3d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, LBk, UBk, & & Nghost, & & EWperiodic(ng), NSperiodic(ng), & & Awrk(:,:,:,Nnew)) ! ! Update integration indices. ! Nsav=Nold Nold=Nnew Nnew=Nsav END DO ! !----------------------------------------------------------------------- ! Integerate vertical diffusion equation implicitly. !----------------------------------------------------------------------- ! DO step=1,NVsteps ! ! Compute diagonal matrix coefficients BC and load right-hand-side ! terms for the diffusion equation into DC. ! DO j=JstrV,Jend DO k=1,N(ng) DO i=Istr,Iend cff=0.5_r8*(Hz(i,j-1,k)+Hz(i,j,k)) BC(i,k)=cff-FC(i,j,k)-FC(i,j,k-1) DC(i,k)=Awrk(i,j,k,Nold)*cff END DO END DO ! ! Solve the tridiagonal system. ! DO i=Istr,Iend cff=1.0_r8/BC(i,1) CF(i,1)=cff*FC(i,j,1) DC(i,1)=cff*DC(i,1) END DO DO k=2,N(ng)-1 DO i=Istr,Iend cff=1.0_r8/(BC(i,k)-FC(i,j,k-1)*CF(i,k-1)) CF(i,k)=cff*FC(i,j,k) DC(i,k)=cff*(DC(i,k)-FC(i,j,k-1)*DC(i,k-1)) END DO END DO ! ! Compute new solution by back substitution. ! DO i=Istr,Iend DC(i,N(ng))=(DC(i,N(ng))- & & FC(i,j,N(ng)-1)*DC(i,N(ng)-1))/ & & (BC(i,N(ng))-FC(i,j,N(ng)-1)*CF(i,N(ng)-1)) Awrk(i,j,N(ng),Nnew)=DC(i,N(ng)) Awrk(i,j,N(ng),Nnew)=Awrk(i,j,N(ng),Nnew)*vmask(i,j) END DO DO k=N(ng)-1,1,-1 DO i=Istr,Iend DC(i,k)=DC(i,k)-CF(i,k)*DC(i,k+1) Awrk(i,j,k,Nnew)=DC(i,k) Awrk(i,j,k,Nnew)=Awrk(i,j,k,Nnew)*vmask(i,j) END DO END DO END DO ! ! Update integration indices. ! Nsav=Nold Nold=Nnew Nnew=Nsav END DO ! !----------------------------------------------------------------------- ! Load convolved solution. !----------------------------------------------------------------------- ! DO k=1,N(ng) DO j=JstrV,Jend DO i=Istr,Iend A(i,j,k)=Awrk(i,j,k,Nold) END DO END DO END DO CALL dabc_v3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, LBk, UBk, & & A) CALL mp_exchange3d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, LBk, UBk, & & Nghost, & & EWperiodic(ng), NSperiodic(ng), & & A) RETURN END SUBROUTINE conv_v3d_tile END MODULE conv_3d_mod