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

Functions/Subroutines

subroutine, public rp_set_massflux (ng, tile, model)
 
subroutine rp_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)
 

Function/Subroutine Documentation

◆ rp_set_massflux()

subroutine, public rp_set_massflux_mod::rp_set_massflux ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model )

Definition at line 30 of file rp_set_massflux.F.

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 rp_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
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, dimension(:), allocatable nrhs
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3

References mod_grid::grid, mod_stepping::nrhs, mod_ocean::ocean, rp_set_massflux_tile(), wclock_off(), and wclock_on().

Referenced by rp_initial(), and rp_main3d().

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

◆ rp_set_massflux_tile()

subroutine rp_set_massflux_mod::rp_set_massflux_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) nrhs,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) u,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) v,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) tl_u,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) tl_v,
real(r8), dimension(lbi:,lbj:,:), intent(in) hz,
real(r8), dimension(lbi:,lbj:,:), intent(in) tl_hz,
real(r8), dimension(lbi:,lbj:), intent(in) om_v,
real(r8), dimension(lbi:,lbj:), intent(in) on_u,
real(r8), dimension(lbi:,lbj:,:), intent(out) tl_huon,
real(r8), dimension(lbi:,lbj:,:), intent(out) tl_hvom )
private

Definition at line 80 of file rp_set_massflux.F.

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 TL_IOMS
172 & 0.5_r8*on_u(i,j)* &
173 & (hz(i,j,k)+hz(i-1,j,k))*u(i,j,k,nrhs)
174# endif
175# ifdef WEC_MELLOR
176!^ Huon(i,j,k)=Huon(i,j,k)+ &
177!^ & 0.5_r8*(Hz(i,j,k)+Hz(i-1,j,k))* &
178!^ & u_stokes(i,j,k)*on_u(i,j)
179!^
180 tl_huon(i,j,k)=tl_huon(i,j,k)+ &
181 & 0.5_r8*on_u(i,j)* &
182 & ((hz(i,j,k)+hz(i-1,j,k))* &
183 & tl_u_stokes(i,j,k)+ &
184 & (tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
185 & u_stokes(i,j,k))- &
186# ifdef TL_IOMS
187 & 0.5_r8*on_u(i,j)* &
188 & (hz(i,j,k)+hz(i-1,j,k))*u_stokes(i,j,k)
189# endif
190# endif
191 END DO
192 END DO
193 DO j=jstrp,jendt
194 DO i=istrt,iendt
195!^ Hvom(i,j,k)=0.5_r8*(Hz(i,j,k)+Hz(i,j-1,k))*v(i,j,k,nrhs)* &
196!^ & om_v(i,j)
197!^
198 tl_hvom(i,j,k)=0.5_r8*om_v(i,j)* &
199 & ((hz(i,j,k)+hz(i,j-1,k))* &
200 & tl_v(i,j,k,nrhs)+ &
201 & (tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
202 & v(i,j,k,nrhs))- &
203# ifdef TL_IOMS
204 & 0.5_r8*om_v(i,j)* &
205 & (hz(i,j,k)+hz(i,j-1,k))*v(i,j,k,nrhs)
206# endif
207# ifdef WEC_MELLOR
208!^ Hvom(i,j,k)=Hvom(i,j,k)+ &
209!^ & 0.5_r8*(Hz(i,j,k)+Hz(i,j-1,k))* &
210!^ & v_stokes(i,j,k)*om_v(i,j)
211!^
212 tl_hvom(i,j,k)=tl_hvom(i,j,k)+ &
213 & 0.5_r8*om_v(i,j)* &
214 & ((hz(i,j,k)+hz(i,j-1,k))* &
215 & tl_v_stokes(i,j,k)+ &
216 & (tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
217 & v_stokes(i,j,k))- &
218# ifdef TL_IOMS
219 & 0.5_r8*om_v(i,j)* &
220 & (hz(i,j,k)+hz(i,j-1,k))*v_stokes(i,j,k)
221# endif
222# endif
223 END DO
224 END DO
225 END DO
226!
227! Exchange boundary information.
228!
229 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
230!^ CALL exchange_u3d_tile (ng, tile, &
231!^ & LBi, UBi, LBj, UBj, 1, N(ng), &
232!^ & Huon)
233!^
234 CALL exchange_u3d_tile (ng, tile, &
235 & lbi, ubi, lbj, ubj, 1, n(ng), &
236 & tl_huon)
237!^ CALL exchange_v3d_tile (ng, tile, &
238!^ & LBi, UBi, LBj, UBj, 1, N(ng), &
239!^ & Hvom)
240!^
241 CALL exchange_v3d_tile (ng, tile, &
242 & lbi, ubi, lbj, ubj, 1, n(ng), &
243 & tl_hvom)
244 END IF
245
246# ifdef DISTRIBUTE
247!^ CALL mp_exchange3d (ng, tile, model, 2, &
248!^ & LBi, UBi, LBj, UBj, 1, N(ng), &
249!^ & NghostPoints, &
250!^ & EWperiodic(ng), NSperiodic(ng), &
251!^ & Huon, Hvom)
252!^
253 CALL mp_exchange3d (ng, tile, model, 2, &
254 & lbi, ubi, lbj, ubj, 1, n(ng), &
255 & nghostpoints, &
256 & ewperiodic(ng), nsperiodic(ng), &
257 & tl_huon, tl_hvom)
258# endif
259!
260 RETURN
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)
integer, dimension(:), allocatable n
Definition mod_param.F:479
integer nghostpoints
Definition mod_param.F:710
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
subroutine mp_exchange3d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, lbk, ubk, nghost, ew_periodic, ns_periodic, a, b, c, d)

References mod_scalars::ewperiodic, exchange_3d_mod::exchange_u3d_tile(), exchange_3d_mod::exchange_v3d_tile(), mp_exchange_mod::mp_exchange3d(), mod_param::nghostpoints, and mod_scalars::nsperiodic.

Referenced by rp_set_massflux().

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