MODULE mod_coupling ! !git $Id$ !svn $Id: mod_coupling.F 1188 2023-08-03 19:26:47Z arango $ !================================================== Hernan G. Arango === ! Copyright (c) 2002-2023 The ROMS/TOMS Group ! ! Licensed under a MIT/X style license ! ! See License_ROMS.md ! !======================================================================= ! ! ! DU_avg1 Time averaged U-flux for 2D equations (m3/s). ! ! DU_avg2 Time averaged U-flux for 3D equations coupling (m3/s). ! ! DV_avg1 Time averaged V-flux for 2D equations (m3/s). ! ! DV_avg2 Time averaged V-flux for 3D equations coupling (m3/s). ! ! Zt_avg1 Free-surface averaged over all short time-steps (m). ! ! rhoA Normalized vertical averaged density. ! ! rhoS Normalized vertical averaged density perturbation. ! ! rufrc Right-hand-side forcing term for 2D U-momentum (m4/s2) ! ! rvfrc Right-hand-side forcing term for 2D V-momentum (m4/s2) ! ! ! !======================================================================= ! USE mod_kinds ! implicit none ! PUBLIC :: allocate_coupling PUBLIC :: deallocate_coupling PUBLIC :: initialize_coupling ! !----------------------------------------------------------------------- ! Define T_COUPLING structure. !----------------------------------------------------------------------- ! TYPE T_COUPLING ! ! Nonlinear model state. ! real(r8), pointer :: DU_avg1(:,:) real(r8), pointer :: DU_avg2(:,:) real(r8), pointer :: DV_avg1(:,:) real(r8), pointer :: DV_avg2(:,:) real(r8), pointer :: Zt_avg1(:,:) real(r8), pointer :: rufrc(:,:) real(r8), pointer :: rvfrc(:,:) real(r8), pointer :: rhoA(:,:) real(r8), pointer :: rhoS(:,:) ! ! Tangent linear model state. ! real(r8), pointer :: tl_DU_avg1(:,:) real(r8), pointer :: tl_DU_avg2(:,:) real(r8), pointer :: tl_DV_avg1(:,:) real(r8), pointer :: tl_DV_avg2(:,:) real(r8), pointer :: tl_Zt_avg1(:,:) real(r8), pointer :: tl_rufrc(:,:) real(r8), pointer :: tl_rvfrc(:,:) ! ! Adjoint model state. ! real(r8), pointer :: ad_DU_avg1(:,:) real(r8), pointer :: ad_DU_avg2(:,:) real(r8), pointer :: ad_DV_avg1(:,:) real(r8), pointer :: ad_DV_avg2(:,:) real(r8), pointer :: ad_Zt_avg1(:,:) real(r8), pointer :: ad_rufrc(:,:) real(r8), pointer :: ad_rvfrc(:,:) ! ! Latest two records of the nonlinear trajectory used to interpolate ! the background state in the tangent linear and adjoint models. ! real(r8), pointer :: DU_avg1G(:,:,:) real(r8), pointer :: DU_avg2G(:,:,:) real(r8), pointer :: DV_avg1G(:,:,:) real(r8), pointer :: DV_avg2G(:,:,:) real(r8), pointer :: rufrcG(:,:,:) real(r8), pointer :: rvfrcG(:,:,:) END TYPE T_COUPLING ! TYPE (T_COUPLING), allocatable :: COUPLING(:) ! CONTAINS ! SUBROUTINE allocate_coupling (ng, LBi, UBi, LBj, UBj) ! !======================================================================= ! ! ! This routine allocates all variables in the module for all nested ! ! grids. ! ! ! !======================================================================= ! USE mod_param ! ! Imported variable declarations. ! integer, intent(in) :: ng, LBi, UBi, LBj, UBj ! ! Local variable declarations. ! real(r8) :: size2d ! !----------------------------------------------------------------------- ! Initialize module variables. !----------------------------------------------------------------------- ! IF (ng.eq.1) allocate ( COUPLING(Ngrids) ) ! ! Set horizontal array size. ! size2d=REAL((UBi-LBi+1)*(UBj-LBj+1),r8) ! ! Nonlinear model state. ! allocate ( COUPLING(ng) % DU_avg1(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % DU_avg2(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % DV_avg1(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % DV_avg2(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % Zt_avg1(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % rufrc(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % rvfrc(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % rhoA(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % rhoS(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d ! ! Tangent linear model state. ! allocate ( COUPLING(ng) % tl_DU_avg1(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % tl_DU_avg2(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % tl_DV_avg1(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % tl_DV_avg2(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % tl_Zt_avg1(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % tl_rufrc(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % tl_rvfrc(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d ! ! Adjoint model state. ! allocate ( COUPLING(ng) % ad_DU_avg1(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % ad_DU_avg2(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % ad_DV_avg1(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % ad_DV_avg2(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % ad_Zt_avg1(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % ad_rufrc(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( COUPLING(ng) % ad_rvfrc(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d ! ! Latest two records of the nonlinear trajectory used to interpolate ! the background state in the tangent linear and adjoint models. ! allocate ( COUPLING(ng) % DU_avg1G(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( COUPLING(ng) % DU_avg2G(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( COUPLING(ng) % DV_avg1G(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( COUPLING(ng) % DV_avg2G(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( COUPLING(ng) % rufrcG(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( COUPLING(ng) % rvfrcG(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d ! RETURN END SUBROUTINE allocate_coupling ! SUBROUTINE deallocate_coupling (ng) ! !======================================================================= ! ! ! This routine deallocates all variables in the module for all nested ! ! grids. ! ! ! !======================================================================= ! USE mod_param, ONLY : Ngrids ! ! Imported variable declarations. ! integer, intent(in) :: ng ! ! Local variable declarations. ! character (len=*), parameter :: MyFile = & & "ROMS/Modules/mod_coupling.F"//", deallocate_coupling" ! !----------------------------------------------------------------------- ! Deallocate derived-type COUPLING structure. !----------------------------------------------------------------------- ! IF (ng.eq.Ngrids) THEN IF (allocated(COUPLING)) deallocate ( COUPLING ) END IF ! RETURN END SUBROUTINE deallocate_coupling ! SUBROUTINE initialize_coupling (ng, tile, model) ! !======================================================================= ! ! ! This routine initialize all variables in the module using first ! ! touch distribution policy. In shared-memory configuration, this ! ! operation actually performs propagation of the "shared arrays" ! ! across the cluster, unless another policy is specified to ! ! override the default. ! ! ! !======================================================================= ! USE mod_param ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model ! ! Local variable declarations. ! integer :: Imin, Imax, Jmin, Jmax integer :: i, j real(r8), parameter :: IniVal = 0.0_r8 ! !----------------------------------------------------------------------- ! 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 ! ! Set array initialization range. ! Imin=BOUNDS(ng)%LBi(tile) Imax=BOUNDS(ng)%UBi(tile) Jmin=BOUNDS(ng)%LBj(tile) Jmax=BOUNDS(ng)%UBj(tile) ! !----------------------------------------------------------------------- ! Initialize module variables. !----------------------------------------------------------------------- ! ! Nonlinear model state. ! IF ((model.eq.0).or.(model.eq.iNLM)) THEN DO j=Jmin,Jmax DO i=Imin,Imax COUPLING(ng) % DU_avg1(i,j) = IniVal COUPLING(ng) % DU_avg2(i,j) = IniVal COUPLING(ng) % DV_avg1(i,j) = IniVal COUPLING(ng) % DV_avg2(i,j) = IniVal COUPLING(ng) % Zt_avg1(i,j) = IniVal COUPLING(ng) % rufrc(i,j) = IniVal COUPLING(ng) % rvfrc(i,j) = IniVal COUPLING(ng) % rhoA(i,j) = IniVal COUPLING(ng) % rhoS(i,j) = IniVal END DO END DO END IF ! ! Tangent linear model state. ! IF ((model.eq.0).or.(model.eq.iTLM).or.(model.eq.iRPM)) THEN DO j=Jmin,Jmax DO i=Imin,Imax COUPLING(ng) % tl_DU_avg1(i,j) = IniVal COUPLING(ng) % tl_DU_avg2(i,j) = IniVal COUPLING(ng) % tl_DV_avg1(i,j) = IniVal COUPLING(ng) % tl_DV_avg2(i,j) = IniVal COUPLING(ng) % tl_Zt_avg1(i,j) = IniVal COUPLING(ng) % tl_rufrc(i,j) = IniVal COUPLING(ng) % tl_rvfrc(i,j) = IniVal END DO END DO END IF ! ! Adjoint model state. ! IF ((model.eq.0).or.(model.eq.iADM)) THEN DO j=Jmin,Jmax DO i=Imin,Imax COUPLING(ng) % ad_DU_avg1(i,j) = IniVal COUPLING(ng) % ad_DU_avg2(i,j) = IniVal COUPLING(ng) % ad_DV_avg1(i,j) = IniVal COUPLING(ng) % ad_DV_avg2(i,j) = IniVal COUPLING(ng) % ad_Zt_avg1(i,j) = IniVal COUPLING(ng) % ad_rufrc(i,j) = IniVal COUPLING(ng) % ad_rvfrc(i,j) = IniVal END DO END DO END IF ! ! Latest two records of the nonlinear trajectory used to interpolate ! the background state in the tangent linear and adjoint models. ! IF (model.eq.0) THEN DO j=Jmin,Jmax DO i=Imin,Imax COUPLING(ng) % DU_avg1G(i,j,1) = IniVal COUPLING(ng) % DU_avg1G(i,j,2) = IniVal COUPLING(ng) % DU_avg2G(i,j,1) = IniVal COUPLING(ng) % DU_avg2G(i,j,2) = IniVal COUPLING(ng) % DV_avg1G(i,j,1) = IniVal COUPLING(ng) % DV_avg1G(i,j,2) = IniVal COUPLING(ng) % DV_avg2G(i,j,1) = IniVal COUPLING(ng) % DV_avg2G(i,j,2) = IniVal COUPLING(ng) % rufrcG(i,j,1) = IniVal COUPLING(ng) % rufrcG(i,j,2) = IniVal COUPLING(ng) % rvfrcG(i,j,1) = IniVal COUPLING(ng) % rvfrcG(i,j,2) = IniVal END DO END DO END IF ! RETURN END SUBROUTINE initialize_coupling END MODULE mod_coupling