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

Functions/Subroutines

subroutine, public ad_u2dbc (ng, tile, kout)
 
subroutine, public ad_u2dbc_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, krhs, kstp, kout, ubar, vbar, zeta, ad_ubar, ad_vbar, ad_zeta)
 

Function/Subroutine Documentation

◆ ad_u2dbc()

subroutine, public ad_u2dbc_mod::ad_u2dbc ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) kout )

Definition at line 28 of file ad_u2dbc_im.F.

29!***********************************************************************
30!
31 USE mod_param
32 USE mod_ocean
33 USE mod_stepping
34!
35! Imported variable declarations.
36!
37 integer, intent(in) :: ng, tile, kout
38!
39! Local variable declarations.
40!
41# include "tile.h"
42!
43 CALL ad_u2dbc_tile (ng, tile, &
44 & lbi, ubi, lbj, ubj, &
45 & imins, imaxs, jmins, jmaxs, &
46 & krhs(ng), kstp(ng), kout, &
47 & ocean(ng) % ubar, &
48 & ocean(ng) % vbar, &
49 & ocean(ng) % zeta, &
50 & ocean(ng) % ad_ubar, &
51 & ocean(ng) % ad_ubar, &
52 & ocean(ng) % ad_zeta)
53
54 RETURN
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, dimension(:), allocatable kstp
integer, dimension(:), allocatable krhs

References ad_u2dbc_tile(), mod_stepping::krhs, mod_stepping::kstp, and mod_ocean::ocean.

Here is the call graph for this function:

◆ ad_u2dbc_tile()

subroutine, public ad_u2dbc_mod::ad_u2dbc_tile ( integer, intent(in) ng,
integer, intent(in) tile,
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) krhs,
integer, intent(in) kstp,
integer, intent(in) kout,
real(r8), dimension(lbi:,lbj:,:), intent(in) ubar,
real(r8), dimension(lbi:,lbj:,:), intent(in) vbar,
real(r8), dimension(lbi:,lbj:,:), intent(in) zeta,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_ubar,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_vbar,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_zeta )

Definition at line 59 of file ad_u2dbc_im.F.

65!***********************************************************************
66!
67 USE mod_param
68 USE mod_boundary
69 USE mod_clima
70 USE mod_forces
71 USE mod_grid
72 USE mod_ncparam
73 USE mod_scalars
74!
75! Imported variable declarations.
76!
77 integer, intent(in) :: ng, tile
78 integer, intent(in) :: LBi, UBi, LBj, UBj
79 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
80 integer, intent(in) :: krhs, kstp, kout
81!
82# ifdef ASSUMED_SHAPE
83 real(r8), intent(in) :: ubar(LBi:,LBj:,:)
84 real(r8), intent(in) :: vbar(LBi:,LBj:,:)
85 real(r8), intent(in) :: zeta(LBi:,LBj:,:)
86
87 real(r8), intent(inout) :: ad_ubar(LBi:,LBj:,:)
88 real(r8), intent(inout) :: ad_vbar(LBi:,LBj:,:)
89 real(r8), intent(inout) :: ad_zeta(LBi:,LBj:,:)
90# else
91 real(r8), intent(in) :: ubar(LBi:UBi,LBj:UBj,:)
92 real(r8), intent(in) :: vbar(LBi:UBi,LBj:UBj,:)
93 real(r8), intent(in) :: zeta(LBi:UBi,LBj:UBj,:)
94
95 real(r8), intent(inout) :: ad_ubar(LBi:UBi,LBj:UBj,:)
96 real(r8), intent(inout) :: ad_vbar(LBi:UBi,LBj:UBj,:)
97 real(r8), intent(inout) :: ad_zeta(LBi:UBi,LBj:UBj,:)
98# endif
99!
100! Local variable declarations.
101!
102 integer :: Imin, Imax
103 integer :: i, j, know
104
105 real(r8) :: Ce, Cx, Zx
106 real(r8) :: bry_pgr, bry_cor, bry_str, bry_val
107 real(r8) :: cff, cff1, cff2, cff3, dt2d
108 real(r8) :: obc_in, obc_out, tau
109# if defined ATM_PRESS && defined PRESS_COMPENSATE
110 real(r8) :: OneAtm, fac
111# endif
112
113 real(r8) :: ad_Ce, ad_Cx
114 real(r8) :: ad_bry_pgr, ad_bry_cor, ad_bry_str, ad_bry_val, ad_Zx
115 real(r8) :: ad_cff, ad_cff1, ad_cff2, ad_cff3
116 real(r8) :: adfac
117
118 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_grad
119
120# include "set_bounds.h"
121!
122!-----------------------------------------------------------------------
123! Initialize adjoint private variables.
124!-----------------------------------------------------------------------
125!
126 ad_ce=0.0_r8
127 ad_cx=0.0_r8
128 ad_zx=0.0_r8
129 ad_cff=0.0_r8
130 ad_cff1=0.0_r8
131 ad_cff2=0.0_r8
132 ad_cff3=0.0_r8
133 ad_bry_pgr=0.0_r8
134 ad_bry_cor=0.0_r8
135 ad_bry_str=0.0_r8
136 ad_bry_val=0.0_r8
137
138 ad_grad(lbi:ubi,lbj:ubj)=0.0_r8
139!
140!-----------------------------------------------------------------------
141! Set time-indices
142!-----------------------------------------------------------------------
143!
144 IF (first_2d_step) THEN
145 know=krhs
146 dt2d=dtfast(ng)
147 ELSE IF (predictor_2d_step(ng)) THEN
148 know=krhs
149 dt2d=2.0_r8*dtfast(ng)
150 ELSE
151 know=kstp
152 dt2d=dtfast(ng)
153 END IF
154# if defined ATM_PRESS && defined PRESS_COMPENSATE
155 oneatm=1013.25_r8 ! 1 atm = 1013.25 mb
156 fac=100.0_r8/(g*rho0)
157# endif
158
159# if defined WET_DRY_NOT_YET
160!
161!-----------------------------------------------------------------------
162! Impose wetting and drying conditions.
163!
164! HGA: need ADM code here for the NLM code below.
165!-----------------------------------------------------------------------
166!
167 IF (.not.(ewperiodic(ng).or.nsperiodic(ng))) THEN
168 IF (domain(ng)%NorthEast_Corner(tile)) THEN
169 IF ((lbc_apply(ng)%north(iend+1).and. &
170 & lbc_apply(ng)%east (jend+1)).or. &
171 & (lbc(ieast,isubar,ng)%nested.and. &
172 & lbc(inorth,isubar,ng)%nested)) THEN
173!^ cff1=ABS(ABS(GRID(ng)%umask_wet(Iend+1,Jend+1))-1.0_r8)
174!^ cff2=0.5_r8+DSIGN(0.5_r8,ubar(Iend+1,Jend+1,kout))* &
175!^ & GRID(ng)%umask_wet(Iend+1,Jend+1)
176!^ cff=0.5_r8*GRID(ng)%umask_wet(Iend+1,Jend+1)+cff1+ &
177!^ & cff2*(1.0_r8-cff1)
178!^ ubar(Iend+1,Jend+1,kout)=ubar(Iend+1,Jend+1,kout)*cff
179 END IF
180 END IF
181 IF (domain(ng)%NorthWest_Corner(tile)) THEN
182 IF ((lbc_apply(ng)%north(istr ).and. &
183 & lbc_apply(ng)%west (jend+1)).or. &
184 & (lbc(iwest,isubar,ng)%nested.and. &
185 & lbc(inorth,isubar,ng)%nested)) THEN
186!^ cff1=ABS(ABS(GRID(ng)%umask_wet(Istr,Jend+1))-1.0_r8)
187!^ cff2=0.5_r8+DSIGN(0.5_r8,ubar(Istr,Jend+1,kout))* &
188!^ & GRID(ng)%umask_wet(Istr,Jend+1)
189!^ cff=0.5_r8*GRID(ng)%umask_wet(Istr,Jend+1)*cff1+ &
190!^ & cff2*(1.0_r8-cff1)
191!^ ubar(Istr,Jend+1,kout)=ubar(Istr,Jend+1,kout)*cff
192 END IF
193 END IF
194 IF (domain(ng)%SouthEast_Corner(tile)) THEN
195 IF ((lbc_apply(ng)%south(iend+1).and. &
196 & lbc_apply(ng)%east (jstr-1)).or. &
197 & (lbc(ieast,isubar,ng)%nested.and. &
198 & lbc(isouth,isubar,ng)%nested)) THEN
199!^ cff1=ABS(ABS(GRID(ng)%umask_wet(Iend+1,Jstr-1))-1.0_r8)
200!^ cff2=0.5_r8+DSIGN(0.5_r8,ubar(Iend+1,Jstr-1,kout))* &
201!^ & GRID(ng)%umask_wet(Iend+1,Jstr-1)
202!^ cff=0.5_r8*GRID(ng)%umask_wet(Iend+1,Jstr-1)*cff1+ &
203!^ & cff2*(1.0_r8-cff1)
204!^ ubar(Iend+1,Jstr-1,kout)=ubar(Iend+1,Jstr-1,kout)*cff
205 END IF
206 END IF
207 IF (domain(ng)%SouthWest_Corner(tile)) THEN
208 IF ((lbc_apply(ng)%south(istr ).and. &
209 & lbc_apply(ng)%west (jstr-1)).or. &
210 & (lbc(iwest,isubar,ng)%nested.and. &
211 & lbc(isouth,isubar,ng)%nested)) THEN
212!^ cff1=ABS(ABS(GRID(ng)%umask_wet(Istr,Jstr-1))-1.0_r8)
213!^ cff2=0.5_r8+DSIGN(0.5_r8,ubar(Istr,Jstr-1,kout))* &
214!^ & GRID(ng)%umask_wet(Istr,Jstr-1)
215!^ cff=0.5_r8*GRID(ng)%umask_wet(Istr,Jstr-1)*cff1+ &
216!^ & cff2*(1.0_r8-cff1)
217!^ ubar(Istr,Jstr-1,kout)=ubar(Istr,Jstr-1,kout)*cff
218 END IF
219 END IF
220 END IF
221!
222 IF (.not.nsperiodic(ng)) THEN
223 IF (domain(ng)%Northern_Edge(tile)) THEN
224 DO i=istr,iend
225 IF (lbc_apply(ng)%north(i).or. &
226 & lbc(inorth,isubar,ng)%nested) THEN
227!^ cff1=ABS(ABS(GRID(ng)%umask_wet(i,Jend+1))-1.0_r8)
228!^ cff2=0.5_r8+DSIGN(0.5_r8,ubar(i,Jend+1,kout))* &
229!^ & GRID(ng)%umask_wet(i,Jend+1)
230!^ cff=0.5_r8*GRID(ng)%umask_wet(i,Jend+1)*cff1+ &
231!^ & cff2*(1.0_r8-cff1)
232!^ ubar(i,Jend+1,kout)=ubar(i,Jend+1,kout)*cff
233 END IF
234 END DO
235 END IF
236 IF (domain(ng)%Southern_Edge(tile)) THEN
237 DO i=istru,iend
238 IF (lbc_apply(ng)%south(i).or. &
239 & lbc(isouth,isubar,ng)%nested) THEN
240!^ cff1=ABS(ABS(GRID(ng)%umask_wet(i,Jstr-1))-1.0_r8)
241!^ cff2=0.5_r8+DSIGN(0.5_r8,ubar(i,Jstr-1,kout))* &
242!^ & GRID(ng)%umask_wet(i,Jstr-1)
243!^ cff=0.5_r8*GRID(ng)%umask_wet(i,Jstr-1)*cff1+ &
244!^ & cff2*(1.0_r8-cff1)
245!^ ubar(i,Jstr-1,kout)=ubar(i,Jstr-1,kout)*cff
246 END IF
247 END DO
248 END IF
249 END IF
250!
251 IF (.not.ewperiodic(ng)) THEN
252 IF (domain(ng)%Eastern_Edge(tile)) THEN
253 DO j=jstr,jend
254 IF (lbc_apply(ng)%east(j).or. &
255 & lbc(ieast,isubar,ng)%nested) THEN
256!^ cff1=ABS(ABS(GRID(ng)%umask_wet(Iend+1,j))-1.0_r8)
257!^ cff2=0.5_r8+DSIGN(0.5_r8,ubar(Iend+1,j,kout))* &
258!^ & GRID(ng)%umask_wet(Iend+1,j)
259!^ cff=0.5_r8*GRID(ng)%umask_wet(Iend+1,j)*cff1+ &
260!^ & cff2*(1.0_r8-cff1)
261!^ ubar(Iend+1,j,kout)=ubar(Iend+1,j,kout)*cff
262 END IF
263 END DO
264 END IF
265 IF (domain(ng)%Western_Edge(tile)) THEN
266 DO j=jstr,jend
267 IF (lbc_apply(ng)%west(j).or. &
268 & lbc(iwest,isubar,ng)%nested) THEN
269!^ cff1=ABS(ABS(GRID(ng)%umask_wet(Istr,j))-1.0_r8)
270!^ cff2=0.5_r8+DSIGN(0.5_r8,ubar(Istr,j,kout))* &
271!^ & GRID(ng)%umask_wet(Istr,j)
272!^ cff=0.5_r8*GRID(ng)%umask_wet(Istr,j)*cff1+ &
273!^ & cff2*(1.0_r8-cff1)
274!^ ubar(Istr,j,kout)=ubar(Istr,j,kout)*cff
275 END IF
276 END DO
277 END IF
278 END IF
279# endif
280!
281!-----------------------------------------------------------------------
282! Boundary corners.
283!-----------------------------------------------------------------------
284!
285 IF (.not.(ewperiodic(ng).or.nsperiodic(ng))) THEN
286 IF (domain(ng)%NorthEast_Corner(tile)) THEN
287 IF (lbc_apply(ng)%north(iend+1).and. &
288 & lbc_apply(ng)%east (jend+1)) THEN
289!^ tl_ubar(Iend+1,Jend+1,kout)=0.5_r8* &
290!^ & (tl_ubar(Iend+1,Jend ,kout)+ &
291!^ & tl_ubar(Iend ,Jend+1,kout))
292!^
293 adfac=0.5_r8*ad_ubar(iend+1,jend+1,kout)
294 ad_ubar(iend+1,jend ,kout)=ad_ubar(iend+1,jend ,kout)+ &
295 & adfac
296 ad_ubar(iend ,jend+1,kout)=ad_ubar(iend ,jend+1,kout)+ &
297 & adfac
298 ad_ubar(iend+1,jend+1,kout)=0.0_r8
299 END IF
300 END IF
301 IF (domain(ng)%NorthWest_Corner(tile)) THEN
302 IF (lbc_apply(ng)%north(istr ).and. &
303 & lbc_apply(ng)%west (jend+1)) THEN
304!^ tl_ubar(Istr,Jend+1,kout)=0.5_r8* &
305!^ & (tl_ubar(Istr ,Jend ,kout)+ &
306!^ & tl_ubar(Istr+1,Jend+1,kout))
307!^
308 adfac=0.5_r8*ad_ubar(istr,jend+1,kout)
309 ad_ubar(istr ,jend ,kout)=ad_ubar(istr ,jend ,kout)+ &
310 & adfac
311 ad_ubar(istr+1,jend+1,kout)=ad_ubar(istr+1,jend+1,kout)+ &
312 & adfac
313 ad_ubar(istr ,jend+1,kout)=0.0_r8
314 END IF
315 END IF
316 IF (domain(ng)%SouthEast_Corner(tile)) THEN
317 IF (lbc_apply(ng)%south(iend+1).and. &
318 & lbc_apply(ng)%east (jstr-1)) THEN
319!^ tl_ubar(Iend+1,Jstr-1,kout)=0.5_r8* &
320!^ & (tl_ubar(Iend ,Jstr-1,kout)+ &
321!^ & tl_ubar(Iend+1,Jstr ,kout))
322!^
323 adfac=0.5_r8*ad_ubar(iend+1,jstr-1,kout)
324 ad_ubar(iend ,jstr-1,kout)=ad_ubar(iend ,jstr-1,kout)+ &
325 & adfac
326 ad_ubar(iend+1,jstr ,kout)=ad_ubar(iend+1,jstr ,kout)+ &
327 & adfac
328 ad_ubar(iend+1,jstr-1,kout)=0.0_r8
329 END IF
330 END IF
331 IF (domain(ng)%SouthWest_Corner(tile)) THEN
332 IF (lbc_apply(ng)%south(istr ).and. &
333 & lbc_apply(ng)%west (jstr-1)) THEN
334!^ ubar(Istr,Jstr-1,kout)=0.5_r8* &
335!^ & (ubar(Istr+1,Jstr-1,kout)+ &
336!^ & ubar(Istr ,Jstr ,kout))
337!^
338 adfac=0.5_r8*ad_ubar(istr,jstr-1,kout)
339 ad_ubar(istr+1,jstr-1,kout)=ad_ubar(istr+1,jstr-1,kout)+ &
340 & adfac
341 ad_ubar(istr ,jstr ,kout)=ad_ubar(istr ,jstr ,kout)+ &
342 & adfac
343 ad_ubar(istr ,jstr-1,kout)=0.0_r8
344 END IF
345 END IF
346 END IF
347!
348!-----------------------------------------------------------------------
349! Lateral boundary conditions at the northern edge.
350!-----------------------------------------------------------------------
351!
352 IF (domain(ng)%Northern_Edge(tile)) THEN
353!
354! Northern edge, implicit upstream radiation condition.
355!
356 IF (ad_lbc(inorth,isubar,ng)%radiation) THEN
357 IF (iic(ng).ne.0) THEN
358 DO i=istru,iend
359 IF (lbc_apply(ng)%north(i)) THEN
360# if defined CELERITY_READ && defined FORWARD_READ
361 IF (ad_lbc(inorth,isubar,ng)%nudging) THEN
362 IF (lnudgem2clm(ng)) THEN
363 obc_out=0.5_r8* &
364 & (clima(ng)%M2nudgcof(i-1,jend+1)+ &
365 & clima(ng)%M2nudgcof(i ,jend+1))
366 obc_in =obcfac(ng)*obc_out
367 ELSE
368 obc_out=m2obc_out(ng,inorth)
369 obc_in =m2obc_in(ng,inorth)
370 END IF
371 IF (boundary(ng)%ubar_north_Ce(i).lt.0.0_r8) THEN
372 tau=obc_in
373 ELSE
374 tau=obc_out
375 END IF
376 tau=tau*dt2d
377 END IF
378# ifdef RADIATION_2D
379 cx=boundary(ng)%ubar_north_Cx(i)
380# else
381 cx=0.0_r8
382# endif
383 ce=boundary(ng)%ubar_north_Ce(i)
384 cff=boundary(ng)%ubar_north_C2(i)
385# endif
386# ifdef MASKING
387!^ tl_ubar(i,Jend+1,kout)=tl_ubar(i,Jend+1,kout)* &
388!^ & GRID(ng)%umask(i,Jend+1)
389!^
390 ad_ubar(i,jend+1,kout)=ad_ubar(i,jend+1,kout)* &
391 & grid(ng)%umask(i,jend+1)
392# endif
393 IF (ad_lbc(inorth,isubar,ng)%nudging) THEN
394!^ tl_ubar(i,Jend+1,kout)=tl_ubar(i,Jend+1,kout)- &
395!^ & tau*tl_ubar(i,Jend+1,know)
396!^
397 ad_ubar(i,jend+1 ,know)=ad_ubar(i,jend+1 ,know)- &
398 & tau*ad_ubar(i,jend+1,kout)
399 END IF
400!^ tl_ubar(i,Jend+1,kout)=(cff*tl_ubar(i,Jend+1,know)+ &
401!^ & Ce *tl_ubar(i,Jend ,kout)- &
402!^ & MAX(Cx,0.0_r8)* &
403!^ & tl_grad(i-1,Jend+1)- &
404!^ & MIN(Cx,0.0_r8)* &
405!^ & tl_grad(i ,Jend+1))/ &
406!^ & (cff+Ce)
407!^
408 adfac=ad_ubar(i,jend+1,kout)/(cff+ce)
409 ad_grad(i-1,jend+1)=ad_grad(i-1,jend+1)- &
410 & max(cx,0.0_r8)*adfac
411 ad_grad(i ,jend+1)=ad_grad(i ,jend+1)- &
412 & min(cx,0.0_r8)*adfac
413 ad_ubar(i,jend ,kout)=ad_ubar(i,jend ,kout)+ce *adfac
414 ad_ubar(i,jend+1,know)=ad_ubar(i,jend+1,know)+cff*adfac
415 ad_ubar(i,jend+1,kout)=0.0_r8
416 END IF
417 END DO
418 END IF
419!
420! Northern edge, Chapman boundary condition.
421!
422 ELSE IF (ad_lbc(inorth,isubar,ng)%Flather.or. &
423 & ad_lbc(inorth,isubar,ng)%Shchepetkin.or. &
424 & ad_lbc(inorth,isubar,ng)%reduced) THEN
425 DO i=istru,iend
426 IF (lbc_apply(ng)%north(i)) THEN
427 cff=dt2d*0.5_r8*(grid(ng)%pn(i-1,jend)+ &
428 & grid(ng)%pn(i ,jend))
429 cff1=sqrt(g*0.5_r8*(grid(ng)%h(i-1,jend)+ &
430 & zeta(i-1,jend,know)+ &
431 & grid(ng)%h(i ,jend)+ &
432 & zeta(i ,jend,know)))
433 ce=cff*cff1
434 cff2=1.0_r8/(1.0_r8+ce)
435# ifdef MASKING
436!^ tl_ubar(i,Jend+1,kout)=tl_ubar(i,Jend+1,kout)* &
437!^ & GRID(ng)%umask(i,Jend+1)
438!^
439 ad_ubar(i,jend+1,kout)=ad_ubar(i,jend+1,kout)* &
440 & grid(ng)%umask(i,jend+1)
441# endif
442!^ tl_ubar(i,Jend+1,kout)=tl_cff2*(ubar(i,Jend+1,know)+ &
443!^ & Ce*ubar(i,Jend,kout))+ &
444!^ & cff2*(tl_ubar(i,Jend+1,know)+ &
445!^ & tl_Ce*ubar(i,Jend,kout)+ &
446!^ & Ce*tl_ubar(i,Jend,kout))
447!^
448 adfac=cff2*ad_ubar(i,jend+1,kout)
449 ad_ubar(i,jend ,kout)=ad_ubar(i,jend ,kout)+ce*adfac
450 ad_ubar(i,jend+1,know)=ad_ubar(i,jend+1,know)+adfac
451 ad_ce=ad_ce+ubar(i,jend,kout)*adfac
452 ad_cff2=ad_cff2+ &
453 & (ubar(i,jend+1,know)+ &
454 & ce*ubar(i,jend,kout))*ad_ubar(i,jend+1,kout)
455 ad_ubar(i,jend+1,kout)=0.0_r8
456!^ tl_cff2=-cff2*cff2*tl_Ce
457!^
458 ad_ce=ad_ce-cff2*cff2*ad_cff2
459 ad_cff2=0.0_r8
460!^ tl_Ce=cff*tl_cff1
461!^
462 ad_cff1=ad_cff1+cff*ad_ce
463 ad_ce=0.0_r8
464!^ tl_cff1=0.25_r8*g*(GRID(ng)%tl_h(i-1,Jend)+ &
465!^ & tl_zeta(i-1,Jend,know)+ &
466!^ & GRID(ng)%tl_h(i ,Jend)+ &
467!^ & tl_zeta(i ,Jend,know))/cff1
468!^
469 adfac=0.25_r8*g*ad_cff1/cff1
470 grid(ng)%ad_h(i-1,jend)=grid(ng)%ad_h(i-1,jend)+adfac
471 grid(ng)%ad_h(i ,jend)=grid(ng)%ad_h(i ,jend)+adfac
472 ad_zeta(i-1,jend,know)=ad_zeta(i-1,jend,know)+adfac
473 ad_zeta(i ,jend,know)=ad_zeta(i ,jend,know)+adfac
474 ad_cff1=0.0_r8
475 END IF
476 END DO
477!
478! Northern edge, clamped boundary condition.
479!
480 ELSE IF (ad_lbc(inorth,isubar,ng)%clamped) THEN
481 DO i=istru,iend
482 IF (lbc_apply(ng)%north(i)) THEN
483# ifdef MASKING
484!^ tl_ubar(i,Jend+1,kout)=tl_ubar(i,Jend+1,kout)* &
485!^ & GRID(ng)%umask(i,Jend+1)
486!^
487 ad_ubar(i,jend+1,kout)=ad_ubar(i,jend+1,kout)* &
488 & grid(ng)%umask(i,jend+1)
489# endif
490# ifdef ADJUST_BOUNDARY
491 IF (lobc(inorth,isubar,ng)) THEN
492!^ tl_ubar(i,Jend+1,kout)=BOUNDARY(ng)%tl_ubar_north(i)
493!^
494 boundary(ng)%ad_ubar_north(i)=boundary(ng)% &
495 & ad_ubar_north(i)+ &
496 & ad_ubar(i,jend+1,kout)
497 ad_ubar(i,jend+1,kout)=0.0_r8
498 ELSE
499!^ tl_ubar(i,Jend+1,kout)=0.0_r8
500!^
501 ad_ubar(i,jend+1,kout)=0.0_r8
502 END IF
503# else
504!^ tl_ubar(i,Jend+1,kout)=0.0_r8
505!^
506 ad_ubar(i,jend+1,kout)=0.0_r8
507# endif
508 END IF
509 END DO
510!
511! Northern edge, gradient boundary condition.
512!
513 ELSE IF (ad_lbc(inorth,isubar,ng)%gradient) THEN
514 DO i=istru,iend
515 IF (lbc_apply(ng)%north(i)) THEN
516# ifdef MASKING
517!^ tl_ubar(i,Jend+1,kout)=tl_ubar(i,Jend+1,kout)* &
518!^ & GRID(ng)%umask(i,Jend+1)
519!^
520 ad_ubar(i,jend+1,kout)=ad_ubar(i,jend+1,kout)* &
521 & grid(ng)%umask(i,jend+1)
522# endif
523!^ tl_ubar(i,Jend+1,kout)=tl_ubar(i,Jend,kout)
524!^
525 ad_ubar(i,jend,kout)=ad_ubar(i,jend,kout)+ &
526 & ad_ubar(i,jend+1,kout)
527 ad_ubar(i,jend+1,kout)=0.0_r8
528 END IF
529 END DO
530!
531! Northern edge, closed boundary condition: free slip (gamma2=1) or
532! no slip (gamma2=-1).
533!
534 ELSE IF (ad_lbc(inorth,isubar,ng)%closed) THEN
535 IF (ewperiodic(ng)) THEN
536 imin=istru
537 imax=iend
538 ELSE
539 imin=istr
540 imax=iendr
541 END IF
542 DO i=imin,imax
543 IF (lbc_apply(ng)%north(i)) THEN
544# ifdef MASKING
545!^ tl_ubar(i,Jend+1,kout)=tl_ubar(i,Jend+1,kout)* &
546!^ & GRID(ng)%umask(i,Jend+1)
547!^
548 ad_ubar(i,jend+1,kout)=ad_ubar(i,jend+1,kout)* &
549 & grid(ng)%umask(i,jend+1)
550# endif
551!^ tl_ubar(i,Jend+1,kout)=gamma2(ng)*tl_ubar(i,Jend,kout)
552!^
553 ad_ubar(i,jend ,kout)=ad_ubar(i,jend,kout)+ &
554 & gamma2(ng)*ad_ubar(i,jend+1,kout)
555 ad_ubar(i,jend+1,kout)=0.0_r8
556 END IF
557 END DO
558 END IF
559 END IF
560!
561!-----------------------------------------------------------------------
562! Lateral boundary conditions at the southern edge.
563!-----------------------------------------------------------------------
564!
565 IF (domain(ng)%Southern_Edge(tile)) THEN
566!
567! Southern edge, implicit upstream radiation condition.
568!
569 IF (ad_lbc(isouth,isubar,ng)%radiation) THEN
570 IF (iic(ng).ne.0) THEN
571 DO i=istru,iend
572 IF (lbc_apply(ng)%south(i)) THEN
573# if defined CELERITY_READ && defined FORWARD_READ
574 IF (ad_lbc(isouth,isubar,ng)%nudging) THEN
575 IF (lnudgem2clm(ng)) THEN
576 obc_out=0.5_r8* &
577 & (clima(ng)%M2nudgcof(i-1,jstr-1)+ &
578 & clima(ng)%M2nudgcof(i ,jstr-1))
579 obc_in =obcfac(ng)*obc_out
580 ELSE
581 obc_out=m2obc_out(ng,isouth)
582 obc_in =m2obc_in(ng,isouth)
583 END IF
584 IF (boundary(ng)%ubar_south_Ce(i).lt.0.0_r8) THEN
585 tau=obc_in
586 ELSE
587 tau=obc_out
588 END IF
589 tau=tau*dt2d
590 END IF
591# ifdef RADIATION_2D
592 cx=boundary(ng)%ubar_south_Cx(i)
593# else
594 cx=0.0_r8
595# endif
596 ce=boundary(ng)%ubar_south_Ce(i)
597 cff=boundary(ng)%ubar_south_C2(i)
598# endif
599# ifdef MASKING
600!^ tl_ubar(i,Jstr-1,kout)=tl_ubar(i,Jstr-1,kout)* &
601!^ & GRID(ng)%umask(i,Jstr-1)
602!^
603 ad_ubar(i,jstr-1,kout)=ad_ubar(i,jstr-1,kout)* &
604 & grid(ng)%umask(i,jstr-1)
605# endif
606 IF (ad_lbc(isouth,isubar,ng)%nudging) THEN
607!^ tl_ubar(i,Jstr-1,kout)=tl_ubar(i,Jstr-1,kout)- &
608!^ & tau*tl_ubar(i,Jstr-1,know)
609!^
610 ad_ubar(i,jstr-1,know)=ad_ubar(i,jstr-1,know)- &
611 & tau*ad_ubar(i,jstr-1,kout)
612 END IF
613!^ tl_ubar(i,Jstr-1,kout)=(cff*tl_ubar(i,Jstr-1,know)+ &
614!^ & Ce *tl_ubar(i,Jstr ,kout)- &
615!^ & MAX(Cx,0.0_r8)* &
616!^ & tl_grad(i-1,Jstr-1)- &
617!^ & MIN(Cx,0.0_r8)* &
618!^ & tl_grad(i ,Jstr-1))/ &
619!^ & (cff+Ce)
620!^
621 adfac=ad_ubar(i,jstr-1,kout)/(cff+ce)
622 ad_grad(i-1,jstr-1)=ad_grad(i-1,jstr-1)- &
623 & max(cx,0.0_r8)*adfac
624 ad_grad(i ,jstr-1)=ad_grad(i ,jstr-1)- &
625 & min(cx,0.0_r8)*adfac
626 ad_ubar(i,jstr-1,know)=ad_ubar(i,jstr-1,know)+cff*adfac
627 ad_ubar(i,jstr ,kout)=ad_ubar(i,jstr ,kout)+ce *adfac
628 ad_ubar(i,jstr-1,kout)=0.0_r8
629 END IF
630 END DO
631 END IF
632!
633! Southern edge, Chapman boundary condition.
634!
635 ELSE IF (ad_lbc(isouth,isubar,ng)%Flather.or. &
636 & ad_lbc(isouth,isubar,ng)%Shchepetkin.or. &
637 & ad_lbc(isouth,isubar,ng)%reduced) THEN
638 DO i=istru,iend
639 IF (lbc_apply(ng)%south(i)) THEN
640 cff=dt2d*0.5_r8*(grid(ng)%pn(i-1,jstr)+ &
641 & grid(ng)%pn(i ,jstr))
642 cff1=sqrt(g*0.5_r8*(grid(ng)%h(i-1,jstr)+ &
643 & zeta(i-1,jstr,know)+ &
644 & grid(ng)%h(i ,jstr)+ &
645 & zeta(i ,jstr,know)))
646 ce=cff*cff1
647 cff2=1.0_r8/(1.0_r8+ce)
648# ifdef MASKING
649!^ tl_ubar(i,Jstr-1,kout)=tl_ubar(i,Jstr-1,kout)* &
650!^ & GRID(ng)%umask(i,Jstr-1)
651!^
652 ad_ubar(i,jstr-1,kout)=ad_ubar(i,jstr-1,kout)* &
653 & grid(ng)%umask(i,jstr-1)
654# endif
655!^ tl_ubar(i,Jstr-1,kout)=tl_cff2*(ubar(i,Jstr-1,know)+ &
656!^ & Ce*ubar(i,Jstr,kout))+ &
657!^ & cff2*(tl_ubar(i,Jstr-1,know)+ &
658!^ & tl_Ce*ubar(i,Jstr,kout)+ &
659!^ & Ce*tl_ubar(i,Jstr,kout))
660!^
661 adfac=cff2*ad_ubar(i,jstr-1,kout)
662 ad_ubar(i,jstr-1,know)=ad_ubar(i,jstr-1,know)+adfac
663 ad_ubar(i,jstr ,kout)=ad_ubar(i,jstr ,kout)+ce*adfac
664 ad_ce=ad_ce+ubar(i,jstr,kout)*adfac
665 ad_cff2=ad_cff2+ &
666 & (ubar(i,jstr-1,know)+ &
667 & ce*ubar(i,jstr,kout))*ad_ubar(i,jstr-1,kout)
668 ad_ubar(i,jstr-1,kout)=0.0_r8
669!^ tl_cff2=-cff2*cff2*tl_Ce
670!^
671 ad_ce=ad_ce-cff2*cff2*ad_cff2
672 ad_cff2=0.0_r8
673!^ tl_Ce=cff*tl_cff1
674!^
675 ad_cff1=ad_cff1+cff*ad_ce
676 ad_ce=0.0_r8
677!^ tl_cff1=0.25_r8*g*(GRID(ng)%tl_h(i-1,Jstr)+ &
678!^ & tl_zeta(i-1,Jstr,know)+ &
679!^ & GRID(ng)%tl_h(i ,Jstr)+ &
680!^ & tl_zeta(i ,Jstr,know))/cff1
681!^
682 adfac=0.25_r8*g*ad_cff1/cff1
683 grid(ng)%ad_h(i-1,jstr)=grid(ng)%ad_h(i-1,jstr)+adfac
684 grid(ng)%ad_h(i ,jstr)=grid(ng)%ad_h(i ,jstr)+adfac
685 ad_zeta(i-1,jstr,know)=ad_zeta(i-1,jstr,know)+adfac
686 ad_zeta(i ,jstr,know)=ad_zeta(i ,jstr,know)+adfac
687 ad_cff1=0.0_r8
688 END IF
689 END DO
690!
691! Southern edge, clamped boundary condition.
692!
693 ELSE IF (ad_lbc(isouth,isubar,ng)%clamped) THEN
694 DO i=istru,iend
695 IF (lbc_apply(ng)%south(i)) THEN
696# ifdef MASKING
697!^ tl_ubar(i,Jstr-1,kout)=tl_ubar(i,Jstr-1,kout)* &
698!^ & GRID(ng)%umask(i,Jstr-1)
699!^
700 ad_ubar(i,jstr-1,kout)=ad_ubar(i,jstr-1,kout)* &
701 & grid(ng)%umask(i,jstr-1)
702# endif
703# ifdef ADJUST_BOUNDARY
704 IF (lobc(isouth,isubar,ng)) THEN
705!^ tl_ubar(i,Jstr-1,kout)=BOUNDARY(ng)%tl_ubar_south(i)
706!^
707 boundary(ng)%ad_ubar_south(i)=boundary(ng)% &
708 & ad_ubar_south(i)+ &
709 & ad_ubar(i,jstr-1,kout)
710 ad_ubar(i,jstr-1,kout)=0.0_r8
711 ELSE
712!^ tl_ubar(i,Jstr-1,kout)=0.0_r8
713!^
714 ad_ubar(i,jstr-1,kout)=0.0_r8
715 END IF
716# else
717!^ tl_ubar(i,Jstr-1,kout)=0.0_r8
718!^
719 ad_ubar(i,jstr-1,kout)=0.0_r8
720# endif
721 END IF
722 END DO
723!
724! Southern edge, gradient boundary condition.
725!
726 ELSE IF (ad_lbc(isouth,isubar,ng)%gradient) THEN
727 DO i=istru,iend
728 IF (lbc_apply(ng)%south(i)) THEN
729# ifdef MASKING
730!^ tl_ubar(i,Jstr-1,kout)=tl_ubar(i,Jstr-1,kout)* &
731!^ & GRID(ng)%umask(i,Jstr-1)
732!^
733 ad_ubar(i,jstr-1,kout)=ad_ubar(i,jstr-1,kout)* &
734 & grid(ng)%umask(i,jstr-1)
735# endif
736!^ tl_ubar(i,Jstr-1,kout)=tl_ubar(i,Jstr,kout)
737!^
738 ad_ubar(i,jstr ,kout)=ad_ubar(i,jstr ,kout)+ &
739 & ad_ubar(i,jstr-1,kout)
740 ad_ubar(i,jstr-1,kout)=0.0_r8
741 END IF
742 END DO
743!
744! Southern edge, closed boundary condition: free slip (gamma2=1) or
745! no slip (gamma2=-1).
746!
747 ELSE IF (ad_lbc(isouth,isubar,ng)%closed) THEN
748 IF (ewperiodic(ng)) THEN
749 imin=istru
750 imax=iend
751 ELSE
752 imin=istr
753 imax=iendr
754 END IF
755 DO i=imin,imax
756 IF (lbc_apply(ng)%south(i)) THEN
757# ifdef MASKING
758!^ tl_ubar(i,Jstr-1,kout)=tl_ubar(i,Jstr-1,kout)* &
759!^ & GRID(ng)%umask(i,Jstr-1)
760!^
761 ad_ubar(i,jstr-1,kout)=ad_ubar(i,jstr-1,kout)* &
762 & grid(ng)%umask(i,jstr-1)
763# endif
764!^ tl_ubar(i,Jstr-1,kout)=gamma2(ng)*tl_ubar(i,Jstr,kout)
765!^
766 ad_ubar(i,jstr ,kout)=ad_ubar(i,jstr,kout)+ &
767 & gamma2(ng)*ad_ubar(i,jstr-1,kout)
768 ad_ubar(i,jstr-1,kout)=0.0_r8
769 END IF
770 END DO
771 END IF
772 END IF
773!
774!-----------------------------------------------------------------------
775! Lateral boundary conditions at the eastern edge.
776!-----------------------------------------------------------------------
777!
778 IF (domain(ng)%Eastern_Edge(tile)) THEN
779!
780! Eastern edge, implicit upstream radiation condition.
781!
782 IF (ad_lbc(ieast,isubar,ng)%radiation) THEN
783 IF (iic(ng).ne.0) THEN
784 DO j=jstr,jend
785 IF (lbc_apply(ng)%east(j)) THEN
786# if defined CELERITY_READ && defined FORWARD_READ
787 IF (ad_lbc(ieast,isubar,ng)%nudging) THEN
788 IF (lnudgem2clm(ng)) THEN
789 obc_out=0.5_r8* &
790 & (clima(ng)%M2nudgcof(iend ,j)+ &
791 & clima(ng)%M2nudgcof(iend+1,j))
792 obc_in =obcfac(ng)*obc_out
793 ELSE
794 obc_out=m2obc_out(ng,ieast)
795 obc_in =m2obc_in(ng,ieast)
796 END IF
797 IF (boundary(ng)%ubar_east_Cx(j).lt.0.0_r8) THEN
798 tau=obc_in
799 ELSE
800 tau=obc_out
801 END IF
802 tau=tau*dt2d
803 END IF
804 cx=boundary(ng)%ubar_east_Cx(j)
805# ifdef RADIATION_2D
806 ce=boundary(ng)%ubar_east_Ce(j)
807# else
808 ce=0.0_r8
809# endif
810 cff=boundary(ng)%ubar_east_C2(j)
811# endif
812# ifdef MASKING
813!^ tl_ubar(Iend+1,j,kout)=tl_ubar(Iend+1,j,kout)* &
814!^ & GRID(ng)%umask(Iend+1,j)
815!^
816 ad_ubar(iend+1,j,kout)=ad_ubar(iend+1,j,kout)* &
817 & grid(ng)%umask(iend+1,j)
818# endif
819 IF (ad_lbc(ieast,isubar,ng)%nudging) THEN
820!^ tl_ubar(Iend+1,j,kout)=tl_ubar(Iend+1,j,kout)- &
821!^ & tau*tl_ubar(Iend+1,j,know)
822!^
823 ad_ubar(iend+1,j,know)=ad_ubar(iend+1,j,know)- &
824 & tau*ad_ubar(iend+1,j,kout)
825 END IF
826!^ tl_ubar(Iend+1,j,kout)=(cff*tl_ubar(Iend+1,j,know)+ &
827!^ & Cx *tl_ubar(Iend ,j,kout)- &
828!^ & MAX(Ce,0.0_r8)* &
829!^ & tl_grad(Iend+1,j )- &
830!^ & MIN(Ce,0.0_r8)* &
831!^ & tl_grad(Iend+1,j+1))/ &
832!^ & (cff+Cx)
833!^
834 adfac=ad_ubar(iend+1,j,kout)/(cff+cx)
835 ad_grad(iend+1,j )=ad_grad(iend+1,j )- &
836 & max(ce,0.0_r8)*adfac
837 ad_grad(iend+1,j+1)=ad_grad(iend+1,j+1)- &
838 & min(ce,0.0_r8)*adfac
839 ad_ubar(iend ,j,kout)=ad_ubar(iend ,j,kout)+cx *adfac
840 ad_ubar(iend+1,j,know)=ad_ubar(iend+1,j,know)+cff*adfac
841 ad_ubar(iend+1,j,kout)=0.0_r8
842 END IF
843 END DO
844 END IF
845!
846! Eastern edge, Flather boundary condition.
847!
848 ELSE IF (ad_lbc(ieast,isubar,ng)%Flather) THEN
849 DO j=jstr,jend
850 IF (lbc_apply(ng)%east(j)) THEN
851 cff=1.0_r8/(0.5_r8*(grid(ng)%h(iend ,j)+ &
852 & zeta(iend ,j,know)+ &
853 & grid(ng)%h(iend+1,j)+ &
854 & zeta(iend+1,j,know)))
855 cx=sqrt(g*cff)
856# ifdef MASKING
857!^ tl_ubar(Iend+1,j,kout)=tl_ubar(Iend+1,j,kout)* &
858!^ & GRID(ng)%umask(Iend+1,j)
859!^
860 ad_ubar(iend+1,j,kout)=ad_ubar(iend+1,j,kout)* &
861 & grid(ng)%umask(iend+1,j)
862# endif
863# ifdef ADJUST_BOUNDARY
864 IF (lobc(ieast,isubar,ng)) THEN
865!^ tl_ubar(Iend+1,j,kout)=tl_ubar(Iend+1,j,kout)- &
866!^ & Cx*BOUNDARY(ng)%tl_zeta_east(j)
867!^
868 boundary(ng)%ad_zeta_east(j)= &
869 & boundary(ng)%ad_zeta_east(j)- &
870 & cx*ad_ubar(iend+1,j,kout)
871 END IF
872# endif
873# if defined ATM_PRESS && defined PRESS_COMPENSATE
874!^ tl_ubar(Iend+1,j,kout)=tl_bry_val+ &
875!^ & tl_Cx* &
876!^ & (0.5_r8* &
877!^ & (zeta(Iend ,j,know)+ &
878!^ & zeta(Iend+1,j,know)+ &
879!^ & fac*(FORCES(ng)%Pair(Iend ,j)+ &
880!^ & FORCES(ng)%Pair(Iend+1,j)- &
881!^ & 2.0_r8*OneAtm))- &
882!^ & BOUNDARY(ng)%zeta_east(j))+ &
883!^ & Cx* &
884!^ & (0.5_r8*(tl_zeta(Iend ,j,know)+ &
885!^ & tl_zeta(Iend+1,j,know)))
886!^
887 adfac=cx*0.5_r8*ad_ubar(iend+1,j,kout)
888 ad_zeta(iend ,j,know)=ad_zeta(iend ,j,know)+adfac
889 ad_zeta(iend+1,j,know)=ad_zeta(iend+1,j,know)+adfac
890 ad_cx=ad_cx+ &
891 & (0.5_r8*(zeta(iend ,j,know)+ &
892 & zeta(iend+1,j,know)+ &
893 & fac*(forces(ng)%Pair(iend ,j)+ &
894 & forces(ng)%Pair(iend+1,j)- &
895 & 2.0_r8*oneatm))- &
896 & boundary(ng)%zeta_east(j))*ad_ubar(iend+1,j,kout)
897 ad_bry_val=ad_bry_val+ad_ubar(iend+1,j,kout)
898 ad_ubar(iend+1,j,kout)=0.0_r8
899# else
900!^ tl_ubar(Iend+1,j,kout)=tl_bry_val+ &
901!^ & tl_Cx* &
902!^ & (0.5_r8*(zeta(Iend ,j,know)+ &
903!^ & zeta(Iend+1,j,know))- &
904!^ & BOUNDARY(ng)%zeta_east(j))+ &
905!^ & Cx* &
906!^ (0.5_r8*(tl_zeta(Iend ,j,know)+ &
907!^ & tl_zeta(Iend+1,j,know)))
908!^
909 adfac=cx*0.5_r8*ad_ubar(iend+1,j,kout)
910 ad_zeta(iend ,j,know)=ad_zeta(iend ,j,know)+adfac
911 ad_zeta(iend+1,j,know)=ad_zeta(iend+1,j,know)+adfac
912 ad_cx=ad_cx+ &
913 & (0.5_r8*(zeta(iend ,j,know)+ &
914 & zeta(iend+1,j,know))- &
915 & boundary(ng)%zeta_east(j))*ad_ubar(iend+1,j,kout)
916 ad_bry_val=ad_bry_val+ad_ubar(iend+1,j,kout)
917 ad_ubar(iend+1,j,kout)=0.0_r8
918# endif
919!^ tl_Cx=0.5_r8*g*tl_cff/Cx
920!^
921 ad_cff=ad_cff+0.5_r8*g*ad_cx/cx
922 ad_cx=0.0_r8
923!^ tl_cff=-cff*cff*(0.5_r8*(GRID(ng)%tl_h(Iend ,j)+ &
924!^ & tl_zeta(Iend ,j,know)+ &
925!^ & GRID(ng)%tl_h(Iend+1,j)+ &
926!^ & tl_zeta(Iend+1,j,know)))
927!^
928 adfac=-cff*cff*0.5_r8*ad_cff
929 grid(ng)%ad_h(iend ,j)=grid(ng)%ad_h(iend ,j)+adfac
930 grid(ng)%ad_h(iend+1,j)=grid(ng)%ad_h(iend+1,j)+adfac
931 ad_zeta(iend ,j,know)=ad_zeta(iend ,j,know)+adfac
932 ad_zeta(iend+1,j,know)=ad_zeta(iend+1,j,know)+adfac
933 ad_cff=0.0_r8
934
935# if defined SSH_TIDES_NOT_YET && !defined UV_TIDES_NOT_YET
936 IF (ad_lbc(ieast,isfsur,ng)%acquire) THEN
937 bry_pgr=-g*(boundary(ng)%zeta_east(j)- &
938 & zeta(iend,j,know))* &
939 & 0.5_r8*grid(ng)%pm(iend,j)
940 ELSE
941 bry_pgr=-g*(zeta(iend+1,j,know)- &
942 & zeta(iend ,j,know))* &
943 & 0.5_r8*(grid(ng)%pm(iend ,j)+ &
944 & grid(ng)%pm(iend+1,j))
945 END IF
946# ifdef UV_COR
947 bry_cor=0.125_r8*(vbar(iend ,j ,know)+ &
948 & vbar(iend ,j+1,know)+ &
949 & vbar(iend+1,j ,know)+ &
950 & vbar(iend+1,j+1,know))* &
951 & (grid(ng)%f(iend ,j)+ &
952 & grid(ng)%f(iend+1,j))
953# else
954 bry_cor=0.0_r8
955# endif
956 cff1=1.0_r8/(0.5_r8*(grid(ng)%h(iend ,j)+ &
957 & zeta(iend ,j,know)+ &
958 & grid(ng)%h(iend+1,j)+ &
959 & zeta(iend+1,j,know)))
960 bry_str=cff1*(forces(ng)%sustr(iend+1,j)- &
961 & forces(ng)%bustr(iend+1,j))
962 cx=1.0_r8/sqrt(g*0.5_r8*(grid(ng)%h(iend+1,j)+ &
963 & zeta(iend+1,j,know)+ &
964 & grid(ng)%h(iend ,j)+ &
965 & zeta(iend ,j,know)))
966 cff2=grid(ng)%om_u(iend+1,j)*cx
967!^ tl_bry_val=tl_ubar(Iend,j,know)+ &
968!^ & tl_cff2*(bry_pgr+ &
969!^ & bry_cor+ &
970!^ & bry_str)+ &
971!^ & cff2*(tl_bry_pgr+ &
972!^ & tl_bry_cor+ &
973!^ & tl_bry_str)
974!^
975 adfac=cff2*ad_bry_val
976 ad_bry_pgr=ad_bry_pgr+adfac
977 ad_bry_cor=ad_bry_cor+adfac
978 ad_bry_str=ad_bry_str+adfac
979 ad_cff2=ad_cff2+(bry_pgr+ &
980 & bry_cor+ &
981 & bry_str)*ad_bry_val
982 ad_ubar(iend,j,know)=ad_ubar(iend,j,know)+ad_bry_val
983 ad_bry_val=0.0_r8
984!^ tl_cff2=GRID(ng)%om_u(Iend+1,j)*tl_Cx
985!^
986 ad_cx=ad_cx+grid(ng)%om_u(iend+1,j)*ad_cff2
987 ad_cff2=0.0_r8
988!^ tl_Cx=-Cx*Cx*Cx*0.25_r8*g*(GRID(ng)%tl_h(Iend+1,j)+ &
989!^ & tl_zeta(Iend+1,j,know)+ &
990!^ & GRID(ng)%tl_h(Iend ,j)+ &
991!^ & tl_zeta(Iend ,j,know))
992!^
993 adfac=-cx*cx*cx*0.25_r8*g*ad_cx
994 ad_zeta(iend ,j,know)=ad_zeta(iend ,j,know)+adfac
995 ad_zeta(iend+1,j,know)=ad_zeta(iend+1,j,know)+adfac
996 grid(ng)%ad_h(iend ,j)=grid(ng)%ad_h(iend ,j)+adfac
997 grid(ng)%ad_h(iend+1,j)=grid(ng)%ad_h(iend+1,j)+adfac
998 ad_cx=0.0_r8
999!^ tl_bry_str=tl_cff1*(FORCES(ng)%sustr(Iend+1,j)- &
1000!^ & FORCES(ng)%bustr(Iend+1,j))+ &
1001!^ & cff1*(FORCES(ng)%tl_sustr(Iend+1,j)- &
1002!^ & FORCES(ng)%tl_bustr(Iend+1,j))
1003!^
1004 adfac=cff1*ad_bry_str
1005 forces(ng)%ad_sustr(iend+1,j)= &
1006 & forces(ng)%ad_sustr(iend+1,j)+ &
1007 & adfac
1008 forces(ng)%tl_bustr(iend+1,j)= &
1009 & forces(ng)%tl_bustr(iend+1,j)- &
1010 & adfac
1011 ad_cff1=ad_cff1+(forces(ng)%sustr(iend+1,j)- &
1012 & forces(ng)%bustr(iend+1,j))*ad_bry_str
1013 ad_bry_str=0.0_r8
1014!^ tl_cff1=-cff1*cff1*0.5_r8*(GRID(ng)%tl_h(Iend ,j)+ &
1015!^ & tl_zeta(Iend ,j,know)+ &
1016!^ & GRID(ng)%tl_h(Iend+1,j)+ &
1017!^ & tl_zeta(Iend+1,j,know))
1018!^
1019 adfac=-cff1*cff1*0.5_r8*ad_cff1
1020 ad_zeta(iend ,j,know)=ad_zeta(iend ,j,know)+adfac
1021 ad_zeta(iend+1,j,know)=ad_zeta(iend+1,j,know)+adfac
1022 grid(ng)%ad_h(iend ,j)=grid(ng)%ad_h(iend ,j)+adfac
1023 grid(ng)%ad_h(iend+1,j)=grid(ng)%ad_h(iend+1,j)+adfac
1024 ad_cff1=0.0_r8
1025# ifdef UV_COR
1026!^ tl_bry_cor=0.125_r8*(tl_vbar(Iend ,j ,know)+ &
1027!^ & tl_vbar(Iend ,j+1,know)+ &
1028!^ & tl_vbar(Iend+1,j ,know)+ &
1029!^ & tl_vbar(Iend+1,j+1,know))* &
1030!^ & (GRID(ng)%f(Iend ,j)+ &
1031!^ & GRID(ng)%f(Iend+1,j))
1032!^
1033 adfac=0.125_r8*(grid(ng)%f(iend ,j)+ &
1034 & grid(ng)%f(iend+1,j))*ad_bry_cor
1035 ad_vbar(iend ,j ,know)=ad_vbar(iend ,j ,know)+adfac
1036 ad_vbar(iend+1,j ,know)=ad_vbar(iend+1,j ,know)+adfac
1037 ad_vbar(iend ,j+1,know)=ad_vbar(iend ,j+1,know)+adfac
1038 ad_vbar(iend+1,j+1,know)=ad_vbar(iend+1,j+1,know)+adfac
1039 ad_bry_cor=0.0_r8
1040# else
1041!^ tl_bry_cor=0.0_r8
1042!^
1043 ad_bry_cor=0.0_r8
1044# endif
1045 IF (ad_lbc(ieast,isfsur,ng)%acquire) THEN
1046# ifdef ADJUST_BOUNDARY
1047 IF (lobc(ieast,isubar,ng)) THEN
1048!^ tl_bry_pgr=tl_bry_pgr- &
1049!^ & g*BOUNDARY(ng)%tl_zeta_east(j)* &
1050!^ & 0.5_r8*GRID(ng)%pm(Iend,j)
1051!^
1052 boundary(ng)%ad_zeta_east(j)=boundary(ng)% &
1053 & ad_zeta_east(j)- &
1054 & g*0.5_r8* &
1055 & grid(ng)%pm(iend,j)* &
1056 & ad_bry_pgr
1057 END IF
1058# endif
1059!^ tl_bry_pgr=g*tl_zeta(Iend,j,know)* &
1060!^ & 0.5_r8*GRID(ng)%pm(Iend,j)
1061!^
1062 ad_zeta(iend,j,know)=ad_zeta(iend,j,know)+ &
1063 & g*0.5_r8*grid(ng)%pm(iend,j)* &
1064 & ad_bry_pgr
1065 ad_bry_pgr=0.0_r8
1066 ELSE
1067!^ tl_bry_pgr=-g*(tl_zeta(Iend+1,j,know)- &
1068!^ & tl_zeta(Iend ,j,know))* &
1069!^ & 0.5_r8*(GRID(ng)%pm(Iend ,j)+ &
1070!^ & GRID(ng)%pm(Iend+1,j))
1071!^
1072 adfac=-g*0.5_r8*(grid(ng)%pm(iend ,j)+ &
1073 & grid(ng)%pm(iend+1,j))*ad_bry_pgr
1074 ad_zeta(iend ,j,know)=ad_zeta(iend ,j,know)-adfac
1075 ad_zeta(iend+1,j,know)=ad_zeta(iend+1,j,know)+adfac
1076 ad_bry_pgr=0.0_r8
1077 END IF
1078# else
1079# ifdef ADJUST_BOUNDARY
1080 IF (lobc(ieast,isubar,ng)) THEN
1081!^ tl_bry_val=BOUNDARY(ng)%tl_ubar_east(j)
1082!^
1083 boundary(ng)%ad_ubar_east(j)= &
1084 & boundary(ng)%ad_ubar_east(j)+ &
1085 & ad_bry_val
1086 ad_bry_val=0.0_r8
1087 ELSE
1088!^ tl_bry_val=0.0_r8
1089!^
1090 ad_bry_val=0.0_r8
1091 END IF
1092# else
1093!^ tl_bry_val=0.0_r8
1094!^
1095 ad_bry_val=0.0_r8
1096# endif
1097# endif
1098 END IF
1099 END DO
1100!
1101! Eastern edge, Shchepetkin boundary condition (Maison et al., 2010).
1102!
1103 ELSE IF (ad_lbc(ieast,isubar,ng)%Shchepetkin) THEN
1104 DO j=jstr,jend
1105 IF (lbc_apply(ng)%east(j)) THEN
1106# ifdef WET_DRY_NOT_YET
1107 cff=0.5_r8*(grid(ng)%h(iend ,j)+ &
1108 & zeta(iend ,j,know)+ &
1109 & grid(ng)%h(iend+1,j)+ &
1110 & zeta(iend+1,j,know))
1111# else
1112 cff=0.5_r8*(grid(ng)%h(iend ,j)+ &
1113 & grid(ng)%h(iend+1,j))
1114# endif
1115 cff1=sqrt(g/cff)
1116 cx=dt2d*cff1*cff*0.5_r8*(grid(ng)%pm(iend ,j)+ &
1117 & grid(ng)%pm(iend+1,j))
1118 zx=(0.5_r8+cx)*zeta(iend ,j,know)+ &
1119 & (0.5_r8-cx)*zeta(iend+1,j,know)
1120 IF (cx.gt.co) THEN
1121 cff2=(1.0_r8-co/cx)**2
1122 cff3=zeta(iend,j,kout)+ &
1123 & cx*zeta(iend+1,j,know)- &
1124 & (1.0_r8+cx)*zeta(iend,j,know)
1125 zx=zx+cff2*cff3
1126 END IF
1127# ifdef MASKING
1128!^ tl_ubar(Iend+1,j,kout)=tl_ubar(Iend+1,j,kout)* &
1129!^ & GRID(ng)%umask(Iend+1,j)
1130!^
1131 ad_ubar(iend+1,j,kout)=ad_ubar(iend+1,j,kout)* &
1132 & grid(ng)%umask(iend+1,j)
1133# endif
1134# ifdef ADJUST_BOUNDARY
1135 IF (lobc(ieast,isubar,ng)) THEN
1136!^ tl_ubar(Iend+1,j,kout)=tl_ubar(Iend+1,j,kout)- &
1137!^ & 0.5_r8*cff1* &
1138!^ & BOUNDARY(ng)%tl_zeta_east(j)
1139!^
1140 boundary(ng)%ad_zeta_east(j)= &
1141 & boundary(ng)%ad_zeta_east(j)- &
1142 & 0.5_r8*cff1* &
1143 & ad_ubar(iend+1,j,kout)
1144 END IF
1145# endif
1146!^ tl_ubar(Iend+1,j,kout)=0.5_r8* &
1147!^ & ((1.0_r8-Cx)* &
1148!^ & tl_ubar(Iend+1,j,know)+ &
1149!^ & tl_Cx*(ubar(Iend ,j,know)- &
1150!^ & ubar(Iend+1,j,know))+ &
1151!^ & Cx*tl_ubar(Iend,j,know)+ &
1152!^ & tl_bry_val+ &
1153!^ & tl_cff1* &
1154!^ & (Zx-BOUNDARY(ng)%zeta_east(j))- &
1155!^ & cff1*tl_Zx)
1156!^
1157 adfac=0.5_r8*ad_ubar(iend+1,j,kout)
1158 ad_ubar(iend+1,j,know)=ad_ubar(iend+1,j,know)+ &
1159 & (1.0_r8-cx)*adfac
1160 ad_ubar(iend ,j,know)=ad_ubar(iend ,j,know)+ &
1161 & cx*adfac
1162 ad_cx=ad_cx+ &
1163 & (ubar(iend ,j,know)- &
1164 & ubar(iend+1,j,know))*adfac
1165 ad_bry_val=ad_bry_val+adfac
1166 ad_cff1=ad_cff1+ &
1167 & (zx-boundary(ng)%zeta_east(j))*adfac
1168 ad_zx=ad_zx-cff1*adfac
1169 ad_ubar(iend+1,j,kout)=0.0_r8
1170 IF (cx.gt.co) THEN
1171!^ tl_Zx=tl_Zx+cff2*tl_cff3+ &
1172!^ & tl_cff2*cff3
1173!^
1174 ad_cff2=ad_cff2+cff3*ad_zx
1175 ad_cff3=ad_cff3+cff2*ad_zx
1176!^ tl_cff3=tl_zeta(Iend,j,kout)+ &
1177!^ & Cx*tl_zeta(Iend+1,j,know)+ &
1178!^ & tl_Cx*(zeta(Iend ,j,know)+ &
1179!^ & zeta(Iend+1,j,know))- &
1180!^ & (1.0_r8+Cx)*tl_zeta(Iend,j,know)
1181!^
1182 ad_zeta(iend ,j,kout)=ad_zeta(iend ,j,kout)+ &
1183 & ad_cff3
1184 ad_zeta(iend ,j,know)=ad_zeta(iend ,j,know)- &
1185 & (1.0_r8+cx)*ad_cff3
1186 ad_zeta(iend+1,j,know)=ad_zeta(iend+1,j,know)+ &
1187 & cx*ad_cff3
1188 ad_cx=ad_cx+ &
1189 & (zeta(iend ,j,know)+ &
1190 & zeta(iend+1,j,know))*ad_cff3
1191 ad_cff3=0.0_r8
1192!^ tl_cff2=2.0_r8*cff2*Co*tl_Cx/(Cx*Cx)
1193!^
1194 ad_cx=ad_cx+ &
1195 & 2.0_r8*cff2*co*ad_cff2/(cx*cx)
1196 ad_cff2=0.0_r8
1197 END IF
1198!^ tl_Zx=(0.5_r8+Cx)*tl_zeta(Iend ,j,know)+ &
1199!^ & (0.5_r8-Cx)*tl_zeta(Iend+1,j,know)+ &
1200!^ & tl_Cx*(zeta(Iend ,j,know)- &
1201!^ & zeta(Iend+1,j,know))
1202!^
1203 ad_zeta(iend ,j,know)=ad_zeta(iend ,j,know)+ &
1204 & (0.5_r8+cx)*ad_zx
1205 ad_zeta(iend+1,j,know)=ad_zeta(iend+1,j,know)+ &
1206 & (0.5_r8-cx)*ad_zx
1207 ad_cx=ad_cx+ &
1208 & (zeta(iend ,j,know)- &
1209 & zeta(iend+1,j,know))*ad_zx
1210 ad_zx=0.0_r8
1211!^ tl_Cx=dt2d*0.5_r8*(GRID(ng)%pm(Iend ,j)+ &
1212!^ & GRID(ng)%pm(Iend+1,j))* &
1213!^ & (cff1*tl_cff+ &
1214!^ & tl_cff1*cff)
1215!^
1216 adfac=dt2d*0.5_r8*(grid(ng)%pm(iend ,j)+ &
1217 & grid(ng)%pm(iend+1,j))*ad_cx
1218 ad_cff=ad_cff+cff1*adfac
1219 ad_cff1=ad_cff1+cff*adfac
1220 ad_cx=0.0_r8
1221!^ tl_cff1=-0.5_r8*cff1*tl_cff/cff
1222!^
1223 ad_cff=ad_cff- &
1224 & 0.5_r8*cff1*ad_cff1/cff
1225 ad_cff1=0.0_r8
1226# ifdef WET_DRY_NOT_YET
1227!^ tl_cff=0.5_r8*(GRID(ng)%tl_h(Iend ,j)+ &
1228!^ & tl_zeta(Iend ,j,know)+ &
1229!^ & GRID(ng)%tl_h(Iend+1,j)+ &
1230!^ & tl_zeta(Iend+1,j,know))
1231!^
1232 adfac=0.5_r8*ad_cff
1233 grid(ng)%ad_h(iend ,j)=grid(ng)%ad_h(iend ,j)+adfac
1234 grid(ng)%ad_h(iend+1,j)=grid(ng)%ad_h(iend+1,j)+adfac
1235 ad_zeta(iend ,j,know)=ad_zeta(iend ,j,know)+adfac
1236 ad_zeta(iend+1,j,know)=tl_zeta(iend+1,j,know)+adfac
1237 ad_cff=0.0_r8
1238# else
1239!^ tl_cff=0.5_r8*(GRID(ng)%tl_h(Iend ,j)+ &
1240!^ & GRID(ng)%tl_h(Iend+1,j))
1241!^
1242 adfac=0.5_r8*ad_cff
1243 grid(ng)%ad_h(iend ,j)=grid(ng)%ad_h(iend ,j)+adfac
1244 grid(ng)%ad_h(iend+1,j)=grid(ng)%ad_h(iend+1,j)+adfac
1245 ad_cff=0.0_r8
1246# endif
1247# if defined SSH_TIDES_NOT_YET && !defined UV_TIDES_NOT_YET
1248 IF (ad_lbc(ieast,isfsur,ng)%acquire) THEN
1249 bry_pgr=-g*(boundary(ng)%zeta_east(j)- &
1250 & zeta(iend,j,know))* &
1251 & 0.5_r8*grid(ng)%pm(iend,j)
1252 ELSE
1253 bry_pgr=-g*(zeta(iend+1,j,know)- &
1254 & zeta(iend ,j,know))* &
1255 & 0.5_r8*(grid(ng)%pm(iend ,j)+ &
1256 & grid(ng)%pm(iend+1,j))
1257 END IF
1258# ifdef UV_COR
1259 bry_cor=0.125_r8*(vbar(iend ,j ,know)+ &
1260 & vbar(iend ,j+1,know)+ &
1261 & vbar(iend+1,j ,know)+ &
1262 & vbar(iend+1,j+1,know))* &
1263 & (grid(ng)%f(iend ,j)+ &
1264 & grid(ng)%f(iend+1,j))
1265# else
1266 bry_cor=0.0_r8
1267# endif
1268 cff1=1.0_r8/(0.5_r8*(grid(ng)%h(iend ,j)+ &
1269 & zeta(iend ,j,know)+ &
1270 & grid(ng)%h(iend+1,j)+ &
1271 & zeta(iend+1,j,know)))
1272 bry_str=cff1*(forces(ng)%sustr(iend+1,j)- &
1273 & forces(ng)%bustr(iend+1,j))
1274 cx=1.0_r8/sqrt(g*0.5_r8*(grid(ng)%h(iend+1,j)+ &
1275 & zeta(iend+1,j,know)+ &
1276 & grid(ng)%h(iend ,j)+ &
1277 & zeta(iend ,j,know)))
1278 cff2=grid(ng)%om_u(iend+1,j)*cx
1279!^ tl_bry_val=tl_ubar(Iend,j,know)+ &
1280!^ & tl_cff2*(bry_pgr+ &
1281!^ & bry_cor+ &
1282!^ & bry_str)+ &
1283!^ & cff2*(tl_bry_pgr+ &
1284!^ & tl_bry_cor+ &
1285!^ & tl_bry_str)
1286!^
1287 adfac=cff2*ad_bry_val
1288 ad_bry_pgr=ad_bry_pgr+adfac
1289 ad_bry_cor=ad_bry_cor+adfac
1290 ad_bry_str=ad_bry_str+adfac
1291 ad_cff2=ad_cff2+(bry_pgr+ &
1292 & bry_cor+ &
1293 & bry_str)*ad_bry_val
1294 ad_ubar(iend,j,know)=ad_ubar(iend,j,know)+ad_bry_val
1295 ad_bry_val=0.0_r8
1296!^ tl_cff2=GRID(ng)%om_u(Iend+1,j)*tl_Cx
1297!^
1298 ad_cx=ad_cx+grid(ng)%om_u(iend+1,j)*ad_cff2
1299 ad_cff2=0.0_r8
1300!^ tl_Cx=-Cx*Cx*Cx*0.25_r8*g*(GRID(ng)%tl_h(Iend+1,j)+ &
1301!^ & tl_zeta(Iend+1,j,know)+ &
1302!^ & GRID(ng)%tl_h(Iend ,j)+ &
1303!^ & tl_zeta(Iend ,j,know))
1304!^
1305 adfac=-cx*cx*cx*0.25_r8*g*ad_cx
1306 ad_zeta(iend ,j,know)=ad_zeta(iend ,j,know)+adfac
1307 ad_zeta(iend+1,j,know)=ad_zeta(iend+1,j,know)+adfac
1308 grid(ng)%ad_h(iend ,j)=grid(ng)%ad_h(iend ,j)+adfac
1309 grid(ng)%ad_h(iend+1,j)=grid(ng)%ad_h(iend+1,j)+adfac
1310 ad_cx=0.0_r8
1311!^ tl_bry_str=tl_cff1*(FORCES(ng)%sustr(Iend+1,j)- &
1312!^ & FORCES(ng)%bustr(Iend+1,j))+ &
1313!^ & cff1*(FORCES(ng)%tl_sustr(Iend+1,j)- &
1314!^ & FORCES(ng)%tl_bustr(Iend+1,j))
1315!^
1316 adfac=cff1*ad_bry_str
1317 forces(ng)%ad_sustr(iend+1,j)= &
1318 & forces(ng)%ad_sustr(iend+1,j)+ &
1319 & adfac
1320 forces(ng)%tl_bustr(iend+1,j)= &
1321 & forces(ng)%tl_bustr(iend+1,j)- &
1322 & adfac
1323 ad_cff1=ad_cff1+(forces(ng)%sustr(iend+1,j)- &
1324 & forces(ng)%bustr(iend+1,j))*ad_bry_str
1325 ad_bry_str=0.0_r8
1326!^ tl_cff1=-cff1*cff1*0.5_r8*(GRID(ng)%tl_h(Iend ,j)+ &
1327!^ & tl_zeta(Iend ,j,know)+ &
1328!^ & GRID(ng)%tl_h(Iend+1,j)+ &
1329!^ & tl_zeta(Iend+1,j,know))
1330!^
1331 adfac=-cff1*cff1*0.5_r8*ad_cff1
1332 ad_zeta(iend ,j,know)=ad_zeta(iend ,j,know)+adfac
1333 ad_zeta(iend+1,j,know)=ad_zeta(iend+1,j,know)+adfac
1334 grid(ng)%ad_h(iend ,j)=grid(ng)%ad_h(iend ,j)+adfac
1335 grid(ng)%ad_h(iend+1,j)=grid(ng)%ad_h(iend+1,j)+adfac
1336 ad_cff1=0.0_r8
1337# ifdef UV_COR
1338!^ tl_bry_cor=0.125_r8*(tl_vbar(Iend ,j ,know)+ &
1339!^ & tl_vbar(Iend ,j+1,know)+ &
1340!^ & tl_vbar(Iend+1,j ,know)+ &
1341!^ & tl_vbar(Iend+1,j+1,know))* &
1342!^ & (GRID(ng)%f(Iend ,j)+ &
1343!^ & GRID(ng)%f(Iend+1,j))
1344!^
1345 adfac=0.125_r8*(grid(ng)%f(iend ,j)+ &
1346 & grid(ng)%f(iend+1,j))*ad_bry_cor
1347 ad_vbar(iend ,j ,know)=ad_vbar(iend ,j ,know)+adfac
1348 ad_vbar(iend+1,j ,know)=ad_vbar(iend+1,j ,know)+adfac
1349 ad_vbar(iend ,j+1,know)=ad_vbar(iend ,j+1,know)+adfac
1350 ad_vbar(iend+1,j+1,know)=ad_vbar(iend+1,j+1,know)+adfac
1351 ad_bry_cor=0.0_r8
1352# else
1353!^ tl_bry_cor=0.0_r8
1354!^
1355 ad_bry_cor=0.0_r8
1356# endif
1357 IF (ad_lbc(ieast,isfsur,ng)%acquire) THEN
1358# ifdef ADJUST_BOUNDARY
1359 IF (lobc(ieast,isubar,ng)) THEN
1360!^ tl_bry_pgr=tl_bry_pgr- &
1361!^ & g*BOUNDARY(ng)%tl_zeta_east(j)* &
1362!^ & 0.5_r8*GRID(ng)%pm(Iend,j)
1363!^
1364 boundary(ng)%ad_zeta_east(j)=boundary(ng)% &
1365 & ad_zeta_east(j)- &
1366 & g*0.5_r8* &
1367 & grid(ng)%pm(iend,j)* &
1368 & ad_bry_pgr
1369 END IF
1370# endif
1371!^ tl_bry_pgr=g*tl_zeta(Iend,j,know)* &
1372!^ & 0.5_r8*GRID(ng)%pm(Iend,j)
1373!^
1374 ad_zeta(iend,j,know)=ad_zeta(iend,j,know)+ &
1375 & g*0.5_r8*grid(ng)%pm(iend,j)* &
1376 & ad_bry_pgr
1377 ad_bry_pgr=0.0_r8
1378 ELSE
1379!^ tl_bry_pgr=-g*(tl_zeta(Iend+1,j,know)- &
1380!^ & tl_zeta(Iend ,j,know))* &
1381!^ & 0.5_r8*(GRID(ng)%pm(Iend ,j)+ &
1382!^ & GRID(ng)%pm(Iend+1,j))
1383!^
1384 adfac=-g*0.5_r8*(grid(ng)%pm(iend ,j)+ &
1385 & grid(ng)%pm(iend+1,j))*ad_bry_pgr
1386 ad_zeta(iend ,j,know)=ad_zeta(iend ,j,know)-adfac
1387 ad_zeta(iend+1,j,know)=ad_zeta(iend+1,j,know)+adfac
1388 ad_bry_pgr=0.0_r8
1389 END IF
1390# else
1391# ifdef ADJUST_BOUNDARY
1392 IF (lobc(ieast,isubar,ng)) THEN
1393!^ tl_bry_val=BOUNDARY(ng)%tl_ubar_east(j)
1394!^
1395 boundary(ng)%ad_ubar_east(j)= &
1396 & boundary(ng)%ad_ubar_east(j)+ &
1397 & ad_bry_val
1398 ad_bry_val=0.0_r8
1399 ELSE
1400!^ tl_bry_val=0.0_r8
1401!^
1402 ad_bry_val=0.0_r8
1403 END IF
1404# else
1405!^ tl_bry_val=0.0_r8
1406!^
1407 ad_bry_val=0.0_r8
1408# endif
1409# endif
1410 END IF
1411 END DO
1412!
1413! Eastern edge, clamped boundary condition.
1414!
1415 ELSE IF (ad_lbc(ieast,isubar,ng)%clamped) THEN
1416 DO j=jstr,jend
1417 IF (lbc_apply(ng)%east(j)) THEN
1418# ifdef MASKING
1419!^ tl_ubar(Iend+1,j,kout)=tl_ubar(Iend+1,j,kout)* &
1420!^ & GRID(ng)%umask(Iend+1,j)
1421!^
1422 ad_ubar(iend+1,j,kout)=ad_ubar(iend+1,j,kout)* &
1423 & grid(ng)%umask(iend+1,j)
1424# endif
1425# ifdef ADJUST_BOUNDARY
1426 IF (lobc(ieast,isubar,ng)) THEN
1427!^ tl_ubar(Iend+1,j,kout)=BOUNDARY(ng)%tl_ubar_east(j)
1428!^
1429 boundary(ng)%ad_ubar_east(j)= &
1430 & boundary(ng)%ad_ubar_east(j)+ &
1431 & ad_ubar(iend+1,j,kout)
1432 ad_ubar(iend+1,j,kout)=0.0_r8
1433 ELSE
1434!^ tl_ubar(Iend+1,j,kout)=0.0_r8
1435!^
1436 ad_ubar(iend+1,j,kout)=0.0_r8
1437 END IF
1438# else
1439!^ tl_ubar(Iend+1,j,kout)=0.0_r8
1440!^
1441 ad_ubar(iend+1,j,kout)=0.0_r8
1442# endif
1443 END IF
1444 END DO
1445!
1446! Eastern edge, gradient boundary condition.
1447!
1448 ELSE IF (ad_lbc(ieast,isubar,ng)%gradient) THEN
1449 DO j=jstr,jend
1450 IF (lbc_apply(ng)%east(j)) THEN
1451# ifdef MASKING
1452!^ tl_ubar(Iend+1,j,kout)=tl_ubar(Iend+1,j,kout)* &
1453!^ & GRID(ng)%umask(Iend+1,j)
1454!^
1455 ad_ubar(iend+1,j,kout)=ad_ubar(iend+1,j,kout)* &
1456 & grid(ng)%umask(iend+1,j)
1457# endif
1458!^ tl_ubar(Iend+1,j,kout)=tl_ubar(Iend,j,kout)
1459!^
1460 ad_ubar(iend,j,kout)=ad_ubar(iend,j,kout)+ &
1461 & ad_ubar(iend+1,j,kout)
1462 ad_ubar(iend+1,j,kout)=0.0_r8
1463 END IF
1464 END DO
1465!
1466! Eastern edge, reduced-physics boundary condition.
1467!
1468 ELSE IF (ad_lbc(ieast,isubar,ng)%reduced) THEN
1469 DO j=jstr,jend
1470 IF (lbc_apply(ng)%east(j)) THEN
1471 cff=1.0_r8/(0.5_r8*(grid(ng)%h(iend ,j)+ &
1472 & zeta(iend ,j,know)+ &
1473 & grid(ng)%h(iend+1,j)+ &
1474 & zeta(iend+1,j,know)))
1475# ifdef MASKING
1476!^ tl_ubar(Iend+1,j,kout)=tl_ubar(Iend+1,j,kout)* &
1477!^ & GRID(ng)%umask(Iend+1,j)
1478!^
1479 ad_ubar(iend+1,j,kout)=ad_ubar(iend+1,j,kout)* &
1480 & grid(ng)%umask(iend+1,j)
1481# endif
1482!^ tl_ubar(Iend+1,j,kout)=tl_ubar(Iend+1,j,know)+ &
1483!^ & dt2d*(tl_bry_pgr+ &
1484!^ & tl_bry_cor+ &
1485!^ & tl_bry_str)
1486!^
1487 adfac=dt2d*ad_ubar(iend+1,j,kout)
1488 ad_bry_pgr=ad_bry_pgr+adfac
1489 ad_bry_cor=ad_bry_cor+adfac
1490 ad_bry_str=ad_bry_str+adfac
1491 ad_ubar(iend+1,j,know)=ad_ubar(iend+1,j,know)+ &
1492 & ad_ubar(iend+1,j,kout)
1493 ad_ubar(iend+1,j,kout)=0.0_r8
1494!^ tl_bry_str=tl_cff*(FORCES(ng)%sustr(Iend+1,j)- &
1495!^ & FORCES(ng)%bustr(Iend+1,j))+ &
1496!^ & cff*(FORCES(ng)%tl_sustr(Iend+1,j)- &
1497!^ & FORCES(ng)%tl_bustr(Iend+1,j))
1498!^
1499 adfac=cff*ad_bry_str
1500 forces(ng)%ad_sustr(iend+1,j)= &
1501 & forces(ng)%ad_sustr(iend+1,j)+ &
1502 & adfac
1503 forces(ng)%ad_bustr(iend+1,j)= &
1504 & forces(ng)%ad_bustr(iend+1,j)- &
1505 & adfac
1506 ad_cff=ad_cff+(forces(ng)%sustr(iend+1,j)- &
1507 & forces(ng)%bustr(iend+1,j))*ad_bry_str
1508 ad_bry_str=0.0_r8
1509# ifdef UV_COR
1510!^ tl_bry_cor=0.125_r8*(tl_vbar(Iend, j ,know)+ &
1511!^ & tl_vbar(Iend ,j+1,know)+ &
1512!^ & tl_vbar(Iend+1,j ,know)+ &
1513!^ & tl_vbar(Iend+1,j+1,know))* &
1514!^ & (GRID(ng)%f(Iend ,j)+ &
1515!^ & GRID(ng)%f(Iend+1,j))
1516!^
1517 adfac=0.125_r8*(grid(ng)%f(iend,j)+ &
1518 & grid(ng)%f(iend+1,j))*ad_bry_cor
1519 ad_vbar(iend, j ,know)=ad_vbar(iend ,j ,know)+adfac
1520 ad_vbar(iend+1,j ,know)=ad_vbar(iend+1,j ,know)+adfac
1521 ad_vbar(iend ,j+1,know)=ad_vbar(iend ,j+1,know)+adfac
1522 ad_vbar(iend+1,j+1,know)=ad_vbar(iend+1,j+1,know)+adfac
1523 ad_bry_cor=0.0_r8
1524# else
1525!^ tl_bry_cor=0.0_r8
1526!^
1527 ad_bry_cor=0.0_r8
1528# endif
1529 IF (ad_lbc(ieast,isfsur,ng)%acquire) THEN
1530# ifdef ADJUST_BOUNDARY
1531 IF (lobc(ieast,isubar,ng)) THEN
1532!^ tl_bry_pgr=tl_bry_pgr- &
1533!^ & g*BOUNDARY(ng)%tl_zeta_east(j)* &
1534!^ & 0.5_r8*GRID(ng)%pm(Iend,j)
1535!^
1536 boundary(ng)%ad_zeta_east(j)=boundary(ng)% &
1537 & ad_zeta_east(j)- &
1538 & g*0.5_r8* &
1539 & grid(ng)%pm(iend,j)* &
1540 & ad_bry_pgr
1541 END IF
1542# endif
1543!^ tl_bry_pgr=g*tl_zeta(Iend,j,know)* &
1544!^ & 0.5_r8*GRID(ng)%pm(Iend,j)
1545!^
1546 ad_zeta(iend,j,know)=ad_zeta(iend,j,know)+ &
1547 & g*0.5_r8*grid(ng)%pm(iend,j)* &
1548 & ad_bry_pgr
1549 ad_bry_pgr=0.0_r8
1550 ELSE
1551!^ tl_bry_pgr=-g*(tl_zeta(Iend+1,j,know)- &
1552!^ & tl_zeta(Iend ,j,know))* &
1553!^ & 0.5_r8*(GRID(ng)%pm(Iend ,j)+ &
1554!^ & GRID(ng)%pm(Iend+1,j))
1555!^
1556 adfac=-g*0.5_r8*(grid(ng)%pm(iend ,j)+ &
1557 & grid(ng)%pm(iend+1,j))*ad_bry_pgr
1558 ad_zeta(iend ,j,know)=ad_zeta(iend ,j,know)-adfac
1559 ad_zeta(iend+1,j,know)=ad_zeta(iend+1,j,know)+adfac
1560 ad_bry_pgr=0.0_r8
1561 END IF
1562 END IF
1563 END DO
1564!
1565! Eastern edge, closed boundary condition.
1566!
1567 ELSE IF (ad_lbc(ieast,isubar,ng)%closed) THEN
1568 DO j=jstr,jend
1569 IF (lbc_apply(ng)%east(j)) THEN
1570!^ tl_ubar(Iend+1,j,kout)=0.0_r8
1571!^
1572 ad_ubar(iend+1,j,kout)=0.0_r8
1573 END IF
1574 END DO
1575 END IF
1576 END IF
1577!
1578!-----------------------------------------------------------------------
1579! Lateral boundary conditions at the western edge.
1580!-----------------------------------------------------------------------
1581!
1582 IF (domain(ng)%Western_Edge(tile)) THEN
1583!
1584! Western edge, implicit upstream radiation condition.
1585!
1586 IF (ad_lbc(iwest,isubar,ng)%radiation) THEN
1587 IF (iic(ng).ne.0) THEN
1588 DO j=jstr,jend
1589 IF (lbc_apply(ng)%west(j)) THEN
1590# if defined CELERITY_READ && defined FORWARD_READ
1591 IF (ad_lbc(iwest,isubar,ng)%nudging) THEN
1592 IF (lnudgem2clm(ng)) THEN
1593 obc_out=0.5_r8* &
1594 & (clima(ng)%M2nudgcof(istr-1,j)+ &
1595 & clima(ng)%M2nudgcof(istr ,j))
1596 obc_in =obcfac(ng)*obc_out
1597 ELSE
1598 obc_out=m2obc_out(ng,iwest)
1599 obc_in =m2obc_in(ng,iwest)
1600 END IF
1601 IF (boundary(ng)%ubar_west_Cx(j).lt.0.0_r8) THEN
1602 tau=obc_in
1603 ELSE
1604 tau=obc_out
1605 END IF
1606 tau=tau*dt2d
1607 END IF
1608 cx=boundary(ng)%ubar_west_Cx(j)
1609# ifdef RADIATION_2D
1610 ce=boundary(ng)%ubar_west_Ce(j)
1611# else
1612 ce=0.0_r8
1613# endif
1614 cff=boundary(ng)%ubar_west_C2(j)
1615# endif
1616# ifdef MASKING
1617!^ tl_ubar(Istr,j,kout)=tl_ubar(Istr,j,kout)* &
1618!^ & GRID(ng)%umask(Istr,j)
1619!^
1620 ad_ubar(istr,j,kout)=ad_ubar(istr,j,kout)* &
1621 & grid(ng)%umask(istr,j)
1622# endif
1623 IF (ad_lbc(iwest,isubar,ng)%nudging) THEN
1624!^ tl_ubar(Istr,j,kout)=tl_ubar(Istr,j,kout)- &
1625!^ & tau*tl_ubar(Istr,j,know)
1626!^
1627 ad_ubar(istr,j,know)=ad_ubar(istr,j,know)- &
1628 & tau*ad_ubar(istr,j,kout)
1629 END IF
1630!^ tl_ubar(Istr,j,kout)=(cff*tl_ubar(Istr ,j,know)+ &
1631!^ & Cx *tl_ubar(Istr+1,j,kout)- &
1632!^ & MAX(Ce,0.0_r8)* &
1633!^ & tl_grad(Istr,j )- &
1634!^ & MIN(Ce,0.0_r8)* &
1635!^ & tl_grad(Istr,j+1))/ &
1636!^ & (cff+Cx)
1637!^
1638 adfac=ad_ubar(istr,j,kout)/(cff+cx)
1639 ad_grad(istr,j )=ad_grad(istr,j )-max(ce,0.0_r8)*adfac
1640 ad_grad(istr,j+1)=ad_grad(istr,j+1)-min(ce,0.0_r8)*adfac
1641 ad_ubar(istr ,j,know)=ad_ubar(istr ,j,know)+cff*adfac
1642 ad_ubar(istr+1,j,kout)=ad_ubar(istr+1,j,kout)+cx *adfac
1643 ad_ubar(istr ,j,kout)=0.0_r8
1644 END IF
1645 END DO
1646 END IF
1647!
1648! Western edge, Flather boundary condition.
1649!
1650 ELSE IF (ad_lbc(iwest,isubar,ng)%Flather) THEN
1651 DO j=jstr,jend
1652 IF (lbc_apply(ng)%west(j)) THEN
1653 cff=1.0_r8/(0.5_r8*(grid(ng)%h(istr-1,j)+ &
1654 & zeta(istr-1,j,know)+ &
1655 & grid(ng)%h(istr ,j)+ &
1656 & zeta(istr ,j,know)))
1657 cx=sqrt(g*cff)
1658# ifdef MASKING
1659!^ tl_ubar(Istr,j,kout)=tl_ubar(Istr,j,kout)* &
1660!^ & GRID(ng)%umask(Istr,j)
1661!^
1662 ad_ubar(istr,j,kout)=ad_ubar(istr,j,kout)* &
1663 & grid(ng)%umask(istr,j)
1664# endif
1665# ifdef ADJUST_BOUNDARY
1666 IF (lobc(iwest,isubar,ng)) THEN
1667!^ tl_ubar(Istr,j,kout)=tl_ubar(Istr,j,kout)+ &
1668!^ & Cx*BOUNDARY(ng)%tl_zeta_west(j)
1669!^
1670 boundary(ng)%ad_zeta_west(j)= &
1671 & boundary(ng)%ad_zeta_west(j)+ &
1672 & cx*ad_ubar(istr,j,kout)
1673 END IF
1674# endif
1675# if defined ATM_PRESS && defined PRESS_COMPENSATE
1676!^ tl_ubar(Istr,j,kout)=tl_bry_val- &
1677!^ & tl_Cx* &
1678!^ & (0.5_r8* &
1679!^ & (zeta(Istr-1,j,know)+ &
1680!^ & zeta(Istr ,j,know)+ &
1681!^ & fac*(FORCES(ng)%Pair(Istr-1,j)+ &
1682!^ & FORCES(ng)%Pair(Istr ,j)- &
1683!^ & 2.0_r8*OneAtm))- &
1684!^ & BOUNDARY(ng)%zeta_west(j))- &
1685!^ & Cx* &
1686!^ & (0.5_r8* &
1687!^ & (tl_zeta(Istr-1,j,know)+ &
1688!^ & tl_zeta(Istr ,j,know)))
1689!^
1690 adfac=cx*0.5_r8*ad_ubar(istr,j,kout)
1691 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)-adfac
1692 ad_zeta(istr ,j,know)=ad_zeta(istr ,j,know)-adfac
1693 ad_cx=ad_cx- &
1694 & (0.5_r8*(zeta(istr-1,j,know)+ &
1695 & zeta(istr ,j,know)+ &
1696 & fac*(forces(ng)%Pair(istr-1,j)+ &
1697 & forces(ng)%Pair(istr ,j)- &
1698 & 2.0_r8*oneatm))- &
1699 & boundary(ng)%zeta_west(j))*ad_ubar(istr,j,kout)
1700 ad_bry_val=ad_bry_val+ad_ubar(istr,j,kout)
1701 ad_ubar(istr,j,kout)=0.0_r8
1702# else
1703!^ tl_ubar(Istr,j,kout)=tl_bry_val- &
1704!^ & tl_Cx* &
1705!^ & (0.5_r8*(zeta(Istr-1,j,know)+ &
1706!^ & zeta(Istr ,j,know))- &
1707!^ & BOUNDARY(ng)%zeta_west(j))- &
1708!^ & Cx* &
1709!^ & (0.5_r8*(tl_zeta(Istr-1,j,know)+ &
1710!^ & tl_zeta(Istr ,j,know)))
1711!^
1712 adfac=cx*0.5_r8*ad_ubar(istr,j,kout)
1713 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)-adfac
1714 ad_zeta(istr ,j,know)=ad_zeta(istr ,j,know)-adfac
1715 ad_cx=ad_cx- &
1716 & (0.5_r8*(zeta(istr-1,j,know)+ &
1717 & zeta(istr ,j,know))- &
1718 & boundary(ng)%zeta_west(j))*ad_ubar(istr,j,kout)
1719 ad_bry_val=ad_bry_val+ad_ubar(istr,j,kout)
1720 ad_ubar(istr,j,kout)=0.0_r8
1721# endif
1722!^ tl_Cx=0.5_r8*g*tl_cff/Cx
1723!^
1724 ad_cff=ad_cff+0.5_r8*g*ad_cx/cx
1725 ad_cx=0.0_r8
1726!^ tl_cff=-cff*cff*(0.5_r8*(GRID(ng)%tl_h(Istr-1,j)+ &
1727!^ & tl_zeta(Istr-1,j,know)+ &
1728!^ & GRID(ng)%tl_h(Istr ,j)+ &
1729!^ & tl_zeta(Istr ,j,know)))
1730!^
1731 adfac=-cff*cff*0.5_r8*ad_cff
1732 grid(ng)%ad_h(istr-1,j)=grid(ng)%ad_h(istr-1,j)+adfac
1733 grid(ng)%ad_h(istr ,j)=grid(ng)%ad_h(istr ,j)+adfac
1734 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)+adfac
1735 ad_zeta(istr ,j,know)=ad_zeta(istr ,j,know)+adfac
1736 ad_cff=0.0_r8
1737
1738# if defined SSH_TIDES_NOT_YET && !defined UV_TIDES_NOT_YET
1739 IF (ad_lbc(iwest,isfsur,ng)%acquire) THEN
1740 bry_pgr=-g*(zeta(istr,j,know)- &
1741 & boundary(ng)%zeta_west(j))* &
1742 & 0.5_r8*grid(ng)%pm(istr,j)
1743 ELSE
1744 bry_pgr=-g*(zeta(istr ,j,know)- &
1745 & zeta(istr-1,j,know))* &
1746 & 0.5_r8*(grid(ng)%pm(istr-1,j)+ &
1747 & grid(ng)%pm(istr,j))
1748 END IF
1749# ifdef UV_COR
1750 bry_cor=0.125_r8*(vbar(istr-1,j ,know)+ &
1751 & vbar(istr-1,j+1,know)+ &
1752 & vbar(istr ,j ,know)+ &
1753 & vbar(istr ,j+1,know))* &
1754 & (grid(ng)%f(istr-1,j)+ &
1755 & grid(ng)%f(istr,j))
1756# else
1757 bry_cor=0.0_r8
1758# endif
1759 cff1=1.0_r8/(0.5_r8*(grid(ng)%h(istr-1,j)+ &
1760 & zeta(istr-1,j,know)+ &
1761 & grid(ng)%h(istr ,j)+ &
1762 & zeta(istr ,j,know)))
1763 bry_str=cff1*(forces(ng)%sustr(istr,j)- &
1764 & forces(ng)%bustr(istr,j))
1765 cx=1.0_r8/sqrt(g*0.5_r8*(grid(ng)%h(istr-1,j)+ &
1766 & zeta(istr-1,j,know)+ &
1767 & grid(ng)%h(istr ,j)+ &
1768 & zeta(istr ,j,know)))
1769 cff2=grid(ng)%om_u(istr,j)*cx
1770!^ tl_bry_val=tl_ubar(Istr+1,j,know)+ &
1771!^ & tl_cff2*(bry_pgr+ &
1772!^ & bry_cor+ &
1773!^ & bry_str)+ &
1774!^ & cff2*(tl_bry_pgr+ &
1775!^ & tl_bry_cor+ &
1776!^ & tl_bry_str)
1777!^
1778 adfac=cff2*ad_bry_val
1779 ad_bry_pgr=ad_bry_pgr+adfac
1780 ad_bry_cor=ad_bry_cor+adfac
1781 ad_bry_str=ad_bry_str+adfac
1782 ad_cff2=ad_cff2+(bry_pgr+ &
1783 & bry_cor+ &
1784 & bry_str)*ad_bry_val
1785 ad_ubar(istr+1,j,know)=ad_ubar(istr+1,j,know)+ad_bry_val
1786 ad_bry_val=0.0_r8
1787!^ tl_cff2=GRID(ng)%om_u(Istr,j)*tl_Cx
1788!^
1789 ad_cx=ad_cx+grid(ng)%om_u(istr,j)*ad_cff2
1790 ad_cff2=0.0_r8
1791!^ tl_Cx=-Cx*Cx*Cx*0.25_r8*g*(GRID(ng)%tl_h(Istr-1,j)+ &
1792!^ & tl_zeta(Istr-1,j,know)+ &
1793!^ & GRID(ng)%tl_h(Istr ,j)+ &
1794!^ & tl_zeta(Istr ,j,know))
1795!^
1796 adfac=-cx*cx*cx*0.25_r8*g*ad_cx
1797 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)+adfac
1798 ad_zeta(istr ,j,know)=ad_zeta(istr ,j,know)+adfac
1799 grid(ng)%ad_h(istr-1,j)=grid(ng)%ad_h(istr-1,j)+adfac
1800 grid(ng)%ad_h(istr ,j)=grid(ng)%ad_h(istr ,j)+adfac
1801 ad_cx=0.0_r8
1802!^ tl_bry_str=tl_cff1*(FORCES(ng)%sustr(Istr,j)- &
1803!^ & FORCES(ng)%bustr(Istr,j))+ &
1804!^ & cff1*(FORCES(ng)%tl_sustr(Istr,j)- &
1805!^ & FORCES(ng)%tl_bustr(Istr,j))
1806!^
1807 adfac=cff1*ad_bry_str
1808 forces(ng)%ad_sustr(istr,j)=forces(ng)%ad_sustr(istr,j)+ &
1809 & adfac
1810 forces(ng)%ad_bustr(istr,j)=forces(ng)%ad_bustr(istr,j)- &
1811 & adfac
1812 ad_cff1=ad_cff1+(forces(ng)%sustr(istr,j)- &
1813 & forces(ng)%bustr(istr,j))*ad_bry_str
1814 ad_bry_str=0.0_r8
1815!^ tl_cff1=-cff1*cff1*0.5_r8*(GRID(ng)%tl_h(Istr-1,j)+ &
1816!^ & tl_zeta(Istr-1,j,know)+ &
1817!^ & GRID(ng)%tl_h(Istr ,j)+ &
1818!^ & tl_zeta(Istr ,j,know))
1819!^
1820 adfac=-cff1*cff1*0.5_r8*ad_cff1
1821 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)+adfac
1822 ad_zeta(istr ,j,know)=ad_zeta(istr ,j,know)+adfac
1823 grid(ng)%ad_h(istr-1,j)=grid(ng)%ad_h(istr-1,j)+adfac
1824 grid(ng)%ad_h(istr ,j)=grid(ng)%ad_h(istr ,j)+adfac
1825 ad_cff1=0.0_r8
1826# ifdef UV_COR
1827!^ tl_bry_cor=0.125_r8*(tl_vbar(Istr-1,j ,know)+ &
1828!^ & tl_vbar(Istr-1,j+1,know)+ &
1829!^ & tl_vbar(Istr ,j ,know)+ &
1830!^ & tl_vbar(Istr ,j+1,know))* &
1831!^ & (GRID(ng)%f(Istr-1,j)+ &
1832!^ & GRID(ng)%f(Istr ,j))
1833!^
1834 adfac=0.125_r8*(grid(ng)%f(istr-1,j)+ &
1835 & grid(ng)%f(istr ,j))*ad_bry_cor
1836 ad_vbar(istr-1,j ,know)=ad_vbar(istr-1,j ,know)+adfac
1837 ad_vbar(istr ,j ,know)=ad_vbar(istr ,j ,know)+adfac
1838 ad_vbar(istr-1,j+1,know)=ad_vbar(istr-1,j+1,know)+adfac
1839 ad_vbar(istr ,j+1,know)=ad_vbar(istr ,j+1,know)+adfac
1840 ad_bry_cor=0.0_r8
1841# else
1842!^ tl_bry_cor=0.0_r8
1843!^
1844 ad_bry_cor=0.0_r8
1845# endif
1846 IF (ad_lbc(iwest,isfsur,ng)%acquire) THEN
1847# ifdef ADJUST_BOUNDARY
1848 IF (lobc(iwest,isubar,ng)) THEN
1849!^ tl_bry_pgr=tl_bry_pgr+ &
1850!^ & g*BOUNDARY(ng)%tl_zeta_west(j)* &
1851!^ & 0.5_r8*GRID(ng)%pm(Istr,j)
1852!^
1853 boundary(ng)%ad_zeta_west(j)=boundary(ng)% &
1854 & ad_zeta_west(j)+ &
1855 & g*0.5_r8* &
1856 & grid(ng)%pm(istr,j)* &
1857 & ad_bry_pgr
1858 END IF
1859# endif
1860!^ tl_bry_pgr=-g*tl_zeta(Istr,j,know)* &
1861!^ & 0.5_r8*GRID(ng)%pm(Istr,j)
1862!^
1863 ad_zeta(istr,j,know)=ad_zeta(istr,j,know)- &
1864 & g*0.5_r8*grid(ng)%pm(istr,j)* &
1865 & ad_bry_pgr
1866 ad_bry_pgr=0.0_r8
1867 ELSE
1868!^ tl_bry_pgr=-g*(tl_zeta(Istr ,j,know)- &
1869!^ & tl_zeta(Istr-1,j,know))* &
1870!^ & 0.5_r8*(GRID(ng)%pm(Istr-1,j)+ &
1871!^ & GRID(ng)%pm(Istr ,j))
1872!^
1873 adfac=-g*0.5_r8*(grid(ng)%pm(istr-1,j)+ &
1874 & grid(ng)%pm(istr ,j))*ad_bry_pgr
1875 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)-adfac
1876 ad_zeta(istr ,j,know)=ad_zeta(istr ,j,know)+adfac
1877 ad_bry_pgr=0.0_r8
1878 END IF
1879# else
1880# ifdef ADJUST_BOUNDARY
1881 IF (lobc(iwest,isubar,ng)) THEN
1882!^ tl_bry_val=BOUNDARY(ng)%tl_ubar_west(j)
1883!^
1884 boundary(ng)%ad_ubar_west(j)= &
1885 & boundary(ng)%ad_ubar_west(j)+ &
1886 & ad_bry_val
1887 ad_bry_val=0.0_r8
1888 ELSE
1889!^ tl_bry_val=0.0_r8
1890!^
1891 ad_bry_val=0.0_r8
1892 END IF
1893# else
1894!^ tl_bry_val=0.0_r8
1895!^
1896 ad_bry_val=0.0_r8
1897# endif
1898# endif
1899 END IF
1900 END DO
1901!
1902! Western edge, Shchepetkin boundary condition (Maison et al., 2010).
1903!
1904 ELSE IF (ad_lbc(iwest,isubar,ng)%Shchepetkin) THEN
1905 DO j=jstr,jend
1906 IF (lbc_apply(ng)%west(j)) THEN
1907# ifdef WET_DRY_NOT_YET
1908 cff=0.5_r8*(grid(ng)%h(istr-1,j)+ &
1909 & zeta(istr-1,j,know)+ &
1910 & grid(ng)%h(istr ,j)+ &
1911 & zeta(istr ,j,know))
1912# else
1913 cff=0.5_r8*(grid(ng)%h(istr-1,j)+ &
1914 & grid(ng)%h(istr ,j))
1915# endif
1916 cff1=sqrt(g/cff)
1917 cx=dt2d*cff1*cff*0.5_r8*(grid(ng)%pm(istr-1,j)+ &
1918 & grid(ng)%pm(istr ,j))
1919 zx=(0.5_r8+cx)*zeta(istr ,j,know)+ &
1920 & (0.5_r8-cx)*zeta(istr-1,j,know)
1921 IF (cx.gt.co) THEN
1922 cff2=(1.0_r8-co/cx)**2
1923 cff3=zeta(istr,j,kout)+ &
1924 & cx*zeta(istr-1,j,know)- &
1925 & (1.0_r8+cx)*zeta(istr,j,know)
1926 zx=zx+cff2*cff3
1927 END IF
1928# ifdef MASKING
1929!^ tl_ubar(Istr,j,kout)=tl_ubar(Istr,j,kout)* &
1930!^ & GRID(ng)%umask(Istr,j)
1931!^
1932 ad_ubar(istr,j,kout)=ad_ubar(istr,j,kout)* &
1933 & grid(ng)%umask(istr,j)
1934# endif
1935# ifdef ADJUST_BOUNDARY
1936 IF (lobc(iwest,isubar,ng)) THEN
1937!^ tl_ubar(Istr,j,kout)=tl_ubar(Istr,j,kout)+ &
1938!^ & Cx*BOUNDARY(ng)%tl_zeta_west(j)
1939!^
1940 boundary(ng)%ad_zeta_west(j)= &
1941 & boundary(ng)%ad_zeta_west(j)+ &
1942 & cx*ad_ubar(istr,j,kout)
1943 END IF
1944# endif
1945!^ tl_ubar(Istr,j,kout)=0.5_r8* &
1946!^ & ((1.0_r8-Cx)* &
1947!^ & tl_ubar(Istr,j,know)- &
1948!^ & tl_Cx*(ubar(Istr ,j,know)- &
1949!^ & ubar(Istr+1,j,know))+ &
1950!^ & Cx*tl_ubar(Istr+1,j,know)+ &
1951!^ & tl_bry_val- &
1952!^ & tl_cff1* &
1953!^ & (Zx-BOUNDARY(ng)%zeta_west(j))- &
1954!^ & cff1*tl_Zx)
1955!^
1956 adfac=0.5_r8*ad_ubar(istr,j,kout)
1957 ad_ubar(istr, j,know)=ad_ubar(istr ,j,know)+ &
1958 & (1.0_r8-cx)*adfac
1959 ad_ubar(istr+1,j,know)=ad_ubar(istr+1,j,know)+ &
1960 & cx*adfac
1961 ad_cx=ad_cx- &
1962 (ubar(istr ,j,know)- &
1963 & ubar(istr+1,j,know))*adfac
1964 ad_bry_val=ad_bry_val+adfac
1965 ad_cff1=ad_cff1- &
1966 & (zx-boundary(ng)%zeta_west(j))*adfac
1967 ad_zx=ad_zx-cff1*adfac
1968 ad_ubar(istr,j,kout)=0.0_r8
1969 IF (cx.gt.co) THEN
1970!^ tl_Zx=tl_Zx+cff2*tl_cff3+ &
1971!^ & tl_cff2*cff3
1972!^
1973 ad_cff2=ad_cff2+cff3*ad_zx
1974 ad_cff3=ad_cff3+cff2*ad_zx
1975!^ tl_cff3=tl_zeta(Istr,j,kout)+ &
1976!^ & Cx*tl_zeta(Istr-1,j,know)+ &
1977!^ & tl_Cx*(zeta(Istr-1,j,know)+ &
1978!^ & zeta(Istr ,j,know))- &
1979!^ & (1.0_r8+Cx)*tl_zeta(Istr,j,know)
1980!^
1981 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)+ &
1982 & cx*ad_cff3
1983 ad_zeta(istr ,j,kout)=ad_zeta(istr ,j,kout)+ &
1984 & ad_cff3
1985 ad_zeta(istr ,j,know)=ad_zeta(istr ,j,know)- &
1986 & (1.0_r8+cx)*ad_cff3
1987 ad_cx=ad_cx+ &
1988 & (zeta(istr-1,j,know)+ &
1989 & zeta(istr ,j,know))*ad_cff3
1990 ad_cff3=0.0_r8
1991!^ tl_cff2=2.0_r8*cff2*Co*tl_Cx/(Cx*Cx)
1992!^
1993 ad_cx=ad_cx+ &
1994 & 2.0_r8*cff2*co*ad_cff2/(cx*cx)
1995 ad_cff2=0.0_r8
1996 END IF
1997!^ tl_Zx=(0.5_r8+Cx)*tl_zeta(Istr ,j,know)+ &
1998!^ & (0.5_r8-Cx)*tl_zeta(Istr-1,j,know)+ &
1999!^ & tl_Cx*(zeta(Istr ,j,know)- &
2000!^ & zeta(Istr-1,j,know))
2001!^
2002 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)+ &
2003 & (0.5_r8-cx)*ad_zx
2004 ad_zeta(istr ,j,know)=ad_zeta(istr ,j,know)+ &
2005 & (0.5_r8+cx)*ad_zx
2006 ad_cx=ad_cx+ &
2007 & (zeta(istr ,j,know)- &
2008 & zeta(istr-1,j,know))*ad_zx
2009 ad_zx=0.0_r8
2010!^ tl_Cx=dt2d*0.5_r8*(GRID(ng)%pm(Istr-1,j)+ &
2011!^ & GRID(ng)%pm(Istr ,j))* &
2012!^ & (cff1*tl_cff+ &
2013!^ & tl_cff1*cff)
2014!^
2015 adfac=dt2d*0.5_r8*(grid(ng)%pm(istr-1,j)+ &
2016 & grid(ng)%pm(istr ,j))*ad_cx
2017 ad_cff=ad_cff+cff1*adfac
2018 ad_cff1=ad_cff1+cff*adfac
2019 ad_cx=0.0_r8
2020!^ tl_cff1=-0.5_r8*cff1*tl_cff/cff
2021!^
2022 ad_cff=ad_cff- &
2023 & 0.5_r8*cff1*ad_cff1/cff
2024 ad_cff1=0.0_r8
2025# ifdef WET_DRY_NOT_YET
2026!^ tl_cff=0.5_r8*(GRID(ng)%tl_h(Istr-1,j)+ &
2027!^ & tl_zeta(Istr-1,j,know)+ &
2028!^ & GRID(ng)%tl_h(Istr ,j)+ &
2029!^ & tl_zeta(Istr ,j,know))
2030!^
2031 adfac=0.5_r8*ad_cff
2032 grid(ng)%ad_h(istr-1,j)=grid(ng)%ad_h(istr-1,j)+adfac
2033 grid(ng)%ad_h(istr ,j)=grid(ng)%ad_h(istr ,j)+adfac
2034 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)+adfac
2035 ad_zeta(istr ,j,know)=ad_zeta(istr ,j,know)+adfac
2036 ad_cff=0.0_r8
2037# else
2038!^ tl_cff=0.5_r8*(GRID(ng)%tl_h(Istr-1,j)+ &
2039!^ & GRID(ng)%tl_h(Istr ,j))
2040!^
2041 adfac=0.5_r8*ad_cff
2042 grid(ng)%ad_h(istr-1,j)=grid(ng)%ad_h(istr-1,j)+adfac
2043 grid(ng)%ad_h(istr ,j)=grid(ng)%ad_h(istr ,j)+adfac
2044 ad_cff=0.0_r8
2045# endif
2046# if defined SSH_TIDES_NOT_YET && !defined UV_TIDES_NOT_YET
2047 IF (ad_lbc(iwest,isfsur,ng)%acquire) THEN
2048 bry_pgr=-g*(zeta(istr,j,know)- &
2049 & boundary(ng)%zeta_west(j))* &
2050 & 0.5_r8*grid(ng)%pm(istr,j)
2051 ELSE
2052 bry_pgr=-g*(zeta(istr ,j,know)- &
2053 & zeta(istr-1,j,know))* &
2054 & 0.5_r8*(grid(ng)%pm(istr-1,j)+ &
2055 & grid(ng)%pm(istr ,j))
2056 END IF
2057# ifdef UV_COR
2058 bry_cor=0.125_r8*(vbar(istr-1,j ,know)+ &
2059 & vbar(istr-1,j+1,know)+ &
2060 & vbar(istr ,j ,know)+ &
2061 & vbar(istr ,j+1,know))* &
2062 & (grid(ng)%f(istr-1,j)+ &
2063 & grid(ng)%f(istr ,j))
2064# else
2065 bry_cor=0.0_r8
2066# endif
2067 cff1=1.0_r8/(0.5_r8*(grid(ng)%h(istr-1,j)+ &
2068 & zeta(istr-1,j,know)+ &
2069 & grid(ng)%h(istr ,j)+ &
2070 & zeta(istr ,j,know)))
2071 bry_str=cff1*(forces(ng)%sustr(istr,j)- &
2072 & forces(ng)%bustr(istr,j))
2073 cx=1.0_r8/sqrt(g*0.5_r8*(grid(ng)%h(istr-1,j)+ &
2074 & zeta(istr-1,j,know)+ &
2075 & grid(ng)%h(istr ,j)+ &
2076 & zeta(istr ,j,know)))
2077 cff2=grid(ng)%om_u(istr,j)*cx
2078!^ tl_bry_val=tl_ubar(Istr+1,j,know)+ &
2079!^ & tl_cff2*(bry_pgr+ &
2080!^ & bry_cor+ &
2081!^ & bry_str)+ &
2082!^ & cff2*(tl_bry_pgr+ &
2083!^ & tl_bry_cor+ &
2084!^ & tl_bry_str)
2085!^
2086 adfac=cff2*ad_bry_val
2087 ad_bry_pgr=ad_bry_pgr+adfac
2088 ad_bry_cor=ad_bry_cor+adfac
2089 ad_bry_str=ad_bry_str+adfac
2090 ad_cff2=ad_cff2+ &
2091 & (bry_pgr+ &
2092 & bry_cor+ &
2093 & bry_str)*ad_bry_val
2094 ad_ubar(istr+1,j,know)=ad_ubar(istr+1,j,know)+ &
2095 & ad_bry_val
2096 ad_bry_val=0.0_r8
2097!^ tl_cff2=GRID(ng)%om_u(Istr,j)*tl_Cx
2098!^
2099 ad_cx=ad_cx+grid(ng)%om_u(istr,j)*ad_cff2
2100 ad_cff2=0.0_r8
2101!^ tl_Cx=-Cx*Cx*Cx*0.25_r8*g*(GRID(ng)%tl_h(Istr-1,j)+ &
2102!^ & tl_zeta(Istr-1,j,know)+ &
2103!^ & GRID(ng)%tl_h(Istr ,j)+ &
2104!^ & tl_zeta(Istr ,j,know))
2105!^
2106 adfac=-cx*cx*cx*0.25_r8*g*ad_cx
2107 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)+adfac
2108 ad_zeta(istr ,j,know)=ad_zeta(istr ,j,know)+adfac
2109 grid(ng)%ad_h(istr-1,j)=grid(ng)%ad_h(istr-1,j)+adfac
2110 grid(ng)%ad_h(istr ,j)=grid(ng)%ad_h(istr ,j)+adfac
2111 ad_cx=0.0_r8
2112!^ tl_bry_str=tl_cff1*(FORCES(ng)%sustr(Istr,j)- &
2113!^ & FORCES(ng)%bustr(Istr,j))+ &
2114!^ & cff1*(FORCES(ng)%tl_sustr(Istr,j)- &
2115!^ & FORCES(ng)%tl_bustr(Istr,j))
2116!^
2117 adfac=cff1*ad_bry_str
2118 forces(ng)%ad_sustr(istr,j)=forces(ng)%ad_sustr(istr,j)+ &
2119 & adfac
2120 forces(ng)%ad_bustr(istr,j)=forces(ng)%ad_bustr(istr,j)- &
2121 & adfac
2122 ad_cff1=ad_cff1+(forces(ng)%sustr(istr,j)- &
2123 & forces(ng)%bustr(istr,j))*ad_bry_str
2124 ad_bry_str=0.0_r8
2125!^ tl_cff1=-cff1*cff1*0.5_r8*(GRID(ng)%tl_h(Istr-1,j)+ &
2126!^ & tl_zeta(Istr-1,j,know)+ &
2127!^ & GRID(ng)%tl_h(Istr ,j)+ &
2128!^ & tl_zeta(Istr ,j,know))
2129!^
2130 adfac=-cff1*cff1*0.5_r8*ad_cff1
2131 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)+adfac
2132 ad_zeta(istr ,j,know)=ad_zeta(istr ,j,know)+adfac
2133 grid(ng)%ad_h(istr-1,j)=grid(ng)%ad_h(istr-1,j)+adfac
2134 grid(ng)%ad_h(istr ,j)=grid(ng)%ad_h(istr ,j)+adfac
2135 ad_cff1=0.0_r8
2136# ifdef UV_COR
2137!^ tl_bry_cor=0.125_r8*(tl_vbar(Istr-1,j ,know)+ &
2138!^ & tl_vbar(Istr-1,j+1,know)+ &
2139!^ & tl_vbar(Istr ,j ,know)+ &
2140!^ & tl_vbar(Istr ,j+1,know))* &
2141!^ & (GRID(ng)%f(Istr-1,j)+ &
2142!^ & GRID(ng)%f(Istr ,j))
2143!^
2144 adfac=0.125_r8*(grid(ng)%f(istr-1,j)+ &
2145 & grid(ng)%f(istr ,j))*ad_bry_cor
2146 ad_vbar(istr-1,j ,know)=ad_vbar(istr-1,j ,know)+adfac
2147 ad_vbar(istr ,j ,know)=ad_vbar(istr ,j ,know)+adfac
2148 ad_vbar(istr-1,j+1,know)=ad_vbar(istr-1,j+1,know)+adfac
2149 ad_vbar(istr ,j+1,know)=ad_vbar(istr ,j+1,know)+adfac
2150 ad_bry_cor=0.0_r8
2151# else
2152!^ tl_bry_cor=0.0_r8
2153!^
2154 ad_bry_cor=0.0_r8
2155# endif
2156 IF (ad_lbc(iwest,isfsur,ng)%acquire) THEN
2157# ifdef ADJUST_BOUNDARY
2158 IF (lobc(iwest,isubar,ng)) THEN
2159!^ tl_bry_pgr=tl_bry_pgr+ &
2160!^ & g*BOUNDARY(ng)%tl_zeta_west(j)* &
2161!^ & 0.5_r8*GRID(ng)%pm(Istr,j)
2162!^
2163 boundary(ng)%ad_zeta_west(j)=boundary(ng)% &
2164 & ad_zeta_west(j)+ &
2165 & g*0.5_r8* &
2166 & grid(ng)%pm(istr,j)* &
2167 & ad_bry_pgr
2168 END IF
2169# endif
2170!^ tl_bry_pgr=-g*tl_zeta(Istr,j,know)* &
2171!^ & 0.5_r8*GRID(ng)%pm(Istr,j)
2172!^
2173 ad_zeta(istr,j,know)=ad_zeta(istr,j,know)- &
2174 & g*0.5_r8*grid(ng)%pm(istr,j)* &
2175 & ad_bry_pgr
2176 ad_bry_pgr=0.0_r8
2177 ELSE
2178!^ tl_bry_pgr=-g*(tl_zeta(Istr ,j,know)- &
2179!^ & tl_zeta(Istr-1,j,know))* &
2180!^ & 0.5_r8*(GRID(ng)%pm(Istr-1,j)+ &
2181!^ & GRID(ng)%pm(Istr ,j))
2182!^
2183 adfac=-g*0.5_r8*(grid(ng)%pm(istr-1,j)+ &
2184 & grid(ng)%pm(istr ,j))*ad_bry_pgr
2185 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)-adfac
2186 ad_zeta(istr ,j,know)=ad_zeta(istr ,j,know)+adfac
2187 ad_bry_pgr=0.0_r8
2188 END IF
2189# else
2190# ifdef ADJUST_BOUNDARY
2191 IF (lobc(iwest,isubar,ng)) THEN
2192!^ tl_bry_val=BOUNDARY(ng)%tl_ubar_west(j)
2193!^
2194 boundary(ng)%ad_ubar_west(j)=boundary(ng)% &
2195 & ad_ubar_west(j)+ &
2196 & ad_bry_val
2197 ad_bry_val=0.0_r8
2198 ELSE
2199!^ tl_bry_val=0.0_r8
2200!^
2201 ad_bry_val=0.0_r8
2202 END IF
2203# else
2204!^ tl_bry_val=0.0_r8
2205!^
2206 ad_bry_val=0.0_r8
2207# endif
2208# endif
2209 END IF
2210 END DO
2211!
2212! Western edge, clamped boundary condition.
2213!
2214 ELSE IF (ad_lbc(iwest,isubar,ng)%clamped) THEN
2215 DO j=jstr,jend
2216 IF (lbc_apply(ng)%west(j)) THEN
2217# ifdef MASKING
2218!^ tl_ubar(Istr,j,kout)=tl_ubar(Istr,j,kout)* &
2219!^ & GRID(ng)%umask(Istr,j)
2220!^
2221 ad_ubar(istr,j,kout)=ad_ubar(istr,j,kout)* &
2222 & grid(ng)%umask(istr,j)
2223# endif
2224# ifdef ADJUST_BOUNDARY
2225 IF (lobc(iwest,isubar,ng)) THEN
2226!^ tl_ubar(Istr,j,kout)=BOUNDARY(ng)%tl_ubar_west(j)
2227!^
2228 boundary(ng)%ad_ubar_west(j)= &
2229 & boundary(ng)%ad_ubar_west(j)+ &
2230 & ad_ubar(istr,j,kout)
2231 ad_ubar(istr,j,kout)=0.0_r8
2232 ELSE
2233!^ tl_ubar(Istr,j,kout)=0.0_r8
2234!^
2235 ad_ubar(istr,j,kout)=0.0_r8
2236 END IF
2237# else
2238!^ tl_ubar(Istr,j,kout)=0.0_r8
2239!^
2240 ad_ubar(istr,j,kout)=0.0_r8
2241# endif
2242 END IF
2243 END DO
2244!
2245! Western edge, gradient boundary condition.
2246!
2247 ELSE IF (ad_lbc(iwest,isubar,ng)%gradient) THEN
2248 DO j=jstr,jend
2249 IF (lbc_apply(ng)%west(j)) THEN
2250# ifdef MASKING
2251!^ tl_ubar(Istr,j,kout)=tl_ubar(Istr,j,kout)* &
2252!^ & GRID(ng)%umask(Istr,j)
2253!^
2254 ad_ubar(istr,j,kout)=ad_ubar(istr,j,kout)* &
2255 & grid(ng)%umask(istr,j)
2256# endif
2257!^ tl_ubar(Istr ,j,kout)=tl_ubar(Istr+1,j,kout)
2258!^
2259 ad_ubar(istr+1,j,kout)=ad_ubar(istr+1,j,kout)+ &
2260 & ad_ubar(istr,j,kout)
2261 ad_ubar(istr ,j,kout)=0.0_r8
2262 END IF
2263 END DO
2264!
2265! Western edge, reduced-physics boundary condition.
2266!
2267 ELSE IF (ad_lbc(iwest,isubar,ng)%reduced) THEN
2268 DO j=jstr,jend
2269 IF (lbc_apply(ng)%west(j)) THEN
2270 cff=1.0_r8/(0.5_r8*(grid(ng)%h(iend ,j)+ &
2271 & zeta(iend ,j,know)+ &
2272 & grid(ng)%h(iend+1,j)+ &
2273 & zeta(iend+1,j,know)))
2274# ifdef MASKING
2275!^ tl_ubar(Istr,j,kout)=tl_ubar(Istr,j,kout)* &
2276!^ & GRID(ng)%umask(Istr,j)
2277!^
2278 ad_ubar(istr,j,kout)=ad_ubar(istr,j,kout)* &
2279 & grid(ng)%umask(istr,j)
2280# endif
2281!^ tl_ubar(Istr,j,kout)=tl_ubar(Istr,j,know)+ &
2282!^ & dt2d*(tl_bry_pgr+ &
2283!^ & tl_bry_cor+ &
2284!^ & tl_bry_str)
2285!^
2286 adfac=dt2d*ad_ubar(istr,j,kout)
2287 ad_bry_pgr=ad_bry_pgr+adfac
2288 ad_bry_cor=ad_bry_cor+adfac
2289 ad_bry_str=ad_bry_str+adfac
2290 ad_ubar(istr,j,know)=ad_ubar(istr,j,know)+ &
2291 & ad_ubar(istr,j,kout)
2292 ad_ubar(istr,j,kout)=0.0_r8
2293!^ tl_bry_str=tl_cff*(FORCES(ng)%sustr(Istr,j)- &
2294!^ & FORCES(ng)%bustr(Istr,j))+ &
2295!^ & cff*(FORCES(ng)%tl_sustr(Istr,j)- &
2296!^ & FORCES(ng)%tl_bustr(Istr,j))
2297!^
2298 adfac=cff*ad_bry_str
2299 forces(ng)%ad_sustr(istr,j)=forces(ng)%ad_sustr(istr,j)+ &
2300 & adfac
2301 forces(ng)%ad_bustr(istr,j)=forces(ng)%ad_bustr(istr,j)- &
2302 & adfac
2303 ad_cff=ad_cff+(forces(ng)%sustr(istr,j)- &
2304 & forces(ng)%bustr(istr,j))*ad_bry_str
2305 ad_bry_str=0.0_r8
2306!^ tl_cff=-cff*cff*0.5_r8*(GRID(ng)%tl_h(Istr-1,j)+ &
2307!^ & tl_zeta(Istr-1,j,know)+ &
2308!^ & GRID(ng)%tl_h(Istr ,j)+ &
2309!^ & tl_zeta(Istr ,j,know))
2310!^
2311 adfac=-cff*cff*0.5_r8*ad_cff
2312 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)+adfac
2313 ad_zeta(istr ,j,know)=ad_zeta(istr ,j,know)+adfac
2314 grid(ng)%ad_h(istr-1,j)=grid(ng)%ad_h(istr-1,j)+adfac
2315 grid(ng)%ad_h(istr ,j)=grid(ng)%ad_h(istr ,j)+adfac
2316 ad_cff=0.0_r8
2317# ifdef UV_COR
2318!^ tl_bry_cor=0.125_r8*(tl_vbar(Istr-1,j ,know)+ &
2319!^ & tl_vbar(Istr-1,j+1,know)+ &
2320!^ & tl_vbar(Istr ,j ,know)+ &
2321!^ & tl_vbar(Istr ,j+1,know))* &
2322!^ & (GRID(ng)%f(Istr-1,j)+ &
2323!^ & GRID(ng)%f(Istr ,j))
2324!^
2325 adfac=0.125_r8*(grid(ng)%f(istr-1,j)+ &
2326 & grid(ng)%f(istr ,j))*ad_bry_cor
2327 ad_vbar(istr-1,j ,know)=ad_vbar(istr-1,j ,know)+adfac
2328 ad_vbar(istr ,j ,know)=ad_vbar(istr ,j ,know)+adfac
2329 ad_vbar(istr-1,j+1,know)=ad_vbar(istr-1,j+1,know)+adfac
2330 ad_vbar(istr ,j+1,know)=ad_vbar(istr ,j+1,know)+adfac
2331 ad_bry_cor=0.0_r8
2332# else
2333!^ tl_bry_cor=0.0_r8
2334!^
2335 ad_bry_cor=0.0_r8
2336# endif
2337 IF (ad_lbc(iwest,isfsur,ng)%acquire) THEN
2338# ifdef ADJUST_BOUNDARY
2339 IF (lobc(iwest,isubar,ng)) THEN
2340!^ tl_bry_pgr=tl_bry_pgr+ &
2341!^ & g*BOUNDARY(ng)%tl_zeta_west(j)* &
2342!^ & 0.5_r8*GRID(ng)%pm(Istr,j)
2343!^
2344 boundary(ng)%ad_zeta_west(j)=boundary(ng)% &
2345 & ad_zeta_west(j)+ &
2346 & g*0.5_r8* &
2347 & grid(ng)%pm(istr,j)* &
2348 & ad_bry_pgr
2349 END IF
2350# endif
2351!^ tl_bry_pgr=-g*(zeta(Istr,j,know)* &
2352!^ & 0.5_r8*GRID(ng)%pm(Istr,j)
2353!^
2354 ad_zeta(istr,j,know)=ad_zeta(istr,j,know)- &
2355 & g*0.5_r8*grid(ng)%pm(istr,j)* &
2356 & ad_bry_pgr
2357 ad_bry_pgr=0.0_r8
2358 ELSE
2359!^ tl_bry_pgr=-g*(tl_zeta(Istr ,j,know)- &
2360!^ & tl_zeta(Istr-1,j,know))* &
2361!^ & 0.5_r8*(GRID(ng)%pm(Istr-1,j)+ &
2362!^ & GRID(ng)%pm(Istr ,j))
2363!^
2364 adfac=-g*0.5_r8*(grid(ng)%pm(istr-1,j)+ &
2365 & grid(ng)%pm(istr ,j))*ad_bry_pgr
2366 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)-adfac
2367 ad_zeta(istr ,j,know)=ad_zeta(istr ,j,know)+adfac
2368 ad_bry_pgr=0.0_r8
2369 END IF
2370 END IF
2371 END DO
2372!
2373! Western edge, closed boundary condition.
2374!
2375 ELSE IF (ad_lbc(iwest,isubar,ng)%closed) THEN
2376 DO j=jstr,jend
2377 IF (lbc_apply(ng)%west(j)) THEN
2378!^ tl_ubar(Istr,j,kout)=0.0_r8
2379!^
2380 ad_ubar(istr,j,kout)=0.0_r8
2381 END IF
2382 END DO
2383 END IF
2384 END IF
2385
2386 RETURN
type(t_boundary), dimension(:), allocatable boundary
type(t_apply), dimension(:), allocatable lbc_apply
type(t_clima), dimension(:), allocatable clima
Definition mod_clima.F:153
type(t_forces), dimension(:), allocatable forces
Definition mod_forces.F:554
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
integer isfsur
integer isubar
type(t_lbc), dimension(:,:,:), allocatable ad_lbc
Definition mod_param.F:378
type(t_lbc), dimension(:,:,:), allocatable lbc
Definition mod_param.F:375
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
real(r8) co
logical, dimension(:), allocatable lnudgem2clm
integer, dimension(:), allocatable iic
logical, dimension(:,:,:), allocatable lobc
logical, dimension(:), allocatable ewperiodic
integer, parameter iwest
logical, dimension(:), allocatable nsperiodic
logical, dimension(:), allocatable predictor_2d_step
real(dp), dimension(:), allocatable obcfac
real(r8), dimension(:), allocatable gamma2
integer, parameter isouth
real(dp), dimension(:), allocatable dtfast
real(dp), dimension(:,:), allocatable m2obc_out
integer, parameter ieast
real(dp) g
real(dp) rho0
integer, parameter inorth
real(dp), dimension(:,:), allocatable m2obc_in

References mod_param::ad_lbc, mod_boundary::boundary, mod_clima::clima, mod_scalars::co, mod_param::domain, mod_scalars::dtfast, mod_scalars::ewperiodic, mod_forces::forces, mod_scalars::g, mod_scalars::gamma2, mod_grid::grid, mod_scalars::ieast, mod_scalars::iic, mod_scalars::inorth, mod_ncparam::isfsur, mod_scalars::isouth, mod_ncparam::isubar, mod_scalars::iwest, mod_param::lbc, mod_boundary::lbc_apply, mod_scalars::lnudgem2clm, mod_scalars::lobc, mod_scalars::m2obc_in, mod_scalars::m2obc_out, mod_scalars::nsperiodic, mod_scalars::obcfac, mod_scalars::predictor_2d_step, and mod_scalars::rho0.

Referenced by ad_ini_fields_mod::ad_ini_fields_tile(), ad_ini_fields_mod::ad_out_fields_tile(), ad_step2d_mod::ad_step2d_tile(), ad_step2d_mod::ad_step2d_tile(), ad_step2d_mod::ad_step2d_tile(), and ad_u2dbc().

Here is the caller graph for this function: