SUBROUTINE get_data (ng) ! !git $Id$ !svn $Id: get_data.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 ! !======================================================================= ! ! ! This routine reads in forcing, climatology and other data from ! ! NetCDF files. If there is more than one time-record, data is ! ! loaded into global two-time record arrays. The interpolation ! ! is carried elsewhere. ! ! ! ! Currently, this routine is only executed in serial mode by the ! ! main thread. ! ! ! !======================================================================= ! USE mod_param USE mod_boundary USE mod_clima USE mod_forces USE mod_grid USE mod_iounits USE mod_ncparam USE mod_ocean USE mod_scalars USE mod_sources USE mod_stepping ! USE strings_mod, ONLY : FoundError ! implicit none ! ! Imported variable declarations. ! integer, intent(in) :: ng ! ! Local variable declarations. ! logical :: Lprocess logical, save :: recordless = .FALSE. logical, dimension(3) :: update = (/ .FALSE., .FALSE., .FALSE. /) ! integer :: ILB, IUB, JLB, JUB integer :: LBi, UBi, LBj, UBj integer :: i, ic, my_tile ! character (len=*), parameter :: MyFile = & & "ROMS/Nonlinear/get_data.F" ! ! 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) ! ! Lower and upper bounds for tiled arrays. ! LBi=LBOUND(GRID(ng)%h,DIM=1) UBi=UBOUND(GRID(ng)%h,DIM=1) LBj=LBOUND(GRID(ng)%h,DIM=2) UBj=UBOUND(GRID(ng)%h,DIM=2) ! !----------------------------------------------------------------------- ! Turn on input data time wall clock. !----------------------------------------------------------------------- ! CALL wclock_on (ng, iNLM, 3, 91, MyFile) ! !======================================================================= ! Point Sources/Sinks time dependent data. !======================================================================= ! ! Point Source/Sink vertically integrated mass transport. ! IF (LuvSrc(ng).or.LwSrc(ng)) THEN CALL get_ngfld (ng, iNLM, idRtra, SSF(ng)%ncid, & & 1, SSF(ng), recordless, update(1), & & 1, Nsrc(ng), 1, 2, 1, Nsrc(ng), 1, & & SOURCES(ng) % QbarG) IF (FoundError(exit_flag, NoError, 110, MyFile)) RETURN END IF ! ! Tracer Sources/Sinks. ! DO i=1,NT(ng) IF (LtracerSrc(i,ng)) THEN CALL get_ngfld (ng, iNLM, idRtrc(i), SSF(ng)%ncid, & & 1, SSF(ng), recordless, update(1), & & 1, Nsrc(ng), N(ng), 2, 1, Nsrc(ng), N(ng), & & SOURCES(ng) % TsrcG(:,:,:,i)) IF (FoundError(exit_flag, NoError, 126, MyFile)) RETURN END IF END DO ! !======================================================================= ! Read in forcing data from FORCING NetCDF file(s). !======================================================================= ! ! Set switch to process surface atmospheric fields. ! ! In 4D-Var data assimilation applications, the user have the option ! to fix the prior (background phase) surface fluxes in the successive ! outer loops (Nouter>1) or the final analysis phase. In such case, the ! fluxes are read from the background trajector ! IF (Nrun.eq.1) THEN Lprocess=.TRUE. ELSE Lprocess=.FALSE. END IF ! ! Surface wind components. ! IF (Lprocess) THEN CALL get_2dfld (ng, iNLM, idUair, FRCncid(idUair,ng), & & nFfiles(ng), FRC(1,ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % rmask, & & FORCES(ng) % UwindG) IF (FoundError(exit_flag, NoError, 170, MyFile)) RETURN ! CALL get_2dfld (ng , iNLM, idVair, FRCncid(idVair,ng), & & nFfiles(ng), FRC(1,ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % rmask, & & FORCES(ng) % VwindG) IF (FoundError(exit_flag, NoError, 182, MyFile)) RETURN END IF ! ! Surface air pressure. ! IF (Lprocess) THEN CALL get_2dfld (ng, iNLM, idPair, FRCncid(idPair,ng), & & nFfiles(ng), FRC(1,ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % rmask, & & FORCES(ng) % PairG) IF (FoundError(exit_flag, NoError, 233, MyFile)) RETURN END IF ! ! Surface solar shortwave radiation. ! IF (Lprocess) THEN CALL get_2dfld (ng, iNLM, idSrad, FRCncid(idSrad,ng), & & nFfiles(ng), FRC(1,ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % rmask, & & FORCES(ng) % srflxG) IF (FoundError(exit_flag, NoError, 437, MyFile)) RETURN END IF ! ! Surface downwelling longwave radiation. ! IF (Lprocess) THEN CALL get_2dfld (ng, iNLM, idLdwn, FRCncid(idLdwn,ng), & & nFfiles(ng), FRC(1,ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % rmask, & & FORCES(ng) % lrflxG) IF (FoundError(exit_flag, NoError, 494, MyFile)) RETURN END IF ! ! Surface air temperature. ! IF (Lprocess) THEN CALL get_2dfld (ng, iNLM, idTair, FRCncid(idTair,ng), & & nFfiles(ng), FRC(1,ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % rmask, & & FORCES(ng) % TairG) IF (FoundError(exit_flag, NoError, 515, MyFile)) RETURN END IF ! ! Surface air humidity. ! IF (Lprocess) THEN CALL get_2dfld (ng, iNLM, idQair, FRCncid(idQair,ng), & & nFfiles(ng), FRC(1,ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % rmask, & & FORCES(ng) % HairG) IF (FoundError(exit_flag, NoError, 535, MyFile)) RETURN END IF ! ! Rain fall rate. ! IF (Lprocess) THEN CALL get_2dfld (ng, iNLM, idrain, FRCncid(idrain,ng), & & nFfiles(ng), FRC(1,ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % rmask, & & FORCES(ng) % rainG) IF (FoundError(exit_flag, NoError, 554, MyFile)) RETURN END IF ! !======================================================================= ! In 4D-Var data assimilation algorithms, the user has the option to ! impose the prior (background phase) surface fluxes in the successive ! outer loops (Nouter>1) or the final analysis phase. Such fluxes were ! computed by routine "bulk_fluxes" and stored in the NLM background ! trajectory, BLK structure. Notice that "bulk_fluxes" is called in ! "main3d" only during the prior trajectory computation. ! ! The strategy is to save the 2D surface fluxes in the quicksave ! NetCDF file to allow frequent data records to resolve the daily ! cycle while maintaining its size small. !======================================================================= ! IF (.not.Lprocess) THEN ! ! Get prior surface wind stress components. ! CALL get_2dfld (ng, iNLM, idUsms, BLK(ng)%ncid, & & 1, BLK(ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % umask, & & FORCES(ng) % sustrG) IF (FoundError(exit_flag, NoError, 813, MyFile)) RETURN ! CALL get_2dfld (ng, iNLM, idVsms, BLK(ng)%ncid, & & 1, BLK(ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % vmask, & & FORCES(ng) % svstrG) IF (FoundError(exit_flag, NoError, 825, MyFile)) RETURN ! ! Get prior shortwave radiation flux. For consistency, we need to ! process the same values using in the computation of the net heat ! flux since there is a time interpolation from snapshots. ! CALL get_2dfld (ng, iNLM, idSrad, BLK(ng)%ncid, & & 1, BLK(ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % rmask, & & FORCES(ng) % srflxG) IF (FoundError(exit_flag, NoError, 861, MyFile)) RETURN ! ! Get prior surface net heat flux. ! CALL get_2dfld (ng, iNLM, idTsur(itemp), BLK(ng)%ncid, & & 1, BLK(ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % rmask, & & FORCES(ng) % stfluxG(:,:,:,itemp)) IF (FoundError(exit_flag, NoError, 876, MyFile)) RETURN ! ! Get prior net freshwater flux (E-P). The kinematic net freshwater ! will be computed in "set_vbc" by multiplying with surface salinity. ! CALL get_2dfld (ng, iNLM, idEmPf, BLK(ng)%ncid, & & 1, BLK(ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % rmask, & & FORCES(ng) % stfluxG(:,:,:,isalt)) IF (FoundError(exit_flag, NoError, 893, MyFile)) RETURN END IF ! !======================================================================= ! Read in open boundary conditions from BOUNDARY NetCDF file. !======================================================================= ! ! Free-surface open boundary conditions. ! IF (LprocessOBC(ng)) THEN IF (LBC(iwest,isFsur,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idZbry(iwest), & & BRYncid(idZbry(iwest),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & JLB, JUB, 1, 2, 0, Mm(ng)+1, 1, & & BOUNDARY(ng) % zetaG_west) IF (FoundError(exit_flag, NoError, 918, MyFile)) RETURN END IF ! IF (LBC(ieast,isFsur,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idZbry(ieast), & & BRYncid(idZbry(ieast),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & JLB, JUB, 1, 2, 0, Mm(ng)+1, 1, & & BOUNDARY(ng) % zetaG_east) IF (FoundError(exit_flag, NoError, 931, MyFile)) RETURN END IF ! IF (LBC(isouth,isFsur,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idZbry(isouth), & & BRYncid(idZbry(isouth),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & ILB, IUB, 1, 2, 0, Lm(ng)+1, 1, & & BOUNDARY(ng) % zetaG_south) IF (FoundError(exit_flag, NoError, 944, MyFile)) RETURN END IF ! IF (LBC(inorth,isFsur,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idZbry(inorth), & & BRYncid(idZbry(inorth),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & ILB, IUB, 1, 2, 0, Lm(ng)+1, 1, & & BOUNDARY(ng) % zetaG_north) IF (FoundError(exit_flag, NoError, 957, MyFile)) RETURN END IF END IF ! ! 2D momentum components open boundary conditions. ! IF (LprocessOBC(ng)) THEN IF (LBC(iwest,isUbar,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idU2bc(iwest), & & BRYncid(idU2bc(iwest),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & JLB, JUB, 1, 2, 0, Mm(ng)+1, 1, & & BOUNDARY(ng) % ubarG_west) IF (FoundError(exit_flag, NoError, 977, MyFile)) RETURN END IF ! IF (LBC(iwest,isVbar,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idV2bc(iwest), & & BRYncid(idV2bc(iwest),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & JLB, JUB, 1, 2, 1, Mm(ng)+1, 1, & & BOUNDARY(ng) % vbarG_west) IF (FoundError(exit_flag, NoError, 990, MyFile)) RETURN END IF ! IF (LBC(ieast,isUbar,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idU2bc(ieast), & & BRYncid(idU2bc(ieast),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & JLB, JUB, 1, 2, 0, Mm(ng)+1, 1, & & BOUNDARY(ng) % ubarG_east) IF (FoundError(exit_flag, NoError, 1003, MyFile)) RETURN END IF ! IF (LBC(ieast,isVbar,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idV2bc(ieast), & & BRYncid(idV2bc(ieast),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & JLB, JUB, 1, 2, 1, Mm(ng)+1, 1, & & BOUNDARY(ng) % vbarG_east) IF (FoundError(exit_flag, NoError, 1016, MyFile)) RETURN END IF ! IF (LBC(isouth,isUbar,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idU2bc(isouth), & & BRYncid(idU2bc(isouth),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & ILB, IUB, 1, 2, 1, Lm(ng)+1, 1, & & BOUNDARY(ng) % ubarG_south) IF (FoundError(exit_flag, NoError, 1029, MyFile)) RETURN END IF ! IF (LBC(isouth,isVbar,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idV2bc(isouth), & & BRYncid(idV2bc(isouth),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & ILB, IUB, 1, 2, 0, Lm(ng)+1, 1, & & BOUNDARY(ng) % vbarG_south) IF (FoundError(exit_flag, NoError, 1042, MyFile)) RETURN END IF ! IF (LBC(inorth,isUbar,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idU2bc(inorth), & & BRYncid(idU2bc(inorth),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & ILB, IUB, 1, 2, 1, Lm(ng)+1, 1, & & BOUNDARY(ng) % ubarG_north) IF (FoundError(exit_flag, NoError, 1055, MyFile)) RETURN END IF ! IF (LBC(inorth,isVbar,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idV2bc(inorth), & & BRYncid(idV2bc(inorth),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & ILB, IUB, 1, 2, 0, Lm(ng)+1, 1, & & BOUNDARY(ng) % vbarG_north) IF (FoundError(exit_flag, NoError, 1068, MyFile)) RETURN END IF END IF ! ! 3D momentum components open boundary conditions. ! IF (LprocessOBC(ng)) THEN IF (LBC(iwest,isUvel,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idU3bc(iwest), & & BRYncid(idU3bc(iwest),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & JLB, JUB, N(ng), 2, 0, Mm(ng)+1, N(ng), & & BOUNDARY(ng) % uG_west) IF (FoundError(exit_flag, NoError, 1089, MyFile)) RETURN END IF ! IF (LBC(iwest,isVvel,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idV3bc(iwest), & & BRYncid(idV3bc(iwest),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & JLB, JUB, N(ng), 2, 1, Mm(ng)+1, N(ng), & & BOUNDARY(ng) % vG_west) IF (FoundError(exit_flag, NoError, 1102, MyFile)) RETURN END IF ! IF (LBC(ieast,isUvel,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idU3bc(ieast), & & BRYncid(idU3bc(ieast),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & JLB, JUB, N(ng), 2, 0, Mm(ng)+1, N(ng), & & BOUNDARY(ng) % uG_east) IF (FoundError(exit_flag, NoError, 1115, MyFile)) RETURN END IF ! IF (LBC(ieast,isVvel,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idV3bc(ieast), & & BRYncid(idV3bc(ieast),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & JLB, JUB, N(ng), 2, 1, Mm(ng)+1, N(ng), & & BOUNDARY(ng) % vG_east) IF (FoundError(exit_flag, NoError, 1128, MyFile)) RETURN END IF ! IF (LBC(isouth,isUvel,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idU3bc(isouth), & & BRYncid(idU3bc(isouth),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & ILB, IUB, N(ng), 2, 1, Lm(ng)+1, N(ng), & & BOUNDARY(ng) % uG_south) IF (FoundError(exit_flag, NoError, 1141, MyFile)) RETURN END IF ! IF (LBC(isouth,isVvel,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idV3bc(isouth), & & BRYncid(idV3bc(isouth),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & ILB, IUB, N(ng), 2, 0, Lm(ng)+1, N(ng), & & BOUNDARY(ng) % vG_south) IF (FoundError(exit_flag, NoError, 1154, MyFile)) RETURN END IF ! IF (LBC(inorth,isUvel,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idU3bc(inorth), & & BRYncid(idU3bc(inorth),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & ILB, IUB, N(ng), 2, 1, Lm(ng)+1, N(ng), & & BOUNDARY(ng) % uG_north) IF (FoundError(exit_flag, NoError, 1167, MyFile)) RETURN END IF ! IF (LBC(inorth,isVvel,ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idV3bc(inorth), & & BRYncid(idV3bc(inorth),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & ILB, IUB, N(ng), 2, 0, Lm(ng)+1, N(ng), & & BOUNDARY(ng) % vG_north) IF (FoundError(exit_flag, NoError, 1180, MyFile)) RETURN END IF END IF ! ! Tracer variables open boundary conditions. ! IF (LprocessOBC(ng)) THEN DO i=1,NT(ng) IF (LBC(iwest,isTvar(i),ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idTbry(iwest,i), & & BRYncid(idTbry(iwest,i),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & JLB, JUB, N(ng), 2, 0, Mm(ng)+1, N(ng), & & BOUNDARY(ng) % tG_west(:,:,:,i)) IF (FoundError(exit_flag, NoError, 1201, MyFile)) RETURN END IF END DO ! DO i=1,NT(ng) IF (LBC(ieast,isTvar(i),ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idTbry(ieast,i), & & BRYncid(idTbry(ieast,i),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & JLB, JUB, N(ng), 2, 0, Mm(ng)+1, N(ng), & & BOUNDARY(ng) % tG_east(:,:,:,i)) IF (FoundError(exit_flag, NoError, 1216, MyFile)) RETURN END IF END DO ! DO i=1,NT(ng) IF (LBC(isouth,isTvar(i),ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idTbry(isouth,i), & & BRYncid(idTbry(isouth,i),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & ILB, IUB, N(ng), 2, 0, Lm(ng)+1, N(ng), & & BOUNDARY(ng) % tG_south(:,:,:,i)) IF (FoundError(exit_flag, NoError, 1231, MyFile)) RETURN END IF END DO ! DO i=1,NT(ng) IF (LBC(inorth,isTvar(i),ng)%acquire) THEN CALL get_ngfld (ng, iNLM, idTbry(inorth,i), & & BRYncid(idTbry(inorth,i),ng), & & nBCfiles(ng), BRY(1,ng), & & recordless, update(1), & & ILB, IUB, N(ng), 2, 0, Lm(ng)+1, N(ng), & & BOUNDARY(ng) % tG_north(:,:,:,i)) IF (FoundError(exit_flag, NoError, 1246, MyFile)) RETURN END IF END DO END IF ! !======================================================================= ! Read in data from Climatology NetCDF file. !======================================================================= ! ! Free-surface climatology. ! IF (LsshCLM(ng)) THEN CALL get_2dfld (ng, iNLM, idSSHc, CLMncid(idSSHc,ng), & & nCLMfiles(ng), CLM(1,ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % rmask, & & CLIMA(ng) % sshG) IF (FoundError(exit_flag, NoError, 1272, MyFile)) RETURN END IF ! ! 2D momentum components climatology. ! IF (Lm2CLM(ng)) THEN CALL get_2dfld (ng, iNLM, idUbcl, CLMncid(idUbcl,ng), & & nCLMfiles(ng), CLM(1,ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % umask, & & CLIMA(ng) % ubarclmG) IF (FoundError(exit_flag, NoError, 1290, MyFile)) RETURN ! CALL get_2dfld (ng, iNLM, idVbcl, CLMncid(idVbcl,ng), & & nCLMfiles(ng), CLM(1,ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % vmask, & & CLIMA(ng) % vbarclmG) IF (FoundError(exit_flag, NoError, 1302, MyFile)) RETURN END IF ! ! 3D momentum components climatology. ! IF (Lm3CLM(ng)) THEN CALL get_3dfld (ng, iNLM, idUclm, CLMncid(idUclm,ng), & & nCLMfiles(ng), CLM(1,ng), update(1), & & LBi, UBi, LBj, UBj, 1, N(ng), 2, 1, & & GRID(ng) % umask, & & CLIMA(ng) % uclmG) IF (FoundError(exit_flag, NoError, 1321, MyFile)) RETURN ! CALL get_3dfld (ng, iNLM, idVclm, CLMncid(idVclm,ng), & & nCLMfiles(ng), CLM(1,ng), update(1), & & LBi, UBi, LBj, UBj, 1, N(ng), 2, 1, & & GRID(ng) % vmask, & & CLIMA(ng) % vclmG) IF (FoundError(exit_flag, NoError, 1333, MyFile)) RETURN END IF ! ! Tracers variables climatology. ! ic=0 DO i=1,NT(ng) IF (LtracerCLM(i,ng)) THEN ic=ic+1 CALL get_3dfld (ng, iNLM, idTclm(i), & & CLMncid(idTclm(i),ng), & & nCLMfiles(ng), CLM(1,ng), update(1), & & LBi, UBi, LBj, UBj, 1, N(ng), 2, 1, & & GRID(ng) % rmask, & & CLIMA(ng) % tclmG(:,:,:,:,ic)) IF (FoundError(exit_flag, NoError, 1355, MyFile)) RETURN END IF END DO ! !======================================================================= ! Read weak contraint forcing snapshots. Notice that the forward ! basic state snapshops arrays are reused here. !======================================================================= ! IF (FrequentImpulse(ng)) THEN ! ! Read in free-surface. ! CALL get_2dfld (ng, iNLM, idFsur, TLF(ng)%ncid, & & 1, TLF(ng), update(1), & & LBi, UBi, LBj, UBj, 2, 1, & & GRID(ng) % rmask, & & OCEAN(ng) % zetaG) IF (FoundError(exit_flag, NoError, 1484, MyFile)) RETURN ! ! Read in 3D momentum. ! CALL get_3dfld (ng, iNLM, idUvel, TLF(ng)%ncid, & & 1, TLF(ng), update(1), & & LBi, UBi, LBj, UBj, 1, N(ng), 2, 1, & & GRID(ng) % umask, & & OCEAN(ng) % uG) IF (FoundError(exit_flag, NoError, 1527, MyFile)) RETURN CALL get_3dfld (ng, iNLM, idVvel, TLF(ng)%ncid, & & 1, TLF(ng), update(1), & & LBi, UBi, LBj, UBj, 1, N(ng), 2, 1, & & GRID(ng) % vmask, & & OCEAN(ng) % vG) IF (FoundError(exit_flag, NoError, 1539, MyFile)) RETURN ! ! Read in 3D tracers. ! DO i=1,NT(ng) CALL get_3dfld (ng, iNLM, idTvar(i), TLF(ng)%ncid, & & 1, TLF(ng), update(1), & & LBi, UBi, LBj, UBj, 1, N(ng), 2, 1, & & GRID(ng) % rmask, & & OCEAN(ng) % tG(:,:,:,:,i)) IF (FoundError(exit_flag, NoError, 1554, MyFile)) RETURN END DO END IF ! !----------------------------------------------------------------------- ! Turn off input data time wall clock. !----------------------------------------------------------------------- ! CALL wclock_off (ng, iNLM, 3, 1566, MyFile) ! RETURN END SUBROUTINE get_data