ROMS
Loading...
Searching...
No Matches
tl_uv3dmix2_s.h
Go to the documentation of this file.
1 MODULE tl_uv3dmix2_mod
2!
3!git $Id$
4!================================================== Hernan G. Arango ===
5! Copyright (c) 2002-2025 The ROMS Group Andrew M. Moore !
6! Licensed under a MIT/X style license !
7! See License_ROMS.md !
8!=======================================================================
9! !
10! This routine computes tangent linear harmonic mixing of momentum, !
11! along constant S-surfaces, from the horizontal divergence of the !
12! stress tensor. A transverse isotropy is assumed so the stress !
13! tensor is splitted into vertical and horizontal subtensors. !
14! !
15! Reference: !
16! !
17! Wajsowicz, R.C, 1993: A consistent formulation of the !
18! anisotropic stress tensor for use in models of the !
19! large-scale ocean circulation, JCP, 105, 333-338. !
20! !
21! Sadourny, R. and K. Maynard, 1997: Formulations of !
22! lateral diffusion in geophysical fluid dynamics !
23! models, In Numerical Methods of Atmospheric and !
24! Oceanic Modelling. Lin, Laprise, and Ritchie, !
25! Eds., NRC Research Press, 547-556. !
26! !
27! Griffies, S.M. and R.W. Hallberg, 2000: Biharmonic !
28! friction with a Smagorinsky-like viscosity for !
29! use in large-scale eddy-permitting ocean models, !
30! Monthly Weather Rev., 128, 8, 2935-2946. !
31! !
32! BASIC STATE variables needed: visc2, u, v, Hz !
33! !
34!=======================================================================
35!
36 implicit none
37!
38 PRIVATE
39 PUBLIC tl_uv3dmix2
40!
41 CONTAINS
42!
43!***********************************************************************
44 SUBROUTINE tl_uv3dmix2 (ng, tile)
45!***********************************************************************
46!
47 USE mod_param
48 USE mod_coupling
49!!#ifdef DIAGNOSTICS_UV
50!! USE mod_diags
51!!#endif
52 USE mod_grid
53 USE mod_mixing
54 USE mod_ocean
55 USE mod_stepping
56!
57! Imported variable declarations.
58!
59 integer, intent(in) :: ng, tile
60!
61! Local variable declarations.
62!
63 character (len=*), parameter :: MyFile = &
64 & __FILE__
65!
66#include "tile.h"
67!
68#ifdef PROFILE
69 CALL wclock_on (ng, itlm, 30, __line__, myfile)
70#endif
71 CALL tl_uv3dmix2_s_tile (ng, tile, &
72 & lbi, ubi, lbj, ubj, &
73 & imins, imaxs, jmins, jmaxs, &
74 & nrhs(ng), nnew(ng), &
75#ifdef MASKING
76 & grid(ng) % pmask, &
77#endif
78 & grid(ng) % Hz, &
79 & grid(ng) % tl_Hz, &
80 & grid(ng) % om_p, &
81 & grid(ng) % om_r, &
82 & grid(ng) % on_p, &
83 & grid(ng) % on_r, &
84 & grid(ng) % pm, &
85 & grid(ng) % pmon_p, &
86 & grid(ng) % pmon_r, &
87 & grid(ng) % pn, &
88 & grid(ng) % pnom_p, &
89 & grid(ng) % pnom_r, &
90 & mixing(ng) % visc2_p, &
91 & mixing(ng) % visc2_r, &
92!!#ifdef DIAGNOSTICS_UV
93!! & DIAGS(ng) % DiaRUfrc, &
94!! & DIAGS(ng) % DiaRVfrc, &
95!! & DIAGS(ng) % DiaU3wrk, &
96!! & DIAGS(ng) % DiaV3wrk, &
97!!#endif
98 & ocean(ng) % u, &
99 & ocean(ng) % v, &
100 & coupling(ng) % tl_rufrc, &
101 & coupling(ng) % tl_rvfrc, &
102 & ocean(ng) % tl_u, &
103 & ocean(ng) % tl_v)
104#ifdef PROFILE
105 CALL wclock_off (ng, itlm, 30, __line__, myfile)
106#endif
107!
108 RETURN
109 END SUBROUTINE tl_uv3dmix2
110
111!
112!***********************************************************************
113 SUBROUTINE tl_uv3dmix2_s_tile (ng, tile, &
114 & LBi, UBi, LBj, UBj, &
115 & IminS, ImaxS, JminS, JmaxS, &
116 & nrhs, nnew, &
117#ifdef MASKING
118 & pmask, &
119#endif
120 & Hz, tl_Hz, &
121 & om_p, om_r, on_p, on_r, &
122 & pm, pmon_p, pmon_r, &
123 & pn, pnom_p, pnom_r, &
124 & visc2_p, visc2_r, &
125!!#ifdef DIAGNOSTICS_UV
126!! & DiaRUfrc, DiaRVfrc, &
127!! & DiaU3wrk, DiaV3wrk, &
128!!#endif
129 & u, v, &
130 & tl_rufrc, tl_rvfrc, &
131 & tl_u, tl_v)
132!***********************************************************************
133!
134 USE mod_param
135 USE mod_scalars
136!
137! Imported variable declarations.
138!
139 integer, intent(in) :: ng, tile
140 integer, intent(in) :: LBi, UBi, LBj, UBj
141 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
142 integer, intent(in) :: nrhs, nnew
143
144#ifdef ASSUMED_SHAPE
145# ifdef MASKING
146 real(r8), intent(in) :: pmask(LBi:,LBj:)
147# endif
148 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
149 real(r8), intent(in) :: tl_Hz(LBi:,LBj:,:)
150 real(r8), intent(in) :: om_p(LBi:,LBj:)
151 real(r8), intent(in) :: om_r(LBi:,LBj:)
152 real(r8), intent(in) :: on_p(LBi:,LBj:)
153 real(r8), intent(in) :: on_r(LBi:,LBj:)
154 real(r8), intent(in) :: pm(LBi:,LBj:)
155 real(r8), intent(in) :: pmon_p(LBi:,LBj:)
156 real(r8), intent(in) :: pmon_r(LBi:,LBj:)
157 real(r8), intent(in) :: pn(LBi:,LBj:)
158 real(r8), intent(in) :: pnom_p(LBi:,LBj:)
159 real(r8), intent(in) :: pnom_r(LBi:,LBj:)
160 real(r8), intent(in) :: visc2_p(LBi:,LBj:)
161 real(r8), intent(in) :: visc2_r(LBi:,LBj:)
162
163 real(r8), intent(in) :: u(LBi:,LBj:,:,:)
164 real(r8), intent(in) :: v(LBi:,LBj:,:,:)
165
166 real(r8), intent(inout) :: tl_rufrc(LBi:,LBj:)
167 real(r8), intent(inout) :: tl_rvfrc(LBi:,LBj:)
168 real(r8), intent(inout) :: tl_u(LBi:,LBj:,:,:)
169 real(r8), intent(inout) :: tl_v(LBi:,LBj:,:,:)
170#else
171# ifdef MASKING
172 real(r8), intent(in) :: pmask(LBi:UBi,LBj:UBj)
173# endif
174 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
175 real(r8), intent(in) :: tl_Hz(LBi:UBi,LBj:UBj,N(ng))
176 real(r8), intent(in) :: om_p(LBi:UBi,LBj:UBj)
177 real(r8), intent(in) :: om_r(LBi:UBi,LBj:UBj)
178 real(r8), intent(in) :: on_p(LBi:UBi,LBj:UBj)
179 real(r8), intent(in) :: on_r(LBi:UBi,LBj:UBj)
180 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
181 real(r8), intent(in) :: pmon_p(LBi:UBi,LBj:UBj)
182 real(r8), intent(in) :: pmon_r(LBi:UBi,LBj:UBj)
183 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
184 real(r8), intent(in) :: pnom_p(LBi:UBi,LBj:UBj)
185 real(r8), intent(in) :: pnom_r(LBi:UBi,LBj:UBj)
186 real(r8), intent(in) :: visc2_p(LBi:UBi,LBj:UBj)
187 real(r8), intent(in) :: visc2_r(LBi:UBi,LBj:UBj)
188
189 real(r8), intent(in) :: u(LBi:UBi,LBj:UBj,N(ng),2)
190 real(r8), intent(in) :: v(LBi:UBi,LBj:UBj,N(ng),2)
191
192 real(r8), intent(inout) :: tl_rufrc(LBi:UBi,LBj:UBj)
193 real(r8), intent(inout) :: tl_rvfrc(LBi:UBi,LBj:UBj)
194 real(r8), intent(inout) :: tl_u(LBi:UBi,LBj:UBj,N(ng),2)
195 real(r8), intent(inout) :: tl_v(LBi:UBi,LBj:UBj,N(ng),2)
196#endif
197!
198! Local variable declarations.
199!
200 integer :: i, j, k
201
202 real(r8) :: cff, tl_cff, tl_cff1, tl_cff2
203
204 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_UFe
205 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_VFe
206 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_UFx
207 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_VFx
208
209#include "set_bounds.h"
210!
211
212!-----------------------------------------------------------------------
213! Compute horizontal harmonic viscosity along constant S-surfaces.
214!-----------------------------------------------------------------------
215!
216 k_loop : DO k=1,n(ng)
217!
218! Compute flux-components of the horizontal divergence of the stress
219! tensor (m5/s2) in XI- and ETA-directions.
220!
221 DO j=jstrv-1,jend
222 DO i=istru-1,iend
223!^ cff=visc2_r(i,j)*Hz(i,j,k)*0.5_r8* &
224!^ & (pmon_r(i,j)* &
225!^ & ((pn(i ,j)+pn(i+1,j))*u(i+1,j,k,nrhs)- &
226!^ & (pn(i-1,j)+pn(i ,j))*u(i ,j,k,nrhs))- &
227!^ & pnom_r(i,j)* &
228!^ & ((pm(i,j )+pm(i,j+1))*v(i,j+1,k,nrhs)- &
229!^ & (pm(i,j-1)+pm(i,j ))*v(i,j ,k,nrhs)))
230!^
231 tl_cff=visc2_r(i,j)*0.5_r8* &
232 & (tl_hz(i,j,k)* &
233 & (pmon_r(i,j)* &
234 & ((pn(i ,j)+pn(i+1,j))*u(i+1,j,k,nrhs)- &
235 & (pn(i-1,j)+pn(i ,j))*u(i ,j,k,nrhs))- &
236 & pnom_r(i,j)* &
237 & ((pm(i,j )+pm(i,j+1))*v(i,j+1,k,nrhs)- &
238 & (pm(i,j-1)+pm(i,j ))*v(i,j ,k,nrhs)))+ &
239 & hz(i,j,k)* &
240 & (pmon_r(i,j)* &
241 & ((pn(i ,j)+pn(i+1,j))*tl_u(i+1,j,k,nrhs)- &
242 & (pn(i-1,j)+pn(i ,j))*tl_u(i ,j,k,nrhs))- &
243 & pnom_r(i,j)* &
244 & ((pm(i,j )+pm(i,j+1))*tl_v(i,j+1,k,nrhs)- &
245 & (pm(i,j-1)+pm(i,j ))*tl_v(i,j ,k,nrhs))))
246!^ UFx(i,j)=on_r(i,j)*on_r(i,j)*cff
247!^
248 tl_ufx(i,j)=on_r(i,j)*on_r(i,j)*tl_cff
249!^ VFe(i,j)=om_r(i,j)*om_r(i,j)*cff
250!^
251 tl_vfe(i,j)=om_r(i,j)*om_r(i,j)*tl_cff
252 END DO
253 END DO
254 DO j=jstr,jend+1
255 DO i=istr,iend+1
256!^ cff=visc2_p(i,j)*0.125_r8*(Hz(i-1,j ,k)+Hz(i,j ,k)+ &
257!^ & Hz(i-1,j-1,k)+Hz(i,j-1,k))* &
258!^ & (pmon_p(i,j)* &
259!^ & ((pn(i ,j-1)+pn(i ,j))*v(i ,j,k,nrhs)- &
260!^ & (pn(i-1,j-1)+pn(i-1,j))*v(i-1,j,k,nrhs))+ &
261!^ & pnom_p(i,j)* &
262!^ & ((pm(i-1,j )+pm(i,j ))*u(i,j ,k,nrhs)- &
263!^ & (pm(i-1,j-1)+pm(i,j-1))*u(i,j-1,k,nrhs)))
264!^
265 tl_cff=visc2_p(i,j)*0.125_r8* &
266 & ((tl_hz(i-1,j ,k)+tl_hz(i,j ,k)+ &
267 & tl_hz(i-1,j-1,k)+tl_hz(i,j-1,k))* &
268 & (pmon_p(i,j)* &
269 & ((pn(i ,j-1)+pn(i ,j))*v(i ,j,k,nrhs)- &
270 & (pn(i-1,j-1)+pn(i-1,j))*v(i-1,j,k,nrhs))+ &
271 & pnom_p(i,j)* &
272 & ((pm(i-1,j )+pm(i,j ))*u(i,j ,k,nrhs)- &
273 & (pm(i-1,j-1)+pm(i,j-1))*u(i,j-1,k,nrhs)))+ &
274 & (hz(i-1,j ,k)+hz(i,j ,k)+ &
275 & hz(i-1,j-1,k)+hz(i,j-1,k))* &
276 & (pmon_p(i,j)* &
277 & ((pn(i ,j-1)+pn(i ,j))*tl_v(i ,j,k,nrhs)- &
278 & (pn(i-1,j-1)+pn(i-1,j))*tl_v(i-1,j,k,nrhs))+ &
279 & pnom_p(i,j)* &
280 & ((pm(i-1,j )+pm(i,j ))*tl_u(i,j ,k,nrhs)- &
281 & (pm(i-1,j-1)+pm(i,j-1))*tl_u(i,j-1,k,nrhs))))
282#ifdef MASKING
283!^ cff=cff*pmask(i,j)
284!^
285 tl_cff=tl_cff*pmask(i,j)
286#endif
287!^ UFe(i,j)=om_p(i,j)*om_p(i,j)*cff
288!^
289 tl_ufe(i,j)=om_p(i,j)*om_p(i,j)*tl_cff
290!^ VFx(i,j)=on_p(i,j)*on_p(i,j)*cff
291!^
292 tl_vfx(i,j)=on_p(i,j)*on_p(i,j)*tl_cff
293 END DO
294 END DO
295!
296! Time-step harmonic, S-surfaces viscosity term. Notice that momentum
297! at this stage is HzU and HzV and has m2/s units. Add contribution for
298! barotropic forcing terms.
299!
300 DO j=jstr,jend
301 DO i=istru,iend
302 cff=0.25_r8*(pm(i-1,j)+pm(i,j))*(pn(i-1,j)+pn(i,j))
303!^ cff1=0.5_r8*((pn(i-1,j)+pn(i,j))* &
304!^ & (UFx(i,j )-UFx(i-1,j))+ &
305!^ & (pm(i-1,j)+pm(i,j))* &
306!^ & (UFe(i,j+1)-UFe(i ,j)))
307!^
308 tl_cff1=0.5_r8*((pn(i-1,j)+pn(i,j))* &
309 & (tl_ufx(i,j )-tl_ufx(i-1,j))+ &
310 & (pm(i-1,j)+pm(i,j))* &
311 & (tl_ufe(i,j+1)-tl_ufe(i ,j)))
312!^ cff2=dt(ng)*cff*cff1
313!^
314 tl_cff2=dt(ng)*cff*tl_cff1
315!^ rufrc(i,j)=rufrc(i,j)+cff1
316!^
317 tl_rufrc(i,j)=tl_rufrc(i,j)+tl_cff1
318
319!^ u(i,j,k,nnew)=u(i,j,k,nnew)+cff2
320!!#ifdef DIAGNOSTICS_UV
321!! DiaRUfrc(i,j,3,M2hvis)=DiaRUfrc(i,j,3,M2hvis)+cff1
322!! DiaU3wrk(i,j,k,M3hvis)=cff2
323!!#endif
324!^
325 tl_u(i,j,k,nnew)=tl_u(i,j,k,nnew)+tl_cff2
326 END DO
327 END DO
328 DO j=jstrv,jend
329 DO i=istr,iend
330 cff=0.25_r8*(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1))
331!^ cff1=0.5_r8*((pn(i,j-1)+pn(i,j))* &
332!^ & (VFx(i+1,j)-VFx(i,j ))- &
333!^ & (pm(i,j-1)+pm(i,j))* &
334!^ & (VFe(i ,j)-VFe(i,j-1)))
335!^
336 tl_cff1=0.5_r8*((pn(i,j-1)+pn(i,j))* &
337 & (tl_vfx(i+1,j)-tl_vfx(i,j ))- &
338 & (pm(i,j-1)+pm(i,j))* &
339 & (tl_vfe(i ,j)-tl_vfe(i,j-1)))
340!^ cff2=dt(ng)*cff*cff1
341!^
342 tl_cff2=dt(ng)*cff*tl_cff1
343!^ rvfrc(i,j)=rvfrc(i,j)+cff1
344!^
345 tl_rvfrc(i,j)=tl_rvfrc(i,j)+tl_cff1
346!^ v(i,j,k,nnew)=v(i,j,k,nnew)+cff2
347!!#ifdef DIAGNOSTICS_UV
348!! DiaRVfrc(i,j,3,M2hvis)=DiaRVfrc(i,j,3,M2hvis)+cff1
349!! DiaV3wrk(i,j,k,M3hvis)=cff2
350!!#endif
351!^
352 tl_v(i,j,k,nnew)=tl_v(i,j,k,nnew)+tl_cff2
353 END DO
354 END DO
355 END DO k_loop
356!
357 RETURN
358 END SUBROUTINE tl_uv3dmix2_s_tile
359
360 END MODULE tl_uv3dmix2_mod
type(t_coupling), dimension(:), allocatable coupling
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
integer, parameter itlm
Definition mod_param.F:663
real(dp), dimension(:), allocatable dt
integer, dimension(:), allocatable nrhs
integer, dimension(:), allocatable nnew
subroutine tl_uv3dmix2_s_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, nnew, pmask, hz, tl_hz, om_p, om_r, on_p, on_r, pm, pmon_p, pmon_r, pn, pnom_p, pnom_r, visc2_p, visc2_r, u, v, tl_rufrc, tl_rvfrc, tl_u, tl_v)
subroutine, public tl_uv3dmix2(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