MODULE mod_ocean ! !git $Id$ !svn $Id: mod_ocean.F 1178 2023-07-11 17:50:57Z arango $ !================================================== Hernan G. Arango === ! Copyright (c) 2002-2023 The ROMS/TOMS Group ! ! Licensed under a MIT/X style license ! ! See License_ROMS.md ! !======================================================================= ! ! ! 2D Primitive Variables. ! ! ! ! rubar Right-hand-side of 2D U-momentum equation (m4/s2). ! ! rvbar Right-hand-side of 2D V-momentum equation (m4/s2). ! ! rzeta Right-hand-side of free surface equation (m3/s). ! ! ubar Vertically integrated U-momentum component (m/s). ! ! vbar Vertically integrated V-momentum component (m/s). ! ! zeta Free surface (m). ! ! ! ! 3D Primitive Variables. ! ! ! ! pden Potential Density anomaly (kg/m3). ! ! rho Density anomaly (kg/m3). ! ! ru Right-hand-side of 3D U-momentum equation (m4/s2). ! ! rv Right hand side of 3D V-momentum equation (m4/s2). ! ! t Tracer type variables (active and passive). ! ! u 3D U-momentum component (m/s). ! ! v 3D V-momentum component (m/s). ! ! W S-coordinate (omega*Hz/mn) vertical velocity (m3/s). ! ! ! !======================================================================= ! USE mod_kinds ! implicit none ! PUBLIC :: allocate_ocean PUBLIC :: deallocate_ocean PUBLIC :: initialize_ocean ! !----------------------------------------------------------------------- ! Define T_OCEAN structure. !----------------------------------------------------------------------- ! TYPE T_OCEAN ! ! Nonlinear model state. ! real(r8), pointer :: rubar(:,:,:) real(r8), pointer :: rvbar(:,:,:) real(r8), pointer :: rzeta(:,:,:) real(r8), pointer :: ubar(:,:,:) real(r8), pointer :: vbar(:,:,:) real(r8), pointer :: zeta(:,:,:) real(r8), pointer :: pden(:,:,:) real(r8), pointer :: rho(:,:,:) real(r8), pointer :: ru(:,:,:,:) real(r8), pointer :: rv(:,:,:,:) real(r8), pointer :: t(:,:,:,:,:) real(r8), pointer :: u(:,:,:,:) real(r8), pointer :: v(:,:,:,:) real(r8), pointer :: W(:,:,:) real(r8), pointer :: wvel(:,:,:) ! ! Tangent linear model state. ! real(r8), pointer :: tl_rubar(:,:,:) real(r8), pointer :: tl_rvbar(:,:,:) real(r8), pointer :: tl_rzeta(:,:,:) real(r8), pointer :: tl_ubar(:,:,:) real(r8), pointer :: tl_vbar(:,:,:) real(r8), pointer :: tl_zeta(:,:,:) real(r8), pointer :: tl_pden(:,:,:) real(r8), pointer :: tl_rho(:,:,:) real(r8), pointer :: tl_ru(:,:,:,:) real(r8), pointer :: tl_rv(:,:,:,:) real(r8), pointer :: tl_t(:,:,:,:,:) real(r8), pointer :: tl_u(:,:,:,:) real(r8), pointer :: tl_v(:,:,:,:) real(r8), pointer :: tl_W(:,:,:) ! ! Adjoint model state. ! real(r8), pointer :: ad_rubar(:,:,:) real(r8), pointer :: ad_rvbar(:,:,:) real(r8), pointer :: ad_rzeta(:,:,:) real(r8), pointer :: ad_ubar(:,:,:) real(r8), pointer :: ad_vbar(:,:,:) real(r8), pointer :: ad_zeta(:,:,:) real(r8), pointer :: ad_ubar_sol(:,:) real(r8), pointer :: ad_vbar_sol(:,:) real(r8), pointer :: ad_zeta_sol(:,:) real(r8), pointer :: ad_pden(:,:,:) real(r8), pointer :: ad_rho(:,:,:) real(r8), pointer :: ad_ru(:,:,:,:) real(r8), pointer :: ad_rv(:,:,:,:) real(r8), pointer :: ad_t(:,:,:,:,:) real(r8), pointer :: ad_u(:,:,:,:) real(r8), pointer :: ad_v(:,:,:,:) real(r8), pointer :: ad_W(:,:,:) real(r8), pointer :: ad_wvel(:,:,:) real(r8), pointer :: ad_t_sol(:,:,:,:) real(r8), pointer :: ad_u_sol(:,:,:) real(r8), pointer :: ad_v_sol(:,:,:) real(r8), pointer :: ad_W_sol(:,:,:) ! ! Working arrays to store adjoint impulse forcing, error covariance, ! standard deviations, or descent conjugate vectors (directions). ! real(r8), pointer :: b_ubar(:,:,:) real(r8), pointer :: b_vbar(:,:,:) real(r8), pointer :: b_zeta(:,:,:) real(r8), pointer :: b_t(:,:,:,:,:) real(r8), pointer :: b_u(:,:,:,:) real(r8), pointer :: b_v(:,:,:,:) real(r8), pointer :: d_ubar(:,:) real(r8), pointer :: d_vbar(:,:) real(r8), pointer :: d_zeta(:,:) real(r8), pointer :: d_t(:,:,:,:) real(r8), pointer :: d_u(:,:,:) real(r8), pointer :: d_v(:,:,:) real(r8), pointer :: e_ubar(:,:,:) real(r8), pointer :: e_vbar(:,:,:) real(r8), pointer :: e_zeta(:,:,:) real(r8), pointer :: e_t(:,:,:,:,:) real(r8), pointer :: e_u(:,:,:,:) real(r8), pointer :: e_v(:,:,:,:) real(r8), pointer :: f_ubar(:,:) real(r8), pointer :: f_vbar(:,:) real(r8), pointer :: f_zeta(:,:) real(r8), pointer :: f_t(:,:,:,:) real(r8), pointer :: f_u(:,:,:) real(r8), pointer :: f_v(:,:,:) ! ! Latest two records of the nonlinear trajectory used to interpolate ! the background state in the tangent linear and adjoint models. ! real(r8), pointer :: ubarG(:,:,:) real(r8), pointer :: vbarG(:,:,:) real(r8), pointer :: zetaG(:,:,:) real(r8), pointer :: tG(:,:,:,:,:) real(r8), pointer :: uG(:,:,:,:) real(r8), pointer :: vG(:,:,:,:) real(r8), pointer :: f_zetaG(:,:,:) real(r8), pointer :: f_tG(:,:,:,:,:) real(r8), pointer :: f_uG(:,:,:,:) real(r8), pointer :: f_vG(:,:,:,:) real(r8), pointer :: f_ubarG(:,:,:) real(r8), pointer :: f_vbarG(:,:,:) END TYPE T_OCEAN ! TYPE (T_OCEAN), allocatable :: OCEAN(:) ! CONTAINS ! SUBROUTINE allocate_ocean (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 ! !----------------------------------------------------------------------- ! Allocate and initialize module variables. !----------------------------------------------------------------------- ! IF (ng.eq.1) allocate ( OCEAN(Ngrids) ) ! ! Set horizontal array size. ! size2d=REAL((UBi-LBi+1)*(UBj-LBj+1),r8) ! ! Nonlinear model state. ! allocate ( OCEAN(ng) % rubar(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( OCEAN(ng) % rvbar(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( OCEAN(ng) % rzeta(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( OCEAN(ng) % ubar(LBi:UBi,LBj:UBj,3) ) Dmem(ng)=Dmem(ng)+3.0_r8*size2d allocate ( OCEAN(ng) % vbar(LBi:UBi,LBj:UBj,3) ) Dmem(ng)=Dmem(ng)+3.0_r8*size2d allocate ( OCEAN(ng) % zeta(LBi:UBi,LBj:UBj,3) ) Dmem(ng)=Dmem(ng)+3.0_r8*size2d allocate ( OCEAN(ng) % pden(LBi:UBi,LBj:UBj,N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % rho(LBi:UBi,LBj:UBj,N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % ru(LBi:UBi,LBj:UBj,0:N(ng),2) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng)+1,r8)*size2d allocate ( OCEAN(ng) % rv(LBi:UBi,LBj:UBj,0:N(ng),2) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng)+1,r8)*size2d allocate ( OCEAN(ng) % t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng)) ) Dmem(ng)=Dmem(ng)+3.0_r8*REAL(N(ng)*NT(ng),r8)*size2d allocate ( OCEAN(ng) % u(LBi:UBi,LBj:UBj,N(ng),2) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % v(LBi:UBi,LBj:UBj,N(ng),2) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % W(LBi:UBi,LBj:UBj,0:N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng)+1,r8)*size2d allocate ( OCEAN(ng) % wvel(LBi:UBi,LBj:UBj,0:N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng)+1,r8)*size2d ! ! Tangent linear model state. ! allocate ( OCEAN(ng) % tl_rubar(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( OCEAN(ng) % tl_rvbar(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( OCEAN(ng) % tl_rzeta(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( OCEAN(ng) % tl_ubar(LBi:UBi,LBj:UBj,3) ) Dmem(ng)=Dmem(ng)+3.0_r8*size2d allocate ( OCEAN(ng) % tl_vbar(LBi:UBi,LBj:UBj,3) ) Dmem(ng)=Dmem(ng)+3.0_r8*size2d allocate ( OCEAN(ng) % tl_zeta(LBi:UBi,LBj:UBj,3) ) Dmem(ng)=Dmem(ng)+3.0_r8*size2d allocate ( OCEAN(ng) % tl_pden(LBi:UBi,LBj:UBj,N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % tl_rho(LBi:UBi,LBj:UBj,N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % tl_ru(LBi:UBi,LBj:UBj,0:N(ng),2) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng)+1,r8)*size2d allocate ( OCEAN(ng) % tl_rv(LBi:UBi,LBj:UBj,0:N(ng),2) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng)+1,r8)*size2d allocate ( OCEAN(ng) % tl_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng)) ) Dmem(ng)=Dmem(ng)+3.0_r8*REAL(N(ng)*NT(ng),r8)*size2d allocate ( OCEAN(ng) % tl_u(LBi:UBi,LBj:UBj,N(ng),2) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % tl_v(LBi:UBi,LBj:UBj,N(ng),2) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % tl_W(LBi:UBi,LBj:UBj,0:N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng)+1,r8)*size2d ! ! Adjoint model state. ! allocate ( OCEAN(ng) % ad_rubar(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( OCEAN(ng) % ad_rvbar(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( OCEAN(ng) % ad_rzeta(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( OCEAN(ng) % ad_ubar(LBi:UBi,LBj:UBj,3) ) Dmem(ng)=Dmem(ng)+3.0_r8*size2d allocate ( OCEAN(ng) % ad_vbar(LBi:UBi,LBj:UBj,3) ) Dmem(ng)=Dmem(ng)+3.0_r8*size2d allocate ( OCEAN(ng) % ad_zeta(LBi:UBi,LBj:UBj,3) ) Dmem(ng)=Dmem(ng)+3.0_r8*size2d allocate ( OCEAN(ng) % ad_ubar_sol(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( OCEAN(ng) % ad_vbar_sol(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( OCEAN(ng) % ad_zeta_sol(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( OCEAN(ng) % ad_pden(LBi:UBi,LBj:UBj,N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % ad_rho(LBi:UBi,LBj:UBj,N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % ad_ru(LBi:UBi,LBj:UBj,0:N(ng),2) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng)+1,r8)*size2d allocate ( OCEAN(ng) % ad_rv(LBi:UBi,LBj:UBj,0:N(ng),2) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng)+1,r8)*size2d allocate ( OCEAN(ng) % ad_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng)) ) Dmem(ng)=Dmem(ng)+3.0_r8*REAL(N(ng)*NT(ng),r8)*size2d allocate ( OCEAN(ng) % ad_u(LBi:UBi,LBj:UBj,N(ng),2) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % ad_v(LBi:UBi,LBj:UBj,N(ng),2) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % ad_W(LBi:UBi,LBj:UBj,0:N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng)+1,r8)*size2d allocate ( OCEAN(ng) % ad_wvel(LBi:UBi,LBj:UBj,0:N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng)+1,r8)*size2d allocate ( OCEAN(ng) % ad_t_sol(LBi:UBi,LBj:UBj,N(ng),NT(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng)*NT(ng),r8)*size2d allocate ( OCEAN(ng) % ad_u_sol(LBi:UBi,LBj:UBj,N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % ad_v_sol(LBi:UBi,LBj:UBj,N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % ad_W_sol(LBi:UBi,LBj:UBj,0:N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng)+1,r8)*size2d ! ! Working arrays to store adjoint impulse forcing, background error ! covariance, background-error standard deviations, or descent ! conjugate vectors (directions). ! allocate ( OCEAN(ng) % b_ubar(LBi:UBi,LBj:UBj,NSA) ) Dmem(ng)=Dmem(ng)+REAL(NSA,r8)*size2d allocate ( OCEAN(ng) % b_vbar(LBi:UBi,LBj:UBj,NSA) ) Dmem(ng)=Dmem(ng)+REAL(NSA,r8)*size2d allocate ( OCEAN(ng) % b_zeta(LBi:UBi,LBj:UBj,NSA) ) Dmem(ng)=Dmem(ng)+REAL(NSA,r8)*size2d allocate ( OCEAN(ng) % b_t(LBi:UBi,LBj:UBj,N(ng),NSA,NT(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng)*NT(ng)*NSA,r8)*size2d allocate ( OCEAN(ng) % b_u(LBi:UBi,LBj:UBj,N(ng),NSA) ) Dmem(ng)=Dmem(ng)+REAL(N(ng)*NSA,r8)*size2d allocate ( OCEAN(ng) % b_v(LBi:UBi,LBj:UBj,N(ng),NSA) ) Dmem(ng)=Dmem(ng)+REAL(N(ng)*NSA,r8)*size2d allocate ( OCEAN(ng) % d_ubar(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( OCEAN(ng) % d_vbar(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( OCEAN(ng) % d_zeta(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( OCEAN(ng) % d_t(LBi:UBi,LBj:UBj,N(ng),NT(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng)*NT(ng),r8)*size2d allocate ( OCEAN(ng) % d_u(LBi:UBi,LBj:UBj,N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % d_v(LBi:UBi,LBj:UBj,N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % e_ubar(LBi:UBi,LBj:UBj,NSA) ) Dmem(ng)=Dmem(ng)+REAL(NSA,r8)*size2d allocate ( OCEAN(ng) % e_vbar(LBi:UBi,LBj:UBj,NSA) ) Dmem(ng)=Dmem(ng)+REAL(NSA,r8)*size2d allocate ( OCEAN(ng) % e_zeta(LBi:UBi,LBj:UBj,NSA) ) Dmem(ng)=Dmem(ng)+REAL(NSA,r8)*size2d allocate ( OCEAN(ng) % e_t(LBi:UBi,LBj:UBj,N(ng),NSA,NT(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng)*NT(ng)*NSA,r8)*size2d allocate ( OCEAN(ng) % e_u(LBi:UBi,LBj:UBj,N(ng),NSA) ) Dmem(ng)=Dmem(ng)+REAL(N(ng)*NSA,r8)*size2d allocate ( OCEAN(ng) % e_v(LBi:UBi,LBj:UBj,N(ng),NSA) ) Dmem(ng)=Dmem(ng)+REAL(N(ng)*NSA,r8)*size2d allocate ( OCEAN(ng) % f_t(LBi:UBi,LBj:UBj,N(ng),NT(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng)*NT(ng),r8)*size2d allocate ( OCEAN(ng) % f_u(LBi:UBi,LBj:UBj,N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % f_v(LBi:UBi,LBj:UBj,N(ng)) ) Dmem(ng)=Dmem(ng)+REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % f_ubar(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( OCEAN(ng) % f_vbar(LBi:UBi,LBj:UBj) ) Dmem(ng)=Dmem(ng)+size2d allocate ( OCEAN(ng) % f_zeta(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 ( OCEAN(ng) % ubarG(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( OCEAN(ng) % vbarG(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( OCEAN(ng) % zetaG(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( OCEAN(ng) % tG(LBi:UBi,LBj:UBj,N(ng),2,NT(ng)) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng)*NT(ng),r8)*size2d allocate ( OCEAN(ng) % uG(LBi:UBi,LBj:UBj,N(ng),2) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % vG(LBi:UBi,LBj:UBj,N(ng),2) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % f_tG(LBi:UBi,LBj:UBj,N(ng),2,NT(ng)) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng)*NT(ng),r8)*size2d allocate ( OCEAN(ng) % f_uG(LBi:UBi,LBj:UBj,N(ng),2) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % f_vG(LBi:UBi,LBj:UBj,N(ng),2) ) Dmem(ng)=Dmem(ng)+2.0_r8*REAL(N(ng),r8)*size2d allocate ( OCEAN(ng) % f_ubarG(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( OCEAN(ng) % f_vbarG(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d allocate ( OCEAN(ng) % f_zetaG(LBi:UBi,LBj:UBj,2) ) Dmem(ng)=Dmem(ng)+2.0_r8*size2d ! RETURN END SUBROUTINE allocate_ocean ! SUBROUTINE deallocate_ocean (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_ocean.F"//", deallocate_ocean" ! !----------------------------------------------------------------------- ! Deallocate derived-type OCEAN structure. !----------------------------------------------------------------------- ! IF (ng.eq.Ngrids) THEN IF (allocated(OCEAN)) deallocate ( OCEAN ) END IF ! RETURN END SUBROUTINE deallocate_ocean ! SUBROUTINE initialize_ocean (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, rec integer :: itrc, k 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 OCEAN(ng) % rubar(i,j,1) = IniVal OCEAN(ng) % rubar(i,j,2) = IniVal OCEAN(ng) % rvbar(i,j,1) = IniVal OCEAN(ng) % rvbar(i,j,2) = IniVal OCEAN(ng) % rzeta(i,j,1) = IniVal OCEAN(ng) % rzeta(i,j,2) = IniVal OCEAN(ng) % ubar(i,j,1) = IniVal OCEAN(ng) % ubar(i,j,2) = IniVal OCEAN(ng) % ubar(i,j,3) = IniVal OCEAN(ng) % vbar(i,j,1) = IniVal OCEAN(ng) % vbar(i,j,2) = IniVal OCEAN(ng) % vbar(i,j,3) = IniVal OCEAN(ng) % zeta(i,j,1) = IniVal OCEAN(ng) % zeta(i,j,2) = IniVal OCEAN(ng) % zeta(i,j,3) = IniVal END DO DO k=1,N(ng) DO i=Imin,Imax OCEAN(ng) % pden(i,j,k) = IniVal OCEAN(ng) % rho(i,j,k) = IniVal OCEAN(ng) % u(i,j,k,1) = IniVal OCEAN(ng) % u(i,j,k,2) = IniVal OCEAN(ng) % v(i,j,k,1) = IniVal OCEAN(ng) % v(i,j,k,2) = IniVal END DO END DO DO k=0,N(ng) DO i=Imin,Imax OCEAN(ng) % ru(i,j,k,1) = IniVal OCEAN(ng) % ru(i,j,k,2) = IniVal OCEAN(ng) % rv(i,j,k,1) = IniVal OCEAN(ng) % rv(i,j,k,2) = IniVal OCEAN(ng) % W(i,j,k) = IniVal OCEAN(ng) % wvel(i,j,k) = IniVal END DO END DO DO itrc=1,NT(ng) DO k=1,N(ng) DO i=Imin,Imax OCEAN(ng) % t(i,j,k,1,itrc) = IniVal OCEAN(ng) % t(i,j,k,2,itrc) = IniVal OCEAN(ng) % t(i,j,k,3,itrc) = IniVal END DO END DO 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 OCEAN(ng) % tl_rubar(i,j,1) = IniVal OCEAN(ng) % tl_rubar(i,j,2) = IniVal OCEAN(ng) % tl_rvbar(i,j,1) = IniVal OCEAN(ng) % tl_rvbar(i,j,2) = IniVal OCEAN(ng) % tl_rzeta(i,j,1) = IniVal OCEAN(ng) % tl_rzeta(i,j,2) = IniVal OCEAN(ng) % tl_ubar(i,j,1) = IniVal OCEAN(ng) % tl_ubar(i,j,2) = IniVal OCEAN(ng) % tl_ubar(i,j,3) = IniVal OCEAN(ng) % tl_vbar(i,j,1) = IniVal OCEAN(ng) % tl_vbar(i,j,2) = IniVal OCEAN(ng) % tl_vbar(i,j,3) = IniVal OCEAN(ng) % tl_zeta(i,j,1) = IniVal OCEAN(ng) % tl_zeta(i,j,2) = IniVal OCEAN(ng) % tl_zeta(i,j,3) = IniVal END DO DO k=1,N(ng) DO i=Imin,Imax OCEAN(ng) % tl_pden(i,j,k) = IniVal OCEAN(ng) % tl_rho(i,j,k) = IniVal OCEAN(ng) % tl_u(i,j,k,1) = IniVal OCEAN(ng) % tl_u(i,j,k,2) = IniVal OCEAN(ng) % tl_v(i,j,k,1) = IniVal OCEAN(ng) % tl_v(i,j,k,2) = IniVal END DO END DO DO k=0,N(ng) DO i=Imin,Imax OCEAN(ng) % tl_ru(i,j,k,1) = IniVal OCEAN(ng) % tl_ru(i,j,k,2) = IniVal OCEAN(ng) % tl_rv(i,j,k,1) = IniVal OCEAN(ng) % tl_rv(i,j,k,2) = IniVal OCEAN(ng) % tl_W(i,j,k) = IniVal END DO END DO DO itrc=1,NT(ng) DO k=1,N(ng) DO i=Imin,Imax OCEAN(ng) % tl_t(i,j,k,1,itrc) = IniVal OCEAN(ng) % tl_t(i,j,k,2,itrc) = IniVal OCEAN(ng) % tl_t(i,j,k,3,itrc) = IniVal END DO END DO 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 OCEAN(ng) % ad_rubar(i,j,1) = IniVal OCEAN(ng) % ad_rubar(i,j,2) = IniVal OCEAN(ng) % ad_rvbar(i,j,1) = IniVal OCEAN(ng) % ad_rvbar(i,j,2) = IniVal OCEAN(ng) % ad_rzeta(i,j,1) = IniVal OCEAN(ng) % ad_rzeta(i,j,2) = IniVal OCEAN(ng) % ad_ubar(i,j,1) = IniVal OCEAN(ng) % ad_ubar(i,j,2) = IniVal OCEAN(ng) % ad_ubar(i,j,3) = IniVal OCEAN(ng) % ad_vbar(i,j,1) = IniVal OCEAN(ng) % ad_vbar(i,j,2) = IniVal OCEAN(ng) % ad_vbar(i,j,3) = IniVal OCEAN(ng) % ad_zeta(i,j,1) = IniVal OCEAN(ng) % ad_zeta(i,j,2) = IniVal OCEAN(ng) % ad_zeta(i,j,3) = IniVal OCEAN(ng) % ad_ubar_sol(i,j) = IniVal OCEAN(ng) % ad_vbar_sol(i,j) = IniVal OCEAN(ng) % ad_zeta_sol(i,j) = IniVal END DO DO k=1,N(ng) DO i=Imin,Imax OCEAN(ng) % ad_pden(i,j,k) = IniVal OCEAN(ng) % ad_rho(i,j,k) = IniVal OCEAN(ng) % ad_u(i,j,k,1) = IniVal OCEAN(ng) % ad_u(i,j,k,2) = IniVal OCEAN(ng) % ad_v(i,j,k,1) = IniVal OCEAN(ng) % ad_v(i,j,k,2) = IniVal OCEAN(ng) % ad_u_sol(i,j,k) = IniVal OCEAN(ng) % ad_v_sol(i,j,k) = IniVal END DO END DO DO k=0,N(ng) DO i=Imin,Imax OCEAN(ng) % ad_ru(i,j,k,1) = IniVal OCEAN(ng) % ad_ru(i,j,k,2) = IniVal OCEAN(ng) % ad_rv(i,j,k,1) = IniVal OCEAN(ng) % ad_rv(i,j,k,2) = IniVal OCEAN(ng) % ad_W(i,j,k) = IniVal OCEAN(ng) % ad_W_sol(i,j,k) = IniVal OCEAN(ng) % ad_wvel(i,j,k) = IniVal END DO END DO DO itrc=1,NT(ng) DO k=1,N(ng) DO i=Imin,Imax OCEAN(ng) % ad_t(i,j,k,1,itrc) = IniVal OCEAN(ng) % ad_t(i,j,k,2,itrc) = IniVal OCEAN(ng) % ad_t(i,j,k,3,itrc) = IniVal OCEAN(ng) % ad_t_sol(i,j,k,itrc) = IniVal END DO END DO END DO END DO END IF ! ! Working arrays to store adjoint impulse forcing, background error ! covariance, background-error standard deviations, or descent ! conjugate vectors (directions). ! IF (model.eq.0) THEN DO j=Jmin,Jmax DO rec=1,NSA DO i=Imin,Imax OCEAN(ng) % b_ubar(i,j,rec) = IniVal OCEAN(ng) % b_vbar(i,j,rec) = IniVal OCEAN(ng) % b_zeta(i,j,rec) = IniVal OCEAN(ng) % e_ubar(i,j,rec) = IniVal OCEAN(ng) % e_vbar(i,j,rec) = IniVal OCEAN(ng) % e_zeta(i,j,rec) = IniVal END DO END DO DO i=Imin,Imax OCEAN(ng) % d_ubar(i,j) = IniVal OCEAN(ng) % d_vbar(i,j) = IniVal OCEAN(ng) % d_zeta(i,j) = IniVal OCEAN(ng) % f_ubar(i,j) = IniVal OCEAN(ng) % f_vbar(i,j) = IniVal OCEAN(ng) % f_zeta(i,j) = IniVal END DO DO rec=1,NSA DO k=1,N(ng) DO i=Imin,Imax OCEAN(ng) % b_u(i,j,k,rec) = IniVal OCEAN(ng) % b_v(i,j,k,rec) = IniVal OCEAN(ng) % e_u(i,j,k,rec) = IniVal OCEAN(ng) % e_v(i,j,k,rec) = IniVal END DO END DO END DO DO k=1,N(ng) DO i=Imin,Imax OCEAN(ng) % d_u(i,j,k) = IniVal OCEAN(ng) % d_v(i,j,k) = IniVal OCEAN(ng) % f_u(i,j,k) = IniVal OCEAN(ng) % f_v(i,j,k) = IniVal END DO END DO DO itrc=1,NT(ng) DO rec=1,NSA DO k=1,N(ng) DO i=Imin,Imax OCEAN(ng) % b_t(i,j,k,rec,itrc) = IniVal OCEAN(ng) % e_t(i,j,k,rec,itrc) = IniVal END DO END DO END DO DO k=1,N(ng) DO i=Imin,Imax OCEAN(ng) % d_t(i,j,k,itrc) = IniVal OCEAN(ng) % f_t(i,j,k,itrc) = IniVal END DO END DO 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 OCEAN(ng) % ubarG(i,j,1) = IniVal OCEAN(ng) % ubarG(i,j,2) = IniVal OCEAN(ng) % vbarG(i,j,1) = IniVal OCEAN(ng) % vbarG(i,j,2) = IniVal OCEAN(ng) % zetaG(i,j,1) = IniVal OCEAN(ng) % zetaG(i,j,2) = IniVal OCEAN(ng) % f_zetaG(i,j,1) = IniVal OCEAN(ng) % f_zetaG(i,j,2) = IniVal OCEAN(ng) % f_ubarG(i,j,1) = IniVal OCEAN(ng) % f_ubarG(i,j,2) = IniVal OCEAN(ng) % f_vbarG(i,j,1) = IniVal OCEAN(ng) % f_vbarG(i,j,2) = IniVal END DO DO k=1,N(ng) DO i=Imin,Imax OCEAN(ng) % uG(i,j,k,1) = IniVal OCEAN(ng) % uG(i,j,k,2) = IniVal OCEAN(ng) % vG(i,j,k,1) = IniVal OCEAN(ng) % vG(i,j,k,2) = IniVal OCEAN(ng) % f_uG(i,j,k,1) = IniVal OCEAN(ng) % f_uG(i,j,k,2) = IniVal OCEAN(ng) % f_vG(i,j,k,1) = IniVal OCEAN(ng) % f_vG(i,j,k,2) = IniVal END DO END DO DO itrc=1,NT(ng) DO k=1,N(ng) DO i=Imin,Imax OCEAN(ng) % tG(i,j,k,1,itrc) = IniVal OCEAN(ng) % tG(i,j,k,2,itrc) = IniVal OCEAN(ng) % f_tG(i,j,k,1,itrc) = IniVal OCEAN(ng) % f_tG(i,j,k,2,itrc) = IniVal END DO END DO END DO END DO END IF RETURN END SUBROUTINE initialize_ocean END MODULE mod_ocean