SUBROUTINE main3d (ng) ! !svn $Id: main3d.F 409 2009-11-18 20:39:05Z arango $ !======================================================================= ! Copyright (c) 2002-2009 The ROMS/TOMS Group ! ! Licensed under a MIT/X style license ! ! See License_ROMS.txt Hernan G. Arango ! !========================================== Alexander F. Shchepetkin === ! ! ! This subroutine is the main driver for nonlinear ROMS/TOMS when ! ! configurated as a full 3D baroclinic ocean model. It advances ! ! forward the primitive equations for a single time step. ! ! ! !======================================================================= ! USE mod_param USE mod_parallel USE mod_iounits USE mod_scalars USE mod_stepping ! USE bbl_mod, ONLY : bblm USE bulk_flux_mod, ONLY : bulk_flux USE diag_mod, ONLY : diag USE gls_corstep_mod, ONLY : gls_corstep USE gls_prestep_mod, ONLY : gls_prestep USE ini_fields_mod, ONLY : ini_fields, ini_zeta USE omega_mod, ONLY : omega USE rho_eos_mod, ONLY : rho_eos USE rhs3d_mod, ONLY : rhs3d USE sediment_mod, ONLY : sediment USE set_depth_mod, ONLY : set_depth USE set_massflux_mod, ONLY : set_massflux USE set_tides_mod, ONLY : set_tides USE set_vbc_mod, ONLY : set_vbc USE set_zeta_mod, ONLY : set_zeta USE step2d_mod, ONLY : step2d USE step3d_t_mod, ONLY : step3d_t USE step3d_uv_mod, ONLY : step3d_uv USE step_floats_mod, ONLY : step_floats USE wvelocity_mod, ONLY : wvelocity ! implicit none ! ! Imported variable declarations. ! integer, intent(in) :: ng ! ! Local variable declarations. ! integer :: my_iif, next_indx1, subs, tile, thread integer :: Lend, Lstr, chunk_size ! !======================================================================= ! Time-step nonlinear 3D primitive equations. !======================================================================= ! ! Set time indices and time clock. ! nstp(ng)=1+MOD(iic(ng)-ntstart(ng),2) nnew(ng)=3-nstp(ng) nrhs(ng)=nstp(ng) time(ng)=time(ng)+dt(ng) tdays(ng)=time(ng)*sec2day CALL time_string (time(ng), time_code(ng)) ! !----------------------------------------------------------------------- ! Read in required data, if any, from input NetCDF files. !----------------------------------------------------------------------- ! CALL get_data (ng) IF (exit_flag.ne.NoError) RETURN ! !----------------------------------------------------------------------- ! If applicable, process input data: time interpolate between data ! snapshots. !----------------------------------------------------------------------- ! DO thread=0,numthreads-1 subs=NtileX(ng)*NtileE(ng)/numthreads DO tile=subs*thread,subs*(thread+1)-1,+1 CALL set_data (ng, MyRank) END DO END DO IF (exit_flag.ne.NoError) RETURN ! !----------------------------------------------------------------------- ! Initialize all time levels and compute other initial fields. !----------------------------------------------------------------------- ! IF (iic(ng).eq.ntstart(ng)) THEN ! ! Initialize free-surface. ! DO thread=0,numthreads-1 subs=NtileX(ng)*NtileE(ng)/numthreads DO tile=subs*thread,subs*(thread+1)-1,+1 CALL ini_zeta (ng, MyRank, iNLM) END DO END DO ! ! Initialize other state variables. ! DO thread=0,numthreads-1 subs=NtileX(ng)*NtileE(ng)/numthreads DO tile=subs*(thread+1)-1,subs*thread,-1 CALL ini_fields (ng, MyRank, iNLM) END DO END DO END IF ! !----------------------------------------------------------------------- ! Compute horizontal mass fluxes (Hz*u/n and Hz*v/m), density related ! quatities and report global diagnostics. !----------------------------------------------------------------------- ! DO thread=0,numthreads-1 subs=NtileX(ng)*NtileE(ng)/numthreads DO tile=subs*thread,subs*(thread+1)-1,+1 CALL set_massflux (ng, MyRank) CALL rho_eos (ng, MyRank) CALL diag (ng, MyRank) END DO END DO IF (exit_flag.ne.NoError) RETURN ! !----------------------------------------------------------------------- ! Set fields for vertical boundary conditions. Process tidal forcing, ! if any. !----------------------------------------------------------------------- ! DO thread=0,numthreads-1 subs=NtileX(ng)*NtileE(ng)/numthreads DO tile=subs*thread,subs*(thread+1)-1,+1 CALL bulk_flux (ng, MyRank) CALL bblm (ng, MyRank) CALL set_vbc (ng, MyRank) CALL set_tides (ng, MyRank) END DO END DO ! !----------------------------------------------------------------------- ! Compute time-dependent vertical/horizontal mixing coefficients for ! momentum and tracers. Compute S-coordinate vertical velocity, ! diagnostically from horizontal mass divergence. !----------------------------------------------------------------------- ! DO thread=0,numthreads-1 subs=NtileX(ng)*NtileE(ng)/numthreads DO tile=subs*(thread+1)-1,subs*thread,-1 CALL omega (ng, MyRank) CALL wvelocity (ng, MyRank, nstp(ng)) END DO END DO ! !----------------------------------------------------------------------- ! Set free-surface to it time-averaged value. If applicable, ! accumulate time-averaged output data which needs a irreversible ! loop in shared-memory jobs. !----------------------------------------------------------------------- ! DO thread=0,numthreads-1 subs=NtileX(ng)*NtileE(ng)/numthreads DO tile=subs*thread,subs*(thread+1)-1,+1 ! irreversible loop CALL set_zeta (ng, MyRank) END DO END DO ! !----------------------------------------------------------------------- ! If appropriate, write out fields into output NetCDF files. Notice ! that IO data is written in delayed and serial mode. Exit if last ! time step. !----------------------------------------------------------------------- ! CALL output (ng) IF ((exit_flag.ne.NoError).or.(iic(ng).eq.(ntend(ng)+1))) RETURN ! !----------------------------------------------------------------------- ! Compute right-hand-side terms for 3D equations. !----------------------------------------------------------------------- ! DO thread=0,numthreads-1 subs=NtileX(ng)*NtileE(ng)/numthreads DO tile=subs*(thread+1)-1,subs*thread,-1 CALL rhs3d (ng, MyRank) CALL gls_prestep (ng, MyRank) END DO END DO ! !----------------------------------------------------------------------- ! Solve the vertically integrated primitive equations for the ! free-surface and barotropic momentum components. !----------------------------------------------------------------------- ! DO my_iif=1,nfast(ng)+1 ! ! Set time indices for predictor step. The PREDICTOR_2D_STEP switch ! it is assumed to be false before the first time-step. ! next_indx1=3-indx1(ng) IF (.not.PREDICTOR_2D_STEP(ng)) THEN PREDICTOR_2D_STEP(ng)=.TRUE. iif(ng)=my_iif IF (iif(ng).eq.1) THEN kstp(ng)=indx1(ng) ELSE kstp(ng)=3-indx1(ng) END IF knew(ng)=3 krhs(ng)=indx1(ng) END IF ! ! Predictor step - Advance barotropic equations using 2D time-step ! ============== predictor scheme. No actual time-stepping is ! performed during the auxiliary (nfast+1) time-step. It is needed ! to finalize the fast-time averaging of 2D fields, if any, and ! compute the new time-evolving depths. ! DO thread=0,numthreads-1 subs=NtileX(ng)*NtileE(ng)/numthreads DO tile=subs*(thread+1)-1,subs*thread,-1 CALL step2d (ng, MyRank) END DO END DO ! ! Set time indices for corrector step. ! IF (PREDICTOR_2D_STEP(ng)) THEN PREDICTOR_2D_STEP(ng)=.FALSE. knew(ng)=next_indx1 kstp(ng)=3-knew(ng) krhs(ng)=3 IF (iif(ng).lt.(nfast(ng)+1)) indx1(ng)=next_indx1 END IF ! ! Corrector step - Apply 2D time-step corrector scheme. Notice that ! ============== there is not need for a corrector step during the ! auxiliary (nfast+1) time-step. ! IF (iif(ng).lt.(nfast(ng)+1)) THEN DO thread=0,numthreads-1 subs=NtileX(ng)*NtileE(ng)/numthreads DO tile=subs*thread,subs*(thread+1)-1,+1 CALL step2d (ng, MyRank) END DO END DO END IF END DO ! !----------------------------------------------------------------------- ! Recompute depths and thicknesses using the new time filtered ! free-surface. This call was moved from "step2d" to here. !----------------------------------------------------------------------- ! DO thread=0,numthreads-1 subs=NtileX(ng)*NtileE(ng)/numthreads DO tile=subs*(thread+1)-1,subs*thread,-1 CALL set_depth (ng, MyRank) END DO END DO ! !----------------------------------------------------------------------- ! Time-step 3D momentum equations. !----------------------------------------------------------------------- ! ! Time-step 3D momentum equations and couple with vertically ! integrated equations. ! DO thread=0,numthreads-1 subs=NtileX(ng)*NtileE(ng)/numthreads DO tile=subs*(thread+1)-1,subs*thread,-1 CALL step3d_uv (ng, MyRank) END DO END DO ! !----------------------------------------------------------------------- ! Time-step vertical mixing turbulent equations and passive tracer ! source and sink terms, if applicable. !----------------------------------------------------------------------- ! DO thread=0,numthreads-1 subs=NtileX(ng)*NtileE(ng)/numthreads DO tile=subs*thread,subs*(thread+1)-1,+1 CALL omega (ng, MyRank) CALL gls_corstep (ng, MyRank) CALL sediment (ng, MyRank) END DO END DO ! !----------------------------------------------------------------------- ! Time-step tracer equations. !----------------------------------------------------------------------- ! DO thread=0,numthreads-1 subs=NtileX(ng)*NtileE(ng)/numthreads DO tile=subs*(thread+1)-1,subs*thread,-1 CALL step3d_t (ng, MyRank) END DO END DO ! !----------------------------------------------------------------------- ! Compute Lagrangian drifters trajectories. !----------------------------------------------------------------------- ! IF (Lfloats(Ng)) THEN DO thread=0,numthreads-1 chunk_size=(Nfloats(ng)+numthreads-1)/numthreads Lstr=1+thread*chunk_size Lend=MIN(Nfloats(ng),Lstr+chunk_size-1) CALL step_floats (ng, Lstr, Lend) END DO ! ! Shift floats time indices. ! nfp1(ng)=MOD(nfp1(ng)+1,NFT+1) nf(ng) =MOD(nf(ng) +1,NFT+1) nfm1(ng)=MOD(nfm1(ng)+1,NFT+1) nfm2(ng)=MOD(nfm2(ng)+1,NFT+1) nfm3(ng)=MOD(nfm3(ng)+1,NFT+1) END IF RETURN END SUBROUTINE main3d