ROMS
Loading...
Searching...
No Matches
step2d_mod Module Reference

Functions/Subroutines

subroutine, public step2d (ng, tile)
 
subroutine step2d_tile (ng, tile, lbi, ubi, lbj, ubj, ubk, imins, imaxs, jmins, jmaxs, krhs, kstp, knew, nstp, nnew, pmask, rmask, umask, vmask, pmask_wet, pmask_full, rmask_wet, rmask_full, umask_wet, umask_full, vmask_wet, vmask_full, rmask_wet_avg, fomn, h, om_u, om_v, on_u, on_v, pm, pn, dndx, dmde, rdrag, rdrag2, pmon_r, pnom_r, pmon_p, pnom_p, om_r, on_r, om_p, on_p, visc2_p, visc2_r, visc4_p, visc4_r, bed_thick, eq_tide, sustr, svstr,
 
subroutine step2d_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, kstp, knew, nstp, nnew, pmask, rmask, umask, vmask, pmask_wet, pmask_full, rmask_wet, rmask_full, umask_wet, umask_full, vmask_wet, vmask_full, rmask_wet_avg, fomn, h, om_u, om_v, on_u, on_v, omn, pm, pn, dndx, dmde, pmon_r, pnom_r, pmon_p, pnom_p, om_r, on_r, om_p, on_p, visc2_p, visc2_r, bed_thick, eq_tide, sustr, svstr, bustr, bvstr, pair, rhoa, rhos, du_avg1, du_avg2, dv_avg1, dv_avg2, zt_avg1, rufrc, rvfrc, rufrc_bak, rvfrc_bak, diau2wrk, diav2wrk, diarubar, diarvbar, diau2int, diav2int, diarufrc, diarvfrc, du_flux, dv_flux, ubar, vbar, zeta)
 
subroutine step2d_tile (ng, tile, lbi, ubi, lbj, ubj, ubk, imins, imaxs, jmins, jmaxs, krhs, kstp, knew, nstp, nnew, pmask, rmask, umask, vmask, pmask_wet, pmask_full, rmask_wet, rmask_full, umask_wet, umask_full, vmask_wet, vmask_full, rmask_wet_avg, fomn, h, om_u, om_v, on_u, on_v, omn, pm, pn, dndx, dmde, pmon_r, pnom_r, pmon_p, pnom_p, om_r, on_r, om_p, on_p, visc2_p, visc2_r, visc4_p, visc4_r, bed_thick, rurol2d, rvrol2d, rubst2d, rvbst2d, russt2d, rvsst2d, rubrk2d, rvbrk2d, rukvf2d, rvkvf2d, bh, qsp, zetaw, ubar_stokes, vbar_stokes, eq_tide, sustr, svstr, bustr, bvstr, pair, rhoa, rhos, du_avg1, du_avg2, dv_avg1, dv_avg2, zt_avg1, rufrc, rvfrc, ru, rv, diau2wrk, diav2wrk, diarubar, diarvbar, diau2int, diav2int, diarufrc, diarvfrc, du_flux, dv_flux, rubar, rvbar, rzeta, ubar, vbar, zeta)
 

Function/Subroutine Documentation

◆ step2d()

subroutine public step2d_mod::step2d ( integer, intent(in) ng,
integer, intent(in) tile )

Definition at line 74 of file step2d_FB.h.

75!***********************************************************************
76!
77! Imported variable declarations.
78!
79 integer, intent(in) :: ng, tile
80!
81! Local variable declarations.
82!
83 character (len=*), parameter :: MyFile = &
84 & __FILE__
85!
86#include "tile.h"
87!
88#ifdef PROFILE
89 CALL wclock_on (ng, inlm, 9, __line__, myfile)
90#endif
91 CALL step2d_tile (ng, tile, &
92 & lbi, ubi, lbj, ubj, n(ng), &
93 & imins, imaxs, jmins, jmaxs, &
94 & krhs(ng), kstp(ng), knew(ng), &
95#ifdef SOLVE3D
96 & nstp(ng), nnew(ng), &
97#endif
98#ifdef MASKING
99 & grid(ng) % pmask, grid(ng) % rmask, &
100 & grid(ng) % umask, grid(ng) % vmask, &
101#endif
102#ifdef WET_DRY
103 & grid(ng) % pmask_wet, grid(ng) % pmask_full, &
104 & grid(ng) % rmask_wet, grid(ng) % rmask_full, &
105 & grid(ng) % umask_wet, grid(ng) % umask_full, &
106 & grid(ng) % vmask_wet, grid(ng) % vmask_full, &
107# ifdef SOLVE3D
108 & grid(ng) % rmask_wet_avg, &
109# endif
110#endif
111#if (defined UV_COR && !defined SOLVE3D) || defined STEP2D_CORIOLIS
112 & grid(ng) % fomn, &
113#endif
114 & grid(ng) % h, &
115 & grid(ng) % om_u, grid(ng) % om_v, &
116 & grid(ng) % on_u, grid(ng) % on_v, &
117 & grid(ng) % pm, grid(ng) % pn, &
118#if defined CURVGRID && defined UV_ADV && !defined SOLVE3D
119 & grid(ng) % dndx, grid(ng) % dmde, &
120#endif
121 & grid(ng) % rdrag, &
122#if defined UV_QDRAG && !defined SOLVE3D
123 & grid(ng) % rdrag2, &
124#endif
125#if (defined UV_VIS2 || defined UV_VIS4) && !defined SOLVE3D
126 & grid(ng) % pmon_r, grid(ng) % pnom_r, &
127 & grid(ng) % pmon_p, grid(ng) % pnom_p, &
128 & grid(ng) % om_r, grid(ng) % on_r, &
129 & grid(ng) % om_p, grid(ng) % on_p, &
130# ifdef UV_VIS2
131 & mixing(ng) % visc2_p, mixing(ng) % visc2_r, &
132# endif
133# ifdef UV_VIS4
134 & mixing(ng) % visc4_p, mixing(ng) % visc4_r, &
135# endif
136#endif
137#if defined TIDE_GENERATING_FORCES && !defined SOLVE3D
138 & ocean(ng) % eq_tide, &
139#endif
140#ifndef SOLVE3D
141 & forces(ng) % sustr, forces(ng) % svstr, &
142# ifdef ATM_PRESS
143 & forces(ng) % Pair, &
144# endif
145#else
146# ifdef VAR_RHO_2D
147 & coupling(ng) % rhoA, coupling(ng) % rhoS, &
148# endif
149 & coupling(ng) % DU_avg1, coupling(ng) % DU_avg2, &
150 & coupling(ng) % DV_avg1, coupling(ng) % DV_avg2, &
151 & coupling(ng) % Zt_avg1, &
152 & coupling(ng) % rufrc, &
153 & coupling(ng) % rvfrc, &
154 & coupling(ng) % rufrc_bak, &
155 & coupling(ng) % rvfrc_bak, &
156#endif
157#if defined NESTING && !defined SOLVE3D
158 & ocean(ng) % DU_flux, ocean(ng) % DV_flux, &
159#endif
160 & ocean(ng) % ubar, ocean(ng) % vbar, &
161 & ocean(ng) % zeta)
162#ifdef PROFILE
163 CALL wclock_off (ng, inlm, 9, __line__, myfile)
164#endif
165!
166 RETURN
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3

References mod_coupling::coupling, mod_forces::forces, mod_grid::grid, mod_param::inlm, mod_stepping::knew, mod_stepping::krhs, mod_stepping::kstp, mod_mixing::mixing, mod_param::n, mod_stepping::nnew, mod_stepping::nstp, mod_ocean::ocean, step2d_tile(), wclock_off(), and wclock_on().

Referenced by main3d().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ step2d_tile() [1/3]

subroutine step2d_mod::step2d_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) kstp,
integer, intent(in) knew,
integer, intent(in) nstp,
integer, intent(in) nnew,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pmask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rmask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) umask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) vmask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) pmask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) pmask_full,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) rmask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) rmask_full,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) umask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) umask_full,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) vmask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) vmask_full,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) rmask_wet_avg,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) fomn,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) h,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_u,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_v,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_u,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_v,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) omn,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pm,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pn,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) dndx,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) dmde,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pmon_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pnom_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pmon_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pnom_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) visc2_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) visc2_r,
real(r8), dimension(lbi:ubi,lbj:ubj,1:3), intent(in) bed_thick,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) eq_tide,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) sustr,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) svstr,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) bustr,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) bvstr,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pair,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rhoa,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rhos,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) du_avg1,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) du_avg2,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) dv_avg1,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) dv_avg2,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) zt_avg1,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) rufrc,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) rvfrc,
real(r8), dimension(lbi:ubi,lbj:ubj,2), intent(inout) rufrc_bak,
real(r8), dimension(lbi:ubi,lbj:ubj,2), intent(inout) rvfrc_bak,
real(r8), dimension(lbi:ubi,lbj:ubj,ndm2d), intent(inout) diau2wrk,
real(r8), dimension(lbi:ubi,lbj:ubj,ndm2d), intent(inout) diav2wrk,
real(r8), dimension(lbi:ubi,lbj:ubj,2,ndm2d-1), intent(inout) diarubar,
real(r8), dimension(lbi:ubi,lbj:ubj,2,ndm2d-1), intent(inout) diarvbar,
real(r8), dimension(lbi:ubi,lbj:ubj,ndm2d), intent(inout) diau2int,
real(r8), dimension(lbi:ubi,lbj:ubj,ndm2d), intent(inout) diav2int,
real(r8), dimension(lbi:ubi,lbj:ubj,3,ndm2d-1), intent(inout) diarufrc,
real(r8), dimension(lbi:ubi,lbj:ubj,3,ndm2d-1), intent(inout) diarvfrc,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(out) du_flux,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(out) dv_flux,
real(r8), dimension(lbi:ubi,lbj:ubj,:), intent(inout) ubar,
real(r8), dimension(lbi:ubi,lbj:ubj,:), intent(inout) vbar,
real(r8), dimension(lbi:ubi,lbj:ubj,:), intent(inout) zeta )
private

Definition at line 171 of file step2d_FB_LF_AM3.h.

236!***********************************************************************
237!
238! Imported variable declarations.
239!
240 integer, intent(in ) :: ng, tile
241 integer, intent(in ) :: LBi, UBi, LBj, UBj
242 integer, intent(in ) :: IminS, ImaxS, JminS, JmaxS
243 integer, intent(in ) :: kstp, knew
244#ifdef SOLVE3D
245 integer, intent(in ) :: nstp, nnew
246#endif
247!
248#ifdef ASSUMED_SHAPE
249# ifdef MASKING
250 real(r8), intent(in ) :: pmask(LBi:,LBj:)
251 real(r8), intent(in ) :: rmask(LBi:,LBj:)
252 real(r8), intent(in ) :: umask(LBi:,LBj:)
253 real(r8), intent(in ) :: vmask(LBi:,LBj:)
254# endif
255# if (defined UV_COR && !defined SOLVE3D) || defined STEP2D_CORIOLIS
256 real(r8), intent(in ) :: fomn(LBi:,LBj:)
257# endif
258# if defined SEDIMENT && defined SED_MORPH
259 real(r8), intent(inout) :: h(LBi:,LBj:)
260# else
261 real(r8), intent(in ) :: h(LBi:,LBj:)
262# endif
263 real(r8), intent(in ) :: om_u(LBi:,LBj:)
264 real(r8), intent(in ) :: om_v(LBi:,LBj:)
265 real(r8), intent(in ) :: on_u(LBi:,LBj:)
266 real(r8), intent(in ) :: on_v(LBi:,LBj:)
267 real(r8), intent(in ) :: omn(LBi:,LBj:)
268 real(r8), intent(in ) :: pm(LBi:,LBj:)
269 real(r8), intent(in ) :: pn(LBi:,LBj:)
270# if defined CURVGRID && defined UV_ADV && !defined SOLVE3D
271 real(r8), intent(in ) :: dndx(LBi:,LBj:)
272 real(r8), intent(in ) :: dmde(LBi:,LBj:)
273# endif
274# if defined UV_VIS2 && !defined SOLVE3D
275 real(r8), intent(in ) :: pmon_r(LBi:,LBj:)
276 real(r8), intent(in ) :: pnom_r(LBi:,LBj:)
277 real(r8), intent(in ) :: pmon_p(LBi:,LBj:)
278 real(r8), intent(in ) :: pnom_p(LBi:,LBj:)
279 real(r8), intent(in ) :: om_r(LBi:,LBj:)
280 real(r8), intent(in ) :: on_r(LBi:,LBj:)
281 real(r8), intent(in ) :: om_p(LBi:,LBj:)
282 real(r8), intent(in ) :: on_p(LBi:,LBj:)
283 real(r8), intent(in ) :: visc2_p(LBi:,LBj:)
284 real(r8), intent(in ) :: visc2_r(LBi:,LBj:)
285# endif
286# if defined SEDIMENT && defined SED_MORPH
287 real(r8), intent(in ) :: bed_thick(LBi:,LBj:,:)
288# endif
289# if defined TIDE_GENERATING_FORCES && !defined SOLVE3D
290 real(r8), intent(in ) :: eq_tide(LBi:,LBj:)
291# endif
292# ifndef SOLVE3D
293 real(r8), intent(in ) :: sustr(LBi:,LBj:)
294 real(r8), intent(in ) :: svstr(LBi:,LBj:)
295 real(r8), intent(in ) :: bustr(LBi:,LBj:)
296 real(r8), intent(in ) :: bvstr(LBi:,LBj:)
297# ifdef ATM_PRESS
298 real(r8), intent(in ) :: Pair(LBi:,LBj:)
299# endif
300# else
301# ifdef VAR_RHO_2D
302 real(r8), intent(in ) :: rhoA(LBi:,LBj:)
303 real(r8), intent(in ) :: rhoS(LBi:,LBj:)
304# endif
305 real(r8), intent(inout) :: DU_avg1(LBi:,LBj:)
306 real(r8), intent(inout) :: DU_avg2(LBi:,LBj:)
307 real(r8), intent(inout) :: DV_avg1(LBi:,LBj:)
308 real(r8), intent(inout) :: DV_avg2(LBi:,LBj:)
309 real(r8), intent(inout) :: Zt_avg1(LBi:,LBj:)
310 real(r8), intent(inout) :: rufrc(LBi:,LBj:)
311 real(r8), intent(inout) :: rvfrc(LBi:,LBj:)
312 real(r8), intent(inout) :: rufrc_bak(LBi:,LBj:,:)
313 real(r8), intent(inout) :: rvfrc_bak(LBi:,LBj:,:)
314# endif
315# ifdef WET_DRY
316 real(r8), intent(inout) :: pmask_full(LBi:,LBj:)
317 real(r8), intent(inout) :: rmask_full(LBi:,LBj:)
318 real(r8), intent(inout) :: umask_full(LBi:,LBj:)
319 real(r8), intent(inout) :: vmask_full(LBi:,LBj:)
320
321 real(r8), intent(inout) :: pmask_wet(LBi:,LBj:)
322 real(r8), intent(inout) :: rmask_wet(LBi:,LBj:)
323 real(r8), intent(inout) :: umask_wet(LBi:,LBj:)
324 real(r8), intent(inout) :: vmask_wet(LBi:,LBj:)
325# ifdef SOLVE3D
326 real(r8), intent(inout) :: rmask_wet_avg(LBi:,LBj:)
327# endif
328# endif
329# ifdef DIAGNOSTICS_UV
330 real(r8), intent(inout) :: DiaU2wrk(LBi:,LBj:,:)
331 real(r8), intent(inout) :: DiaV2wrk(LBi:,LBj:,:)
332 real(r8), intent(inout) :: DiaRUbar(LBi:,LBj:,:,:)
333 real(r8), intent(inout) :: DiaRVbar(LBi:,LBj:,:,:)
334# ifdef SOLVE3D
335 real(r8), intent(inout) :: DiaU2int(LBi:,LBj:,:)
336 real(r8), intent(inout) :: DiaV2int(LBi:,LBj:,:)
337 real(r8), intent(inout) :: DiaRUfrc(LBi:,LBj:,:,:)
338 real(r8), intent(inout) :: DiaRVfrc(LBi:,LBj:,:,:)
339# endif
340# endif
341 real(r8), intent(inout) :: ubar(LBi:,LBj:,:)
342 real(r8), intent(inout) :: vbar(LBi:,LBj:,:)
343 real(r8), intent(inout) :: zeta(LBi:,LBj:,:)
344# if defined NESTING && !defined SOLVE3D
345 real(r8), intent(out ) :: DU_flux(LBi:,LBj:)
346 real(r8), intent(out ) :: DV_flux(LBi:,LBj:)
347# endif
348
349#else
350
351# ifdef MASKING
352 real(r8), intent(in ) :: pmask(LBi:UBi,LBj:UBj)
353 real(r8), intent(in ) :: rmask(LBi:UBi,LBj:UBj)
354 real(r8), intent(in ) :: umask(LBi:UBi,LBj:UBj)
355 real(r8), intent(in ) :: vmask(LBi:UBi,LBj:UBj)
356# endif
357# if (defined UV_COR && !defined SOLVE3D) || defined STEP2D_CORIOLIS
358 real(r8), intent(in ) :: fomn(LBi:UBi,LBj:UBj)
359# endif
360# if defined SEDIMENT && defined SED_MORPH
361 real(r8), intent(inout) :: h(LBi:UBi,LBj:UBj)
362# else
363 real(r8), intent(in ) :: h(LBi:UBi,LBj:UBj)
364# endif
365 real(r8), intent(in ) :: om_u(LBi:UBi,LBj:UBj)
366 real(r8), intent(in ) :: om_v(LBi:UBi,LBj:UBj)
367 real(r8), intent(in ) :: on_u(LBi:UBi,LBj:UBj)
368 real(r8), intent(in ) :: on_v(LBi:UBi,LBj:UBj)
369 real(r8), intent(in ) :: omn(LBi:UBi,LBj:UBj)
370 real(r8), intent(in ) :: pm(LBi:UBi,LBj:UBj)
371 real(r8), intent(in ) :: pn(LBi:UBi,LBj:UBj)
372# if defined CURVGRID && defined UV_ADV && !defined SOLVE3D
373 real(r8), intent(in ) :: dndx(LBi:UBi,LBj:UBj)
374 real(r8), intent(in ) :: dmde(LBi:UBi,LBj:UBj)
375# endif
376# if defined UV_VIS2 && !defined SOLVE3D
377 real(r8), intent(in ) :: pmon_r(LBi:UBi,LBj:UBj)
378 real(r8), intent(in ) :: pnom_r(LBi:UBi,LBj:UBj)
379 real(r8), intent(in ) :: pmon_p(LBi:UBi,LBj:UBj)
380 real(r8), intent(in ) :: pnom_p(LBi:UBi,LBj:UBj)
381 real(r8), intent(in ) :: om_r(LBi:UBi,LBj:UBj)
382 real(r8), intent(in ) :: on_r(LBi:UBi,LBj:UBj)
383 real(r8), intent(in ) :: om_p(LBi:UBi,LBj:UBj)
384 real(r8), intent(in ) :: on_p(LBi:UBi,LBj:UBj)
385 real(r8), intent(in ) :: visc2_p(LBi:UBi,LBj:UBj)
386 real(r8), intent(in ) :: visc2_r(LBi:UBi,LBj:UBj)
387# endif
388# if defined SEDIMENT && defined SED_MORPH
389 real(r8), intent(in ) :: bed_thick(LBi:UBi,LBj:UBj,1:3)
390# endif
391# if defined TIDE_GENERATING_FORCES && !defined SOLVE3D
392 real(r8), intent(in ) :: eq_tide(LBi:UBi,LBj:UBj)
393# endif
394# ifndef SOLVE3D
395 real(r8), intent(in ) :: sustr(LBi:UBi,LBj:UBj)
396 real(r8), intent(in ) :: svstr(LBi:UBi,LBj:UBj)
397 real(r8), intent(in ) :: bustr(LBi:UBi,LBj:UBj)
398 real(r8), intent(in ) :: bvstr(LBi:UBi,LBj:UBj)
399# ifdef ATM_PRESS
400 real(r8), intent(in ) :: Pair(LBi:UBi,LBj:UBj)
401# endif
402# else
403# ifdef VAR_RHO_2D
404 real(r8), intent(in ) :: rhoA(LBi:UBi,LBj:UBj)
405 real(r8), intent(in ) :: rhoS(LBi:UBi,LBj:UBj)
406# endif
407 real(r8), intent(inout) :: DU_avg1(LBi:UBi,LBj:UBj)
408 real(r8), intent(inout) :: DU_avg2(LBi:UBi,LBj:UBj)
409 real(r8), intent(inout) :: DV_avg1(LBi:UBi,LBj:UBj)
410 real(r8), intent(inout) :: DV_avg2(LBi:UBi,LBj:UBj)
411 real(r8), intent(inout) :: Zt_avg1(LBi:UBi,LBj:UBj)
412 real(r8), intent(inout) :: rufrc(LBi:UBi,LBj:UBj)
413 real(r8), intent(inout) :: rvfrc(LBi:UBi,LBj:UBj)
414 real(r8), intent(inout) :: rufrc_bak(LBi:UBi,LBj:UBj,2)
415 real(r8), intent(inout) :: rvfrc_bak(LBi:UBi,LBj:UBj,2)
416# endif
417# ifdef WET_DRY
418 real(r8), intent(inout) :: pmask_full(LBi:UBi,LBj:UBj)
419 real(r8), intent(inout) :: rmask_full(LBi:UBi,LBj:UBj)
420 real(r8), intent(inout) :: umask_full(LBi:UBi,LBj:UBj)
421 real(r8), intent(inout) :: vmask_full(LBi:UBi,LBj:UBj)
422
423 real(r8), intent(inout) :: pmask_wet(LBi:UBi,LBj:UBj)
424 real(r8), intent(inout) :: rmask_wet(LBi:UBi,LBj:UBj)
425 real(r8), intent(inout) :: umask_wet(LBi:UBi,LBj:UBj)
426 real(r8), intent(inout) :: vmask_wet(LBi:UBi,LBj:UBj)
427# ifdef SOLVE3D
428 real(r8), intent(inout) :: rmask_wet_avg(LBi:UBi,LBj:UBj)
429# endif
430# endif
431# ifdef DIAGNOSTICS_UV
432 real(r8), intent(inout) :: DiaU2wrk(LBi:UBi,LBj:UBj,NDM2d)
433 real(r8), intent(inout) :: DiaV2wrk(LBi:UBi,LBj:UBj,NDM2d)
434 real(r8), intent(inout) :: DiaRUbar(LBi:UBi,LBj:UBj,2,NDM2d-1)
435 real(r8), intent(inout) :: DiaRVbar(LBi:UBi,LBj:UBj,2,NDM2d-1)
436# ifdef SOLVE3D
437 real(r8), intent(inout) :: DiaU2int(LBi:UBi,LBj:UBj,NDM2d)
438 real(r8), intent(inout) :: DiaV2int(LBi:UBi,LBj:UBj,NDM2d)
439 real(r8), intent(inout) :: DiaRUfrc(LBi:UBi,LBj:UBj,3,NDM2d-1)
440 real(r8), intent(inout) :: DiaRVfrc(LBi:UBi,LBj:UBj,3,NDM2d-1)
441# endif
442# endif
443 real(r8), intent(inout) :: ubar(LBi:UBi,LBj:UBj,:)
444 real(r8), intent(inout) :: vbar(LBi:UBi,LBj:UBj,:)
445 real(r8), intent(inout) :: zeta(LBi:UBi,LBj:UBj,:)
446# if defined NESTING && !defined SOLVE3D
447 real(r8), intent(out ) :: DU_flux(LBi:UBi,LBj:UBj)
448 real(r8), intent(out ) :: DV_flux(LBi:UBi,LBj:UBj)
449# endif
450#endif
451!
452! Local variable declarations.
453!
454 integer :: i, is, j
455 integer :: krhs, kbak
456#ifdef DIAGNOSTICS_UV
457 integer :: idiag
458#endif
459!
460 real(r8) :: cff, cff1, cff2, cff3, cff4
461#ifdef WET_DRY
462 real(r8) :: cff5, cff6, cff7
463#endif
464 real(r8) :: fac, fac1, fac2
465!
466#if defined UV_C4ADVECTION && !defined SOLVE3D
467 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Dgrad
468#endif
469 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Dnew
470 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Drhs
471#if defined UV_VIS2 && !defined SOLVE3D
472 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Drhs_p
473#endif
474 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Dstp
475 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: DUon
476 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: DVom
477#if defined STEP2D_CORIOLIS || !defined SOLVE3D
478 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: UFx
479 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: VFe
480#endif
481#if !defined SOLVE3D
482 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: UFe
483 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: VFx
484#endif
485#if defined UV_C4ADVECTION && !defined SOLVE3D
486 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: grad
487#endif
488 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rubar
489 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rvbar
490 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rzeta
491 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rzeta2
492#if defined VAR_RHO_2D && defined SOLVE3D
493 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rzetaSA
494#endif
495 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: urhs
496 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: vrhs
497 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: zwrk
498#ifdef DIAGNOSTICS_UV
499 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Uwrk
500 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Vwrk
501 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,NDM2d-1) :: DiaU2rhs
502 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,NDM2d-1) :: DiaV2rhs
503#endif
504!
505 real(r8), allocatable :: zeta_new(:,:)
506!
507! The following stability limits are obtained empirically using 3/4
508! degree Atlantic model configuration. In all these cases barotropic
509! mode time step is about 180...250 seconds, which is much less than
510! the inertial period. The maximum stability coefficients turned out
511! to be slightly different than predicted by linear theory, although
512! all theoretical tendencies agree with the practice. Note the nearly
513! 70% gain in stability compared with LF-TR for appropriate
514! coefficients (linear theory predicts beta=0.166, epsil=0.84).
515!
516!! real(r8), parameter :: gamma=0.0_r8, &
517!! & beta =0.0_r8, epsil=0.0_r8 !--> Cu=0.818
518!! real(r8), parameter :: gamma=1.0_r8/12.0_r8, &
519!! & beta =0.0_r8, epsil=0.0_r8 !--> Cu=0.878
520!! real(r8), parameter :: gamma=1./12., &
521!! beta =0.1_r8, epsil=0.6_r8 !--> Cu=1.050
522 real(r8), parameter :: gamma=0.0_r8, &
523 & beta =0.14_r8, epsil=0.74_r8 !==> Cu=1.341
524
525#include "set_bounds.h"
526!
527!-----------------------------------------------------------------------
528! Timestep vertically integrated (barotropic) equations.
529!-----------------------------------------------------------------------
530!
531! In the code below it is assumed that variables with time index "krhs"
532! are time-centered at step "n" in barotropic time during predictor
533! sub-step and "n+1/2" during corrector.
534!
535 IF (predictor_2d_step) THEN
536 krhs=kstp
537 ELSE
538 krhs=3
539 END IF
540 IF (first_2d_step) THEN
541 kbak=kstp ! "kbak" is used as "from"
542 ELSE ! time index for LF timestep
543 kbak=3-kstp
544 END IF
545!
546!-----------------------------------------------------------------------
547! Preliminary steps.
548!-----------------------------------------------------------------------
549!
550! Compute total depth of the water column and vertically integrated
551! mass fluxes, which are used in computation of free-surface elevation
552! time tendency and advection terms for the barotropic momentum
553! equations.
554!
555#if defined DISTRIBUTE && !defined NESTING
556# define IR_RANGE IstrUm2-1,Iendp2
557# define JR_RANGE JstrVm2-1,Jendp2
558# define IU_RANGE IstrUm1-1,Iendp2
559# define JU_RANGE Jstrm1-1,Jendp2
560# define IV_RANGE Istrm1-1,Iendp2
561# define JV_RANGE JstrVm1-1,Jendp2
562#else
563# define IR_RANGE IstrUm2-1,Iendp2
564# define JR_RANGE JstrVm2-1,Jendp2
565# define IU_RANGE IstrUm2,Iendp2
566# define JU_RANGE JstrVm2-1,Jendp2
567# define IV_RANGE IstrUm2-1,Iendp2
568# define JV_RANGE JstrVm2,Jendp2
569#endif
570
571 DO j=jr_range
572 DO i=ir_range
573 drhs(i,j)=zeta(i,j,krhs)+h(i,j)
574 END DO
575 END DO
576 DO j=ju_range
577 DO i=iu_range
578 cff=0.5_r8*on_u(i,j)
579 cff1=cff*(drhs(i,j)+drhs(i-1,j))
580 duon(i,j)=ubar(i,j,krhs)*cff1
581 END DO
582 END DO
583 DO j=jv_range
584 DO i=iv_range
585 cff=0.5_r8*om_v(i,j)
586 cff1=cff*(drhs(i,j)+drhs(i,j-1))
587 dvom(i,j)=vbar(i,j,krhs)*cff1
588 END DO
589 END DO
590
591#undef IR_RANGE
592#undef IU_RANGE
593#undef IV_RANGE
594#undef JR_RANGE
595#undef JU_RANGE
596#undef JV_RANGE
597
598#if defined DISTRIBUTE && \
599 defined uv_adv && defined uv_c4advection && !defined SOLVE3D
600!
601! In distributed-memory, the I- and J-ranges are different and a
602! special exchange is done here to avoid having three ghost points
603! for high-order numerical stencils. Notice that a private array is
604! passed below to the exchange routine. It also applies periodic
605! boundary conditions, if appropriate and no partitions in I- or
606! J-directions.
607!
608 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
609 CALL exchange_u2d_tile (ng, tile, &
610 & imins, imaxs, jmins, jmaxs, &
611 & duon)
612 CALL exchange_v2d_tile (ng, tile, &
613 & imins, imaxs, jmins, jmaxs, &
614 & dvom)
615 END IF
616 CALL mp_exchange2d (ng, tile, inlm, 2, &
617 & imins, imaxs, jmins, jmaxs, &
618 & nghostpoints, &
619 & ewperiodic(ng), nsperiodic(ng), &
620 & duon, dvom)
621#endif
622!
623! Set vertically integrated mass fluxes DUon and DVom along the open
624! boundaries in such a way that the integral volume is conserved.
625!
626 IF (any(volcons(:,ng))) THEN
627 CALL set_duv_bc_tile (ng, tile, &
628 & lbi, ubi, lbj, ubj, &
629 & imins, imaxs, jmins, jmaxs, &
630 & krhs, &
631#ifdef MASKING
632 & umask, vmask, &
633#endif
634 & om_v, on_u, &
635 & ubar, vbar, &
636 & drhs, duon, dvom)
637 END IF
638
639#ifdef SOLVE3D
640!
641!-----------------------------------------------------------------------
642! Fields averaged over all barotropic time steps.
643!-----------------------------------------------------------------------
644!
645! Notice that the index ranges here are designed to include physical
646! boundaries only. Periodic ghost points and internal mpi computational
647! margins are NOT included.
648!
649! Reset all barotropic mode time-averaged arrays during the first
650! predictor step. At all subsequent time steps, accumulate averages
651! of the first kind using the DELAYED way. For example, "Zt_avg1" is
652! not summed immediately after the corrector step when computed but
653! during the subsequent predictor substep. It allows saving operations
654! because "DUon" and "DVom" are calculated anyway. The last time step
655! has a special code to add all three barotropic variables after the
656! last corrector substep.
657!
658 IF (predictor_2d_step) THEN ! PREDICTOR STEP
659 IF (first_2d_step) THEN
660 DO j=jstrr,jendr
661 DO i=istrr,iendr
662 zt_avg1(i,j)=0.0_r8
663 du_avg1(i,j)=0.0_r8
664 dv_avg1(i,j)=0.0_r8
665 du_avg2(i,j)=0.0_r8
666 dv_avg2(i,j)=0.0_r8
667 END DO
668 END DO
669 ELSE
670 cff=weight(1,iif(ng)-1,ng)
671 DO j=jstrr,jendr
672 DO i=istrr,iendr
673 zt_avg1(i,j)=zt_avg1(i,j)+cff*zeta(i,j,krhs)
674 IF (i.ge.istr) THEN
675 du_avg1(i,j)=du_avg1(i,j)+cff*duon(i,j)
676 END IF
677 IF (j.ge.jstr) THEN
678 dv_avg1(i,j)=dv_avg1(i,j)+cff*dvom(i,j)
679 END IF
680 END DO
681 END DO
682 END IF
683 ELSE ! CORRECTOR STEP
684 cff=weight(2,iif(ng),ng)
685 DO j=jstrr,jendr
686 DO i=istrr,iendr
687 IF (i.ge.istr) THEN
688 du_avg2(i,j)=du_avg2(i,j)+cff*duon(i,j)
689 END IF
690 IF (j.ge.jstr) THEN
691 dv_avg2(i,j)=dv_avg2(i,j)+cff*dvom(i,j)
692 END IF
693 END DO
694 END DO
695 END IF
696#endif
697!
698!-----------------------------------------------------------------------
699! Advance free-surface.
700!-----------------------------------------------------------------------
701!
702! Notice that the new local free-surface is allocated so it can be
703! passed as an argumment to "zetabc" to avoid memory issues.
704!
705 allocate ( zeta_new(imins:imaxs,jmins:jmaxs) )
706 zeta_new = 0.0_r8
707!
708! Compute "zeta_new" at the new time step and interpolate backward for
709! the subsequent computation of barotropic pressure-gradient terms.
710! Notice that during the predictor of the first 2D step in 3D mode,
711! the pressure gradient terms are computed using just zeta(:,:,kstp),
712! i.e., like in the Forward Euler step, rather than the more accurate
713! predictor of generalized RK2. This is to keep it consistent with the
714! computation of pressure gradient in 3D mode, which uses precisely
715! the initial value of "zeta" rather than the value changed by the
716! first barotropic predictor step. Later in this code, just after
717! "rufrc, rvfrc" are finalized, a correction term based on the
718! difference zeta_new(:,:)-zeta(:,:,kstp) to "rubar, rvbar" to make
719! them consistent with generalized RK2 stepping for pressure gradient
720! terms.
721!
722 IF (predictor_2d_step) THEN
723 IF (first_2d_step) THEN ! Modified RK2 time step (with
724 cff=dtfast(ng) ! Forward-Backward feedback with
725#ifdef SOLVE3D
726 cff1=0.0_r8 !==> Forward Euler
727 cff2=1.0_r8
728#else
729 cff1=0.333333333333_r8 ! optimally chosen beta=1/3 and
730 cff2=0.666666666667_r8 ! epsilon=2/3, see below) is used
731#endif
732 cff3=0.0_r8 ! here for the start up.
733 ELSE
734 cff=2.0_r8*dtfast(ng) ! In the code below "zwrk" is
735 cff1=beta ! time-centered at time step "n"
736 cff2=1.0_r8-2.0_r8*beta ! in the case of LF (for all but
737 cff3=beta ! the first time step)
738 END IF
739!
740 DO j=jstrv-1,jend
741 DO i=istru-1,iend
742 fac=cff*pm(i,j)*pn(i,j)
743 zeta_new(i,j)=zeta(i,j,kbak)+ &
744 & fac*(duon(i,j)-duon(i+1,j)+ &
745 & dvom(i,j)-dvom(i,j+1))
746#ifdef MASKING
747 zeta_new(i,j)=zeta_new(i,j)*rmask(i,j)
748# ifdef WET_DRY
749 zeta_new(i,j)=zeta_new(i,j)+ &
750 & (dcrit(ng)-h(i,j))*(1.0_r8-rmask(i,j))
751# endif
752#endif
753 dnew(i,j)=zeta_new(i,j)+h(i,j)
754
755 zwrk(i,j)=cff1*zeta_new(i,j)+ &
756 & cff2*zeta(i,j,kstp)+ &
757 & cff3*zeta(i,j,kbak)
758#if defined VAR_RHO_2D && defined SOLVE3D
759 rzeta(i,j)=(1.0_r8+rhos(i,j))*zwrk(i,j)
760 rzeta2(i,j)=rzeta(i,j)*zwrk(i,j)
761 rzetasa(i,j)=zwrk(i,j)*(rhos(i,j)-rhoa(i,j))
762#else
763 rzeta(i,j)=zwrk(i,j)
764 rzeta2(i,j)=zwrk(i,j)*zwrk(i,j)
765#endif
766 END DO
767 END DO
768 ELSE !--> CORRECTOR STEP
769 IF (first_2d_step) THEN
770 cff =0.333333333333_r8 ! Modified RK2 weighting:
771 cff1=0.333333333333_r8 ! here "zwrk" is time-
772 cff2=0.333333333333_r8 ! centered at "n+1/2".
773 cff3=0.0_r8
774 ELSE
775 cff =1.0_r8-epsil ! zwrk is always time-
776 cff1=(0.5_r8-gamma)*epsil ! centered at n+1/2
777 cff2=(0.5_r8+2.0_r8*gamma)*epsil ! during corrector sub-
778 cff3=-gamma *epsil ! step.
779 END IF
780!
781 DO j=jstrv-1,jend
782 DO i=istru-1,iend
783 fac=dtfast(ng)*pm(i,j)*pn(i,j)
784 zeta_new(i,j)=zeta(i,j,kstp)+ &
785 & fac*(duon(i,j)-duon(i+1,j)+ &
786 & dvom(i,j)-dvom(i,j+1))
787#ifdef MASKING
788 zeta_new(i,j)=zeta_new(i,j)*rmask(i,j)
789#endif
790 dnew(i,j)=zeta_new(i,j)+h(i,j)
791
792 zwrk(i,j)=cff *zeta(i,j,krhs)+ &
793 & cff1*zeta_new(i,j)+ &
794 & cff2*zeta(i,j,kstp)+ &
795 & cff3*zeta(i,j,kbak)
796#if defined VAR_RHO_2D && defined SOLVE3D
797 rzeta(i,j)=(1.0_r8+rhos(i,j))*zwrk(i,j)
798 rzeta2(i,j)=rzeta(i,j)*zwrk(i,j)
799 rzetasa(i,j)=zwrk(i,j)*(rhos(i,j)-rhoa(i,j))
800#else
801 rzeta(i,j)=zwrk(i,j)
802 rzeta2(i,j)=zwrk(i,j)*zwrk(i,j)
803#endif
804 END DO
805 END DO
806 END IF
807!
808! Apply mass point sources (volume vertical influx), if any.
809!
810! Dsrc(is) = 2, flow across grid cell w-face (positive or negative)
811!
812 IF (lwsrc(ng)) THEN
813 DO is=1,nsrc(ng)
814 IF (int(sources(ng)%Dsrc(is)).eq.2) THEN
815 i=sources(ng)%Isrc(is)
816 j=sources(ng)%Jsrc(is)
817 IF (((istrr.le.i).and.(i.le.iendr)).and. &
818 & ((jstrr.le.j).and.(j.le.jendr))) THEN
819 zeta_new(i,j)=zeta_new(i,j)+ &
820 & sources(ng)%Qbar(is)* &
821 & pm(i,j)*pn(i,j)*dtfast(ng)
822 END IF
823 END IF
824 END DO
825 END IF
826!
827! Apply boundary conditions to newly computed free-surface "zeta_new"
828! and load into global state array. Notice that "zeta_new" is always
829! centered at time step "m+1", while zeta(:,:,knew) should be centered
830! either at "m+1/2" after predictor step and at "m+1" after corrector.
831! Chosing it to be this way makes it possible avoid storing RHS for
832! zeta, ubar, and vbar between predictor and corrector sub-steps.
833!
834! Here, we use the local "zetabc" since the private array "zeta_new"
835! is passed as an argument to allow computing the lateral boundary
836! conditions on the range IstrU-1:Iend and JstrV-1:Jend, so parallel
837! tile exchanges are avoided.
838!
839 CALL zetabc_local (ng, tile, &
840 & lbi, ubi, lbj, ubj, &
841 & imins, imaxs, jmins, jmaxs, &
842 & kstp, &
843 & zeta, &
844 & zeta_new)
845!
846 IF (predictor_2d_step) THEN
847 IF (first_2d_step) THEN
848 cff1=0.5_r8
849 cff2=0.5_r8
850 cff3=0.0_r8
851 ELSE
852 cff1=0.5_r8-gamma
853 cff2=0.5_r8+2.0_r8*gamma
854 cff3=-gamma
855 END IF
856 DO j=jstrr,jendr
857 DO i=istrr,iendr
858 zeta(i,j,knew)=cff1*zeta_new(i,j)+ &
859 & cff2*zeta(i,j,kstp)+ &
860 & cff3*zeta(i,j,kbak)
861 END DO
862 END DO
863 ELSE
864 DO j=jstrr,jendr
865 DO i=istrr,iendr
866 zeta(i,j,knew)=zeta_new(i,j)
867 END DO
868 END DO
869 END IF
870!
871!=======================================================================
872! Compute right-hand-side for the 2D momentum equations.
873!=======================================================================
874#ifdef SOLVE3D
875!
876! Notice that we are suppressing the computation of momentum advection,
877! Coriolis, and lateral viscosity terms in 3D Applications because
878! these terms are already included in the baroclinic-to-barotropic
879! forcing arrays "rufrc" and "rvfrc". It does not mean we are entirely
880! omitting them, but it is a choice between recomputing them at every
881! barotropic step or keeping them "frozen" during the fast-time
882! stepping.
883# ifdef STEP2D_CORIOLIS
884! However, in some coarse grid applications with larger baroclinic
885! timestep (say, DT around 20 minutes or larger), adding the Coriolis
886! term in the barotropic equations is useful since f*DT is no longer
887! small.
888# endif
889#endif
890!
891!-----------------------------------------------------------------------
892! Compute pressure-gradient terms.
893!-----------------------------------------------------------------------
894!
895! Notice that "rubar" and "rvbar" are computed within the same to allow
896! shared references to array elements (i,j), which increases the
897! computational density by almost a factor of 1.5 resulting in overall
898! more efficient code.
899!
900 cff1=0.5*g
901 cff2=0.333333333333_r8
902#if !defined SOLVE3D && defined ATM_PRESS
903 fac=0.5_r8*100.0_r8/rho0
904#endif
905 DO j=jstr,jend
906 DO i=istr,iend
907 IF (i.ge.istru) THEN
908 rubar(i,j)=cff1*on_u(i,j)* &
909 & ((h(i-1,j)+ &
910 & h(i ,j))* &
911 & (rzeta(i-1,j)- &
912 & rzeta(i ,j))+ &
913#if defined VAR_RHO_2D && defined SOLVE3D
914 & (h(i-1,j)- &
915 & h(i ,j))* &
916 & (rzetasa(i-1,j)+ &
917 & rzetasa(i ,j)+ &
918 & cff2*(rhoa(i-1,j)- &
919 & rhoa(i ,j))* &
920 & (zwrk(i-1,j)- &
921 & zwrk(i ,j)))+ &
922#endif
923 & (rzeta2(i-1,j)- &
924 & rzeta2(i ,j)))
925#if defined ATM_PRESS && !defined SOLVE3D
926 rubar(i,j)=rubar(i,j)- &
927 & fac*on_u(i,j)* &
928 & (h(i-1,j)+h(i,j)+ &
929 & rzeta(i-1,j)+rzeta(i,j))* &
930 & (pair(i,j)-pair(i-1,j))
931#endif
932#if defined TIDE_GENERATING_FORCES && !defined SOLVE3D
933 rubar(i,j)=rubar(i,j)- &
934 & cff1*on_u(i,j)* &
935 & (h(i-1,j)+h(i,j)+ &
936 & rzeta(i-1,j)+rzeta(i,j))* &
937 & (eq_tide(i,j)-eq_tide(i-1,j))
938#endif
939#ifdef DIAGNOSTICS_UV
940 diau2rhs(i,j,m2pgrd)=rubar(i,j)
941#endif
942 END IF
943!
944 IF (j.ge.jstrv) THEN
945 rvbar(i,j)=cff1*om_v(i,j)* &
946 & ((h(i,j-1)+ &
947 & h(i,j ))* &
948 & (rzeta(i,j-1)- &
949 & rzeta(i,j ))+ &
950#if defined VAR_RHO_2D && defined SOLVE3D
951 & (h(i,j-1)- &
952 & h(i,j ))* &
953 & (rzetasa(i,j-1)+ &
954 & rzetasa(i,j )+ &
955 & cff2*(rhoa(i,j-1)- &
956 & rhoa(i,j ))* &
957 & (zwrk(i,j-1)- &
958 & zwrk(i,j )))+ &
959#endif
960 & (rzeta2(i,j-1)- &
961 & rzeta2(i,j )))
962#if defined ATM_PRESS && !defined SOLVE3D
963 rvbar(i,j)=rvbar(i,j)- &
964 & fac*om_v(i,j)* &
965 & (h(i,j-1)+h(i,j)+ &
966 & rzeta(i,j-1)+rzeta(i,j))* &
967 & (pair(i,j)-pair(i,j-1))
968#endif
969#if defined TIDE_GENERATING_FORCES && !defined SOLVE3D
970 rvbar(i,j)=rvbar(i,j)- &
971 & cff1*om_v(i,j)* &
972 & (h(i,j-1)+h(i,j)+ &
973 & rzeta(i,j-1)+rzeta(i,j))* &
974 & (eq_tide(i,j)-eq_tide(i,j-1))
975#endif
976#ifdef DIAGNOSTICS_UV
977 diav2rhs(i,j,m2pgrd)=rvbar(i,j)
978#endif
979 END IF
980 END DO
981 END DO
982
983#if defined UV_ADV && !defined SOLVE3D
984!
985!-----------------------------------------------------------------------
986! Add in horizontal advection of momentum.
987!-----------------------------------------------------------------------
988
989# ifdef UV_C2ADVECTION
990!
991! Second-order, centered differences advection fluxes.
992!
993 DO j=jstr,jend
994 DO i=istru-1,iend
995 ufx(i,j)=0.25_r8* &
996 & (duon(i,j)+duon(i+1,j))* &
997 & (ubar(i ,j,krhs)+ &
998 & ubar(i+1,j,krhs))
999 END DO
1000 END DO
1001!
1002 DO j=jstr,jend+1
1003 DO i=istru,iend
1004 ufe(i,j)=0.25_r8* &
1005 & (dvom(i,j)+dvom(i-1,j))* &
1006 & (ubar(i,j ,krhs)+ &
1007 & ubar(i,j-1,krhs))
1008 END DO
1009 END DO
1010!
1011 DO j=jstrv,jend
1012 DO i=istr,iend+1
1013 vfx(i,j)=0.25_r8* &
1014 & (duon(i,j)+duon(i,j-1))* &
1015 & (vbar(i ,j,krhs)+ &
1016 & vbar(i-1,j,krhs))
1017 END DO
1018 END DO
1019!
1020 DO j=jstrv-1,jend
1021 DO i=istr,iend
1022 vfe(i,j)=0.25_r8* &
1023 & (dvom(i,j)+dvom(i,j+1))* &
1024 & (vbar(i,j ,krhs)+ &
1025 & vbar(i,j+1,krhs))
1026 END DO
1027 END DO
1028
1029# elif defined UV_C4ADVECTION
1030!
1031! Fourth-order, centered differences u-momentum advection fluxes.
1032!
1033 DO j=jstr,jend
1034 DO i=istrum1,iendp1
1035 grad(i,j)=ubar(i-1,j,krhs)-2.0_r8*ubar(i,j,krhs)+ &
1036 & ubar(i+1,j,krhs)
1037 dgrad(i,j)=duon(i-1,j)-2.0_r8*duon(i,j)+duon(i+1,j)
1038 END DO
1039 END DO
1040 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
1041 IF (domain(ng)%Western_Edge(tile)) THEN
1042 DO j=jstr,jend
1043 grad(istr,j)=grad(istr+1,j)
1044 dgrad(istr,j)=dgrad(istr+1,j)
1045 END DO
1046 END IF
1047 END IF
1048 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
1049 IF (domain(ng)%Eastern_Edge(tile)) THEN
1050 DO j=jstr,jend
1051 grad(iend+1,j)=grad(iend,j)
1052 dgrad(iend+1,j)=dgrad(iend,j)
1053 END DO
1054 END IF
1055 END IF
1056! d/dx(Duu/n)
1057 cff=1.0_r8/6.0_r8
1058 DO j=jstr,jend
1059 DO i=istru-1,iend
1060 ufx(i,j)=0.25_r8*(ubar(i ,j,krhs)+ &
1061 & ubar(i+1,j,krhs)- &
1062 & cff*(grad(i,j)+grad(i+1,j)))* &
1063 & (duon(i,j)+duon(i+1,j)- &
1064 & cff*(dgrad(i,j)+dgrad(i+1,j)))
1065 END DO
1066 END DO
1067!
1068 DO j=jstrm1,jendp1
1069 DO i=istru,iend
1070 grad(i,j)=ubar(i,j-1,krhs)-2.0_r8*ubar(i,j,krhs)+ &
1071 & ubar(i,j+1,krhs)
1072 END DO
1073 END DO
1074!
1075 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
1076 IF (domain(ng)%Southern_Edge(tile)) THEN
1077 DO i=istru,iend
1078 grad(i,jstr-1)=grad(i,jstr)
1079 END DO
1080 END IF
1081 END IF
1082 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
1083 IF (domain(ng)%Northern_Edge(tile)) THEN
1084 DO i=istru,iend
1085 grad(i,jend+1)=grad(i,jend)
1086 END DO
1087 END IF
1088 END IF
1089 DO j=jstr,jend+1
1090 DO i=istru-1,iend
1091 dgrad(i,j)=dvom(i-1,j)-2.0_r8*dvom(i,j)+dvom(i+1,j)
1092 END DO
1093 END DO
1094! d/dy(Duv/m)
1095 cff=1.0_r8/6.0_r8
1096 DO j=jstr,jend+1
1097 DO i=istru,iend
1098 ufe(i,j)=0.25_r8*(ubar(i,j ,krhs)+ &
1099 & ubar(i,j-1,krhs)- &
1100 & cff*(grad(i,j)+grad(i,j-1)))* &
1101 & (dvom(i,j)+dvom(i-1,j)- &
1102 & cff*(dgrad(i,j)+dgrad(i-1,j)))
1103 END DO
1104 END DO
1105!
1106! Fourth-order, centered differences v-momentum advection fluxes.
1107!
1108 DO j=jstrv,jend
1109 DO i=istrm1,iendp1
1110 grad(i,j)=vbar(i-1,j,krhs)-2.0_r8*vbar(i,j,krhs)+ &
1111 & vbar(i+1,j,krhs)
1112 END DO
1113 END DO
1114!
1115 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
1116 IF (domain(ng)%Western_Edge(tile)) THEN
1117 DO j=jstrv,jend
1118 grad(istr-1,j)=grad(istr,j)
1119 END DO
1120 END IF
1121 END IF
1122 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
1123 IF (domain(ng)%Eastern_Edge(tile)) THEN
1124 DO j=jstrv,jend
1125 grad(iend+1,j)=grad(iend,j)
1126 END DO
1127 END IF
1128 END IF
1129 DO j=jstrv-1,jend
1130 DO i=istr,iend+1
1131 dgrad(i,j)=duon(i,j-1)-2.0_r8*duon(i,j)+duon(i,j+1)
1132 END DO
1133 END DO
1134! d/dx(Duv/n)
1135 cff=1.0_r8/6.0_r8
1136 DO j=jstrv,jend
1137 DO i=istr,iend+1
1138 vfx(i,j)=0.25_r8*(vbar(i ,j,krhs)+ &
1139 & vbar(i-1,j,krhs)- &
1140 & cff*(grad(i,j)+grad(i-1,j)))* &
1141 & (duon(i,j)+duon(i,j-1)- &
1142 & cff*(dgrad(i,j)+dgrad(i,j-1)))
1143 END DO
1144 END DO
1145!
1146 DO j=jstrvm1,jendp1
1147 DO i=istr,iend
1148 grad(i,j)=vbar(i,j-1,krhs)-2.0_r8*vbar(i,j,krhs)+ &
1149 & vbar(i,j+1,krhs)
1150 dgrad(i,j)=dvom(i,j-1)-2.0_r8*dvom(i,j)+dvom(i,j+1)
1151 END DO
1152 END DO
1153 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
1154 IF (domain(ng)%Southern_Edge(tile)) THEN
1155 DO i=istr,iend
1156 grad(i,jstr)=grad(i,jstr+1)
1157 dgrad(i,jstr)=dgrad(i,jstr+1)
1158 END DO
1159 END IF
1160 END IF
1161 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
1162 IF (domain(ng)%Northern_Edge(tile)) THEN
1163 DO i=istr,iend
1164 grad(i,jend+1)=grad(i,jend)
1165 dgrad(i,jend+1)=dgrad(i,jend)
1166 END DO
1167 END IF
1168 END IF
1169! d/dy(Dvv/m)
1170 cff=1.0_r8/6.0_r8
1171 DO j=jstrv-1,jend
1172 DO i=istr,iend
1173 vfe(i,j)=0.25_r8*(vbar(i,j ,krhs)+ &
1174 & vbar(i,j+1,krhs)- &
1175 & cff*(grad(i,j)+grad(i,j+1)))* &
1176 & (dvom(i,j)+dvom(i,j+1)- &
1177 & cff*(dgrad(i,j)+dgrad(i,j+1)))
1178 END DO
1179 END DO
1180# endif
1181!
1182! Add advection to RHS terms.
1183!
1184 DO j=jstr,jend
1185 DO i=istr,iend
1186 IF (i.ge.istru) THEN
1187 cff1=ufx(i,j)-ufx(i-1,j)
1188 cff2=ufe(i,j+1)-ufe(i,j)
1189 fac=cff1+cff2
1190 rubar(i,j)=rubar(i,j)-fac
1191# if defined DIAGNOSTICS_UV
1192 diau2rhs(i,j,m2xadv)=-cff1
1193 diau2rhs(i,j,m2yadv)=-cff2
1194 diau2rhs(i,j,m2hadv)=-fac
1195# endif
1196 END IF
1197!
1198 IF (j.ge.jstrv) THEN
1199 cff1=vfx(i+1,j)-vfx(i,j)
1200 cff2=vfe(i,j)-vfe(i,j-1)
1201 fac=cff1+cff2
1202 rvbar(i,j)=rvbar(i,j)-fac
1203# if defined DIAGNOSTICS_UV
1204 diav2rhs(i,j,m2xadv)=-cff1
1205 diav2rhs(i,j,m2yadv)=-cff2
1206 diav2rhs(i,j,m2hadv)=-fac
1207# endif
1208 END IF
1209 END DO
1210 END DO
1211#endif
1212
1213#if (defined UV_COR & !defined SOLVE3D) || defined STEP2D_CORIOLIS
1214!
1215!-----------------------------------------------------------------------
1216! Add in Coriolis term.
1217!-----------------------------------------------------------------------
1218!
1219 DO j=jstrv-1,jend
1220 DO i=istru-1,iend
1221 cff=0.5_r8*drhs(i,j)*fomn(i,j)
1222 ufx(i,j)=cff*(vbar(i,j ,krhs)+ &
1223 & vbar(i,j+1,krhs))
1224 vfe(i,j)=cff*(ubar(i ,j,krhs)+ &
1225 & ubar(i+1,j,krhs))
1226 END DO
1227 END DO
1228!
1229 DO j=jstr,jend
1230 DO i=istr,iend
1231 IF (i.ge.istru) THEN
1232 fac1=0.5_r8*(ufx(i,j)+ufx(i-1,j))
1233 rubar(i,j)=rubar(i,j)+fac1
1234# if defined DIAGNOSTICS_UV
1235 diau2rhs(i,j,m2fcor)=fac1
1236# endif
1237 END IF
1238!
1239 IF (j.ge.jstrv) THEN
1240 fac2=0.5_r8*(vfe(i,j)+vfe(i,j-1))
1241 rvbar(i,j)=rvbar(i,j)-fac2
1242# if defined DIAGNOSTICS_UV
1243 diav2rhs(i,j,m2fcor)=-fac2
1244# endif
1245 END IF
1246 END DO
1247 END DO
1248#endif
1249
1250#if (defined CURVGRID && defined UV_ADV) && !defined SOLVE3D
1251!
1252!-----------------------------------------------------------------------
1253! Add in curvilinear transformation terms.
1254!-----------------------------------------------------------------------
1255!
1256 DO j=jstrv-1,jend
1257 DO i=istru-1,iend
1258 cff1=0.5_r8*(vbar(i,j ,krhs)+ &
1259 & vbar(i,j+1,krhs))
1260 cff2=0.5_r8*(ubar(i ,j,krhs)+ &
1261 & ubar(i+1,j,krhs))
1262 cff3=cff1*dndx(i,j)
1263 cff4=cff2*dmde(i,j)
1264 cff=drhs(i,j)*(cff3-cff4)
1265 ufx(i,j)=cff*cff1
1266 vfe(i,j)=cff*cff2
1267# if defined DIAGNOSTICS_UV
1268 cff=drhs(i,j)*cff4
1269 uwrk(i,j)=-cff*cff1 ! ubar equation, ETA-term
1270 vwrk(i,j)=-cff*cff2 ! vbar equation, ETA-term
1271# endif
1272 END DO
1273 END DO
1274!
1275 DO j=jstr,jend
1276 DO i=istr,iend
1277 IF (i.ge.istru) THEN
1278 fac1=0.5_r8*(ufx(i,j)+ufx(i-1,j))
1279 rubar(i,j)=rubar(i,j)+fac1
1280# if defined DIAGNOSTICS_UV
1281 fac2=0.5_r8*(uwrk(i,j)+uwrk(i-1,j))
1282 diau2rhs(i,j,m2xadv)=diau2rhs(i,j,m2xadv)+fac1-fac2
1283 diau2rhs(i,j,m2yadv)=diau2rhs(i,j,m2yadv)+fac2
1284 diau2rhs(i,j,m2hadv)=diau2rhs(i,j,m2hadv)+fac1
1285# endif
1286 END IF
1287!
1288 IF (j.ge.jstrv) THEN
1289 fac1=0.5_r8*(vfe(i,j)+vfe(i,j-1))
1290 rvbar(i,j)=rvbar(i,j)-fac1
1291# if defined DIAGNOSTICS_UV
1292 fac2=0.5_r8*(vwrk(i,j)+vwrk(i,j-1))
1293 diav2rhs(i,j,m2xadv)=diav2rhs(i,j,m2xadv)-fac1+fac2
1294 diav2rhs(i,j,m2yadv)=diav2rhs(i,j,m2yadv)-fac2
1295 diav2rhs(i,j,m2hadv)=diav2rhs(i,j,m2hadv)-fac1
1296# endif
1297 END IF
1298 END DO
1299 END DO
1300#endif
1301
1302#if defined UV_VIS2 && !defined SOLVE3D
1303!
1304!-----------------------------------------------------------------------
1305! Add in horizontal harmonic viscosity.
1306!-----------------------------------------------------------------------
1307!
1308! Compute total depth at PSI-points.
1309!
1310 DO j=jstr,jend+1
1311 DO i=istr,iend+1
1312 drhs_p(i,j)=0.25_r8*(drhs(i,j )+drhs(i-1,j )+ &
1313 & drhs(i,j-1)+drhs(i-1,j-1))
1314 END DO
1315 END DO
1316!
1317! Compute flux-components of the horizontal divergence of the stress
1318! tensor (m5/s2) in XI- and ETA-directions.
1319!
1320 DO j=jstrv-1,jend
1321 DO i=istru-1,iend
1322 cff=visc2_r(i,j)*drhs(i,j)*0.5_r8* &
1323 & (pmon_r(i,j)* &
1324 & ((pn(i ,j)+pn(i+1,j))*ubar(i+1,j,krhs)- &
1325 & (pn(i-1,j)+pn(i ,j))*ubar(i ,j,krhs))- &
1326 & pnom_r(i,j)* &
1327 & ((pm(i,j )+pm(i,j+1))*vbar(i,j+1,krhs)- &
1328 & (pm(i,j-1)+pm(i,j ))*vbar(i,j ,krhs)))
1329 ufx(i,j)=on_r(i,j)*on_r(i,j)*cff
1330 vfe(i,j)=om_r(i,j)*om_r(i,j)*cff
1331 END DO
1332 END DO
1333!
1334 DO j=jstr,jend+1
1335 DO i=istr,iend+1
1336 cff=visc2_p(i,j)*drhs_p(i,j)*0.5_r8* &
1337 & (pmon_p(i,j)* &
1338 & ((pn(i ,j-1)+pn(i ,j))*vbar(i ,j,krhs)- &
1339 & (pn(i-1,j-1)+pn(i-1,j))*vbar(i-1,j,krhs))+ &
1340 & pnom_p(i,j)* &
1341 & ((pm(i-1,j )+pm(i,j ))*ubar(i,j ,krhs)- &
1342 & (pm(i-1,j-1)+pm(i,j-1))*ubar(i,j-1,krhs)))
1343# ifdef MASKING
1344 cff=cff*pmask(i,j)
1345# endif
1346# ifdef WET_DRY
1347 cff=cff*pmask_wet(i,j)
1348# endif
1349 ufe(i,j)=om_p(i,j)*om_p(i,j)*cff
1350 vfx(i,j)=on_p(i,j)*on_p(i,j)*cff
1351 END DO
1352 END DO
1353!
1354! Add in harmonic viscosity.
1355!
1356 DO j=jstr,jend
1357 DO i=istr,iend
1358 IF (i.ge.istru) THEN
1359 cff1=0.5_r8*(pn(i-1,j)+pn(i,j))*(ufx(i,j )-ufx(i-1,j))
1360 cff2=0.5_r8*(pm(i-1,j)+pm(i,j))*(ufe(i,j+1)-ufe(i ,j))
1361 fac=cff1+cff2
1362 rubar(i,j)=rubar(i,j)+fac
1363# if defined DIAGNOSTICS_UV
1364 diau2rhs(i,j,m2hvis)=fac
1365 diau2rhs(i,j,m2xvis)=cff1
1366 diau2rhs(i,j,m2yvis)=cff2
1367# endif
1368 END IF
1369!
1370 IF (j.ge.jstrv) THEN
1371 cff1=0.5_r8*(pn(i,j-1)+pn(i,j))*(vfx(i+1,j)-vfx(i,j ))
1372 cff2=0.5_r8*(pm(i,j-1)+pm(i,j))*(vfe(i ,j)-vfe(i,j-1))
1373 fac=cff1-cff2
1374 rvbar(i,j)=rvbar(i,j)+fac
1375# if defined DIAGNOSTICS_UV
1376 diav2rhs(i,j,m2hvis)=fac
1377 diav2rhs(i,j,m2xvis)= cff1
1378 diav2rhs(i,j,m2yvis)=-cff2
1379# endif
1380 END IF
1381 END DO
1382 END DO
1383#endif
1384
1385#ifndef SOLVE3D
1386!
1387!-----------------------------------------------------------------------
1388! Add in bottom stress.
1389!-----------------------------------------------------------------------
1390!
1391 DO j=jstr,jend
1392 DO i=istru,iend
1393 fac=bustr(i,j)*om_u(i,j)*on_u(i,j)
1394 rubar(i,j)=rubar(i,j)-fac
1395# ifdef DIAGNOSTICS_UV
1396 diau2rhs(i,j,m2bstr)=-fac
1397# endif
1398 END DO
1399 END DO
1400 DO j=jstrv,jend
1401 DO i=istr,iend
1402 fac=bvstr(i,j)*om_v(i,j)*on_v(i,j)
1403 rvbar(i,j)=rvbar(i,j)-fac
1404# ifdef DIAGNOSTICS_UV
1405 diav2rhs(i,j,m2bstr)=-fac
1406# endif
1407 END DO
1408 END DO
1409#else
1410# ifdef DIAGNOSTICS_UV
1411!
1412! Initialize the stress term if no bottom friction is defined.
1413!
1414 DO j=jstr,jend
1415 DO i=istru,iend
1416 diau2rhs(i,j,m2bstr)=0.0_r8
1417 END DO
1418 END DO
1419 DO j=jstrv,jend
1420 DO i=istr,iend
1421 diav2rhs(i,j,m2bstr)=0.0_r8
1422 END DO
1423 END DO
1424# endif
1425#endif
1426
1427#ifdef SOLVE3D
1428!
1429!-----------------------------------------------------------------------
1430! Coupling between 2D and 3D equations.
1431!-----------------------------------------------------------------------
1432!
1433! Before the predictor step of the first barotropic time step, arrays
1434! "rufrc" and "rvfrc" contain vertical integrals of the 3D RHS terms
1435! for the momentum equations (including surface and bottom stresses,
1436! if so prescribed). During the first barotropic time step, convert
1437! them into forcing terms by subtracting the fast-time "rubar" and
1438! "rvbar" from them.
1439!
1440! These forcing terms are then extrapolated forward in time using
1441! optimized Adams-Bashforth weights, so that the resultant "rufrc"
1442! and "rvfrc" are centered effectively at time n+1/2 in baroclinic
1443! time.
1444!
1445! From now on, these newly computed forcing terms remain unchanged
1446! during the fast time stepping and will be added to "rubar" and
1447! "rvbar" during all subsequent barotropic time steps.
1448!
1449! Thus, the algorithm below is designed for coupling during the 3D
1450! predictor sub-step. The forcing terms "rufrc" and "rvfrc" are
1451! computed as instantaneous values at 3D time index "nstp" first and
1452! then extrapolated half-step forward using AM3-like weights optimized
1453! for maximum stability (with particular care for startup).
1454!
1455 IF (first_2d_step.and.predictor_2d_step) THEN
1456 IF (first_time_step) THEN
1457 cff3=0.0_r8
1458 cff2=0.0_r8
1459 cff1=1.0_r8
1460 ELSE IF (first_time_step+1) THEN
1461 cff3=0.0_r8
1462 cff2=-0.5_r8
1463 cff1=1.5_r8
1464 ELSE
1465 cff3=0.281105_r8
1466 cff2=-0.5_r8-2.0_r8*cff3
1467 cff1=1.5_r8+cff3
1468 END IF
1469!
1470 DO j=jstr,jend
1471 DO i=istru,iend
1472 cff=rufrc(i,j)-rubar(i,j)
1473 rufrc(i,j)=cff1*cff+ &
1474 & cff2*rufrc_bak(i,j,3-nstp)+ &
1475 & cff3*rufrc_bak(i,j,nstp )
1476 rufrc_bak(i,j,nstp)=cff
1477 END DO
1478 END DO
1479 DO j=jstrv,jend
1480 DO i=istr,iend
1481 cff=rvfrc(i,j)-rvbar(i,j)
1482 rvfrc(i,j)=cff1*cff+ &
1483 & cff2*rvfrc_bak(i,j,3-nstp)+ &
1484 & cff3*rvfrc_bak(i,j,nstp )
1485 rvfrc_bak(i,j,nstp)=cff
1486 END DO
1487 END DO
1488!
1489! Since coupling requires that the pressure gradient term is computed
1490! using zeta(:,:,kstp) instead of 1/3 toward zeta_new(:,:) as needed
1491! by generalized RK2 scheme, apply compensation to shift pressure
1492! gradient terms from "kstp" to 1/3 toward "knew".
1493!
1494 cff1=0.5_r8*g
1495 cff2=0.333333333333_r8
1496 cff3=1.666666666666_r8
1497
1498 DO j=jstrv-1,jend
1499 DO i=istru-1,iend
1500 zwrk(i,j)=cff2*(zeta_new(i,j)-zeta(i,j,kstp))
1501# if defined VAR_RHO_2D && defined SOLVE3D
1502 rzeta(i,j)=(1.0_r8+rhos(i,j))*zwrk(i,j)
1503 rzeta2(i,j)=rzeta(i,j)* &
1504 & (cff2*zeta_new(i,j)+ &
1505 & cff3*zeta(i,j,kstp))
1506 rzetasa(i,j)=zwrk(i,j)*(rhos(i,j)-rhoa(i,j))
1507# else
1508 rzeta(i,j)=zwrk(i,j)
1509 rzeta2(i,j)=zwrk(i,j)* &
1510 & (cff2*zeta_new(i,j)+ &
1511 & cff3*zeta(i,j,kstp))
1512# endif
1513 END DO
1514 END DO
1515!
1516 DO j=jstr,jend
1517 DO i=istr,iend
1518 IF (i.ge.istru) THEN
1519 rubar(i,j)=rubar(i,j)+ &
1520 & cff1*on_u(i,j)* &
1521 & ((h(i-1,j)+ &
1522 & h(i ,j))* &
1523 & (rzeta(i-1,j)- &
1524 & rzeta(i ,j))+ &
1525# if defined VAR_RHO_2D && defined SOLVE3D
1526 & (h(i-1,j)- &
1527 & h(i ,j))* &
1528 & (rzetasa(i-1,j)+ &
1529 & rzetasa(i ,j)+ &
1530 & cff2*(rhoa(i-1,j)- &
1531 & rhoa(i ,j))* &
1532 & (zwrk(i-1,j)- &
1533 & zwrk(i ,j)))+ &
1534# endif
1535 & (rzeta2(i-1,j)- &
1536 & rzeta2(i ,j)))
1537# ifdef DIAGNOSTICS_UV
1538 diau2rhs(i,j,m2pgrd)=diau2rhs(i,j,m2pgrd)+ &
1539 & rubar(i,j)
1540# endif
1541 END IF
1542!
1543 IF (j.ge.jstrv) THEN
1544 rvbar(i,j)=rvbar(i,j)+ &
1545 & cff1*om_v(i,j)* &
1546 & ((h(i,j-1)+ &
1547 & h(i,j ))* &
1548 & (rzeta(i,j-1)- &
1549 & rzeta(i,j ))+ &
1550# if defined VAR_RHO_2D && defined SOLVE3D
1551 & (h(i,j-1)- &
1552 & h(i,j ))* &
1553 & (rzetasa(i,j-1)+ &
1554 & rzetasa(i,j )+ &
1555 & cff2*(rhoa(i,j-1)- &
1556 & rhoa(i,j ))* &
1557 & (zwrk(i,j-1)- &
1558 & zwrk(i,j )))+ &
1559# endif
1560 & (rzeta2(i,j-1)- &
1561 & rzeta2(i,j )))
1562# ifdef DIAGNOSTICS_UV
1563 diav2rhs(i,j,m2pgrd)=diav2rhs(i,j,m2pgrd)+ &
1564 & rvbar(i,j)
1565# endif
1566 END IF
1567 END DO
1568 END DO
1569 END IF
1570#endif
1571!
1572!=======================================================================
1573! Time step 2D momentum equations.
1574!=======================================================================
1575!
1576! Compute total water column depth.
1577!
1578 IF (first_2d_step.or.(.not.predictor_2d_step)) THEN
1579 DO j=jstrv-1,jend
1580 DO i=istru-1,iend
1581 dstp(i,j)=h(i,j)+zeta(i,j,kstp)
1582 END DO
1583 END DO
1584 ELSE
1585 DO j=jstrv-1,jend
1586 DO i=istru-1,iend
1587 dstp(i,j)=h(i,j)+zeta(i,j,kbak)
1588 END DO
1589 END DO
1590 END IF
1591!
1592! During the predictor sub-step, once newly computed "ubar" and "vbar"
1593! become available, interpolate them half-step backward in barotropic
1594! time (i.e., they end up time-centered at n+1/2) in order to use it
1595! during subsequent corrector sub-step.
1596!
1597 IF (predictor_2d_step) THEN
1598 IF (first_2d_step) THEN
1599 cff1=0.5_r8*dtfast(ng)
1600 cff2=0.5_r8
1601 cff3=0.5_r8
1602 cff4=0.0_r8
1603 ELSE
1604 cff1=dtfast(ng)
1605 cff2=0.5_r8-gamma
1606 cff3=0.5_r8+2.0_r8*gamma
1607 cff4=-gamma
1608 ENDIF
1609!
1610 DO j=jstr,jend
1611 DO i=istru,iend
1612 cff=cff1*(pm(i,j)+pm(i-1,j))*(pn(i,j)+pn(i-1,j))
1613 fac1=1.0_r8/(dnew(i,j)+dnew(i-1,j))
1614 ubar(i,j,knew)=fac1* &
1615 & (ubar(i,j,kbak)* &
1616 & (dstp(i,j)+dstp(i-1,j))+ &
1617#ifdef SOLVE3D
1618 & cff*(rubar(i,j)+rufrc(i,j)))
1619#else
1620 & cff*rubar(i,j)+4.0_r8*cff1*sustr(i,j))
1621#endif
1622#ifdef MASKING
1623 ubar(i,j,knew)=ubar(i,j,knew)*umask(i,j)
1624#endif
1625 ubar(i,j,knew)=cff2*ubar(i,j,knew)+ &
1626 & cff3*ubar(i,j,kstp)+ &
1627 & cff4*ubar(i,j,kbak)
1628#ifdef WET_DRY
1629 cff5=abs(abs(umask_wet(i,j))-1.0_r8)
1630 cff6=0.5_r8+dsign(0.5_r8,ubar(i,j,knew))*umask_wet(i,j)
1631 cff7=0.5_r8*umask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
1632 ubar(i,j,knew)=ubar(i,j,knew)*cff7
1633#endif
1634#if defined NESTING && !defined SOLVE3D
1635 du_flux(i,j)=0.5_r8*on_u(i,j)* &
1636 & (dnew(i,j)+dnew(i-1,j))*ubar(i,j,knew)
1637#endif
1638 END DO
1639 END DO
1640!
1641 DO j=jstrv,jend
1642 DO i=istr,iend
1643 cff=cff1*(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1))
1644 fac2=1.0_r8/(dnew(i,j)+dnew(i,j-1))
1645 vbar(i,j,knew)=fac2* &
1646 & (vbar(i,j,kbak)* &
1647 & (dstp(i,j)+dstp(i,j-1))+ &
1648#ifdef SOLVE3D
1649 & cff*(rvbar(i,j)+rvfrc(i,j)))
1650#else
1651 & cff*rvbar(i,j)+4.0_r8*cff1*svstr(i,j))
1652#endif
1653#ifdef MASKING
1654 vbar(i,j,knew)=vbar(i,j,knew)*vmask(i,j)
1655#endif
1656 vbar(i,j,knew)=cff2*vbar(i,j,knew)+ &
1657 & cff3*vbar(i,j,kstp)+ &
1658 & cff4*vbar(i,j,kbak)
1659#ifdef WET_DRY
1660 cff5=abs(abs(vmask_wet(i,j))-1.0_r8)
1661 cff6=0.5_r8+dsign(0.5_r8,vbar(i,j,knew))*vmask_wet(i,j)
1662 cff7=0.5_r8*vmask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
1663 vbar(i,j,knew)=vbar(i,j,knew)*cff7
1664#endif
1665#if defined NESTING && !defined SOLVE3D
1666 dv_flux(i,j)=0.5_r8*om_v(i,j)* &
1667 & (dnew(i,j)+dnew(i,j-1))*vbar(i,j,knew)
1668#endif
1669 END DO
1670 END DO
1671
1672 ELSE !--> CORRECTOR_2D_STEP
1673
1674 cff1=0.5_r8*dtfast(ng)
1675 DO j=jstr,jend
1676 DO i=istru,iend
1677 cff=cff1*(pm(i,j)+pm(i-1,j))*(pn(i,j)+pn(i-1,j))
1678 fac1=1.0_r8/(dnew(i,j)+dnew(i-1,j))
1679 ubar(i,j,knew)=fac1* &
1680 & (ubar(i,j,kstp)* &
1681 & (dstp(i,j)+dstp(i-1,j))+ &
1682#ifdef SOLVE3D
1683 & cff*(rubar(i,j)+rufrc(i,j)))
1684#else
1685 & cff*rubar(i,j)+4.0_r8*cff1*sustr(i,j))
1686#endif
1687#ifdef MASKING
1688 ubar(i,j,knew)=ubar(i,j,knew)*umask(i,j)
1689#endif
1690#ifdef WET_DRY
1691 cff5=abs(abs(umask_wet(i,j))-1.0_r8)
1692 cff6=0.5_r8+dsign(0.5_r8,ubar(i,j,knew))*umask_wet(i,j)
1693 cff7=0.5_r8*umask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
1694 ubar(i,j,knew)=ubar(i,j,knew)*cff7
1695#endif
1696#if defined NESTING && !defined SOLVE3D
1697 du_flux(i,j)=0.5_r8*on_u(i,j)* &
1698 & (dnew(i,j)+dnew(i-1,j))*ubar(i,j,knew)
1699#endif
1700 END DO
1701 END DO
1702!
1703 DO j=jstrv,jend
1704 DO i=istr,iend
1705 cff=cff1*(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1))
1706 fac2=1.0_r8/(dnew(i,j)+dnew(i,j-1))
1707 vbar(i,j,knew)=fac2* &
1708 & (vbar(i,j,kstp)* &
1709 & (dstp(i,j)+dstp(i,j-1))+ &
1710#ifdef SOLVE3D
1711 & cff*(rvbar(i,j)+rvfrc(i,j)))
1712#else
1713 & cff*rvbar(i,j)+4.0_r8*cff1*svstr(i,j))
1714#endif
1715#ifdef MASKING
1716 vbar(i,j,knew)=vbar(i,j,knew)*vmask(i,j)
1717#endif
1718#ifdef WET_DRY
1719 cff5=abs(abs(vmask_wet(i,j))-1.0_r8)
1720 cff6=0.5_r8+dsign(0.5_r8,vbar(i,j,knew))*vmask_wet(i,j)
1721 cff7=0.5_r8*vmask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
1722 vbar(i,j,knew)=vbar(i,j,knew)*cff7
1723#endif
1724#if defined NESTING && !defined SOLVE3D
1725 dv_flux(i,j)=0.5_r8*om_v(i,j)* &
1726 & (dnew(i,j)+dnew(i,j-1))*vbar(i,j,knew)
1727#endif
1728 END DO
1729 END DO
1730 END IF
1731!
1732! Apply lateral boundary conditions.
1733!
1734 CALL u2dbc_tile (ng, tile, &
1735 & lbi, ubi, lbj, ubj, &
1736 & imins, imaxs, jmins, jmaxs, &
1737 & krhs, kstp, knew, &
1738 & ubar, vbar, zeta)
1739 CALL v2dbc_tile (ng, tile, &
1740 & lbi, ubi, lbj, ubj, &
1741 & imins, imaxs, jmins, jmaxs, &
1742 & krhs, kstp, knew, &
1743 & ubar, vbar, zeta)
1744!
1745! Compute integral mass flux across open boundaries and adjust
1746! for volume conservation.
1747!
1748 IF (any(volcons(:,ng))) THEN
1749 CALL obc_flux_tile (ng, tile, &
1750 & lbi, ubi, lbj, ubj, &
1751 & imins, imaxs, jmins, jmaxs, &
1752 & knew, &
1753#ifdef MASKING
1754 & umask, vmask, &
1755#endif
1756 & h, om_v, on_u, &
1757 & ubar, vbar, zeta)
1758 END IF
1759
1760#if defined NESTING && !defined SOLVE3D
1761!
1762! Set barotropic fluxes along physical boundaries.
1763!
1764 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
1765 IF (domain(ng)%Western_Edge(tile)) THEN
1766 DO j=jstr-1,jendr
1767 dnew(istr-1,j)=h(istr-1,j)+zeta_new(istr-1,j)
1768 END DO
1769 END IF
1770 END IF
1771 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
1772 IF (domain(ng)%Eastern_Edge(tile)) THEN
1773 DO j=jstr-1,jendr
1774 dnew(iend+1,j)=h(iend+1,j)+zeta_new(iend+1,j)
1775 END DO
1776 END IF
1777 END IF
1778 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
1779 IF (domain(ng)%Southern_Edge(tile)) THEN
1780 DO i=istr-1,iendr
1781 dnew(i,jstr-1)=h(i,jstr-1)+zeta_new(i,jstr-1)
1782 END DO
1783 END IF
1784 END IF
1785 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
1786 IF (domain(ng)%Northern_Edge(tile)) THEN
1787 DO i=istr-1,iendr
1788 dnew(i,jend+1)=h(i,jend+1)+zeta_new(i,jend+1)
1789 END DO
1790 END IF
1791 END IF
1792!
1793 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
1794 IF (domain(ng)%Western_Edge(tile)) THEN
1795 DO j=jstrr,jendr
1796 du_flux(istru-1,j)=0.5_r8*on_u(istru-1,j)* &
1797 & (dnew(istru-1,j)+dnew(istru-2,j))* &
1798 & ubar(istru-1,j,knew)
1799 END DO
1800 DO j=jstrv,jend
1801 dv_flux(istr-1,j)=0.5_r8*om_v(istr-1,j)* &
1802 & (dnew(istr-1,j)+dnew(istr-1,j-1))* &
1803 & vbar(istr-1,j,knew)
1804 END DO
1805 END IF
1806 END IF
1807 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
1808 IF (domain(ng)%Eastern_Edge(tile)) THEN
1809 DO j=jstrr,jendr
1810 du_flux(iend+1,j)=0.5_r8*on_u(iend+1,j)* &
1811 & (dnew(iend+1,j)+dnew(iend,j))* &
1812 & ubar(iend+1,j,knew)
1813 END DO
1814 DO j=jstrv,jend
1815 dv_flux(iend+1,j)=0.5_r8*om_v(iend+1,j)* &
1816 & (dnew(iend+1,j)+dnew(iend+1,j-1))* &
1817 & vbar(iend+1,j,knew)
1818 END DO
1819 END IF
1820 END IF
1821 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
1822 IF (domain(ng)%Southern_Edge(tile)) THEN
1823 DO i=istru,iend
1824 du_flux(i,jstr-1)=0.5_r8*on_u(i,jstr-1)* &
1825 & (dnew(i,jstr-1)+dnew(i-1,jstr-1))* &
1826 & ubar(i,jstr-1,knew)
1827 END DO
1828 DO i=istrr,iendr
1829 dv_flux(i,jstrv-1)=0.5_r8*om_v(i,jstrv-1)* &
1830 & (dnew(i,jstrv-1)+dnew(i,jstrv-2))* &
1831 & vbar(i,jstrv-1,knew)
1832 END DO
1833 END IF
1834 END IF
1835 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
1836 IF (domain(ng)%Northern_Edge(tile)) THEN
1837 DO i=istru,iend
1838 du_flux(i,jend+1)=0.5_r8*on_u(i,jend+1)* &
1839 & (dnew(i,jend+1)+dnew(i-1,jend+1))* &
1840 & ubar(i,jend+1,knew)
1841 END DO
1842 DO i=istrr,iendr
1843 dv_flux(i,jend+1)=0.5_r8*om_v(i,jend+1)* &
1844 & (dnew(i,jend+1)+dnew(i,jend))* &
1845 & vbar(i,jend+1,knew)
1846 END DO
1847 END IF
1848 END IF
1849#endif
1850!
1851! Apply momentum transport point sources (like river runoff), if any.
1852!
1853! Dsrc(is) = 0, flow across grid cell u-face (positive or negative)
1854! Dsrc(is) = 1, flow across grid cell v-face (positive or negative)
1855!
1856 IF (luvsrc(ng)) THEN
1857 DO is=1,nsrc(ng)
1858 i=sources(ng)%Isrc(is)
1859 j=sources(ng)%Jsrc(is)
1860 IF (((istrr.le.i).and.(i.le.iendr)).and. &
1861 & ((jstrr.le.j).and.(j.le.jendr))) THEN
1862 IF (int(sources(ng)%Dsrc(is)).eq.0) THEN
1863 cff=1.0_r8/(on_u(i,j)* &
1864 & 0.5_r8*(dnew(i-1,j)+dnew(i,j)))
1865 ubar(i,j,knew)=sources(ng)%Qbar(is)*cff
1866#ifdef SOLVE3D
1867 du_avg1(i,j)=sources(ng)%Qbar(is)
1868#endif
1869#if defined NESTING && !defined SOLVE3D
1870 du_flux(i,j)=sources(ng)%Qbar(is)
1871#endif
1872 ELSE IF (int(sources(ng)%Dsrc(is)).eq.1) THEN
1873 cff=1.0_r8/(om_v(i,j)* &
1874 & 0.5_r8*(dnew(i,j-1)+dnew(i,j)))
1875 vbar(i,j,knew)=sources(ng)%Qbar(is)*cff
1876#ifdef SOLVE3D
1877 dv_avg1(i,j)=sources(ng)%Qbar(is)
1878#endif
1879#if defined NESTING && !defined SOLVE3D
1880 dv_flux(i,j)=sources(ng)%Qbar(is)
1881#endif
1882 END IF
1883 END IF
1884 END DO
1885 END IF
1886
1887#ifdef SOLVE3D
1888!
1889!-----------------------------------------------------------------------
1890! Finalize computation of barotropic mode averages.
1891!-----------------------------------------------------------------------
1892!
1893! This procedure starts with filling in boundary rows of total depths
1894! at the new time step, which is needed to be done only during the
1895! last barotropic time step, Normally, the computation of averages
1896! occurs at the beginning of the next predictor step because "DUon"
1897! and "DVom" are being computed anyway. Strictly speaking, the filling
1898! the boundaries are necessary only in the case of open boundaries,
1899! otherwise, the associated fluxes are all zeros.
1900!
1901 IF ((iif(ng).eq.nfast(ng)).and.(knew.lt.3)) THEN
1902 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
1903 IF (domain(ng)%Western_Edge(tile)) THEN
1904 DO j=jstr-1,jendr
1905 dnew(istr-1,j)=h(istr-1,j)+zeta_new(istr-1,j)
1906 END DO
1907 END IF
1908 END IF
1909 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
1910 IF (domain(ng)%Eastern_Edge(tile)) THEN
1911 DO j=jstr-1,jendr
1912 dnew(iend+1,j)=h(iend+1,j)+zeta_new(iend+1,j)
1913 END DO
1914 END IF
1915 END IF
1916 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
1917 IF (domain(ng)%Southern_Edge(tile)) THEN
1918 DO i=istr-1,iendr
1919 dnew(i,jstr-1)=h(i,jstr-1)+zeta_new(i,jstr-1)
1920 END DO
1921 END IF
1922 END IF
1923 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
1924 IF (domain(ng)%Northern_Edge(tile)) THEN
1925 DO i=istr-1,iendr
1926 dnew(i,jend+1)=h(i,jend+1)+zeta_new(i,jend+1)
1927 END DO
1928 END IF
1929 END IF
1930!
1931! At the end of the last 2D time step replace the new free-surface
1932! zeta(:,:,knew) with it fast time-averaged value, Zt_avg1. Recall
1933! this is state variable is the one that communicates with the 3D
1934! kernel. Then, compute time-dependent depths.
1935!
1936 cff=weight(1,iif(ng),ng)
1937 cff1=0.5*cff
1938 DO j=jstrr,jendr
1939 DO i=istrr,iendr
1940 zt_avg1(i,j)=zt_avg1(i,j)+ &
1941 & cff*zeta(i,j,knew)
1942 IF (i.ge.istr) THEN
1943 du_avg1(i,j)=du_avg1(i,j)+ &
1944 & cff1*on_u(i,j)* &
1945 & (dnew(i,j)+dnew(i-1,j))*ubar(i,j,knew)
1946 END IF
1947 IF (j.ge.jstr) THEN
1948 dv_avg1(i,j)=dv_avg1(i,j)+ &
1949 & cff1*om_v(i,j)* &
1950 & (dnew(i,j)+dnew(i,j-1))*vbar(i,j,knew)
1951 END IF
1952 zeta(i,j,knew)=zt_avg1(i,j)
1953 END DO
1954 END DO
1955 CALL set_depth (ng, tile, inlm)
1956
1957# ifdef NESTING
1958!
1959! After all fast time steps are completed, apply boundary conditions
1960! to time averaged fields.
1961!
1962! In nesting applications with refinement grids, we need to exchange
1963! the DU_avg2 and DV_avg2 fluxes boundary information for the case
1964! that a contact point is at a tile partition. Notice that in such
1965! cases, we need i+1 and j+1 values for spatial/temporal interpolation.
1966!
1967 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
1968 CALL exchange_r2d_tile (ng, tile, &
1969 & lbi, ubi, lbj, ubj, &
1970 & zt_avg1)
1971 CALL exchange_u2d_tile (ng, tile, &
1972 & lbi, ubi, lbj, ubj, &
1973 & du_avg1)
1974 CALL exchange_v2d_tile (ng, tile, &
1975 & lbi, ubi, lbj, ubj, &
1976 & dv_avg1)
1977 CALL exchange_u2d_tile (ng, tile, &
1978 & lbi, ubi, lbj, ubj, &
1979 & du_avg2)
1980 CALL exchange_v2d_tile (ng, tile, &
1981 & lbi, ubi, lbj, ubj, &
1982 & dv_avg2)
1983 END IF
1984
1985# ifdef DISTRIBUTE
1986 CALL mp_exchange2d (ng, tile, inlm, 3, &
1987 & lbi, ubi, lbj, ubj, &
1988 & nghostpoints, &
1989 & ewperiodic(ng), nsperiodic(ng), &
1990 & zt_avg1, du_avg1, dv_avg1)
1991 CALL mp_exchange2d (ng, tile, inlm, 2, &
1992 & lbi, ubi, lbj, ubj, &
1993 & nghostpoints, &
1994 & ewperiodic(ng), nsperiodic(ng), &
1995 & du_avg2, dv_avg2)
1996# endif
1997# endif
1998 END IF
1999#endif
2000#if defined NESTING && !defined SOLVE3D
2001!
2002! In nesting applications with refinement grids, we need to exchange
2003! the DU_flux and DV_flux fluxes boundary information for the case
2004! that a contact point is at a tile partition. Notice that in such
2005! cases, we need i+1 and j+1 values for spatial/temporal interpolation.
2006!
2007 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
2008 CALL exchange_u2d_tile (ng, tile, &
2009 & lbi, ubi, lbj, ubj, &
2010 & du_flux)
2011 CALL exchange_v2d_tile (ng, tile, &
2012 & lbi, ubi, lbj, ubj, &
2013 & dv_flux)
2014 END IF
2015
2016# ifdef DISTRIBUTE
2017!
2018 CALL mp_exchange2d (ng, tile, inlm, 2, &
2019 & lbi, ubi, lbj, ubj, &
2020 & nghostpoints, &
2021 & ewperiodic(ng), nsperiodic(ng), &
2022 & du_flux, dv_flux)
2023# endif
2024#endif
2025!
2026! Deallocate local new free-surface.
2027!
2028 deallocate ( zeta_new )
2029
2030#ifdef WET_DRY
2031!
2032!-----------------------------------------------------------------------
2033! Compute new wet/dry masks.
2034!-----------------------------------------------------------------------
2035!
2036 CALL wetdry_tile (ng, tile, &
2037 & lbi, ubi, lbj, ubj, &
2038 & imins, imaxs, jmins, jmaxs, &
2039# ifdef MASKING
2040 & pmask, rmask, umask, vmask, &
2041# endif
2042 & h, zeta(:,:,knew), &
2043# ifdef SOLVE3D
2044 & du_avg1, dv_avg1, &
2045 & rmask_wet_avg, &
2046# endif
2047 & pmask_wet, pmask_full, &
2048 & rmask_wet, rmask_full, &
2049 & umask_wet, umask_full, &
2050 & vmask_wet, vmask_full)
2051#endif
2052!
2053!-----------------------------------------------------------------------
2054! Exchange boundary information.
2055!-----------------------------------------------------------------------
2056!
2057 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
2058 CALL exchange_r2d_tile (ng, tile, &
2059 & lbi, ubi, lbj, ubj, &
2060 & zeta(:,:,knew))
2061 CALL exchange_u2d_tile (ng, tile, &
2062 & lbi, ubi, lbj, ubj, &
2063 & ubar(:,:,knew))
2064 CALL exchange_v2d_tile (ng, tile, &
2065 & lbi, ubi, lbj, ubj, &
2066 & vbar(:,:,knew))
2067 END IF
2068
2069#ifdef DISTRIBUTE
2070!
2071 CALL mp_exchange2d (ng, tile, inlm, 3, &
2072 & lbi, ubi, lbj, ubj, &
2073 & nghostpoints, &
2074 & ewperiodic(ng), nsperiodic(ng), &
2075 & zeta(:,:,knew), &
2076 & ubar(:,:,knew), &
2077 & vbar(:,:,knew))
2078#endif
2079!
2080 RETURN

References mod_scalars::compositegrid, mod_scalars::dcrit, mod_param::domain, mod_scalars::dtfast, mod_scalars::ewperiodic, exchange_2d_mod::exchange_r2d_tile(), exchange_2d_mod::exchange_u2d_tile(), exchange_2d_mod::exchange_v2d_tile(), mod_scalars::g, mod_scalars::ieast, mod_scalars::iif, mod_param::inlm, mod_scalars::inorth, mod_scalars::isouth, mod_scalars::iwest, mod_scalars::luvsrc, mod_scalars::lwsrc, mod_scalars::m2bstr, mod_scalars::m2fcor, mod_scalars::m2hadv, mod_scalars::m2hvis, mod_scalars::m2pgrd, mod_scalars::m2xadv, mod_scalars::m2xvis, mod_scalars::m2yadv, mod_scalars::m2yvis, mp_exchange_mod::mp_exchange2d(), mod_scalars::nfast, mod_param::nghostpoints, mod_scalars::nsperiodic, mod_sources::nsrc, obc_volcons_mod::obc_flux_tile(), mod_scalars::predictor_2d_step, mod_scalars::rho0, set_depth_mod::set_depth(), obc_volcons_mod::set_duv_bc_tile(), mod_sources::sources, u2dbc_mod::u2dbc_tile(), v2dbc_mod::v2dbc_tile(), mod_scalars::volcons, mod_scalars::weight, and wetdry_mod::wetdry_tile().

Here is the call graph for this function:

◆ step2d_tile() [2/3]

subroutine step2d_mod::step2d_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) ubk,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) krhs,
integer, intent(in) kstp,
integer, intent(in) knew,
integer, intent(in) nstp,
integer, intent(in) nnew,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pmask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rmask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) umask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) vmask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) pmask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) pmask_full,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) rmask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) rmask_full,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) umask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) umask_full,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) vmask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) vmask_full,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) rmask_wet_avg,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) fomn,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) h,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_u,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_v,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_u,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_v,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) omn,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pm,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pn,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) dndx,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) dmde,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pmon_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pnom_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pmon_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pnom_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) visc2_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) visc2_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) visc4_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) visc4_r,
real(r8), dimension(lbi:ubi,lbj:ubj,1:3), intent(in) bed_thick,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rurol2d,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rvrol2d,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rubst2d,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rvbst2d,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) russt2d,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rvsst2d,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rubrk2d,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rvbrk2d,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rukvf2d,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rvkvf2d,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) bh,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) qsp,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) zetaw,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) ubar_stokes,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) vbar_stokes,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) eq_tide,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) sustr,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) svstr,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) bustr,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) bvstr,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pair,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rhoa,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rhos,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) du_avg1,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) du_avg2,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) dv_avg1,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) dv_avg2,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) zt_avg1,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) rufrc,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) rvfrc,
real(r8), dimension(lbi:ubi,lbj:ubj,0:ubk,2), intent(inout) ru,
real(r8), dimension(lbi:ubi,lbj:ubj,0:ubk,2), intent(inout) rv,
real(r8), dimension(lbi:ubi,lbj:ubj,ndm2d), intent(inout) diau2wrk,
real(r8), dimension(lbi:ubi,lbj:ubj,ndm2d), intent(inout) diav2wrk,
real(r8), dimension(lbi:ubi,lbj:ubj,2,ndm2d-1), intent(inout) diarubar,
real(r8), dimension(lbi:ubi,lbj:ubj,2,ndm2d-1), intent(inout) diarvbar,
real(r8), dimension(lbi:ubi,lbj:ubj,ndm2d), intent(inout) diau2int,
real(r8), dimension(lbi:ubi,lbj:ubj,ndm2d), intent(inout) diav2int,
real(r8), dimension(lbi:ubi,lbj:ubj,3,ndm2d-1), intent(inout) diarufrc,
real(r8), dimension(lbi:ubi,lbj:ubj,3,ndm2d-1), intent(inout) diarvfrc,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(out) du_flux,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(out) dv_flux,
real(r8), dimension(lbi:ubi,lbj:ubj,2), intent(inout) rubar,
real(r8), dimension(lbi:ubi,lbj:ubj,2), intent(inout) rvbar,
real(r8), dimension(lbi:ubi,lbj:ubj,2), intent(inout) rzeta,
real(r8), dimension(lbi:ubi,lbj:ubj,:), intent(inout) ubar,
real(r8), dimension(lbi:ubi,lbj:ubj,:), intent(inout) vbar,
real(r8), dimension(lbi:ubi,lbj:ubj,:), intent(inout) zeta )
private

Definition at line 163 of file step2d_LF_AM3.h.

247!***********************************************************************
248!
249 USE mod_param
250 USE mod_clima
251 USE mod_ncparam
252 USE mod_scalars
253#if defined SEDIMENT && defined SED_MORPH
254 USE mod_sedbed
255#endif
256 USE mod_sources
257!
259#ifdef DISTRIBUTE
261#endif
263 USE u2dbc_mod, ONLY : u2dbc_tile
264 USE v2dbc_mod, ONLY : v2dbc_tile
265 USE zetabc_mod, ONLY : zetabc_tile
266#ifdef WET_DRY
267 USE wetdry_mod, ONLY : wetdry_tile
268#endif
269!
270! Imported variable declarations.
271!
272 integer, intent(in ) :: ng, tile
273 integer, intent(in ) :: LBi, UBi, LBj, UBj, UBk
274 integer, intent(in ) :: IminS, ImaxS, JminS, JmaxS
275 integer, intent(in ) :: krhs, kstp, knew
276#ifdef SOLVE3D
277 integer, intent(in ) :: nstp, nnew
278#endif
279!
280#ifdef ASSUMED_SHAPE
281# ifdef MASKING
282 real(r8), intent(in ) :: pmask(LBi:,LBj:)
283 real(r8), intent(in ) :: rmask(LBi:,LBj:)
284 real(r8), intent(in ) :: umask(LBi:,LBj:)
285 real(r8), intent(in ) :: vmask(LBi:,LBj:)
286# endif
287 real(r8), intent(in ) :: fomn(LBi:,LBj:)
288# if defined SEDIMENT && defined SED_MORPH
289 real(r8), intent(inout) :: h(LBi:,LBj:)
290# else
291 real(r8), intent(in ) :: h(LBi:,LBj:)
292# endif
293 real(r8), intent(in ) :: om_u(LBi:,LBj:)
294 real(r8), intent(in ) :: om_v(LBi:,LBj:)
295 real(r8), intent(in ) :: on_u(LBi:,LBj:)
296 real(r8), intent(in ) :: on_v(LBi:,LBj:)
297 real(r8), intent(in ) :: omn(LBi:,LBj:)
298 real(r8), intent(in ) :: pm(LBi:,LBj:)
299 real(r8), intent(in ) :: pn(LBi:,LBj:)
300# if defined CURVGRID && defined UV_ADV
301 real(r8), intent(in ) :: dndx(LBi:,LBj:)
302 real(r8), intent(in ) :: dmde(LBi:,LBj:)
303# endif
304# if defined UV_VIS2 || defined UV_VIS4
305 real(r8), intent(in ) :: pmon_r(LBi:,LBj:)
306 real(r8), intent(in ) :: pnom_r(LBi:,LBj:)
307 real(r8), intent(in ) :: pmon_p(LBi:,LBj:)
308 real(r8), intent(in ) :: pnom_p(LBi:,LBj:)
309 real(r8), intent(in ) :: om_r(LBi:,LBj:)
310 real(r8), intent(in ) :: on_r(LBi:,LBj:)
311 real(r8), intent(in ) :: om_p(LBi:,LBj:)
312 real(r8), intent(in ) :: on_p(LBi:,LBj:)
313# ifdef UV_VIS2
314 real(r8), intent(in ) :: visc2_p(LBi:,LBj:)
315 real(r8), intent(in ) :: visc2_r(LBi:,LBj:)
316# endif
317# ifdef UV_VIS4
318 real(r8), intent(in ) :: visc4_p(LBi:,LBj:)
319 real(r8), intent(in ) :: visc4_r(LBi:,LBj:)
320# endif
321# endif
322# if defined SEDIMENT && defined SED_MORPH
323 real(r8), intent(in ) :: bed_thick(LBi:,LBj:,:)
324# endif
325# ifdef WEC
326# ifdef WEC_VF
327# ifdef WEC_ROLLER
328 real(r8), intent(in) :: rurol2d(LBi:,LBj:)
329 real(r8), intent(in) :: rvrol2d(LBi:,LBj:)
330# endif
331# ifdef BOTTOM_STREAMING
332 real(r8), intent(in) :: rubst2d(LBi:,LBj:)
333 real(r8), intent(in) :: rvbst2d(LBi:,LBj:)
334# endif
335# ifdef SURFACE_STREAMING
336 real(r8), intent(in) :: russt2d(LBi:,LBj:)
337 real(r8), intent(in) :: rvsst2d(LBi:,LBj:)
338# endif
339 real(r8), intent(in) :: rubrk2d(LBi:,LBj:)
340 real(r8), intent(in) :: rvbrk2d(LBi:,LBj:)
341 real(r8), intent(in) :: rukvf2d(LBi:,LBj:)
342 real(r8), intent(in) :: rvkvf2d(LBi:,LBj:)
343 real(r8), intent(in) :: bh(LBi:,LBj:)
344 real(r8), intent(in) :: qsp(LBi:,LBj:)
345 real(r8), intent(in) :: zetaw(LBi:,LBj:)
346# endif
347 real(r8), intent(in) :: ubar_stokes(LBi:,LBj:)
348 real(r8), intent(in) :: vbar_stokes(LBi:,LBj:)
349# endif
350# if defined TIDE_GENERATING_FORCES && !defined SOLVE3D
351 real(r8), intent(in ) :: eq_tide(LBi:,LBj:)
352# endif
353# ifndef SOLVE3D
354 real(r8), intent(in ) :: sustr(LBi:,LBj:)
355 real(r8), intent(in ) :: svstr(LBi:,LBj:)
356 real(r8), intent(in ) :: bustr(LBi:,LBj:)
357 real(r8), intent(in ) :: bvstr(LBi:,LBj:)
358# ifdef ATM_PRESS
359 real(r8), intent(in ) :: Pair(LBi:,LBj:)
360# endif
361# else
362# ifdef VAR_RHO_2D
363 real(r8), intent(in ) :: rhoA(LBi:,LBj:)
364 real(r8), intent(in ) :: rhoS(LBi:,LBj:)
365# endif
366 real(r8), intent(inout) :: DU_avg1(LBi:,LBj:)
367 real(r8), intent(inout) :: DU_avg2(LBi:,LBj:)
368 real(r8), intent(inout) :: DV_avg1(LBi:,LBj:)
369 real(r8), intent(inout) :: DV_avg2(LBi:,LBj:)
370 real(r8), intent(inout) :: Zt_avg1(LBi:,LBj:)
371 real(r8), intent(inout) :: rufrc(LBi:,LBj:)
372 real(r8), intent(inout) :: rvfrc(LBi:,LBj:)
373 real(r8), intent(inout) :: ru(LBi:,LBj:,0:,:)
374 real(r8), intent(inout) :: rv(LBi:,LBj:,0:,:)
375# endif
376# ifdef WET_DRY
377 real(r8), intent(inout) :: pmask_full(LBi:,LBj:)
378 real(r8), intent(inout) :: rmask_full(LBi:,LBj:)
379 real(r8), intent(inout) :: umask_full(LBi:,LBj:)
380 real(r8), intent(inout) :: vmask_full(LBi:,LBj:)
381
382 real(r8), intent(inout) :: pmask_wet(LBi:,LBj:)
383 real(r8), intent(inout) :: rmask_wet(LBi:,LBj:)
384 real(r8), intent(inout) :: umask_wet(LBi:,LBj:)
385 real(r8), intent(inout) :: vmask_wet(LBi:,LBj:)
386# ifdef SOLVE3D
387 real(r8), intent(inout) :: rmask_wet_avg(LBi:,LBj:)
388# endif
389# endif
390# ifdef DIAGNOSTICS_UV
391 real(r8), intent(inout) :: DiaU2wrk(LBi:,LBj:,:)
392 real(r8), intent(inout) :: DiaV2wrk(LBi:,LBj:,:)
393 real(r8), intent(inout) :: DiaRUbar(LBi:,LBj:,:,:)
394 real(r8), intent(inout) :: DiaRVbar(LBi:,LBj:,:,:)
395# ifdef SOLVE3D
396 real(r8), intent(inout) :: DiaU2int(LBi:,LBj:,:)
397 real(r8), intent(inout) :: DiaV2int(LBi:,LBj:,:)
398 real(r8), intent(inout) :: DiaRUfrc(LBi:,LBj:,:,:)
399 real(r8), intent(inout) :: DiaRVfrc(LBi:,LBj:,:,:)
400# endif
401# endif
402 real(r8), intent(inout) :: rubar(LBi:,LBj:,:)
403 real(r8), intent(inout) :: rvbar(LBi:,LBj:,:)
404 real(r8), intent(inout) :: rzeta(LBi:,LBj:,:)
405 real(r8), intent(inout) :: ubar(LBi:,LBj:,:)
406 real(r8), intent(inout) :: vbar(LBi:,LBj:,:)
407 real(r8), intent(inout) :: zeta(LBi:,LBj:,:)
408# if defined NESTING && !defined SOLVE3D
409 real(r8), intent(out ) :: DU_flux(LBi:,LBj:)
410 real(r8), intent(out ) :: DV_flux(LBi:,LBj:)
411# endif
412
413#else
414
415# ifdef MASKING
416 real(r8), intent(in ) :: pmask(LBi:UBi,LBj:UBj)
417 real(r8), intent(in ) :: rmask(LBi:UBi,LBj:UBj)
418 real(r8), intent(in ) :: umask(LBi:UBi,LBj:UBj)
419 real(r8), intent(in ) :: vmask(LBi:UBi,LBj:UBj)
420# endif
421 real(r8), intent(in ) :: fomn(LBi:UBi,LBj:UBj)
422# if defined SEDIMENT && defined SED_MORPH
423 real(r8), intent(inout) :: h(LBi:UBi,LBj:UBj)
424# else
425 real(r8), intent(in ) :: h(LBi:UBi,LBj:UBj)
426# endif
427 real(r8), intent(in ) :: om_u(LBi:UBi,LBj:UBj)
428 real(r8), intent(in ) :: om_v(LBi:UBi,LBj:UBj)
429 real(r8), intent(in ) :: on_u(LBi:UBi,LBj:UBj)
430 real(r8), intent(in ) :: on_v(LBi:UBi,LBj:UBj)
431 real(r8), intent(in ) :: omn(LBi:UBi,LBj:UBj)
432 real(r8), intent(in ) :: pm(LBi:UBi,LBj:UBj)
433 real(r8), intent(in ) :: pn(LBi:UBi,LBj:UBj)
434# if defined CURVGRID && defined UV_ADV
435 real(r8), intent(in ) :: dndx(LBi:UBi,LBj:UBj)
436 real(r8), intent(in ) :: dmde(LBi:UBi,LBj:UBj)
437# endif
438# if defined UV_VIS2 || defined UV_VIS4
439 real(r8), intent(in ) :: pmon_r(LBi:UBi,LBj:UBj)
440 real(r8), intent(in ) :: pnom_r(LBi:UBi,LBj:UBj)
441 real(r8), intent(in ) :: pmon_p(LBi:UBi,LBj:UBj)
442 real(r8), intent(in ) :: pnom_p(LBi:UBi,LBj:UBj)
443 real(r8), intent(in ) :: om_r(LBi:UBi,LBj:UBj)
444 real(r8), intent(in ) :: on_r(LBi:UBi,LBj:UBj)
445 real(r8), intent(in ) :: om_p(LBi:UBi,LBj:UBj)
446 real(r8), intent(in ) :: on_p(LBi:UBi,LBj:UBj)
447# ifdef UV_VIS2
448 real(r8), intent(in ) :: visc2_p(LBi:UBi,LBj:UBj)
449 real(r8), intent(in ) :: visc2_r(LBi:UBi,LBj:UBj)
450# endif
451# ifdef UV_VIS4
452 real(r8), intent(in ) :: visc4_p(LBi:UBi,LBj:UBj)
453 real(r8), intent(in ) :: visc4_r(LBi:UBi,LBj:UBj)
454# endif
455# endif
456# if defined SEDIMENT && defined SED_MORPH
457 real(r8), intent(in ) :: bed_thick(LBi:UBi,LBj:UBj,1:3)
458# endif
459# ifdef WEC
460# ifdef WEC_VF
461# ifdef WEC_ROLLER
462 real(r8), intent(in) :: rurol2d(LBi:UBi,LBj:UBj)
463 real(r8), intent(in) :: rvrol2d(LBi:UBi,LBj:UBj)
464# endif
465# ifdef BOTTOM_STREAMING
466 real(r8), intent(in) :: rubst2d(LBi:UBi,LBj:UBj)
467 real(r8), intent(in) :: rvbst2d(LBi:UBi,LBj:UBj)
468# endif
469# ifdef SURFACE_STREAMING
470 real(r8), intent(in) :: russt2d(LBi:UBi,LBj:UBj)
471 real(r8), intent(in) :: rvsst2d(LBi:UBi,LBj:UBj)
472# endif
473 real(r8), intent(in) :: rubrk2d(LBi:UBi,LBj:UBj)
474 real(r8), intent(in) :: rvbrk2d(LBi:UBi,LBj:UBj)
475 real(r8), intent(in) :: rukvf2d(LBi:UBi,LBj:UBj)
476 real(r8), intent(in) :: rvkvf2d(LBi:UBi,LBj:UBj)
477 real(r8), intent(in) :: bh(LBi:UBi,LBj:UBj)
478 real(r8), intent(in) :: qsp(LBi:UBi,LBj:UBj)
479 real(r8), intent(in) :: zetaw(LBi:UBi,LBj:UBj)
480# endif
481 real(r8), intent(in) :: ubar_stokes(LBi:UBi,LBj:UBj)
482 real(r8), intent(in) :: vbar_stokes(LBi:UBi,LBj:UBj)
483# endif
484# if defined TIDE_GENERATING_FORCES && !defined SOLVE3D
485 real(r8), intent(in ) :: eq_tide(LBi:UBi,LBj:UBj)
486# endif
487# ifndef SOLVE3D
488 real(r8), intent(in ) :: sustr(LBi:UBi,LBj:UBj)
489 real(r8), intent(in ) :: svstr(LBi:UBi,LBj:UBj)
490 real(r8), intent(in ) :: bustr(LBi:UBi,LBj:UBj)
491 real(r8), intent(in ) :: bvstr(LBi:UBi,LBj:UBj)
492# ifdef ATM_PRESS
493 real(r8), intent(in ) :: Pair(LBi:UBi,LBj:UBj)
494# endif
495# else
496# ifdef VAR_RHO_2D
497 real(r8), intent(in ) :: rhoA(LBi:UBi,LBj:UBj)
498 real(r8), intent(in ) :: rhoS(LBi:UBi,LBj:UBj)
499# endif
500 real(r8), intent(inout) :: DU_avg1(LBi:UBi,LBj:UBj)
501 real(r8), intent(inout) :: DU_avg2(LBi:UBi,LBj:UBj)
502 real(r8), intent(inout) :: DV_avg1(LBi:UBi,LBj:UBj)
503 real(r8), intent(inout) :: DV_avg2(LBi:UBi,LBj:UBj)
504 real(r8), intent(inout) :: Zt_avg1(LBi:UBi,LBj:UBj)
505 real(r8), intent(inout) :: rufrc(LBi:UBi,LBj:UBj)
506 real(r8), intent(inout) :: rvfrc(LBi:UBi,LBj:UBj)
507 real(r8), intent(inout) :: ru(LBi:UBi,LBj:UBj,0:UBk,2)
508 real(r8), intent(inout) :: rv(LBi:UBi,LBj:UBj,0:UBk,2)
509# endif
510# ifdef WET_DRY
511 real(r8), intent(inout) :: pmask_full(LBi:UBi,LBj:UBj)
512 real(r8), intent(inout) :: rmask_full(LBi:UBi,LBj:UBj)
513 real(r8), intent(inout) :: umask_full(LBi:UBi,LBj:UBj)
514 real(r8), intent(inout) :: vmask_full(LBi:UBi,LBj:UBj)
515
516 real(r8), intent(inout) :: pmask_wet(LBi:UBi,LBj:UBj)
517 real(r8), intent(inout) :: rmask_wet(LBi:UBi,LBj:UBj)
518 real(r8), intent(inout) :: umask_wet(LBi:UBi,LBj:UBj)
519 real(r8), intent(inout) :: vmask_wet(LBi:UBi,LBj:UBj)
520# ifdef SOLVE3D
521 real(r8), intent(inout) :: rmask_wet_avg(LBi:UBi,LBj:UBj)
522# endif
523# endif
524# ifdef DIAGNOSTICS_UV
525 real(r8), intent(inout) :: DiaU2wrk(LBi:UBi,LBj:UBj,NDM2d)
526 real(r8), intent(inout) :: DiaV2wrk(LBi:UBi,LBj:UBj,NDM2d)
527 real(r8), intent(inout) :: DiaRUbar(LBi:UBi,LBj:UBj,2,NDM2d-1)
528 real(r8), intent(inout) :: DiaRVbar(LBi:UBi,LBj:UBj,2,NDM2d-1)
529# ifdef SOLVE3D
530 real(r8), intent(inout) :: DiaU2int(LBi:UBi,LBj:UBj,NDM2d)
531 real(r8), intent(inout) :: DiaV2int(LBi:UBi,LBj:UBj,NDM2d)
532 real(r8), intent(inout) :: DiaRUfrc(LBi:UBi,LBj:UBj,3,NDM2d-1)
533 real(r8), intent(inout) :: DiaRVfrc(LBi:UBi,LBj:UBj,3,NDM2d-1)
534# endif
535# endif
536 real(r8), intent(inout) :: rubar(LBi:UBi,LBj:UBj,2)
537 real(r8), intent(inout) :: rvbar(LBi:UBi,LBj:UBj,2)
538 real(r8), intent(inout) :: rzeta(LBi:UBi,LBj:UBj,2)
539 real(r8), intent(inout) :: ubar(LBi:UBi,LBj:UBj,:)
540 real(r8), intent(inout) :: vbar(LBi:UBi,LBj:UBj,:)
541 real(r8), intent(inout) :: zeta(LBi:UBi,LBj:UBj,:)
542# if defined NESTING && !defined SOLVE3D
543 real(r8), intent(out ) :: DU_flux(LBi:UBi,LBj:UBj)
544 real(r8), intent(out ) :: DV_flux(LBi:UBi,LBj:UBj)
545# endif
546#endif
547!
548! Local variable declarations.
549!
550 logical :: CORRECTOR_2D_STEP
551!
552 integer :: i, is, j, ptsk
553#ifdef DIAGNOSTICS_UV
554 integer :: idiag
555#endif
556!
557 real(r8) :: cff, cff1, cff2, cff3, cff4, cff5, cff6, cff7
558 real(r8) :: fac, fac1, fac2, fac3
559!
560 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Dgrad
561 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Dnew
562 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Drhs
563 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Drhs_p
564 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Dstp
565 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: DUon
566 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: DVom
567#ifdef WEC
568 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: DUSon
569 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: DVSom
570#endif
571#ifdef WEC_VF
572 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: DUSom
573 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: DVSon
574#endif
575#ifdef UV_VIS4
576 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: LapU
577 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: LapV
578#endif
579 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: UFe
580 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: UFx
581 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: VFe
582 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: VFx
583 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: grad
584 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: gzeta
585 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: gzeta2
586#if defined VAR_RHO_2D && defined SOLVE3D
587 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: gzetaSA
588#endif
589 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rhs_ubar
590 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rhs_vbar
591 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rhs_zeta
592 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: zeta_new
593 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: zwrk
594#ifdef WET_DRY
595 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: wetdry
596#endif
597#ifdef DIAGNOSTICS_UV
598 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Uwrk
599 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Vwrk
600 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,NDM2d-1) :: DiaU2rhs
601 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,NDM2d-1) :: DiaV2rhs
602#endif
603
604#include "set_bounds.h"
605!
606 ptsk=3-kstp
607 corrector_2d_step=.not.predictor_2d_step(ng)
608!
609!-----------------------------------------------------------------------
610! Compute total depth (m) and vertically integrated mass fluxes.
611!-----------------------------------------------------------------------
612!
613#if defined DISTRIBUTE && !defined NESTING
614
615! In distributed-memory, the I- and J-ranges are different and a
616! special exchange is done to avoid having three ghost points for
617! high order numerical stencils. Notice that a private array is
618! passed below to the exchange routine. It also applies periodic
619! boundary conditions, if appropriate and no partitions in I- or
620! J-directions.
621!
622 DO j=jstrv-2,jendp2
623 DO i=istru-2,iendp2
624 drhs(i,j)=zeta(i,j,krhs)+h(i,j)
625 END DO
626 END DO
627 DO j=jstrv-2,jendp2
628 DO i=istru-1,iendp2
629 cff=0.5_r8*on_u(i,j)
630 cff1=cff*(drhs(i,j)+drhs(i-1,j))
631 duon(i,j)=ubar(i,j,krhs)*cff1
632# ifdef WEC
633# ifdef WET_DRY
634 cff5=abs(abs(umask_wet(i,j))-1.0_r8)
635 cff6=0.5_r8+dsign(0.5_r8,ubar_stokes(i,j))*umask_wet(i,j)
636 cff7=0.5_r8*umask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
637 cff1=cff1*cff7
638# endif
639 duson(i,j)=ubar_stokes(i,j)*cff1
640 duon(i,j)=duon(i,j)+duson(i,j)
641# endif
642 END DO
643 END DO
644 DO j=jstrv-1,jendp2
645 DO i=istru-2,iendp2
646 cff=0.5_r8*om_v(i,j)
647 cff1=cff*(drhs(i,j)+drhs(i,j-1))
648 dvom(i,j)=vbar(i,j,krhs)*cff1
649# ifdef WEC
650# ifdef WET_DRY
651 cff5=abs(abs(vmask_wet(i,j))-1.0_r8)
652 cff6=0.5_r8+dsign(0.5_r8,vbar_stokes(i,j))*vmask_wet(i,j)
653 cff7=0.5_r8*vmask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
654 cff1=cff1*cff7
655# endif
656 dvsom(i,j)=vbar_stokes(i,j)*cff1
657 dvom(i,j)=dvom(i,j)+dvsom(i,j)
658# endif
659 END DO
660 END DO
661
662#else
663
664 DO j=jstrvm2-1,jendp2
665 DO i=istrum2-1,iendp2
666 drhs(i,j)=zeta(i,j,krhs)+h(i,j)
667 END DO
668 END DO
669 DO j=jstrvm2-1,jendp2
670 DO i=istrum2,iendp2
671 cff=0.5_r8*on_u(i,j)
672 cff1=cff*(drhs(i,j)+drhs(i-1,j))
673 duon(i,j)=ubar(i,j,krhs)*cff1
674# ifdef WEC
675# ifdef WET_DRY
676 cff5=abs(abs(umask_wet(i,j))-1.0_r8)
677 cff6=0.5_r8+dsign(0.5_r8,ubar_stokes(i,j))*umask_wet(i,j)
678 cff7=0.5_r8*umask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
679 cff1=cff1*cff7
680# endif
681 duson(i,j)=ubar_stokes(i,j)*cff1
682 duon(i,j)=duon(i,j)+duson(i,j)
683# endif
684 END DO
685 END DO
686 DO j=jstrvm2,jendp2
687 DO i=istrum2-1,iendp2
688 cff=0.5_r8*om_v(i,j)
689 cff1=cff*(drhs(i,j)+drhs(i,j-1))
690 dvom(i,j)=vbar(i,j,krhs)*cff1
691# ifdef WEC
692# ifdef WET_DRY
693 cff5=abs(abs(vmask_wet(i,j))-1.0_r8)
694 cff6=0.5_r8+dsign(0.5_r8,vbar_stokes(i,j))*vmask_wet(i,j)
695 cff7=0.5_r8*vmask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
696 cff1=cff1*cff7
697# endif
698 dvsom(i,j)=vbar_stokes(i,j)*cff1
699 dvom(i,j)=dvom(i,j)+dvsom(i,j)
700# endif
701 END DO
702 END DO
703#endif
704#ifdef DISTRIBUTE
705!
706 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
707 CALL exchange_u2d_tile (ng, tile, &
708 & imins, imaxs, jmins, jmaxs, &
709 & duon)
710 CALL exchange_v2d_tile (ng, tile, &
711 & imins, imaxs, jmins, jmaxs, &
712 & dvom)
713 END IF
714 CALL mp_exchange2d (ng, tile, inlm, 2, &
715 & imins, imaxs, jmins, jmaxs, &
716 & nghostpoints, &
717 & ewperiodic(ng), nsperiodic(ng), &
718 & duon, dvom)
719#endif
720!
721! Set vertically integrated mass fluxes DUon and DVom along the open
722! boundaries in such a way that the integral volume is conserved.
723!
724 IF (any(volcons(:,ng))) THEN
725 CALL set_duv_bc_tile (ng, tile, &
726 & lbi, ubi, lbj, ubj, &
727 & imins, imaxs, jmins, jmaxs, &
728 & krhs, &
729#ifdef MASKING
730 & umask, vmask, &
731#endif
732 & om_v, on_u, &
733 & ubar, vbar, &
734 & drhs, duon, dvom)
735 END IF
736#ifdef SOLVE3D
737!
738!-----------------------------------------------------------------------
739! Compute time averaged fields over all short time-steps.
740!-----------------------------------------------------------------------
741!
742 IF (predictor_2d_step(ng)) THEN
743 IF (first_2d_step) THEN
744!
745! Reset arrays for 2D fields averaged within the short time-steps.
746!
747 cff2=(-1.0_r8/12.0_r8)*weight(2,iif(ng)+1,ng)
748 DO j=jstrr,jendr
749 DO i=istrr,iendr
750 zt_avg1(i,j)=0.0_r8
751 END DO
752 DO i=istr,iendr
753 du_avg1(i,j)=0.0_r8
754 du_avg2(i,j)=cff2*duon(i,j)
755 END DO
756 END DO
757 DO j=jstr,jendr
758 DO i=istrr,iendr
759 dv_avg1(i,j)=0.0_r8
760 dv_avg2(i,j)=cff2*dvom(i,j)
761 END DO
762 END DO
763 ELSE
764!
765! Accumulate field averages of previous time-step after they are
766! computed in the previous corrector step, updated their boundaries,
767! and synchronized.
768!
769 cff1=weight(1,iif(ng)-1,ng)
770 cff2=(8.0_r8/12.0_r8)*weight(2,iif(ng) ,ng)- &
771 & (1.0_r8/12.0_r8)*weight(2,iif(ng)+1,ng)
772 DO j=jstrr,jendr
773 DO i=istrr,iendr
774 zt_avg1(i,j)=zt_avg1(i,j)+cff1*zeta(i,j,krhs)
775 END DO
776 DO i=istr,iendr
777 du_avg1(i,j)=du_avg1(i,j)+cff1*duon(i,j)
778# ifdef WEC
779 du_avg1(i,j)=du_avg1(i,j)-cff1*duson(i,j)
780# endif
781 du_avg2(i,j)=du_avg2(i,j)+cff2*duon(i,j)
782 END DO
783 END DO
784 DO j=jstr,jendr
785 DO i=istrr,iendr
786 dv_avg1(i,j)=dv_avg1(i,j)+cff1*dvom(i,j)
787# ifdef WEC
788 dv_avg1(i,j)=dv_avg1(i,j)-cff1*dvsom(i,j)
789# endif
790 dv_avg2(i,j)=dv_avg2(i,j)+cff2*dvom(i,j)
791 END DO
792 END DO
793 END IF
794 ELSE
795 IF (first_2d_step) THEN
796 cff2=weight(2,iif(ng),ng)
797 ELSE
798 cff2=(5.0_r8/12.0_r8)*weight(2,iif(ng),ng)
799 END IF
800 DO j=jstrr,jendr
801 DO i=istr,iendr
802 du_avg2(i,j)=du_avg2(i,j)+cff2*duon(i,j)
803 END DO
804 END DO
805 DO j=jstr,jendr
806 DO i=istrr,iendr
807 dv_avg2(i,j)=dv_avg2(i,j)+cff2*dvom(i,j)
808 END DO
809 END DO
810 END IF
811!
812! After all fast time steps are completed, apply boundary conditions
813! to time averaged fields.
814# ifdef NESTING
815! In nesting applications with refinement grids, we need to exchange
816! the DU_avg2 and DV_avg2 fluxes boundary information for the case
817! that a contact point is at a tile partition. Notice that in such
818! cases, we need i+1 and j+1 values for spatial/temporal interpolation.
819# endif
820!
821 IF ((iif(ng).eq.(nfast(ng)+1)).and.predictor_2d_step(ng)) THEN
822 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
823 CALL exchange_r2d_tile (ng, tile, &
824 & lbi, ubi, lbj, ubj, &
825 & zt_avg1)
826 CALL exchange_u2d_tile (ng, tile, &
827 & lbi, ubi, lbj, ubj, &
828 & du_avg1)
829 CALL exchange_v2d_tile (ng, tile, &
830 & lbi, ubi, lbj, ubj, &
831 & dv_avg1)
832# ifdef NESTING
833 CALL exchange_u2d_tile (ng, tile, &
834 & lbi, ubi, lbj, ubj, &
835 & du_avg2)
836 CALL exchange_v2d_tile (ng, tile, &
837 & lbi, ubi, lbj, ubj, &
838 & dv_avg2)
839# endif
840 END IF
841# ifdef DISTRIBUTE
842 CALL mp_exchange2d (ng, tile, inlm, 3, &
843 & lbi, ubi, lbj, ubj, &
844 & nghostpoints, &
845 & ewperiodic(ng), nsperiodic(ng), &
846 & zt_avg1, du_avg1, dv_avg1)
847# ifdef NESTING
848 CALL mp_exchange2d (ng, tile, inlm, 2, &
849 & lbi, ubi, lbj, ubj, &
850 & nghostpoints, &
851 & ewperiodic(ng), nsperiodic(ng), &
852 & du_avg2, dv_avg2)
853# endif
854# endif
855 END IF
856#endif
857#ifdef WET_DRY
858!
859!-----------------------------------------------------------------------
860! Compute new wet/dry masks.
861!-----------------------------------------------------------------------
862!
863 CALL wetdry_tile (ng, tile, &
864 & lbi, ubi, lbj, ubj, &
865 & imins, imaxs, jmins, jmaxs, &
866# ifdef MASKING
867 & pmask, rmask, umask, vmask, &
868# endif
869 & h, zeta(:,:,kstp), &
870# ifdef SOLVE3D
871 & du_avg1, dv_avg1, &
872 & rmask_wet_avg, &
873# endif
874 & pmask_wet, pmask_full, &
875 & rmask_wet, rmask_full, &
876 & umask_wet, umask_full, &
877 & vmask_wet, vmask_full)
878#endif
879!
880! Do not perform the actual time stepping during the auxiliary
881! (nfast(ng)+1) time step.
882!
883 IF (iif(ng).gt.nfast(ng)) RETURN
884!
885!=======================================================================
886! Time step free-surface equation.
887!=======================================================================
888!
889! During the first time-step, the predictor step is Forward-Euler
890! and the corrector step is Backward-Euler. Otherwise, the predictor
891! step is Leap-frog and the corrector step is Adams-Moulton.
892#if defined VAR_RHO_2D && defined SOLVE3D
893! Recall that the vertical averaged density (rhoA) and density
894! pertubation (rhoS) are nondimensional quantities.
895!
896 fac=1000.0_r8/rho0 ! nondimensional
897#endif
898!
899 IF (first_2d_step) THEN
900 cff1=dtfast(ng)
901 DO j=jstrv-1,jend
902 DO i=istru-1,iend
903 rhs_zeta(i,j)=(duon(i,j)-duon(i+1,j))+ &
904 & (dvom(i,j)-dvom(i,j+1))
905 zeta_new(i,j)=zeta(i,j,kstp)+ &
906 & pm(i,j)*pn(i,j)*cff1*rhs_zeta(i,j)
907#ifdef MASKING
908 zeta_new(i,j)=zeta_new(i,j)*rmask(i,j)
909#endif
910 dnew(i,j)=zeta_new(i,j)+h(i,j)
911!
912 zwrk(i,j)=0.5_r8*(zeta(i,j,kstp)+zeta_new(i,j))
913#if defined VAR_RHO_2D && defined SOLVE3D
914 gzeta(i,j)=(fac+rhos(i,j))*zwrk(i,j)
915 gzeta2(i,j)=gzeta(i,j)*zwrk(i,j)
916 gzetasa(i,j)=zwrk(i,j)*(rhos(i,j)-rhoa(i,j))
917#else
918 gzeta(i,j)=zwrk(i,j)
919 gzeta2(i,j)=zwrk(i,j)*zwrk(i,j)
920#endif
921 END DO
922 END DO
923 ELSE IF (predictor_2d_step(ng)) THEN
924 cff1=2.0_r8*dtfast(ng)
925 cff4=4.0_r8/25.0_r8
926 cff5=1.0_r8-2.0_r8*cff4
927 DO j=jstrv-1,jend
928 DO i=istru-1,iend
929 rhs_zeta(i,j)=(duon(i,j)-duon(i+1,j))+ &
930 & (dvom(i,j)-dvom(i,j+1))
931 zeta_new(i,j)=zeta(i,j,kstp)+ &
932 & pm(i,j)*pn(i,j)*cff1*rhs_zeta(i,j)
933#ifdef MASKING
934 zeta_new(i,j)=zeta_new(i,j)*rmask(i,j)
935#endif
936 dnew(i,j)=zeta_new(i,j)+h(i,j)
937!
938 zwrk(i,j)=cff5*zeta(i,j,krhs)+ &
939 & cff4*(zeta(i,j,kstp)+zeta_new(i,j))
940#if defined VAR_RHO_2D && defined SOLVE3D
941 gzeta(i,j)=(fac+rhos(i,j))*zwrk(i,j)
942 gzeta2(i,j)=gzeta(i,j)*zwrk(i,j)
943 gzetasa(i,j)=zwrk(i,j)*(rhos(i,j)-rhoa(i,j))
944#else
945 gzeta(i,j)=zwrk(i,j)
946 gzeta2(i,j)=zwrk(i,j)*zwrk(i,j)
947#endif
948 END DO
949 END DO
950 ELSE IF (corrector_2d_step) THEN
951 cff1=dtfast(ng)*5.0_r8/12.0_r8
952 cff2=dtfast(ng)*8.0_r8/12.0_r8
953 cff3=dtfast(ng)*1.0_r8/12.0_r8
954 cff4=2.0_r8/5.0_r8
955 cff5=1.0_r8-cff4
956 DO j=jstrv-1,jend
957 DO i=istru-1,iend
958 cff=cff1*((duon(i,j)-duon(i+1,j))+ &
959 & (dvom(i,j)-dvom(i,j+1)))
960 zeta_new(i,j)=zeta(i,j,kstp)+ &
961 & pm(i,j)*pn(i,j)*(cff+ &
962 & cff2*rzeta(i,j,kstp)- &
963 & cff3*rzeta(i,j,ptsk))
964#ifdef MASKING
965 zeta_new(i,j)=zeta_new(i,j)*rmask(i,j)
966#endif
967 dnew(i,j)=zeta_new(i,j)+h(i,j)
968!
969 zwrk(i,j)=cff5*zeta_new(i,j)+cff4*zeta(i,j,krhs)
970#if defined VAR_RHO_2D && defined SOLVE3D
971 gzeta(i,j)=(fac+rhos(i,j))*zwrk(i,j)
972 gzeta2(i,j)=gzeta(i,j)*zwrk(i,j)
973 gzetasa(i,j)=zwrk(i,j)*(rhos(i,j)-rhoa(i,j))
974#else
975 gzeta(i,j)=zwrk(i,j)
976 gzeta2(i,j)=zwrk(i,j)*zwrk(i,j)
977#endif
978 END DO
979 END DO
980 END IF
981!
982! Load new free-surface values into shared array at both predictor
983! and corrector steps.
984#ifdef WET_DRY
985! Modify new free-surface to Ensure that depth is > Dcrit for masked
986! cells.
987#endif
988!
989 DO j=jstr,jend
990 DO i=istr,iend
991 zeta(i,j,knew)=zeta_new(i,j)
992#if defined WET_DRY && defined MASKING
993 zeta(i,j,knew)=zeta(i,j,knew)+ &
994 & (dcrit(ng)-h(i,j))*(1.0_r8-rmask(i,j))
995#endif
996 END DO
997 END DO
998!
999! If predictor step, load right-side-term into shared array.
1000!
1001 IF (predictor_2d_step(ng)) THEN
1002 DO j=jstr,jend
1003 DO i=istr,iend
1004 rzeta(i,j,krhs)=rhs_zeta(i,j)
1005 END DO
1006 END DO
1007 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
1008 CALL exchange_r2d_tile (ng, tile, &
1009 & lbi, ubi, lbj, ubj, &
1010 & rzeta(:,:,krhs))
1011 END IF
1012#ifdef DISTRIBUTE
1013 CALL mp_exchange2d (ng, tile, inlm, 1, &
1014 & lbi, ubi, lbj, ubj, &
1015 & nghostpoints, &
1016 & ewperiodic(ng), nsperiodic(ng), &
1017 & rzeta(:,:,krhs))
1018#endif
1019 END IF
1020!
1021! Apply mass point sources (volume vertical influx), if any.
1022!
1023! Dsrc(is) = 2, flow across grid cell w-face (positive or negative)
1024!
1025 IF (lwsrc(ng)) THEN
1026 DO is=1,nsrc(ng)
1027 IF (int(sources(ng)%Dsrc(is)).eq.2) THEN
1028 i=sources(ng)%Isrc(is)
1029 j=sources(ng)%Jsrc(is)
1030 IF (((istrr.le.i).and.(i.le.iendr)).and. &
1031 & ((jstrr.le.j).and.(j.le.jendr))) THEN
1032 zeta(i,j,knew)=zeta(i,j,knew)+ &
1033 & sources(ng)%Qbar(is)* &
1034 & pm(i,j)*pn(i,j)*dtfast(ng)
1035 END IF
1036 END IF
1037 END DO
1038 END IF
1039
1040#if defined SEDIMENT && defined SED_MORPH
1041!
1042! Scale the bed change with the fast time stepping. The half is
1043! becasue we do predictor and corrector. The "ndtfast/nfast" is
1044! becasue we do "nfast" steps to here.
1045!
1046 fac=0.5_r8*dtfast(ng)*ndtfast(ng)/(nfast(ng)*dt(ng))
1047 DO j=jstr,jend
1048 DO i=istr,iend
1049 cff=fac*(bed_thick(i,j,nstp)-bed_thick(i,j,nnew))
1050 h(i,j)=h(i,j)-cff
1051 END DO
1052 END DO
1053#endif
1054!
1055! Set free-surface lateral boundary conditions.
1056!
1057 CALL zetabc_tile (ng, tile, &
1058 & lbi, ubi, lbj, ubj, &
1059 & imins, imaxs, jmins, jmaxs, &
1060 & krhs, kstp, knew, &
1061 & zeta)
1062 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
1063 CALL exchange_r2d_tile (ng, tile, &
1064 & lbi, ubi, lbj, ubj, &
1065 & zeta(:,:,knew))
1066 END IF
1067#ifdef DISTRIBUTE
1068 CALL mp_exchange2d (ng, tile, inlm, 1, &
1069 & lbi, ubi, lbj, ubj, &
1070 & nghostpoints, &
1071 & ewperiodic(ng), nsperiodic(ng), &
1072 & zeta(:,:,knew))
1073#endif
1074!
1075!=======================================================================
1076! Compute right-hand-side for the 2D momentum equations.
1077!=======================================================================
1078!
1079!-----------------------------------------------------------------------
1080! Compute pressure gradient terms.
1081!-----------------------------------------------------------------------
1082!
1083 cff1=0.5_r8*g
1084 cff2=1.0_r8/3.0_r8
1085#if !defined SOLVE3D && defined ATM_PRESS
1086 fac3=0.5_r8*100.0_r8/rho0
1087#endif
1088 DO j=jstr,jend
1089 DO i=istru,iend
1090 rhs_ubar(i,j)=cff1*on_u(i,j)* &
1091 & ((h(i-1,j)+ &
1092 & h(i ,j))* &
1093 & (gzeta(i-1,j)- &
1094 & gzeta(i ,j))+ &
1095#if defined VAR_RHO_2D && defined SOLVE3D
1096 & (h(i-1,j)- &
1097 & h(i ,j))* &
1098 & (gzetasa(i-1,j)+ &
1099 & gzetasa(i ,j)+ &
1100 & cff2*(rhoa(i-1,j)- &
1101 & rhoa(i ,j))* &
1102 & (zwrk(i-1,j)- &
1103 & zwrk(i ,j)))+ &
1104#endif
1105 & (gzeta2(i-1,j)- &
1106 & gzeta2(i ,j)))
1107#if defined ATM_PRESS && !defined SOLVE3D
1108 rhs_ubar(i,j)=rhs_ubar(i,j)- &
1109 & fac3*on_u(i,j)* &
1110 & (h(i-1,j)+h(i,j)+ &
1111 & gzeta(i-1,j)+gzeta(i,j))* &
1112 & (pair(i,j)-pair(i-1,j))
1113#endif
1114#if defined TIDE_GENERATING_FORCES && !defined SOLVE3D
1115 rhs_ubar(i,j)=rhs_ubar(i,j)- &
1116 & cff1*on_u(i,j)* &
1117 & (h(i-1,j)+h(i,j)+ &
1118 & gzeta(i-1,j)+gzeta(i,j))* &
1119 & (eq_tide(i,j)-eq_tide(i-1,j))
1120#endif
1121#ifdef DIAGNOSTICS_UV
1122 diau2rhs(i,j,m2pgrd)=rhs_ubar(i,j)
1123#endif
1124#if defined WEC_VF
1125 cff3=0.5_r8*on_u(i,j)* &
1126 & (h(i-1,j)+h(i,j)+ &
1127 & gzeta(i-1,j)+gzeta(i,j))
1128 cff4=cff3*g*(zetaw(i-1,j)-zetaw(i,j))
1129 cff5=cff3*g*(qsp(i-1,j)-qsp(i,j))
1130 cff6=cff3*(bh(i-1,j)-bh(i,j))
1131 cff7=rukvf2d(i,j)
1132 rhs_ubar(i,j)=rhs_ubar(i,j)-cff4-cff5+cff6+cff7
1133# ifdef DIAGNOSTICS_UV
1134 diau2rhs(i,j,m2zeta)=diau2rhs(i,j,m2pgrd)
1135 diau2rhs(i,j,m2pgrd)=diau2rhs(i,j,m2pgrd)-cff4-cff5+cff6
1136 diau2rhs(i,j,m2zetw)=-cff4
1137 diau2rhs(i,j,m2zqsp)=-cff5
1138 diau2rhs(i,j,m2zbeh)=cff6
1139 diau2rhs(i,j,m2kvrf)=cff7
1140# ifndef UV_ADV
1141 diau2rhs(i,j,m2hjvf)=0.0_r8
1142# endif
1143# endif
1144#endif
1145 END DO
1146 IF (j.ge.jstrv) THEN
1147 DO i=istr,iend
1148 rhs_vbar(i,j)=cff1*om_v(i,j)* &
1149 & ((h(i,j-1)+ &
1150 & h(i,j ))* &
1151 & (gzeta(i,j-1)- &
1152 & gzeta(i,j ))+ &
1153#if defined VAR_RHO_2D && defined SOLVE3D
1154 & (h(i,j-1)- &
1155 & h(i,j ))* &
1156 & (gzetasa(i,j-1)+ &
1157 & gzetasa(i,j )+ &
1158 & cff2*(rhoa(i,j-1)- &
1159 & rhoa(i,j ))* &
1160 & (zwrk(i,j-1)- &
1161 & zwrk(i,j )))+ &
1162#endif
1163 & (gzeta2(i,j-1)- &
1164 & gzeta2(i,j )))
1165#if defined ATM_PRESS && !defined SOLVE3D
1166 rhs_vbar(i,j)=rhs_vbar(i,j)- &
1167 & fac3*om_v(i,j)* &
1168 & (h(i,j-1)+h(i,j)+ &
1169 & gzeta(i,j-1)+gzeta(i,j))* &
1170 & (pair(i,j)-pair(i,j-1))
1171#endif
1172#if defined TIDE_GENERATING_FORCES && !defined SOLVE3D
1173 rhs_vbar(i,j)=rhs_vbar(i,j)- &
1174 & cff1*om_v(i,j)* &
1175 & (h(i,j-1)+h(i,j)+ &
1176 & gzeta(i,j-1)+gzeta(i,j))* &
1177 & (eq_tide(i,j)-eq_tide(i,j-1))
1178#endif
1179#ifdef DIAGNOSTICS_UV
1180 diav2rhs(i,j,m2pgrd)=rhs_vbar(i,j)
1181#endif
1182#if defined WEC_VF
1183 cff3=0.5_r8*om_v(i,j)* &
1184 & (h(i,j-1)+h(i,j)+ &
1185 & gzeta(i,j-1)+gzeta(i,j))
1186 cff4=cff3*g*(zetaw(i,j-1)-zetaw(i,j))
1187 cff5=cff3*g*(qsp(i,j-1)-qsp(i,j))
1188 cff6=cff3*(bh(i,j-1)-bh(i,j))
1189 cff7=rvkvf2d(i,j)
1190 rhs_vbar(i,j)=rhs_vbar(i,j)-cff4-cff5+cff6+cff7
1191# ifdef DIAGNOSTICS_UV
1192 diav2rhs(i,j,m2zeta)=diav2rhs(i,j,m2pgrd)
1193 diav2rhs(i,j,m2pgrd)=diav2rhs(i,j,m2pgrd)-cff4-cff5+cff6
1194 diav2rhs(i,j,m2zetw)=-cff4
1195 diav2rhs(i,j,m2zqsp)=-cff5
1196 diav2rhs(i,j,m2zbeh)=cff6
1197 diav2rhs(i,j,m2kvrf)=cff7
1198# ifndef UV_ADV
1199 diav2rhs(i,j,m2hjvf)=0.0_r8
1200# endif
1201# endif
1202#endif
1203 END DO
1204 END IF
1205 END DO
1206#ifdef UV_ADV
1207!
1208!-----------------------------------------------------------------------
1209! Add in horizontal advection of momentum.
1210!-----------------------------------------------------------------------
1211
1212# ifdef UV_C2ADVECTION
1213!
1214! Second-order, centered differences advection fluxes.
1215!
1216 DO j=jstr,jend
1217 DO i=istru-1,iend
1218 ufx(i,j)=0.25_r8*(duon(i,j)+duon(i+1,j))* &
1219 & (ubar(i ,j,krhs)+ &
1220 & ubar(i+1,j,krhs))
1221 END DO
1222 END DO
1223!
1224 DO j=jstr,jend+1
1225 DO i=istru,iend
1226 ufe(i,j)=0.25_r8*(dvom(i,j)+dvom(i-1,j))* &
1227 & (ubar(i,j ,krhs)+ &
1228 & ubar(i,j-1,krhs))
1229 END DO
1230 END DO
1231!
1232 DO j=jstrv,jend
1233 DO i=istr,iend+1
1234 vfx(i,j)=0.25_r8*(duon(i,j)+duon(i,j-1))* &
1235 & (vbar(i ,j,krhs)+ &
1236 & vbar(i-1,j,krhs))
1237 END DO
1238 END DO
1239!
1240 DO j=jstrv-1,jend
1241 DO i=istr,iend
1242 vfe(i,j)=0.25_r8*(dvom(i,j)+dvom(i,j+1))* &
1243 & (vbar(i,j ,krhs)+ &
1244 & vbar(i,j+1,krhs))
1245 END DO
1246 END DO
1247# else
1248!
1249! Fourth-order, centered differences advection fluxes.
1250!
1251 DO j=jstr,jend
1252 DO i=istrum1,iendp1
1253 grad(i,j)=ubar(i-1,j,krhs)-2.0_r8*ubar(i,j,krhs)+ &
1254 & ubar(i+1,j,krhs)
1255 dgrad(i,j)=duon(i-1,j)-2.0_r8*duon(i,j)+duon(i+1,j)
1256 END DO
1257 END DO
1258 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
1259 IF (domain(ng)%Western_Edge(tile)) THEN
1260 DO j=jstr,jend
1261 grad(istr,j)=grad(istr+1,j)
1262 dgrad(istr,j)=dgrad(istr+1,j)
1263 END DO
1264 END IF
1265 END IF
1266 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
1267 IF (domain(ng)%Eastern_Edge(tile)) THEN
1268 DO j=jstr,jend
1269 grad(iend+1,j)=grad(iend,j)
1270 dgrad(iend+1,j)=dgrad(iend,j)
1271 END DO
1272 END IF
1273 END IF
1274
1275 cff=1.0_r8/6.0_r8
1276 DO j=jstr,jend
1277 DO i=istru-1,iend
1278 ufx(i,j)=0.25_r8*(ubar(i ,j,krhs)+ &
1279 & ubar(i+1,j,krhs)- &
1280 & cff*(grad(i,j)+grad(i+1,j)))* &
1281 & (duon(i,j)+duon(i+1,j)- &
1282 & cff*(dgrad(i,j)+dgrad(i+1,j)))
1283 END DO
1284 END DO
1285!
1286 DO j=jstrm1,jendp1
1287 DO i=istru,iend
1288 grad(i,j)=ubar(i,j-1,krhs)-2.0_r8*ubar(i,j,krhs)+ &
1289 & ubar(i,j+1,krhs)
1290 END DO
1291 END DO
1292 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
1293 IF (domain(ng)%Southern_Edge(tile)) THEN
1294 DO i=istru,iend
1295 grad(i,jstr-1)=grad(i,jstr)
1296 END DO
1297 END IF
1298 END IF
1299 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
1300 IF (domain(ng)%Northern_Edge(tile)) THEN
1301 DO i=istru,iend
1302 grad(i,jend+1)=grad(i,jend)
1303 END DO
1304 END IF
1305 END IF
1306 DO j=jstr,jend+1
1307 DO i=istru-1,iend
1308 dgrad(i,j)=dvom(i-1,j)-2.0_r8*dvom(i,j)+dvom(i+1,j)
1309 END DO
1310 END DO
1311
1312 cff=1.0_r8/6.0_r8
1313 DO j=jstr,jend+1
1314 DO i=istru,iend
1315 ufe(i,j)=0.25_r8*(ubar(i,j ,krhs)+ &
1316 & ubar(i,j-1,krhs)- &
1317 & cff*(grad(i,j)+grad(i,j-1)))* &
1318 & (dvom(i,j)+dvom(i-1,j)- &
1319 & cff*(dgrad(i,j)+dgrad(i-1,j)))
1320 END DO
1321 END DO
1322!
1323 DO j=jstrv,jend
1324 DO i=istrm1,iendp1
1325 grad(i,j)=vbar(i-1,j,krhs)-2.0_r8*vbar(i,j,krhs)+ &
1326 & vbar(i+1,j,krhs)
1327 END DO
1328 END DO
1329 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
1330 IF (domain(ng)%Western_Edge(tile)) THEN
1331 DO j=jstrv,jend
1332 grad(istr-1,j)=grad(istr,j)
1333 END DO
1334 END IF
1335 END IF
1336 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
1337 IF (domain(ng)%Eastern_Edge(tile)) THEN
1338 DO j=jstrv,jend
1339 grad(iend+1,j)=grad(iend,j)
1340 END DO
1341 END IF
1342 END IF
1343 DO j=jstrv-1,jend
1344 DO i=istr,iend+1
1345 dgrad(i,j)=duon(i,j-1)-2.0_r8*duon(i,j)+duon(i,j+1)
1346 END DO
1347 END DO
1348
1349 cff=1.0_r8/6.0_r8
1350 DO j=jstrv,jend
1351 DO i=istr,iend+1
1352 vfx(i,j)=0.25_r8*(vbar(i ,j,krhs)+ &
1353 & vbar(i-1,j,krhs)- &
1354 & cff*(grad(i,j)+grad(i-1,j)))* &
1355 & (duon(i,j)+duon(i,j-1)- &
1356 & cff*(dgrad(i,j)+dgrad(i,j-1)))
1357 END DO
1358 END DO
1359!
1360 DO j=jstrvm1,jendp1
1361 DO i=istr,iend
1362 grad(i,j)=vbar(i,j-1,krhs)-2.0_r8*vbar(i,j,krhs)+ &
1363 & vbar(i,j+1,krhs)
1364 dgrad(i,j)=dvom(i,j-1)-2.0_r8*dvom(i,j)+dvom(i,j+1)
1365 END DO
1366 END DO
1367 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
1368 IF (domain(ng)%Southern_Edge(tile)) THEN
1369 DO i=istr,iend
1370 grad(i,jstr)=grad(i,jstr+1)
1371 dgrad(i,jstr)=dgrad(i,jstr+1)
1372 END DO
1373 END IF
1374 END IF
1375 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
1376 IF (domain(ng)%Northern_Edge(tile)) THEN
1377 DO i=istr,iend
1378 grad(i,jend+1)=grad(i,jend)
1379 dgrad(i,jend+1)=dgrad(i,jend)
1380 END DO
1381 END IF
1382 END IF
1383
1384 cff=1.0_r8/6.0_r8
1385 DO j=jstrv-1,jend
1386 DO i=istr,iend
1387 vfe(i,j)=0.25_r8*(vbar(i,j ,krhs)+ &
1388 & vbar(i,j+1,krhs)- &
1389 & cff*(grad(i,j)+grad(i,j+1)))* &
1390 & (dvom(i,j)+dvom(i,j+1)- &
1391 & cff*(dgrad(i,j)+dgrad(i,j+1)))
1392 END DO
1393 END DO
1394# endif
1395!
1396! Add advection to RHS terms.
1397!
1398 DO j=jstr,jend
1399 DO i=istru,iend
1400 cff1=ufx(i,j)-ufx(i-1,j)
1401 cff2=ufe(i,j+1)-ufe(i,j)
1402 fac=cff1+cff2
1403 rhs_ubar(i,j)=rhs_ubar(i,j)-fac
1404# if defined DIAGNOSTICS_UV
1405 diau2rhs(i,j,m2xadv)=-cff1
1406 diau2rhs(i,j,m2yadv)=-cff2
1407 diau2rhs(i,j,m2hadv)=-fac
1408# endif
1409 END DO
1410 END DO
1411 DO j=jstrv,jend
1412 DO i=istr,iend
1413 cff1=vfx(i+1,j)-vfx(i,j)
1414 cff2=vfe(i,j)-vfe(i,j-1)
1415 fac=cff1+cff2
1416 rhs_vbar(i,j)=rhs_vbar(i,j)-fac
1417# if defined DIAGNOSTICS_UV
1418 diav2rhs(i,j,m2xadv)=-cff1
1419 diav2rhs(i,j,m2yadv)=-cff2
1420 diav2rhs(i,j,m2hadv)=-fac
1421# endif
1422 END DO
1423 END DO
1424#endif
1425
1426#ifdef UV_COR
1427!
1428!-----------------------------------------------------------------------
1429! Add in Coriolis term.
1430!-----------------------------------------------------------------------
1431!
1432 DO j=jstrv-1,jend
1433 DO i=istru-1,iend
1434 cff=0.5_r8*drhs(i,j)*fomn(i,j)
1435 ufx(i,j)=cff*(vbar(i,j ,krhs)+ &
1436 & vbar(i,j+1,krhs))
1437 vfe(i,j)=cff*(ubar(i ,j,krhs)+ &
1438 & ubar(i+1,j,krhs))
1439 END DO
1440 END DO
1441 DO j=jstr,jend
1442 DO i=istru,iend
1443 fac1=0.5_r8*(ufx(i,j)+ufx(i-1,j))
1444 rhs_ubar(i,j)=rhs_ubar(i,j)+fac1
1445# if defined DIAGNOSTICS_UV
1446 diau2rhs(i,j,m2fcor)=fac1
1447# endif
1448 END DO
1449 END DO
1450 DO j=jstrv,jend
1451 DO i=istr,iend
1452 fac1=0.5_r8*(vfe(i,j)+vfe(i,j-1))
1453 rhs_vbar(i,j)=rhs_vbar(i,j)-fac1
1454# if defined DIAGNOSTICS_UV
1455 diav2rhs(i,j,m2fcor)=-fac1
1456# endif
1457 END DO
1458 END DO
1459!
1460# ifdef WEC
1461 DO j=jstrv-1,jend
1462 DO i=istru-1,iend
1463 cff=0.5_r8*drhs(i,j)*fomn(i,j)
1464 ufx(i,j)=cff*(vbar_stokes(i,j )+ &
1465 & vbar_stokes(i,j+1))
1466 vfe(i,j)=cff*(ubar_stokes(i ,j)+ &
1467 & ubar_stokes(i+1,j))
1468 END DO
1469 END DO
1470 DO j=jstr,jend
1471 DO i=istru,iend
1472 fac1=0.5_r8*(ufx(i,j)+ufx(i-1,j))
1473 rhs_ubar(i,j)=rhs_ubar(i,j)+fac1
1474# if defined DIAGNOSTICS_UV
1475 diau2rhs(i,j,m2fsco)=fac1
1476# endif
1477 END DO
1478 END DO
1479 DO j=jstrv,jend
1480 DO i=istr,iend
1481 fac1=0.5_r8*(vfe(i,j)+vfe(i,j-1))
1482 rhs_vbar(i,j)=rhs_vbar(i,j)-fac1
1483# if defined DIAGNOSTICS_UV
1484 diav2rhs(i,j,m2fsco)=-fac1
1485# endif
1486 END DO
1487 END DO
1488# endif
1489#endif
1490
1491#if defined CURVGRID && defined UV_ADV
1492!
1493!-----------------------------------------------------------------------
1494! Add in curvilinear transformation terms.
1495!-----------------------------------------------------------------------
1496!
1497 DO j=jstrv-1,jend
1498 DO i=istru-1,iend
1499 cff1=0.5_r8*(vbar(i,j ,krhs)+ &
1500# ifdef WEC
1501 & vbar_stokes(i,j )+ &
1502 & vbar_stokes(i,j+1)+ &
1503# endif
1504 & vbar(i,j+1,krhs))
1505 cff2=0.5_r8*(ubar(i ,j,krhs)+ &
1506# ifdef WEC
1507 & ubar_stokes(i ,j)+ &
1508 & ubar_stokes(i+1,j)+ &
1509# endif
1510 & ubar(i+1,j,krhs))
1511 cff3=cff1*dndx(i,j)
1512 cff4=cff2*dmde(i,j)
1513# ifdef WEC_VF
1514 cff5=0.5_r8*(vbar_stokes(i,j )+ &
1515 & vbar_stokes(i,j+1))
1516 cff6=0.5_r8*(ubar_stokes(i ,j)+ &
1517 & ubar_stokes(i+1,j))
1518 cff7=cff5*dndx(i,j)
1519 cff8=cff6*dmde(i,j)
1520# endif
1521 cff=drhs(i,j)*(cff3-cff4)
1522 ufx(i,j)=cff*cff1
1523 vfe(i,j)=cff*cff2
1524# ifdef WEC_VF
1525 ufx(i,j)=ufx(i,j)-(cff5*drhs(i,j)*(cff7-cff8))
1526 vfe(i,j)=vfe(i,j)-(cff6*drhs(i,j)*(cff7-cff8))
1527# endif
1528# if defined DIAGNOSTICS_UV
1529 cff=drhs(i,j)*cff4
1530 uwrk(i,j)=-cff*cff1 ! ubar equation, ETA-term
1531 vwrk(i,j)=-cff*cff2 ! vbar equation, ETA-term
1532# ifdef WEC_VF
1533 uwrk(i,j)=uwrk(i,j)+drhs(i,j)*cff5*cff8
1534 vwrk(i,j)=vwrk(i,j)-drhs(i,j)*cff6*cff8
1535# endif
1536# endif
1537 END DO
1538 END DO
1539 DO j=jstr,jend
1540 DO i=istru,iend
1541 fac1=0.5_r8*(ufx(i,j)+ufx(i-1,j))
1542 rhs_ubar(i,j)=rhs_ubar(i,j)+fac1
1543# if defined DIAGNOSTICS_UV
1544 fac2=0.5_r8*(uwrk(i,j)+uwrk(i-1,j))
1545 diau2rhs(i,j,m2xadv)=diau2rhs(i,j,m2xadv)+fac1-fac2
1546 diau2rhs(i,j,m2yadv)=diau2rhs(i,j,m2yadv)+fac2
1547 diau2rhs(i,j,m2hadv)=diau2rhs(i,j,m2hadv)+fac1
1548# endif
1549 END DO
1550 END DO
1551 DO j=jstrv,jend
1552 DO i=istr,iend
1553 fac1=0.5_r8*(vfe(i,j)+vfe(i,j-1))
1554 rhs_vbar(i,j)=rhs_vbar(i,j)-fac1
1555# if defined DIAGNOSTICS_UV
1556 fac2=0.5_r8*(vwrk(i,j)+vwrk(i,j-1))
1557 diav2rhs(i,j,m2xadv)=diav2rhs(i,j,m2xadv)-fac1+fac2
1558 diav2rhs(i,j,m2yadv)=diav2rhs(i,j,m2yadv)-fac2
1559 diav2rhs(i,j,m2hadv)=diav2rhs(i,j,m2hadv)-fac1
1560# endif
1561 END DO
1562 END DO
1563#endif
1564#if defined UV_VIS2 || defined UV_VIS4
1565!
1566!-----------------------------------------------------------------------
1567! If horizontal mixing, compute total depth at PSI-points.
1568!-----------------------------------------------------------------------
1569!
1570# ifdef UV_VIS4
1571 DO j=jstrm1,jendp2
1572 DO i=istrm1,iendp2
1573# else
1574 DO j=jstr,jend+1
1575 DO i=istr,iend+1
1576# endif
1577 drhs_p(i,j)=0.25_r8*(drhs(i,j )+drhs(i-1,j )+ &
1578 & drhs(i,j-1)+drhs(i-1,j-1))
1579 END DO
1580 END DO
1581#endif
1582#ifdef UV_VIS2
1583!
1584!-----------------------------------------------------------------------
1585! Add in horizontal harmonic viscosity.
1586!-----------------------------------------------------------------------
1587!
1588! Compute flux-components of the horizontal divergence of the stress
1589! tensor (m5/s2) in XI- and ETA-directions.
1590!
1591 DO j=jstrv-1,jend
1592 DO i=istru-1,iend
1593 cff=visc2_r(i,j)*drhs(i,j)*0.5_r8* &
1594 & (pmon_r(i,j)* &
1595 & ((pn(i ,j)+pn(i+1,j))*ubar(i+1,j,krhs)- &
1596 & (pn(i-1,j)+pn(i ,j))*ubar(i ,j,krhs))- &
1597 & pnom_r(i,j)* &
1598 & ((pm(i,j )+pm(i,j+1))*vbar(i,j+1,krhs)- &
1599 & (pm(i,j-1)+pm(i,j ))*vbar(i,j ,krhs)))
1600 ufx(i,j)=on_r(i,j)*on_r(i,j)*cff
1601 vfe(i,j)=om_r(i,j)*om_r(i,j)*cff
1602 END DO
1603 END DO
1604 DO j=jstr,jend+1
1605 DO i=istr,iend+1
1606 cff=visc2_p(i,j)*drhs_p(i,j)*0.5_r8* &
1607 & (pmon_p(i,j)* &
1608 & ((pn(i ,j-1)+pn(i ,j))*vbar(i ,j,krhs)- &
1609 & (pn(i-1,j-1)+pn(i-1,j))*vbar(i-1,j,krhs))+ &
1610 & pnom_p(i,j)* &
1611 & ((pm(i-1,j )+pm(i,j ))*ubar(i,j ,krhs)- &
1612 & (pm(i-1,j-1)+pm(i,j-1))*ubar(i,j-1,krhs)))
1613# ifdef MASKING
1614 cff=cff*pmask(i,j)
1615# endif
1616# ifdef WET_DRY
1617 cff=cff*pmask_wet(i,j)
1618# endif
1619 ufe(i,j)=om_p(i,j)*om_p(i,j)*cff
1620 vfx(i,j)=on_p(i,j)*on_p(i,j)*cff
1621 END DO
1622 END DO
1623!
1624! Add in harmonic viscosity.
1625!
1626 DO j=jstr,jend
1627 DO i=istru,iend
1628 cff1=0.5_r8*(pn(i-1,j)+pn(i,j))*(ufx(i,j )-ufx(i-1,j))
1629 cff2=0.5_r8*(pm(i-1,j)+pm(i,j))*(ufe(i,j+1)-ufe(i ,j))
1630 fac=cff1+cff2
1631 rhs_ubar(i,j)=rhs_ubar(i,j)+fac
1632# if defined DIAGNOSTICS_UV
1633 diau2rhs(i,j,m2hvis)=fac
1634 diau2rhs(i,j,m2xvis)=cff1
1635 diau2rhs(i,j,m2yvis)=cff2
1636# endif
1637 END DO
1638 END DO
1639 DO j=jstrv,jend
1640 DO i=istr,iend
1641 cff1=0.5_r8*(pn(i,j-1)+pn(i,j))*(vfx(i+1,j)-vfx(i,j ))
1642 cff2=0.5_r8*(pm(i,j-1)+pm(i,j))*(vfe(i ,j)-vfe(i,j-1))
1643 fac=cff1-cff2
1644 rhs_vbar(i,j)=rhs_vbar(i,j)+fac
1645# if defined DIAGNOSTICS_UV
1646 diav2rhs(i,j,m2hvis)=fac
1647 diav2rhs(i,j,m2xvis)= cff1
1648 diav2rhs(i,j,m2yvis)=-cff2
1649# endif
1650 END DO
1651 END DO
1652#endif
1653#ifdef UV_VIS4
1654!
1655!-----------------------------------------------------------------------
1656! Add in horizontal biharmonic viscosity. The biharmonic operator
1657! is computed by applying the harmonic operator twice.
1658!-----------------------------------------------------------------------
1659!
1660! Compute flux-components of the horizontal divergence of the stress
1661! tensor (m4 s^-3/2) in XI- and ETA-directions. It is assumed here
1662! that "visc4_r" and "visc4_p" are the squared root of the biharmonic
1663! viscosity coefficient. For momentum balance purposes, the total
1664! thickness "D" appears only when computing the second harmonic
1665! operator.
1666!
1667 DO j=jstrvm2,jendp1
1668 DO i=istrum2,iendp1
1669 cff=visc4_r(i,j)*0.5_r8* &
1670 & (pmon_r(i,j)* &
1671 & ((pn(i ,j)+pn(i+1,j))*ubar(i+1,j,krhs)- &
1672 & (pn(i-1,j)+pn(i ,j))*ubar(i ,j,krhs))- &
1673 & pnom_r(i,j)* &
1674 & ((pm(i,j )+pm(i,j+1))*vbar(i,j+1,krhs)- &
1675 & (pm(i,j-1)+pm(i,j ))*vbar(i,j ,krhs)))
1676 ufx(i,j)=on_r(i,j)*on_r(i,j)*cff
1677 vfe(i,j)=om_r(i,j)*om_r(i,j)*cff
1678 END DO
1679 END DO
1680 DO j=jstrm1,jendp2
1681 DO i=istrm1,iendp2
1682 cff=visc4_p(i,j)*0.5_r8* &
1683 & (pmon_p(i,j)* &
1684 & ((pn(i ,j-1)+pn(i ,j))*vbar(i ,j,krhs)- &
1685 & (pn(i-1,j-1)+pn(i-1,j))*vbar(i-1,j,krhs))+ &
1686 & pnom_p(i,j)* &
1687 & ((pm(i-1,j )+pm(i,j ))*ubar(i,j ,krhs)- &
1688 & (pm(i-1,j-1)+pm(i,j-1))*ubar(i,j-1,krhs)))
1689# ifdef MASKING
1690 cff=cff*pmask(i,j)
1691# endif
1692# ifdef WET_DRY
1693 cff=cff*pmask_wet(i,j)
1694# endif
1695 ufe(i,j)=om_p(i,j)*om_p(i,j)*cff
1696 vfx(i,j)=on_p(i,j)*on_p(i,j)*cff
1697 END DO
1698 END DO
1699!
1700! Compute first harmonic operator (m s^-3/2).
1701!
1702 DO j=jstrm1,jendp1
1703 DO i=istrum1,iendp1
1704 lapu(i,j)=0.125_r8* &
1705 & (pm(i-1,j)+pm(i,j))*(pn(i-1,j)+pn(i,j))* &
1706 & ((pn(i-1,j)+pn(i,j))* &
1707 & (ufx(i,j )-ufx(i-1,j))+ &
1708 & (pm(i-1,j)+pm(i,j))* &
1709 & (ufe(i,j+1)-ufe(i ,j)))
1710 END DO
1711 END DO
1712 DO j=jstrvm1,jendp1
1713 DO i=istrm1,iendp1
1714 lapv(i,j)=0.125_r8* &
1715 & (pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1))* &
1716 & ((pn(i,j-1)+pn(i,j))* &
1717 & (vfx(i+1,j)-vfx(i,j ))- &
1718 & (pm(i,j-1)+pm(i,j))* &
1719 & (vfe(i ,j)-vfe(i,j-1)))
1720 END DO
1721 END DO
1722!
1723! Apply boundary conditions (other than periodic) to the first
1724! harmonic operator. These are gradient or closed (free slip or
1725! no slip) boundary conditions.
1726!
1727 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
1728 IF (domain(ng)%Western_Edge(tile)) THEN
1729 IF (lbc(iwest,isubar,ng)%closed) THEN
1730 DO j=jstrm1,jendp1
1731 lapu(istru-1,j)=0.0_r8
1732 END DO
1733 ELSE
1734 DO j=jstrm1,jendp1
1735 lapu(istru-1,j)=lapu(istru,j)
1736 END DO
1737 END IF
1738 IF (lbc(iwest,isvbar,ng)%closed) THEN
1739 DO j=jstrvm1,jendp1
1740 lapv(istr-1,j)=gamma2(ng)*lapv(istr,j)
1741 END DO
1742 ELSE
1743 DO j=jstrvm1,jendp1
1744 lapv(istr-1,j)=0.0_r8
1745 END DO
1746 END IF
1747 END IF
1748 END IF
1749!
1750 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
1751 IF (domain(ng)%Eastern_Edge(tile)) THEN
1752 IF (lbc(ieast,isubar,ng)%closed) THEN
1753 DO j=jstrm1,jendp1
1754 lapu(iend+1,j)=0.0_r8
1755 END DO
1756 ELSE
1757 DO j=jstrm1,jendp1
1758 lapu(iend+1,j)=lapu(iend,j)
1759 END DO
1760 END IF
1761 IF (lbc(ieast,isvbar,ng)%closed) THEN
1762 DO j=jstrvm1,jendp1
1763 lapv(iend+1,j)=gamma2(ng)*lapv(iend,j)
1764 END DO
1765 ELSE
1766 DO j=jstrvm1,jendp1
1767 lapv(iend+1,j)=0.0_r8
1768 END DO
1769 END IF
1770 END IF
1771 END IF
1772!
1773 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
1774 IF (domain(ng)%Southern_Edge(tile)) THEN
1775 IF (lbc(isouth,isubar,ng)%closed) THEN
1776 DO i=istrum1,iendp1
1777 lapu(i,jstr-1)=gamma2(ng)*lapu(i,jstr)
1778 END DO
1779 ELSE
1780 DO i=istrum1,iendp1
1781 lapu(i,jstr-1)=0.0_r8
1782 END DO
1783 END IF
1784 IF (lbc(isouth,isvbar,ng)%closed) THEN
1785 DO i=istrm1,iendp1
1786 lapv(i,jstrv-1)=0.0_r8
1787 END DO
1788 ELSE
1789 DO i=istrm1,iendp1
1790 lapv(i,jstrv-1)=lapv(i,jstrv)
1791 END DO
1792 END IF
1793 END IF
1794 END IF
1795!
1796 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
1797 IF (domain(ng)%Northern_Edge(tile)) THEN
1798 IF (lbc(inorth,isubar,ng)%closed) THEN
1799 DO i=istrum1,iendp1
1800 lapu(i,jend+1)=gamma2(ng)*lapu(i,jend)
1801 END DO
1802 ELSE
1803 DO i=istrum1,iendp1
1804 lapu(i,jend+1)=0.0_r8
1805 END DO
1806 END IF
1807 IF (lbc(inorth,isvbar,ng)%closed) THEN
1808 DO i=istrm1,iendp1
1809 lapv(i,jend+1)=0.0_r8
1810 END DO
1811 ELSE
1812 DO i=istrm1,iendp1
1813 lapv(i,jend+1)=lapv(i,jend)
1814 END DO
1815 END IF
1816 END IF
1817 END IF
1818!
1819 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng).or. &
1820 & compositegrid(iwest ,ng).or.ewperiodic(ng))) THEN
1821 IF (domain(ng)%SouthWest_Corner(tile)) THEN
1822 lapu(istr ,jstr-1)=0.5_r8*(lapu(istr+1,jstr-1)+ &
1823 & lapu(istr ,jstr ))
1824 lapv(istr-1,jstr )=0.5_r8*(lapv(istr-1,jstr+1)+ &
1825 & lapv(istr ,jstr ))
1826 END IF
1827 END IF
1828
1829 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng).or. &
1830 & compositegrid(ieast ,ng).or.ewperiodic(ng))) THEN
1831 IF (domain(ng)%SouthEast_Corner(tile)) THEN
1832 lapu(iend+1,jstr-1)=0.5_r8*(lapu(iend ,jstr-1)+ &
1833 & lapu(iend+1,jstr ))
1834 lapv(iend+1,jstr )=0.5_r8*(lapv(iend ,jstr )+ &
1835 & lapv(iend+1,jstr+1))
1836 END IF
1837 END IF
1838
1839 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng).or. &
1840 & compositegrid(iwest ,ng).or.ewperiodic(ng))) THEN
1841 IF (domain(ng)%NorthWest_Corner(tile)) THEN
1842 lapu(istr ,jend+1)=0.5_r8*(lapu(istr+1,jend+1)+ &
1843 & lapu(istr ,jend ))
1844 lapv(istr-1,jend+1)=0.5_r8*(lapv(istr ,jend+1)+ &
1845 & lapv(istr-1,jend ))
1846 END IF
1847 END IF
1848
1849 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng).or. &
1850 & compositegrid(ieast ,ng).or.ewperiodic(ng))) THEN
1851 IF (domain(ng)%NorthEast_Corner(tile)) THEN
1852 lapu(iend+1,jend+1)=0.5_r8*(lapu(iend ,jend+1)+ &
1853 & lapu(iend+1,jend ))
1854 lapv(iend+1,jend+1)=0.5_r8*(lapv(iend ,jend+1)+ &
1855 & lapv(iend+1,jend ))
1856 END IF
1857 END IF
1858!
1859! Compute flux-components of the horizontal divergence of the
1860! biharmonic stress tensor (m4/s2) in XI- and ETA-directions.
1861!
1862 DO j=jstrv-1,jend
1863 DO i=istru-1,iend
1864 cff=visc4_r(i,j)*drhs(i,j)*0.5_r8* &
1865 & (pmon_r(i,j)* &
1866 & ((pn(i ,j)+pn(i+1,j))*lapu(i+1,j)- &
1867 & (pn(i-1,j)+pn(i ,j))*lapu(i ,j))- &
1868 & pnom_r(i,j)* &
1869 & ((pm(i,j )+pm(i,j+1))*lapv(i,j+1)- &
1870 & (pm(i,j-1)+pm(i,j ))*lapv(i,j )))
1871 ufx(i,j)=on_r(i,j)*on_r(i,j)*cff
1872 vfe(i,j)=om_r(i,j)*om_r(i,j)*cff
1873 END DO
1874 END DO
1875 DO j=jstr,jend+1
1876 DO i=istr,iend+1
1877 cff=visc4_p(i,j)*drhs_p(i,j)*0.5_r8* &
1878 & (pmon_p(i,j)* &
1879 & ((pn(i ,j-1)+pn(i ,j))*lapv(i ,j)- &
1880 & (pn(i-1,j-1)+pn(i-1,j))*lapv(i-1,j))+ &
1881 & pnom_p(i,j)* &
1882 & ((pm(i-1,j )+pm(i,j ))*lapu(i,j )- &
1883 & (pm(i-1,j-1)+pm(i,j-1))*lapu(i,j-1)))
1884# ifdef MASKING
1885 cff=cff*pmask(i,j)
1886# endif
1887# ifdef WET_DRY
1888 cff=cff*pmask_wet(i,j)
1889# endif
1890 ufe(i,j)=om_p(i,j)*om_p(i,j)*cff
1891 vfx(i,j)=on_p(i,j)*on_p(i,j)*cff
1892 END DO
1893 END DO
1894!
1895! Add in biharmonic viscosity.
1896!
1897 DO j=jstr,jend
1898 DO i=istru,iend
1899 cff1=0.5_r8*(pn(i-1,j)+pn(i,j))*(ufx(i,j )-ufx(i-1,j))
1900 cff2=0.5_r8*(pm(i-1,j)+pm(i,j))*(ufe(i,j+1)-ufe(i ,j))
1901 fac=cff1+cff2
1902 rhs_ubar(i,j)=rhs_ubar(i,j)-fac
1903# if defined DIAGNOSTICS_UV
1904 diau2rhs(i,j,m2hvis)=-fac
1905 diau2rhs(i,j,m2xvis)=-cff1
1906 diau2rhs(i,j,m2yvis)=-cff2
1907# endif
1908 END DO
1909 END DO
1910 DO j=jstrv,jend
1911 DO i=istr,iend
1912 cff1=0.5_r8*(pn(i,j-1)+pn(i,j))*(vfx(i+1,j)-vfx(i,j ))
1913 cff2=0.5_r8*(pm(i,j-1)+pm(i,j))*(vfe(i ,j)-vfe(i,j-1))
1914 fac=cff1-cff2
1915 rhs_vbar(i,j)=rhs_vbar(i,j)-fac
1916# if defined DIAGNOSTICS_UV
1917 diav2rhs(i,j,m2hvis)=-fac
1918 diav2rhs(i,j,m2xvis)=-cff1
1919 diav2rhs(i,j,m2yvis)= cff2
1920# endif
1921 END DO
1922 END DO
1923#endif
1924#if defined WEC_VF && defined SOLVE3D
1925!
1926!-----------------------------------------------------------------------
1927! Add in non-conservative roller terms.
1928!-----------------------------------------------------------------------
1929!
1930 DO j=jstr,jend
1931 DO i=istru,iend
1932 cff1=rubrk2d(i,j)
1933 rhs_ubar(i,j)=rhs_ubar(i,j)+cff1
1934# ifdef DIAGNOSTICS_UV
1935 diau2rhs(i,j,m2wbrk)=cff1
1936# endif
1937# ifdef WEC_ROLLER
1938 cff1=rurol2d(i,j)
1939 rhs_ubar(i,j)=rhs_ubar(i,j)+cff1
1940# ifdef DIAGNOSTICS_UV
1941 diau2rhs(i,j,m2wrol)=cff1
1942# endif
1943# else
1944# ifdef DIAGNOSTICS_UV
1945 diau2rhs(i,j,m2wrol)=0.0_r8
1946# endif
1947# endif
1948# ifdef BOTTOM_STREAMING
1949 cff1=rubst2d(i,j)
1950 rhs_ubar(i,j)=rhs_ubar(i,j)+cff1
1951# ifdef DIAGNOSTICS_UV
1952 diau2rhs(i,j,m2bstm)=cff1
1953# endif
1954# endif
1955# ifdef SURFACE_STREAMING
1956 cff1=russt2d(i,j)
1957 rhs_ubar(i,j)=rhs_ubar(i,j)+cff1
1958# ifdef DIAGNOSTICS_UV
1959 diau2rhs(i,j,m2sstm)=cff1
1960# endif
1961# endif
1962 END DO
1963 IF (j.ge.jstrv) THEN
1964 DO i=istr,iend
1965 cff1=rvbrk2d(i,j)
1966 rhs_vbar(i,j)=rhs_vbar(i,j)+cff1
1967# ifdef DIAGNOSTICS_UV
1968 diav2rhs(i,j,m2wbrk)=cff1
1969# endif
1970# ifdef WEC_ROLLER
1971 cff1=rvrol2d(i,j)
1972 rhs_vbar(i,j)=rhs_vbar(i,j)+cff1
1973# ifdef DIAGNOSTICS_UV
1974 diav2rhs(i,j,m2wrol)=cff1
1975# endif
1976# else
1977# ifdef DIAGNOSTICS_UV
1978 diav2rhs(i,j,m2wrol)=0.0_r8
1979# endif
1980# endif
1981# ifdef BOTTOM_STREAMING
1982 cff1=rvbst2d(i,j)
1983 rhs_vbar(i,j)=rhs_vbar(i,j)+cff1
1984# ifdef DIAGNOSTICS_UV
1985 diav2rhs(i,j,m2bstm)=cff1
1986# endif
1987# endif
1988# ifdef SURFACE_STREAMING
1989 cff1=rvsst2d(i,j)
1990 rhs_vbar(i,j)=rhs_vbar(i,j)+cff1
1991# ifdef DIAGNOSTICS_UV
1992 diav2rhs(i,j,m2sstm)=cff1
1993# endif
1994# endif
1995 END DO
1996 END IF
1997 END DO
1998# ifdef UV_ADV
1999# ifdef DIAGNOSTICS_UV
2000!
2001!---------------------------------------------------------------------------
2002! To obtain the full horizotal 'J' vortex force term:
2003! Compute term for diagnostics only. Subtract from hadv and add to vorf.
2004!---------------------------------------------------------------------------
2005!
2006 DO j=jstr,jend
2007 DO i=istru,iend
2008 cff=0.5_r8*(drhs(i-1,j)+drhs(i,j))
2009 dvsom(i,j)=0.25_r8*cff*om_u(i,j)* &
2010 & (vbar_stokes(i ,j )+ &
2011 & vbar_stokes(i ,j+1)+ &
2012 & vbar_stokes(i-1,j )+ &
2013 & vbar_stokes(i-1,j+1))
2014 END DO
2015 END DO
2016 DO j=jstr,jend+1
2017 DO i=istru,iend
2018 ufx(i,j)=0.5_r8*(ubar(i ,j-1,krhs)+ &
2019 ubar(i ,j ,krhs))
2020 END DO
2021 END DO
2022 DO j=jstr,jend
2023 DO i=istru,iend
2024 cff1=ufx(i,j+1)-ufx(i,j)
2025 cff=cff1*dvsom(i,j)
2026 diau2rhs(i,j,m2xadv)=diau2rhs(i,j,m2xadv)+cff
2027 diau2rhs(i,j,m2hadv)=diau2rhs(i,j,m2hadv)+cff
2028 diau2rhs(i,j,m2hjvf)=-cff
2029 END DO
2030 END DO
2031 DO j=jstrv,jend
2032 DO i=istr,iend
2033 cff=0.5_r8*(drhs(i,j)+drhs(i,j-1))
2034 duson(i,j)=cff*0.25_r8*on_v(i,j)* &
2035 & (ubar_stokes(i ,j )+ &
2036 & ubar_stokes(i+1,j )+ &
2037 & ubar_stokes(i ,j-1)+ &
2038 & ubar_stokes(i+1,j-1))
2039 END DO
2040 END DO
2041 DO j=jstrv,jend
2042 DO i=istr,iend+1
2043 vfe(i,j)=0.5_r8*(vbar(i-1,j ,krhs)+ &
2044 & vbar(i ,j ,krhs))
2045 END DO
2046 END DO
2047 DO i=istr,iend
2048 DO j=jstrv,jend
2049 cff2=vfe(i+1,j)-vfe(i,j)
2050 cff=cff2*duson(i,j)
2051 diav2rhs(i,j,m2yadv)=diav2rhs(i,j,m2yadv)+cff
2052 diav2rhs(i,j,m2hadv)=diav2rhs(i,j,m2hadv)+cff
2053 diav2rhs(i,j,m2hjvf)=-cff
2054 END DO
2055 END DO
2056# endif
2057!
2058!---------------------------------------------------------------------------
2059! Contribution of a term corresponding to product of
2060! Stokes and Eulerian Velocity Eqn. 26 and 27.
2061! This removes terms that were unneccessarily added in flux form.
2062!---------------------------------------------------------------------------
2063!
2064 DO j=jstr,jend
2065 DO i=istru,iend
2066 cff=0.5_r8*(drhs(i-1,j)+drhs(i,j))
2067 duson(i,j)=cff*on_u(i,j)*ubar_stokes(i,j)
2068 dvson(i,j)=0.25_r8*cff*on_u(i,j)* &
2069 & (vbar_stokes(i ,j )+ &
2070 & vbar_stokes(i ,j+1)+ &
2071 & vbar_stokes(i-1,j )+ &
2072 & vbar_stokes(i-1,j+1))
2073 END DO
2074 DO i=istru-1,iend
2075 ufx(i,j)=0.5_r8*(ubar(i ,j ,krhs)+ &
2076 ubar(i+1,j ,krhs))
2077 vfx(i,j)=0.5_r8*(vbar(i ,j ,krhs)+ &
2078 & vbar(i ,j+1,krhs))
2079 END DO
2080 END DO
2081 DO j=jstrv,jend
2082 DO i=istr,iend
2083 cff=0.5_r8*(drhs(i,j)+drhs(i,j-1))
2084 dusom(i,j)=cff*0.25_r8*om_v(i,j)* &
2085 & (ubar_stokes(i ,j )+ &
2086 & ubar_stokes(i+1,j )+ &
2087 & ubar_stokes(i ,j-1)+ &
2088 & ubar_stokes(i+1,j-1))
2089 dvsom(i,j)=cff*om_v(i,j)*vbar_stokes(i,j)
2090 END DO
2091 END DO
2092 DO j=jstrv-1,jend
2093 DO i=istr,iend
2094 cff=0.5_r8*(drhs(i,j)+drhs(i,j-1))
2095 ufe(i,j)=0.5_r8*(ubar(i+1,j ,krhs)+ &
2096 & ubar(i ,j ,krhs))
2097 vfe(i,j)=0.5_r8*(vbar(i ,j ,krhs)+ &
2098 & vbar(i ,j+1,krhs))
2099 END DO
2100 END DO
2101 DO j=jstr,jend
2102 DO i=istru,iend
2103 cff1=ufx(i,j)-ufx(i-1,j)
2104 cff2=vfx(i,j)-vfx(i-1,j)
2105 cff3=duson(i,j)*cff1
2106 cff4=dvson(i,j)*cff2
2107 rhs_ubar(i,j)=rhs_ubar(i,j)+cff3+cff4
2108! rustr2d(i,j)=rustr2d(i,j)-cff3-cff4
2109# ifdef DIAGNOSTICS_UV
2110 diau2rhs(i,j,m2xadv)=diau2rhs(i,j,m2xadv)+cff3
2111 diau2rhs(i,j,m2hadv)=diau2rhs(i,j,m2hadv)+cff3
2112 diau2rhs(i,j,m2hjvf)=diau2rhs(i,j,m2hjvf)+cff4
2113# endif
2114 END DO
2115 END DO
2116 DO i=istr,iend
2117 DO j=jstrv,jend
2118 cff1=ufe(i,j)-ufe(i,j-1)
2119 cff2=vfe(i,j)-vfe(i,j-1)
2120 cff3=dusom(i,j)*cff1
2121 cff4=dvsom(i,j)*cff2
2122 rhs_vbar(i,j)=rhs_vbar(i,j)+cff3+cff4
2123! rvstr2d(i,j)=rvstr2d(i,j,k)-cff3-cff4
2124# ifdef DIAGNOSTICS_UV
2125 diav2rhs(i,j,m2yadv)=diav2rhs(i,j,m2yadv)+cff4
2126 diav2rhs(i,j,m2hadv)=diav2rhs(i,j,m2hadv)+cff4
2127 diav2rhs(i,j,m2hjvf)=diav2rhs(i,j,m2hjvf)+cff3
2128# endif
2129 END DO
2130 END DO
2131# endif
2132#endif
2133#ifndef SOLVE3D
2134!
2135!-----------------------------------------------------------------------
2136! Add in bottom stress.
2137!-----------------------------------------------------------------------
2138!
2139 DO j=jstr,jend
2140 DO i=istru,iend
2141 fac=bustr(i,j)*om_u(i,j)*on_u(i,j)
2142 rhs_ubar(i,j)=rhs_ubar(i,j)-fac
2143# ifdef DIAGNOSTICS_UV
2144 diau2rhs(i,j,m2bstr)=-fac
2145# endif
2146 END DO
2147 END DO
2148 DO j=jstrv,jend
2149 DO i=istr,iend
2150 fac=bvstr(i,j)*om_v(i,j)*on_v(i,j)
2151 rhs_vbar(i,j)=rhs_vbar(i,j)-fac
2152# ifdef DIAGNOSTICS_UV
2153 diav2rhs(i,j,m2bstr)=-fac
2154# endif
2155 END DO
2156 END DO
2157#else
2158# ifdef DIAGNOSTICS_UV
2159!
2160! Initialize the stress term if no bottom friction is defined.
2161!
2162 DO j=jstr,jend
2163 DO i=istru,iend
2164 diau2rhs(i,j,m2bstr)=0.0_r8
2165 END DO
2166 END DO
2167 DO j=jstrv,jend
2168 DO i=istr,iend
2169 diav2rhs(i,j,m2bstr)=0.0_r8
2170 END DO
2171 END DO
2172# endif
2173#endif
2174!
2175!-----------------------------------------------------------------------
2176! Add in nudging of 2D momentum climatology.
2177!-----------------------------------------------------------------------
2178!
2179 IF (lnudgem2clm(ng)) THEN
2180 DO j=jstr,jend
2181 DO i=istru,iend
2182 cff=0.25_r8*(clima(ng)%M2nudgcof(i-1,j)+ &
2183 & clima(ng)%M2nudgcof(i ,j))* &
2184 & om_u(i,j)*on_u(i,j)
2185 rhs_ubar(i,j)=rhs_ubar(i,j)+ &
2186 & cff*(drhs(i-1,j)+drhs(i,j))* &
2187 & (clima(ng)%ubarclm(i,j)- &
2188 & ubar(i,j,krhs))
2189 END DO
2190 END DO
2191 DO j=jstrv,jend
2192 DO i=istr,iend
2193 cff=0.25_r8*(clima(ng)%M2nudgcof(i,j-1)+ &
2194 & clima(ng)%M2nudgcof(i,j ))* &
2195 & om_v(i,j)*on_v(i,j)
2196 rhs_vbar(i,j)=rhs_vbar(i,j)+ &
2197 & cff*(drhs(i,j-1)+drhs(i,j))* &
2198 & (clima(ng)%vbarclm(i,j)- &
2199 & vbar(i,j,krhs))
2200 END DO
2201 END DO
2202 END IF
2203
2204#ifdef SOLVE3D
2205# ifdef WET_DRY
2206 DO j=jstr,jend
2207 DO i=istru,iend
2208 cff5=abs(abs(umask_wet(i,j))-1.0_r8)
2209 cff6=0.5_r8+dsign(0.5_r8,rhs_ubar(i,j))*umask_wet(i,j)
2210 cff7=0.5_r8*umask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
2211 rhs_ubar(i,j)=rhs_ubar(i,j)*cff7
2212 END DO
2213 END DO
2214 DO j=jstrv,jend
2215 DO i=istr,iend
2216 cff5=abs(abs(vmask_wet(i,j))-1.0_r8)
2217 cff6=0.5_r8+dsign(0.5_r8,rhs_vbar(i,j))*vmask_wet(i,j)
2218 cff7=0.5_r8*vmask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
2219 rhs_vbar(i,j)=rhs_vbar(i,j)*cff7
2220 END DO
2221 END DO
2222# endif
2223!
2224!-----------------------------------------------------------------------
2225! Coupling between 2D and 3D equations.
2226!-----------------------------------------------------------------------
2227!
2228! Before the predictor step of the first barotropic time-step,
2229! arrays "rufrc" and "rvfrc" contain the vertical integrals of
2230! the 3D right-hand-side terms for momentum equations (including
2231! surface and bottom stresses, if so prescribed).
2232!
2233! Convert them into forcing terms by subtracting the fast time
2234! "rhs_ubar" and "rhs_vbar" from them; Also, immediately apply
2235! these forcing terms "rhs_ubar" and "rhs_vbar".
2236!
2237! From now on, these newly computed forcing terms will remain
2238! constant during the fast time stepping and will added to
2239! "rhs_ubar" and "rhs_vbar" during all subsequent time steps.
2240!
2241 IF (first_2d_step.and.predictor_2d_step(ng)) THEN
2242 IF (iic(ng).eq.ntfirst(ng)) THEN
2243 DO j=jstr,jend
2244 DO i=istru,iend
2245 rufrc(i,j)=rufrc(i,j)-rhs_ubar(i,j)
2246 rhs_ubar(i,j)=rhs_ubar(i,j)+rufrc(i,j)
2247 ru(i,j,0,nstp)=rufrc(i,j)
2248# ifdef DIAGNOSTICS_UV
2249 DO idiag=1,m2pgrd
2250 diarufrc(i,j,3,idiag)=diarufrc(i,j,3,idiag)- &
2251 & diau2rhs(i,j,idiag)
2252 diau2rhs(i,j,idiag)=diau2rhs(i,j,idiag)+ &
2253 & diarufrc(i,j,3,idiag)
2254 diarufrc(i,j,nstp,idiag)=diarufrc(i,j,3,idiag)
2255 END DO
2256 diau2rhs(i,j,m2sstr)=diarufrc(i,j,3,m2sstr)
2257 diarufrc(i,j,nstp,m2sstr)=diarufrc(i,j,3,m2sstr)
2258 diau2rhs(i,j,m2bstr)=diarufrc(i,j,3,m2bstr)
2259 diarufrc(i,j,nstp,m2bstr)=diarufrc(i,j,3,m2bstr)
2260# ifdef WEC_VF
2261 diau2rhs(i,j,m2zeta)=diau2rhs(i,j,m2zeta)+ &
2262 & diarufrc(i,j,3,m2pgrd)
2263# endif
2264# endif
2265 END DO
2266 END DO
2267 DO j=jstrv,jend
2268 DO i=istr,iend
2269 rvfrc(i,j)=rvfrc(i,j)-rhs_vbar(i,j)
2270 rhs_vbar(i,j)=rhs_vbar(i,j)+rvfrc(i,j)
2271 rv(i,j,0,nstp)=rvfrc(i,j)
2272# ifdef DIAGNOSTICS_UV
2273 DO idiag=1,m2pgrd
2274 diarvfrc(i,j,3,idiag)=diarvfrc(i,j,3,idiag)- &
2275 & diav2rhs(i,j,idiag)
2276 diav2rhs(i,j,idiag)=diav2rhs(i,j,idiag)+ &
2277 & diarvfrc(i,j,3,idiag)
2278 diarvfrc(i,j,nstp,idiag)=diarvfrc(i,j,3,idiag)
2279 END DO
2280 diav2rhs(i,j,m2sstr)=diarvfrc(i,j,3,m2sstr)
2281 diarvfrc(i,j,nstp,m2sstr)=diarvfrc(i,j,3,m2sstr)
2282 diav2rhs(i,j,m2bstr)=diarvfrc(i,j,3,m2bstr)
2283 diarvfrc(i,j,nstp,m2bstr)=diarvfrc(i,j,3,m2bstr)
2284# ifdef WEC_VF
2285 diav2rhs(i,j,m2zeta)=diav2rhs(i,j,m2zeta)+ &
2286 & diarvfrc(i,j,3,m2pgrd)
2287# endif
2288# endif
2289 END DO
2290 END DO
2291 ELSE IF (iic(ng).eq.(ntfirst(ng)+1)) THEN
2292 DO j=jstr,jend
2293 DO i=istru,iend
2294 rufrc(i,j)=rufrc(i,j)-rhs_ubar(i,j)
2295 rhs_ubar(i,j)=rhs_ubar(i,j)+ &
2296 & 1.5_r8*rufrc(i,j)-0.5_r8*ru(i,j,0,nnew)
2297 ru(i,j,0,nstp)=rufrc(i,j)
2298# ifdef DIAGNOSTICS_UV
2299 DO idiag=1,m2pgrd
2300 diarufrc(i,j,3,idiag)=diarufrc(i,j,3,idiag)- &
2301 & diau2rhs(i,j,idiag)
2302 diau2rhs(i,j,idiag)=diau2rhs(i,j,idiag)+ &
2303 & 1.5_r8*diarufrc(i,j,3,idiag)- &
2304 & 0.5_r8*diarufrc(i,j,nnew,idiag)
2305 diarufrc(i,j,nstp,idiag)=diarufrc(i,j,3,idiag)
2306 END DO
2307 diau2rhs(i,j,m2sstr)=1.5_r8*diarufrc(i,j,3,m2sstr)- &
2308 & 0.5_r8*diarufrc(i,j,nnew,m2sstr)
2309 diarufrc(i,j,nstp,m2sstr)=diarufrc(i,j,3,m2sstr)
2310 diau2rhs(i,j,m2bstr)=1.5_r8*diarufrc(i,j,3,m2bstr)- &
2311 & 0.5_r8*diarufrc(i,j,nnew,m2bstr)
2312 diarufrc(i,j,nstp,m2bstr)=diarufrc(i,j,3,m2bstr)
2313# ifdef WEC_VF
2314 diau2rhs(i,j,m2zeta)=diau2rhs(i,j,m2zeta)+ &
2315 & 1.5_r8*diarufrc(i,j,3,m2pgrd)- &
2316 & 0.5_r8*diarufrc(i,j,nnew,m2pgrd)
2317# endif
2318# endif
2319 END DO
2320 END DO
2321 DO j=jstrv,jend
2322 DO i=istr,iend
2323 rvfrc(i,j)=rvfrc(i,j)-rhs_vbar(i,j)
2324 rhs_vbar(i,j)=rhs_vbar(i,j)+ &
2325 & 1.5_r8*rvfrc(i,j)-0.5_r8*rv(i,j,0,nnew)
2326 rv(i,j,0,nstp)=rvfrc(i,j)
2327# ifdef DIAGNOSTICS_UV
2328 DO idiag=1,m2pgrd
2329 diarvfrc(i,j,3,idiag)=diarvfrc(i,j,3,idiag)- &
2330 & diav2rhs(i,j,idiag)
2331 diav2rhs(i,j,idiag)=diav2rhs(i,j,idiag)+ &
2332 & 1.5_r8*diarvfrc(i,j,3,idiag)- &
2333 & 0.5_r8*diarvfrc(i,j,nnew,idiag)
2334 diarvfrc(i,j,nstp,idiag)=diarvfrc(i,j,3,idiag)
2335 END DO
2336 diav2rhs(i,j,m2sstr)=1.5_r8*diarvfrc(i,j,3,m2sstr)- &
2337 & 0.5_r8*diarvfrc(i,j,nnew,m2sstr)
2338 diarvfrc(i,j,nstp,m2sstr)=diarvfrc(i,j,3,m2sstr)
2339 diav2rhs(i,j,m2bstr)=1.5_r8*diarvfrc(i,j,3,m2bstr)- &
2340 & 0.5_r8*diarvfrc(i,j,nnew,m2bstr)
2341 diarvfrc(i,j,nstp,m2bstr)=diarvfrc(i,j,3,m2bstr)
2342# ifdef WEC_VF
2343 diav2rhs(i,j,m2zeta)=diav2rhs(i,j,m2zeta)+ &
2344 & 1.5_r8*diarvfrc(i,j,3,m2pgrd)- &
2345 & 0.5_r8*diarvfrc(i,j,nnew,m2pgrd)
2346# endif
2347# endif
2348 END DO
2349 END DO
2350 ELSE
2351 cff1=23.0_r8/12.0_r8
2352 cff2=16.0_r8/12.0_r8
2353 cff3= 5.0_r8/12.0_r8
2354 DO j=jstr,jend
2355 DO i=istru,iend
2356 rufrc(i,j)=rufrc(i,j)-rhs_ubar(i,j)
2357 rhs_ubar(i,j)=rhs_ubar(i,j)+ &
2358 & cff1*rufrc(i,j)- &
2359 & cff2*ru(i,j,0,nnew)+ &
2360 & cff3*ru(i,j,0,nstp)
2361 ru(i,j,0,nstp)=rufrc(i,j)
2362# ifdef DIAGNOSTICS_UV
2363 DO idiag=1,m2pgrd
2364 diarufrc(i,j,3,idiag)=diarufrc(i,j,3,idiag)- &
2365 & diau2rhs(i,j,idiag)
2366 diau2rhs(i,j,idiag)=diau2rhs(i,j,idiag)+ &
2367 & cff1*diarufrc(i,j,3,idiag)- &
2368 & cff2*diarufrc(i,j,nnew,idiag)+ &
2369 & cff3*diarufrc(i,j,nstp,idiag)
2370 diarufrc(i,j,nstp,idiag)=diarufrc(i,j,3,idiag)
2371 END DO
2372 diau2rhs(i,j,m2sstr)=cff1*diarufrc(i,j,3,m2sstr)- &
2373 & cff2*diarufrc(i,j,nnew,m2sstr)+ &
2374 & cff3*diarufrc(i,j,nstp,m2sstr)
2375 diarufrc(i,j,nstp,m2sstr)=diarufrc(i,j,3,m2sstr)
2376 diau2rhs(i,j,m2bstr)=cff1*diarufrc(i,j,3,m2bstr)- &
2377 & cff2*diarufrc(i,j,nnew,m2bstr)+ &
2378 & cff3*diarufrc(i,j,nstp,m2bstr)
2379 diarufrc(i,j,nstp,m2bstr)=diarufrc(i,j,3,m2bstr)
2380# ifdef WEC_VF
2381 diau2rhs(i,j,m2zeta)=diau2rhs(i,j,m2zeta)+ &
2382 & cff1*diarufrc(i,j,3,m2pgrd)- &
2383 & cff2*diarufrc(i,j,nnew,m2pgrd)+ &
2384 & cff3*diarufrc(i,j,nstp,m2pgrd)
2385# endif
2386# endif
2387 END DO
2388 END DO
2389 DO j=jstrv,jend
2390 DO i=istr,iend
2391 rvfrc(i,j)=rvfrc(i,j)-rhs_vbar(i,j)
2392 rhs_vbar(i,j)=rhs_vbar(i,j)+ &
2393 & cff1*rvfrc(i,j)- &
2394 & cff2*rv(i,j,0,nnew)+ &
2395 & cff3*rv(i,j,0,nstp)
2396 rv(i,j,0,nstp)=rvfrc(i,j)
2397# ifdef DIAGNOSTICS_UV
2398 DO idiag=1,m2pgrd
2399 diarvfrc(i,j,3,idiag)=diarvfrc(i,j,3,idiag)- &
2400 & diav2rhs(i,j,idiag)
2401 diav2rhs(i,j,idiag)=diav2rhs(i,j,idiag)+ &
2402 & cff1*diarvfrc(i,j,3,idiag)- &
2403 & cff2*diarvfrc(i,j,nnew,idiag)+ &
2404 & cff3*diarvfrc(i,j,nstp,idiag)
2405 diarvfrc(i,j,nstp,idiag)=diarvfrc(i,j,3,idiag)
2406 END DO
2407 diav2rhs(i,j,m2sstr)=cff1*diarvfrc(i,j,3,m2sstr)- &
2408 & cff2*diarvfrc(i,j,nnew,m2sstr)+ &
2409 & cff3*diarvfrc(i,j,nstp,m2sstr)
2410 diarvfrc(i,j,nstp,m2sstr)=diarvfrc(i,j,3,m2sstr)
2411 diav2rhs(i,j,m2bstr)=cff1*diarvfrc(i,j,3,m2bstr)- &
2412 & cff2*diarvfrc(i,j,nnew,m2bstr)+ &
2413 & cff3*diarvfrc(i,j,nstp,m2bstr)
2414 diarvfrc(i,j,nstp,m2bstr)=diarvfrc(i,j,3,m2bstr)
2415# ifdef WEC_VF
2416 diav2rhs(i,j,m2zeta)=diav2rhs(i,j,m2zeta)+ &
2417 & cff1*diarvfrc(i,j,3,m2pgrd)- &
2418 & cff2*diarvfrc(i,j,nnew,m2pgrd)+ &
2419 & cff3*diarvfrc(i,j,nstp,m2pgrd)
2420# endif
2421# endif
2422 END DO
2423 END DO
2424 END IF
2425 ELSE
2426 DO j=jstr,jend
2427 DO i=istru,iend
2428 rhs_ubar(i,j)=rhs_ubar(i,j)+rufrc(i,j)
2429# ifdef DIAGNOSTICS_UV
2430 DO idiag=1,m2pgrd
2431 diau2rhs(i,j,idiag)=diau2rhs(i,j,idiag)+ &
2432 & diarufrc(i,j,3,idiag)
2433 END DO
2434 diau2rhs(i,j,m2sstr)=diarufrc(i,j,3,m2sstr)
2435 diau2rhs(i,j,m2bstr)=diarufrc(i,j,3,m2bstr)
2436# ifdef WEC_VF
2437 diau2rhs(i,j,m2zeta)=diau2rhs(i,j,m2zeta)+ &
2438 & diarufrc(i,j,3,m2pgrd)
2439# endif
2440# endif
2441 END DO
2442 END DO
2443 DO j=jstrv,jend
2444 DO i=istr,iend
2445 rhs_vbar(i,j)=rhs_vbar(i,j)+rvfrc(i,j)
2446# ifdef DIAGNOSTICS_UV
2447 DO idiag=1,m2pgrd
2448 diav2rhs(i,j,idiag)=diav2rhs(i,j,idiag)+ &
2449 & diarvfrc(i,j,3,idiag)
2450 END DO
2451 diav2rhs(i,j,m2sstr)=diarvfrc(i,j,3,m2sstr)
2452 diav2rhs(i,j,m2bstr)=diarvfrc(i,j,3,m2bstr)
2453# ifdef WEC_VF
2454 diav2rhs(i,j,m2zeta)=diav2rhs(i,j,m2zeta)+ &
2455 & diarvfrc(i,j,3,m2pgrd)
2456# endif
2457# endif
2458 END DO
2459 END DO
2460 END IF
2461#else
2462!
2463!-----------------------------------------------------------------------
2464! Add in surface momentum stress.
2465!-----------------------------------------------------------------------
2466!
2467 DO j=jstr,jend
2468 DO i=istru,iend
2469 fac=sustr(i,j)*om_u(i,j)*on_u(i,j)
2470 rhs_ubar(i,j)=rhs_ubar(i,j)+fac
2471# ifdef DIAGNOSTICS_UV
2472 diau2rhs(i,j,m2sstr)=fac
2473# endif
2474 END DO
2475 END DO
2476 DO j=jstrv,jend
2477 DO i=istr,iend
2478 fac=svstr(i,j)*om_v(i,j)*on_v(i,j)
2479 rhs_vbar(i,j)=rhs_vbar(i,j)+fac
2480# ifdef DIAGNOSTICS_UV
2481 diav2rhs(i,j,m2sstr)=fac
2482# endif
2483 END DO
2484 END DO
2485#endif
2486!
2487!=======================================================================
2488! Time step 2D momentum equations.
2489!=======================================================================
2490!
2491! Compute total water column depth.
2492!
2493 DO j=jstrv-1,jend
2494 DO i=istru-1,iend
2495 dstp(i,j)=zeta(i,j,kstp)+h(i,j)
2496 END DO
2497 END DO
2498!
2499! During the first time-step, the predictor step is Forward-Euler
2500! and the corrector step is Backward-Euler. Otherwise, the predictor
2501! step is Leap-frog and the corrector step is Adams-Moulton.
2502!
2503 IF (first_2d_step) THEN
2504 cff1=0.5_r8*dtfast(ng)
2505#ifdef WET_DRY
2506 cff2=1.0_r8/cff1
2507#endif
2508 DO j=jstr,jend
2509 DO i=istru,iend
2510 cff=(pm(i,j)+pm(i-1,j))*(pn(i,j)+pn(i-1,j))
2511 fac=1.0_r8/(dnew(i,j)+dnew(i-1,j))
2512 ubar(i,j,knew)=(ubar(i,j,kstp)* &
2513 & (dstp(i,j)+dstp(i-1,j))+ &
2514 & cff*cff1*rhs_ubar(i,j))*fac
2515#ifdef MASKING
2516 ubar(i,j,knew)=ubar(i,j,knew)*umask(i,j)
2517#endif
2518#ifdef WET_DRY
2519 cff5=abs(abs(umask_wet(i,j))-1.0_r8)
2520 cff6=0.5_r8+dsign(0.5_r8,ubar(i,j,knew))*umask_wet(i,j)
2521 cff7=0.5_r8*umask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
2522 ubar(i,j,knew)=ubar(i,j,knew)*cff7
2523 rhs_ubar(i,j)=rhs_ubar(i,j)*cff7
2524# ifdef SOLVE3D
2525 IF (predictor_2d_step(ng)) THEN
2526 rufrc(i,j)=rufrc(i,j)*cff7
2527 ru(i,j,0,nstp)=rufrc(i,j)
2528 END IF
2529# endif
2530#endif
2531#if defined NESTING && !defined SOLVE3D
2532 du_flux(i,j)=ubar(i,j,knew)* &
2533 & 0.5_r8*(dnew(i,j)+dnew(i-1,j))*on_u(i,j)
2534#endif
2535 END DO
2536 END DO
2537 DO j=jstrv,jend
2538 DO i=istr,iend
2539 cff=(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1))
2540 fac=1.0_r8/(dnew(i,j)+dnew(i,j-1))
2541 vbar(i,j,knew)=(vbar(i,j,kstp)* &
2542 & (dstp(i,j)+dstp(i,j-1))+ &
2543 & cff*cff1*rhs_vbar(i,j))*fac
2544#ifdef MASKING
2545 vbar(i,j,knew)=vbar(i,j,knew)*vmask(i,j)
2546#endif
2547#ifdef WET_DRY
2548 cff5=abs(abs(vmask_wet(i,j))-1.0_r8)
2549 cff6=0.5_r8+dsign(0.5_r8,vbar(i,j,knew))*vmask_wet(i,j)
2550 cff7=0.5_r8*vmask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
2551 vbar(i,j,knew)=vbar(i,j,knew)*cff7
2552 rhs_vbar(i,j)=rhs_vbar(i,j)*cff7
2553# ifdef SOLVE3D
2554 IF (predictor_2d_step(ng)) THEN
2555 rvfrc(i,j)=rvfrc(i,j)*cff7
2556 rv(i,j,0,nstp)=rvfrc(i,j)
2557 END IF
2558# endif
2559#endif
2560#if defined NESTING && !defined SOLVE3D
2561 dv_flux(i,j)=vbar(i,j,knew)* &
2562 & 0.5_r8*(dnew(i,j)+dnew(i,j-1))*om_v(i,j)
2563#endif
2564 END DO
2565 END DO
2566 ELSE IF (predictor_2d_step(ng)) THEN
2567 cff1=dtfast(ng)
2568#ifdef WET_DRY
2569 cff2=1.0_r8/cff1
2570#endif
2571 DO j=jstr,jend
2572 DO i=istru,iend
2573 cff=(pm(i,j)+pm(i-1,j))*(pn(i,j)+pn(i-1,j))
2574 fac=1.0_r8/(dnew(i,j)+dnew(i-1,j))
2575 ubar(i,j,knew)=(ubar(i,j,kstp)* &
2576 & (dstp(i,j)+dstp(i-1,j))+ &
2577 & cff*cff1*rhs_ubar(i,j))*fac
2578#ifdef MASKING
2579 ubar(i,j,knew)=ubar(i,j,knew)*umask(i,j)
2580#endif
2581#ifdef WET_DRY
2582 cff5=abs(abs(umask_wet(i,j))-1.0_r8)
2583 cff6=0.5_r8+dsign(0.5_r8,ubar(i,j,knew))*umask_wet(i,j)
2584 cff7=0.5_r8*umask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
2585 ubar(i,j,knew)=ubar(i,j,knew)*cff7
2586 rhs_ubar(i,j)=rhs_ubar(i,j)*cff7
2587#endif
2588#if defined NESTING && !defined SOLVE3D
2589 du_flux(i,j)=ubar(i,j,knew)* &
2590 & 0.5_r8*(dnew(i,j)+dnew(i-1,j))*on_u(i,j)
2591#endif
2592 END DO
2593 END DO
2594 DO j=jstrv,jend
2595 DO i=istr,iend
2596 cff=(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1))
2597 fac=1.0_r8/(dnew(i,j)+dnew(i,j-1))
2598 vbar(i,j,knew)=(vbar(i,j,kstp)* &
2599 & (dstp(i,j)+dstp(i,j-1))+ &
2600 & cff*cff1*rhs_vbar(i,j))*fac
2601#ifdef MASKING
2602 vbar(i,j,knew)=vbar(i,j,knew)*vmask(i,j)
2603#endif
2604#ifdef WET_DRY
2605 cff5=abs(abs(vmask_wet(i,j))-1.0_r8)
2606 cff6=0.5_r8+dsign(0.5_r8,vbar(i,j,knew))*vmask_wet(i,j)
2607 cff7=0.5_r8*vmask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
2608 vbar(i,j,knew)=vbar(i,j,knew)*cff7
2609 rhs_vbar(i,j)=rhs_vbar(i,j)*cff7
2610#endif
2611#if defined NESTING && !defined SOLVE3D
2612 dv_flux(i,j)=vbar(i,j,knew)* &
2613 & 0.5_r8*(dnew(i,j)+dnew(i,j-1))*om_v(i,j)
2614#endif
2615 END DO
2616 END DO
2617 ELSE IF (corrector_2d_step) THEN
2618 cff1=0.5_r8*dtfast(ng)*5.0_r8/12.0_r8
2619 cff2=0.5_r8*dtfast(ng)*8.0_r8/12.0_r8
2620 cff3=0.5_r8*dtfast(ng)*1.0_r8/12.0_r8
2621#ifdef WET_DRY
2622 cff4=1.0_r8/cff1
2623#endif
2624 DO j=jstr,jend
2625 DO i=istru,iend
2626 cff=(pm(i,j)+pm(i-1,j))*(pn(i,j)+pn(i-1,j))
2627 fac=1.0_r8/(dnew(i,j)+dnew(i-1,j))
2628 ubar(i,j,knew)=(ubar(i,j,kstp)* &
2629 & (dstp(i,j)+dstp(i-1,j))+ &
2630 & cff*(cff1*rhs_ubar(i,j)+ &
2631 & cff2*rubar(i,j,kstp)- &
2632 & cff3*rubar(i,j,ptsk)))*fac
2633#ifdef MASKING
2634 ubar(i,j,knew)=ubar(i,j,knew)*umask(i,j)
2635#endif
2636#ifdef WET_DRY
2637 cff5=abs(abs(umask_wet(i,j))-1.0_r8)
2638 cff6=0.5_r8+dsign(0.5_r8,ubar(i,j,knew))*umask_wet(i,j)
2639 cff7=0.5_r8*umask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
2640 ubar(i,j,knew)=ubar(i,j,knew)*cff7
2641 rhs_ubar(i,j)=rhs_ubar(i,j)*cff7
2642#endif
2643#if defined NESTING && !defined SOLVE3D
2644 du_flux(i,j)=ubar(i,j,knew)* &
2645 & 0.5_r8*(dnew(i,j)+dnew(i-1,j))*on_u(i,j)
2646#endif
2647 END DO
2648 END DO
2649 DO j=jstrv,jend
2650 DO i=istr,iend
2651 cff=(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1))
2652 fac=1.0_r8/(dnew(i,j)+dnew(i,j-1))
2653 vbar(i,j,knew)=(vbar(i,j,kstp)* &
2654 & (dstp(i,j)+dstp(i,j-1))+ &
2655 & cff*(cff1*rhs_vbar(i,j)+ &
2656 & cff2*rvbar(i,j,kstp)- &
2657 & cff3*rvbar(i,j,ptsk)))*fac
2658#ifdef MASKING
2659 vbar(i,j,knew)=vbar(i,j,knew)*vmask(i,j)
2660#endif
2661#ifdef WET_DRY
2662 cff5=abs(abs(vmask_wet(i,j))-1.0_r8)
2663 cff6=0.5_r8+dsign(0.5_r8,vbar(i,j,knew))*vmask_wet(i,j)
2664 cff7=0.5_r8*vmask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
2665 vbar(i,j,knew)=vbar(i,j,knew)*cff7
2666 rhs_vbar(i,j)=rhs_vbar(i,j)*cff7
2667#endif
2668#if defined NESTING && !defined SOLVE3D
2669 dv_flux(i,j)=vbar(i,j,knew)* &
2670 & 0.5_r8*(dnew(i,j)+dnew(i,j-1))*om_v(i,j)
2671#endif
2672 END DO
2673 END DO
2674 END IF
2675
2676#ifdef DIAGNOSTICS_UV
2677!
2678!-----------------------------------------------------------------------
2679! Time step 2D momentum diagnostic terms.
2680!-----------------------------------------------------------------------
2681
2682# ifdef MASKING
2683!
2684! Apply land/sea mask.
2685!
2686 DO idiag=1,ndm2d-1
2687 DO j=jstr,jend
2688 DO i=istru,iend
2689 diau2rhs(i,j,idiag)=diau2rhs(i,j,idiag)*umask(i,j)
2690 END DO
2691 END DO
2692 DO j=jstrv,jend
2693 DO i=istr,iend
2694 diav2rhs(i,j,idiag)=diav2rhs(i,j,idiag)*vmask(i,j)
2695 END DO
2696 END DO
2697 END DO
2698# endif
2699# ifdef SOLVE3D
2700!
2701! The arrays "DiaU2rhs" and "DiaV2rhs" contain the contributions of
2702! each of the 2D right-hand-side terms for the momentum equations.
2703!
2704! These values are integrated, time-stepped and converted to mass flux
2705! units (m3 s-1) for coupling with the 3D diagnostic terms.
2706!
2707 fac=weight(1,iif(ng),ng)
2708 IF (first_2d_step.and.corrector_2d_step) THEN
2709 cff1=0.5_r8*dtfast(ng)
2710 DO idiag=1,ndm2d-1
2711 DO j=jstr,jend
2712 DO i=istru,iend
2713 diau2int(i,j,idiag)=cff1*diau2rhs(i,j,idiag)
2714 diau2wrk(i,j,idiag)=diau2int(i,j,idiag)* &
2715 & (pm(i-1,j)+pm(i,j))*fac
2716 END DO
2717 END DO
2718 DO j=jstrv,jend
2719 DO i=istr,iend
2720 diav2int(i,j,idiag)=cff1*diav2rhs(i,j,idiag)
2721 diav2wrk(i,j,idiag)=diav2int(i,j,idiag)* &
2722 & (pn(i,j)+pn(i,j-1))*fac
2723 END DO
2724 END DO
2725 END DO
2726 ELSE IF (corrector_2d_step) THEN
2727 cff1=0.5_r8*dtfast(ng)*5.0_r8/12.0_r8
2728 cff2=0.5_r8*dtfast(ng)*8.0_r8/12.0_r8
2729 cff3=0.5_r8*dtfast(ng)*1.0_r8/12.0_r8
2730 DO idiag=1,ndm2d-1
2731 DO j=jstr,jend
2732 DO i=istru,iend
2733 diau2int(i,j,idiag)=diau2int(i,j,idiag)+ &
2734 & (cff1*diau2rhs(i,j,idiag)+ &
2735 & cff2*diarubar(i,j,kstp,idiag)- &
2736 & cff3*diarubar(i,j,ptsk,idiag))
2737 diau2wrk(i,j,idiag)=diau2wrk(i,j,idiag)+ &
2738 & diau2int(i,j,idiag)* &
2739 & (pm(i-1,j)+pm(i,j))*fac
2740 END DO
2741 END DO
2742 DO j=jstrv,jend
2743 DO i=istr,iend
2744 diav2int(i,j,idiag)=diav2int(i,j,idiag)+ &
2745 & (cff1*diav2rhs(i,j,idiag)+ &
2746 & cff2*diarvbar(i,j,kstp,idiag)- &
2747 & cff3*diarvbar(i,j,ptsk,idiag))
2748 diav2wrk(i,j,idiag)=diav2wrk(i,j,idiag)+ &
2749 & diav2int(i,j,idiag)* &
2750 & (pn(i,j)+pn(i,j-1))*fac
2751 END DO
2752 END DO
2753 END DO
2754 END IF
2755# else
2756!
2757! Time-step the diagnostic terms.
2758!
2759 IF (first_2d_step.and.corrector_2d_step) THEN
2760 cff1=0.5_r8*dtfast(ng)
2761 DO idiag=1,ndm2d-1
2762 DO j=jstr,jend
2763 DO i=istru,iend
2764 cff=(pm(i,j)+pm(i-1,j))*(pn(i,j)+pn(i-1,j))
2765 fac=1.0_r8/(dnew(i,j)+dnew(i-1,j))
2766 diau2wrk(i,j,idiag)=cff*cff1*diau2rhs(i,j,idiag)*fac
2767 END DO
2768 END DO
2769 DO j=jstrv,jend
2770 DO i=istr,iend
2771 cff=(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1))
2772 fac=1.0_r8/(dnew(i,j)+dnew(i,j-1))
2773 diav2wrk(i,j,idiag)=cff*cff1*diav2rhs(i,j,idiag)*fac
2774 END DO
2775 END DO
2776 END DO
2777 DO j=jstr,jend
2778 DO i=istru,iend
2779 fac=1.0_r8/(dnew(i,j)+dnew(i-1,j))
2780 diau2wrk(i,j,m2rate)=ubar(i,j,knew)-ubar(i,j,kstp)* &
2781 & (dstp(i,j)+dstp(i-1,j))*fac
2782 END DO
2783 END DO
2784 DO j=jstrv,jend
2785 DO i=istr,iend
2786 fac=1.0_r8/(dnew(i,j)+dnew(i,j-1))
2787 diav2wrk(i,j,m2rate)=vbar(i,j,knew)-vbar(i,j,kstp)* &
2788 & (dstp(i,j)+dstp(i,j-1))*fac
2789 END DO
2790 END DO
2791 ELSE IF (corrector_2d_step) THEN
2792 cff1=0.5_r8*dtfast(ng)*5.0_r8/12.0_r8
2793 cff2=0.5_r8*dtfast(ng)*8.0_r8/12.0_r8
2794 cff3=0.5_r8*dtfast(ng)*1.0_r8/12.0_r8
2795 DO idiag=1,ndm2d-1
2796 DO j=jstr,jend
2797 DO i=istru,iend
2798 cff=(pm(i,j)+pm(i-1,j))*(pn(i,j)+pn(i-1,j))
2799 fac=1.0_r8/(dnew(i,j)+dnew(i-1,j))
2800 diau2wrk(i,j,idiag)=cff*(cff1*diau2rhs(i,j,idiag)+ &
2801 & cff2*diarubar(i,j,kstp,idiag)- &
2802 & cff3*diarubar(i,j,ptsk,idiag))* &
2803 & fac
2804 END DO
2805 END DO
2806 DO j=jstrv,jend
2807 DO i=istr,iend
2808 cff=(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1))
2809 fac=1.0_r8/(dnew(i,j)+dnew(i,j-1))
2810 diav2wrk(i,j,idiag)=cff*(cff1*diav2rhs(i,j,idiag)+ &
2811 & cff2*diarvbar(i,j,kstp,idiag)- &
2812 & cff3*diarvbar(i,j,ptsk,idiag))* &
2813 & fac
2814 END DO
2815 END DO
2816 END DO
2817 DO j=jstr,jend
2818 DO i=istru,iend
2819 fac=1.0_r8/(dnew(i,j)+dnew(i-1,j))
2820 diau2wrk(i,j,m2rate)=ubar(i,j,knew)- &
2821 & ubar(i,j,kstp)* &
2822 & (dstp(i,j)+dstp(i-1,j))*fac
2823 END DO
2824 END DO
2825 DO j=jstrv,jend
2826 DO i=istr,iend
2827 fac=1.0_r8/(dnew(i,j)+dnew(i,j-1))
2828 diav2wrk(i,j,m2rate)=vbar(i,j,knew)- &
2829 & vbar(i,j,kstp)* &
2830 & (dstp(i,j)+dstp(i,j-1))*fac
2831 END DO
2832 END DO
2833 END IF
2834# endif
2835#endif
2836!
2837! If predictor step, load right-side-term into shared arrays for
2838! future use during the subsequent corrector step.
2839!
2840 IF (predictor_2d_step(ng)) THEN
2841 DO j=jstr,jend
2842 DO i=istru,iend
2843 rubar(i,j,krhs)=rhs_ubar(i,j)
2844 END DO
2845 END DO
2846 DO j=jstrv,jend
2847 DO i=istr,iend
2848 rvbar(i,j,krhs)=rhs_vbar(i,j)
2849 END DO
2850 END DO
2851#ifdef DIAGNOSTICS_UV
2852 DO idiag=1,ndm2d-1
2853 DO j=jstr,jend
2854 DO i=istru,iend
2855 diarubar(i,j,krhs,idiag)=diau2rhs(i,j,idiag)
2856 END DO
2857 END DO
2858 DO j=jstrv,jend
2859 DO i=istr,iend
2860 diarvbar(i,j,krhs,idiag)=diav2rhs(i,j,idiag)
2861 END DO
2862 END DO
2863 END DO
2864#endif
2865 END IF
2866!
2867!-----------------------------------------------------------------------
2868! Apply lateral boundary conditions.
2869!-----------------------------------------------------------------------
2870!
2871 CALL u2dbc_tile (ng, tile, &
2872 & lbi, ubi, lbj, ubj, &
2873 & imins, imaxs, jmins, jmaxs, &
2874 & krhs, kstp, knew, &
2875 & ubar, vbar, zeta)
2876 CALL v2dbc_tile (ng, tile, &
2877 & lbi, ubi, lbj, ubj, &
2878 & imins, imaxs, jmins, jmaxs, &
2879 & krhs, kstp, knew, &
2880 & ubar, vbar, zeta)
2881!
2882! Compute integral mass flux across open boundaries and adjust
2883! for volume conservation.
2884!
2885 IF (any(volcons(:,ng))) THEN
2886 CALL obc_flux_tile (ng, tile, &
2887 & lbi, ubi, lbj, ubj, &
2888 & imins, imaxs, jmins, jmaxs, &
2889 & knew, &
2890#ifdef MASKING
2891 & umask, vmask, &
2892#endif
2893 & h, om_v, on_u, &
2894 & ubar, vbar, zeta)
2895 END IF
2896
2897#if defined NESTING && !defined SOLVE3D
2898!
2899!-----------------------------------------------------------------------
2900! Set barotropic fluxes along physical boundaries.
2901!-----------------------------------------------------------------------
2902!
2903 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
2904 IF (domain(ng)%Western_Edge(tile)) THEN
2905 DO j=jstr-1,jendr
2906 dnew(istr-1,j)=h(istr-1,j)+zeta_new(istr-1,j)
2907 END DO
2908 END IF
2909 END IF
2910 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
2911 IF (domain(ng)%Eastern_Edge(tile)) THEN
2912 DO j=jstr-1,jendr
2913 dnew(iend+1,j)=h(iend+1,j)+zeta_new(iend+1,j)
2914 END DO
2915 END IF
2916 END IF
2917 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
2918 IF (domain(ng)%Southern_Edge(tile)) THEN
2919 DO i=istr-1,iendr
2920 dnew(i,jstr-1)=h(i,jstr-1)+zeta_new(i,jstr-1)
2921 END DO
2922 END IF
2923 END IF
2924 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
2925 IF (domain(ng)%Northern_Edge(tile)) THEN
2926 DO i=istr-1,iendr
2927 dnew(i,jend+1)=h(i,jend+1)+zeta_new(i,jend+1)
2928 END DO
2929 END IF
2930 END IF
2931!
2932 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
2933 IF (domain(ng)%Western_Edge(tile)) THEN
2934 DO j=jstrr,jendr
2935 du_flux(istru-1,j)=ubar(istru-1,j,knew)*on_u(istru-1,j)* &
2936 & 0.5_r8*(dnew(istru-1,j)+dnew(istru-2,j))
2937 END DO
2938 DO j=jstrv,jend
2939 dv_flux(istr-1,j)=vbar(istr-1,j,knew)*om_v(istr-1,j)* &
2940 & 0.5_r8*(dnew(istr-1,j)+dnew(istr-1,j-1))
2941 END DO
2942 END IF
2943 END IF
2944 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
2945 IF (domain(ng)%Eastern_Edge(tile)) THEN
2946 DO j=jstrr,jendr
2947 du_flux(iend+1,j)=ubar(iend+1,j,knew)*on_u(iend+1,j)* &
2948 & 0.5_r8*(dnew(iend+1,j)+dnew(iend,j))
2949 END DO
2950 DO j=jstrv,jend
2951 dv_flux(iend+1,j)=vbar(iend+1,j,knew)*om_v(iend+1,j)* &
2952 & 0.5_r8*(dnew(iend+1,j)+dnew(iend+1,j-1))
2953 END DO
2954 END IF
2955 END IF
2956 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
2957 IF (domain(ng)%Southern_Edge(tile)) THEN
2958 DO i=istru,iend
2959 du_flux(i,jstr-1)=ubar(i,jstr-1,knew)*on_u(i,jstr-1)* &
2960 & 0.5_r8*(dnew(i,jstr-1)+dnew(i-1,jstr-1))
2961 END DO
2962 DO i=istrr,iendr
2963 dv_flux(i,jstrv-1)=vbar(i,jstrv-1,knew)*om_v(i,jstrv-1)* &
2964 & 0.5_r8*(dnew(i,jstrv-1)+dnew(i,jstrv-2))
2965 END DO
2966 END IF
2967 END IF
2968 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
2969 IF (domain(ng)%Northern_Edge(tile)) THEN
2970 DO i=istru,iend
2971 du_flux(i,jend+1)=ubar(i,jend+1,knew)*on_u(i,jend+1)* &
2972 & 0.5_r8*(dnew(i,jend+1)+dnew(i-1,jend+1))
2973 END DO
2974 DO i=istrr,iendr
2975 dv_flux(i,jend+1)=vbar(i,jend+1,knew)*om_v(i,jend+1)* &
2976 & 0.5_r8*(dnew(i,jend+1)+dnew(i,jend))
2977 END DO
2978 END IF
2979 END IF
2980#endif
2981!
2982!-----------------------------------------------------------------------
2983! Apply momentum transport point sources (like river runoff), if any.
2984!
2985! Dsrc(is) = 0, flow across grid cell u-face (positive or negative)
2986! Dsrc(is) = 1, flow across grid cell v-face (positive or negative)
2987!-----------------------------------------------------------------------
2988!
2989 IF (luvsrc(ng)) THEN
2990 DO is=1,nsrc(ng)
2991 i=sources(ng)%Isrc(is)
2992 j=sources(ng)%Jsrc(is)
2993 IF (((istrr.le.i).and.(i.le.iendr)).and. &
2994 & ((jstrr.le.j).and.(j.le.jendr))) THEN
2995 IF (int(sources(ng)%Dsrc(is)).eq.0) THEN
2996 cff=1.0_r8/(on_u(i,j)* &
2997 & 0.5_r8*(zeta(i-1,j,knew)+h(i-1,j)+ &
2998 & zeta(i ,j,knew)+h(i ,j)))
2999 ubar(i,j,knew)=sources(ng)%Qbar(is)*cff
3000#if defined NESTING && !defined SOLVE3D
3001 du_flux(i,j)=sources(ng)%Qbar(is)
3002#endif
3003 ELSE IF (int(sources(ng)%Dsrc(is)).eq.1) THEN
3004 cff=1.0_r8/(om_v(i,j)* &
3005 & 0.5_r8*(zeta(i,j-1,knew)+h(i,j-1)+ &
3006 & zeta(i,j ,knew)+h(i,j )))
3007 vbar(i,j,knew)=sources(ng)%Qbar(is)*cff
3008#if defined NESTING && !defined SOLVE3D
3009 dv_flux(i,j)=sources(ng)%Qbar(is)
3010#endif
3011 END IF
3012 END IF
3013 END DO
3014 END IF
3015!
3016!-----------------------------------------------------------------------
3017! Exchange boundary information.
3018!-----------------------------------------------------------------------
3019!
3020 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
3021 CALL exchange_u2d_tile (ng, tile, &
3022 & lbi, ubi, lbj, ubj, &
3023 & ubar(:,:,knew))
3024 CALL exchange_v2d_tile (ng, tile, &
3025 & lbi, ubi, lbj, ubj, &
3026 & vbar(:,:,knew))
3027
3028#if defined NESTING && !defined SOLVE3D
3029 CALL exchange_u2d_tile (ng, tile, &
3030 & lbi, ubi, lbj, ubj, &
3031 & du_flux)
3032 CALL exchange_v2d_tile (ng, tile, &
3033 & lbi, ubi, lbj, ubj, &
3034 & dv_flux)
3035#endif
3036 END IF
3037
3038#ifdef DISTRIBUTE
3039!
3040# if defined NESTING && !defined SOLVE3D
3041 CALL mp_exchange2d (ng, tile, inlm, 4, &
3042# else
3043 CALL mp_exchange2d (ng, tile, inlm, 2, &
3044# endif
3045 & lbi, ubi, lbj, ubj, &
3046 & nghostpoints, &
3047 & ewperiodic(ng), nsperiodic(ng), &
3048# if defined NESTING && !defined SOLVE3D
3049 & du_flux, dv_flux, &
3050# endif
3051 & ubar(:,:,knew), &
3052 & vbar(:,:,knew))
3053#endif
3054!
3055 RETURN
subroutine exchange_r2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
subroutine exchange_u2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
subroutine exchange_v2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
type(t_clima), dimension(:), allocatable clima
Definition mod_clima.F:153
integer isvbar
integer isubar
integer, parameter inlm
Definition mod_param.F:662
integer nghostpoints
Definition mod_param.F:710
type(t_lbc), dimension(:,:,:), allocatable lbc
Definition mod_param.F:375
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer ndm2d
Definition mod_param.F:578
logical, dimension(:), allocatable luvsrc
integer m2fcor
logical, dimension(:), allocatable lnudgem2clm
integer, dimension(:), allocatable iic
real(dp), dimension(:), allocatable dt
real(r8), dimension(:), allocatable dcrit
integer m2pgrd
logical, dimension(:), allocatable ewperiodic
integer, parameter iwest
logical, dimension(:), allocatable nsperiodic
integer m2xadv
integer m2yadv
logical, dimension(:), allocatable predictor_2d_step
real(dp), dimension(:,:,:), allocatable weight
logical, dimension(:,:), allocatable volcons
integer, dimension(:), allocatable nfast
integer, dimension(:), allocatable ndtfast
logical, dimension(:), allocatable lwsrc
real(r8), dimension(:), allocatable gamma2
integer m2hvis
logical, dimension(:,:), allocatable compositegrid
integer, parameter isouth
integer m2rate
integer m2yvis
real(dp), dimension(:), allocatable dtfast
integer, dimension(:), allocatable ntfirst
integer m2sstr
integer, parameter ieast
real(dp) g
integer m2hadv
real(dp) rho0
integer, parameter inorth
integer m2xvis
integer m2bstr
integer, dimension(:), allocatable iif
type(t_sources), dimension(:), allocatable sources
Definition mod_sources.F:90
integer, dimension(:), allocatable nsrc
Definition mod_sources.F:97
subroutine mp_exchange2d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, nghost, ew_periodic, ns_periodic, a, b, c, d)
subroutine, public set_duv_bc_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, kinp, umask, vmask, om_v, on_u, ubar, vbar, drhs, duon, dvom)
subroutine, public obc_flux_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, kinp, umask, vmask, h, om_v, on_u, ubar, vbar, zeta)
Definition obc_volcons.F:69
subroutine, public u2dbc_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, krhs, kstp, kout, ubar, vbar, zeta)
Definition u2dbc_im.F:56
subroutine, public v2dbc_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, krhs, kstp, kout, ubar, vbar, zeta)
Definition v2dbc_im.F:57
subroutine wetdry_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, pmask, rmask, umask, vmask, h, zeta, du_avg1, dv_avg1, rmask_wet_avg, pmask_wet, pmask_full, rmask_wet, rmask_full, umask_wet, umask_full, vmask_wet, vmask_full)
Definition wetdry.F:108
subroutine, public zetabc_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, krhs, kstp, kout, zeta)
Definition zetabc.F:65

References mod_clima::clima, mod_scalars::compositegrid, mod_scalars::dcrit, mod_param::domain, mod_scalars::dt, mod_scalars::dtfast, mod_scalars::ewperiodic, exchange_2d_mod::exchange_r2d_tile(), exchange_2d_mod::exchange_u2d_tile(), exchange_2d_mod::exchange_v2d_tile(), mod_scalars::g, mod_scalars::gamma2, mod_scalars::ieast, mod_scalars::iic, mod_scalars::iif, mod_param::inlm, mod_scalars::inorth, mod_scalars::isouth, mod_ncparam::isubar, mod_ncparam::isvbar, mod_scalars::iwest, mod_param::lbc, mod_scalars::lnudgem2clm, mod_scalars::luvsrc, mod_scalars::lwsrc, mod_scalars::m2bstr, mod_scalars::m2fcor, mod_scalars::m2hadv, mod_scalars::m2hvis, mod_scalars::m2pgrd, mod_scalars::m2rate, mod_scalars::m2sstr, mod_scalars::m2xadv, mod_scalars::m2xvis, mod_scalars::m2yadv, mod_scalars::m2yvis, mp_exchange_mod::mp_exchange2d(), mod_scalars::ndtfast, mod_scalars::nfast, mod_param::nghostpoints, mod_scalars::nsperiodic, mod_sources::nsrc, mod_scalars::ntfirst, obc_volcons_mod::obc_flux_tile(), mod_scalars::predictor_2d_step, mod_scalars::rho0, obc_volcons_mod::set_duv_bc_tile(), mod_sources::sources, u2dbc_mod::u2dbc_tile(), v2dbc_mod::v2dbc_tile(), mod_scalars::volcons, mod_scalars::weight, wetdry_mod::wetdry_tile(), and zetabc_mod::zetabc_tile().

Here is the call graph for this function:

◆ step2d_tile() [3/3]

subroutine step2d_mod::step2d_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) ubk,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) krhs,
integer, intent(in) kstp,
integer, intent(in) knew,
integer, intent(in) nstp,
integer, intent(in) nnew,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pmask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rmask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) umask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) vmask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) pmask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) pmask_full,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) rmask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) rmask_full,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) umask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) umask_full,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) vmask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) vmask_full,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) rmask_wet_avg,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) fomn,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) h,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_u,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_v,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_u,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_v,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pm,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pn,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) dndx,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) dmde,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rdrag,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rdrag2,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pmon_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pnom_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pmon_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pnom_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) visc2_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) visc2_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) visc4_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) visc4_r,
real(r8), dimension(lbi:ubi,lbj:ubj,1:3), intent(in) bed_thick,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) eq_tide,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) sustr,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) svstr )
private

Definition at line 170 of file step2d_FB.h.

219# ifdef ATM_PRESS
220 & pair, &
221# endif
222#else
223# ifdef VAR_RHO_2D
224 & rhoa, rhos, &
225# endif
226 & du_avg1, du_avg2, &
227 & dv_avg1, dv_avg2, &
228 & zt_avg1, &
229 & rufrc, rvfrc, &
230 & rufrc_bak, rvfrc_bak, &
231#endif
232#ifdef DIAGNOSTICS_UV
233 & diau2wrk, diav2wrk, &
234 & diarubar, diarvbar, &
235# ifdef SOLVE3D
236 & diau2int, diav2int, &
237 & diarufrc, diarvfrc, &
238# endif
239#endif
240#if defined NESTING && !defined SOLVE3D
241 & du_flux, dv_flux, &
242#endif
243 & ubar, vbar, zeta)
244!***********************************************************************
245!
246! Imported variable declarations.
247!
248 integer, intent(in ) :: ng, tile
249 integer, intent(in ) :: LBi, UBi, LBj, UBj, UBk
250 integer, intent(in ) :: IminS, ImaxS, JminS, JmaxS
251 integer, intent(in ) :: krhs, kstp, knew
252#ifdef SOLVE3D
253 integer, intent(in ) :: nstp, nnew
254#endif
255!
256#ifdef ASSUMED_SHAPE
257# ifdef MASKING
258 real(r8), intent(in ) :: pmask(LBi:,LBj:)
259 real(r8), intent(in ) :: rmask(LBi:,LBj:)
260 real(r8), intent(in ) :: umask(LBi:,LBj:)
261 real(r8), intent(in ) :: vmask(LBi:,LBj:)
262# endif
263# if (defined UV_COR && !defined SOLVE3D) || defined STEP2D_CORIOLIS
264 real(r8), intent(in ) :: fomn(LBi:,LBj:)
265# endif
266# if defined SEDIMENT && defined SED_MORPH
267 real(r8), intent(inout) :: h(LBi:,LBj:)
268# else
269 real(r8), intent(in ) :: h(LBi:,LBj:)
270# endif
271 real(r8), intent(in ) :: om_u(LBi:,LBj:)
272 real(r8), intent(in ) :: om_v(LBi:,LBj:)
273 real(r8), intent(in ) :: on_u(LBi:,LBj:)
274 real(r8), intent(in ) :: on_v(LBi:,LBj:)
275 real(r8), intent(in ) :: pm(LBi:,LBj:)
276 real(r8), intent(in ) :: pn(LBi:,LBj:)
277# if defined CURVGRID && defined UV_ADV && !defined SOLVE3D
278 real(r8), intent(in ) :: dndx(LBi:,LBj:)
279 real(r8), intent(in ) :: dmde(LBi:,LBj:)
280# endif
281 real(r8), intent(in ) :: rdrag(LBi:,LBj:)
282# if defined UV_QDRAG && !defined SOLVE3D
283 real(r8), intent(in ) :: rdrag2(LBi:,LBj:)
284# endif
285# if (defined UV_VIS2 || defined UV_VIS4) && !defined SOLVE3D
286 real(r8), intent(in ) :: pmon_r(LBi:,LBj:)
287 real(r8), intent(in ) :: pnom_r(LBi:,LBj:)
288 real(r8), intent(in ) :: pmon_p(LBi:,LBj:)
289 real(r8), intent(in ) :: pnom_p(LBi:,LBj:)
290 real(r8), intent(in ) :: om_r(LBi:,LBj:)
291 real(r8), intent(in ) :: on_r(LBi:,LBj:)
292 real(r8), intent(in ) :: om_p(LBi:,LBj:)
293 real(r8), intent(in ) :: on_p(LBi:,LBj:)
294# ifdef UV_VIS2
295 real(r8), intent(in ) :: visc2_p(LBi:,LBj:)
296 real(r8), intent(in ) :: visc2_r(LBi:,LBj:)
297# endif
298# ifdef UV_VIS4
299 real(r8), intent(in ) :: visc4_p(LBi:,LBj:)
300 real(r8), intent(in ) :: visc4_r(LBi:,LBj:)
301# endif
302# endif
303# if defined SEDIMENT && defined SED_MORPH
304 real(r8), intent(in ) :: bed_thick(LBi:,LBj:,:)
305# endif
306# if defined TIDE_GENERATING_FORCES && !defined SOLVE3D
307 real(r8), intent(in ) :: eq_tide(LBi:,LBj:)
308# endif
309# ifndef SOLVE3D
310 real(r8), intent(in ) :: sustr(LBi:,LBj:)
311 real(r8), intent(in ) :: svstr(LBi:,LBj:)
312# ifdef ATM_PRESS
313 real(r8), intent(in ) :: Pair(LBi:,LBj:)
314# endif
315# else
316# ifdef VAR_RHO_2D
317 real(r8), intent(in ) :: rhoA(LBi:,LBj:)
318 real(r8), intent(in ) :: rhoS(LBi:,LBj:)
319# endif
320 real(r8), intent(inout) :: DU_avg1(LBi:,LBj:)
321 real(r8), intent(inout) :: DU_avg2(LBi:,LBj:)
322 real(r8), intent(inout) :: DV_avg1(LBi:,LBj:)
323 real(r8), intent(inout) :: DV_avg2(LBi:,LBj:)
324 real(r8), intent(inout) :: Zt_avg1(LBi:,LBj:)
325 real(r8), intent(inout) :: rufrc(LBi:,LBj:)
326 real(r8), intent(inout) :: rvfrc(LBi:,LBj:)
327 real(r8), intent(inout) :: rufrc_bak(LBi:,LBj:,:)
328 real(r8), intent(inout) :: rvfrc_bak(LBi:,LBj:,:)
329# endif
330# ifdef WET_DRY
331 real(r8), intent(inout) :: pmask_full(LBi:,LBj:)
332 real(r8), intent(inout) :: rmask_full(LBi:,LBj:)
333 real(r8), intent(inout) :: umask_full(LBi:,LBj:)
334 real(r8), intent(inout) :: vmask_full(LBi:,LBj:)
335
336 real(r8), intent(inout) :: pmask_wet(LBi:,LBj:)
337 real(r8), intent(inout) :: rmask_wet(LBi:,LBj:)
338 real(r8), intent(inout) :: umask_wet(LBi:,LBj:)
339 real(r8), intent(inout) :: vmask_wet(LBi:,LBj:)
340# ifdef SOLVE3D
341 real(r8), intent(inout) :: rmask_wet_avg(LBi:,LBj:)
342# endif
343# endif
344 real(r8), intent(inout) :: ubar(LBi:,LBj:,:)
345 real(r8), intent(inout) :: vbar(LBi:,LBj:,:)
346 real(r8), intent(inout) :: zeta(LBi:,LBj:,:)
347# if defined NESTING && !defined SOLVE3D
348 real(r8), intent(out ) :: DU_flux(LBi:,LBj:)
349 real(r8), intent(out ) :: DV_flux(LBi:,LBj:)
350# endif
351
352#else
353
354# ifdef MASKING
355 real(r8), intent(in ) :: pmask(LBi:UBi,LBj:UBj)
356 real(r8), intent(in ) :: rmask(LBi:UBi,LBj:UBj)
357 real(r8), intent(in ) :: umask(LBi:UBi,LBj:UBj)
358 real(r8), intent(in ) :: vmask(LBi:UBi,LBj:UBj)
359# endif
360# if (defined UV_COR && !defined SOLVE3D) || defined STEP2D_CORIOLIS
361 real(r8), intent(in ) :: fomn(LBi:UBi,LBj:UBj)
362# endif
363# if defined SEDIMENT && defined SED_MORPH
364 real(r8), intent(inout) :: h(LBi:UBi,LBj:UBj)
365# else
366 real(r8), intent(in ) :: h(LBi:UBi,LBj:UBj)
367# endif
368 real(r8), intent(in ) :: om_u(LBi:UBi,LBj:UBj)
369 real(r8), intent(in ) :: om_v(LBi:UBi,LBj:UBj)
370 real(r8), intent(in ) :: on_u(LBi:UBi,LBj:UBj)
371 real(r8), intent(in ) :: on_v(LBi:UBi,LBj:UBj)
372 real(r8), intent(in ) :: pm(LBi:UBi,LBj:UBj)
373 real(r8), intent(in ) :: pn(LBi:UBi,LBj:UBj)
374# if defined CURVGRID && defined UV_ADV && !defined SOLVE3D
375 real(r8), intent(in ) :: dndx(LBi:UBi,LBj:UBj)
376 real(r8), intent(in ) :: dmde(LBi:UBi,LBj:UBj)
377# endif
378 real(r8), intent(in ) :: rdrag(LBi:UBi,LBj:UBj)
379# if defined UV_QDRAG && !defined SOLVE3D
380 real(r8), intent(in ) :: rdrag2(LBi:UBi,LBj:UBj)
381# endif
382# if (defined UV_VIS2 || defined UV_VIS4) && !defined SOLVE3D
383 real(r8), intent(in ) :: pmon_r(LBi:UBi,LBj:UBj)
384 real(r8), intent(in ) :: pnom_r(LBi:UBi,LBj:UBj)
385 real(r8), intent(in ) :: pmon_p(LBi:UBi,LBj:UBj)
386 real(r8), intent(in ) :: pnom_p(LBi:UBi,LBj:UBj)
387 real(r8), intent(in ) :: om_r(LBi:UBi,LBj:UBj)
388 real(r8), intent(in ) :: on_r(LBi:UBi,LBj:UBj)
389 real(r8), intent(in ) :: om_p(LBi:UBi,LBj:UBj)
390 real(r8), intent(in ) :: on_p(LBi:UBi,LBj:UBj)
391# ifdef UV_VIS2
392 real(r8), intent(in ) :: visc2_p(LBi:UBi,LBj:UBj)
393 real(r8), intent(in ) :: visc2_r(LBi:UBi,LBj:UBj)
394# endif
395# ifdef UV_VIS4
396 real(r8), intent(in ) :: visc4_p(LBi:UBi,LBj:UBj)
397 real(r8), intent(in ) :: visc4_r(LBi:UBi,LBj:UBj)
398# endif
399# endif
400# if defined SEDIMENT && defined SED_MORPH
401 real(r8), intent(in ) :: bed_thick(LBi:UBi,LBj:UBj,1:3)
402# endif
403# if defined TIDE_GENERATING_FORCES && !defined SOLVE3D
404 real(r8), intent(in ) :: eq_tide(LBi:UBi,LBj:UBj)
405# endif
406# ifndef SOLVE3D
407 real(r8), intent(in ) :: sustr(LBi:UBi,LBj:UBj)
408 real(r8), intent(in ) :: svstr(LBi:UBi,LBj:UBj)
409# ifdef ATM_PRESS
410 real(r8), intent(in ) :: Pair(LBi:UBi,LBj:UBj)
411# endif
412# else
413# ifdef VAR_RHO_2D
414 real(r8), intent(in ) :: rhoA(LBi:UBi,LBj:UBj)
415 real(r8), intent(in ) :: rhoS(LBi:UBi,LBj:UBj)
416# endif
417 real(r8), intent(inout) :: DU_avg1(LBi:UBi,LBj:UBj)
418 real(r8), intent(inout) :: DU_avg2(LBi:UBi,LBj:UBj)
419 real(r8), intent(inout) :: DV_avg1(LBi:UBi,LBj:UBj)
420 real(r8), intent(inout) :: DV_avg2(LBi:UBi,LBj:UBj)
421 real(r8), intent(inout) :: Zt_avg1(LBi:UBi,LBj:UBj)
422 real(r8), intent(inout) :: rufrc(LBi:UBi,LBj:UBj)
423 real(r8), intent(inout) :: rvfrc(LBi:UBi,LBj:UBj)
424 real(r8), intent(inout) :: rufrc_bak(LBi:UBi,LBj:UBj,2)
425 real(r8), intent(inout) :: rvfrc_bak(LBi:UBi,LBj:UBj,2)
426# endif
427# ifdef WET_DRY
428 real(r8), intent(inout) :: pmask_full(LBi:UBi,LBj:UBj)
429 real(r8), intent(inout) :: rmask_full(LBi:UBi,LBj:UBj)
430 real(r8), intent(inout) :: umask_full(LBi:UBi,LBj:UBj)
431 real(r8), intent(inout) :: vmask_full(LBi:UBi,LBj:UBj)
432
433 real(r8), intent(inout) :: pmask_wet(LBi:UBi,LBj:UBj)
434 real(r8), intent(inout) :: rmask_wet(LBi:UBi,LBj:UBj)
435 real(r8), intent(inout) :: umask_wet(LBi:UBi,LBj:UBj)
436 real(r8), intent(inout) :: vmask_wet(LBi:UBi,LBj:UBj)
437# ifdef SOLVE3D
438 real(r8), intent(inout) :: rmask_wet_avg(LBi:UBi,LBj:UBj)
439# endif
440# endif
441 real(r8), intent(inout) :: ubar(LBi:UBi,LBj:UBj,:)
442 real(r8), intent(inout) :: vbar(LBi:UBi,LBj:UBj,:)
443 real(r8), intent(inout) :: zeta(LBi:UBi,LBj:UBj,:)
444# if defined NESTING && !defined SOLVE3D
445 real(r8), intent(out ) :: DU_flux(LBi:UBi,LBj:UBj)
446 real(r8), intent(out ) :: DV_flux(LBi:UBi,LBj:UBj)
447# endif
448#endif
449!
450! Local variable declarations.
451!
452 integer :: i, is, j
453 integer :: kbak, kold
454#ifdef DIAGNOSTICS_UV
455 integer :: idiag
456#endif
457!
458 real(r8) :: bkw0, bkw1, bkw2, bkw_new
459 real(r8) :: fwd0, fwd1, fwd2
460#ifdef SOLVE3D
461 real(r8) :: cfwd0, cfwd1, cfwd2
462#endif
463 real(r8) :: cff, cff1, cff2, cff3, cff4
464#ifdef WET_DRY
465 real(r8) :: cff5, cff6, cff7
466#endif
467 real(r8) :: fac, fac1, fac2
468#ifdef DEBUG
469 real(r8), parameter :: IniVal = 0.0_r8
470#endif
471!
472#if defined UV_C4ADVECTION && !defined SOLVE3D
473 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Dgrad
474#endif
475 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Dnew
476 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Dnew_rd
477 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Drhs
478#if (defined UV_VIS2 || defined UV_VIS4) && !defined SOLVE3D
479 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Drhs_p
480#endif
481 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Dstp
482 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: DUon
483 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: DVom
484#if defined STEP2D_CORIOLIS || !defined SOLVE3D
485 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: UFx
486 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: VFe
487#endif
488#if !defined SOLVE3D
489 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: UFe
490 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: VFx
491#endif
492#if defined UV_C4ADVECTION && !defined SOLVE3D
493 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: grad
494#endif
495 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rzeta2
496#if defined VAR_RHO_2D && defined SOLVE3D
497 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rzetaSA
498#endif
499 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rubar
500 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rvbar
501 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rzeta
502 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: urhs
503 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: vrhs
504 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: zwrk
505#ifdef WET_DRY
506 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: wetdry
507#endif
508#ifdef DIAGNOSTICS_UV
509 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Uwrk
510 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Vwrk
511 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,NDM2d-1) :: DiaU2rhs
512 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,NDM2d-1) :: DiaV2rhs
513#endif
514!
515 real(r8), allocatable :: zeta_new(:,:)
516
517#include "set_bounds.h"
518
519#ifdef DEBUG
520!
521!-----------------------------------------------------------------------
522! Initialize private arrays for debugging.
523!-----------------------------------------------------------------------
524!
525# if defined UV_C4ADVECTION && !defined SOLVE3D
526 dgrad=inival
527# endif
528 dnew=inival
529 dnew_rd=inival
530 drhs=inival
531# if (defined UV_VIS2 || defined UV_VIS4) && !defined SOLVE3D
532 drhs_p=inival
533# endif
534 dstp=inival
535 duon=inival
536 dvom=inival
537# if defined STEP2D_CORIOLIS || !defined SOLVE3D
538 ufx=inival
539 vfe=inival
540# endif
541# if !defined SOLVE3D
542 ufe=inival
543 vfx=inival
544# endif
545# if defined UV_C4ADVECTION && !defined SOLVE3D
546 grad=inival
547# endif
548 rzeta2=inival
549# if defined VAR_RHO_2D && defined SOLVE3D
550 rzetasa=inival
551# endif
552 rzeta=inival
553 rubar=inival
554 rvbar=inival
555 urhs=inival
556 vrhs=inival
557 zwrk=inival
558# ifdef WET_DRY
559 wetdry=inival
560# endif
561# ifdef DIAGNOSTICS_UV
562 uwrk=inival
563 vwrk=inival
564 diau2rhs=inival
565 diav2rhs=inival
566# endif
567#endif
568!
569!-----------------------------------------------------------------------
570! Set coefficients for AB3-AM4 forward-backward algorithm.
571!-----------------------------------------------------------------------
572!
573! Because the Forward Euler step is used to update "zeta" during the
574! first barotropic step, the pressure-gradient term in the momentum
575! equation must be computed via the Backward step to keep it
576! numerically stable. However, this interferes with the computation
577! of forcing terms "rufrc" and "rvfrc" because the free surface in
578! pressure gradient computation in 3D is exactly at the time
579! corresponding to baroclinic step "nstp" (rather than ahead by one
580! barotropic step after it updated by a normal forward-backward step).
581! To resolve this conflict, the pressure gradient term is computed in
582! two stages during the first barotropic step. It uses zeta(:,:,kstp)
583! at first to ensure exact consistency with 3D model. Then, after
584! vertical integrals of 3D right-hand-side "rufrc" and "rvfrc" are
585! converted into forcing terms, add correction based on the difference
586! zeta_new(:,:)-zeta(:,:,kstp) to "rubar" and "rvbar" to make them
587! consistent with the Backward step for pressure gradient.
588! For pressure gradient terms, search for the label PGF_FB_CORRECTION
589! below.
590!
591 IF (first_2d_step) THEN ! Meaning of time indices
592 kbak=kstp !------------------------
593 kold=kstp ! m-2 m-1 m m+1
594 fwd0=1.0_r8 ! kold kbak kstp knew
595 fwd1=0.0_r8 ! fwd2 fwd1 fwd0
596 fwd2=0.0_r8 ! bkw2 bkw1 bkw0 bkw_new
597#ifdef SOLVE3D
598 bkw_new=0.0_r8
599 bkw0=1.0_r8
600#else
601 bkw_new=1.0_r8
602 bkw0=0.0_r8
603#endif
604 bkw1=0.0_r8
605 bkw2=0.0_r8
606 ELSE IF (first_2d_step+1) THEN
607 kbak=kstp-1
608 IF (kbak.lt.1) kbak=4
609 kold=kbak
610 fwd0=1.0_r8 ! Logically AB2-AM3 forward-
611 fwd1=0.0_r8 ! backward scheme with maximum
612 fwd2=0.0_r8 ! stability coefficients while
613 bkw_new=1.0833333333333_r8 ! maintaining third-order
614 bkw0=-0.1666666666666_r8 ! accuracy, alpha_max=1.73
615 bkw1= 0.0833333333333_r8
616 bkw2= 0.0_r8
617 ELSE
618 kbak=kstp-1
619 IF (kbak.lt.1) kbak=4
620 kold=kbak-1
621 IF (kold.lt.1) kold=4
622 fwd0=1.781105_r8
623 fwd1=-1.06221_r8
624 fwd2=0.281105_r8
625 bkw_new=0.614_r8
626 bkw0=0.285_r8
627 bkw1=0.0880_r8
628 bkw2=0.013_r8
629 END IF
630
631#ifdef DEBUG
632!
633 IF (master) THEN
634 WRITE (20,10) iic(ng), iif(ng), kold, kbak, kstp, knew
635 10 FORMAT (' iic = ',i5.5,' iif = ',i3.3, &
636 & ' kold = ',i1,' kbak = ',i1,' kstp = ',i1,' knew = ',i1)
637 END IF
638#endif
639!
640!-----------------------------------------------------------------------
641! Preliminary steps.
642!-----------------------------------------------------------------------
643!
644! Compute total depth of water column and vertically integrated fluxes
645! needed for computing horizontal divergence to advance free surface
646! and for nonlinear advection terms for the barotropic momentum
647! equations.
648!
649#if defined DISTRIBUTE && !defined NESTING
650# define IR_RANGE IstrUm2-1,Iendp2
651# define JR_RANGE JstrVm2-1,Jendp2
652# define IU_RANGE IstrUm1-1,Iendp2
653# define JU_RANGE Jstrm1-1,Jendp2
654# define IV_RANGE Istrm1-1,Iendp2
655# define JV_RANGE JstrVm1-1,Jendp2
656#else
657# define IR_RANGE IstrUm2-1,Iendp2
658# define JR_RANGE JstrVm2-1,Jendp2
659# define IU_RANGE IstrUm2,Iendp2
660# define JU_RANGE JstrVm2-1,Jendp2
661# define IV_RANGE IstrUm2-1,Iendp2
662# define JV_RANGE JstrVm2,Jendp2
663#endif
664
665 DO j=jr_range
666 DO i=ir_range
667 drhs(i,j)=h(i,j)+fwd0*zeta(i,j,kstp)+ &
668 & fwd1*zeta(i,j,kbak)+ &
669 & fwd2*zeta(i,j,kold)
670 END DO
671 END DO
672!
673 DO j=ju_range
674 DO i=iu_range
675 cff=0.5_r8*on_u(i,j)
676 cff1=cff*(drhs(i,j)+drhs(i-1,j))
677 urhs(i,j)=fwd0*ubar(i,j,kstp)+ &
678 & fwd1*ubar(i,j,kbak)+ &
679 & fwd2*ubar(i,j,kold)
680 duon(i,j)=urhs(i,j)*cff1
681 END DO
682 END DO
683!
684 DO j=jv_range
685 DO i=iv_range
686 cff=0.5_r8*om_v(i,j)
687 cff1=cff*(drhs(i,j)+drhs(i,j-1))
688 vrhs(i,j)=fwd0*vbar(i,j,kstp)+ &
689 & fwd1*vbar(i,j,kbak)+ &
690 & fwd2*vbar(i,j,kold)
691 dvom(i,j)=vrhs(i,j)*cff1
692 END DO
693 END DO
694
695#undef IR_RANGE
696#undef IU_RANGE
697#undef IV_RANGE
698#undef JR_RANGE
699#undef JU_RANGE
700#undef JV_RANGE
701
702#if defined DISTRIBUTE && \
703 defined uv_adv && defined uv_c4advection && !defined SOLVE3D
704!
705! In distributed-memory, the I- and J-ranges are different and a
706! special exchange is done here to avoid having three ghost points
707! for high-order numerical stencils. Notice that a private array is
708! passed below to the exchange routine. It also applies periodic
709! boundary conditions, if appropriate and no partitions in I- or
710! J-directions.
711!
712 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
713 CALL exchange_u2d_tile (ng, tile, &
714 & imins, imaxs, jmins, jmaxs, &
715 & duon)
716 CALL exchange_v2d_tile (ng, tile, &
717 & imins, imaxs, jmins, jmaxs, &
718 & dvom)
719 END IF
720 CALL mp_exchange2d (ng, tile, inlm, 2, &
721 & imins, imaxs, jmins, jmaxs, &
722 & nghostpoints, &
723 & ewperiodic(ng), nsperiodic(ng), &
724 & duon, dvom)
725#endif
726!
727! Set vertically integrated mass fluxes DUon and DVom along the open
728! boundaries in such a way that the integral volume is conserved.
729! HGA: Need to resolve 'krhs' index here.
730!
731 IF (any(volcons(:,ng))) THEN
732 CALL set_duv_bc_tile (ng, tile, &
733 & lbi, ubi, lbj, ubj, &
734 & imins, imaxs, jmins, jmaxs, &
735 & krhs, &
736#ifdef MASKING
737 & umask, vmask, &
738#endif
739 & om_v, on_u, &
740 & ubar, vbar, &
741 & drhs, duon, dvom)
742 END IF
743!
744!-----------------------------------------------------------------------
745! Advance free-surface.
746!-----------------------------------------------------------------------
747!
748! Notice that the new local free-surface is allocated so it can be
749! passed as an argumment to "zetabc_local". An automatic array cannot
750! be used here because of weird memory problems.
751!
752 allocate ( zeta_new(imins:imaxs,jmins:jmaxs) )
753 zeta_new = 0.0_r8
754!
755! Compute "zeta_new" at new time step and interpolate it half-step
756! backward, "zwrk" for the subsequent computation of the tangent
757! linear barotropic pressure gradient. Here, we use the BASIC STATE
758! values. Thus, the nonlinear correction to the pressure-gradient
759! term from "kstp" to "knew" is not needed for Forward-Euler to
760! Forward-Backward steps (PGF_FB_CORRECTION method).
761!
762 DO j=jstrv-1,jend
763 DO i=istru-1,iend
764 fac=dtfast(ng)*pm(i,j)*pn(i,j)
765 zeta_new(i,j)=zeta(i,j,kstp)+ &
766 & fac*(duon(i,j)-duon(i+1,j)+ &
767 & dvom(i,j)-dvom(i,j+1))
768#ifdef MASKING
769 zeta_new(i,j)=zeta_new(i,j)*rmask(i,j)
770# ifdef WET_DRY
771 zeta_new(i,j)=zeta_new(i,j)+ &
772 & (dcrit(ng)-h(i,j))*(1.0_r8-rmask(i,j))
773# endif
774#endif
775 zwrk(i,j)=bkw_new*zeta_new(i,j)+ &
776 & bkw0*zeta(i,j,kstp)+ &
777 & bkw1*zeta(i,j,kbak)+ &
778 & bkw2*zeta(i,j,kold)
779
780#if defined VAR_RHO_2D && defined SOLVE3D
781 rzeta(i,j)=(1.0_r8+rhos(i,j))*zwrk(i,j)
782 rzeta2(i,j)=rzeta(i,j)*zwrk(i,j)
783 rzetasa(i,j)=zwrk(i,j)*(rhos(i,j)-rhoa(i,j))
784#else
785 rzeta(i,j)=zwrk(i,j)
786 rzeta2(i,j)=zwrk(i,j)*zwrk(i,j)
787#endif
788 END DO
789 END DO
790!
791! Apply mass point sources (volume vertical influx), if any.
792!
793! Dsrc(is) = 2, flow across grid cell w-face (positive or negative)
794!
795 IF (lwsrc(ng)) THEN
796 DO is=1,nsrc(ng)
797 IF (int(sources(ng)%Dsrc(is)).eq.2) THEN
798 i=sources(ng)%Isrc(is)
799 j=sources(ng)%Jsrc(is)
800 IF (((istrr.le.i).and.(i.le.iendr)).and. &
801 & ((jstrr.le.j).and.(j.le.jendr))) THEN
802 zeta_new(i,j)=zeta_new(i,j)+ &
803 & sources(ng)%Qbar(is)* &
804 & pm(i,j)*pn(i,j)*dtfast(ng)
805 END IF
806 END IF
807 END DO
808 END IF
809!
810! Apply boundary conditions to newly computed free-surface "zeta_new".
811! Notice that we are using the local routine, which passes the private
812! "zeta_new" array as argument.
813!
814! Here, we use the local "zetabc" since the private array "zeta_new"
815! is passed as an argument to allow computing the lateral boundary
816! conditions on the range IstrU-1:Iend and JstrV-1:Jend, so parallel
817! tile exchanges are avoided.
818!
819 CALL zetabc_local (ng, tile, &
820 & lbi, ubi, lbj, ubj, &
821 & imins, imaxs, jmins, jmaxs, &
822 & kstp, &
823 & zeta, &
824 & zeta_new)
825!
826! Load new computed free-surface into global state array.
827!
828 DO j=jstrr,jendr
829 DO i=istrr,iendr
830 zeta(i,j,knew)=zeta_new(i,j)
831 END DO
832 END DO
833
834#ifdef SOLVE3D
835!
836!-----------------------------------------------------------------------
837! Compute fast-time-averaged fields over all barotropic time steps.
838!-----------------------------------------------------------------------
839!
840! Reset/initialize arrays for averaged fields during the first
841! barotropic time step. Then, accumulate it time average. Include
842! physical boundary points, but not periodic ghost points or
843! computation distributed-memory computational margins.
844!
845 cff1=weight(1,iif(ng),ng)
846 cff2=weight(2,iif(ng),ng)
847!
848 IF (first_2d_step) THEN
849 DO j=jstrr,jendr
850 DO i=istrr,iendr
851 zt_avg1(i,j)=cff1*zeta(i,j,knew)
852 IF (i.ge.istr) THEN
853 du_avg1(i,j)=0.0_r8
854 du_avg2(i,j)=cff2*duon(i,j)
855 END IF
856 IF (j.ge.jstr) THEN
857 dv_avg1(i,j)=0.0_r8
858 dv_avg2(i,j)=cff2*dvom(i,j)
859 END IF
860 END DO
861 END DO
862 ELSE
863 DO j=jstrr,jendr
864 DO i=istrr,iendr
865 zt_avg1(i,j)=zt_avg1(i,j)+cff1*zeta(i,j,knew)
866 IF (i.ge.istr) THEN
867 du_avg2(i,j)=du_avg2(i,j)+cff2*duon(i,j)
868 END IF
869 IF (j.ge.jstr) THEN
870 dv_avg2(i,j)=dv_avg2(i,j)+cff2*dvom(i,j)
871 END IF
872 END DO
873 END DO
874 END IF
875#endif
876!
877!=======================================================================
878! Compute right-hand-side for the 2D momentum equations.
879!=======================================================================
880#ifdef SOLVE3D
881!
882! Notice that we are suppressing the computation of momentum advection,
883! Coriolis, and lateral viscosity terms in 3D Applications because
884! these terms are already included in the baroclinic-to-barotropic
885! forcing arrays "rufrc" and "rvfrc". It does not mean we are entirely
886! omitting them, but it is a choice between recomputing them at every
887! barotropic step or keeping them "frozen" during the fast-time
888! stepping.
889# ifdef STEP2D_CORIOLIS
890! However, in some coarse grid applications with larger baroclinic
891! timestep (say, DT around 20 minutes or larger), adding the Coriolis
892! term in the barotropic equations is useful since f*DT is no longer
893! small.
894# endif
895#endif
896!
897!-----------------------------------------------------------------------
898! Compute pressure-gradient terms.
899!-----------------------------------------------------------------------
900!
901 cff1=0.5_r8*g
902#if defined VAR_RHO_2D && defined SOLVE3D
903 cff2=0.333333333333_r8
904#endif
905#if defined ATM_PRESS && !defined SOLVE3D
906 cff3=0.5_r8*100.0_r8/rho0
907#endif
908 DO j=jstr,jend
909 DO i=istr,iend
910 IF (i.ge.istru) THEN
911 rubar(i,j)=cff1*on_u(i,j)* &
912 & ((h(i-1,j)+ &
913 & h(i ,j))* &
914 & (rzeta(i-1,j)- &
915 & rzeta(i ,j))+ &
916#if defined VAR_RHO_2D && defined SOLVE3D
917 & (h(i-1,j)- &
918 & h(i ,j))* &
919 & (rzetasa(i-1,j)+ &
920 & rzetasa(i ,j)+ &
921 & cff2*(rhoa(i-1,j)- &
922 & rhoa(i ,j))* &
923 & (zwrk(i-1,j)- &
924 & zwrk(i,j)))+ &
925#endif
926 & (rzeta2(i-1,j)- &
927 & rzeta2(i ,j)))
928#if defined ATM_PRESS && !defined SOLVE3D
929 rubar(i,j)=rubar(i,j)- &
930 & cff3*on_u(i,j)* &
931 & (h(i-1,j)+h(i,j)+ &
932 & rzeta(i-1,j)+rzeta(i,j))* &
933 & (pair(i,j)-pair(i-1,j))
934#endif
935#if defined TIDE_GENERATING_FORCES && !defined SOLVE3D
936 rubar(i,j)=rubar(i,j)- &
937 & cff1*on_u(i,j)* &
938 & (h(i-1,j)+h(i,j)+ &
939 & rzeta(i-1,j)+rzeta(i,j))* &
940 & (eq_tide(i,j)-eq_tide(i-1,j))
941#endif
942#ifdef DIAGNOSTICS_UV
943 diau2rhs(i,j,m2pgrd)=rubar(i,j)
944#endif
945 END IF
946!
947 IF (j.ge.jstrv) THEN
948 rvbar(i,j)=cff1*om_v(i,j)* &
949 & ((h(i,j-1)+ &
950 & h(i,j ))* &
951 & (rzeta(i,j-1)- &
952 & rzeta(i,j ))+ &
953#if defined VAR_RHO_2D && defined SOLVE3D
954 & (h(i,j-1)- &
955 & h(i,j ))* &
956 & (rzetasa(i,j-1)+ &
957 & rzetasa(i,j )+ &
958 & cff2*(rhoa(i,j-1)- &
959 & rhoa(i,j ))* &
960 & (zwrk(i,j-1)- &
961 & zwrk(i,j )))+ &
962#endif
963 & (rzeta2(i,j-1)- &
964 & rzeta2(i,j )))
965#if defined ATM_PRESS && !defined SOLVE3D
966 rvbar(i,j)=rvbar(i,j)- &
967 & cff3*om_v(i,j)* &
968 & (h(i,j-1)+h(i,j)+ &
969 & rzeta(i,j-1)+rzeta(i,j))* &
970 & (pair(i,j)-pair(i,j-1))
971#endif
972#if defined TIDE_GENERATING_FORCES && !defined SOLVE3D
973 rvbar(i,j)=rvbar(i,j)- &
974 & cff1*om_v(i,j)* &
975 & (h(i,j-1)+h(i,j)+ &
976 & rzeta(i,j-1)+rzeta(i,j))* &
977 & (eq_tide(i,j)-eq_tide(i,j-1))
978#endif
979#ifdef DIAGNOSTICS_UV
980 diav2rhs(i,j,m2pgrd)=rvbar(i,j)
981#endif
982 END IF
983 END DO
984 END DO
985
986#if defined UV_ADV && !defined SOLVE3D
987!
988!-----------------------------------------------------------------------
989! Add in horizontal advection of momentum.
990!-----------------------------------------------------------------------
991
992# ifdef UV_C2ADVECTION
993!
994! Second-order, centered differences advection fluxes.
995!
996 DO j=jstr,jend
997 DO i=istr-1,iend
998 IF (i.ge.istru-1) THEN
999 ufx(i,j)=0.25_r8* &
1000 & (duon(i,j)+duon(i+1,j))* &
1001 & (urhs(i ,j)+ &
1002 & urhs(i+1,j))
1003 END IF
1004!
1005 vfx(i+1,j)=0.25_r8* &
1006# ifdef MASKING
1007 & pmask(i+1,j)* &
1008# endif
1009 & (duon(i+1,j)+duon(i+1,j-1))* &
1010 & (vrhs(i+1,j)+ &
1011 & vrhs(i ,j))
1012 END DO
1013 END DO
1014!
1015 DO j=jstr-1,jend
1016 DO i=istr,iend
1017 IF (j.ge.jstrv-1) THEN
1018 vfe(i,j)=0.25_r8* &
1019 & (dvom(i,j)+dvom(i,j+1))* &
1020 & (vrhs(i,j )+ &
1021 & vrhs(i,j+1))
1022 END IF
1023!
1024 ufe(i,j+1)=0.25_r8* &
1025# ifdef MASKING
1026 & pmask(i,j+1)* &
1027# endif
1028 & (dvom(i,j+1)+dvom(i-1,j+1))* &
1029 & (urhs(i,j+1)+ &
1030 & urhs(i,j ))
1031 END DO
1032 END DO
1033
1034# elif defined UV_C4ADVECTION
1035!
1036! Fourth-order, centered differences u-momentum advection fluxes.
1037!
1038 DO j=jstr,jend
1039 DO i=istrum1,iendp1
1040 grad(i,j)=urhs(i-1,j)-2.0_r8*urhs(i,j)+ &
1041 & urhs(i+1,j)
1042 dgrad(i,j)=duon(i-1,j)-2.0_r8*duon(i,j)+duon(i+1,j)
1043 END DO
1044 END DO
1045 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
1046 IF (domain(ng)%Western_Edge(tile)) THEN
1047 DO j=jstr,jend
1048 grad(istr,j)=grad(istr+1,j)
1049 dgrad(istr,j)=dgrad(istr+1,j)
1050 END DO
1051 END IF
1052 END IF
1053 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
1054 IF (domain(ng)%Eastern_Edge(tile)) THEN
1055 DO j=jstr,jend
1056 grad(iend+1,j)=grad(iend,j)
1057 dgrad(iend+1,j)=dgrad(iend,j)
1058 END DO
1059 END IF
1060 END IF
1061! d/dx(Duu/n)
1062 cff=1.0_r8/6.0_r8
1063 DO j=jstr,jend
1064 DO i=istru-1,iend
1065 ufx(i,j)=0.25_r8*(urhs(i ,j)+ &
1066 & urhs(i+1,j)- &
1067 & cff*(grad(i,j)+grad(i+1,j)))* &
1068 & (duon(i,j)+duon(i+1,j)- &
1069 & cff*(dgrad(i,j)+dgrad(i+1,j)))
1070 END DO
1071 END DO
1072!
1073 DO j=jstrm1,jendp1
1074 DO i=istru,iend
1075 grad(i,j)=urhs(i,j-1)-2.0_r8*urhs(i,j)+ &
1076 & urhs(i,j+1)
1077 END DO
1078 END DO
1079 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
1080 IF (domain(ng)%Southern_Edge(tile)) THEN
1081 DO i=istru,iend
1082 grad(i,jstr-1)=grad(i,jstr)
1083 END DO
1084 END IF
1085 END IF
1086 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
1087 IF (domain(ng)%Northern_Edge(tile)) THEN
1088 DO i=istru,iend
1089 grad(i,jend+1)=grad(i,jend)
1090 END DO
1091 END IF
1092 END IF
1093 DO j=jstr,jend+1
1094 DO i=istru-1,iend
1095 dgrad(i,j)=dvom(i-1,j)-2.0_r8*dvom(i,j)+dvom(i+1,j)
1096 END DO
1097 END DO
1098! d/dy(Duv/m)
1099 cff=1.0_r8/6.0_r8
1100 DO j=jstr,jend+1
1101 DO i=istru,iend
1102 ufe(i,j)=0.25_r8*(urhs(i,j )+ &
1103 & urhs(i,j-1)- &
1104 & cff*(grad(i,j)+grad(i,j-1)))* &
1105 & (dvom(i,j)+dvom(i-1,j)- &
1106 & cff*(dgrad(i,j)+dgrad(i-1,j)))
1107 END DO
1108 END DO
1109!
1110! Fourth-order, centered differences v-momentum advection fluxes.
1111!
1112 DO j=jstrv,jend
1113 DO i=istrm1,iendp1
1114 grad(i,j)=vrhs(i-1,j)-2.0_r8*vrhs(i,j)+ &
1115 & vrhs(i+1,j)
1116 END DO
1117 END DO
1118 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
1119 IF (domain(ng)%Western_Edge(tile)) THEN
1120 DO j=jstrv,jend
1121 grad(istr-1,j)=grad(istr,j)
1122 END DO
1123 END IF
1124 END IF
1125 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
1126 IF (domain(ng)%Eastern_Edge(tile)) THEN
1127 DO j=jstrv,jend
1128 grad(iend+1,j)=grad(iend,j)
1129 END DO
1130 END IF
1131 END IF
1132 DO j=jstrv-1,jend
1133 DO i=istr,iend+1
1134 dgrad(i,j)=duon(i,j-1)-2.0_r8*duon(i,j)+duon(i,j+1)
1135 END DO
1136 END DO
1137! d/dx(Duv/n)
1138 cff=1.0_r8/6.0_r8
1139 DO j=jstrv,jend
1140 DO i=istr,iend+1
1141 vfx(i,j)=0.25_r8*(vrhs(i ,j)+ &
1142 & vrhs(i-1,j)- &
1143 & cff*(grad(i,j)+grad(i-1,j)))* &
1144 & (duon(i,j)+duon(i,j-1)- &
1145 & cff*(dgrad(i,j)+dgrad(i,j-1)))
1146 END DO
1147 END DO
1148!
1149 DO j=jstrvm1,jendp1
1150 DO i=istr,iend
1151 grad(i,j)=vrhs(i,j-1)-2.0_r8*vrhs(i,j)+ &
1152 & vrhs(i,j+1)
1153 dgrad(i,j)=dvom(i,j-1)-2.0_r8*dvom(i,j)+dvom(i,j+1)
1154 END DO
1155 END DO
1156 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
1157 IF (domain(ng)%Southern_Edge(tile)) THEN
1158 DO i=istr,iend
1159 grad(i,jstr)=grad(i,jstr+1)
1160 dgrad(i,jstr)=dgrad(i,jstr+1)
1161 END DO
1162 END IF
1163 END IF
1164 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
1165 IF (domain(ng)%Northern_Edge(tile)) THEN
1166 DO i=istr,iend
1167 grad(i,jend+1)=grad(i,jend)
1168 dgrad(i,jend+1)=dgrad(i,jend)
1169 END DO
1170 END IF
1171 END IF
1172! d/dy(Dvv/m)
1173 cff=1.0_r8/6.0_r8
1174 DO j=jstrv-1,jend
1175 DO i=istr,iend
1176 vfe(i,j)=0.25_r8*(vrhs(i,j )+ &
1177 & vrhs(i,j+1)- &
1178 & cff*(grad(i,j)+grad(i,j+1)))* &
1179 & (dvom(i,j)+dvom(i,j+1)- &
1180 & cff*(dgrad(i,j)+dgrad(i,j+1)))
1181 END DO
1182 END DO
1183# endif
1184!
1185! Add advection to RHS terms.
1186!
1187 DO j=jstr,jend
1188 DO i=istr,iend
1189 IF (i.ge.istru) THEN
1190 cff1=ufx(i,j)-ufx(i-1,j)
1191 cff2=ufe(i,j+1)-ufe(i,j)
1192 fac1=cff1+cff2
1193 rubar(i,j)=rubar(i,j)-fac1
1194# if defined DIAGNOSTICS_UV
1195 diau2rhs(i,j,m2xadv)=-cff1
1196 diau2rhs(i,j,m2yadv)=-cff2
1197 diau2rhs(i,j,m2hadv)=-fac1
1198# endif
1199 END IF
1200!
1201 IF (j.ge.jstrv) THEN
1202 cff3=vfx(i+1,j)-vfx(i,j)
1203 cff4=vfe(i,j)-vfe(i,j-1)
1204 fac2=cff3+cff4
1205 rvbar(i,j)=rvbar(i,j)-fac2
1206# if defined DIAGNOSTICS_UV
1207 diav2rhs(i,j,m2xadv)=-cff3
1208 diav2rhs(i,j,m2yadv)=-cff4
1209 diav2rhs(i,j,m2hadv)=-fac2
1210# endif
1211 END IF
1212 END DO
1213 END DO
1214#endif
1215
1216#if (defined UV_COR && !defined SOLVE3D) || defined STEP2D_CORIOLIS
1217!
1218!-----------------------------------------------------------------------
1219! Add in Coriolis term.
1220!-----------------------------------------------------------------------
1221!
1222 DO j=jstrv-1,jend
1223 DO i=istru-1,iend
1224 cff=0.5_r8*drhs(i,j)*fomn(i,j)
1225 ufx(i,j)=cff*(vrhs(i,j )+ &
1226 & vrhs(i,j+1))
1227 vfe(i,j)=cff*(urhs(i ,j)+ &
1228 & urhs(i+1,j))
1229 END DO
1230 END DO
1231!
1232 DO j=jstr,jend
1233 DO i=istr,iend
1234 IF (i.ge.istru) THEN
1235 fac1=0.5_r8*(ufx(i,j)+ufx(i-1,j))
1236 rubar(i,j)=rubar(i,j)+fac1
1237# if defined DIAGNOSTICS_UV
1238 diau2rhs(i,j,m2fcor)=fac1
1239# endif
1240
1241 END IF
1242!
1243 IF (j.ge.jstrv) THEN
1244 fac2=0.5_r8*(vfe(i,j)+vfe(i,j-1))
1245 rvbar(i,j)=rvbar(i,j)-fac2
1246# if defined DIAGNOSTICS_UV
1247 diav2rhs(i,j,m2fcor)=-fac2
1248# endif
1249 END IF
1250 END DO
1251 END DO
1252#endif
1253
1254#if (defined CURVGRID && defined UV_ADV) && !defined SOLVE3D
1255!
1256!-----------------------------------------------------------------------
1257! Add in curvilinear transformation terms.
1258!-----------------------------------------------------------------------
1259!
1260 DO j=jstrv-1,jend
1261 DO i=istru-1,iend
1262 cff1=0.5_r8*(vrhs(i,j )+ &
1263 & vrhs(i,j+1))
1264 cff2=0.5_r8*(urhs(i ,j)+
1265 & urhs(i+1,j))
1266 cff3=cff1*dndx(i,j)
1267 cff4=cff2*dmde(i,j)
1268 cff=drhs(i,j)*(cff3-cff4)
1269 ufx(i,j)=cff*cff1
1270 vfe(i,j)=cff*cff2
1271# if defined DIAGNOSTICS_UV
1272 cff=drhs(i,j)*cff4
1273 uwrk(i,j)=-cff*cff1 ! ubar equation, ETA-term
1274 vwrk(i,j)=-cff*cff2 ! vbar equation, ETA-term
1275# endif
1276 END DO
1277 END DO
1278!
1279 DO j=jstr,jend
1280 DO i=istr,iend
1281 IF (i.ge.istru) THEN
1282 fac1=0.5_r8*(ufx(i,j)+ufx(i-1,j))
1283 rubar(i,j)=rubar(i,j)+fac1
1284# if defined DIAGNOSTICS_UV
1285 fac2=0.5_r8*(uwrk(i,j)+uwrk(i-1,j))
1286 diau2rhs(i,j,m2xadv)=diau2rhs(i,j,m2xadv)+fac1-fac2
1287 diau2rhs(i,j,m2yadv)=diau2rhs(i,j,m2yadv)+fac2
1288 diau2rhs(i,j,m2hadv)=diau2rhs(i,j,m2hadv)+fac1
1289# endif
1290 END IF
1291!
1292 IF (j.ge.jstrv) THEN
1293 fac1=0.5_r8*(vfe(i,j)+vfe(i,j-1))
1294 rvbar(i,j)=rvbar(i,j)-fac1
1295# if defined DIAGNOSTICS_UV
1296 fac2=0.5_r8*(vwrk(i,j)+vwrk(i,j-1))
1297 diav2rhs(i,j,m2xadv)=diav2rhs(i,j,m2xadv)-fac1+fac2
1298 diav2rhs(i,j,m2yadv)=diav2rhs(i,j,m2yadv)-fac2
1299 diav2rhs(i,j,m2hadv)=diav2rhs(i,j,m2hadv)-fac1
1300# endif
1301 END IF
1302 END DO
1303 END DO
1304#endif
1305
1306#if defined UV_VIS2 && !defined SOLVE3D
1307!
1308!-----------------------------------------------------------------------
1309! Add in horizontal harmonic viscosity.
1310!-----------------------------------------------------------------------
1311!
1312! Compute total depth at PSI-points.
1313!
1314 DO j=jstr,jend+1
1315 DO i=istr,iend+1
1316 drhs_p(i,j)=0.25_r8*(drhs(i,j )+drhs(i-1,j )+ &
1317 & drhs(i,j-1)+drhs(i-1,j-1))
1318 END DO
1319 END DO
1320!
1321! Compute flux-components of the horizontal divergence of the stress
1322! tensor (m5/s2) in XI- and ETA-directions.
1323!
1324 DO j=jstrv-1,jend
1325 DO i=istru-1,iend
1326 cff=visc2_r(i,j)*drhs(i,j)*0.5_r8* &
1327 & (pmon_r(i,j)* &
1328 & ((pn(i ,j)+pn(i+1,j))*ubar(i+1,j,kstp)- &
1329 & (pn(i-1,j)+pn(i ,j))*ubar(i ,j,kstp))- &
1330 & pnom_r(i,j)* &
1331 & ((pm(i,j )+pm(i,j+1))*vbar(i,j+1,kstp)- &
1332 & (pm(i,j-1)+pm(i,j ))*vbar(i,j ,kstp)))
1333 ufx(i,j)=on_r(i,j)*on_r(i,j)*cff
1334 vfe(i,j)=om_r(i,j)*om_r(i,j)*cff
1335 END DO
1336 END DO
1337!
1338 DO j=jstr,jend+1
1339 DO i=istr,iend+1
1340 cff=visc2_p(i,j)*drhs_p(i,j)*0.5_r8* &
1341 & (pmon_p(i,j)* &
1342 & ((pn(i ,j-1)+pn(i ,j))*vbar(i ,j,kstp)- &
1343 & (pn(i-1,j-1)+pn(i-1,j))*vbar(i-1,j,kstp))+ &
1344 & pnom_p(i,j)* &
1345 & ((pm(i-1,j )+pm(i,j ))*ubar(i,j ,kstp)- &
1346 & (pm(i-1,j-1)+pm(i,j-1))*ubar(i,j-1,kstp)))
1347# ifdef MASKING
1348 cff=cff*pmask(i,j)
1349# endif
1350# ifdef WET_DRY
1351 cff=cff*pmask_wet(i,j)
1352# endif
1353 ufe(i,j)=om_p(i,j)*om_p(i,j)*cff
1354 vfx(i,j)=on_p(i,j)*on_p(i,j)*cff
1355 END DO
1356 END DO
1357!
1358! Add in harmonic viscosity.
1359!
1360 DO j=jstr,jend
1361 DO i=istr,iend
1362 IF (i.ge.istru) THEN
1363 cff1=0.5_r8*(pn(i-1,j)+pn(i,j))*(ufx(i,j )-ufx(i-1,j))
1364 cff2=0.5_r8*(pm(i-1,j)+pm(i,j))*(ufe(i,j+1)-ufe(i ,j))
1365 fac1=cff1+cff2
1366 rubar(i,j)=rubar(i,j)+fac1
1367# if defined DIAGNOSTICS_UV
1368 diau2rhs(i,j,m2hvis)=fac1
1369 diau2rhs(i,j,m2xvis)=cff1
1370 diau2rhs(i,j,m2yvis)=cff2
1371# endif
1372 END IF
1373!
1374 IF (j.ge.jstrv) THEN
1375 cff1=0.5_r8*(pn(i,j-1)+pn(i,j))*(vfx(i+1,j)-vfx(i,j ))
1376 cff2=0.5_r8*(pm(i,j-1)+pm(i,j))*(vfe(i ,j)-vfe(i,j-1))
1377 fac1=cff1-cff2
1378 rvbar(i,j)=rvbar(i,j)+fac1
1379# if defined DIAGNOSTICS_UV
1380 diav2rhs(i,j,m2hvis)=fac1
1381 diav2rhs(i,j,m2xvis)= cff1
1382 diav2rhs(i,j,m2yvis)=-cff2
1383# endif
1384 END IF
1385 END DO
1386 END DO
1387#endif
1388
1389#ifdef SOLVE3D
1390!
1391!-----------------------------------------------------------------------
1392! Coupling between 2D and 3D equations.
1393!-----------------------------------------------------------------------
1394!
1395! Before the first barotropic time step, arrays "rufrc" and "rvfrc"
1396! contain vertical integrals of the 3D right-hand-side terms for the
1397! momentum equations (including surface and bottom stresses). During
1398! the first barotropic time step, convert them into forcing terms by
1399! subtracting the fast-time "rubar" and "rvbar" from them.
1400!
1401! In the predictor-coupled mode, the resultant forcing terms "rufrc"
1402! and "rvfrc" are extrapolated forward in time, so they become
1403! centered effectively at time n+1/2. This is done using optimized
1404! Adams-Bashforth weights. In the code below, rufrc_bak(:,:,nstp) is
1405! at (n-1)time step, while rufrc_bak(:,:,3-nstp) is at (n-2). After
1406! its use as input, the latter is overwritten by the value at time
1407! step "nstp" (mathematically "n") during the next step.
1408!
1409! From now on, the computed forcing terms "rufrc" and "rvfrc" will
1410! remain constant during the fast-time stepping and will be added
1411! to "rubar" and "rvbar" during all subsequent barotropic steps.
1412!
1413 coupled_step : IF (first_2d_step) THEN
1414!
1415! Predictor coupled barotropic mode: Set coefficients for AB3-like
1416! forward-in-time extrapolation of 3D to 2D forcing terms "rufrc" and
1417! "rvfrc".
1418!
1419 IF (iic(ng).eq.ntstart(ng)) THEN
1420 cfwd0=1.0_r8
1421 cfwd1=0.0_r8
1422 cfwd2=0.0_r8
1423 ELSE IF (iic(ng).eq.ntstart(ng)+1) THEN
1424 cfwd0=1.5_r8
1425 cfwd1=-0.5_r8
1426 cfwd2=0.0_r8
1427 ELSE
1428 cfwd2=0.281105_r8
1429 cfwd1=-0.5_r8-2.0_r8*cfwd2
1430 cfwd0=1.5_r8+cfwd2
1431 END IF
1432!
1433 DO j=jstr,jend
1434 DO i=istr,iend
1435!
1436! Compensate for (cancel) bottom drag terms: at input into step2d
1437! "rufrc" and "rvfrc" contain bottom drag terms computed by 3D mode.
1438! However, there are no 2D counterparts in "rubar" and "rvbar" because
1439! 2D bottom drag will be computed implicitly during the final stage of
1440! updating ubar(:,:,knew) and vbar(:,:,knew) below. Note that unlike
1441! the other terms, the bottom drag should not be extrapolated forward,
1442! if "rufrc" and "rvfrc" are, so this cancelation needs to be done
1443! right now rather than at the bottom of this loop.
1444!
1445 IF (i.ge.istru) THEN
1446 rufrc(i,j)=rufrc(i,j)+ &
1447 & 0.5_r8*(rdrag(i,j)+rdrag(i-1,j))* &
1448 & om_u(i,j)*on_u(i,j)*ubar(i,j,kstp)
1449 END IF
1450!
1451 IF (j.ge.jstrv) THEN
1452 rvfrc(i,j)=rvfrc(i,j)+ &
1453 & 0.5_r8*(rdrag(i,j)+rdrag(i,j-1))* &
1454 & om_v(i,j)*on_v(i,j)*vbar(i,j,kstp)
1455 END IF
1456!
1457! Barotropic mode running predictor stage: forward extrapolation.
1458!
1459 IF (i.ge.istru) THEN
1460 cff1=rufrc(i,j)-rubar(i,j)
1461 rufrc(i,j)=cfwd0*cff1+ &
1462 & cfwd1*rufrc_bak(i,j, nstp)+ &
1463 & cfwd2*rufrc_bak(i,j,3-nstp)
1464 rufrc_bak(i,j,3-nstp)=cff1
1465 END IF
1466!
1467 IF (j.ge.jstrv) THEN
1468 cff2=rvfrc(i,j)-rvbar(i,j)
1469 rvfrc(i,j)=cfwd0*cff2+ &
1470 & cfwd1*rvfrc_bak(i,j, nstp)+ &
1471 & cfwd2*rvfrc_bak(i,j,3-nstp)
1472 rvfrc_bak(i,j,3-nstp)=cff2
1473 END IF
1474 END DO
1475 END DO
1476!
1477! Add correction term to shift pressure-gradient terms from "kstp" to
1478! "knew". That is, it converts the first 2D step from Forward-Euler
1479! to Forward-Backward (this is PGF_FB_CORRECTION mentioned above).
1480!
1481 DO j=jstrv-1,jend
1482 DO i=istru-1,iend
1483 zwrk(i,j)=zeta_new(i,j)-zeta(i,j,kstp)
1484# if defined VAR_RHO_2D && defined SOLVE3D
1485 rzeta(i,j)=(1.0_r8+rhos(i,j))*zwrk(i,j)
1486 rzeta2(i,j)=rzeta(i,j)*(zeta_new(i,j)+zeta(i,j,kstp))
1487 rzetasa(i,j)=zwrk(i,j)*(rhos(i,j)-rhoa(i,j))
1488# else
1489 rzeta(i,j)=zwrk(i,j)
1490 rzeta2(i,j)=zwrk(i,j)*(zeta_new(i,j)+zeta(i,j,kstp))
1491# endif
1492 END DO
1493 END DO
1494!
1495 cff1=0.5*g
1496# if defined VAR_RHO_2D && defined SOLVE3D
1497 cff2=0.333333333333_r8
1498# endif
1499 DO j=jstr,jend
1500 DO i=istr,iend
1501 IF (i.ge.istru) THEN
1502 rubar(i,j)=rubar(i,j)+ &
1503 & cff1*on_u(i,j)* &
1504 & ((h(i-1,j)+ &
1505 & h(i ,j))* &
1506 & (rzeta(i-1,j)- &
1507 & rzeta(i ,j))+ &
1508# if defined VAR_RHO_2D && defined SOLVE3D
1509 & (h(i-1,j)- &
1510 & h(i ,j))* &
1511 & (rzetasa(i-1,j)+ &
1512 & rzetasa(i ,j)+ &
1513 & cff2*(rhoa(i-1,j)- &
1514 & rhoa(i ,j))* &
1515 & (zwrk(i-1,j)- &
1516 & zwrk(i ,j)))+ &
1517# endif
1518 & (rzeta2(i-1,j)- &
1519 & rzeta2(i ,j)))
1520# ifdef DIAGNOSTICS_UV
1521 diau2rhs(i,j,m2pgrd)=diau2rhs(i,j,m2pgrd)+ &
1522 & rubar(i,j)
1523# endif
1524 END IF
1525!
1526 IF (j.ge.jstrv) THEN
1527 rvbar(i,j)=rvbar(i,j)+ &
1528 & cff1*om_v(i,j)* &
1529 & ((h(i,j-1)+ &
1530 & h(i,j ))* &
1531 & (rzeta(i,j-1)- &
1532 & rzeta(i,j ))+ &
1533# if defined VAR_RHO_2D && defined SOLVE3D
1534 & (h(i,j-1)- &
1535 & h(i,j ))* &
1536 & (rzetasa(i,j-1)+ &
1537 & rzetasa(i,j )+ &
1538 & cff2*(rhoa(i,j-1)- &
1539 & rhoa(i,j ))* &
1540 & (zwrk(i,j-1)- &
1541 & zwrk(i,j )))+ &
1542# endif
1543 & (rzeta2(i,j-1)- &
1544 & rzeta2(i,j )))
1545# ifdef DIAGNOSTICS_UV
1546 diav2rhs(i,j,m2pgrd)=diav2rhs(i,j,m2pgrd)+ &
1547 & rvbar(i,j)
1548# endif
1549 END IF
1550 END DO
1551 END DO
1552 END IF coupled_step
1553#endif
1554!
1555!-----------------------------------------------------------------------
1556! Time step 2D momentum equations.
1557!-----------------------------------------------------------------------
1558!
1559! Advance 2D momentum components while simultaneously adding them to
1560! accumulate fast-time averages to compute barotropic fluxes. Doing so
1561! straight away yields a more computationally dense code. However, the
1562! fast-time averaged fluxes (DU_avg1 and DV_avg1) are needed both at
1563! the interior and physical boundary points. Thus, we need separate
1564! loops along the domain boundaries after setting "ubar" and "vbar"
1565! lateral boundary conditions. Also, note that bottom drag is treated
1566! implicitly:
1567!
1568! Dnew*ubar(:,:,m+1) = Dold*ubar(:,:,m) +
1569! dtfast(ng)*rhs2D(:,:) -
1570! dtfast(ng)*rdrag(:,:)*ubar(:,:,m+1)
1571! hence
1572!
1573! ubar(:,:,m+1)=[Dold * ubar(..,m) + dtfast(ng) * rhs2D(:,:)] /
1574! [Dnew + dtfast(ng) * rdrag(:,:)]
1575!
1576! DU_avg1 = DU_avg1 +
1577! weight(m+1) * Dnew * ubar(:,:,m+1) * on_u(:,:)
1578!
1579! where it should be noted that Dnew .ne. Dnew + dtfast * rdrag
1580!
1581 DO j=jstrv-1,jend
1582 DO i=istru-1,iend
1583 dnew(i,j)=h(i,j)+zeta_new(i,j)
1584 dnew_rd(i,j)=dnew(i,j)+dtfast(ng)*rdrag(i,j)
1585 dstp(i,j)=h(i,j)+zeta(i,j,kstp)
1586 END DO
1587 END DO
1588
1589#if defined UV_QDRAG && !defined SOLVE3D
1590!
1591! Add quadratic drag term associated in shallow-water applications.
1592!
1593! Here, the SQRT(3) is due to a linear interpolation with second order
1594! accuaracy that ensures positive and negative values of the velocity
1595! components:
1596!
1597! u^2(i+1/2) = (1/3)*[u(i)*u(i) + u(i)*u(i+1) + u(i+1)*u(i+1)]
1598!
1599! If u(i)=1 and u(i+1)=-1, then u^2(i+1/2)=1/3 as it should be.
1600!
1601 cff=dtfast(ng)/sqrt(3.0_r8)
1602 DO j=jstrv-1,jend
1603 DO i=istru-1,iend
1604 cff1=ubar(i ,j,kstp)**2+ &
1605 & ubar(i+1,j,kstp)**2+ &
1606 & ubar(i ,j,kstp)*ubar(i+1,j,kstp)+ &
1607 & vbar(i,j ,kstp)**2+ &
1608 & vbar(i,j+1,kstp)**2+ &
1609 & vbar(i,j ,kstp)*vbar(i,j+1,kstp)
1610 cff2=sqrt(cff1)
1611 dnew_rd(i,j)=dnew_rd(i,j)+ &
1612 & cff*rdrag2(i,j)*cff2
1613 END DO
1614 END DO
1615#endif
1616!
1617! Step 2D momentum equations.
1618!
1619 cff=0.5_r8*dtfast(ng)
1620#ifdef SOLVE3D
1621 cff1=0.5_r8*weight(1,iif(ng),ng)
1622#else
1623 cff2=2.0_r8*dtfast(ng)
1624#endif
1625 DO j=jstr,jend
1626 DO i=istru,iend
1627 cff3=cff*(pm(i,j)+pm(i-1,j))*(pn(i,j)+pn(i-1,j))
1628 fac1=1.0_r8/(dnew_rd(i,j)+dnew_rd(i-1,j))
1629 ubar(i,j,knew)=fac1*((dstp(i,j)+dstp(i-1,j))*ubar(i,j,kstp)+ &
1630#ifdef SOLVE3D
1631 & cff3*(rubar(i,j)+rufrc(i,j)))
1632#else
1633 & cff3*rubar(i,j)+cff2*sustr(i,j))
1634#endif
1635#ifdef MASKING
1636 ubar(i,j,knew)=ubar(i,j,knew)*umask(i,j)
1637#endif
1638#ifdef WET_DRY
1639 cff5=abs(abs(umask_wet(i,j))-1.0_r8)
1640 cff6=0.5_r8+dsign(0.5_r8,ubar(i,j,knew))*umask_wet(i,j)
1641 cff7=0.5_r8*umask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
1642 ubar(i,j,knew)=ubar(i,j,knew)*cff7
1643#endif
1644#ifdef SOLVE3D
1645 du_avg1(i,j)=du_avg1(i,j)+ &
1646 & cff1*on_u(i,j)* &
1647 & (dnew(i,j)+dnew(i-1,j))*ubar(i,j,knew)
1648#endif
1649#if defined NESTING && !defined SOLVE3D
1650 du_flux(i,j)=0.5_r8*on_u(i,j)* &
1651 & (dnew(i,j)+dnew(i-1,j))*ubar(i,j,knew)
1652#endif
1653 END DO
1654 END DO
1655!
1656 DO j=jstrv,jend
1657 DO i=istr,iend
1658 cff3=cff*(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1))
1659 fac2=1.0_r8/(dnew_rd(i,j)+dnew_rd(i,j-1))
1660 vbar(i,j,knew)=fac2*((dstp(i,j)+dstp(i,j-1))*vbar(i,j,kstp)+ &
1661#ifdef SOLVE3D
1662 & cff3*(rvbar(i,j)+rvfrc(i,j)))
1663#else
1664 & cff3*rvbar(i,j)+cff2*svstr(i,j))
1665#endif
1666#ifdef MASKING
1667 vbar(i,j,knew)=vbar(i,j,knew)*vmask(i,j)
1668#endif
1669#ifdef WET_DRY
1670 cff5=abs(abs(vmask_wet(i,j))-1.0_r8)
1671 cff6=0.5_r8+dsign(0.5_r8,vbar(i,j,knew))*vmask_wet(i,j)
1672 cff7=0.5_r8*vmask_wet(i,j)*cff5+cff6*(1.0_r8-cff5)
1673 vbar(i,j,knew)=vbar(i,j,knew)*cff7
1674#endif
1675#ifdef SOLVE3D
1676 dv_avg1(i,j)=dv_avg1(i,j)+ &
1677 & cff1*om_v(i,j)* &
1678 & (dnew(i,j)+dnew(i,j-1))*vbar(i,j,knew)
1679#endif
1680#if defined NESTING && !defined SOLVE3D
1681 dv_flux(i,j)=0.5_r8*om_v(i,j)* &
1682 & (dnew(i,j)+dnew(i,j-1))*vbar(i,j,knew)
1683#endif
1684 END DO
1685 END DO
1686!
1687! Apply lateral boundary conditions.
1688!
1689 CALL u2dbc_tile (ng, tile, &
1690 & lbi, ubi, lbj, ubj, &
1691 & imins, imaxs, jmins, jmaxs, &
1692 & krhs, kstp, knew, &
1693 & ubar, vbar, zeta)
1694 CALL v2dbc_tile (ng, tile, &
1695 & lbi, ubi, lbj, ubj, &
1696 & imins, imaxs, jmins, jmaxs, &
1697 & krhs, kstp, knew, &
1698 & ubar, vbar, zeta)
1699!
1700! Compute integral mass flux across open boundaries and adjust
1701! for volume conservation.
1702!
1703 IF (any(volcons(:,ng))) THEN
1704 CALL obc_flux_tile (ng, tile, &
1705 & lbi, ubi, lbj, ubj, &
1706 & imins, imaxs, jmins, jmaxs, &
1707 & knew, &
1708#ifdef MASKING
1709 & umask, vmask, &
1710#endif
1711 & h, om_v, on_u, &
1712 & ubar, vbar, zeta)
1713 END IF
1714
1715#if defined SOLVE3D || (defined NESTING && !defined SOLVE3D)
1716!
1717! Set barotropic fluxes along physical boundaries.
1718!
1719 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
1720 IF (domain(ng)%Western_Edge(tile)) THEN
1721 DO j=jstr-1,jendr
1722 dnew(istr-1,j)=h(istr-1,j)+zeta_new(istr-1,j)
1723 END DO
1724 END IF
1725 END IF
1726 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
1727 IF (domain(ng)%Eastern_Edge(tile)) THEN
1728 DO j=jstr-1,jendr
1729 dnew(iend+1,j)=h(iend+1,j)+zeta_new(iend+1,j)
1730 END DO
1731 END IF
1732 END IF
1733 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
1734 IF (domain(ng)%Southern_Edge(tile)) THEN
1735 DO i=istr-1,iendr
1736 dnew(i,jstr-1)=h(i,jstr-1)+zeta_new(i,jstr-1)
1737 END DO
1738 END IF
1739 END IF
1740 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
1741 IF (domain(ng)%Northern_Edge(tile)) THEN
1742 DO i=istr-1,iendr
1743 dnew(i,jend+1)=h(i,jend+1)+zeta_new(i,jend+1)
1744 END DO
1745 END IF
1746 END IF
1747
1748# ifdef SOLVE3D
1749!
1750 cff1=0.5*weight(1,iif(ng),ng)
1751# endif
1752!
1753 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
1754 IF (domain(ng)%Western_Edge(tile)) THEN
1755 DO j=jstrr,jendr
1756# if defined NESTING && !defined SOLVE3D
1757 du_flux(istru-1,j)=0.5_r8*on_u(istru-1,j)* &
1758 & (dnew(istru-1,j)+dnew(istru-2,j))* &
1759 & ubar(istru-1,j,knew)
1760# else
1761 du_avg1(istru-1,j)=du_avg1(istru-1,j)+ &
1762 & cff1*on_u(istru-1,j)* &
1763 & (dnew(istru-1,j)+dnew(istru-2,j))* &
1764 & ubar(istru-1,j,knew)
1765# endif
1766 END DO
1767 DO j=jstrv,jend
1768# if defined NESTING && !defined SOLVE3D
1769 dv_flux(istr-1,j)=0.5_r8*om_v(istr-1,j)* &
1770 & (dnew(istr-1,j)+dnew(istr-1,j-1))* &
1771 & vbar(istr-1,j,knew)
1772# else
1773 dv_avg1(istr-1,j)=dv_avg1(istr-1,j)+ &
1774 & cff1*om_v(istr-1,j)* &
1775 & (dnew(istr-1,j)+dnew(istr-1,j-1))* &
1776 & vbar(istr-1,j,knew)
1777# endif
1778 END DO
1779 END IF
1780 END IF
1781 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
1782 IF (domain(ng)%Eastern_Edge(tile)) THEN
1783 DO j=jstrr,jendr
1784# if defined NESTING && !defined SOLVE3D
1785 du_flux(iend+1,j)=0.5_r8*on_u(iend+1,j)* &
1786 & (dnew(iend+1,j)+dnew(iend,j))* &
1787 & ubar(iend+1,j,knew)
1788# else
1789 du_avg1(iend+1,j)=du_avg1(iend+1,j)+ &
1790 & cff1*on_u(iend+1,j)* &
1791 & (dnew(iend+1,j)+dnew(iend,j))* &
1792 & ubar(iend+1,j,knew)
1793# endif
1794 END DO
1795 DO j=jstrv,jend
1796# if defined NESTING && !defined SOLVE3D
1797 dv_flux(iend+1,j)=0.5_r8*om_v(iend+1,j)* &
1798 & (dnew(iend+1,j)+dnew(iend+1,j-1))* &
1799 & vbar(iend+1,j,knew)
1800# else
1801 dv_avg1(iend+1,j)=dv_avg1(iend+1,j)+ &
1802 & cff1*om_v(iend+1,j)* &
1803 & (dnew(iend+1,j)+dnew(iend+1,j-1))* &
1804 & vbar(iend+1,j,knew)
1805# endif
1806 END DO
1807 END IF
1808 END IF
1809 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
1810 IF (domain(ng)%Southern_Edge(tile)) THEN
1811 DO i=istru,iend
1812# if defined NESTING && !defined SOLVE3D
1813 du_flux(i,jstr-1)=0.5_r8*on_u(i,jstr-1)* &
1814 & (dnew(i,jstr-1)+dnew(i-1,jstr-1))* &
1815 & ubar(i,jstr-1,knew)
1816# else
1817 du_avg1(i,jstr-1)=du_avg1(i,jstr-1)+ &
1818 & cff1*on_u(i,jstr-1)* &
1819 & (dnew(i,jstr-1)+dnew(i-1,jstr-1))* &
1820 & ubar(i,jstr-1,knew)
1821# endif
1822 END DO
1823 DO i=istrr,iendr
1824# if defined NESTING && !defined SOLVE3D
1825 dv_flux(i,jstrv-1)=0.5_r8*om_v(i,jstrv-1)* &
1826 & (dnew(i,jstrv-1)+dnew(i,jstrv-2))* &
1827 & vbar(i,jstrv-1,knew)
1828# else
1829 dv_avg1(i,jstrv-1)=dv_avg1(i,jstrv-1)+ &
1830 & cff1*om_v(i,jstrv-1)* &
1831 & (dnew(i,jstrv-1)+dnew(i,jstrv-2))* &
1832 & vbar(i,jstrv-1,knew)
1833# endif
1834 END DO
1835 END IF
1836 END IF
1837 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
1838 IF (domain(ng)%Northern_Edge(tile)) THEN
1839 DO i=istru,iend
1840# if defined NESTING && !defined SOLVE3D
1841 du_flux(i,jend+1)=0.5_r8*on_u(i,jend+1)* &
1842 & (dnew(i,jend+1)+dnew(i-1,jend+1))* &
1843 & ubar(i,jend+1,knew)
1844# else
1845 du_avg1(i,jend+1)=du_avg1(i,jend+1)+ &
1846 & cff1*on_u(i,jend+1)* &
1847 & (dnew(i,jend+1)+dnew(i-1,jend+1))* &
1848 & ubar(i,jend+1,knew)
1849# endif
1850 END DO
1851 DO i=istrr,iendr
1852# if defined NESTING && !defined SOLVE3D
1853 dv_flux(i,jend+1)=0.5_r8*om_v(i,jend+1)* &
1854 & (dnew(i,jend+1)+dnew(i,jend))* &
1855 & vbar(i,jend+1,knew)
1856# else
1857 dv_avg1(i,jend+1)=dv_avg1(i,jend+1)+ &
1858 & cff1*om_v(i,jend+1)* &
1859 & (dnew(i,jend+1)+dnew(i,jend))* &
1860 & vbar(i,jend+1,knew)
1861# endif
1862 END DO
1863 END IF
1864 END IF
1865#endif
1866!
1867!-----------------------------------------------------------------------
1868! Apply momentum transport point sources (like river runoff), if any.
1869!
1870! Dsrc(is) = 0, flow across grid cell u-face (positive or negative)
1871! Dsrc(is) = 1, flow across grid cell v-face (positive or negative)
1872!-----------------------------------------------------------------------
1873!
1874 IF (luvsrc(ng)) THEN
1875 DO is=1,nsrc(ng)
1876 i=sources(ng)%Isrc(is)
1877 j=sources(ng)%Jsrc(is)
1878 IF (((istrr.le.i).and.(i.le.iendr)).and. &
1879 & ((jstrr.le.j).and.(j.le.jendr))) THEN
1880 IF (int(sources(ng)%Dsrc(is)).eq.0) THEN
1881 cff=1.0_r8/(on_u(i,j)* &
1882 & 0.5_r8*(dnew(i-1,j)+dnew(i,j)))
1883 ubar(i,j,knew)=sources(ng)%Qbar(is)*cff
1884#ifdef SOLVE3D
1885 du_avg1(i,j)=sources(ng)%Qbar(is)
1886#endif
1887#if defined NESTING && !defined SOLVE3D
1888 du_flux(i,j)=sources(ng)%Qbar(is)
1889#endif
1890 ELSE IF (int(sources(ng)%Dsrc(is)).eq.1) THEN
1891 cff=1.0_r8/(om_v(i,j)* &
1892 & 0.5_r8*(dnew(i,j-1)+dnew(i,j)))
1893 vbar(i,j,knew)=sources(ng)%Qbar(is)*cff
1894#ifdef SOLVE3D
1895 dv_avg1(i,j)=sources(ng)%Qbar(is)
1896#endif
1897#if defined NESTING && !defined SOLVE3D
1898 dv_flux(i,j)=sources(ng)%Qbar(is)
1899#endif
1900 END IF
1901 END IF
1902 END DO
1903 END IF
1904!
1905! Deallocate local new free-surface.
1906!
1907 deallocate ( zeta_new )
1908
1909#ifdef WET_DRY
1910!
1911!-----------------------------------------------------------------------
1912! Compute new wet/dry masks.
1913!-----------------------------------------------------------------------
1914!
1915 CALL wetdry_tile (ng, tile, &
1916 & lbi, ubi, lbj, ubj, &
1917 & imins, imaxs, jmins, jmaxs, &
1918# ifdef MASKING
1919 & pmask, rmask, umask, vmask, &
1920# endif
1921 & h, zeta(:,:,knew), &
1922# ifdef SOLVE3D
1923 & du_avg1, dv_avg1, &
1924 & rmask_wet_avg, &
1925# endif
1926 & pmask_wet, pmask_full, &
1927 & rmask_wet, rmask_full, &
1928 & umask_wet, umask_full, &
1929 & vmask_wet, vmask_full)
1930#endif
1931
1932#ifdef SOLVE3D
1933!
1934!-----------------------------------------------------------------------
1935! At the end of the last 2D time step replace the new free-surface
1936! zeta(:,:,knew) with it fast time-averaged value, Zt_avg1. Recall
1937! this is state variable is the one that communicates with the 3D
1938! kernel. Then, compute time-dependent depths.
1939!-----------------------------------------------------------------------
1940!
1941 IF (iif(ng).eq.nfast(ng)) THEN
1942 DO j=jstrr,jendr
1943 DO i=istrr,iendr
1944 zeta(i,j,knew)=zt_avg1(i,j)
1945 END DO
1946 END DO
1947 CALL set_depth (ng, tile, inlm)
1948 END IF
1949#endif
1950
1951#ifdef NESTING
1952# ifdef SOLVE3D
1953!
1954!-----------------------------------------------------------------------
1955! If nesting and after all fast time steps are completed, exchange
1956! halo information to time averaged fields.
1957!-----------------------------------------------------------------------
1958!
1959 IF (iif(ng).eq.nfast(ng)) THEN
1960!
1961! In nesting applications with refinement grids, we need to exchange
1962! the DU_avg2 and DV_avg2 fluxes boundary information for the case
1963! that a contact point is at a tile partition. Notice that in such
1964! cases, we need i+1 and j+1 values for spatial/temporal interpolation.
1965!
1966 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
1967 CALL exchange_r2d_tile (ng, tile, &
1968 & lbi, ubi, lbj, ubj, &
1969 & zt_avg1)
1970 CALL exchange_u2d_tile (ng, tile, &
1971 & lbi, ubi, lbj, ubj, &
1972 & du_avg1)
1973 CALL exchange_v2d_tile (ng, tile, &
1974 & lbi, ubi, lbj, ubj, &
1975 & dv_avg1)
1976 CALL exchange_u2d_tile (ng, tile, &
1977 & lbi, ubi, lbj, ubj, &
1978 & du_avg2)
1979 CALL exchange_v2d_tile (ng, tile, &
1980 & lbi, ubi, lbj, ubj, &
1981 & dv_avg2)
1982 END IF
1983
1984# ifdef DISTRIBUTE
1985!
1986 CALL mp_exchange2d (ng, tile, inlm, 3, &
1987 & lbi, ubi, lbj, ubj, &
1988 & nghostpoints, &
1989 & ewperiodic(ng), nsperiodic(ng), &
1990 & zt_avg1, du_avg1, dv_avg1)
1991 CALL mp_exchange2d (ng, tile, inlm, 2, &
1992 & lbi, ubi, lbj, ubj, &
1993 & nghostpoints, &
1994 & ewperiodic(ng), nsperiodic(ng), &
1995 & du_avg2, dv_avg2)
1996# endif
1997 END IF
1998# else
1999!
2000! In nesting applications with refinement grids, we need to exchange
2001! the DU_flux and DV_flux fluxes boundary information for the case
2002! that a contact point is at a tile partition. Notice that in such
2003! cases, we need i+1 and j+1 values for spatial/temporal interpolation.
2004!
2005 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
2006 CALL exchange_u2d_tile (ng, tile, &
2007 & lbi, ubi, lbj, ubj, &
2008 & du_flux)
2009 CALL exchange_v2d_tile (ng, tile, &
2010 & lbi, ubi, lbj, ubj, &
2011 & dv_flux)
2012 END IF
2013
2014# ifdef DISTRIBUTE
2015!
2016 CALL mp_exchange2d (ng, tile, inlm, 2, &
2017 & lbi, ubi, lbj, ubj, &
2018 & nghostpoints, &
2019 & ewperiodic(ng), nsperiodic(ng), &
2020 & du_flux, dv_flux)
2021# endif
2022# endif
2023#endif
2024!
2025!-----------------------------------------------------------------------
2026! Exchange halo tile information.
2027!-----------------------------------------------------------------------
2028!
2029 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
2030 CALL exchange_r2d_tile (ng, tile, &
2031 & lbi, ubi, lbj, ubj, &
2032 & zeta(:,:,knew))
2033 CALL exchange_u2d_tile (ng, tile, &
2034 & lbi, ubi, lbj, ubj, &
2035 & ubar(:,:,knew))
2036 CALL exchange_v2d_tile (ng, tile, &
2037 & lbi, ubi, lbj, ubj, &
2038 & vbar(:,:,knew))
2039 END IF
2040
2041#ifdef DISTRIBUTE
2042!
2043 CALL mp_exchange2d (ng, tile, inlm, 3, &
2044 & lbi, ubi, lbj, ubj, &
2045 & nghostpoints, &
2046 & ewperiodic(ng), nsperiodic(ng), &
2047 & zeta(:,:,knew), &
2048 & ubar(:,:,knew), &
2049 & vbar(:,:,knew))
2050#endif
2051!
2052 RETURN

References mod_scalars::compositegrid, mod_scalars::dcrit, mod_param::domain, mod_scalars::dtfast, mod_scalars::ewperiodic, exchange_2d_mod::exchange_r2d_tile(), exchange_2d_mod::exchange_u2d_tile(), exchange_2d_mod::exchange_v2d_tile(), mod_scalars::g, mod_scalars::ieast, mod_scalars::iic, mod_scalars::iif, mod_param::inlm, mod_scalars::inorth, mod_scalars::isouth, mod_scalars::iwest, mod_scalars::luvsrc, mod_scalars::lwsrc, mod_scalars::m2fcor, mod_scalars::m2hadv, mod_scalars::m2hvis, mod_scalars::m2pgrd, mod_scalars::m2xadv, mod_scalars::m2xvis, mod_scalars::m2yadv, mod_scalars::m2yvis, mod_parallel::master, mp_exchange_mod::mp_exchange2d(), mod_scalars::nfast, mod_param::nghostpoints, mod_scalars::nsperiodic, mod_sources::nsrc, mod_scalars::ntstart, obc_volcons_mod::obc_flux_tile(), mod_scalars::rho0, set_depth_mod::set_depth(), obc_volcons_mod::set_duv_bc_tile(), mod_sources::sources, u2dbc_mod::u2dbc_tile(), v2dbc_mod::v2dbc_tile(), mod_scalars::volcons, mod_scalars::weight, and wetdry_mod::wetdry_tile().

Referenced by step2d().

Here is the call graph for this function:
Here is the caller graph for this function: