ROMS
Loading...
Searching...
No Matches
tl_set_massflux.F
Go to the documentation of this file.
1#include "cppdefs.h"
3
4#if defined TANGENT && defined SOLVE3D
5!
6!git $Id$
7!================================================== Hernan G. Arango ===
8! Copyright (c) 2002-2025 The ROMS Group Andrew M. Moore !
9! Licensed under a MIT/X style license !
10! See License_ROMS.md !
11!=======================================================================
12! !
13! This routine computes tangent linear horizontal mass fluxes, !
14! Hz*u/n and Hz*v/m. !
15! !
16! BASIC STATE variables required: Hz, u, v !
17! Dependend variables: tl_Huon, tl_Hvom !
18! Independend variables: tl_Hz, tl_u, tl_v !
19! !
20!=======================================================================
21!
22 implicit none
23!
24 PRIVATE
25 PUBLIC :: tl_set_massflux
26!
27 CONTAINS
28!
29!***********************************************************************
30 SUBROUTINE tl_set_massflux (ng, tile, model)
31!***********************************************************************
32!
33 USE mod_param
34 USE mod_grid
35 USE mod_ocean
36 USE mod_stepping
37!
38! Imported variable declarations.
39!
40 integer, intent(in) :: ng, tile, model
41!
42! Local variable declarations.
43!
44 character (len=*), parameter :: myfile = &
45 & __FILE__
46!
47# include "tile.h"
48!
49# ifdef PROFILE
50 CALL wclock_on (ng, model, 12, __line__, myfile)
51# endif
52 CALL tl_set_massflux_tile (ng, tile, model, &
53 & lbi, ubi, lbj, ubj, &
54 & imins, imaxs, jmins, jmaxs, &
55 & nrhs(ng), &
56 & ocean(ng) % u, &
57 & ocean(ng) % v, &
58 & ocean(ng) % tl_u, &
59 & ocean(ng) % tl_v, &
60# ifdef WEC_MELLOR
61 & ocean(ng) % u_stokes, &
62 & ocean(ng) % v_stokes, &
63 & ocean(ng) % tl_u_stokes, &
64 & ocean(ng) % tl_v_stokes, &
65# endif
66 & grid(ng) % Hz, &
67 & grid(ng) % tl_Hz, &
68 & grid(ng) % om_v, &
69 & grid(ng) % on_u, &
70 & grid(ng) % tl_Huon, &
71 & grid(ng) % tl_Hvom)
72# ifdef PROFILE
73 CALL wclock_off (ng, model, 12, __line__, myfile)
74# endif
75!
76 RETURN
77 END SUBROUTINE tl_set_massflux
78!
79!***********************************************************************
80 SUBROUTINE tl_set_massflux_tile (ng, tile, model, &
81 & LBi, UBi, LBj, UBj, &
82 & IminS, ImaxS, JminS, JmaxS, &
83 & nrhs, &
84 & u, v, &
85 & tl_u, tl_v, &
86# ifdef WEC_MELLOR
87 & u_stokes, v_stokes, &
88 & tl_u_stokes, tl_v_stokes, &
89# endif
90 & Hz, tl_Hz, &
91 & om_v, on_u, &
92 & tl_Huon, tl_Hvom)
93!***********************************************************************
94!
95 USE mod_param
96 USE mod_scalars
97!
99# ifdef DISTRIBUTE
101# endif
102!
103! Imported variable declarations.
104!
105 integer, intent(in) :: ng, tile, model
106 integer, intent(in) :: LBi, UBi, LBj, UBj
107 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
108 integer, intent(in) :: nrhs
109!
110# ifdef ASSUMED_SHAPE
111 real(r8), intent(in) :: u(LBi:,LBj:,:,:)
112 real(r8), intent(in) :: v(LBi:,LBj:,:,:)
113 real(r8), intent(in) :: tl_u(LBi:,LBj:,:,:)
114 real(r8), intent(in) :: tl_v(LBi:,LBj:,:,:)
115# ifdef WEC_MELLOR
116 real(r8), intent(in) :: u_stokes(LBi:,LBj:,:)
117 real(r8), intent(in) :: v_stokes(LBi:,LBj:,:)
118 real(r8), intent(in) :: tl_u_stokes(LBi:,LBj:,:)
119 real(r8), intent(in) :: tl_v_stokes(LBi:,LBj:,:)
120# endif
121 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
122 real(r8), intent(in) :: tl_Hz(LBi:,LBj:,:)
123 real(r8), intent(in) :: om_v(LBi:,LBj:)
124 real(r8), intent(in) :: on_u(LBi:,LBj:)
125
126 real(r8), intent(out) :: tl_Huon(LBi:,LBj:,:)
127 real(r8), intent(out) :: tl_Hvom(LBi:,LBj:,:)
128# else
129 real(r8), intent(in) :: u(LBi:UBi,LBj:UBj,N(ng),2)
130 real(r8), intent(in) :: v(LBi:UBi,LBj:UBj,N(ng),2)
131 real(r8), intent(in) :: tl_u(LBi:UBi,LBj:UBj,N(ng),2)
132 real(r8), intent(in) :: tl_v(LBi:UBi,LBj:UBj,N(ng),2)
133# ifdef WEC_MELLOR
134 real(r8), intent(in) :: u_stokes(LBi:UBi,LBj:UBj,N(ng))
135 real(r8), intent(in) :: v_stokes(LBi:UBi,LBj:UBj,N(ng))
136 real(r8), intent(in) :: tl_u_stokes(LBi:UBi,LBj:UBj,N(ng))
137 real(r8), intent(in) :: tl_v_stokes(LBi:UBi,LBj:UBj,N(ng))
138# endif
139 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
140 real(r8), intent(in) :: tl_Hz(LBi:UBi,LBj:UBj,N(ng))
141 real(r8), intent(in) :: om_v(LBi:UBi,LBj:UBj)
142 real(r8), intent(in) :: on_u(LBi:UBi,LBj:UBj)
143
144 real(r8), intent(out) :: tl_Huon(LBi:UBi,LBj:UBj,N(ng))
145 real(r8), intent(out) :: tl_Hvom(LBi:UBi,LBj:UBj,N(ng))
146# endif
147!
148! Local variable declarations.
149!
150 integer :: i, j, k
151
152# include "set_bounds.h"
153!
154!-----------------------------------------------------------------------
155! Compute horizontal mass fluxes, Hz*u/n and Hz*v/m.
156!-----------------------------------------------------------------------
157!
158! Compute horizontal mass fluxes.
159!
160 DO k=1,n(ng)
161 DO j=jstrt,jendt
162 DO i=istrp,iendt
163!^ Huon(i,j,k)=0.5_r8*(Hz(i,j,k)+Hz(i-1,j,k))*u(i,j,k,nrhs)* &
164!^ & on_u(i,j)
165!^
166 tl_huon(i,j,k)=0.5_r8*on_u(i,j)* &
167 & ((hz(i,j,k)+hz(i-1,j,k))* &
168 & tl_u(i,j,k,nrhs)+ &
169 & (tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
170 & u(i,j,k,nrhs))
171# ifdef WEC_MELLOR
172!^ Huon(i,j,k)=Huon(i,j,k)+ &
173!^ & 0.5_r8*(Hz(i,j,k)+Hz(i-1,j,k))* &
174!^ & u_stokes(i,j,k)*on_u(i,j)
175!^
176 tl_huon(i,j,k)=tl_huon(i,j,k)+ &
177 & 0.5_r8*on_u(i,j)* &
178 & ((hz(i,j,k)+hz(i-1,j,k))* &
179 & tl_u_stokes(i,j,k)+ &
180 & (tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
181 & u_stokes(i,j,k))
182# endif
183 END DO
184 END DO
185 DO j=jstrp,jendt
186 DO i=istrt,iendt
187!^ Hvom(i,j,k)=0.5_r8*(Hz(i,j,k)+Hz(i,j-1,k))*v(i,j,k,nrhs)* &
188!^ & om_v(i,j)
189!^
190 tl_hvom(i,j,k)=0.5_r8*om_v(i,j)* &
191 & ((hz(i,j,k)+hz(i,j-1,k))* &
192 & tl_v(i,j,k,nrhs)+ &
193 & (tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
194 & v(i,j,k,nrhs))
195# ifdef WEC_MELLOR
196!^ Hvom(i,j,k)=Hvom(i,j,k)+ &
197!^ & 0.5_r8*(Hz(i,j,k)+Hz(i,j-1,k))* &
198!^ & v_stokes(i,j,k)*om_v(i,j)
199!^
200 tl_hvom(i,j,k)=tl_hvom(i,j,k)+ &
201 & 0.5_r8*om_v(i,j)* &
202 & ((hz(i,j,k)+hz(i,j-1,k))* &
203 & tl_v_stokes(i,j,k)+ &
204 & (tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
205 & v_stokes(i,j,k))
206# endif
207 END DO
208 END DO
209 END DO
210!
211! Exchange boundary information.
212!
213 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
214!^ CALL exchange_u3d_tile (ng, tile, &
215!^ & LBi, UBi, LBj, UBj, 1, N(ng), &
216!^ & Huon)
217!^
218 CALL exchange_u3d_tile (ng, tile, &
219 & lbi, ubi, lbj, ubj, 1, n(ng), &
220 & tl_huon)
221!^ CALL exchange_v3d_tile (ng, tile, &
222!^ & LBi, UBi, LBj, UBj, 1, N(ng), &
223!^ & Hvom)
224!^
225 CALL exchange_v3d_tile (ng, tile, &
226 & lbi, ubi, lbj, ubj, 1, n(ng), &
227 & tl_hvom)
228 END IF
229
230# ifdef DISTRIBUTE
231!^ CALL mp_exchange3d (ng, tile, model, 2, &
232!^ & LBi, UBi, LBj, UBj, 1, N(ng), &
233!^ & NghostPoints, &
234!^ & EWperiodic(ng), NSperiodic(ng), &
235!^ & Huon, Hvom)
236!^
237 CALL mp_exchange3d (ng, tile, model, 2, &
238 & lbi, ubi, lbj, ubj, 1, n(ng), &
239 & nghostpoints, &
240 & ewperiodic(ng), nsperiodic(ng), &
241 & tl_huon, tl_hvom)
242# endif
243!
244 RETURN
245 END SUBROUTINE tl_set_massflux_tile
246#endif
247 END MODULE tl_set_massflux_mod
248
subroutine exchange_u3d_tile(ng, tile, lbi, ubi, lbj, ubj, lbk, ubk, a)
subroutine exchange_v3d_tile(ng, tile, lbi, ubi, lbj, ubj, lbk, ubk, a)
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer nghostpoints
Definition mod_param.F:710
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
integer, dimension(:), allocatable nrhs
subroutine mp_exchange3d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, lbk, ubk, nghost, ew_periodic, ns_periodic, a, b, c, d)
subroutine tl_set_massflux_tile(ng, tile, model, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, u, v, tl_u, tl_v, hz, tl_hz, om_v, on_u, tl_huon, tl_hvom)
subroutine, public tl_set_massflux(ng, tile, model)
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