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