SUBROUTINE tl_set_data (ng, tile) ! !git $Id$ !svn $Id: tl_set_data.F 1151 2023-02-09 03:08:53Z arango $ !================================================== Hernan G. Arango === ! Copyright (c) 2002-2023 The ROMS/TOMS Group ! ! Licensed under a MIT/X style license ! ! See License_ROMS.md ! !======================================================================= ! ! ! This subroutine processes forcing, boundary, climatology, and ! ! other input data. It time-interpolates between snapshots. ! ! ! !======================================================================= ! USE mod_param ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile ! ! Local variable declarations. ! character (len=*), parameter :: MyFile = & & "ROMS/Tangent/tl_set_data.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, iTLM, 4, 32, MyFile) CALL tl_set_data_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS) CALL wclock_off (ng, iTLM, 4, 38, MyFile) ! RETURN END SUBROUTINE tl_set_data ! !*********************************************************************** SUBROUTINE tl_set_data_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS) !*********************************************************************** ! USE mod_param USE mod_boundary USE mod_clima USE mod_coupling USE mod_forces USE mod_grid USE mod_mixing USE mod_ncparam USE mod_ocean USE mod_stepping USE mod_scalars USE mod_sources ! USE analytical_mod USE exchange_2d_mod USE set_2dfld_mod USE set_3dfld_mod USE mp_exchange_mod, ONLY : mp_exchange2d USE mp_exchange_mod, ONLY : mp_exchange3d USE strings_mod, ONLY : FoundError ! implicit none ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile integer, intent(in) :: LBi, UBi, LBj, UBj integer, intent(in) :: IminS, ImaxS, JminS, JmaxS ! ! Local variable declarations. ! logical :: SetBC logical :: update = .FALSE. ! integer :: ILB, IUB, JLB, JUB integer :: i, ic, itrc, j, k, my_tile ! real(r8) :: cff, cff1, cff2 ! character (len=*), parameter :: MyFile = & & "ROMS/Tangent/tl_set_data.F"//", tl_set_data_tile" ! !----------------------------------------------------------------------- ! 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 ! ! Lower and upper bounds for nontiled (global values) boundary arrays. ! my_tile=-1 ! for global values ILB=BOUNDS(ng)%LBi(my_tile) IUB=BOUNDS(ng)%UBi(my_tile) JLB=BOUNDS(ng)%LBj(my_tile) JUB=BOUNDS(ng)%UBj(my_tile) ! !----------------------------------------------------------------------- ! Set kinematic surface solar shortwave radiation flux (degC m/s). !----------------------------------------------------------------------- ! CALL set_2dfld_tile (ng, tile, iTLM, idSrad, & & LBi, UBi, LBj, UBj, & & FORCES(ng)%srflxG, & & FORCES(ng)%srflx, & & update) IF (FoundError(exit_flag, NoError, 194, MyFile)) RETURN ! !----------------------------------------------------------------------- ! Set kinematic surface net heat flux (degC m/s). !----------------------------------------------------------------------- ! CALL set_2dfld_tile (ng, tile, iTLM, idTsur(itemp), & & LBi, UBi, LBj, UBj, & & FORCES(ng)%stfluxG(:,:,:,itemp), & & FORCES(ng)%stflux (:,:,itemp), & & update) IF (FoundError(exit_flag, NoError, 333, MyFile)) RETURN ! !----------------------------------------------------------------------- ! Set kinematic bottom net heat flux (degC m/s). !----------------------------------------------------------------------- ! CALL ana_btflux (ng, tile, iTLM, itemp) ! !----------------------------------------------------------------------- ! Set kinematic surface freshwater (E-P) flux (m/s). !----------------------------------------------------------------------- ! CALL set_2dfld_tile (ng, tile, iTLM, idEmPf, & & LBi, UBi, LBj, UBj, & & FORCES(ng)%stfluxG(:,:,:,isalt), & & FORCES(ng)%stflux (:,:,isalt), & & update) IF (FoundError(exit_flag, NoError, 407, MyFile)) RETURN ! !----------------------------------------------------------------------- ! Set kinematic bottom salt flux (m/s). !----------------------------------------------------------------------- ! CALL ana_btflux (ng, tile, iTLM, isalt) ! !----------------------------------------------------------------------- ! Set kinematic surface momentum flux (m2/s2). !----------------------------------------------------------------------- ! CALL set_2dfld_tile (ng, tile, iTLM, idUsms, & & LBi, UBi, LBj, UBj, & & FORCES(ng)%sustrG, & & FORCES(ng)%sustr, & & update) IF (FoundError(exit_flag, NoError, 492, MyFile)) RETURN ! CALL set_2dfld_tile (ng, tile, iTLM, idVsms, & & LBi, UBi, LBj, UBj, & & FORCES(ng)%svstrG, & & FORCES(ng)%svstr, & & update) IF (FoundError(exit_flag, NoError, 499, MyFile)) RETURN ! ! If input point wind stress, rotate to curvilinear grid. Notice ! that rotation is done at RHO-points. It does not matter. ! IF (.not.Linfo(1,idUsms,ng).or. & & (Iinfo(5,idUsms,ng).ne.Lm(ng)+1).or. & & (Iinfo(6,idUsms,ng).ne.Mm(ng)+2)) THEN DO j=JstrR,JendR DO i=IstrR,IendR cff1=FORCES(ng)%sustr(i,j)*GRID(ng)%CosAngler(i,j)+ & & FORCES(ng)%svstr(i,j)*GRID(ng)%SinAngler(i,j) cff2=FORCES(ng)%svstr(i,j)*GRID(ng)%CosAngler(i,j)- & & FORCES(ng)%sustr(i,j)*GRID(ng)%SinAngler(i,j) FORCES(ng)%sustr(i,j)=cff1 FORCES(ng)%svstr(i,j)=cff2 END DO END DO IF (EWperiodic(ng).or.NSperiodic(ng)) THEN CALL exchange_u2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & FORCES(ng)%sustr) CALL exchange_v2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & FORCES(ng)%svstr) END IF CALL mp_exchange2d (ng, tile, iTLM, 2, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & FORCES(ng)%sustr, & & FORCES(ng)%svstr) END IF ! !----------------------------------------------------------------------- ! Set point Sources/Sinks (river runoff). !----------------------------------------------------------------------- ! IF (DOMAIN(ng)%SouthWest_Test(tile)) THEN IF (LuvSrc(ng).or.LwSrc(ng)) THEN CALL set_ngfld (ng, iTLM, idRtra, 1, Nsrc(ng), 1, & & 1, Nsrc(ng), 1, & & SOURCES(ng) % QbarG, & & SOURCES(ng) % Qbar, & & update) IF (FoundError(exit_flag, NoError, 719, MyFile)) RETURN DO k=1,N(ng) DO i=1,Nsrc(ng) SOURCES(ng)%Qsrc(i,k)=SOURCES(ng)%Qbar(i)* & & SOURCES(ng)%Qshape(i,k) END DO END DO END IF DO itrc=1,NT(ng) IF (LtracerSrc(itrc,ng)) THEN CALL set_ngfld (ng, iTLM, idRtrc(itrc), 1, Nsrc(ng), N(ng), & & 1, Nsrc(ng), N(ng), & & SOURCES(ng) % TsrcG(:,:,:,itrc), & & SOURCES(ng) % Tsrc(:,:,itrc), & & update) IF (FoundError(exit_flag, NoError, 739, MyFile)) RETURN END IF END DO END IF ! !----------------------------------------------------------------------- ! Set open boundary conditions fields. !----------------------------------------------------------------------- ! ! Free-surface. ! IF (LprocessOBC(ng)) THEN IF (DOMAIN(ng)%SouthWest_Test(tile)) THEN IF (tl_LBC(iwest,isFsur,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idZbry(iwest), JLB, JUB, 1, & & 0, Mm(ng)+1, 1, & & BOUNDARY(ng) % zetaG_west, & & BOUNDARY(ng) % zeta_west, & & update) IF (FoundError(exit_flag, NoError, 763, MyFile)) RETURN END IF IF (tl_LBC(ieast,isFsur,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idZbry(ieast), JLB, JUB, 1, & & 0, Mm(ng)+1, 1, & & BOUNDARY(ng) % zetaG_east, & & BOUNDARY(ng) % zeta_east, & & update) IF (FoundError(exit_flag, NoError, 772, MyFile)) RETURN END IF IF (tl_LBC(isouth,isFsur,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idZbry(isouth), ILB, IUB, 1, & & 0, Lm(ng)+1 ,1, & & BOUNDARY(ng) % zetaG_south, & & BOUNDARY(ng) % zeta_south, & & update) IF (FoundError(exit_flag, NoError, 781, MyFile)) RETURN END IF IF (tl_LBC(inorth,isFsur,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idZbry(inorth), ILB, IUB, 1, & & 0, Lm(ng)+1, 1, & & BOUNDARY(ng) % zetaG_north, & & BOUNDARY(ng) % zeta_north, & & update) IF (FoundError(exit_flag, NoError, 790, MyFile)) RETURN END IF END IF END IF ! ! 2D momentum. ! IF (LprocessOBC(ng)) THEN IF (DOMAIN(ng)%SouthWest_Test(tile)) THEN IF (tl_LBC(iwest,isUbar,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idU2bc(iwest), JLB, JUB, 1, & & 0, Mm(ng)+1, 1, & & BOUNDARY(ng) % ubarG_west, & & BOUNDARY(ng) % ubar_west, & & update) IF (FoundError(exit_flag, NoError, 886, MyFile)) RETURN END IF IF (tl_LBC(iwest,isVbar,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idV2bc(iwest), JLB, JUB, 1, & & 1, Mm(ng)+1, 1, & & BOUNDARY(ng) % vbarG_west, & & BOUNDARY(ng) % vbar_west, & & update) IF (FoundError(exit_flag, NoError, 895, MyFile)) RETURN END IF IF (tl_LBC(ieast,isUbar,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idU2bc(ieast), JLB, JUB, 1, & & 0, Mm(ng)+1, 1, & & BOUNDARY(ng) % ubarG_east, & & BOUNDARY(ng) % ubar_east, & & update) IF (FoundError(exit_flag, NoError, 904, MyFile)) RETURN END IF IF (tl_LBC(ieast,isVbar,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idV2bc(ieast), JLB, JUB, 1, & & 1, Mm(ng)+1, 1, & & BOUNDARY(ng) % vbarG_east, & & BOUNDARY(ng) % vbar_east, & & update) IF (FoundError(exit_flag, NoError, 913, MyFile)) RETURN END IF IF (tl_LBC(isouth,isUbar,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idU2bc(isouth), ILB, IUB, 1, & & 1, Lm(ng)+1, 1, & & BOUNDARY(ng) % ubarG_south, & & BOUNDARY(ng) % ubar_south, & & update) IF (FoundError(exit_flag, NoError, 922, MyFile)) RETURN END IF IF (tl_LBC(isouth,isVbar,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idV2bc(isouth), ILB, IUB, 1, & & 0, Lm(ng)+1, 1, & & BOUNDARY(ng) % vbarG_south, & & BOUNDARY(ng) % vbar_south, & & update) IF (FoundError(exit_flag, NoError, 931, MyFile)) RETURN END IF IF (tl_LBC(inorth,isUbar,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idU2bc(inorth), ILB, IUB, 1, & & 1, Lm(ng)+1, 1, & & BOUNDARY(ng) % ubarG_north, & & BOUNDARY(ng) % ubar_north, & & update) IF (FoundError(exit_flag, NoError, 940, MyFile)) RETURN END IF IF (tl_LBC(inorth,isVbar,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idV2bc(inorth), ILB, IUB, 1, & & 0, Lm(ng)+1, 1, & & BOUNDARY(ng) % vbarG_north, & & BOUNDARY(ng) % vbar_north, & & update) IF (FoundError(exit_flag, NoError, 949, MyFile)) RETURN END IF END IF END IF ! ! 3D momentum. ! IF (LprocessOBC(ng)) THEN IF (DOMAIN(ng)%SouthWest_Test(tile)) THEN IF (tl_LBC(iwest,isUvel,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idU3bc(iwest), JLB, JUB, N(ng), & & 0, Mm(ng)+1, N(ng), & & BOUNDARY(ng) % uG_west, & & BOUNDARY(ng) % u_west, & & update) IF (FoundError(exit_flag, NoError, 970, MyFile)) RETURN END IF IF (tl_LBC(iwest,isVvel,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idV3bc(iwest), JLB, JUB, N(ng), & & 1, Mm(ng)+1, N(ng), & & BOUNDARY(ng) % vG_west, & & BOUNDARY(ng) % v_west, & & update) IF (FoundError(exit_flag, NoError, 979, MyFile)) RETURN END IF IF (tl_LBC(ieast,isUvel,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idU3bc(ieast), JLB, JUB, N(ng), & & 0, Mm(ng)+1, N(ng), & & BOUNDARY(ng) % uG_east, & & BOUNDARY(ng) % u_east, & & update) IF (FoundError(exit_flag, NoError, 988, MyFile)) RETURN END IF IF (tl_LBC(ieast,isVvel,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idV3bc(ieast), JLB, JUB, N(ng), & & 1, Mm(ng)+1, N(ng), & & BOUNDARY(ng) % vG_east, & & BOUNDARY(ng) % v_east, & & update) IF (FoundError(exit_flag, NoError, 997, MyFile)) RETURN END IF IF (tl_LBC(isouth,isUvel,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idU3bc(isouth), ILB, IUB, N(ng), & & 1, Lm(ng)+1, N(ng), & & BOUNDARY(ng) % uG_south, & & BOUNDARY(ng) % u_south, & & update) IF (FoundError(exit_flag, NoError, 1006, MyFile)) RETURN END IF IF (tl_LBC(isouth,isVvel,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idV3bc(isouth), ILB, IUB, N(ng), & & 0, Lm(ng)+1, N(ng), & & BOUNDARY(ng) % vG_south, & & BOUNDARY(ng) % v_south, & & update) IF (FoundError(exit_flag, NoError, 1015, MyFile)) RETURN END IF IF (tl_LBC(inorth,isUvel,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idU3bc(inorth), ILB, IUB, N(ng), & & 1, Lm(ng)+1, N(ng), & & BOUNDARY(ng) % uG_north, & & BOUNDARY(ng) % u_north, & & update) IF (FoundError(exit_flag, NoError, 1024, MyFile)) RETURN END IF IF (tl_LBC(inorth,isVvel,ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idV3bc(inorth), ILB, IUB, N(ng), & & 0, Lm(ng)+1, N(ng), & & BOUNDARY(ng) % vG_north, & & BOUNDARY(ng) % v_north, & & update) IF (FoundError(exit_flag, NoError, 1033, MyFile)) RETURN END IF END IF END IF ! ! Tracers. ! IF (LprocessOBC(ng)) THEN IF (DOMAIN(ng)%SouthWest_Test(tile)) THEN DO itrc=1,NT(ng) IF (tl_LBC(iwest,isTvar(itrc),ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idTbry(iwest,itrc), & & JLB, JUB, N(ng), 0, Mm(ng)+1, N(ng), & & BOUNDARY(ng) % tG_west(:,:,:,itrc), & & BOUNDARY(ng) % t_west(:,:,itrc), & & update) IF (FoundError(exit_flag, NoError, & & 1054, MyFile)) RETURN END IF IF (tl_LBC(ieast,isTvar(itrc),ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idTbry(ieast,itrc), & & JLB, JUB, N(ng), 0, Mm(ng)+1, N(ng), & & BOUNDARY(ng) % tG_east(:,:,:,itrc), & & BOUNDARY(ng) % t_east(:,:,itrc), & & update) IF (FoundError(exit_flag, NoError, & & 1064, MyFile)) RETURN END IF IF (tl_LBC(isouth,isTvar(itrc),ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idTbry(isouth,itrc), & & ILB, IUB, N(ng), 0, Lm(ng)+1, N(ng), & & BOUNDARY(ng) % tG_south(:,:,:,itrc), & & BOUNDARY(ng) % t_south(:,:,itrc), & & update) IF (FoundError(exit_flag, NoError, & & 1074, MyFile)) RETURN END IF IF (tl_LBC(inorth,isTvar(itrc),ng)%acquire) THEN CALL set_ngfld (ng, iTLM, idTbry(inorth,itrc), & & ILB, IUB, N(ng), 0, Lm(ng)+1, N(ng), & & BOUNDARY(ng) % tG_north(:,:,:,itrc), & & BOUNDARY(ng) % t_north(:,:,itrc), & & update) IF (FoundError(exit_flag, NoError, & & 1084, MyFile)) RETURN END IF END DO END IF END IF ! !----------------------------------------------------------------------- ! Set sea surface height climatology (m). !----------------------------------------------------------------------- ! IF (LsshCLM(ng)) THEN CALL set_2dfld_tile (ng, tile, iTLM, idSSHc, & & LBi, UBi, LBj, UBj, & & CLIMA(ng)%sshG, & & CLIMA(ng)%ssh, & & update) IF (FoundError(exit_flag, NoError, 1105, MyFile)) RETURN END IF ! !----------------------------------------------------------------------- ! Set 2D momentum climatology (m/s). !----------------------------------------------------------------------- ! IF (Lm2CLM(ng)) THEN CALL set_2dfld_tile (ng, tile, iTLM, idUbcl, & & LBi, UBi, LBj, UBj, & & CLIMA(ng)%ubarclmG, & & CLIMA(ng)%ubarclm, & & update) IF (FoundError(exit_flag, NoError, 1122, MyFile)) RETURN ! CALL set_2dfld_tile (ng, tile, iTLM, idVbcl, & & LBi, UBi, LBj, UBj, & & CLIMA(ng)%vbarclmG, & & CLIMA(ng)%vbarclm, & & update) IF (FoundError(exit_flag, NoError, 1129, MyFile)) RETURN END IF ! !----------------------------------------------------------------------- ! Set 3D momentum climatology (m/s). !----------------------------------------------------------------------- ! IF (Lm3CLM(ng)) THEN CALL set_3dfld_tile (ng, tile, iTLM, idUclm, & & LBi, UBi, LBj, UBj, 1, N(ng), & & CLIMA(ng)%uclmG, & & CLIMA(ng)%uclm, & & update) IF (FoundError(exit_flag, NoError, 1148, MyFile)) RETURN ! CALL set_3dfld_tile (ng, tile, iTLM, idVclm, & & LBi, UBi, LBj, UBj, 1, N(ng), & & CLIMA(ng)%vclmG, & & CLIMA(ng)%vclm, & & update) IF (FoundError(exit_flag, NoError, 1155, MyFile)) RETURN END IF ! !----------------------------------------------------------------------- ! Set tracer climatology. !----------------------------------------------------------------------- ! ic=0 DO itrc=1,NT(ng) IF (LtracerCLM(itrc,ng)) THEN ic=ic+1 CALL set_3dfld_tile (ng, tile, iTLM, idTclm(itrc), & & LBi, UBi, LBj, UBj, 1, N(ng), & & CLIMA(ng)%tclmG(:,:,:,:,ic), & & CLIMA(ng)%tclm (:,:,:,ic), & & update) IF (FoundError(exit_flag, NoError, 1177, MyFile)) RETURN END IF END DO ! !----------------------------------------------------------------------- ! Set forward solution needed by Tangent/Adjoint models. !----------------------------------------------------------------------- ! ! Set forward free-surface. ! DO k=1,3 CALL set_2dfld_tile (ng, tile, iTLM, idFsur, & & LBi, UBi, LBj, UBj, & & OCEAN(ng)%zetaG, & & OCEAN(ng)%zeta(:,:,k), & & update) IF (FoundError(exit_flag, NoError, 1197, MyFile)) RETURN END DO DO j=JstrR,JendR DO i=IstrR,IendR COUPLING(ng)%Zt_avg1(i,j)=OCEAN(ng)%zeta(i,j,1) END DO END DO ! ! Set forward 2D momentum. ! DO k=1,3 CALL set_2dfld_tile (ng, tile, iTLM, idUbar, & & LBi, UBi, LBj, UBj, & & OCEAN(ng)%ubarG, & & OCEAN(ng)%ubar(:,:,k), & & update) IF (FoundError(exit_flag, NoError, 1216, MyFile)) RETURN CALL set_2dfld_tile (ng, tile, iTLM, idVbar, & & LBi, UBi, LBj, UBj, & & OCEAN(ng)%vbarG, & & OCEAN(ng)%vbar(:,:,k), & & update) IF (FoundError(exit_flag, NoError, 1223, MyFile)) RETURN END DO ! ! Set forward time-averaged barotropic fluxes. ! CALL set_2dfld_tile (ng, tile, iTLM, idUfx1, & & LBi, UBi, LBj, UBj, & & COUPLING(ng)%DU_avg1G, & & COUPLING(ng)%DU_avg1, & & update) IF (FoundError(exit_flag, NoError, 1263, MyFile)) RETURN ! CALL set_2dfld_tile (ng, tile, iTLM, idUfx2, & & LBi, UBi, LBj, UBj, & & COUPLING(ng)%DU_avg2G, & & COUPLING(ng)%DU_avg2, & & update) IF (FoundError(exit_flag, NoError, 1270, MyFile)) RETURN ! CALL set_2dfld_tile (ng, tile, iTLM, idVfx1, & & LBi, UBi, LBj, UBj, & & COUPLING(ng)%DV_avg1G, & & COUPLING(ng)%DV_avg1, & & update) IF (FoundError(exit_flag, NoError, 1277, MyFile)) RETURN ! CALL set_2dfld_tile (ng, tile, iTLM, idVfx2, & & LBi, UBi, LBj, UBj, & & COUPLING(ng)%DV_avg2G, & & COUPLING(ng)%DV_avg2, & & update) IF (FoundError(exit_flag, NoError, 1284, MyFile)) RETURN ! ! Set 3D momentum. ! DO k=1,2 CALL set_3dfld_tile (ng, tile, iTLM, idUvel, & & LBi, UBi, LBj, UBj, 1, N(ng), & & OCEAN(ng)%uG, & & OCEAN(ng)%u(:,:,:,k), & & update) IF (FoundError(exit_flag, NoError, 1294, MyFile)) RETURN ! CALL set_3dfld_tile (ng, tile, iTLM, idVvel, & & LBi, UBi, LBj, UBj, 1, N(ng), & & OCEAN(ng)%vG, & & OCEAN(ng)%v(:,:,:,k), & & update) IF (FoundError(exit_flag, NoError, 1301, MyFile)) RETURN END DO ! ! Set 3D tracers. ! DO itrc=1,NT(ng) DO k=1,3 CALL set_3dfld_tile (ng, tile, iTLM, idTvar(itrc), & & LBi, UBi, LBj, UBj, 1, N(ng), & & OCEAN(ng)%tG(:,:,:,:,itrc), & & OCEAN(ng)%t(:,:,:,k,itrc), & & update) IF (FoundError(exit_flag, NoError, 1334, MyFile)) RETURN END DO END DO ! ! Set forward vertical mixing variables. ! DO itrc=1,NAT CALL set_3dfld_tile (ng, tile, iTLM, idDiff(itrc), & & LBi, UBi, LBj, UBj, 0, N(ng), & & MIXING(ng)%AktG(:,:,:,:,itrc), & & MIXING(ng)%Akt(:,:,:,itrc), & & update) IF (FoundError(exit_flag, NoError, 1348, MyFile)) RETURN END DO ! CALL set_3dfld_tile (ng, tile, iTLM, idVvis, & & LBi, UBi, LBj, UBj, 0, N(ng), & & MIXING(ng)%AkvG, & & MIXING(ng)%Akv, & & update) IF (FoundError(exit_flag, NoError, 1356, MyFile)) RETURN ! !----------------------------------------------------------------------- ! Set weak contraint frequent impulse forcing. !----------------------------------------------------------------------- ! IF (FrequentImpulse(ng)) THEN ! ! Set free-surface forcing. ! CALL set_2dfld_tile (ng, tile, iTLM, idZtlf, & & LBi, UBi, LBj, UBj, & & OCEAN(ng)%f_zetaG, & & OCEAN(ng)%f_zeta, & & update) IF (FoundError(exit_flag, NoError, 1468, MyFile)) RETURN ! ! Set 3D momentum forcing. ! CALL set_3dfld_tile (ng, tile, iTLM, idUtlf, & & LBi, UBi, LBj, UBj, 1, N(ng), & & OCEAN(ng)%f_uG, & & OCEAN(ng)%f_u, & & update) IF (FoundError(exit_flag, NoError, 1496, MyFile)) RETURN ! CALL set_3dfld_tile (ng, tile, iTLM, idVtlf, & & LBi, UBi, LBj, UBj, 1, N(ng), & & OCEAN(ng)%f_vG, & & OCEAN(ng)%f_v, & & update) IF (FoundError(exit_flag, NoError, 1503, MyFile)) RETURN ! ! Set 3D tracers forcing. ! DO itrc=1,NT(ng) CALL set_3dfld_tile (ng, tile, iTLM, idTtlf(itrc), & & LBi, UBi, LBj, UBj, 1, N(ng), & & OCEAN(ng)%f_tG(:,:,:,:,itrc), & & OCEAN(ng)%f_t(:,:,:,itrc), & & update) IF (FoundError(exit_flag, NoError, 1513, MyFile)) RETURN END DO END IF ! RETURN END SUBROUTINE tl_set_data_tile