2#if defined NONLINEAR && !defined SOLVE3D
3 SUBROUTINE main2d (RunInterval)
17# if defined STEP2D_FB_LF_AM3
19# elif defined STEP2D_FB_LF_AM3
29# if defined MODEL_COUPLING && defined MCT_LIB
44# ifdef TIDE_GENERATING_FORCES
47# if defined NLM_OUTER || \
49 defined rbl4dvar_ana_sensitivity || \
50 defined rbl4dvar_fct_sensitivity || \
58# if defined ATM_COUPLING && defined MCT_LIB
61# if defined WAV_COUPLING && defined MCT_LIB
70# if defined ADJUST_BOUNDARY
76# if defined SSH_TIDES || defined UV_TIDES
90 real(dp),
intent(in) :: RunInterval
94 logical :: DoNestLayer, Time_Step
96 integer :: Nsteps, Rsteps
97 integer :: ig, il, istep, ng, nl, tile
100 integer :: Lend, Lstr, chunk_size
103 character (len=*),
parameter :: MyFile = &
116 kernel_loop :
DO WHILE (time_step)
128 nest_layer :
DO WHILE (donestlayer)
143 step_loop :
DO istep=1,nsteps
164 & __line__, myfile))
RETURN
182# if defined NLM_OUTER || \
183 defined rbl4dvar || \
184 defined rbl4dvar_ana_sensitivity || \
185 defined rbl4dvar_fct_sensitivity
197# if defined WEAK_NOINTERP
198 IF ((
iic(ng).gt.1).and.(
iic(ng).ne.
ntend(ng)+1).and. &
199 & (mod(
iic(ng)-1,
nadj(ng)).eq.0))
THEN
201 WRITE (
stdout,*)
' FORCING NLM at iic = ',
iic(ng)
265# ifdef TIDE_GENERATING_FORCES
266 CALL equilibrium_tide (ng, tile,
inlm)
277# if defined ATM_COUPLING && defined MCT_LIB
289 CALL ocn2atm_coupling (ng, tile)
296# ifdef ADJUST_BOUNDARY
316# ifdef ADJUST_WSTRESS
335# if defined WAV_COUPLING && defined MCT_LIB
347 CALL ocn2wav_coupling (ng, tile)
363# if defined SSH_TIDES || defined UV_TIDES
426# ifdef STEP2D_FB_AB3_AM4
443 IF (mod(
knew(ng),2).eq.0)
THEN
458# ifdef STEP2D_FB_LF_AM3
535 next_indx1=3-
indx1(ng)
538 IF (first_2d_step)
THEN
609# if defined MASKING && defined WET_DRY
666 IF (do_twoway(
inlm, nl, il, ng, istep))
THEN
705 lend=min(
nfloats(ng),lstr+chunk_size-1)
subroutine, public time_string(mytime, date_string)
subroutine, public diag(ng, tile)
subroutine, public nl_dotproduct(ng, tile, linp)
subroutine, public forcing(ng, tile, kfrc, nfrc)
subroutine, public load_frc(ng, tile, lout)
subroutine, public frc_adjust(ng, tile, linp)
subroutine, public ini_zeta(ng, tile, model)
subroutine, public ini_fields(ng, tile, model)
integer, dimension(:,:), allocatable couplesteps
real(r8), dimension(:), allocatable twowayinterval
logical, dimension(:), allocatable donortofiner
integer, dimension(:), allocatable first_tile
integer, dimension(:), allocatable last_tile
integer, dimension(:), allocatable nfloats
integer, dimension(:,:), allocatable gridnumber
integer, dimension(:), allocatable gridsinlayer
logical, dimension(:), allocatable lfloats
integer, dimension(:), allocatable iic
real(dp), dimension(:), allocatable dt
integer, dimension(:), allocatable next_kstp
logical, dimension(:), allocatable predictor_2d_step
real(dp), dimension(:), allocatable tdays
logical, dimension(:), allocatable frequentimpulse
integer, dimension(:), allocatable nfast
real(dp), parameter sec2day
integer, dimension(:), allocatable ntend
character(len=22), dimension(:), allocatable time_code
integer, dimension(:), allocatable indx1
logical, dimension(:,:), allocatable compositegrid
real(dp), dimension(:), allocatable time
logical, dimension(:), allocatable refinedgrid
integer, dimension(:), allocatable refinescale
integer, dimension(:), allocatable ntstart
integer, dimension(:), allocatable step_counter
integer, dimension(:), allocatable nadj
integer, dimension(:), allocatable iif
integer, dimension(:), allocatable nfm2
integer, dimension(:), allocatable lbout
integer, dimension(:), allocatable kstp
integer, dimension(:), allocatable knew
integer, dimension(:), allocatable nfm1
integer, dimension(:), allocatable lfinp
integer, dimension(:), allocatable lbinp
integer, dimension(:), allocatable nf
integer, dimension(:), allocatable nfm3
integer, dimension(:), allocatable nfp1
integer, dimension(:), allocatable lnew
integer, dimension(:), allocatable krhs
integer, dimension(:), allocatable lfout
integer, dimension(:), allocatable nstp
subroutine, public nesting(ng, model, isection)
subroutine, public obc_adjust(ng, tile, linp)
subroutine, public load_obc(ng, tile, lout)
subroutine, public set_avg(ng, tile)
subroutine, public set_tides(ng, tile)
subroutine, public set_vbc(ng, tile)
subroutine, public step2d(ng, tile)
subroutine, public step_floats(ng, lstr, lend)
logical function, public founderror(flag, noerr, line, routine)
subroutine ntimesteps(model, runinterval, nl, nsteps, rsteps)
subroutine set_data(ng, tile)
subroutine set_diags(ng, tile)