ROMS
Loading...
Searching...
No Matches
tl_rhs3d.F
Go to the documentation of this file.
1#include "cppdefs.h"
3#if defined TANGENT && defined SOLVE3D
4!
5!git $Id$
6!================================================== Hernan G. Arango ===
7! Copyright (c) 2002-2025 The ROMS Group Andrew M. Moore !
8! Licensed under a MIT/X style license !
9! See License_ROMS.md !
10!=======================================================================
11! !
12! This subroutine evaluates tangent linear right-hand-side terms !
13! for 3D momentum and tracers equations !
14! !
15! BASIC STATE variables needed: Hz, Huon, HVom, u, v, W, uclm, vclm, !
16! sustr, svstr, bustr, bvstr !
17! !
18!=======================================================================
19!
20 implicit none
21!
22 PRIVATE
23 PUBLIC :: tl_rhs3d
24!
25 CONTAINS
26!
27!***********************************************************************
28 SUBROUTINE tl_rhs3d (ng, tile)
29!***********************************************************************
30!
31 USE mod_param
32 USE mod_coupling
33# ifdef DIAGNOSTICS_UV
34 USE mod_diags
35# endif
36 USE mod_forces
37# if defined R4DVAR_ANA_SENSITIVITY && defined RPM_RELAXATION
38 USE mod_fourdvar
39# endif
40 USE mod_grid
41# ifdef WEC_MELLOR
42 USE mod_mixing
43# endif
44 USE mod_ocean
45 USE mod_stepping
46!
48 USE tl_prsgrd_mod, ONLY : tl_prsgrd
49# ifndef TS_FIXED
50# ifdef TS_DIF2
51 USE tl_t3dmix2_mod, ONLY : tl_t3dmix2
52# endif
53# ifdef TS_DIF4
54 USE tl_t3dmix4_mod, ONLY : tl_t3dmix4
55# endif
56# endif
57# if defined R4DVAR_ANA_SENSITIVITY && defined RPM_RELAXATION
60# endif
61# ifdef UV_VIS2
63# endif
64# ifdef UV_VIS4
66# endif
67!
68! Imported variable declarations.
69!
70 integer, intent(in) :: ng, tile
71!
72! Local variable declarations.
73!
74 character (len=*), parameter :: myfile = &
75 & __FILE__
76!
77# include "tile.h"
78!
79!-----------------------------------------------------------------------
80! Initialize computations for new time step of the 3D primitive
81! variables.
82!-----------------------------------------------------------------------
83!
84 CALL tl_pre_step3d (ng, tile)
85!
86!-----------------------------------------------------------------------
87! Compute baroclinic pressure gradient.
88!-----------------------------------------------------------------------
89!
90 CALL tl_prsgrd (ng, tile)
91# ifndef TS_FIXED
92# ifdef TS_DIF2
93!
94!-----------------------------------------------------------------------
95! Compute horizontal harmonic mixing of tracer type variables.
96!-----------------------------------------------------------------------
97!
98 CALL tl_t3dmix2 (ng, tile)
99# endif
100# ifdef TS_DIF4
101!
102!-----------------------------------------------------------------------
103! Compute horizontal biharmonic mixing of tracer type variables.
104!-----------------------------------------------------------------------
105!
106 CALL tl_t3dmix4 (ng, tile)
107# endif
108# endif
109# if defined R4DVAR_ANA_SENSITIVITY && defined RPM_RELAXATION
110!
111!-----------------------------------------------------------------------
112! Improve stability and convergence of the tangent linear representer
113! model tracer type variables by a "diffusive relaxation" to previous
114! Picard iteration solution.
115!-----------------------------------------------------------------------
116!
117 IF (lweakrelax(ng)) THEN
118 CALL tl_t3drelax (ng, tile)
119 END IF
120# endif
121!
122!-----------------------------------------------------------------------
123! Compute right-hand-side terms for the 3D momentum equations.
124!-----------------------------------------------------------------------
125!
126# ifdef PROFILE
127 CALL wclock_on (ng, itlm, 21, __line__, myfile)
128# endif
129 CALL tl_rhs3d_tile (ng, tile, &
130 & lbi, ubi, lbj, ubj, &
131 & imins, imaxs, jmins, jmaxs, &
132 & nrhs(ng), &
133 & grid(ng) % Hz, &
134 & grid(ng) % tl_Hz, &
135 & grid(ng) % Huon, &
136 & grid(ng) % tl_Huon, &
137 & grid(ng) % Hvom, &
138 & grid(ng) % tl_Hvom, &
139# if defined CURVGRID && defined UV_ADV
140 & grid(ng) % dmde, &
141 & grid(ng) % dndx, &
142# endif
143 & grid(ng) % fomn, &
144 & grid(ng) % om_u, &
145 & grid(ng) % om_v, &
146 & grid(ng) % on_u, &
147 & grid(ng) % on_v, &
148 & grid(ng) % pm, &
149 & grid(ng) % pn, &
150# ifdef WET_DRY_NOT_YET
151 & grid(ng)%umask_wet, &
152 & grid(ng)%vmask_wet, &
153# endif
154 & forces(ng) % bustr, &
155 & forces(ng) % tl_bustr, &
156 & forces(ng) % bvstr, &
157 & forces(ng) % tl_bvstr, &
158 & forces(ng) % sustr, &
159 & forces(ng) % tl_sustr, &
160 & forces(ng) % svstr, &
161 & forces(ng) % tl_svstr, &
162 & ocean(ng) % u, &
163 & ocean(ng) % tl_u, &
164 & ocean(ng) % v, &
165 & ocean(ng) % tl_v, &
166 & ocean(ng) % W, &
167 & ocean(ng) % tl_W, &
168# ifdef WEC_MELLOR
169 & ocean(ng) % u_stokes, &
170 & ocean(ng) % tl_u_stokes, &
171 & ocean(ng) % v_stokes, &
172 & ocean(ng) % tl_v_stokes, &
173 & ocean(ng) % tl_rulag3d, &
174 & ocean(ng) % tl_rvlag3d, &
175 & mixing(ng) % tl_rustr3d, &
176 & mixing(ng) % tl_rvstr3d, &
177# endif
178 & coupling(ng) % tl_rufrc, &
179 & coupling(ng) % tl_rvfrc, &
180# ifdef DIAGNOSTICS_UV
181!! & DIAGS(ng) % DiaRUfrc, &
182!! & DIAGS(ng) % DiaRVfrc, &
183!! & DIAGS(ng) % DiaRU, &
184!! & DIAGS(ng) % DiaRV, &
185# endif
186 & ocean(ng) % tl_ru, &
187 & ocean(ng) % tl_rv)
188# ifdef PROFILE
189 CALL wclock_off (ng, itlm, 21, __line__, myfile)
190# endif
191# ifdef UV_VIS2
192!
193!-----------------------------------------------------------------------
194! Compute horizontal, harmonic mixing of momentum.
195!-----------------------------------------------------------------------
196!
197 CALL tl_uv3dmix2 (ng, tile)
198# endif
199# ifdef UV_VIS4
200!
201!-----------------------------------------------------------------------
202! Compute horizontal, biharmonic mixing of momentum.
203!-----------------------------------------------------------------------
204!
205 CALL tl_uv3dmix4 (ng, tile)
206# endif
207# if defined R4DVAR_ANA_SENSITIVITY && defined RPM_RELAXATION
208!
209!-----------------------------------------------------------------------
210! Improve stability and convergence of the tangent linear representer
211! model 3D momentum by a "diffusive relaxation" to previous Picard
212! iteration solution.
213!-----------------------------------------------------------------------
214!
215 IF (lweakrelax(ng)) THEN
216 CALL tl_uv3drelax (ng, tile)
217 END IF
218# endif
219
220 RETURN
221 END SUBROUTINE tl_rhs3d
222!
223!***********************************************************************
224 SUBROUTINE tl_rhs3d_tile (ng, tile, &
225 & LBi, UBi, LBj, UBj, &
226 & IminS, ImaxS, JminS, JmaxS, &
227 & nrhs, &
228 & Hz, tl_Hz, &
229 & Huon, tl_Huon, &
230 & Hvom, tl_Hvom, &
231# if defined CURVGRID && defined UV_ADV
232 & dmde, dndx, &
233# endif
234 & fomn, &
235 & om_u, om_v, on_u, on_v, pm, pn, &
236# ifdef WET_DRY_NOT_YET
237 & umask_wet, vmask_wet, &
238# endif
239 & bustr, tl_bustr, &
240 & bvstr, tl_bvstr, &
241 & sustr, tl_sustr, &
242 & svstr, tl_svstr, &
243 & u, tl_u, &
244 & v, tl_v, &
245 & W, tl_W, &
246# ifdef WEC_MELLOR
247 & u_stokes, tl_u_stokes, &
248 & v_stokes, tl_v_stokes, &
249 & tl_rulag3d, tl_rvlag3d, &
250 & tl_rustr3d, tl_rvstr3d, &
251# endif
252 & tl_rufrc, &
253 & tl_rvfrc, &
254# ifdef DIAGNOSTICS_UV
255!! & DiaRUfrc, DiaRVfrc, &
256!! & DiaRU, DiaRV, &
257# endif
258 & tl_ru, tl_rv)
259!***********************************************************************
260!
261 USE mod_param
262 USE mod_clima
263 USE mod_scalars
264!
265! Imported variable declarations.
266!
267 integer, intent(in) :: ng, tile
268 integer, intent(in) :: LBi, UBi, LBj, UBj
269 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
270 integer, intent(in) :: nrhs
271!
272# ifdef ASSUMED_SHAPE
273 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
274 real(r8), intent(in) :: Huon(LBi:,LBj:,:)
275 real(r8), intent(in) :: Hvom(LBi:,LBj:,:)
276# if defined CURVGRID && defined UV_ADV
277 real(r8), intent(in) :: dmde(LBi:,LBj:)
278 real(r8), intent(in) :: dndx(LBi:,LBj:)
279# endif
280 real(r8), intent(in) :: fomn(LBi:,LBj:)
281 real(r8), intent(in) :: om_u(LBi:,LBj:)
282 real(r8), intent(in) :: om_v(LBi:,LBj:)
283 real(r8), intent(in) :: on_u(LBi:,LBj:)
284 real(r8), intent(in) :: on_v(LBi:,LBj:)
285 real(r8), intent(in) :: pm(LBi:,LBj:)
286 real(r8), intent(in) :: pn(LBi:,LBj:)
287# ifdef WET_DRY_NOT_YET
288 real(r8), intent(in) :: umask_wet(LBi:,LBj:)
289 real(r8), intent(in) :: vmask_wet(LBi:,LBj:)
290# endif
291 real(r8), intent(in) :: bustr(LBi:,LBj:)
292 real(r8), intent(in) :: bvstr(LBi:,LBj:)
293 real(r8), intent(in) :: sustr(LBi:,LBj:)
294 real(r8), intent(in) :: svstr(LBi:,LBj:)
295 real(r8), intent(in) :: u(LBi:,LBj:,:,:)
296 real(r8), intent(in) :: v(LBi:,LBj:,:,:)
297 real(r8), intent(in) :: W(LBi:,LBj:,0:)
298
299 real(r8), intent(in) :: tl_Hz(LBi:,LBj:,:)
300 real(r8), intent(in) :: tl_Huon(LBi:,LBj:,:)
301 real(r8), intent(in) :: tl_Hvom(LBi:,LBj:,:)
302 real(r8), intent(in) :: tl_bustr(LBi:,LBj:)
303 real(r8), intent(in) :: tl_bvstr(LBi:,LBj:)
304 real(r8), intent(in) :: tl_sustr(LBi:,LBj:)
305 real(r8), intent(in) :: tl_svstr(LBi:,LBj:)
306 real(r8), intent(in) :: tl_u(LBi:,LBj:,:,:)
307 real(r8), intent(in) :: tl_v(LBi:,LBj:,:,:)
308 real(r8), intent(in) :: tl_W(LBi:,LBj:,0:)
309# ifdef WEC_MELLOR
310 real(r8), intent(in) :: u_stokes(LBi:,LBj:,:)
311 real(r8), intent(in) :: v_stokes(LBi:,LBj:,:)
312 real(r8), intent(in) :: tl_u_stokes(LBi:,LBj:,:)
313 real(r8), intent(in) :: tl_v_stokes(LBi:,LBj:,:)
314 real(r8), intent(in) :: tl_rulag3d(LBi:,LBj:,:)
315 real(r8), intent(in) :: tl_rvlag3d(LBi:,LBj:,:)
316 real(r8), intent(in) :: tl_rustr3d(LBi:,LBj:,:)
317 real(r8), intent(in) :: tl_rvstr3d(LBi:,LBj:,:)
318# endif
319# ifdef DIAGNOSTICS_UV
320!! real(r8), intent(inout) :: DiaRUfrc(LBi:,LBj:,:,:)
321!! real(r8), intent(inout) :: DiaRVfrc(LBi:,LBj:,:,:)
322!! real(r8), intent(inout) :: DiaRU(LBi:,LBj:,:,:,:)
323!! real(r8), intent(inout) :: DiaRV(LBi:,LBj:,:,:,:)
324# endif
325 real(r8), intent(inout) :: tl_ru(LBi:,LBj:,0:,:)
326 real(r8), intent(inout) :: tl_rv(LBi:,LBj:,0:,:)
327
328 real(r8), intent(out) :: tl_rufrc(LBi:,LBj:)
329 real(r8), intent(out) :: tl_rvfrc(LBi:,LBj:)
330# else
331 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
332 real(r8), intent(in) :: Huon(LBi:UBi,LBj:UBj,N(ng))
333 real(r8), intent(in) :: Hvom(LBi:UBi,LBj:UBj,N(ng))
334# if defined CURVGRID && defined UV_ADV
335 real(r8), intent(in) :: dmde(LBi:UBi,LBj:UBj)
336 real(r8), intent(in) :: dndx(LBi:UBi,LBj:UBj)
337# endif
338 real(r8), intent(in) :: fomn(LBi:UBi,LBj:UBj)
339 real(r8), intent(in) :: om_u(LBi:UBi,LBj:UBj)
340 real(r8), intent(in) :: om_v(LBi:UBi,LBj:UBj)
341 real(r8), intent(in) :: on_u(LBi:UBi,LBj:UBj)
342 real(r8), intent(in) :: on_v(LBi:UBi,LBj:UBj)
343 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
344 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
345# ifdef WET_DRY_NOT_YET
346 real(r8), intent(in) :: umask_wet(LBi:UBi,LBj:UBj)
347 real(r8), intent(in) :: vmask_wet(LBi:UBi,LBj:UBj)
348# endif
349 real(r8), intent(in) :: bustr(LBi:UBi,LBj:UBj)
350 real(r8), intent(in) :: bvstr(LBi:UBi,LBj:UBj)
351 real(r8), intent(in) :: sustr(LBi:UBi,LBj:UBj)
352 real(r8), intent(in) :: svstr(LBi:UBi,LBj:UBj)
353 real(r8), intent(in) :: u(LBi:UBi,LBj:UBj,N(ng),2)
354 real(r8), intent(in) :: v(LBi:UBi,LBj:UBj,N(ng),2)
355 real(r8), intent(in) :: W(LBi:UBi,LBj:UBj,0:N(ng))
356
357 real(r8), intent(in) :: tl_Hz(LBi:UBi,LBj:UBj,N(ng))
358 real(r8), intent(in) :: tl_Huon(LBi:UBi,LBj:UBj,N(ng))
359 real(r8), intent(in) :: tl_Hvom(LBi:UBi,LBj:UBj,N(ng))
360 real(r8), intent(in) :: tl_bustr(LBi:UBi,LBj:UBj)
361 real(r8), intent(in) :: tl_bvstr(LBi:UBi,LBj:UBj)
362 real(r8), intent(in) :: tl_sustr(LBi:UBi,LBj:UBj)
363 real(r8), intent(in) :: tl_svstr(LBi:UBi,LBj:UBj)
364 real(r8), intent(in) :: tl_u(LBi:UBi,LBj:UBj,N(ng),2)
365 real(r8), intent(in) :: tl_v(LBi:UBi,LBj:UBj,N(ng),2)
366 real(r8), intent(in) :: tl_W(LBi:UBi,LBj:UBj,0:N(ng))
367# ifdef WEC_MELLOR
368 real(r8), intent(in) :: u_stokes(LBi:UBi,LBj:UBj,N(ng))
369 real(r8), intent(in) :: v_stokes(LBi:UBi,LBj:UBj,N(ng))
370 real(r8), intent(in) :: tl_u_stokes(LBi:UBi,LBj:UBj,N(ng))
371 real(r8), intent(in) :: tl_v_stokes(LBi:UBi,LBj:UBj,N(ng))
372 real(r8), intent(in) :: tl_rulag3d(LBi:UBi,LBj:UBj,N(ng))
373 real(r8), intent(in) :: tl_rvlag3d(LBi:UBi,LBj:UBj,N(ng))
374 real(r8), intent(in) :: tl_rustr3d(LBi:UBi,LBj:UBj,N(ng))
375 real(r8), intent(in) :: tl_rvstr3d(LBi:UBi,LBj:UBj,N(ng))
376# endif
377# ifdef DIAGNOSTICS_UV
378!! real(r8), intent(inout) :: DiaRUfrc(LBi:UBi,LBj:UBj,3,NDM2d-1)
379!! real(r8), intent(inout) :: DiaRVfrc(LBi:UBi,LBj:UBj,3,NDM2d-1)
380!! real(r8), intent(inout) :: DiaRU(LBi:UBi,LBj:UBj,N(ng),2,NDrhs)
381!! real(r8), intent(inout) :: DiaRV(LBi:UBi,LBj:UBj,N(ng),2,NDrhs)
382# endif
383 real(r8), intent(inout) :: tl_ru(LBi:UBi,LBj:UBj,0:N(ng),2)
384 real(r8), intent(inout) :: tl_rv(LBi:UBi,LBj:UBj,0:N(ng),2)
385
386 real(r8), intent(out) :: tl_rufrc(LBi:UBi,LBj:UBj)
387 real(r8), intent(out) :: tl_rvfrc(LBi:UBi,LBj:UBj)
388# endif
389!
390! Local variable declarations.
391!
392 integer :: i, j, k
393
394 real(r8), parameter :: Gadv = -0.25_r8
395
396 real(r8) :: cff, cff1, cff2, cff3, cff4
397 real(r8) :: tl_cff, tl_cff1, tl_cff2, tl_cff3, tl_cff4
398
399 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: CF
400 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: DC
401 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: FC
402
403 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: tl_CF
404 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: tl_DC
405 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: tl_FC
406
407 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Huee
408 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Huxx
409 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Hvee
410 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Hvxx
411 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: UFx
412 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: UFe
413 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Uwrk
414 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: VFx
415 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: VFe
416 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Vwrk
417 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: uee
418 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: uxx
419 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: vee
420 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: vxx
421 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: wrk
422
423 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_Huee
424 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_Huxx
425 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_Hvee
426 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_Hvxx
427 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_UFx
428 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_UFe
429 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_Uwrk
430 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_VFx
431 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_VFe
432 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_Vwrk
433 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_uee
434 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_uxx
435 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_vee
436 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_vxx
437 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_wrk
438
439# include "set_bounds.h"
440
441# ifdef BODYFORCE
442!
443!-----------------------------------------------------------------------
444! Apply surface stress as a bodyforce: determine the thickness (m)
445! of the surface layer; then add in surface stress as a bodyfoce.
446!-----------------------------------------------------------------------
447!
448# ifdef DIAGNOSTICS_UV
449!! DO k=1,N(ng)
450!! DO j=Jstr,Jend
451!! DO i=Istr,Iend
452!! DiaRU(i,j,k,nrhs,M3vvis)=0.0_r8
453!! DiaRV(i,j,k,nrhs,M3vvis)=0.0_r8
454!! END DO
455!! END DO
456!! END DO
457!! DO j=Jstr,Jend
458!! DO i=IstrU,Iend
459!! DiaRUfrc(i,j,3,M2sstr)=0.0_r8
460!! DiaRUfrc(i,j,3,M2bstr)=0.0_r8
461!! END DO
462!! END DO
463!! DO j=JstrV,Jend
464!! DO i=Istr,Iend
465!! DiaRVfrc(i,j,3,M2sstr)=0.0_r8
466!! DiaRVfrc(i,j,3,M2bstr)=0.0_r8
467!! END DO
468!! END DO
469# endif
470 DO j=jstrv-1,jend
471 DO i=istru-1,iend
472 wrk(i,j)=0.0_r8
473 tl_wrk(i,j)=0.0_r8
474 END DO
475 END DO
476 DO k=n(ng),levsfrc(ng),-1
477 DO j=jstrv-1,jend
478 DO i=istru-1,iend
479 wrk(i,j)=wrk(i,j)+hz(i,j,k)
480 tl_wrk(i,j)=tl_wrk(i,j)+tl_hz(i,j,k)
481 END DO
482 END DO
483 END DO
484 DO j=jstr,jend
485 DO i=istru,iend
486 cff=0.25_r8*(pm(i-1,j)+pm(i,j))* &
487 & (pn(i-1,j)+pn(i,j))
488 cff1=1.0_r8/(cff*(wrk(i-1,j)+wrk(i,j)))
489 tl_cff1=-cff1*cff1*cff*(tl_wrk(i-1,j)+tl_wrk(i,j))
490 uwrk(i,j)=sustr(i,j)*cff1
491 tl_uwrk(i,j)=tl_sustr(i,j)*cff1+ &
492 & sustr(i,j)*tl_cff1
493 END DO
494 END DO
495 DO j=jstrv,jend
496 DO i=istr,iend
497 cff=0.25*(pm(i,j-1)+pm(i,j))* &
498 & (pn(i,j-1)+pn(i,j))
499 cff1=1.0_r8/(cff*(wrk(i,j-1)+wrk(i,j)))
500 tl_cff1=-cff1*cff1*cff*(tl_wrk(i,j-1)+tl_wrk(i,j))
501 vwrk(i,j)=svstr(i,j)*cff1
502 tl_vwrk(i,j)=tl_svstr(i,j)*cff1+ &
503 & svstr(i,j)*tl_cff1
504 END DO
505 END DO
506 DO k=levsfrc(ng),n(ng)
507 DO j=jstr,jend
508 DO i=istru,iend
509!^ cff=Uwrk(i,j)*(Hz(i ,j,k)+ &
510!^ & Hz(i-1,j,k))
511!^
512 tl_cff=tl_uwrk(i,j)*(hz(i ,j,k)+ &
513 & hz(i-1,j,k))+ &
514 & uwrk(i,j)*(tl_hz(i ,j,k)+ &
515 & tl_hz(i-1,j,k))
516!^ ru(i,j,k,nrhs)=ru(i,j,k,nrhs)+cff
517!^
518 tl_ru(i,j,k,nrhs)=tl_ru(i,j,k,nrhs)+tl_cff
519# ifdef DIAGNOSTICS_UV
520!! DiaRU(i,j,k,nrhs,M3vvis)=DiaRU(i,j,k,nrhs,M3vvis)+cff
521!! DiaRUfrc(i,j,3,M2sstr)=DiaRUfrc(i,j,3,M2sstr)+cff
522# endif
523 END DO
524 END DO
525 DO j=jstrv,jend
526 DO i=istr,iend
527!^ cff=Vwrk(i,j)*(Hz(i,j ,k)+ &
528!^ & Hz(i,j-1,k))
529!^
530 tl_cff=tl_vwrk(i,j)*(hz(i,j ,k)+ &
531 & hz(i,j-1,k))+ &
532 & vwrk(i,j)*(tl_hz(i,j ,k)+ &
533 & tl_hz(i,j-1,k))
534!^ rv(i,j,k,nrhs)=rv(i,j,k,nrhs)+cff
535!^
536 tl_rv(i,j,k,nrhs)=tl_rv(i,j,k,nrhs)+tl_cff
537# ifdef DIAGNOSTICS_UV
538!! DiaRV(i,j,k,nrhs,M3vvis)=DiaRV(i,j,k,nrhs,M3vvis)+cff
539!! DiaRVfrc(i,j,3,M2sstr)=DiaRVfrc(i,j,3,M2sstr)+cff
540# endif
541 END DO
542 END DO
543 END DO
544!
545! Apply bottom stress as a bodyforce: determine the thickness (m)
546! of the bottom layer; then add in bottom stress as a bodyfoce.
547!
548 DO j=jstrv-1,jend
549 DO i=istru-1,iend
550 wrk(i,j)=0.0_r8
551 tl_wrk(i,j)=0.0_r8
552 END DO
553 END DO
554 DO k=1,levbfrc(ng)
555 DO j=jstrv-1,jend
556 DO i=istru-1,iend
557 wrk(i,j)=wrk(i,j)+hz(i,j,k)
558 tl_wrk(i,j)=tl_wrk(i,j)+tl_hz(i,j,k)
559 END DO
560 END DO
561 END DO
562 DO j=jstr,jend
563 DO i=istru,iend
564 cff=0.25_r8*(pm(i-1,j)+pm(i,j))* &
565 & (pn(i-1,j)+pn(i,j))
566 cff1=1.0_r8/(cff*(wrk(i-1,j)+wrk(i,j)))
567 tl_cff1=-cff1*cff1*cff*(tl_wrk(i-1,j)+tl_wrk(i,j))
568 uwrk(i,j)=bustr(i,j)*cff1
569 tl_uwrk(i,j)=tl_bustr(i,j)*cff1+ &
570 & bustr(i,j)*tl_cff1
571 END DO
572 END DO
573 DO j=jstrv,jend
574 DO i=istr,iend
575 cff=0.25_r8*(pm(i,j-1)+pm(i,j))* &
576 & (pn(i,j-1)+pn(i,j))
577 cff1=1.0_r8/(cff*(wrk(i,j-1)+wrk(i,j)))
578 tl_cff1=-cff1*cff1*cff*(tl_wrk(i,j-1)+tl_wrk(i,j))
579 vwrk(i,j)=bvstr(i,j)*cff1
580 tl_vwrk(i,j)=tl_bvstr(i,j)*cff1+ &
581 & bvstr(i,j)*tl_cff1
582 END DO
583 END DO
584 DO k=1,levbfrc(ng)
585 DO j=jstr,jend
586 DO i=istru,iend
587!^ cff=Uwrk(i,j)*(Hz(i ,j,k)+ &
588!^ & Hz(i-1,j,k))
589!^
590 tl_cff=tl_uwrk(i,j)*(hz(i ,j,k)+ &
591 & hz(i-1,j,k))+ &
592 & uwrk(i,j)*(tl_hz(i ,j,k)+ &
593 & tl_hz(i-1,j,k))
594!^ ru(i,j,k,nrhs)=ru(i,j,k,nrhs)-cff
595!^
596 tl_ru(i,j,k,nrhs)=tl_ru(i,j,k,nrhs)-tl_cff
597# ifdef DIAGNOSTICS_UV
598!! DiaRU(i,j,k,nrhs,M3vvis)=DiaRU(i,j,k,nrhs,M3vvis)-cff
599!! DiaRUfrc(i,j,3,M2bstr)=DiaRUfrc(i,j,3,M2bstr)-cff
600# endif
601 END DO
602 END DO
603 DO j=jstrv,jend
604 DO i=istr,iend
605!^ cff=Vwrk(i,j)*(Hz(i,j ,k)+ &
606!^ & Hz(i,j-1,k))
607!^
608 tl_cff=tl_vwrk(i,j)*(hz(i,j ,k)+ &
609 & hz(i,j-1,k))+ &
610 & vwrk(i,j)*(tl_hz(i,j ,k)+ &
611 & tl_hz(i,j-1,k))
612!^ rv(i,j,k,nrhs)=rv(i,j,k,nrhs)-cff
613!^
614 tl_rv(i,j,k,nrhs)=tl_rv(i,j,k,nrhs)-tl_cff
615# ifdef DIAGNOSTICS_UV
616!! DiaRV(i,j,k,nrhs,M3vvis)=DiaRV(i,j,k,nrhs,M3vvis)-cff
617!! DiaRVfrc(i,j,3,M2bstr)=DiaRVfrc(i,j,3,M2bstr)-cff
618# endif
619 END DO
620 END DO
621 END DO
622# endif
623!
624 k_loop : DO k=1,n(ng)
625
626# ifdef UV_COR
627!
628!-----------------------------------------------------------------------
629! Add in Coriolis terms.
630!-----------------------------------------------------------------------
631!
632 DO j=jstrv-1,jend
633 DO i=istru-1,iend
634 cff=0.5_r8*hz(i,j,k)*fomn(i,j)
635 tl_cff=0.5_r8*tl_hz(i,j,k)*fomn(i,j)
636!^ UFx(i,j)=cff*(v(i,j ,k,nrhs)+ &
637# ifdef WEC_MELLOR
638!^ & v_stokes(i,j ,k)+ &
639!^ & v_stokes(i,j+1,k)+ &
640# endif
641!^ & v(i,j+1,k,nrhs))
642!^
643 tl_ufx(i,j)=tl_cff*(v(i,j ,k,nrhs)+ &
644# ifdef WEC_MELLOR
645 & v_stokes(i,j ,k)+ &
646 & v_stokes(i,j+1,k)+ &
647# endif
648 & v(i,j+1,k,nrhs))+ &
649 & cff*(tl_v(i,j ,k,nrhs)+ &
650# ifdef WEC_MELLOR
651 & tl_v_stokes(i,j ,k)+ &
652 & tl_v_stokes(i,j+1,k)+ &
653# endif
654 & tl_v(i,j+1,k,nrhs))
655!^ VFe(i,j)=cff*(u(i ,j,k,nrhs)+ &
656# ifdef WEC_MELLOR
657!^ & u_stokes(i ,j,k)+ &
658!^ & u_stokes(i+1,j,k)+ &
659# endif
660!^ & u(i+1,j,k,nrhs))
661!^
662 tl_vfe(i,j)=tl_cff*(u(i ,j,k,nrhs)+ &
663# ifdef WEC_MELLOR
664 & u_stokes(i ,j,k)+ &
665 & u_stokes(i+1,j,k)+ &
666# endif
667 & u(i+1,j,k,nrhs))+ &
668 & cff*(tl_u(i ,j,k,nrhs)+ &
669# ifdef WEC_MELLOR
670 & tl_u_stokes(i ,j,k)+ &
671 & tl_u_stokes(i+1,j,k)+ &
672# endif
673 & tl_u(i+1,j,k,nrhs))
674 END DO
675 END DO
676 DO j=jstr,jend
677 DO i=istru,iend
678!^ cff1=0.5_r8*(UFx(i,j)+UFx(i-1,j))
679!^
680 tl_cff1=0.5_r8*(tl_ufx(i,j)+tl_ufx(i-1,j))
681!^ ru(i,j,k,nrhs)=ru(i,j,k,nrhs)+cff1
682!^
683 tl_ru(i,j,k,nrhs)=tl_ru(i,j,k,nrhs)+tl_cff1
684# ifdef DIAGNOSTICS_UV
685!! DiaRU(i,j,k,nrhs,M3fcor)=cff1
686# endif
687 END DO
688 END DO
689 DO j=jstrv,jend
690 DO i=istr,iend
691!^ cff1=0.5_r8*(VFe(i,j)+VFe(i,j-1))
692!^
693 tl_cff1=0.5_r8*(tl_vfe(i,j)+tl_vfe(i,j-1))
694!^ rv(i,j,k,nrhs)=rv(i,j,k,nrhs)-cff1
695!^
696 tl_rv(i,j,k,nrhs)=tl_rv(i,j,k,nrhs)-tl_cff1
697# ifdef DIAGNOSTICS_UV
698!! DiaRV(i,j,k,nrhs,M3fcor)=-cff1
699# endif
700 END DO
701 END DO
702# endif
703# if defined CURVGRID && defined UV_ADV
704!
705!-----------------------------------------------------------------------
706! Add in curvilinear transformation terms.
707!-----------------------------------------------------------------------
708!
709 DO j=jstrv-1,jend
710 DO i=istru-1,iend
711 cff1=0.5_r8*(v(i,j ,k,nrhs)+ &
712# ifdef WEC_MELLOR
713 & v_stokes(i,j ,k)+ &
714 & v_stokes(i,j+1,k)+ &
715# endif
716 & v(i,j+1,k,nrhs))
717 tl_cff1=0.5_r8*(tl_v(i,j ,k,nrhs)+ &
718# ifdef WEC_MELLOR
719 & tl_v_stokes(i,j ,k)+ &
720 & tl_v_stokes(i,j+1,k)+ &
721# endif
722 & tl_v(i,j+1,k,nrhs))
723 cff2=0.5_r8*(u(i ,j,k,nrhs)+ &
724# ifdef WEC_MELLOR
725 & u_stokes(i ,j,k)+ &
726 & u_stokes(i+1,j,k)+ &
727# endif
728 & u(i+1,j,k,nrhs))
729 tl_cff2=0.5_r8*(tl_u(i ,j,k,nrhs)+ &
730# ifdef WEC_MELLOR
731 & tl_u_stokes(i ,j,k)+ &
732 & tl_u_stokes(i+1,j,k)+ &
733# endif
734 & tl_u(i+1,j,k,nrhs))
735 cff3=cff1*dndx(i,j)
736 tl_cff3=tl_cff1*dndx(i,j)
737 cff4=cff2*dmde(i,j)
738 tl_cff4=tl_cff2*dmde(i,j)
739 cff=hz(i,j,k)*(cff3-cff4)
740 tl_cff=tl_hz(i,j,k)*(cff3-cff4)+ &
741 & hz(i,j,k)*(tl_cff3-tl_cff4)
742!^ UFx(i,j)=cff*cff1
743!^
744 tl_ufx(i,j)=tl_cff*cff1+cff*tl_cff1
745!^ VFe(i,j)=cff*cff2
746!^
747 tl_vfe(i,j)=tl_cff*cff2+cff*tl_cff2
748# if defined DIAGNOSTICS_UV
749!! cff=Hz(i,j,k)*cff4
750!! Uwrk(i,j)=-cff*cff1 ! u equation, ETA-term
751!! Vwrk(i,j)=-cff*cff2 ! v equation, ETA-term
752# endif
753 END DO
754 END DO
755 DO j=jstr,jend
756 DO i=istru,iend
757!^ cff1=0.5_r8*(UFx(i,j)+UFx(i-1,j))
758!^
759 tl_cff1=0.5_r8*(tl_ufx(i,j)+tl_ufx(i-1,j))
760!^ ru(i,j,k,nrhs)=ru(i,j,k,nrhs)+cff1
761!^
762 tl_ru(i,j,k,nrhs)=tl_ru(i,j,k,nrhs)+tl_cff1
763# ifdef DIAGNOSTICS_UV
764!! cff2=0.5_r8*(Uwrk(i,j)+Uwrk(i-1,j))
765!! DiaRU(i,j,k,nrhs,M3xadv)=cff1-cff2
766!! DiaRU(i,j,k,nrhs,M3yadv)=cff2
767!! DiaRU(i,j,k,nrhs,M3hadv)=cff1
768# endif
769 END DO
770 END DO
771 DO j=jstrv,jend
772 DO i=istr,iend
773!^ cff1=0.5_r8*(VFe(i,j)+VFe(i,j-1))
774!^
775 tl_cff1=0.5_r8*(tl_vfe(i,j)+tl_vfe(i,j-1))
776!^ rv(i,j,k,nrhs)=rv(i,j,k,nrhs)-cff1
777!^
778 tl_rv(i,j,k,nrhs)=tl_rv(i,j,k,nrhs)-tl_cff1
779# ifdef DIAGNOSTICS_UV
780!! cff2=0.5_r8*(Vwrk(i,j)+Vwrk(i,j-1))
781!! DiaRV(i,j,k,nrhs,M3xadv)=-cff1+cff2
782!! DiaRV(i,j,k,nrhs,M3yadv)=-cff2
783!! DiaRV(i,j,k,nrhs,M3hadv)=-cff1
784# endif
785 END DO
786 END DO
787# endif
788!
789!-----------------------------------------------------------------------
790! Add in nudging of 3D momentum climatology.
791!-----------------------------------------------------------------------
792!
793 IF (lnudgem3clm(ng)) THEN
794 DO j=jstr,jend
795 DO i=istru,iend
796 cff=0.25_r8*(clima(ng)%M3nudgcof(i-1,j,k)+ &
797 & clima(ng)%M3nudgcof(i ,j,k))* &
798 & om_u(i,j)*on_u(i,j)
799!^ ru(i,j,k,nrhs)=ru(i,j,k,nrhs)+ &
800!^ & cff*(Hz(i-1,j,k)+Hz(i,j,k))* &
801!^ & (CLIMA(ng)%uclm(i,j,k)- &
802!^ & u(i,j,k,nrhs))
803!^
804 tl_ru(i,j,k,nrhs)=tl_ru(i,j,k,nrhs)+ &
805 & cff*((hz(i-1,j,k)+hz(i,j,k))* &
806 & (-tl_u(i,j,k,nrhs))+ &
807 & (tl_hz(i-1,j,k)+tl_hz(i,j,k))* &
808 & (clima(ng)%uclm(i,j,k)- &
809 & u(i,j,k,nrhs)))
810 END DO
811 END DO
812 DO j=jstrv,jend
813 DO i=istr,iend
814 cff=0.25_r8*(clima(ng)%M3nudgcof(i,j-1,k)+ &
815 & clima(ng)%M3nudgcof(i,j ,k))* &
816 & om_v(i,j)*on_v(i,j)
817!^ rv(i,j,k,nrhs)=rv(i,j,k,nrhs)+ &
818!^ & cff*(Hz(i,j-1,k)+Hz(i,j,k))* &
819!^ & (CLIMA(ng)%vclm(i,j,k)- &
820!^ & v(i,j,k,nrhs))
821!^
822 tl_rv(i,j,k,nrhs)=tl_rv(i,j,k,nrhs)+ &
823 & cff*((hz(i,j-1,k)+hz(i,j,k))* &
824 & (-tl_v(i,j,k,nrhs))+ &
825 & (tl_hz(i,j-1,k)+tl_hz(i,j,k))* &
826 & (clima(ng)%vclm(i,j,k)- &
827 & v(i,j,k,nrhs)))
828 END DO
829 END DO
830 END IF
831
832# ifdef UV_ADV
833!
834!-----------------------------------------------------------------------
835! Add in horizontal advection of momentum.
836!-----------------------------------------------------------------------
837!
838! Compute diagonal [UFx,VFe] and off-diagonal [UFe,VFx] components
839! of tensor of momentum flux due to horizontal advection.
840!
841# ifdef UV_C2ADVECTION
842!
843! Second-order, centered differences advection.
844!
845 DO j=jstr,jend
846 DO i=istru-1,iend
847!^ UFx(i,j)=0.25_r8*(u(i ,j,k,nrhs)+ &
848# ifdef WEC_MELLOR
849!^ & u_stokes(i ,j,k)+ &
850!^ & u_stokes(i+1,j,k)+ &
851# endif
852!^ & u(i+1,j,k,nrhs))* &
853!^ & (Huon(i ,j,k)+ &
854!^ & Huon(i+1,j,k))
855!^
856 tl_ufx(i,j)=0.25_r8* &
857 & ((tl_u(i ,j,k,nrhs)+ &
858# ifdef WEC_MELLOR
859 & tl_u_stokes(i ,j,k)+ &
860 & tl_u_stokes(i+1,j,k)+ &
861# endif
862 & tl_u(i+1,j,k,nrhs))* &
863 & (huon(i ,j,k)+ &
864 & huon(i+1,j,k))+ &
865 & (u(i ,j,k,nrhs)+ &
866# ifdef WEC_MELLOR
867 & u_stokes(i ,j,k)+ &
868 & u_stokes(i+1,j,k)+ &
869# endif
870 & u(i+1,j,k,nrhs))* &
871 & (tl_huon(i ,j,k)+ &
872 & tl_huon(i+1,j,k)))
873 END DO
874 END DO
875 DO j=jstr,jend+1
876 DO i=istru,iend
877!^ UFe(i,j)=0.25_r8*(u(i,j-1,k,nrhs)+ &
878# ifdef WEC_MELLOR
879!^ & u_stokes(i,j-1,k)+ &
880!^ & u_stokes(i,j ,k)+ &
881# endif
882!^ & u(i,j ,k,nrhs))* &
883!^ & (Hvom(i-1,j,k)+ &
884!^ & Hvom(i ,j,k))
885!^
886 tl_ufe(i,j)=0.25_r8* &
887 & ((tl_u(i,j-1,k,nrhs)+ &
888# ifdef WEC_MELLOR
889 & tl_u_stokes(i,j-1,k)+ &
890 & tl_u_stokes(i,j ,k)+ &
891# endif
892 & tl_u(i,j ,k,nrhs))* &
893 & (hvom(i-1,j,k)+ &
894 & hvom(i ,j,k))+ &
895 & (u(i,j-1,k,nrhs)+ &
896# ifdef WEC_MELLOR
897 & u_stokes(i,j-1,k)+ &
898 & u_stokes(i,j ,k)+ &
899# endif
900 & u(i,j ,k,nrhs))* &
901 & (tl_hvom(i-1,j,k)+ &
902 & tl_hvom(i ,j,k)))
903 END DO
904 END DO
905 DO j=jstrv,jend
906 DO i=istr,iend+1
907!^ VFx(i,j)=0.25_r8*(v(i-1,j,k,nrhs)+ &
908# ifdef WEC_MELLOR
909!^ & v_stokes(i-1,j,k)+ &
910!^ & v_stokes(i ,j,k)+ &
911# endif
912!^ & v(i ,j,k,nrhs))* &
913!^ & (Huon(i,j-1,k)+ &
914!^ Huon(i,j ,k))
915!^
916 tl_vfx(i,j)=0.25_r8* &
917 & ((tl_v(i-1,j,k,nrhs)+ &
918# ifdef WEC_MELLOR
919 & tl_v_stokes(i-1,j,k)+ &
920 & tl_v_stokes(i ,j,k)+ &
921# endif
922 & tl_v(i ,j,k,nrhs))* &
923 & (huon(i,j-1,k)+ &
924 & huon(i,j ,k))+ &
925 & (v(i-1,j,k,nrhs)+ &
926# ifdef WEC_MELLOR
927 & v_stokes(i-1,j,k)+ &
928 & v_stokes(i ,j,k)+ &
929# endif
930 & v(i ,j,k,nrhs))* &
931 & (tl_huon(i,j-1,k)+ &
932 & tl_huon(i,j ,k)))
933 END DO
934 END DO
935 DO j=jstrv-1,jend
936 DO i=istr,iend
937!^ VFe(i,j)=0.25_r8*(v(i,j ,k,nrhs)+ &
938# ifdef WEC_MELLOR
939!^ & v_stokes(i,j ,k)+ &
940!^ & v_stokes(i,j+1,k)+ &
941# endif
942!^ & v(i,j+1,k,nrhs))* &
943!^ & (Hvom(i,j ,k)+ &
944!^ & Hvom(i,j+1,k))
945!^
946 tl_vfe(i,j)=0.25_r8* &
947 & ((tl_v(i,j ,k,nrhs)+ &
948# ifdef WEC_MELLOR
949 & tl_v_stokes(i,j ,k)+ &
950 & tl_v_stokes(i,j+1,k)+ &
951# endif
952 & tl_v(i,j+1,k,nrhs))* &
953 & (hvom(i,j ,k)+ &
954 & hvom(i,j+1,k))+ &
955 & (v(i,j ,k,nrhs)+ &
956# ifdef WEC_MELLOR
957 & v_stokes(i,j ,k)+ &
958 & v_stokes(i,j+1,k)+ &
959# endif
960 & v(i,j+1,k,nrhs))* &
961 & (tl_hvom(i,j ,k)+ &
962 & tl_hvom(i,j+1,k)))
963 END DO
964 END DO
965# else
966 DO j=jstr,jend
967 DO i=istrum1,iendp1
968 uxx(i,j)=u(i-1,j,k,nrhs)-2.0_r8*u(i,j,k,nrhs)+ &
969# ifdef WEC_MELLOR
970 & u_stokes(i-1,j,k)-2.0_r8*u_stokes(i,j,k)+ &
971 & u_stokes(i+1,j,k)+ &
972# endif
973 & u(i+1,j,k,nrhs)
974 tl_uxx(i,j)=tl_u(i-1,j,k,nrhs)-2.0_r8*tl_u(i,j,k,nrhs)+ &
975# ifdef WEC_MELLOR
976 & tl_u_stokes(i-1,j,k)-2.0_r8*tl_u_stokes(i,j,k)+ &
977 & tl_u_stokes(i+1,j,k)+ &
978# endif
979 & tl_u(i+1,j,k,nrhs)
980 huxx(i,j)=huon(i-1,j,k)-2.0_r8*huon(i,j,k)+huon(i+1,j,k)
981 tl_huxx(i,j)=tl_huon(i-1,j,k)-2.0_r8*tl_huon(i,j,k)+ &
982 & tl_huon(i+1,j,k)
983 END DO
984 END DO
985 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
986 IF (domain(ng)%Western_Edge(tile)) THEN
987 DO j=jstr,jend
988 uxx(istr,j)=uxx(istr+1,j)
989 tl_uxx(istr,j)=tl_uxx(istr+1,j)
990 huxx(istr,j)=huxx(istr+1,j)
991 tl_huxx(istr,j)=tl_huxx(istr+1,j)
992 END DO
993 END IF
994 END IF
995 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
996 IF (domain(ng)%Eastern_Edge(tile)) THEN
997 DO j=jstr,jend
998 uxx(iend+1,j)=uxx(iend,j)
999 tl_uxx(iend+1,j)=tl_uxx(iend,j)
1000 huxx(iend+1,j)=huxx(iend,j)
1001 tl_huxx(iend+1,j)=tl_huxx(iend,j)
1002 END DO
1003 END IF
1004 END IF
1005# ifdef UV_C4ADVECTION
1006!
1007! Fourth-order, centered differences u-momentum horizontal advection.
1008!
1009 cff=1.0_r8/6.0_r8
1010 DO j=jstr,jend
1011 DO i=istru-1,iend
1012!^ UFx(i,j)=0.25_r8*(u(i ,j,k,nrhs)+ &
1013# ifdef WEC_MELLOR
1014!^ & u_stokes(i ,j,k)+ &
1015!^ & u_stokes(i+1,j,k)+ &
1016# endif
1017!^ & u(i+1,j,k,nrhs)- &
1018!^ & cff*(uxx (i ,j)+ &
1019!^ & uxx (i+1,j)))* &
1020!^ & (Huon(i ,j,k)+ &
1021!^ & Huon(i+1,j,k)- &
1022!^ & cff*(Huxx(i ,j)+ &
1023!^ & Huxx(i+1,j)))
1024!^
1025 tl_ufx(i,j)=0.25_r8*((tl_u(i ,j,k,nrhs)+ &
1026# ifdef WEC_MELLOR
1027 & tl_u_stokes(i ,j,k)+ &
1028 & tl_u_stokes(i+1,j,k)+ &
1029# endif
1030 & tl_u(i+1,j,k,nrhs)- &
1031 & cff*(tl_uxx(i ,j)+ &
1032 & tl_uxx(i+1,j)))* &
1033 & (huon(i ,j,k)+ &
1034 & huon(i+1,j,k)- &
1035 & cff*(huxx(i ,j)+ &
1036 & huxx(i+1,j)))+ &
1037 & (u(i ,j,k,nrhs)+ &
1038# ifdef WEC_MELLOR
1039 & u_stokes(i ,j,k)+ &
1040 & u_stokes(i+1,j,k)+ &
1041# endif
1042 & u(i+1,j,k,nrhs)- &
1043 & cff*(uxx(i ,j)+ &
1044 & uxx(i+1,j)))* &
1045 & (tl_huon(i ,j,k)+ &
1046 & tl_huon(i+1,j,k)- &
1047 & cff*(tl_huxx(i ,j)+ &
1048 & tl_huxx(i+1,j))))
1049 END DO
1050 END DO
1051# else
1052!
1053! Third-order, upstream bias u-momentum advection with velocity
1054! dependent hyperdiffusion.
1055!
1056 DO j=jstr,jend
1057 DO i=istru-1,iend
1058 cff1=u(i ,j,k,nrhs)+ &
1059# ifdef WEC_MELLOR
1060 & u_stokes(i ,j,k)+ &
1061 & u_stokes(i+1,j,k)+ &
1062# endif
1063 & u(i+1,j,k,nrhs)
1064 tl_cff1=tl_u(i ,j,k,nrhs)+ &
1065# ifdef WEC_MELLOR
1066 & tl_u_stokes(i ,j,k)+ &
1067 & tl_u_stokes(i+1,j,k)+ &
1068# endif
1069 & tl_u(i+1,j,k,nrhs)
1070 IF (cff1.gt.0.0_r8) THEN
1071 cff=uxx(i,j)
1072 tl_cff=tl_uxx(i,j)
1073 ELSE
1074 cff=uxx(i+1,j)
1075 tl_cff=tl_uxx(i+1,j)
1076 END IF
1077!^ UFx(i,j)=0.25_r8*(cff1+Gadv*cff)* &
1078!^ & (Huon(i ,j,k)+ &
1079!^ & Huon(i+1,j,k)+ &
1080!^ & Gadv*0.5_r8*(Huxx(i ,j)+ &
1081!^ & Huxx(i+1,j)))
1082!^
1083 tl_ufx(i,j)=0.25_r8* &
1084 & ((tl_cff1+gadv*tl_cff)* &
1085 & (huon(i ,j,k)+ &
1086 & huon(i+1,j,k)+ &
1087 & gadv*0.5_r8*(huxx(i ,j)+ &
1088 & huxx(i+1,j)))+ &
1089 & (cff1+gadv*cff)* &
1090 & (tl_huon(i ,j,k)+ &
1091 & tl_huon(i+1,j,k)+ &
1092 & gadv*0.5_r8*(tl_huxx(i ,j)+ &
1093 & tl_huxx(i+1,j))))
1094 END DO
1095 END DO
1096# endif
1097 DO j=jstrm1,jendp1
1098 DO i=istru,iend
1099 uee(i,j)=u(i,j-1,k,nrhs)-2.0_r8*u(i,j,k,nrhs)+ &
1100# ifdef WEC_MELLOR
1101 & u_stokes(i,j-1,k)-2.0_r8*u_stokes(i,j,k)+ &
1102 & u_stokes(i,j+1,k)+ &
1103# endif
1104 & u(i,j+1,k,nrhs)
1105 tl_uee(i,j)=tl_u(i,j-1,k,nrhs)-2.0_r8*tl_u(i,j,k,nrhs)+ &
1106# ifdef WEC_MELLOR
1107 & tl_u_stokes(i,j-1,k)-2.0_r8*tl_u_stokes(i,j,k)+ &
1108 & tl_u_stokes(i,j+1,k)+ &
1109# endif
1110 & tl_u(i,j+1,k,nrhs)
1111 END DO
1112 END DO
1113 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
1114 IF (domain(ng)%Southern_Edge(tile)) THEN
1115 DO i=istru,iend
1116 uee(i,jstr-1)=uee(i,jstr)
1117 tl_uee(i,jstr-1)=tl_uee(i,jstr)
1118 END DO
1119 END IF
1120 END IF
1121 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
1122 IF (domain(ng)%Northern_Edge(tile)) THEN
1123 DO i=istru,iend
1124 uee(i,jend+1)=uee(i,jend)
1125 tl_uee(i,jend+1)=tl_uee(i,jend)
1126 END DO
1127 END IF
1128 END IF
1129 DO j=jstr,jend+1
1130 DO i=istru-1,iend
1131 hvxx(i,j)=hvom(i-1,j,k)-2.0_r8*hvom(i,j,k)+hvom(i+1,j,k)
1132 tl_hvxx(i,j)=tl_hvom(i-1,j,k)-2.0_r8*tl_hvom(i,j,k)+ &
1133 & tl_hvom(i+1,j,k)
1134 END DO
1135 END DO
1136# ifdef UV_C4ADVECTION
1137 cff=1.0_r8/6.0_r8
1138 DO j=jstr,jend+1
1139 DO i=istru,iend
1140!^ UFe(i,j)=0.25_r8*(u(i,j ,k,nrhs)+ &
1141# ifdef WEC_MELLOR
1142!^ & u_stokes(i,j ,k)+ &
1143!^ & u_stokes(i,j-1,k)+ &
1144# endif
1145!^ & u(i,j-1,k,nrhs)- &
1146!^ & cff*(uee (i,j )+ &
1147!^ & uee (i,j-1)))* &
1148!^ & (Hvom(i ,j,k)+ &
1149!^ & Hvom(i-1,j,k)- &
1150!^ & cff*(Hvxx(i ,j)+ &
1151!^ & Hvxx(i-1,j)))
1152!^
1153 tl_ufe(i,j)=0.25_r8*((tl_u(i,j ,k,nrhs)+ &
1154# ifdef WEC_MELLOR
1155 & tl_u_stokes(i,j ,k)+ &
1156 & tl_u_stokes(i,j-1,k)+ &
1157# endif
1158 & tl_u(i,j-1,k,nrhs)- &
1159 & cff*(tl_uee(i,j )+ &
1160 & tl_uee(i,j-1)))* &
1161 & (hvom(i ,j,k)+ &
1162 & hvom(i-1,j,k)- &
1163 & cff*(hvxx(i ,j)+ &
1164 & hvxx(i-1,j)))+ &
1165 & (u(i,j ,k,nrhs)+ &
1166# ifdef WEC_MELLOR
1167 & u_stokes(i,j ,k)+ &
1168 & u_stokes(i,j-1,k)+ &
1169# endif
1170 & u(i,j-1,k,nrhs)- &
1171 & cff*(uee(i,j )+ &
1172 & uee(i,j-1)))* &
1173 & (tl_hvom(i ,j,k)+ &
1174 & tl_hvom(i-1,j,k)- &
1175 & cff*(tl_hvxx(i ,j)+ &
1176 & tl_hvxx(i-1,j))))
1177 END DO
1178 END DO
1179# else
1180 DO j=jstr,jend+1
1181 DO i=istru,iend
1182 cff1=u(i,j ,k,nrhs)+ &
1183# ifdef WEC_MELLOR
1184 & u_stokes(i,j ,k)+ &
1185 & u_stokes(i,j-1,k)+ &
1186# endif
1187 & u(i,j-1,k,nrhs)
1188 tl_cff1=tl_u(i,j,k,nrhs)+ &
1189# ifdef WEC_MELLOR
1190 & tl_u_stokes(i,j ,k)+ &
1191 & tl_u_stokes(i,j-1,k)+ &
1192# endif
1193 & tl_u(i,j-1,k,nrhs)
1194 cff2=hvom(i,j,k)+hvom(i-1,j,k)
1195 tl_cff2=tl_hvom(i,j,k)+tl_hvom(i-1,j,k)
1196 IF (cff2.gt.0.0_r8) THEN
1197 cff=uee(i,j-1)
1198 tl_cff=tl_uee(i,j-1)
1199 ELSE
1200 cff=uee(i,j)
1201 tl_cff=tl_uee(i,j)
1202 END IF
1203!^ UFe(i,j)=0.25_r8*(cff1+Gadv*cff)* &
1204!^ & (cff2+Gadv*0.5_r8*(Hvxx(i ,j)+ &
1205!^ & Hvxx(i-1,j)))
1206!^
1207 tl_ufe(i,j)=0.25_r8* &
1208 & ((tl_cff1+gadv*tl_cff)* &
1209 & (cff2+gadv*0.5_r8*(hvxx(i ,j)+ &
1210 & hvxx(i-1,j)))+ &
1211 & (cff1+gadv*cff)* &
1212 & (tl_cff2+gadv*0.5_r8*(tl_hvxx(i ,j)+ &
1213 & tl_hvxx(i-1,j))))
1214 END DO
1215 END DO
1216# endif
1217 DO j=jstrv,jend
1218 DO i=istrm1,iendp1
1219 vxx(i,j)=v(i-1,j,k,nrhs)-2.0_r8*v(i,j,k,nrhs)+ &
1220# ifdef WEC_MELLOR
1221 & v_stokes(i-1,j,k)-2.0_r8*v_stokes(i,j,k)+ &
1222 & v_stokes(i+1,j,k)+ &
1223# endif
1224 & v(i+1,j,k,nrhs)
1225 tl_vxx(i,j)=tl_v(i-1,j,k,nrhs)-2.0_r8*tl_v(i,j,k,nrhs)+ &
1226# ifdef WEC_MELLOR
1227 & tl_v_stokes(i-1,j,k)-2.0_r8*tl_v_stokes(i,j,k)+ &
1228 & tl_v_stokes(i+1,j,k)+ &
1229# endif
1230 & tl_v(i+1,j,k,nrhs)
1231 END DO
1232 END DO
1233 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
1234 IF (domain(ng)%Western_Edge(tile)) THEN
1235 DO j=jstrv,jend
1236 vxx(istr-1,j)=vxx(istr,j)
1237 tl_vxx(istr-1,j)=tl_vxx(istr,j)
1238 END DO
1239 END IF
1240 END IF
1241 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
1242 IF (domain(ng)%Eastern_Edge(tile)) THEN
1243 DO j=jstrv,jend
1244 vxx(iend+1,j)=vxx(iend,j)
1245 tl_vxx(iend+1,j)=tl_vxx(iend,j)
1246 END DO
1247 END IF
1248 END IF
1249 DO j=jstrv-1,jend
1250 DO i=istr,iend+1
1251 huee(i,j)=huon(i,j-1,k)-2.0_r8*huon(i,j,k)+huon(i,j+1,k)
1252 tl_huee(i,j)=tl_huon(i,j-1,k)-2.0_r8*tl_huon(i,j,k)+ &
1253 & tl_huon(i,j+1,k)
1254 END DO
1255 END DO
1256# ifdef UV_C4ADVECTION
1257!
1258! Fourth-order, centered differences v-momentum horizontal advection.
1259!
1260 cff=1.0_r8/6.0_r8
1261 DO j=jstrv,jend
1262 DO i=istr,iend+1
1263!^ VFx(i,j)=0.25_r8*(v(i ,j,k,nrhs)+ &
1264# ifdef WEC_MELLOR
1265!^ & v_stokes(i ,j,k)+ &
1266!^ & v_stokes(i-1,j,k)+ &
1267# endif
1268!^ & v(i-1,j,k,nrhs)- &
1269!^ & cff*(vxx (i ,j)+ &
1270!^ & vxx (i-1,j)))* &
1271!^ & (Huon(i,j ,k)+ &
1272!^ & Huon(i,j-1,k)- &
1273!^ & cff*(Huee(i,j )+ &
1274!^ & Huee(i,j-1)))
1275!^
1276 tl_vfx(i,j)=0.25_r8*((tl_v(i ,j,k,nrhs)+ &
1277# ifdef WEC_MELLOR
1278 & tl_v_stokes(i ,j,k)+ &
1279 & tl_v_stokes(i-1,j,k)+ &
1280# endif
1281 & tl_v(i-1,j,k,nrhs)- &
1282 & cff*(tl_vxx(i ,j)+ &
1283 & tl_vxx(i-1,j)))* &
1284 & (huon(i,j ,k)+ &
1285 & huon(i,j-1,k)- &
1286 & cff*(huee(i,j )+ &
1287 & huee(i,j-1)))+ &
1288 & (v(i ,j,k,nrhs)+ &
1289# ifdef WEC_MELLOR
1290 & v_stokes(i ,j,k)+ &
1291 & v_stokes(i-1,j,k)+ &
1292# endif
1293 & v(i-1,j,k,nrhs)- &
1294 & cff*(vxx(i ,j)+ &
1295 & vxx(i-1,j)))* &
1296 & (tl_huon(i,j ,k)+ &
1297 & tl_huon(i,j-1,k)- &
1298 & cff*(tl_huee(i,j )+ &
1299 & tl_huee(i,j-1))))
1300 END DO
1301 END DO
1302# else
1303!
1304! Third-order, upstream bias v-momentum advection with velocity
1305! dependent hyperdiffusion.
1306!
1307 DO j=jstrv,jend
1308 DO i=istr,iend+1
1309 cff1=v(i ,j,k,nrhs)+ &
1310# ifdef WEC_MELLOR
1311 & v_stokes(i ,j,k)+ &
1312 & v_stokes(i-1,j,k)+ &
1313# endif
1314 & v(i-1,j,k,nrhs)
1315 tl_cff1=tl_v(i ,j,k,nrhs)+ &
1316# ifdef WEC_MELLOR
1317 & tl_v_stokes(i ,j,k)+ &
1318 & tl_v_stokes(i-1,j,k)+ &
1319# endif
1320 & tl_v(i-1,j,k,nrhs)
1321 cff2=huon(i,j,k)+huon(i,j-1,k)
1322 tl_cff2=tl_huon(i,j,k)+tl_huon(i,j-1,k)
1323 IF (cff2.gt.0.0_r8) THEN
1324 cff=vxx(i-1,j)
1325 tl_cff=tl_vxx(i-1,j)
1326 ELSE
1327 cff=vxx(i,j)
1328 tl_cff=tl_vxx(i,j)
1329 END IF
1330!^ VFx(i,j)=0.25_r8*(cff1+Gadv*cff)* &
1331!^ & (cff2+Gadv*0.5_r8*(Huee(i,j )+ &
1332!^ & Huee(i,j-1)))
1333!^
1334 tl_vfx(i,j)=0.25_r8* &
1335 & ((tl_cff1+gadv*tl_cff)* &
1336 & (cff2+gadv*0.5_r8*(huee(i,j )+ &
1337 & huee(i,j-1)))+ &
1338 & (cff1+gadv*cff)* &
1339 & (tl_cff2+gadv*0.5_r8*(tl_huee(i,j )+ &
1340 & tl_huee(i,j-1))))
1341 END DO
1342 END DO
1343# endif
1344 DO j=jstrvm1,jendp1
1345 DO i=istr,iend
1346 vee(i,j)=v(i,j-1,k,nrhs)-2.0_r8*v(i,j,k,nrhs)+ &
1347# ifdef WEC_MELLOR
1348 & v_stokes(i,j-1,k)-2.0_r8*v_stokes(i,j,k)+ &
1349 & v_stokes(i,j+1,k)+ &
1350# endif
1351 & v(i,j+1,k,nrhs)
1352 tl_vee(i,j)=tl_v(i,j-1,k,nrhs)-2.0_r8*tl_v(i,j,k,nrhs)+ &
1353# ifdef WEC_MELLOR
1354 & tl_v_stokes(i,j-1,k)-2.0_r8*tl_v_stokes(i,j,k)+ &
1355 & tl_v_stokes(i,j+1,k)+ &
1356# endif
1357 & tl_v(i,j+1,k,nrhs)
1358 hvee(i,j)=hvom(i,j-1,k)-2.0_r8*hvom(i,j,k)+hvom(i,j+1,k)
1359 tl_hvee(i,j)=tl_hvom(i,j-1,k)-2.0_r8*tl_hvom(i,j,k)+ &
1360 & tl_hvom(i,j+1,k)
1361 END DO
1362 END DO
1363 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
1364 IF (domain(ng)%Southern_Edge(tile)) THEN
1365 DO i=istr,iend
1366 vee(i,jstr)=vee(i,jstr+1)
1367 tl_vee(i,jstr)=tl_vee(i,jstr+1)
1368 hvee(i,jstr)=hvee(i,jstr+1)
1369 tl_hvee(i,jstr)=tl_hvee(i,jstr+1)
1370 END DO
1371 END IF
1372 END IF
1373 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
1374 IF (domain(ng)%Northern_Edge(tile)) THEN
1375 DO i=istr,iend
1376 vee(i,jend+1)=vee(i,jend)
1377 tl_vee(i,jend+1)=tl_vee(i,jend)
1378 hvee(i,jend+1)=hvee(i,jend)
1379 tl_hvee(i,jend+1)=tl_hvee(i,jend)
1380 END DO
1381 END IF
1382 END IF
1383# ifdef UV_C4ADVECTION
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*(v(i,j,k,nrhs)+ &
1388# ifdef WEC_MELLOR
1389!^ & v_stokes(i,j ,k)+ &
1390!^ & v_stokes(i,j+1,k)+ &
1391# endif
1392!^ & v(i,j+1,k,nrhs)- &
1393!^ & cff*(vee (i,j )+ &
1394!^ & vee (i,j+1)))* &
1395!^ & (Hvom(i,j ,k)+ &
1396!^ & Hvom(i,j+1,k)- &
1397!^ & cff*(Hvee(i,j )+ &
1398!^ & Hvee(i,j+1)))
1399!^
1400 tl_vfe(i,j)=0.25_r8*((tl_v(i,j ,k,nrhs)+ &
1401# ifdef WEC_MELLOR
1402 & tl_v_stokes(i,j ,k)+ &
1403 & tl_v_stokes(i,j+1,k)+ &
1404# endif
1405 & tl_v(i,j+1,k,nrhs)- &
1406 & cff*(tl_vee(i,j )+ &
1407 & tl_vee(i,j+1)))* &
1408 & (hvom(i,j ,k)+ &
1409 & hvom(i,j+1,k)- &
1410 & cff*(hvee(i,j )+ &
1411 & hvee(i,j+1)))+ &
1412 & (v(i,j ,k,nrhs)+ &
1413# ifdef WEC_MELLOR
1414 & v_stokes(i,j ,k)+ &
1415 & v_stokes(i,j+1,k)+ &
1416# endif
1417 & v(i,j+1,k,nrhs)- &
1418 & cff*(vee(i,j )+ &
1419 & vee(i,j+1)))* &
1420 & (tl_hvom(i,j ,k)+ &
1421 & tl_hvom(i,j+1,k)- &
1422 & cff*(tl_hvee(i,j )+ &
1423 & tl_hvee(i,j+1))))
1424 END DO
1425 END DO
1426# else
1427 DO j=jstrv-1,jend
1428 DO i=istr,iend
1429 cff1=v(i,j ,k,nrhs)+ &
1430# ifdef WEC_MELLOR
1431 & v_stokes(i,j ,k)+ &
1432 & v_stokes(i,j+1,k)+ &
1433# endif
1434 & v(i,j+1,k,nrhs)
1435 tl_cff1=tl_v(i,j ,k,nrhs)+ &
1436# ifdef WEC_MELLOR
1437 & tl_v_stokes(i,j ,k)+ &
1438 & tl_v_stokes(i,j+1,k)+ &
1439# endif
1440 & tl_v(i,j+1,k,nrhs)
1441 IF (cff1.gt.0.0_r8) THEN
1442 cff=vee(i,j)
1443 tl_cff=tl_vee(i,j)
1444 ELSE
1445 cff=vee(i,j+1)
1446 tl_cff=tl_vee(i,j+1)
1447 END IF
1448!^ VFe(i,j)=0.25_r8*(cff1+Gadv*cff)* &
1449!^ & (Hvom(i,j ,k)+ &
1450!^ & Hvom(i,j+1,k)+ &
1451!^ & Gadv*0.5_r8*(Hvee(i,j )+ &
1452!^ & Hvee(i,j+1)))
1453!^
1454 tl_vfe(i,j)=0.25_r8* &
1455 & ((tl_cff1+gadv*tl_cff)* &
1456 & (hvom(i,j ,k)+ &
1457 & hvom(i,j+1,k)+ &
1458 & gadv*0.5_r8*(hvee(i,j )+ &
1459 & hvee(i,j+1)))+ &
1460 & (cff1+gadv*cff)* &
1461 & (tl_hvom(i,j ,k)+ &
1462 & tl_hvom(i,j+1,k)+ &
1463 & gadv*0.5_r8*(tl_hvee(i,j )+ &
1464 & tl_hvee(i,j+1))))
1465 END DO
1466 END DO
1467# endif
1468# endif
1469!
1470! Add in horizontal advection.
1471!
1472 DO j=jstr,jend
1473 DO i=istru,iend
1474!^ cff1=UFx(i,j)-UFx(i-1,j)
1475!^
1476 tl_cff1=tl_ufx(i,j)-tl_ufx(i-1,j)
1477!^ cff2=UFe(i,j+1)-UFe(i,j)
1478!^
1479 tl_cff2=tl_ufe(i,j+1)-tl_ufe(i,j)
1480!^ cff=cff1+cff2
1481!^
1482 tl_cff=tl_cff1+tl_cff2
1483!^ ru(i,j,k,nrhs)=ru(i,j,k,nrhs)-cff
1484!^
1485 tl_ru(i,j,k,nrhs)=tl_ru(i,j,k,nrhs)-tl_cff
1486# ifdef DIAGNOSTICS_UV
1487# ifdef CURVGRID
1488!! DiaRU(i,j,k,nrhs,M3xadv)=DiaRU(i,j,k,nrhs,M3xadv)-cff1
1489!! DiaRU(i,j,k,nrhs,M3yadv)=DiaRU(i,j,k,nrhs,M3yadv)-cff2
1490!! DiaRU(i,j,k,nrhs,M3hadv)=DiaRU(i,j,k,nrhs,M3hadv)-cff
1491# else
1492!! DiaRU(i,j,k,nrhs,M3xadv)=-cff1
1493!! DiaRU(i,j,k,nrhs,M3yadv)=-cff2
1494!! DiaRU(i,j,k,nrhs,M3hadv)=-cff
1495# endif
1496# endif
1497 END DO
1498 END DO
1499 DO j=jstrv,jend
1500 DO i=istr,iend
1501!^ cff1=VFx(i+1,j)-VFx(i,j)
1502!^
1503 tl_cff1=tl_vfx(i+1,j)-tl_vfx(i,j)
1504!^ cff2=VFe(i,j)-VFe(i,j-1)
1505!^
1506 tl_cff2=tl_vfe(i,j)-tl_vfe(i,j-1)
1507!^ cff=cff1+cff2
1508!^
1509 tl_cff=tl_cff1+tl_cff2
1510!^ rv(i,j,k,nrhs)=rv(i,j,k,nrhs)-cff
1511!^
1512 tl_rv(i,j,k,nrhs)=tl_rv(i,j,k,nrhs)-tl_cff
1513# ifdef DIAGNOSTICS_UV
1514# ifdef CURVGRID
1515!! DiaRV(i,j,k,nrhs,M3xadv)=DiaRV(i,j,k,nrhs,M3xadv)-cff1
1516!! DiaRV(i,j,k,nrhs,M3yadv)=DiaRV(i,j,k,nrhs,M3yadv)-cff2
1517!! DiaRV(i,j,k,nrhs,M3hadv)=DiaRV(i,j,k,nrhs,M3hadv)-cff
1518# else
1519!! DiaRV(i,j,k,nrhs,M3xadv)=-cff1
1520!! DiaRV(i,j,k,nrhs,M3yadv)=-cff2
1521!! DiaRV(i,j,k,nrhs,M3hadv)=-cff
1522# endif
1523# endif
1524 END DO
1525 END DO
1526# endif
1527# ifdef WEC_MELLOR
1528!
1529!-----------------------------------------------------------------------
1530! Add in radiation stress terms. Convert stresses to m4/s2.
1531!-----------------------------------------------------------------------
1532!
1533 DO j=jstr,jend
1534 DO i=istru,iend
1535!^ ru(i,j,k,nrhs)=ru(i,j,k,nrhs)- &
1536!^ & rustr3d(i,j,k)*om_u(i,j)*on_u(i,j)- &
1537!^ & rulag3d(i,j,k)
1538!^
1539 tl_ru(i,j,k,nrhs)=tl_ru(i,j,k,nrhs)- &
1540 & tl_rustr3d(i,j,k)*om_u(i,j)*on_u(i,j)- &
1541 & tl_rulag3d(i,j,k)
1542 END DO
1543 END DO
1544 DO j=jstrv,jend
1545 DO i=istr,iend
1546!^ rv(i,j,k,nrhs)=rv(i,j,k,nrhs)- &
1547!^ & rvstr3d(i,j,k)*om_v(i,j)*on_v(i,j)- &
1548!^ & rvlag3d(i,j,k)
1549!^
1550 tl_rv(i,j,k,nrhs)=tl_rv(i,j,k,nrhs)- &
1551 & tl_rvstr3d(i,j,k)*om_v(i,j)*on_v(i,j)- &
1552 & tl_rvlag3d(i,j,k)
1553 END DO
1554 END DO
1555# endif
1556
1557 END DO k_loop
1558!
1559 j_loop : DO j=jstr,jend
1560# ifdef UV_ADV
1561!
1562!-----------------------------------------------------------------------
1563! Add in vertical advection.
1564!-----------------------------------------------------------------------
1565!
1566# ifdef UV_SADVECTION
1567!
1568! Apply spline code to BASIC STATE u-momentum which should be in
1569! units of m/s. CF will be used by the tangent linear spline code.
1570!
1571 cff1=9.0_r8/16.0_r8
1572 cff2=1.0_r8/16.0_r8
1573 DO k=1,n(ng)
1574 DO i=istru,iend
1575 dc(i,k)=cff1*(hz(i ,j,k)+ &
1576 & hz(i-1,j,k))- &
1577 & cff2*(hz(i+1,j,k)+ &
1578 & hz(i-2,j,k))
1579 END DO
1580 END DO
1581 DO i=istru,iend
1582 fc(i,0)=0.0_r8
1583 cf(i,0)=0.0_r8
1584 END DO
1585 DO k=1,n(ng)-1
1586 DO i=istru,iend
1587 cff=1.0_r8/(2.0_r8*dc(i,k+1)+dc(i,k)*(2.0_r8-fc(i,k-1)))
1588 fc(i,k)=cff*dc(i,k+1)
1589 cf(i,k)=cff*(6.0_r8*(u(i,j,k+1,nrhs)- &
1590# ifdef WEC_MELLOR
1591 & u_stokes(i,j,k )+ &
1592 & u_stokes(i,j,k+1)- &
1593# endif
1594 & u(i,j,k ,nrhs))- &
1595 & dc(i,k)*cf(i,k-1))
1596 END DO
1597 END DO
1598 DO i=istru,iend
1599 cf(i,n(ng))=0.0_r8
1600 END DO
1601 DO k=n(ng)-1,1,-1
1602 DO i=istru,iend
1603 cf(i,k)=cf(i,k)-fc(i,k)*cf(i,k+1)
1604 END DO
1605 END DO
1606!
1607! Construct tangent linear conservative parabolic splines for the
1608! vertical derivatives "tl_CF" of u-momentum.
1609!
1610 cff1=9.0_r8/16.0_r8
1611 cff2=1.0_r8/16.0_r8
1612 DO k=1,n(ng)
1613 DO i=istru,iend
1614 dc(i,k)=cff1*(hz(i ,j,k)+ &
1615 & hz(i-1,j,k))- &
1616 & cff2*(hz(i+1,j,k)+ &
1617 & hz(i-2,j,k))
1618 tl_dc(i,k)=cff1*(tl_hz(i ,j,k)+ &
1619 & tl_hz(i-1,j,k))- &
1620 & cff2*(tl_hz(i+1,j,k)+ &
1621 & tl_hz(i-2,j,k))
1622 END DO
1623 END DO
1624 DO i=istru,iend
1625 fc(i,0)=0.0_r8
1626 tl_cf(i,0)=0.0_r8
1627 END DO
1628 DO k=1,n(ng)-1
1629 DO i=istru,iend
1630 cff=1.0_r8/(2.0_r8*dc(i,k+1)+dc(i,k)*(2.0_r8-fc(i,k-1)))
1631 fc(i,k)=cff*dc(i,k+1)
1632 tl_cf(i,k)=cff*(6.0_r8*(tl_u(i,j,k+1,nrhs)- &
1633# ifdef WEC_MELLOR
1634 & tl_u_stokes(i,j,k )+ &
1635 & tl_u_stokes(i,j,k+1)- &
1636# endif
1637 & tl_u(i,j,k ,nrhs))- &
1638 & (tl_dc(i,k)*cf(i,k-1)+ &
1639 & 2.0_r8*(tl_dc(i,k)+tl_dc(i,k+1))*cf(i,k)+ &
1640 & tl_dc(i,k+1)*cf(i,k+1))- &
1641 & dc(i,k)*tl_cf(i,k-1))
1642 END DO
1643 END DO
1644 DO i=istru,iend
1645 tl_cf(i,n(ng))=0.0_r8
1646 END DO
1647 DO k=n(ng)-1,1,-1
1648 DO i=istru,iend
1649 tl_cf(i,k)=tl_cf(i,k)-fc(i,k)*tl_cf(i,k+1)
1650 END DO
1651 END DO
1652!
1653! Compute spline-interpolated, vertical advective u-momentum flux.
1654!
1655 cff3=1.0_r8/3.0_r8
1656 cff4=1.0_r8/6.0_r8
1657 DO k=1,n(ng)-1
1658 DO i=istru,iend
1659!^ FC(i,k)=(cff1*(W(i ,j,k)+ &
1660!^ & W(i-1,j,k))- &
1661!^ & cff2*(W(i+1,j,k)+ &
1662!^ & W(i-2,j,k)))* &
1663!^ & (u(i,j,k,nrhs)+ &
1664# ifdef WEC_MELLOR
1665!^ & u_stokes(i,j,k)+ &
1666# endif
1667!^ & DC(i,k)*(cff3*CF(i,k )+ &
1668!^ & cff4*CF(i,k-1)))
1669!^
1670 tl_fc(i,k)=(cff1*(tl_w(i ,j,k)+ &
1671 & tl_w(i-1,j,k))- &
1672 & cff2*(tl_w(i+1,j,k)+ &
1673 & tl_w(i-2,j,k)))* &
1674 & (u(i,j,k,nrhs)+ &
1675# ifdef WEC_MELLOR
1676 & u_stokes(i,j,k)+ &
1677# endif
1678 & dc(i,k)*(cff3*cf(i,k )+ &
1679 & cff4*cf(i,k-1)))+ &
1680 & (cff1*(w(i ,j,k)+ &
1681 & w(i-1,j,k))- &
1682 & cff2*(w(i+1,j,k)+ &
1683 & w(i-2,j,k)))* &
1684 & (tl_u(i,j,k,nrhs)+ &
1685# ifdef WEC_MELLOR
1686 & tl_u_stokes(i,j,k)+ &
1687# endif
1688 & dc(i,k)*(cff3*tl_cf(i,k )+ &
1689 & cff4*tl_cf(i,k-1))+ &
1690 & tl_dc(i,k)*(cff3*cf(i,k )+ &
1691 & cff4*cf(i,k-1)))
1692 END DO
1693 END DO
1694 DO i=istru,iend
1695!^ FC(i,N(ng))=0.0_r8
1696!^
1697 tl_fc(i,n(ng))=0.0_r8
1698!^ FC(i,0)=0.0_r8
1699!^
1700 tl_fc(i,0)=0.0_r8
1701 END DO
1702# elif defined UV_C2ADVECTION
1703 DO k=1,n(ng)-1
1704 DO i=istru,iend
1705!^ FC(i,k)=0.25_r8*(u(i,j,k ,nrhs)+ &
1706# ifdef WEC_MELLOR
1707!^ & u_stokes(i,j,k )+ &
1708!^ & u_stokes(i,j,k+1)+ &
1709# endif
1710!^ & u(i,j,k+1,nrhs))* &
1711!^ & (W(i ,j,k)+ &
1712!^ & W(i-1,j,k))
1713!^
1714 tl_fc(i,k)=0.25_r8*((tl_u(i,j,k ,nrhs)+ &
1715# ifdef WEC_MELLOR
1716 & tl_u_stokes(i,j,k )+ &
1717 & tl_u_stokes(i,j,k+1)+ &
1718# endif
1719 & tl_u(i,j,k+1,nrhs))* &
1720 & (w(i ,j,k)+ &
1721 & w(i-1,j,k))+ &
1722 & (u(i,j,k ,nrhs)+ &
1723# ifdef WEC_MELLOR
1724 & u_stokes(i,j,k )+ &
1725 & u_stokes(i,j,k+1)+ &
1726# endif
1727 & u(i,j,k+1,nrhs))* &
1728 & (tl_w(i ,j,k)+ &
1729 & tl_w(i-1,j,k)))
1730 END DO
1731 END DO
1732 DO i=istru,iend
1733!^ FC(i,0)=0.0_r8
1734!^
1735 tl_fc(i,0)=0.0_r8
1736!^ FC(i,N(ng))=0.0_r8
1737!^
1738 tl_fc(i,n(ng))=0.0_r8
1739 END DO
1740# elif defined UV_C4ADVECTION
1741 cff1=9.0_r8/32.0_r8
1742 cff2=1.0_r8/32.0_r8
1743 DO k=2,n(ng)-2
1744 DO i=istru,iend
1745!^ FC(i,k)=(cff1*(u(i,j,k ,nrhs)+ &
1746# ifdef WEC_MELLOR
1747!^ & u_stokes(i,j,k )+ &
1748!^ & u_stokes(i,j,k+1)+ &
1749# endif
1750!^ & u(i,j,k+1,nrhs))- &
1751!^ & cff2*(u(i,j,k-1,nrhs)+ &
1752# ifdef WEC_MELLOR
1753!^ & u_stokes(i,j,k-1)+ &
1754!^ & u_stokes(i,j,k+2)+ &
1755# endif
1756!^ & u(i,j,k+2,nrhs)))* &
1757!^ & (W(i ,j,k)+ &
1758!^ & W(i-1,j,k))
1759!^
1760 tl_fc(i,k)=(cff1*(tl_u(i,j,k ,nrhs)+ &
1761# ifdef WEC_MELLOR
1762 & tl_u_stokes(i,j,k )+ &
1763 & tl_u_stokes(i,j,k+1)+ &
1764# endif
1765 & tl_u(i,j,k+1,nrhs))- &
1766 & cff2*(tl_u(i,j,k-1,nrhs)+ &
1767# ifdef WEC_MELLOR
1768 & tl_u_stokes(i,j,k-1)+ &
1769 & tl_u_stokes(i,j,k+2)+ &
1770# endif
1771 & tl_u(i,j,k+2,nrhs)))* &
1772 & (w(i ,j,k)+ &
1773 & w(i-1,j,k))+ &
1774 & (cff1*(u(i,j,k ,nrhs)+ &
1775# ifdef WEC_MELLOR
1776 & u_stokes(i,j,k )+ &
1777 & u_stokes(i,j,k+1)+ &
1778# endif
1779 & u(i,j,k+1,nrhs))- &
1780 & cff2*(u(i,j,k-1,nrhs)+ &
1781# ifdef WEC_MELLOR
1782 & u_stokes(i,j,k-1)+ &
1783 & u_stokes(i,j,k+2)+ &
1784# endif
1785 & u(i,j,k+2,nrhs)))* &
1786 & (tl_w(i ,j,k)+ &
1787 & tl_w(i-1,j,k))
1788 END DO
1789 END DO
1790 DO i=istru,iend
1791!^ FC(i,N(ng))=0.0_r8
1792!^
1793 tl_fc(i,n(ng))=0.0_r8
1794!^ FC(i,N(ng)-1)=(cff1*(u(i,j,N(ng)-1,nrhs)+ &
1795# ifdef WEC_MELLOR
1796!^ & u_stokes(i,j,N(ng)-1)+ &
1797!^ & u_stokes(i,j,N(ng) )+ &
1798# endif
1799!^ & u(i,j,N(ng) ,nrhs))- &
1800!^ & cff2*(u(i,j,N(ng)-2,nrhs)+ &
1801# ifdef WEC_MELLOR
1802!^ & u_stokes(i,j,N(ng)-2)+ &
1803!^ & u_stokes(i,j,N(ng) )+ &
1804# endif
1805!^ & u(i,j,N(ng) ,nrhs)))* &
1806!^ & (W(i ,j,N(ng)-1)+ &
1807!^ & W(i-1,j,N(ng)-1))
1808!^
1809 tl_fc(i,n(ng)-1)=(cff1*(tl_u(i,j,n(ng)-1,nrhs)+ &
1810# ifdef WEC_MELLOR
1811 & tl_u_stokes(i,j,n(ng)-1)+ &
1812 & tl_u_stokes(i,j,n(ng) )+ &
1813# endif
1814 & tl_u(i,j,n(ng) ,nrhs))- &
1815 & cff2*(tl_u(i,j,n(ng)-2,nrhs)+ &
1816# ifdef WEC_MELLOR
1817 & tl_u_stokes(i,j,n(ng)-2)+ &
1818 & tl_u_stokes(i,j,n(ng) )+ &
1819# endif
1820 & tl_u(i,j,n(ng) ,nrhs)))* &
1821 & (w(i ,j,n(ng)-1)+ &
1822 & w(i-1,j,n(ng)-1))+ &
1823 & (cff1*(u(i,j,n(ng)-1,nrhs)+ &
1824# ifdef WEC_MELLOR
1825 & u_stokes(i,j,n(ng)-1)+ &
1826 & u_stokes(i,j,n(ng) )+ &
1827# endif
1828 & u(i,j,n(ng) ,nrhs))- &
1829 & cff2*(u(i,j,n(ng)-2,nrhs)+ &
1830# ifdef WEC_MELLOR
1831 & u_stokes(i,j,n(ng)-2)+ &
1832 & u_stokes(i,j,n(ng) )+ &
1833# endif
1834 & u(i,j,n(ng) ,nrhs)))* &
1835 & (tl_w(i ,j,n(ng)-1)+ &
1836 & tl_w(i-1,j,n(ng)-1))
1837!^ FC(i,1)=(cff1*(u(i,j,1,nrhs)+ &
1838# ifdef WEC_MELLOR
1839!^ & u_stokes(i,j,1)+ &
1840!^ & u_stokes(i,j,2)+ &
1841# endif
1842!^ & u(i,j,2,nrhs))- &
1843!^ & cff2*(u(i,j,1,nrhs)+ &
1844# ifdef WEC_MELLOR
1845!^ & u_stokes(i,j,1)+ &
1846!^ & u_stokes(i,j,3)+ &
1847# endif
1848!^ & u(i,j,3,nrhs)))* &
1849!^ & (W(i ,j,1)+ &
1850!^ & W(i-1,j,1))
1851!^
1852 tl_fc(i,1)=(cff1*(tl_u(i,j,1,nrhs)+ &
1853# ifdef WEC_MELLOR
1854 & tl_u_stokes(i,j,1)+ &
1855 & tl_u_stokes(i,j,2)+ &
1856# endif
1857 & tl_u(i,j,2,nrhs))- &
1858 & cff2*(tl_u(i,j,1,nrhs)+ &
1859# ifdef WEC_MELLOR
1860 & tl_u_stokes(i,j,1)+ &
1861 & tl_u_stokes(i,j,3)+ &
1862# endif
1863 & tl_u(i,j,3,nrhs)))* &
1864 & (w(i ,j,1)+ &
1865 & w(i-1,j,1))+ &
1866 & (cff1*(u(i,j,1,nrhs)+ &
1867# ifdef WEC_MELLOR
1868 & u_stokes(i,j,1)+ &
1869 & u_stokes(i,j,2)+ &
1870# endif
1871 & u(i,j,2,nrhs))- &
1872 & cff2*(u(i,j,1,nrhs)+ &
1873# ifdef WEC_MELLOR
1874 & u_stokes(i,j,1)+ &
1875 & u_stokes(i,j,3)+ &
1876# endif
1877 & u(i,j,3,nrhs)))* &
1878 & (tl_w(i ,j,1)+ &
1879 & tl_w(i-1,j,1))
1880!^ FC(i,0)=0.0_r8
1881!^
1882 tl_fc(i,0)=0.0_r8
1883 END DO
1884# else
1885 cff1=9.0_r8/16.0_r8
1886 cff2=1.0_r8/16.0_r8
1887 DO k=2,n(ng)-2
1888 DO i=istru,iend
1889!^ FC(i,k)=(cff1*(u(i,j,k ,nrhs)+ &
1890# ifdef WEC_MELLOR
1891!^ & u_stokes(i,j,k )+ &
1892!^ & u_stokes(i,j,k+1)+ &
1893# endif
1894!^ & u(i,j,k+1,nrhs))- &
1895!^ & cff2*(u(i,j,k-1,nrhs)+ &
1896# ifdef WEC_MELLOR
1897!^ & u_stokes(i,j,k-1)+ &
1898!^ & u_stokes(i,j,k+2)+ &
1899# endif
1900!^ & u(i,j,k+2,nrhs)))* &
1901!^ & (cff1*(W(i ,j,k)+ &
1902!^ & W(i-1,j,k))- &
1903!^ & cff2*(W(i+1,j,k)+ &
1904!^ & W(i-2,j,k)))
1905!^
1906 tl_fc(i,k)=(cff1*(tl_u(i,j,k ,nrhs)+ &
1907# ifdef WEC_MELLOR
1908 & tl_u_stokes(i,j,k )+ &
1909 & tl_u_stokes(i,j,k+1)+ &
1910# endif
1911 & tl_u(i,j,k+1,nrhs))- &
1912 & cff2*(tl_u(i,j,k-1,nrhs)+ &
1913# ifdef WEC_MELLOR
1914 & tl_u_stokes(i,j,k-1)+ &
1915 & tl_u_stokes(i,j,k+2)+ &
1916# endif
1917 & tl_u(i,j,k+2,nrhs)))* &
1918 & (cff1*(w(i ,j,k)+ &
1919 & w(i-1,j,k))- &
1920 & cff2*(w(i+1,j,k)+ &
1921 & w(i-2,j,k)))+ &
1922 & (cff1*(u(i,j,k ,nrhs)+ &
1923# ifdef WEC_MELLOR
1924 & u_stokes(i,j,k )+ &
1925 & u_stokes(i,j,k+1)+ &
1926# endif
1927 & u(i,j,k+1,nrhs))- &
1928 & cff2*(u(i,j,k-1,nrhs)+ &
1929# ifdef WEC_MELLOR
1930 & u_stokes(i,j,k-1)+ &
1931 & u_stokes(i,j,k+2)+ &
1932# endif
1933 & u(i,j,k+2,nrhs)))* &
1934 & (cff1*(tl_w(i ,j,k)+ &
1935 & tl_w(i-1,j,k))- &
1936 & cff2*(tl_w(i+1,j,k)+ &
1937 & tl_w(i-2,j,k)))
1938 END DO
1939 END DO
1940 DO i=istru,iend
1941!^ FC(i,N(ng))=0.0_r8
1942!^
1943 tl_fc(i,n(ng))=0.0_r8
1944!^ FC(i,N(ng)-1)=(cff1*(u(i,j,N(ng)-1,nrhs)+ &
1945# ifdef WEC_MELLOR
1946!^ & u_stokes(i,j,N(ng)-1)+ &
1947!^ & u_stokes(i,j,N(ng) )+ &
1948# endif
1949!^ & u(i,j,N(ng) ,nrhs))- &
1950!^ & cff2*(u(i,j,N(ng)-2,nrhs)+ &
1951# ifdef WEC_MELLOR
1952!^ & u_stokes(i,j,N(ng)-2)+ &
1953!^ & u_stokes(i,j,N(ng) )+ &
1954# endif
1955!^ & u(i,j,N(ng) ,nrhs)))* &
1956!^ & (cff1*(W(i ,j,N(ng)-1)+ &
1957!^ & W(i-1,j,N(ng)-1))- &
1958!^ & cff2*(W(i+1,j,N(ng)-1)+ &
1959!^ & W(i-2,j,N(ng)-1)))
1960!^
1961 tl_fc(i,n(ng)-1)=(cff1*(tl_u(i,j,n(ng)-1,nrhs)+ &
1962# ifdef WEC_MELLOR
1963 & tl_u_stokes(i,j,n(ng)-1)+ &
1964 & tl_u_stokes(i,j,n(ng) )+ &
1965# endif
1966 & tl_u(i,j,n(ng) ,nrhs))- &
1967 & cff2*(tl_u(i,j,n(ng)-2,nrhs)+ &
1968# ifdef WEC_MELLOR
1969 & tl_u_stokes(i,j,n(ng)-2)+ &
1970 & tl_u_stokes(i,j,n(ng) )+ &
1971# endif
1972 & tl_u(i,j,n(ng) ,nrhs)))* &
1973 & (cff1*(w(i ,j,n(ng)-1)+ &
1974 & w(i-1,j,n(ng)-1))- &
1975 & cff2*(w(i+1,j,n(ng)-1)+ &
1976 & w(i-2,j,n(ng)-1)))+ &
1977 & (cff1*(u(i,j,n(ng)-1,nrhs)+ &
1978# ifdef WEC_MELLOR
1979 & u_stokes(i,j,n(ng)-1)+ &
1980 & u_stokes(i,j,n(ng) )+ &
1981# endif
1982 & u(i,j,n(ng) ,nrhs))- &
1983 & cff2*(u(i,j,n(ng)-2,nrhs)+ &
1984# ifdef WEC_MELLOR
1985 & u_stokes(i,j,n(ng)-2)+ &
1986 & u_stokes(i,j,n(ng) )+ &
1987# endif
1988 & u(i,j,n(ng) ,nrhs)))* &
1989 & (cff1*(tl_w(i ,j,n(ng)-1)+ &
1990 & tl_w(i-1,j,n(ng)-1))- &
1991 & cff2*(tl_w(i+1,j,n(ng)-1)+ &
1992 & tl_w(i-2,j,n(ng)-1)))
1993!^ FC(i,1)=(cff1*(u(i,j,1,nrhs)+ &
1994# ifdef WEC_MELLOR
1995!^ & u_stokes(i,j,1)+ &
1996!^ & u_stokes(i,j,2)+ &
1997# endif
1998!^ & u(i,j,2,nrhs))- &
1999!^ & cff2*(u(i,j,1,nrhs)+ &
2000# ifdef WEC_MELLOR
2001!^ & u_stokes(i,j,1)+ &
2002!^ & u_stokes(i,j,3)+ &
2003# endif
2004!^ & u(i,j,3,nrhs)))* &
2005!^ & (cff1*(W(i ,j,1)+ &
2006!^ & W(i-1,j,1))- &
2007!^ & cff2*(W(i+1,j,1)+ &
2008!^ & W(i-2,j,1)))
2009!^
2010 tl_fc(i,1)=(cff1*(tl_u(i,j,1,nrhs)+ &
2011# ifdef WEC_MELLOR
2012 & tl_u_stokes(i,j,1)+ &
2013 & tl_u_stokes(i,j,2)+ &
2014# endif
2015 & tl_u(i,j,2,nrhs))- &
2016 & cff2*(tl_u(i,j,1,nrhs)+ &
2017# ifdef WEC_MELLOR
2018 & tl_u_stokes(i,j,1)+ &
2019 & tl_u_stokes(i,j,3)+ &
2020# endif
2021 & tl_u(i,j,3,nrhs)))* &
2022 & (cff1*(w(i ,j,1)+ &
2023 & w(i-1,j,1))- &
2024 & cff2*(w(i+1,j,1)+ &
2025 & w(i-2,j,1)))+ &
2026 & (cff1*(u(i,j,1,nrhs)+ &
2027# ifdef WEC_MELLOR
2028 & u_stokes(i,j,1)+ &
2029 & u_stokes(i,j,2)+ &
2030# endif
2031 & u(i,j,2,nrhs))- &
2032 & cff2*(u(i,j,1,nrhs)+ &
2033# ifdef WEC_MELLOR
2034 & u_stokes(i,j,1)+ &
2035 & u_stokes(i,j,3)+ &
2036# endif
2037 & u(i,j,3,nrhs)))* &
2038 & (cff1*(tl_w(i ,j,1)+ &
2039 & tl_w(i-1,j,1))- &
2040 & cff2*(tl_w(i+1,j,1)+ &
2041 & tl_w(i-2,j,1)))
2042!^ FC(i,0)=0.0_r8
2043!^
2044 tl_fc(i,0)=0.0_r8
2045 END DO
2046# endif
2047 DO k=1,n(ng)
2048 DO i=istru,iend
2049!^ cff=FC(i,k)-FC(i,k-1)
2050!^
2051 tl_cff=tl_fc(i,k)-tl_fc(i,k-1)
2052!^ ru(i,j,k,nrhs)=ru(i,j,k,nrhs)-cff
2053!^
2054 tl_ru(i,j,k,nrhs)=tl_ru(i,j,k,nrhs)-tl_cff
2055# ifdef DIAGNOSTICS_UV
2056!! DiaRU(i,j,k,nrhs,M3vadv)=-cff
2057# endif
2058 END DO
2059 END DO
2060 IF (j.ge.jstrv) THEN
2061# ifdef UV_SADVECTION
2062!
2063! Apply spline code to BASIC STATE v-momentum which should be in
2064! units of m/s. CF will be used by the tangent linear spline code.
2065!
2066 cff1=9.0_r8/16.0_r8
2067 cff2=1.0_r8/16.0_r8
2068 DO k=1,n(ng)
2069 DO i=istr,iend
2070 dc(i,k)=(cff1*(hz(i,j ,k)+ &
2071 & hz(i,j-1,k))- &
2072 & cff2*(hz(i,j+1,k)+ &
2073 & hz(i,j-2,k)))
2074 END DO
2075 END DO
2076 DO i=istr,iend
2077 fc(i,0)=0.0_r8
2078 cf(i,0)=0.0_r8
2079 END DO
2080 DO k=1,n(ng)-1
2081 DO i=istr,iend
2082 cff=1.0_r8/(2.0_r8*dc(i,k+1)+dc(i,k)*(2.0_r8-fc(i,k-1)))
2083 fc(i,k)=cff*dc(i,k+1)
2084 cf(i,k)=cff*(6.0_r8*(v(i,j,k+1,nrhs)- &
2085# ifdef WEC_MELLOR
2086 & v_stokes(i,j,k )+ &
2087 & v_stokes(i,j,k+1)- &
2088# endif
2089 & v(i,j,k ,nrhs))- &
2090 & dc(i,k)*cf(i,k-1))
2091 END DO
2092 END DO
2093 DO i=istr,iend
2094 cf(i,n(ng))=0.0_r8
2095 END DO
2096 DO k=n(ng)-1,1,-1
2097 DO i=istr,iend
2098 cf(i,k)=cf(i,k)-fc(i,k)*cf(i,k+1)
2099 END DO
2100 END DO
2101!
2102! Construct tangent linear conservative parabolic splines for the
2103! vertical derivatives "tl_CF" of v-momentum.
2104!
2105 cff1=9.0_r8/16.0_r8
2106 cff2=1.0_r8/16.0_r8
2107 DO k=1,n(ng)
2108 DO i=istr,iend
2109 dc(i,k)=(cff1*(hz(i,j ,k)+ &
2110 & hz(i,j-1,k))- &
2111 & cff2*(hz(i,j+1,k)+ &
2112 & hz(i,j-2,k)))
2113 tl_dc(i,k)=(cff1*(tl_hz(i,j ,k)+ &
2114 & tl_hz(i,j-1,k))- &
2115 & cff2*(tl_hz(i,j+1,k)+ &
2116 & tl_hz(i,j-2,k)))
2117 END DO
2118 END DO
2119 DO i=istr,iend
2120 fc(i,0)=0.0_r8
2121 tl_cf(i,0)=0.0_r8
2122 END DO
2123 DO k=1,n(ng)-1
2124 DO i=istr,iend
2125 cff=1.0_r8/(2.0_r8*dc(i,k+1)+dc(i,k)*(2.0_r8-fc(i,k-1)))
2126 fc(i,k)=cff*dc(i,k+1)
2127 tl_cf(i,k)=cff*(6.0_r8*(tl_v(i,j,k+1,nrhs)- &
2128# ifdef WEC_MELLOR
2129 & tl_v_stokes(i,j,k )+ &
2130 & tl_v_stokes(i,j,k+1)- &
2131# endif
2132 & tl_v(i,j,k ,nrhs))- &
2133 & (tl_dc(i,k)*cf(i,k-1)+ &
2134 & 2.0_r8*(tl_dc(i,k )+ &
2135 & tl_dc(i,k+1))*cf(i,k)+ &
2136 & tl_dc(i,k+1)*cf(i,k+1))- &
2137 & dc(i,k)*tl_cf(i,k-1))
2138 END DO
2139 END DO
2140 DO i=istr,iend
2141 tl_cf(i,n(ng))=0.0_r8
2142 END DO
2143 DO k=n(ng)-1,1,-1
2144 DO i=istr,iend
2145 tl_cf(i,k)=tl_cf(i,k)-fc(i,k)*tl_cf(i,k+1)
2146 END DO
2147 END DO
2148!
2149! Compute spline-interpolated, vertical advective v-momentum flux.
2150!
2151 cff3=1.0_r8/3.0_r8
2152 cff4=1.0_r8/6.0_r8
2153 DO k=1,n(ng)-1
2154 DO i=istr,iend
2155!^ FC(i,k)=(cff1*(W(i,j ,k)+ &
2156!^ & W(i,j-1,k))- &
2157!^ & cff2*(W(i,j+1,k)+ &
2158!^ & W(i,j-2,k)))* &
2159!^ & (v(i,j,k,nrhs)+ &
2160# ifdef WEC_MELLOR
2161!^ & v_stokes(i,j,k)+ &
2162# endif
2163!^ & DC(i,k)*(cff3*CF(i,k )+ &
2164!^ & cff4*CF(i,k-1)))
2165!^
2166 tl_fc(i,k)=(cff1*(tl_w(i,j ,k)+ &
2167 & tl_w(i,j-1,k))- &
2168 & cff2*(tl_w(i,j+1,k)+ &
2169 & tl_w(i,j-2,k)))* &
2170 & (v(i,j,k,nrhs)+ &
2171# ifdef WEC_MELLOR
2172 & v_stokes(i,j,k)+ &
2173# endif
2174 & dc(i,k)*(cff3*cf(i,k )+ &
2175 & cff4*cf(i,k-1)))+ &
2176 & (cff1*(w(i,j ,k)+ &
2177 & w(i,j-1,k))- &
2178 & cff2*(w(i,j+1,k)+ &
2179 & w(i,j-2,k)))* &
2180 & (tl_v(i,j,k,nrhs)+ &
2181# ifdef WEC_MELLOR
2182 & tl_v_stokes(i,j,k)+ &
2183# endif
2184 & dc(i,k)*(cff3*tl_cf(i,k )+ &
2185 & cff4*tl_cf(i,k-1))+ &
2186 & tl_dc(i,k)*(cff3*cf(i,k )+ &
2187 & cff4*cf(i,k-1)))
2188 END DO
2189 END DO
2190 DO i=istr,iend
2191!^ FC(i,N(ng))=0.0_r8
2192!^
2193 tl_fc(i,n(ng))=0.0_r8
2194!^ FC(i,0)=0.0_r8
2195!^
2196 tl_fc(i,0)=0.0_r8
2197 END DO
2198# elif defined UV_C2ADVECTION
2199!
2200! Second-order, centered differences vertical advection.
2201!
2202 DO k=1,n(ng)-1
2203 DO i=istr,iend
2204!^ FC(i,k)=0.25_r8*(v(i,j,k ,nrhs)+ &
2205# ifdef WEC_MELLOR
2206!^ & v_stokes(i,j,k )+ &
2207!^ & v_stokes(i,j,k+1)+ &
2208# endif
2209!^ & v(i,j,k+1,nrhs))* &
2210!^ & (W(i,j ,k)+ &
2211!^ & W(i,j-1,k))
2212!^
2213 tl_fc(i,k)=0.25_r8*((tl_v(i,j,k ,nrhs)+ &
2214# ifdef WEC_MELLOR
2215 & tl_v_stokes(i,j,k )+ &
2216 & tl_v_stokes(i,j,k+1)+ &
2217# endif
2218 & tl_v(i,j,k+1,nrhs))* &
2219 & (w(i,j ,k)+ &
2220 & w(i,j-1,k))+ &
2221 & (v(i,j,k ,nrhs)+ &
2222# ifdef WEC_MELLOR
2223 & v_stokes(i,j,k )+ &
2224 & v_stokes(i,j,k+1)+ &
2225# endif
2226 & v(i,j,k+1,nrhs))* &
2227 & (tl_w(i,j ,k)+ &
2228 & tl_w(i,j-1,k)))
2229 END DO
2230 END DO
2231 DO i=istr,iend
2232!^ FC(i,0)=0.0_r8
2233!^
2234 tl_fc(i,0)=0.0_r8
2235!^ FC(i,N(ng))=0.0_r8
2236!^
2237 tl_fc(i,n(ng))=0.0_r8
2238 END DO
2239# elif defined UV_C4ADVECTION
2240!
2241! Forth-order, centered differences vertical advection.
2242!
2243 cff1=9.0_r8/32.0_r8
2244 cff2=1.0_r8/32.0_r8
2245 DO k=2,n(ng)-2
2246 DO i=istr,iend
2247!^ FC(i,k)=(cff1*(v(i,j,k ,nrhs)+ &
2248# ifdef WEC_MELLOR
2249!^ & v_stokes(i,j,k )+ &
2250!^ & v_stokes(i,j,k+1)+ &
2251# endif
2252!^ & v(i,j,k+1,nrhs))- &
2253!^ & cff2*(v(i,j,k-1,nrhs)+ &
2254# ifdef WEC_MELLOR
2255!^ & v_stokes(i,j,k-1)+ &
2256!^ & v_stokes(i,j,k+2)+ &
2257# endif
2258!^ & v(i,j,k+2,nrhs)))* &
2259!^ & (W(i,j ,k)+ &
2260!^ & W(i,j-1,k))
2261!^
2262 tl_fc(i,k)=(cff1*(tl_v(i,j,k ,nrhs)+ &
2263# ifdef WEC_MELLOR
2264 & tl_v_stokes(i,j,k )+ &
2265 & tl_v_stokes(i,j,k+1)+ &
2266# endif
2267 & tl_v(i,j,k+1,nrhs))- &
2268 & cff2*(tl_v(i,j,k-1,nrhs)+ &
2269# ifdef WEC_MELLOR
2270 & tl_v_stokes(i,j,k-1)+ &
2271 & tl_v_stokes(i,j,k+2)+ &
2272# endif
2273 & tl_v(i,j,k+2,nrhs)))* &
2274 & (w(i,j ,k)+ &
2275 & w(i,j-1,k))+ &
2276 & (cff1*(v(i,j,k ,nrhs)+ &
2277# ifdef WEC_MELLOR
2278 & v_stokes(i,j,k )+ &
2279 & v_stokes(i,j,k+1)+ &
2280# endif
2281 & v(i,j,k+1,nrhs))- &
2282 & cff2*(v(i,j,k-1,nrhs)+ &
2283# ifdef WEC_MELLOR
2284 & v_stokes(i,j,k-1)+ &
2285 & v_stokes(i,j,k+2)+ &
2286# endif
2287 & v(i,j,k+2,nrhs)))* &
2288 & (tl_w(i,j ,k)+ &
2289 & tl_w(i,j-1,k))
2290 END DO
2291 END DO
2292 DO i=istr,iend
2293!^ FC(i,N(ng))=0.0_r8
2294!^
2295 tl_fc(i,n(ng))=0.0_r8
2296!^ FC(i,N(ng)-1)=(cff1*(v(i,j,N(ng)-1,nrhs)+ &
2297# ifdef WEC_MELLOR
2298!^ & v_stokes(i,j,N(ng)-1)+ &
2299!^ & v_stokes(i,j,N(ng) )+ &
2300# endif
2301!^ & v(i,j,N(ng) ,nrhs))- &
2302!^ & cff2*(v(i,j,N(ng)-2,nrhs)+ &
2303# ifdef WEC_MELLOR
2304!^ & v_stokes(i,j,N(ng)-2)+ &
2305!^ & v_stokes(i,j,N(ng) )+ &
2306# endif
2307!^ & v(i,j,N(ng) ,nrhs)))* &
2308!^ & (W(i,j ,N(ng)-1)+ &
2309!^ & W(i,j-1,N(ng)-1))
2310!^
2311 tl_fc(i,n(ng)-1)=(cff1*(tl_v(i,j,n(ng)-1,nrhs)+ &
2312# ifdef WEC_MELLOR
2313 & tl_v_stokes(i,j,n(ng)-1)+ &
2314 & tl_v_stokes(i,j,n(ng) )+ &
2315# endif
2316 & tl_v(i,j,n(ng) ,nrhs))- &
2317 & cff2*(tl_v(i,j,n(ng)-2,nrhs)+ &
2318# ifdef WEC_MELLOR
2319 & tl_v_stokes(i,j,n(ng)-2)+ &
2320 & tl_v_stokes(i,j,n(ng) )+ &
2321# endif
2322 & tl_v(i,j,n(ng) ,nrhs)))* &
2323 & (w(i,j ,n(ng)-1)+ &
2324 & w(i,j-1,n(ng)-1))+ &
2325 & (cff1*(v(i,j,n(ng)-1,nrhs)+ &
2326# ifdef WEC_MELLOR
2327 & v_stokes(i,j,n(ng)-1)+ &
2328 & v_stokes(i,j,n(ng) )+ &
2329# endif
2330 & v(i,j,n(ng) ,nrhs))- &
2331 & cff2*(v(i,j,n(ng)-2,nrhs)+ &
2332# ifdef WEC_MELLOR
2333 & v_stokes(i,j,n(ng)-2)+ &
2334 & v_stokes(i,j,n(ng) )+ &
2335# endif
2336 & v(i,j,n(ng) ,nrhs)))* &
2337 & (tl_w(i,j ,n(ng)-1)+ &
2338 & tl_w(i,j-1,n(ng)-1))
2339!^ FC(i,1)=(cff1*(v(i,j,1,nrhs)+ &
2340# ifdef WEC_MELLOR
2341!^ & v_stokes(i,j,1)+ &
2342!^ & v_stokes(i,j,2)+ &
2343# endif
2344!^ & v(i,j,2,nrhs))- &
2345!^ & cff2*(v(i,j,1,nrhs)+ &
2346# ifdef WEC_MELLOR
2347!^ & v_stokes(i,j,1)+ &
2348!^ & v_stokes(i,j,3)+ &
2349# endif
2350!^ & v(i,j,3,nrhs)))* &
2351!^ & (W(i,j ,1)+ &
2352!^ & W(i,j-1,1))
2353!^
2354 tl_fc(i,1)=(cff1*(tl_v(i,j,1,nrhs)+ &
2355# ifdef WEC_MELLOR
2356 & tl_v_stokes(i,j,1)+ &
2357 & tl_v_stokes(i,j,2)+ &
2358# endif
2359 & tl_v(i,j,2,nrhs))- &
2360 & cff2*(tl_v(i,j,1,nrhs)+ &
2361# ifdef WEC_MELLOR
2362 & tl_v_stokes(i,j,1)+ &
2363 & tl_v_stokes(i,j,3)+ &
2364# endif
2365 & tl_v(i,j,3,nrhs)))* &
2366 & (w(i,j ,1)+ &
2367 & w(i,j-1,1))+ &
2368 & (cff1*(v(i,j,1,nrhs)+ &
2369# ifdef WEC_MELLOR
2370 & v_stokes(i,j,1)+ &
2371 & v_stokes(i,j,2)+ &
2372# endif
2373 & v(i,j,2,nrhs))- &
2374 & cff2*(v(i,j,1,nrhs)+ &
2375# ifdef WEC_MELLOR
2376 & v_stokes(i,j,1)+ &
2377 & v_stokes(i,j,3)+ &
2378# endif
2379 & v(i,j,3,nrhs)))* &
2380 & (tl_w(i,j ,1)+ &
2381 & tl_w(i,j-1,1))
2382!^ FC(i,0)=0.0_r8
2383!^
2384 tl_fc(i,0)=0.0_r8
2385 END DO
2386# else
2387 cff1=9.0_r8/16.0_r8
2388 cff2=1.0_r8/16.0_r8
2389 DO k=2,n(ng)-2
2390 DO i=istr,iend
2391!^ FC(i,k)=(cff1*(v(i,j,k ,nrhs)+ &
2392# ifdef WEC_MELLOR
2393!^ & v_stokes(i,j,k )+ &
2394!^ & v_stokes(i,j,k+1)+ &
2395# endif
2396!^ & v(i,j,k+1,nrhs))- &
2397!^ & cff2*(v(i,j,k-1,nrhs)+ &
2398# ifdef WEC_MELLOR
2399!^ & v_stokes(i,j,k-1)+ &
2400!^ & v_stokes(i,j,k+2)+ &
2401# endif
2402!^ & v(i,j,k+2,nrhs)))* &
2403!^ & (cff1*(W(i,j ,k)+ &
2404!^ & W(i,j-1,k))- &
2405!^ & cff2*(W(i,j+1,k)+ &
2406!^ & W(i,j-2,k)))
2407!^
2408 tl_fc(i,k)=(cff1*(tl_v(i,j,k ,nrhs)+ &
2409# ifdef WEC_MELLOR
2410 & tl_v_stokes(i,j,k )+ &
2411 & tl_v_stokes(i,j,k+1)+ &
2412# endif
2413 & tl_v(i,j,k+1,nrhs))- &
2414 & cff2*(tl_v(i,j,k-1,nrhs)+ &
2415# ifdef WEC_MELLOR
2416 & tl_v_stokes(i,j,k-1)+ &
2417 & tl_v_stokes(i,j,k+2)+ &
2418# endif
2419 & tl_v(i,j,k+2,nrhs)))* &
2420 & (cff1*(w(i,j ,k)+ &
2421 & w(i,j-1,k))- &
2422 & cff2*(w(i,j+1,k)+ &
2423 & w(i,j-2,k)))+ &
2424 & (cff1*(v(i,j,k ,nrhs)+ &
2425# ifdef WEC_MELLOR
2426 & v_stokes(i,j,k )+ &
2427 & v_stokes(i,j,k+1)+ &
2428# endif
2429 & v(i,j,k+1,nrhs))- &
2430 & cff2*(v(i,j,k-1,nrhs)+ &
2431# ifdef WEC_MELLOR
2432 & v_stokes(i,j,k-1)+ &
2433 & v_stokes(i,j,k+2)+ &
2434# endif
2435 & v(i,j,k+2,nrhs)))* &
2436 & (cff1*(tl_w(i,j ,k)+ &
2437 & tl_w(i,j-1,k))- &
2438 & cff2*(tl_w(i,j+1,k)+ &
2439 & tl_w(i,j-2,k)))
2440 END DO
2441 END DO
2442 DO i=istr,iend
2443!^ FC(i,N(ng))=0.0_r8
2444!^
2445 tl_fc(i,n(ng))=0.0_r8
2446!^ FC(i,N(ng)-1)=(cff1*(v(i,j,N(ng)-1,nrhs)+ &
2447# ifdef WEC_MELLOR
2448!^ & v_stokes(i,j,N(ng)-1)+ &
2449!^ & v_stokes(i,j,N(ng) )+ &
2450# endif
2451!^ & v(i,j,N(ng) ,nrhs))- &
2452!^ & cff2*(v(i,j,N(ng)-2,nrhs)+ &
2453# ifdef WEC_MELLOR
2454!^ & v_stokes(i,j,N(ng)-2)+ &
2455!^ & v_stokes(i,j,N(ng) )+ &
2456# endif
2457!^ & v(i,j,N(ng) ,nrhs)))* &
2458!^ & (cff1*(W(i,j ,N(ng)-1)+ &
2459!^ & W(i,j-1,N(ng)-1))- &
2460!^ & cff2*(W(i,j+1,N(ng)-1)+ &
2461!^ & W(i,j-2,N(ng)-1)))
2462!^
2463 tl_fc(i,n(ng)-1)=(cff1*(tl_v(i,j,n(ng)-1,nrhs)+ &
2464# ifdef WEC_MELLOR
2465 & tl_v_stokes(i,j,n(ng)-1)+ &
2466 & tl_v_stokes(i,j,n(ng) )+ &
2467# endif
2468 & tl_v(i,j,n(ng) ,nrhs))- &
2469 & cff2*(tl_v(i,j,n(ng)-2,nrhs)+ &
2470# ifdef WEC_MELLOR
2471 & tl_v_stokes(i,j,n(ng)-2)+ &
2472 & tl_v_stokes(i,j,n(ng) )+ &
2473# endif
2474 & tl_v(i,j,n(ng) ,nrhs)))* &
2475 & (cff1*(w(i,j ,n(ng)-1)+ &
2476 & w(i,j-1,n(ng)-1))- &
2477 & cff2*(w(i,j+1,n(ng)-1)+ &
2478 & w(i,j-2,n(ng)-1)))+ &
2479 & (cff1*(v(i,j,n(ng)-1,nrhs)+ &
2480# ifdef WEC_MELLOR
2481 & v_stokes(i,j,n(ng)-1)+ &
2482 & v_stokes(i,j,n(ng) )+ &
2483# endif
2484 & v(i,j,n(ng) ,nrhs))- &
2485 & cff2*(v(i,j,n(ng)-2,nrhs)+ &
2486# ifdef WEC_MELLOR
2487 & v_stokes(i,j,n(ng)-2)+ &
2488 & v_stokes(i,j,n(ng) )+ &
2489# endif
2490 & v(i,j,n(ng) ,nrhs)))* &
2491 & (cff1*(tl_w(i,j ,n(ng)-1)+ &
2492 & tl_w(i,j-1,n(ng)-1))- &
2493 & cff2*(tl_w(i,j+1,n(ng)-1)+ &
2494 & tl_w(i,j-2,n(ng)-1)))
2495!^ FC(i,1)=(cff1*(v(i,j,1,nrhs)+ &
2496# ifdef WEC_MELLOR
2497!^ & v_stokes(i,j,1)+ &
2498!^ & v_stokes(i,j,2)+ &
2499# endif
2500!^ & v(i,j,2,nrhs))- &
2501!^ & cff2*(v(i,j,1,nrhs)+ &
2502# ifdef WEC_MELLOR
2503!^ & v_stokes(i,j,1)+ &
2504!^ & v_stokes(i,j,3)+ &
2505# endif
2506!^ & v(i,j,3,nrhs)))* &
2507!^ & (cff1*(W(i,j ,1)+ &
2508!^ & W(i,j-1,1))- &
2509!^ & cff2*(W(i,j+1,1)+ &
2510!^ & W(i,j-2,1)))
2511!^
2512 tl_fc(i,1)=(cff1*(tl_v(i,j,1,nrhs)+ &
2513# ifdef WEC_MELLOR
2514 & tl_v_stokes(i,j,1)+ &
2515 & tl_v_stokes(i,j,2)+ &
2516# endif
2517 & tl_v(i,j,2,nrhs))- &
2518 & cff2*(tl_v(i,j,1,nrhs)+ &
2519# ifdef WEC_MELLOR
2520 & tl_v_stokes(i,j,1)+ &
2521 & tl_v_stokes(i,j,3)+ &
2522# endif
2523 & tl_v(i,j,3,nrhs)))* &
2524 & (cff1*(w(i,j ,1)+ &
2525 & w(i,j-1,1))- &
2526 & cff2*(w(i,j+1,1)+ &
2527 & w(i,j-2,1)))+ &
2528 & (cff1*(v(i,j,1,nrhs)+ &
2529# ifdef WEC_MELLOR
2530 & v_stokes(i,j,1)+ &
2531 & v_stokes(i,j,2)+ &
2532# endif
2533 & v(i,j,2,nrhs))- &
2534 & cff2*(v(i,j,1,nrhs)+ &
2535# ifdef WEC_MELLOR
2536 & v_stokes(i,j,1)+ &
2537 & v_stokes(i,j,3)+ &
2538# endif
2539 & v(i,j,3,nrhs)))* &
2540 & (cff1*(tl_w(i,j ,1)+ &
2541 & tl_w(i,j-1,1))- &
2542 & cff2*(tl_w(i,j+1,1)+ &
2543 & tl_w(i,j-2,1)))
2544!^ FC(i,0)=0.0_r8
2545!^
2546 tl_fc(i,0)=0.0_r8
2547 END DO
2548# endif
2549 DO k=1,n(ng)
2550 DO i=istr,iend
2551!^ cff=FC(i,k)-FC(i,k-1)
2552!^
2553 tl_cff=tl_fc(i,k)-tl_fc(i,k-1)
2554!^ rv(i,j,k,nrhs)=rv(i,j,k,nrhs)-cff
2555!^
2556 tl_rv(i,j,k,nrhs)=tl_rv(i,j,k,nrhs)-tl_cff
2557# ifdef DIAGNOSTICS_UV
2558!! DiaRV(i,j,k,nrhs,M3vadv)=-cff
2559# endif
2560 END DO
2561 END DO
2562 END IF
2563# endif
2564!
2565!-----------------------------------------------------------------------
2566! Compute forcing term for the 2D momentum equations.
2567!-----------------------------------------------------------------------
2568!
2569! Vertically integrate baroclinic right-hand-side terms. If not
2570! body force stresses, add in the difference between surface and
2571! bottom stresses.
2572!
2573 DO i=istru,iend
2574# ifdef WET_DRY_NOT_YET
2575!^ ru(i,j,1,nrhs)=ru(i,j,1,nrhs)*umask_wet(i,j)
2576!^
2577 tl_ru(i,j,1,nrhs)=tl_ru(i,j,1,nrhs)*umask_wet(i,j)
2578# endif
2579!^ rufrc(i,j)=ru(i,j,1,nrhs)
2580!^
2581 tl_rufrc(i,j)=tl_ru(i,j,1,nrhs)
2582# ifdef DIAGNOSTICS_UV
2583!! DiaRUfrc(i,j,3,M2pgrd)=DiaRU(i,j,1,nrhs,M3pgrd)
2584# ifdef UV_COR
2585!! DiaRUfrc(i,j,3,M2fcor)=DiaRU(i,j,1,nrhs,M3fcor)
2586# endif
2587# ifdef UV_ADV
2588!! DiaRUfrc(i,j,3,M2xadv)=DiaRU(i,j,1,nrhs,M3xadv)
2589!! DiaRUfrc(i,j,3,M2yadv)=DiaRU(i,j,1,nrhs,M3yadv)
2590!! DiaRUfrc(i,j,3,M2hadv)=DiaRU(i,j,1,nrhs,M3hadv)
2591# endif
2592# ifdef WEC_MELLOR
2593!! DiaRUfrc(i,j,3,M2hrad)=DiaRU(i,j,1,nrhs,M3hrad)
2594# endif
2595# if defined UV_VIS2 || defined UV_VIS4
2596!! DiaRUfrc(i,j,3,M2xvis)=0.0_r8
2597!! DiaRUfrc(i,j,3,M2yvis)=0.0_r8
2598!! DiaRUfrc(i,j,3,M2hvis)=0.0_r8
2599# endif
2600# ifdef BODYFORCE
2601!! DiaRUfrc(i,j,3,M2strs)=DiaRU(i,j,1,nrhs,M3vvis)
2602# endif
2603# endif
2604 END DO
2605 DO k=2,n(ng)
2606 DO i=istru,iend
2607# ifdef WET_DRY_NOT_YET
2608!^ ru(i,j,k,nrhs)=ru(i,j,k,nrhs)*umask_wet(i,j)
2609!^
2610 tl_ru(i,j,k,nrhs)=tl_ru(i,j,k,nrhs)*umask_wet(i,j)
2611# endif
2612!^ rufrc(i,j)=rufrc(i,j)+ru(i,j,k,nrhs)
2613!^
2614 tl_rufrc(i,j)=tl_rufrc(i,j)+tl_ru(i,j,k,nrhs)
2615# ifdef DIAGNOSTICS_UV
2616!! DiaRUfrc(i,j,3,M2pgrd)=DiaRUfrc(i,j,3,M2pgrd)+ &
2617!! & DiaRU(i,j,k,nrhs,M3pgrd)
2618# ifdef UV_COR
2619!! DiaRUfrc(i,j,3,M2fcor)=DiaRUfrc(i,j,3,M2fcor)+ &
2620!! & DiaRU(i,j,k,nrhs,M3fcor)
2621# endif
2622# ifdef UV_ADV
2623!! DiaRUfrc(i,j,3,M2xadv)=DiaRUfrc(i,j,3,M2xadv)+ &
2624!! & DiaRU(i,j,k,nrhs,M3xadv)
2625!! DiaRUfrc(i,j,3,M2yadv)=DiaRUfrc(i,j,3,M2yadv)+ &
2626!! & DiaRU(i,j,k,nrhs,M3yadv)
2627!! DiaRUfrc(i,j,3,M2hadv)=DiaRUfrc(i,j,3,M2hadv)+ &
2628!! & DiaRU(i,j,k,nrhs,M3hadv)
2629# endif
2630# ifdef WEC_MELLOR
2631!! DiaRUfrc(i,j,3,M2hrad)=DiaRUfrc(i,j,3,M2hrad)+ &
2632!! & DiaRU(i,j,k,nrhs,M3hrad)
2633# endif
2634# ifdef BODYFORCE
2635!! DiaRUfrc(i,j,3,M2strs)=DiaRUfrc(i,j,3,M2strs)+ &
2636!! & DiaRU(i,j,k,nrhs,M3vvis)
2637# endif
2638# endif
2639 END DO
2640 END DO
2641# ifndef BODYFORCE
2642 DO i=istru,iend
2643 cff=om_u(i,j)*on_u(i,j)
2644!^ cff1= sustr(i,j)*cff
2645!^
2646 tl_cff1= tl_sustr(i,j)*cff
2647!^ cff2=-bustr(i,j)*cff
2648!^
2649 tl_cff2=-tl_bustr(i,j)*cff
2650!^ rufrc(i,j)=rufrc(i,j)+cff1+cff2
2651!^
2652 tl_rufrc(i,j)=tl_rufrc(i,j)+tl_cff1+tl_cff2
2653# ifdef WET_DRY_NOT_YET
2654!^ rufrc(i,j)=rufrc(i,j)*umask_wet(i,j)
2655!^
2656 tl_rufrc(i,j)=tl_rufrc(i,j)*umask_wet(i,j)
2657# endif
2658# ifdef DIAGNOSTICS_UV
2659!! DiaRUfrc(i,j,3,M2sstr)=cff1
2660!! DiaRUfrc(i,j,3,M2bstr)=cff2
2661# endif
2662 END DO
2663# endif
2664 IF (j.ge.jstrv) THEN
2665 DO i=istr,iend
2666# ifdef WET_DRY_NOT_YET
2667!^ rv(i,j,1,nrhs)=rv(i,j,1,nrhs)*vmask_wet(i,j)
2668!^
2669 tl_rv(i,j,1,nrhs)=tl_rv(i,j,1,nrhs)*vmask_wet(i,j)
2670# endif
2671!^ rvfrc(i,j)=rv(i,j,1,nrhs)
2672!^
2673 tl_rvfrc(i,j)=tl_rv(i,j,1,nrhs)
2674# ifdef DIAGNOSTICS_UV
2675!! DiaRVfrc(i,j,3,M2pgrd)=DiaRV(i,j,1,nrhs,M3pgrd)
2676# ifdef UV_COR
2677!! DiaRVfrc(i,j,3,M2fcor)=DiaRV(i,j,1,nrhs,M3fcor)
2678# endif
2679# ifdef UV_ADV
2680!! DiaRVfrc(i,j,3,M2xadv)=DiaRV(i,j,1,nrhs,M3xadv)
2681!! DiaRVfrc(i,j,3,M2yadv)=DiaRV(i,j,1,nrhs,M3yadv)
2682!! DiaRVfrc(i,j,3,M2hadv)=DiaRV(i,j,1,nrhs,M3hadv)
2683# endif
2684# ifdef WEC_MELLOR
2685!! DiaRVfrc(i,j,3,M2hrad)=DiaRV(i,j,1,nrhs,M3hrad)
2686# endif
2687# if defined UV_VIS2 || defined UV_VIS4
2688!! DiaRVfrc(i,j,3,M2hvis)=0.0_r8
2689!! DiaRVfrc(i,j,3,M2xvis)=0.0_r8
2690!! DiaRVfrc(i,j,3,M2yvis)=0.0_r8
2691# endif
2692# ifdef BODYFORCE
2693!! DiaRVfrc(i,j,3,M2strs)=DiaRV(i,j,1,nrhs,M3vvis)
2694# endif
2695# endif
2696 END DO
2697 DO k=2,n(ng)
2698 DO i=istr,iend
2699# ifdef WET_DRY_NOT_YET
2700!^ rv(i,j,k,nrhs)=rv(i,j,k,nrhs)*vmask_wet(i,j)
2701!^
2702 tl_rv(i,j,k,nrhs)=tl_rv(i,j,k,nrhs)*vmask_wet(i,j)
2703# endif
2704!^ rvfrc(i,j)=rvfrc(i,j)+rv(i,j,k,nrhs)
2705!^
2706 tl_rvfrc(i,j)=tl_rvfrc(i,j)+tl_rv(i,j,k,nrhs)
2707# ifdef DIAGNOSTICS_UV
2708!! DiaRVfrc(i,j,3,M2pgrd)=DiaRVfrc(i,j,3,M2pgrd)+ &
2709!! & DiaRV(i,j,k,nrhs,M3pgrd)
2710# ifdef UV_COR
2711!! DiaRVfrc(i,j,3,M2fcor)=DiaRVfrc(i,j,3,M2fcor)+ &
2712!! & DiaRV(i,j,k,nrhs,M3fcor)
2713# endif
2714# ifdef UV_ADV
2715!! DiaRVfrc(i,j,3,M2xadv)=DiaRVfrc(i,j,3,M2xadv)+ &
2716!! & DiaRV(i,j,k,nrhs,M3xadv)
2717!! DiaRVfrc(i,j,3,M2yadv)=DiaRVfrc(i,j,3,M2yadv)+ &
2718!! & DiaRV(i,j,k,nrhs,M3yadv)
2719!! DiaRVfrc(i,j,3,M2hadv)=DiaRVfrc(i,j,3,M2hadv)+ &
2720!! & DiaRV(i,j,k,nrhs,M3hadv)
2721# endif
2722# ifdef WEC_MELLOR
2723!! DiaRVfrc(i,j,3,M2hrad)=DiaRVfrc(i,j,3,M2hrad)+ &
2724!! & DiaRV(i,j,k,nrhs,M3hrad)
2725# endif
2726# ifdef BODYFORCE
2727!! DiaRVfrc(i,j,3,M2strs)=DiaRVfrc(i,j,3,M2strs)+ &
2728!! & DiaRV(i,j,k,nrhs,M3vvis)
2729# endif
2730# endif
2731 END DO
2732 END DO
2733# ifndef BODYFORCE
2734 DO i=istr,iend
2735 cff=om_v(i,j)*on_v(i,j)
2736!^ cff1= svstr(i,j)*cff
2737!^
2738 tl_cff1= tl_svstr(i,j)*cff
2739!^ cff2=-bvstr(i,j)*cff
2740!^
2741 tl_cff2=-tl_bvstr(i,j)*cff
2742!^ rvfrc(i,j)=rvfrc(i,j)+cff1+cff2
2743!^
2744 tl_rvfrc(i,j)=tl_rvfrc(i,j)+tl_cff1+tl_cff2
2745# ifdef WET_DRY_NOT_YET
2746!^ rvfrc(i,j)=rvfrc(i,j)*vmask_wet(i,j)
2747!^
2748 tl_rvfrc(i,j)=tl_rvfrc(i,j)*vmask_wet(i,j)
2749# endif
2750# ifdef DIAGNOSTICS_UV
2751!! DiaRVfrc(i,j,3,M2sstr)=cff1
2752!! DiaRVfrc(i,j,3,M2bstr)=cff2
2753# endif
2754 END DO
2755# endif
2756 END IF
2757 END DO j_loop
2758!
2759 RETURN
2760 END SUBROUTINE tl_rhs3d_tile
2761#endif
2762 END MODULE tl_rhs3d_mod
type(t_clima), dimension(:), allocatable clima
Definition mod_clima.F:153
type(t_coupling), dimension(:), allocatable coupling
type(t_forces), dimension(:), allocatable forces
Definition mod_forces.F:554
logical, dimension(:), allocatable lweakrelax
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_mixing), dimension(:), allocatable mixing
Definition mod_mixing.F:399
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, parameter itlm
Definition mod_param.F:663
integer, dimension(:), allocatable levbfrc
logical, dimension(:), allocatable ewperiodic
integer, parameter iwest
logical, dimension(:), allocatable nsperiodic
logical, dimension(:), allocatable lnudgem3clm
logical, dimension(:,:), allocatable compositegrid
integer, parameter isouth
integer, parameter ieast
integer, parameter inorth
integer, dimension(:), allocatable levsfrc
integer, dimension(:), allocatable nrhs
subroutine, public tl_pre_step3d(ng, tile)
subroutine, public tl_prsgrd(ng, tile)
Definition tl_prsgrd31.h:35
subroutine, public tl_rhs3d(ng, tile)
Definition tl_rhs3d.F:29
subroutine tl_rhs3d_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, hz, tl_hz, huon, tl_huon, hvom, tl_hvom, dmde, dndx, fomn, om_u, om_v, on_u, on_v, pm, pn, umask_wet, vmask_wet, bustr, tl_bustr, bvstr, tl_bvstr, sustr, tl_sustr, svstr, tl_svstr, u, tl_u, v, tl_v, w, tl_w, tl_rufrc, tl_rvfrc, tl_ru, tl_rv)
Definition tl_rhs3d.F:259
subroutine, public tl_t3dmix2(ng, tile)
subroutine, public tl_t3dmix4(ng, tile)
subroutine, public tl_t3drelax(ng, tile)
Definition tl_t3drelax.F:27
subroutine, public tl_uv3dmix2(ng, tile)
subroutine, public tl_uv3dmix4(ng, tile)
subroutine, public tl_uv3drelax(ng, tile)
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