ROMS
Loading...
Searching...
No Matches
ad_wvelocity.F
Go to the documentation of this file.
1#include "cppdefs.h"
3#if defined ADJOINT && defined SOLVE3D
4!
5!git $Id$
6!=======================================================================
7! Copyright (c) 2002-2025 The ROMS Group !
8! Licensed under a MIT/X style license !
9! See License_ROMS.md Hernan G. Arango !
10!=================================================== Andrew M. Moore ===
11! !
12! This subroutines computes vertical velocity (m/s) at W-points !
13! from the vertical mass flux (omega*hz/m*n). This computation !
14! is done solely for output purposes. !
15! !
16!=======================================================================
17!
18 implicit none
19
20 PRIVATE
21 PUBLIC :: ad_wvelocity
22
23 CONTAINS
24!
25!***********************************************************************
26 SUBROUTINE ad_wvelocity (ng, tile, Ninp)
27!***********************************************************************
28!
29 USE mod_param
30 USE mod_coupling
31 USE mod_grid
32 USE mod_ocean
33 USE mod_stepping
34!
35! Imported variable declarations.
36!
37 integer, intent(in) :: ng, tile, ninp
38!
39! Local variable declarations.
40!
41# include "tile.h"
42!
43 CALL ad_wvelocity_tile (ng, tile, &
44 & lbi, ubi, lbj, ubj, &
45 & imins, imaxs, jmins, jmaxs, &
46 & ninp, &
47 & grid(ng) % pm, &
48 & grid(ng) % pn, &
49 & grid(ng) % z_r, &
50 & grid(ng) % z_w, &
51 & grid(ng) % ad_z_r, &
52 & grid(ng) % ad_z_w, &
53 & coupling(ng) % DU_avg1, &
54 & coupling(ng) % DV_avg1, &
55 & coupling(ng) % ad_DU_avg1, &
56 & coupling(ng) % ad_DV_avg1, &
57 & ocean(ng) % u, &
58 & ocean(ng) % v, &
59 & ocean(ng) % W, &
60 & ocean(ng) % ad_u, &
61 & ocean(ng) % ad_v, &
62 & ocean(ng) % ad_W, &
63 & ocean(ng) % ad_wvel)
64 RETURN
65 END SUBROUTINE ad_wvelocity
66!
67!***********************************************************************
68 SUBROUTINE ad_wvelocity_tile (ng, tile, &
69 & LBi, UBi, LBj, UBj, &
70 & IminS, ImaxS, JminS, JmaxS, &
71 & Ninp, &
72 & pm, pn, z_r, z_w, &
73 & ad_z_r, ad_z_w, &
74 & DU_avg1, DV_avg1, &
75 & ad_DU_avg1, ad_DV_avg1, &
76 & u, v, W, &
77 & ad_u, ad_v, ad_W, &
78 & ad_wvel)
79!***********************************************************************
80!
81 USE mod_param
82 USE mod_scalars
83!
84 USE bc_3d_mod, ONLY : bc_w3d_tile
87# ifdef DISTRIBUTE
89# endif
90!
91! Imported variable declarations.
92!
93 integer, intent(in) :: ng, tile
94 integer, intent(in) :: LBi, UBi, LBj, UBj
95 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
96 integer, intent(in) :: Ninp
97!
98# ifdef ASSUMED_SHAPE
99 real(r8), intent(in) :: pm(LBi:,LBj:)
100 real(r8), intent(in) :: pn(LBi:,LBj:)
101 real(r8), intent(in) :: z_r(LBi:,LBj:,:)
102 real(r8), intent(in) :: z_w(LBi:,LBj:,0:)
103 real(r8), intent(in) :: DU_avg1(LBi:,LBj:)
104 real(r8), intent(in) :: DV_avg1(LBi:,LBj:)
105 real(r8), intent(in) :: u(LBi:,LBj:,:,:)
106 real(r8), intent(in) :: v(LBi:,LBj:,:,:)
107 real(r8), intent(in) :: W(LBi:,LBj:,0:)
108 real(r8), intent(inout) :: ad_z_r(LBi:,LBj:,:)
109 real(r8), intent(inout) :: ad_z_w(LBi:,LBj:,0:)
110 real(r8), intent(inout) :: ad_DU_avg1(LBi:,LBj:)
111 real(r8), intent(inout) :: ad_DV_avg1(LBi:,LBj:)
112 real(r8), intent(inout) :: ad_u(LBi:,LBj:,:,:)
113 real(r8), intent(inout) :: ad_v(LBi:,LBj:,:,:)
114 real(r8), intent(inout) :: ad_W(LBi:,LBj:,0:)
115 real(r8), intent(inout) :: ad_wvel(LBi:,LBj:,0:)
116# else
117 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
118 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
119 real(r8), intent(in) :: z_r(LBi:UBi,LBj:UBj,N(ng))
120 real(r8), intent(in) :: z_w(LBi:UBi,LBj:UBj,0:N(ng))
121 real(r8), intent(in) :: DU_avg1(LBi:UBi,LBj:UBj)
122 real(r8), intent(in) :: DV_avg1(LBi:UBi,LBj:UBj)
123 real(r8), intent(in) :: u(LBi:UBi,LBj:UBj,N(ng),2)
124 real(r8), intent(in) :: v(LBi:UBi,LBj:UBj,N(ng),2)
125 real(r8), intent(in) :: W(LBi:UBi,LBj:UBj,0:N(ng))
126 real(r8), intent(inout) :: ad_z_r(LBi:UBi,LBj:UBj,N(ng))
127 real(r8), intent(inout) :: ad_z_w(LBi:UBi,LBj:UBj,0:N(ng))
128 real(r8), intent(inout) :: ad_DU_avg1(LBi:UBi,LBj:UBj)
129 real(r8), intent(inout) :: ad_DV_avg1(LBi:UBi,LBj:UBj)
130 real(r8), intent(inout) :: ad_u(LBi:UBi,LBj:UBj,N(ng),2)
131 real(r8), intent(inout) :: ad_v(LBi:UBi,LBj:UBj,N(ng),2)
132 real(r8), intent(inout) :: ad_W(LBi:UBi,LBj:UBj,0:N(ng))
133 real(r8), intent(inout) :: ad_wvel(LBi:UBi,LBj:UBj,0:N(ng))
134# endif
135
136!
137! Local variable declarations.
138!
139 integer :: i, j, k
140
141 real(r8) :: cff1, cff2, cff3, cff4, cff5, slope , ad_slope
142 real(r8) :: adfac, adfac1, adfac2, adfac3
143
144 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,N(ng)) :: vert
145 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,N(ng)) :: ad_vert
146
147 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: wrk
148 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_wrk
149
150# include "set_bounds.h"
151!
152!-----------------------------------------------------------------------
153! Compute adjoint "true" vertical velocity (m/s).
154!-----------------------------------------------------------------------
155!
156# ifdef DISTRIBUTE
157 CALL ad_mp_exchange3d (ng, tile, iadm, 1, &
158 & lbi, ubi, lbj, ubj, 0, n(ng), &
159 & nghostpoints, &
160 & ewperiodic(ng), nsperiodic(ng), &
161 & ad_wvel)
162# endif
163!
164! Set lateral boundary conditions.
165!
166 CALL ad_bc_w3d_tile (ng, tile, &
167 & lbi, ubi, lbj, ubj, 0, n(ng), &
168 & ad_wvel)
169!
170! Initialize local adjoint variables and arrays.
171!
172 ad_slope=0.0_r8
173
174 DO j=jmins,jmaxs
175 DO i=imins,imaxs
176 ad_wrk(i,j)=0.0_r8
177 END DO
178 END DO
179 DO k=1,n(ng)
180 DO j=jmins,jmaxs
181 DO i=imins,imaxs
182 ad_vert(i,j,k)=0.0_r8
183 END DO
184 END DO
185 END DO
186!
187! Compute vert.
188!
189 DO k=1,n(ng)
190 DO j=jstr,jend
191 DO i=istr,iend+1
192 wrk(i,j)=u(i,j,k,ninp)*(z_r(i,j,k)-z_r(i-1,j,k))* &
193 & (pm(i-1,j)+pm(i,j))
194 END DO
195 DO i=istr,iend
196 vert(i,j,k)=0.25_r8*(wrk(i,j)+wrk(i+1,j))
197 END DO
198 END DO
199 DO j=jstr,jend+1
200 DO i=istr,iend
201 wrk(i,j)=v(i,j,k,ninp)*(z_r(i,j,k)-z_r(i,j-1,k))* &
202 & (pn(i,j-1)+pn(i,j))
203 END DO
204 END DO
205 DO j=jstr,jend
206 DO i=istr,iend
207 vert(i,j,k)=vert(i,j,k)+0.25_r8*(wrk(i,j)+wrk(i,j+1))
208 END DO
209 END DO
210 END DO
211!
212 cff1=3.0_r8/8.0_r8
213 cff2=3.0_r8/4.0_r8
214 cff3=1.0_r8/8.0_r8
215 cff4=9.0_r8/16.0_r8
216 cff5=1.0_r8/16.0_r8
217
218 j_loop : DO j=jstr,jend
219 DO i=istr,iend
220 wrk(i,j)=(du_avg1(i,j)-du_avg1(i+1,j)+ &
221 & dv_avg1(i,j)-dv_avg1(i,j+1))/ &
222 & (z_w(i,j,n(ng))-z_w(i,j,0))
223 END DO
224 DO i=istr,iend
225 slope=(z_w(i,j,n(ng))-z_r(i,j,n(ng) ))/ &
226 & (z_r(i,j,n(ng))-z_r(i,j,n(ng)-1)) ! extrapolation slope
227!^ tl_wvel(i,j,N(ng)-1)=pm(i,j)*pn(i,j)* &
228!^ & (tl_W(i,j,N(ng)-1)+ &
229!^ & tl_wrk(i,j)* &
230!^ & (z_w(i,j,N(ng)-1)-z_w(i,j,0))+ &
231!^ & wrk(i,j)* &
232!^ & (tl_z_w(i,j,N(ng)-1)-tl_z_w(i,j,0)))+ &
233!^ & cff1*tl_vert(i,j,N(ng) )+ &
234!^ & cff2*tl_vert(i,j,N(ng)-1)- &
235!^ & cff3*tl_vert(i,j,N(ng)-2)
236!^
237 adfac=pm(i,j)*pn(i,j)*ad_wvel(i,j,n(ng)-1)
238 adfac1=wrk(i,j)*adfac
239 ad_w(i,j,n(ng)-1)=ad_w(i,j,n(ng)-1)+adfac
240 ad_wrk(i,j)=ad_wrk(i,j)+(z_w(i,j,n(ng)-1)-z_w(i,j,0))*adfac
241 ad_z_w(i,j,n(ng)-1)=ad_z_w(i,j,n(ng)-1)+adfac1
242 ad_z_w(i,j,0)=ad_z_w(i,j,0)-adfac1
243 ad_vert(i,j,n(ng) )=ad_vert(i,j,n(ng))+ &
244 & cff1*ad_wvel(i,j,n(ng)-1)
245 ad_vert(i,j,n(ng)-1)=ad_vert(i,j,n(ng)-1)+ &
246 & cff2*ad_wvel(i,j,n(ng)-1)
247 ad_vert(i,j,n(ng)-2)=ad_vert(i,j,n(ng)-2)- &
248 & cff3*ad_wvel(i,j,n(ng)-1)
249 ad_wvel(i,j,n(ng)-1)=0.0_r8
250!^ tl_wvel(i,j,N(ng))=pm(i,j)*pn(i,j)* &
251!^ & (tl_wrk(i,j)* &
252!^ & (z_w(i,j,N(ng))-z_w(i,j,0))+ &
253!^ & wrk(i,j)* &
254!^ & (tl_z_w(i,j,N(ng))-tl_z_w(i,j,0)))+ &
255!^ & cff1*(tl_vert(i,j,N(ng))+ &
256!^ & tl_slope*(vert(i,j,N(ng) )- &
257!^ & vert(i,j,N(ng)-1))+ &
258!^ & slope*(tl_vert(i,j,N(ng) )- &
259!^ & tl_vert(i,j,N(ng)-1)))+ &
260!^ & cff2*tl_vert(i,j,N(ng) )- &
261!^ & cff3*tl_vert(i,j,N(ng)-1)
262!^
263 adfac=pm(i,j)*pn(i,j)*ad_wvel(i,j,n(ng))
264 adfac1=wrk(i,j)*adfac
265 adfac2=cff1*ad_wvel(i,j,n(ng))
266 adfac3=slope*adfac2
267 ad_wrk(i,j)=ad_wrk(i,j)+(z_w(i,j,n(ng))-z_w(i,j,0))*adfac
268 ad_z_w(i,j,n(ng))=ad_z_w(i,j,n(ng))+adfac1
269 ad_z_w(i,j,0 )=ad_z_w(i,j,0 )-adfac1
270 ad_vert(i,j,n(ng))=ad_vert(i,j,n(ng))+adfac2
271 ad_slope=ad_slope+ &
272 & (vert(i,j,n(ng) )-vert(i,j,n(ng)-1))*adfac2
273 ad_vert(i,j,n(ng) )=ad_vert(i,j,n(ng) )+ &
274 & adfac3+cff2*ad_wvel(i,j,n(ng))
275 ad_vert(i,j,n(ng)-1)=ad_vert(i,j,n(ng)-1)- &
276 & adfac3-cff3*ad_wvel(i,j,n(ng))
277 ad_wvel(i,j,n(ng))=0.0_r8
278!^ tl_slope=(tl_z_w(i,j,N(ng))-tl_z_r(i,j,N(ng) ))/ &
279!^ & (z_r(i,j,N(ng))-z_r(i,j,N(ng)-1))- &
280!^ & (tl_z_r(i,j,N(ng))-tl_z_r(i,j,N(ng)-1))*slope/ &
281!^ & (z_r(i,j,N(ng))-z_r(i,j,N(ng)-1))
282!^
283 adfac1=ad_slope/(z_r(i,j,n(ng))-z_r(i,j,n(ng)-1))
284 adfac2=slope*adfac1
285 ad_z_w(i,j,n(ng))=ad_z_w(i,j,n(ng))+adfac1
286 ad_z_r(i,j,n(ng) )=ad_z_r(i,j,n(ng) )-adfac1-adfac2
287 ad_z_r(i,j,n(ng)-1)=ad_z_r(i,j,n(ng)-1)+adfac2
288 ad_slope=0.0_r8
289 END DO
290 DO k=2,n(ng)-2
291 DO i=istr,iend
292!^ tl_wvel(i,j,k)=pm(i,j)*pn(i,j)* &
293!^ & (tl_W(i,j,k)+ &
294!^ & tl_wrk(i,j)*(z_w(i,j,k)-z_w(i,j,0))+ &
295!^ & wrk(i,j)*(tl_z_w(i,j,k)-tl_z_w(i,j,0)))+ &
296!^ & cff4*(tl_vert(i,j,k )+tl_vert(i,j,k+1))- &
297!^ & cff5*(tl_vert(i,j,k-1)+tl_vert(i,j,k+2))
298!^
299 adfac=pm(i,j)*pn(i,j)*ad_wvel(i,j,k)
300 adfac1=wrk(i,j)*adfac
301 adfac2=cff4*ad_wvel(i,j,k)
302 adfac3=cff5*ad_wvel(i,j,k)
303 ad_w(i,j,k)=ad_w(i,j,k)+adfac
304 ad_wrk(i,j)=ad_wrk(i,j)+(z_w(i,j,k)-z_w(i,j,0))*adfac
305 ad_z_w(i,j,k)=ad_z_w(i,j,k)+adfac1
306 ad_z_w(i,j,0)=ad_z_w(i,j,0)-adfac1
307 ad_vert(i,j,k )=ad_vert(i,j,k )+adfac2
308 ad_vert(i,j,k+1)=ad_vert(i,j,k+1)+adfac2
309 ad_vert(i,j,k-1)=ad_vert(i,j,k-1)-adfac3
310 ad_vert(i,j,k+2)=ad_vert(i,j,k+2)-adfac3
311 ad_wvel(i,j,k)=0.0_r8
312 END DO
313 END DO
314 DO i=istr,iend
315 slope=(z_r(i,j,1)-z_w(i,j,0))/ &
316 & (z_r(i,j,2)-z_r(i,j,1)) ! extrapolation slope
317!^ tl_wvel(i,j,1)=pm(i,j)*pn(i,j)* &
318!^ & (tl_W(i,j,1)+ &
319!^ & tl_wrk(i,j)*(z_w(i,j,1)-z_w(i,j,0))+ &
320!^ & wrk(i,j)*(tl_z_w(i,j,1)-tl_z_w(i,j,0)))+ &
321!^ & cff1*tl_vert(i,j,1)+ &
322!^ & cff2*tl_vert(i,j,2)- &
323!^ & cff3*tl_vert(i,j,3)
324!^
325 adfac=pm(i,j)*pn(i,j)*ad_wvel(i,j,1)
326 adfac1=wrk(i,j)*adfac
327 ad_w(i,j,1)=ad_w(i,j,1)+adfac
328 ad_wrk(i,j)=ad_wrk(i,j)+(z_w(i,j,1)-z_w(i,j,0))*adfac
329 ad_z_w(i,j,1)=ad_z_w(i,j,1)+adfac1
330 ad_z_w(i,j,0)=ad_z_w(i,j,0)-adfac1
331 ad_vert(i,j,1)=ad_vert(i,j,1)+cff1*ad_wvel(i,j,1)
332 ad_vert(i,j,2)=ad_vert(i,j,2)+cff2*ad_wvel(i,j,1)
333 ad_vert(i,j,3)=ad_vert(i,j,3)-cff3*ad_wvel(i,j,1)
334 ad_wvel(i,j,1)=0.0_r8
335!^ tl_wvel(i,j,0)=cff1*(tl_vert(i,j,1)- &
336!^ & tl_slope*(vert(i,j,2)- &
337!^ & vert(i,j,1))- &
338!^ & slope*(tl_vert(i,j,2)- &
339!^ & tl_vert(i,j,1)))+ &
340!^ & cff2*tl_vert(i,j,1)- &
341!^ & cff3*tl_vert(i,j,2)
342!^
343 adfac=cff1*ad_wvel(i,j,0)
344 adfac1=slope*adfac
345 ad_vert(i,j,1)=ad_vert(i,j,1)+adfac
346 ad_slope=ad_slope-(vert(i,j,2)-vert(i,j,1))*adfac
347 ad_vert(i,j,2)=ad_vert(i,j,2)-adfac1
348 ad_vert(i,j,1)=ad_vert(i,j,1)+adfac1
349 ad_vert(i,j,1)=ad_vert(i,j,1)+cff2*ad_wvel(i,j,0)
350 ad_vert(i,j,2)=ad_vert(i,j,2)-cff3*ad_wvel(i,j,0)
351 ad_wvel(i,j,0)=0.0_r8
352!^ tl_slope=(tl_z_r(i,j,1)-tl_z_w(i,j,0))/ &
353!^ & (z_r(i,j,2)-z_r(i,j,1))- &
354!^ & (tl_z_r(i,j,2)-tl_z_r(i,j,1))*slope/ &
355!^ & (z_r(i,j,2)-z_r(i,j,1))
356!^
357 adfac=ad_slope/(z_r(i,j,2)-z_r(i,j,1))
358 adfac1=slope*adfac
359 ad_z_r(i,j,1)=ad_z_r(i,j,1)+adfac
360 ad_z_w(i,j,0)=ad_z_w(i,j,0)-adfac
361 ad_z_r(i,j,2)=ad_z_r(i,j,2)-adfac1
362 ad_z_r(i,j,1)=ad_z_r(i,j,1)+adfac1
363 ad_slope=0.0_r8
364 END DO
365 DO i=istr,iend
366!^ tl_wrk(i,j)=(tl_DU_avg1(i,j)-tl_DU_avg1(i+1,j)+ &
367!^ & tl_DV_avg1(i,j)-tl_DV_avg1(i,j+1))/ &
368!^ & (z_w(i,j,N(ng))-z_w(i,j,0))- &
369!^ & (tl_z_w(i,j,N(ng))-tl_z_w(i,j,0))*wrk(i,j)/ &
370!^ & (z_w(i,j,N(ng))-z_w(i,j,0))
371!^
372 adfac=ad_wrk(i,j)/(z_w(i,j,n(ng))-z_w(i,j,0))
373 adfac1=wrk(i,j)*adfac
374 ad_du_avg1(i ,j)=ad_du_avg1(i ,j)+adfac
375 ad_du_avg1(i+1,j)=ad_du_avg1(i+1,j)-adfac
376 ad_dv_avg1(i,j )=ad_dv_avg1(i,j )+adfac
377 ad_dv_avg1(i,j+1)=ad_dv_avg1(i,j+1)-adfac
378 ad_z_w(i,j,n(ng))=ad_z_w(i,j,n(ng))-adfac1
379 ad_z_w(i,j,0)=ad_z_w(i,j,0)+adfac1
380 ad_wrk(i,j)=0.0_r8
381 END DO
382
383 END DO j_loop
384 DO k=1,n(ng)
385 DO j=jstr,jend
386 DO i=istr,iend
387!^ tl_vert(i,j,k)=tl_vert(i,j,k)+ &
388!^ & 0.25_r8*(tl_wrk(i,j)+tl_wrk(i,j+1))
389!^
390 adfac=0.25_r8*ad_vert(i,j,k)
391 ad_wrk(i,j )=ad_wrk(i,j )+adfac
392 ad_wrk(i,j+1)=ad_wrk(i,j+1)+adfac
393 END DO
394 END DO
395 DO j=jstr,jend+1
396 DO i=istr,iend
397!^ tl_wrk(i,j)=(pn(i,j-1)+pn(i,j))* &
398!^ & (tl_v(i,j,k,Ninp)*(z_r(i,j,k)-z_r(i,j-1,k))+ &
399!^ & v(i,j,k,Ninp)*(tl_z_r(i,j,k)-tl_z_r(i,j-1,k)))
400!^
401 adfac=(pn(i,j-1)+pn(i,j))*ad_wrk(i,j)
402 adfac1=v(i,j,k,ninp)*adfac
403 ad_v(i,j,k,ninp)=ad_v(i,j,k,ninp)+ &
404 & (z_r(i,j,k)-z_r(i,j-1,k))*adfac
405 ad_z_r(i,j ,k)=ad_z_r(i,j ,k)+adfac1
406 ad_z_r(i,j-1,k)=ad_z_r(i,j-1,k)-adfac1
407 ad_wrk(i,j)=0.0_r8
408 END DO
409 END DO
410 DO j=jstr,jend
411 DO i=istr,iend
412!^ tl_vert(i,j,k)=0.25_r8*(tl_wrk(i,j)+tl_wrk(i+1,j))
413!^
414 adfac=0.25_r8*ad_vert(i,j,k)
415 ad_wrk(i ,j)=ad_wrk(i ,j)+adfac
416 ad_wrk(i+1,j)=ad_wrk(i+1,j)+adfac
417 ad_vert(i,j,k)=0.0_r8
418 END DO
419 DO i=istr,iend+1
420!^ tl_wrk(i,j)=(pm(i-1,j)+pm(i,j))* &
421!^ & (tl_u(i,j,k,Ninp)*(z_r(i,j,k)-z_r(i-1,j,k))+ &
422!^ & u(i,j,k,Ninp)*(tl_z_r(i,j,k)-tl_z_r(i-1,j,k)))
423!^
424 adfac=(pm(i-1,j)+pm(i,j))*ad_wrk(i,j)
425 adfac1=u(i,j,k,ninp)*adfac
426 ad_u(i,j,k,ninp)=ad_u(i,j,k,ninp)+ &
427 & (z_r(i,j,k)-z_r(i-1,j,k))*adfac
428 ad_z_r(i ,j,k)=ad_z_r(i ,j,k)+adfac1
429 ad_z_r(i-1,j,k)=ad_z_r(i-1,j,k)-adfac1
430 ad_wrk(i,j)=0.0_r8
431 END DO
432 END DO
433 END DO
434!
435! Exchange time-averaged fields.
436!
437# ifdef DISTRIBUTE
438 CALL ad_mp_exchange2d (ng, tile, iadm, 2, &
439 & lbi, ubi, lbj, ubj, &
440 & nghostpoints, &
441 & ewperiodic(ng), nsperiodic(ng), &
442 & ad_du_avg1, ad_dv_avg1)
443# endif
444
445 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
446 CALL ad_exchange_u2d_tile (ng, tile, &
447 & lbi, ubi, lbj, ubj, &
448 & ad_du_avg1)
449 CALL ad_exchange_v2d_tile (ng, tile, &
450 & lbi, ubi, lbj, ubj, &
451 & ad_dv_avg1)
452 END IF
453
454 RETURN
455 END SUBROUTINE ad_wvelocity_tile
456#endif
457 END MODULE ad_wvelocity_mod
subroutine ad_bc_w3d_tile(ng, tile, lbi, ubi, lbj, ubj, lbk, ubk, ad_a)
Definition ad_bc_3d.F:733
subroutine ad_exchange_u2d_tile(ng, tile, lbi, ubi, lbj, ubj, ad_a)
subroutine ad_exchange_v2d_tile(ng, tile, lbi, ubi, lbj, ubj, ad_a)
subroutine, public ad_wvelocity(ng, tile, ninp)
subroutine ad_wvelocity_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, ninp, pm, pn, z_r, z_w, ad_z_r, ad_z_w, du_avg1, dv_avg1, ad_du_avg1, ad_dv_avg1, u, v, w, ad_u, ad_v, ad_w, ad_wvel)
subroutine bc_w3d_tile(ng, tile, lbi, ubi, lbj, ubj, lbk, ubk, a)
Definition bc_3d.F:591
type(t_coupling), dimension(:), allocatable coupling
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
integer, parameter iadm
Definition mod_param.F:665
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
subroutine ad_mp_exchange3d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, lbk, ubk, nghost, ew_periodic, ns_periodic, ad_a, ad_b, ad_c, ad_d)
subroutine ad_mp_exchange2d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, nghost, ew_periodic, ns_periodic, ad_a, ad_b, ad_c, ad_d)