SUBROUTINE set_data (ng, tile) ! !svn $Id: set_data.F 413 2009-12-07 19:06:32Z arango $ !================================================== Hernan G. Arango === ! Copyright (c) 2002-2009 The ROMS/TOMS Group ! ! Licensed under a MIT/X style license ! ! See License_ROMS.txt ! !======================================================================= ! ! ! This subroutine processes forcing, boundary, climatology, and ! ! assimilation input data. It time-interpolates between snapshots. ! ! ! !======================================================================= ! USE mod_param ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile ! ! Local variable declarations. ! 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)- and MAX(I,J)-directions. ! LBij=BOUNDS(ng)%LBij UBij=BOUNDS(ng)%UBij ! CALL wclock_on (ng, iNLM, 4) CALL set_data_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS) CALL wclock_off (ng, iNLM, 4) RETURN END SUBROUTINE set_data ! !*********************************************************************** SUBROUTINE set_data_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS) !*********************************************************************** ! USE mod_param USE mod_boundary 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 distribute_mod, ONLY : mp_boundary USE mp_exchange_mod, ONLY : mp_exchange2d USE mp_exchange_mod, ONLY : mp_exchange3d USE set_2dfld_mod USE set_3dfld_mod ! 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 :: EWperiodic=.FALSE. logical :: NSperiodic=.FALSE. logical :: update = .FALSE. logical :: bry_update integer :: ILB, IUB, JLB, JUB integer :: i, itrc, j, k, order real(r8) :: Zr, cff, cff1, cff2 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: work1 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: work2 ! !----------------------------------------------------------------------- ! 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, IstrR, IstrT, IstrU, Iend, IendR, IendT integer :: Jstr, JstrR, JstrT, JstrV, Jend, JendR, JendT ! Istr =BOUNDS(ng)%Istr (tile) IstrR=BOUNDS(ng)%IstrR(tile) IstrT=BOUNDS(ng)%IstrT(tile) IstrU=BOUNDS(ng)%IstrU(tile) Iend =BOUNDS(ng)%Iend (tile) IendR=BOUNDS(ng)%IendR(tile) IendT=BOUNDS(ng)%IendT(tile) Jstr =BOUNDS(ng)%Jstr (tile) JstrR=BOUNDS(ng)%JstrR(tile) JstrT=BOUNDS(ng)%JstrT(tile) JstrV=BOUNDS(ng)%JstrV(tile) Jend =BOUNDS(ng)%Jend (tile) JendR=BOUNDS(ng)%JendR(tile) JendT=BOUNDS(ng)%JendT(tile) ! !----------------------------------------------------------------------- ! Set cloud fraction (nondimensional). Notice that clouds are ! processed first in case that they are used to adjust shortwave ! radiation. !----------------------------------------------------------------------- ! CALL set_2dfld_tile (ng, tile, iNLM, idCfra, & & LBi, UBi, LBj, UBj, & & FORCES(ng)%cloudG, & & FORCES(ng)%cloud, & & update) ! !----------------------------------------------------------------------- ! Set surface air temperature (degC). !----------------------------------------------------------------------- ! CALL set_2dfld_tile (ng, tile, iNLM, idTair, & & LBi, UBi, LBj, UBj, & & FORCES(ng)%TairG, & & FORCES(ng)%Tair, & & update) ! !----------------------------------------------------------------------- ! Set surface air relative or specific humidity. !----------------------------------------------------------------------- ! CALL set_2dfld_tile (ng, tile, iNLM, idQair, & & LBi, UBi, LBj, UBj, & & FORCES(ng)%HairG, & & FORCES(ng)%Hair, & & update) ! !----------------------------------------------------------------------- ! Set kinematic surface solar shortwave radiation flux (degC m/s). !----------------------------------------------------------------------- ! CALL ana_srflux (ng, tile, iNLM) ! ! Modulate the averaged shortwave radiation flux by the local diurnal ! cycle. ! CALL ana_srflux (ng, tile, iNLM) ! !----------------------------------------------------------------------- ! Set surface winds (m/s). !----------------------------------------------------------------------- ! CALL set_2dfld_tile (ng, tile, iNLM, idUair, & & LBi, UBi, LBj, UBj, & & FORCES(ng)%UwindG, & & FORCES(ng)%Uwind, & & update) CALL set_2dfld_tile (ng, tile, iNLM, idVair, & & LBi, UBi, LBj, UBj, & & FORCES(ng)%VwindG, & & FORCES(ng)%Vwind, & & update) ! ! If input point surface winds or interpolated from coarse data, rotate ! to curvilinear grid. ! IF (.not.Linfo(1,idUair,ng).or. & & (Iinfo(5,idUair,ng).ne.Lm(ng)+2).or. & & (Iinfo(6,idUair,ng).ne.Mm(ng)+2)) THEN DO j=JstrR,JendR DO i=IstrR,IendR cff1=FORCES(ng)%Uwind(i,j)*GRID(ng)%CosAngler(i,j)+ & & FORCES(ng)%Vwind(i,j)*GRID(ng)%SinAngler(i,j) cff2=FORCES(ng)%Vwind(i,j)*GRID(ng)%CosAngler(i,j)- & & FORCES(ng)%Uwind(i,j)*GRID(ng)%SinAngler(i,j) FORCES(ng)%Uwind(i,j)=cff1 FORCES(ng)%Vwind(i,j)=cff2 END DO END DO CALL mp_exchange2d (ng, tile, iNLM, 2, & & LBi, UBi, LBj, UBj, & & NghostPoints, EWperiodic, NSperiodic, & & FORCES(ng)%UWind, & & FORCES(ng)%VWind) END IF ! !----------------------------------------------------------------------- ! Set rain fall rate (kg/m2/s). !----------------------------------------------------------------------- ! CALL ana_rain (ng, tile, iNLM) ! !----------------------------------------------------------------------- ! Set kinematic bottom net heat flux (degC m/s). !----------------------------------------------------------------------- ! CALL ana_btflux (ng, tile, iNLM, itemp) ! !----------------------------------------------------------------------- ! Set kinematic surface freshwater (E-P) flux (m/s). !----------------------------------------------------------------------- ! CALL ana_stflux (ng, tile, iNLM, isalt) ! !----------------------------------------------------------------------- ! Set kinematic bottom salt flux (m/s). !----------------------------------------------------------------------- ! CALL ana_btflux (ng, tile, iNLM, isalt) ! !----------------------------------------------------------------------- ! Set kinematic surface and bottom pasive tracer fluxes (T m/s). !----------------------------------------------------------------------- ! DO itrc=NAT+1,NT(ng) CALL ana_stflux (ng, tile, iNLM, itrc) CALL ana_btflux (ng, tile, iNLM, itrc) END DO ! !----------------------------------------------------------------------- ! Set surface air pressure (mb). !----------------------------------------------------------------------- ! CALL set_2dfld_tile (ng, tile, iNLM, idPair, & & LBi, UBi, LBj, UBj, & & FORCES(ng)%PairG, & & FORCES(ng)%Pair, & & update) ! !----------------------------------------------------------------------- ! Set surface wind-induced wave amplitude, direction and period. !----------------------------------------------------------------------- ! CALL ana_wwave (ng, tile, iNLM) ! !----------------------------------------------------------------------- ! Set point Sources/Sinks (river runoff). !----------------------------------------------------------------------- ! IF (.TRUE.) THEN CALL set_ngfld (ng, iNLM, idRtra, 1, Nsrc(ng), 1, & & 1, Nsrc(ng), 1, & & SOURCES(ng) % QbarG, & & SOURCES(ng) % Qbar, & & update) 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 DO itrc=1,NT(ng) IF (LtracerSrc(itrc,ng)) THEN CALL set_ngfld (ng, iNLM, idRtrc(itrc), 1, Nsrc(ng), N(ng), & & 1, Nsrc(ng), N(ng), & & SOURCES(ng) % TsrcG(1,1,1,itrc), & & SOURCES(ng) % Tsrc(1,1,itrc), & & update) END IF END DO END IF ! !----------------------------------------------------------------------- ! Set open boundary conditions fields. !----------------------------------------------------------------------- ! ! Lower and upper bounds for nontiled arrays. ! ILB=0 IUB=Im(ng)+1 JLB=0 JUB=Jm(ng)+1 CALL ana_fsobc (ng, tile, iNLM) ! ! Ensure that water level on boundary cells is above bed elevation. ! bry_update=.FALSE. IF (Istr.eq.1) THEN DO j=JstrR,JendR cff=Dcrit(ng)-GRID(ng)%h(0,j) IF (BOUNDARY(ng)%zeta_west(j).le.cff) THEN BOUNDARY(ng)%zeta_west(j)=cff END IF END DO bry_update=.TRUE. END IF CALL mp_boundary (ng, iNLM, JstrR, JendR, JLB, JUB, 1, 1, & & bry_update, & & BOUNDARY(ng)%zeta_west) bry_update=.FALSE. IF (Iend.eq.Lm(ng)) THEN DO j=JstrR,JendR cff=Dcrit(ng)-GRID(ng)%h(Lm(ng)+1,j) IF (BOUNDARY(ng)%zeta_east(j).le.cff) THEN BOUNDARY(ng)%zeta_east(j)=cff END IF END DO bry_update=.TRUE. END IF CALL mp_boundary (ng, iNLM, JstrR, JendR, JLB, JUB, 1, 1, & & bry_update, & & BOUNDARY(ng)%zeta_east) bry_update=.FALSE. IF (Jstr.eq.1) THEN DO i=IstrR,IendR cff=Dcrit(ng)-GRID(ng)%h(i,0) IF (BOUNDARY(ng)%zeta_south(i).le.cff) THEN BOUNDARY(ng)%zeta_south(i)=cff END IF END DO bry_update=.TRUE. END IF CALL mp_boundary (ng, iNLM, IstrR, IendR, ILB, IUB, 1, 1, & & bry_update, & & BOUNDARY(ng)%zeta_south) CALL ana_m2obc (ng, tile, iNLM) IF (.TRUE.) THEN END IF IF (.TRUE.) THEN DO itrc=1,NT(ng) END DO END IF RETURN END SUBROUTINE set_data_tile