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

Functions/Subroutines

subroutine, public ad_zetabc (ng, tile, kout)
 
subroutine, public ad_zetabc_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, krhs, kstp, kout, zeta, ad_zeta)
 

Function/Subroutine Documentation

◆ ad_zetabc()

subroutine, public ad_zetabc_mod::ad_zetabc ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) kout )

Definition at line 27 of file ad_zetabc.F.

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

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

Here is the call graph for this function:

◆ ad_zetabc_tile()

subroutine, public ad_zetabc_mod::ad_zetabc_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) zeta,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_zeta )

Definition at line 52 of file ad_zetabc.F.

57!***********************************************************************
58!
59 USE mod_param
60 USE mod_boundary
61 USE mod_grid
62 USE mod_ncparam
63 USE mod_scalars
64!
65! Imported variable declarations.
66!
67 integer, intent(in) :: ng, tile
68 integer, intent(in) :: LBi, UBi, LBj, UBj
69 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
70 integer, intent(in) :: krhs, kstp, kout
71!
72# ifdef ASSUMED_SHAPE
73 real(r8), intent(in) :: zeta(LBi:,LBj:,:)
74
75 real(r8), intent(inout) :: ad_zeta(LBi:,LBj:,:)
76# else
77 real(r8), intent(in) :: zeta(LBi:UBi,LBj:UBj,3)
78
79 real(r8), intent(inout) :: ad_zeta(LBi:UBi,LBj:UBj,3)
80# endif
81!
82! Local variable declarations.
83!
84 integer :: i, j, know
85
86 real(r8) :: Ce, Cx
87 real(r8) :: cff, cff1, cff2, dt2d, tau
88
89 real(r8) :: ad_Ce, ad_Cx
90 real(r8) :: ad_cff1, ad_cff2
91 real(r8) :: adfac
92
93 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_grad
94
95# include "set_bounds.h"
96!
97!-----------------------------------------------------------------------
98! Initialize adjoint private variables.
99!-----------------------------------------------------------------------
100!
101 ad_ce=0.0_r8
102 ad_cx=0.0_r8
103 ad_cff1=0.0_r8
104 ad_cff2=0.0_r8
105
106 ad_grad(lbi:ubi,lbj:ubj)=0.0_r8
107!
108!-----------------------------------------------------------------------
109! Set time-indices
110!-----------------------------------------------------------------------
111!
112 IF (first_2d_step) THEN
113 know=krhs
114 dt2d=dtfast(ng)
115 ELSE IF (predictor_2d_step(ng)) THEN
116 know=krhs
117 dt2d=2.0_r8*dtfast(ng)
118 ELSE
119 know=kstp
120 dt2d=dtfast(ng)
121 END IF
122
123# if defined WET_DRY
124!
125!-----------------------------------------------------------------------
126! Ensure that water level on boundary cells is above bed elevation.
127!-----------------------------------------------------------------------
128!
129 cff=dcrit(ng)-eps
130
131 IF (.not.(ewperiodic(ng).or.nsperiodic(ng))) THEN
132 IF (domain(ng)%NorthEast_Corner(tile)) THEN
133 IF (lbc_apply(ng)%north(iend+1).and. &
134 & lbc_apply(ng)%east (jend+1)) THEN
135 IF (zeta(iend+1,jend+1,kout).le. &
136 & (dcrit(ng)-grid(ng)%h(iend+1,jend+1))) THEN
137!^ tl_zeta(Iend+1,Jend+1,kout)=-GRID(ng)%tl_h(Iend+1,Jend+1)
138!^
139 grid(ng)%ad_h(iend+1,jend+1)=grid(ng)%ad_h(iend+1,jend+1)-&
140 & ad_zeta(iend+1,jend+1,kout)
141 ad_zeta(iend+1,jend+1,kout)=0.0_r8
142 END IF
143 END IF
144 END IF
145 IF (domain(ng)%NorthWest_Corner(tile)) THEN
146 IF (lbc_apply(ng)%north(istr-1).and. &
147 & lbc_apply(ng)%west (jend+1)) THEN
148 IF (zeta(istr-1,jend+1,kout).le. &
149 & (dcrit(ng)-grid(ng)%h(istr-1,jend+1))) THEN
150!^ tl_zeta(Istr-1,Jend+1,kout)=-GRID(ng)%tl_h(Istr-1,Jend+1)
151!^
152 grid(ng)%ad_h(istr-1,jend+1)=grid(ng)%ad_h(istr-1,jend+1)-&
153 & ad_zeta(istr-1,jend+1,kout)
154 ad_zeta(istr-1,jend+1,kout)=0.0_r8
155 END IF
156 END IF
157 END IF
158 IF (domain(ng)%SouthEast_Corner(tile)) THEN
159 IF (lbc_apply(ng)%south(iend+1).and. &
160 & lbc_apply(ng)%east (jstr-1)) THEN
161 IF (zeta(iend+1,jstr-1,kout).le. &
162 & (dcrit(ng)-grid(ng)%h(iend+1,jstr-1))) THEN
163!^ tl_zeta(Iend+1,Jstr-1,kout)=-GRID(ng)%tl_h(Iend+1,Jstr-1)
164!^
165 grid(ng)%ad_h(iend+1,jstr-1)=grid(ng)%ad_h(iend+1,jstr-1)-&
166 & ad_zeta(iend+1,jstr-1,kout)
167 tl_zeta(iend+1,jstr-1,kout)=0.0_r8
168 END IF
169 END IF
170 END IF
171 IF (domain(ng)%SouthWest_Corner(tile)) THEN
172 IF (lbc_apply(ng)%south(istr-1).and. &
173 & lbc_apply(ng)%west (jstr-1)) THEN
174 IF (zeta(istr-1,jstr-1,kout).le. &
175 & (dcrit(ng)-grid(ng)%h(istr-1,jstr-1))) THEN
176!^ tl_zeta(Istr-1,Jstr-1,kout)=-GRID(ng)%tl_h(Istr-1,Jstr-1)
177!^
178 grid(ng)%ad_h(istr-1,jstr-1)=grid(ng)%ad_h(istr-1,jstr-1)-&
179 & ad_zeta(istr-1,jstr-1,kout)
180 ad_zeta(istr-1,jstr-1,kout)=0.0_r8
181 END IF
182 END IF
183 END IF
184 END IF
185
186 IF (.not.nsperiodic(ng)) THEN
187 IF (domain(ng)%Northern_Edge(tile)) THEN
188 DO i=istr,iend
189 IF (lbc_apply(ng)%north(i)) THEN
190 IF (zeta(i,jend+1,kout).le. &
191 & (dcrit(ng)-grid(ng)%h(i,jend+1))) THEN
192!^ tl_zeta(i,Jend+1,kout)=-GRID(ng)%tl_h(i,Jend+1)
193!^
194 grid(ng)%ad_h(i,jend+1)=grid(ng)%ad_h(i,jend+1)- &
195 & ad_zeta(i,jend+1,kout)
196 ad_zeta(i,jend+1,kout)=0.0_r8
197 END IF
198 END IF
199 END DO
200 END IF
201 IF (domain(ng)%Southern_Edge(tile)) THEN
202 DO i=istr,iend
203 IF (lbc_apply(ng)%south(i)) THEN
204 IF (zeta(i,jstr-1,kout).le. &
205 & (dcrit(ng)-grid(ng)%h(i,jstr-1))) THEN
206!^ tl_zeta(i,Jstr-1,kout)=-GRID(ng)%tl_h(i,Jstr-1)
207!^
208 grid(ng)%ad_h(i,jstr-1)=grid(ng)%ad_h(i,jstr-1)- &
209 & ad_zeta(i,jstr-1,kout)
210 ad_zeta(i,jstr-1,kout)=0.0_r8
211 END IF
212 END IF
213 END DO
214 END IF
215 END IF
216
217 IF (.not.ewperiodic(ng)) THEN
218 IF (domain(ng)%Eastern_Edge(tile)) THEN
219 DO j=jstr,jend
220 IF (lbc_apply(ng)%east(j)) THEN
221 IF (zeta(iend+1,j,kout).le. &
222 & (dcrit(ng)-grid(ng)%h(iend+1,j))) THEN
223!^ tl_zeta(Iend+1,j,kout)=-GRID(ng)%tl_h(Iend+1,j)
224!^
225 grid(ng)%ad_h(iend+1,j)=grid(ng)%ad_h(iend+1,j)- &
226 & ad_zeta(iend+1,j,kout)
227 ad_zeta(iend+1,j,kout)=0.0_r8
228 END IF
229 END IF
230 END DO
231 END IF
232 IF (domain(ng)%Western_Edge(tile)) THEN
233 DO j=jstr,jend
234 IF (lbc_apply(ng)%west(j)) THEN
235 IF (zeta(istr-1,j,kout).le. &
236 & (dcrit(ng)-grid(ng)%h(istr-1,j))) THEN
237!^ tl_zeta(Istr-1,j,kout)=-GRID(ng)%tl_h(Istr-1,j)
238!^
239 grid(ng)%ad_h(istr-1,j)=grid(ng)%ad_h(istr-1,j)- &
240 & ad_zeta(istr-1,j,kout)
241 ad_zeta(istr-1,j,kout)=0.0_r8
242 END IF
243 END IF
244 END DO
245 END IF
246 END IF
247# endif
248!
249!-----------------------------------------------------------------------
250! Boundary corners.
251!-----------------------------------------------------------------------
252!
253 IF (.not.(ewperiodic(ng).or.nsperiodic(ng))) THEN
254 IF (domain(ng)%NorthEast_Corner(tile)) THEN
255 IF (lbc_apply(ng)%north(iend+1).and. &
256 & lbc_apply(ng)%east (jend+1)) THEN
257!^ tl_zeta(Iend+1,Jend+1,kout)=0.5_r8* &
258!^ & (tl_zeta(Iend+1,Jend ,kout)+ &
259!^ & tl_zeta(Iend ,Jend+1,kout))
260!^
261 adfac=0.5_r8*ad_zeta(iend+1,jend+1,kout)
262 ad_zeta(iend ,jend+1,kout)=ad_zeta(iend ,jend+1,kout)+ &
263 & adfac
264 ad_zeta(iend+1,jend ,kout)=ad_zeta(iend+1,jend ,kout)+ &
265 & adfac
266 ad_zeta(iend+1,jend+1,kout)=0.0_r8
267 END IF
268 END IF
269 IF (domain(ng)%NorthWest_Corner(tile)) THEN
270 IF (lbc_apply(ng)%north(istr-1).and. &
271 & lbc_apply(ng)%west (jend+1)) THEN
272!^ tl_zeta(Istr-1,Jend+1,kout)=0.5_r8* &
273!^ & (tl_zeta(Istr-1,Jend ,kout)+ &
274!^ & tl_zeta(Istr ,Jend+1,kout))
275!^
276 adfac=0.5_r8*ad_zeta(istr-1,jend+1,kout)
277 ad_zeta(istr-1,jend ,kout)=ad_zeta(istr-1,jend ,kout)+ &
278 & adfac
279 ad_zeta(istr ,jend+1,kout)=ad_zeta(istr ,jend+1,kout)+ &
280 & adfac
281 ad_zeta(istr-1,jend+1,kout)=0.0_r8
282 END IF
283 END IF
284 IF (domain(ng)%SouthEast_Corner(tile)) THEN
285 IF (lbc_apply(ng)%south(iend+1).and. &
286 & lbc_apply(ng)%east (jstr-1)) THEN
287!^ tl_zeta(Iend+1,Jstr-1,kout)=0.5_r8* &
288!^ & (tl_zeta(Iend ,Jstr-1,kout)+ &
289!^ & tl_zeta(Iend+1,Jstr ,kout))
290!^
291 adfac=0.5_r8*ad_zeta(iend+1,jstr-1,kout)
292 ad_zeta(iend ,jstr-1,kout)=ad_zeta(iend ,jstr-1,kout)+ &
293 & adfac
294 ad_zeta(iend+1,jstr ,kout)=ad_zeta(iend+1,jstr ,kout)+ &
295 & adfac
296 ad_zeta(iend+1,jstr-1,kout)=0.0_r8
297 END IF
298 END IF
299 IF (domain(ng)%SouthWest_Corner(tile)) THEN
300 IF (lbc_apply(ng)%south(istr-1).and. &
301 & lbc_apply(ng)%west (jstr-1)) THEN
302!^ tl_zeta(Istr-1,Jstr-1,kout)=0.5_r8* &
303!^ & (tl_zeta(Istr ,Jstr-1,kout)+ &
304!^ & tl_zeta(Istr-1,Jstr ,kout))
305!^
306 adfac=0.5_r8*ad_zeta(istr-1,jstr-1,kout)
307 ad_zeta(istr ,jstr-1,kout)=ad_zeta(istr ,jstr-1,kout)+ &
308 & adfac
309 ad_zeta(istr-1,jstr ,kout)=ad_zeta(istr-1,jstr ,kout)+ &
310 & adfac
311 ad_zeta(istr-1,jstr-1,kout)=0.0_r8
312 END IF
313 END IF
314 END IF
315!
316!-----------------------------------------------------------------------
317! Lateral boundary conditions at the northern edge.
318!-----------------------------------------------------------------------
319!
320 IF (domain(ng)%Northern_Edge(tile)) THEN
321!
322! Northern edge, implicit upstream radiation condition.
323!
324 IF (ad_lbc(inorth,isfsur,ng)%radiation) THEN
325 IF (iic(ng).ne.0) THEN
326 DO i=istr,iend
327 IF (lbc_apply(ng)%north(i)) THEN
328# if defined CELERITY_READ && defined FORWARD_READ
329 IF (ad_lbc(inorth,isfsur,ng)%nudging) THEN
330 IF (boundary(ng)%zeta_north_Ce(i).eq.0.0_r8) THEN
331 tau=fsobc_in(ng,inorth)
332 ELSE
333 tau=fsobc_out(ng,inorth)
334 END IF
335 tau=tau*dt2d
336 END IF
337# ifdef RADIATION_2D
338 cx=boundary(ng)%zeta_north_Cx(i)
339# else
340 cx=0.0_r8
341# endif
342 ce=boundary(ng)%zeta_north_Ce(i)
343 cff=boundary(ng)%zeta_north_C2(i)
344# endif
345# ifdef MASKING
346!^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend+1,kout)* &
347!^ & GRID(ng)%rmask(i,Jend+1)
348!^
349 ad_zeta(i,jend+1,kout)=ad_zeta(i,jend+1,kout)* &
350 & grid(ng)%rmask(i,jend+1)
351# endif
352 IF (ad_lbc(inorth,isfsur,ng)%nudging) THEN
353!^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend+1,kout)- &
354!^ & tau*tl_zeta(i,Jend+1,know)
355!^
356 ad_zeta(i,jend+1,know)=ad_zeta(i,jend+1,know)- &
357 & tau*ad_zeta(i,jend+1,kout)
358 END IF
359!^ tl_zeta(i,Jend+1,kout)=(cff*tl_zeta(i,Jend+1,know)+ &
360!^ & Ce *tl_zeta(i,Jend ,kout)- &
361!^ & MAX(Cx,0.0_r8)* &
362!^ & tl_grad(i ,Jend+1)- &
363!^ & MIN(Cx,0.0_r8)* &
364!^ & tl_grad(i+1,Jend+1))/ &
365!^ & (cff+Ce)
366!^
367 adfac=ad_zeta(i,jend+1,kout)/(cff+ce)
368 ad_grad(i ,jend+1)=ad_grad(i ,jend+1)- &
369 & max(cx,0.0_r8)*adfac
370 ad_grad(i+1,jend+1)=ad_grad(i+1,jend+1)- &
371 & min(cx,0.0_r8)*adfac
372 ad_zeta(i,jend ,kout)=ad_zeta(i,jend ,kout)+ce *adfac
373 ad_zeta(i,jend+1,know)=ad_zeta(i,jend+1,know)+cff*adfac
374 ad_zeta(i,jend+1,kout)=0.0_r8
375 END IF
376 END DO
377 END IF
378!
379! Northern edge, explicit Chapman boundary condition.
380!
381 ELSE IF (ad_lbc(inorth,isfsur,ng)%Chapman_explicit) THEN
382 DO i=istr,iend
383 IF (lbc_apply(ng)%north(i)) THEN
384 cff=dt2d*grid(ng)%pn(i,jend)
385 cff1=sqrt(g*(grid(ng)%h(i,jend)+ &
386 & zeta(i,jend,know)))
387 ce=cff*cff1
388# ifdef MASKING
389!^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend+1,kout)* &
390!^ & GRID(ng)%rmask(i,Jend+1)
391!^
392 ad_zeta(i,jend+1,kout)=ad_zeta(i,jend+1,kout)* &
393 & grid(ng)%rmask(i,jend+1)
394# endif
395!^ tl_zeta(i,Jend+1,kout)=(1.0_r8-Ce)*tl_zeta(i,Jend+1,know)+&
396!^ & tl_Ce*(zeta(i,Jend+1,know)+ &
397!^ & zeta(i,Jend ,know))+ &
398!^ & Ce*tl_zeta(i,Jend,know)
399!^
400 ad_zeta(i,jend+1,know)=ad_zeta(i,jend+1,know)+ &
401 & (1.0_r8-ce)*ad_zeta(i,jend+1,kout)
402 ad_ce=ad_ce+(zeta(i,jend+1,know)+ &
403 & zeta(i,jend ,know))*ad_zeta(i,jend+1,kout)
404 ad_zeta(i,jend,know)=ad_zeta(i,jend,know)+ &
405 & ce*ad_zeta(i,jend+1,kout)
406 ad_zeta(i,jend+1,kout)=0.0_r8
407!^ tl_Ce=cff*tl_cff1
408!^
409 ad_cff1=ad_cff1+cff*ad_ce
410 ad_ce=0.0_r8
411!^ tl_cff1=0.5_r8*g*(GRID(ng)%tl_h(i,Jend)+ &
412!^ & tl_zeta(i,Jend,know))/cff1
413!^
414 adfac=0.5_r8*g*ad_cff1/cff1
415 grid(ng)%ad_h(i,jend)=grid(ng)%ad_h(i,jend)+adfac
416 ad_zeta(i,jend,know)=ad_zeta(i,jend,know)+adfac
417 ad_cff1=0.0_r8
418 END IF
419 END DO
420!
421! Northern edge, implicit Chapman boundary condition.
422!
423 ELSE IF (ad_lbc(inorth,isfsur,ng)%Chapman_implicit) THEN
424 DO i=istr,iend
425 IF (lbc_apply(ng)%north(i)) THEN
426 cff=dt2d*grid(ng)%pn(i,jend)
427 cff1=sqrt(g*(grid(ng)%h(i,jend)+ &
428 & zeta(i,jend,know)))
429 ce=cff*cff1
430 cff2=1.0_r8/(1.0_r8+ce)
431# ifdef MASKING
432!^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend+1,kout)* &
433!^ & GRID(ng)%rmask(i,Jend+1)
434!^
435 ad_zeta(i,jend+1,kout)=ad_zeta(i,jend+1,kout)* &
436 & grid(ng)%rmask(i,jend+1)
437# endif
438!^ tl_zeta(i,Jend+1,kout)=tl_cff2*(zeta(i,Jend+1,know)+ &
439!^ & Ce*zeta(i,Jend,kout))+ &
440!^ & cff2*(tl_zeta(i,Jend+1,know)+ &
441!^ & tl_Ce*zeta(i,Jend,kout)+ &
442!^ & Ce*tl_zeta(i,Jend,kout))
443!^
444 adfac=cff2*ad_zeta(i,jend+1,kout)
445 ad_zeta(i,jend ,kout)=ad_zeta(i,jend ,kout)+ce*adfac
446 ad_zeta(i,jend+1,know)=ad_zeta(i,jend+1,know)+adfac
447 ad_ce=ad_ce+zeta(i,jend,kout)*adfac
448 ad_cff2=ad_cff2+ &
449 & (zeta(i,jend+1,know)+ &
450 & ce*zeta(i,jend,kout))*ad_zeta(i,jend+1,kout)
451 ad_zeta(i,jend+1,kout)=0.0_r8
452!^ tl_cff2=-cff2*cff2*tl_Ce
453!^
454 ad_ce=ad_ce-cff2*cff2*ad_cff2
455 ad_cff2=0.0_r8
456!^ tl_Ce=cff*tl_cff1
457!^
458 ad_cff1=ad_cff1+cff*ad_ce
459 ad_ce=0.0_r8
460!^ tl_cff1=0.5_r8*g*(GRID(ng)%tl_h(i,Jend)+ &
461!^ & tl_zeta(i,Jend,know))/cff1
462!^
463 adfac=0.5_r8*g*ad_cff1/cff1
464 grid(ng)%ad_h(i,jend)=grid(ng)%ad_h(i,jend)+adfac
465 ad_zeta(i,jend,know)=ad_zeta(i,jend,know)+adfac
466 ad_cff1=0.0_r8
467 END IF
468 END DO
469!
470! Northern edge, clamped boundary condition.
471!
472 ELSE IF (ad_lbc(inorth,isfsur,ng)%clamped) THEN
473 DO i=istr,iend
474 IF (lbc_apply(ng)%north(i)) THEN
475# ifdef MASKING
476!^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend+1,kout)* &
477!^ & GRID(ng)%rmask(i,Jend+1)
478!^
479 ad_zeta(i,jend+1,kout)=ad_zeta(i,jend+1,kout)* &
480 & grid(ng)%rmask(i,jend+1)
481# endif
482# ifdef ADJUST_BOUNDARY
483 IF (lobc(inorth,isfsur,ng)) THEN
484!^ tl_zeta(i,Jend+1,kout)=BOUNDARY(ng)%tl_zeta_north(i)
485!^
486 boundary(ng)%ad_zeta_north(i)=boundary(ng)% &
487 & ad_zeta_north(i)+ &
488 & ad_zeta(i,jend+1,kout)
489 ad_zeta(i,jend+1,kout)=0.0_r8
490 ELSE
491!^ tl_zeta(i,Jend+1,kout)=0.0_r8
492!^
493 ad_zeta(i,jend+1,kout)=0.0_r8
494 END IF
495# else
496!^ tl_zeta(i,Jend+1,kout)=0.0_r8
497!^
498 ad_zeta(i,jend+1,kout)=0.0_r8
499# endif
500 END IF
501 END DO
502!
503! Northern edge, gradient boundary condition.
504!
505 ELSE IF (ad_lbc(inorth,isfsur,ng)%gradient) THEN
506 DO i=istr,iend
507 IF (lbc_apply(ng)%north(i)) THEN
508# ifdef MASKING
509!^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend+1,kout)* &
510!^ & GRID(ng)%rmask(i,Jend+1)
511!^
512 ad_zeta(i,jend+1,kout)=ad_zeta(i,jend+1,kout)* &
513 & grid(ng)%rmask(i,jend+1)
514# endif
515!^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend,kout)
516!^
517 ad_zeta(i,jend ,kout)=ad_zeta(i,jend ,kout)+ &
518 & ad_zeta(i,jend+1,kout)
519 ad_zeta(i,jend+1,kout)=0.0_r8
520 END IF
521 END DO
522!
523! Northern edge, closed boundary condition.
524!
525 ELSE IF (ad_lbc(inorth,isfsur,ng)%closed) THEN
526 DO i=istr,iend
527 IF (lbc_apply(ng)%north(i)) THEN
528# ifdef MASKING
529!^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend+1,kout)* &
530!^ & GRID(ng)%rmask(i,Jend+1)
531!^
532 ad_zeta(i,jend+1,kout)=ad_zeta(i,jend+1,kout)* &
533 & grid(ng)%rmask(i,jend+1)
534# endif
535!^ tl_zeta(i,Jend+1,kout)=tl_zeta(i,Jend,kout)
536!^
537 ad_zeta(i,jend ,kout)=ad_zeta(i,jend ,kout)+ &
538 & ad_zeta(i,jend+1,kout)
539 ad_zeta(i,jend+1,kout)=0.0_r8
540 END IF
541 END DO
542 END IF
543 END IF
544!
545!-----------------------------------------------------------------------
546! Lateral boundary conditions at the southern edge.
547!-----------------------------------------------------------------------
548!
549 IF (domain(ng)%Southern_Edge(tile)) THEN
550!
551! Southern edge, implicit upstream radiation condition.
552!
553 IF (ad_lbc(isouth,isfsur,ng)%radiation) THEN
554 IF (iic(ng).ne.0) THEN
555 DO i=istr,iend
556 IF (lbc_apply(ng)%south(i)) THEN
557# if defined CELERITY_READ && defined FORWARD_READ
558 IF (ad_lbc(isouth,isfsur,ng)%nudging) THEN
559 IF (boundary(ng)%zeta_south_Ce(i).eq.0.0_r8) THEN
560 tau=fsobc_in(ng,isouth)
561 ELSE
562 tau=fsobc_out(ng,isouth)
563 END IF
564 tau=tau*dt2d
565 END IF
566# ifdef RADIATION_2D
567 cx=boundary(ng)%zeta_south_Cx(i)
568# else
569 cx=0.0_r8
570# endif
571 ce=boundary(ng)%zeta_south_Ce(i)
572 cff=boundary(ng)%zeta_south_C2(i)
573# endif
574# ifdef MASKING
575!^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr-1,kout)* &
576!^ & GRID(ng)%rmask(i,Jstr-1)
577!^
578 ad_zeta(i,jstr-1,kout)=ad_zeta(i,jstr-1,kout)* &
579 & grid(ng)%rmask(i,jstr-1)
580# endif
581 IF (ad_lbc(isouth,isfsur,ng)%nudging) THEN
582!^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr-1,kout)- &
583!^ & tau*tl_zeta(i,Jstr-1,know)
584!^
585 ad_zeta(i,jstr-1,know)=ad_zeta(i,jstr-1,know)- &
586 & tau*ad_zeta(i,jstr-1,kout)
587 END IF
588!^ tl_zeta(i,Jstr-1,kout)=(cff*tl_zeta(i,Jstr-1,know)+ &
589!^ & Ce *tl_zeta(i,Jstr ,kout)- &
590!^ & MAX(Cx,0.0_r8)* &
591!^ & tl_grad(i ,Jstr-1)- &
592!^ & MIN(Cx,0.0_r8)* &
593!^ & tl_grad(i+1,Jstr-1))/ &
594!^ & (cff+Ce)
595!^
596 adfac=ad_zeta(i,jstr-1,kout)/(cff+ce)
597 ad_grad(i ,jstr-1)=ad_grad(i ,jstr-1)- &
598 & max(cx,0.0_r8)*adfac
599 ad_grad(i+1,jstr-1)=ad_grad(i+1,jstr-1)- &
600 & min(cx,0.0_r8)*adfac
601 ad_zeta(i,jstr-1,know)=ad_zeta(i,jstr-1,know)+cff*adfac
602 ad_zeta(i,jstr ,kout)=ad_zeta(i,jstr ,kout)+ce *adfac
603 ad_zeta(i,jstr-1,kout)=0.0_r8
604 END IF
605 END DO
606 END IF
607!
608! Southern edge, explicit Chapman boundary condition.
609!
610 ELSE IF (ad_lbc(isouth,isfsur,ng)%Chapman_explicit) THEN
611 DO i=istr,iend
612 IF (lbc_apply(ng)%south(i)) THEN
613 cff=dt2d*grid(ng)%pn(i,jstr)
614 cff1=sqrt(g*(grid(ng)%h(i,jstr)+ &
615 & zeta(i,jstr,know)))
616 ce=cff*cff1
617# ifdef MASKING
618!^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr-1,kout)* &
619!^ & GRID(ng)%rmask(i,Jstr-1)
620!^
621 ad_zeta(i,jstr-1,kout)=ad_zeta(i,jstr-1,kout)* &
622 & grid(ng)%rmask(i,jstr-1)
623# endif
624!^ tl_zeta(i,Jstr-1,kout)=(1.0_r8-Ce)*tl_zeta(i,Jstr-1,know)+&
625!^ & tl_Ce*(zeta(i,Jstr-1,know)+ &
626!^ & zeta(i,Jstr ,know))+ &
627!^ & Ce*tl_zeta(i,Jstr,know)
628!^
629 ad_zeta(i,jstr-1,know)=ad_zeta(i,jstr-1,know)+ &
630 & (1.0_r8-ce)*ad_zeta(i,jstr-1,kout)
631 ad_ce=ad_ce+(zeta(i,jstr-1,know)+ &
632 & zeta(i,jstr ,know))*ad_zeta(i,jstr-1,kout)
633 ad_zeta(i,jstr,know)=ad_zeta(i,jstr,know)+ &
634 & ce*ad_zeta(i,jstr-1,kout)
635 ad_zeta(i,jstr-1,kout)=0.0_r8
636!^ tl_Ce=cff*tl_cff1
637!^
638 ad_cff1=ad_cff1+cff*ad_ce
639 ad_ce=0.0_r8
640!^ tl_cff1=0.5_r8*g*(GRID(ng)%tl_h(i,Jstr)+ &
641!^ & tl_zeta(i,Jstr,know))/cff1
642!^
643 adfac=0.5_r8*g*ad_cff1/cff1
644 grid(ng)%ad_h(i,jstr)=grid(ng)%ad_h(i,jstr)+adfac
645 ad_zeta(i,jstr,know)=ad_zeta(i,jstr,know)+adfac
646 ad_cff1=0.0_r8
647 END IF
648 END DO
649!
650! Southern edge, implicit Chapman boundary condition.
651!
652 ELSE IF (ad_lbc(isouth,isfsur,ng)%Chapman_implicit) THEN
653 DO i=istr,iend
654 IF (lbc_apply(ng)%south(i)) THEN
655 cff=dt2d*grid(ng)%pn(i,jstr)
656 cff1=sqrt(g*(grid(ng)%h(i,jstr)+ &
657 & zeta(i,jstr,know)))
658 ce=cff*cff1
659 cff2=1.0_r8/(1.0_r8+ce)
660# ifdef MASKING
661!^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr-1,kout)* &
662!^ & GRID(ng)%rmask(i,Jstr-1)
663!^
664 ad_zeta(i,jstr-1,kout)=ad_zeta(i,jstr-1,kout)* &
665 & grid(ng)%rmask(i,jstr-1)
666# endif
667!^ tl_zeta(i,Jstr-1,kout)=tl_cff2*(zeta(i,Jstr-1,know)+ &
668!^ & Ce*zeta(i,Jstr,kout))+ &
669!^ & cff2*(tl_zeta(i,Jstr-1,know)+ &
670!^ & tl_Ce*zeta(i,Jstr,kout)+ &
671!^ & Ce*tl_zeta(i,Jstr,kout))
672!^
673 adfac=cff2*ad_zeta(i,jstr-1,kout)
674 ad_zeta(i,jstr-1,know)=ad_zeta(i,jstr-1,know)+adfac
675 ad_zeta(i,jstr ,kout)=ad_zeta(i,jstr ,kout)+ce*adfac
676 ad_ce=ad_ce+zeta(i,jstr,kout)*adfac
677 ad_cff2=ad_cff2+ &
678 & (zeta(i,jstr-1,know)+ &
679 & ce*zeta(i,jstr,kout))*ad_zeta(i,jstr-1,kout)
680 ad_zeta(i,jstr-1,kout)=0.0_r8
681!^ tl_cff2=-cff2*cff2*tl_Ce
682!^
683 ad_ce=ad_ce-cff2*cff2*ad_cff2
684 ad_cff2=0.0_r8
685!^ tl_Ce=cff*tl_cff1
686!^
687 ad_cff1=ad_cff1+cff*ad_ce
688 ad_ce=0.0_r8
689!^ tl_cff1=0.5_r8*g*(GRID(ng)%tl_h(i,Jstr)+ &
690!^ & tl_zeta(i,Jstr,know))/cff1
691!^
692 adfac=0.5_r8*g*ad_cff1/cff1
693 grid(ng)%ad_h(i,jstr)=grid(ng)%ad_h(i,jstr)+adfac
694 ad_zeta(i,jstr,know)=ad_zeta(i,jstr,know)+adfac
695 ad_cff1=0.0_r8
696 END IF
697 END DO
698!
699! Southern edge, clamped boundary condition.
700!
701 ELSE IF (ad_lbc(isouth,isfsur,ng)%clamped) THEN
702 DO i=istr,iend
703 IF (lbc_apply(ng)%south(i)) THEN
704# ifdef MASKING
705!^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr-1,kout)* &
706!^ & GRID(ng)%rmask(i,Jstr-1)
707!^
708 ad_zeta(i,jstr-1,kout)=ad_zeta(i,jstr-1,kout)* &
709 & grid(ng)%rmask(i,jstr-1)
710# endif
711# ifdef ADJUST_BOUNDARY
712 IF (lobc(isouth,isfsur,ng)) THEN
713!^ tl_zeta(i,Jstr-1,kout)=BOUNDARY(ng)%tl_zeta_south(i)
714!^
715 boundary(ng)%ad_zeta_south(i)=boundary(ng)% &
716 & ad_zeta_south(i)+ &
717 & ad_zeta(i,jstr-1,kout)
718 ad_zeta(i,jstr-1,kout)=0.0_r8
719 ELSE
720!^ tl_zeta(i,Jstr-1,kout)=0.0_r8
721!^
722 ad_zeta(i,jstr-1,kout)=0.0_r8
723 END IF
724# else
725!^ tl_zeta(i,Jstr-1,kout)=0.0_r8
726!^
727 ad_zeta(i,jstr-1,kout)=0.0_r8
728# endif
729 END IF
730 END DO
731!
732! Southern edge, gradient boundary condition.
733!
734 ELSE IF (ad_lbc(isouth,isfsur,ng)%gradient) THEN
735 DO i=istr,iend
736 IF (lbc_apply(ng)%south(i)) THEN
737# ifdef MASKING
738!^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr-1,kout)* &
739!^ & GRID(ng)%rmask(i,Jstr-1)
740!^
741 ad_zeta(i,jstr-1,kout)=ad_zeta(i,jstr-1,kout)* &
742 & grid(ng)%rmask(i,jstr-1)
743# endif
744!^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr,kout)
745!^
746 ad_zeta(i,jstr ,kout)=ad_zeta(i,jstr ,kout)+ &
747 & ad_zeta(i,jstr-1,kout)
748 ad_zeta(i,jstr-1,kout)=0.0_r8
749 END IF
750 END DO
751!
752! Southern edge, closed boundary condition.
753!
754 ELSE IF (ad_lbc(isouth,isfsur,ng)%closed) THEN
755 DO i=istr,iend
756 IF (lbc_apply(ng)%south(i)) THEN
757# ifdef MASKING
758!^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr-1,kout)* &
759!^ & GRID(ng)%rmask(i,Jstr-1)
760!^
761 ad_zeta(i,jstr-1,kout)=ad_zeta(i,jstr-1,kout)* &
762 & grid(ng)%rmask(i,jstr-1)
763# endif
764!^ tl_zeta(i,Jstr-1,kout)=tl_zeta(i,Jstr,kout)
765!^
766 ad_zeta(i,jstr ,kout)=ad_zeta(i,jstr ,kout)+ &
767 & ad_zeta(i,jstr-1,kout)
768 ad_zeta(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,isfsur,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,isfsur,ng)%nudging) THEN
788 IF (boundary(ng)%zeta_east_Cx(j).eq.0.0_r8) THEN
789 tau=fsobc_in(ieast)
790 ELSE
791 tau=fsobc_out(ieast)
792 END IF
793 tau=tau*dt2d
794 END IF
795 cx=boundary(ng)%zeta_east_Cx(j)
796# ifdef RADIATION_2D
797 ce=boundary(ng)%zeta_east_Ce(j)
798# else
799 ce=0.0_r8
800# endif
801 cff=boundary(ng)%zeta_east_C2(j)
802# endif
803# ifdef MASKING
804!^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend+1,j,kout)* &
805!^ & GRID(ng)%rmask(Iend+1,j)
806!^
807 ad_zeta(iend+1,j,kout)=ad_zeta(iend+1,j,kout)* &
808 & grid(ng)%rmask(iend+1,j)
809# endif
810 IF (ad_lbc(ieast,isfsur,ng)%nudging) THEN
811!^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend+1,j,kout)- &
812!^ & tau*tl_zeta(Iend+1,j,know)
813!^
814 ad_zeta(iend+1,j,know)=ad_zeta(iend+1,j,know)- &
815 & tau*ad_zeta(iend+1,j,kout)
816 END IF
817!^ tl_zeta(Iend+1,j,kout)=(cff*tl_zeta(Iend+1,j,know)+ &
818!^ & Cx *tl_zeta(Iend ,j,kout)- &
819!^ & MAX(Ce,0.0_r8)* &
820!^ & tl_grad(Iend+1,j )- &
821!^ & MIN(Ce,0.0_r8)* &
822!^ & tl_grad(Iend+1,j+1))/ &
823!^ & (cff+Cx)
824!^
825 adfac=ad_zeta(iend+1,j,kout)/(cff+cx)
826 ad_grad(iend+1,j )=ad_grad(iend+1,j )- &
827 & max(ce,0.0_r8)*adfac
828 ad_grad(iend+1,j+1)=ad_grad(iend+1,j+1)- &
829 & min(ce,0.0_r8)*adfac
830 ad_zeta(iend ,j,kout)=ad_zeta(iend ,j,kout)+cx *adfac
831 ad_zeta(iend+1,j,know)=ad_zeta(iend+1,j,know)+cff*adfac
832 ad_zeta(iend+1,j,kout)=0.0_r8
833 END IF
834 END DO
835 END IF
836!
837! Eastern edge, explicit Chapman boundary condition.
838!
839 ELSE IF (ad_lbc(ieast,isfsur,ng)%Chapman_explicit) THEN
840 DO j=jstr,jend
841 IF (lbc_apply(ng)%east(j)) THEN
842 cff=dt2d*grid(ng)%pm(iend,j)
843 cff1=sqrt(g*(grid(ng)%h(iend,j)+ &
844 & zeta(iend,j,know)))
845 cx=cff*cff1
846# ifdef MASKING
847!^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend+1,j,kout)* &
848!^ & GRID(ng)%rmask(Iend+1,j)
849!^
850 ad_zeta(iend+1,j,kout)=ad_zeta(iend+1,j,kout)* &
851 & grid(ng)%rmask(iend+1,j)
852# endif
853!^ tl_zeta(Iend+1,j,kout)=(1.0_r8-Cx)*tl_zeta(Iend+1,j,know)+&
854!^ & tl_Cx*(zeta(Iend+1,j,know)+ &
855!^ & zeta(Iend ,j,know))+ &
856!^ & Cx*tl_zeta(Iend,j,know)
857!^
858 ad_zeta(iend+1,j,know)=ad_zeta(iend+1,j,know)+ &
859 & (1.0_r8-cx)*ad_zeta(iend+1,j,kout)
860 ad_cx=ad_cx+(zeta(iend+1,j,know)+ &
861 & zeta(iend ,j,know))*ad_zeta(iend+1,j,kout)
862 ad_zeta(iend,j,know)=ad_zeta(iend,j,know)+ &
863 & cx*ad_zeta(iend+1,j,kout)
864 ad_zeta(iend+1,j,kout)=0.0_r8
865!^ tl_Cx=cff*tl_cff1
866!^
867 ad_cff1=ad_cff1+cff*ad_cx
868 ad_cx=0.0_r8
869!^ tl_cff1=0.5_r8*g*(GRID(ng)%tl_h(Iend,j)+ &
870!^ & tl_zeta(Iend,j,know))/cff1
871!^
872 adfac=0.5_r8*g*ad_cff1/cff1
873 grid(ng)%ad_h(iend,j)=grid(ng)%ad_h(iend,j)+adfac
874 ad_zeta(iend,j,know)=ad_zeta(iend,j,know)+adfac
875 ad_cff1=0.0_r8
876 END IF
877 END DO
878!
879! Eastern edge, implicit Chapman boundary condition.
880!
881 ELSE IF (ad_lbc(ieast,isfsur,ng)%Chapman_implicit) THEN
882 DO j=jstr,jend
883 IF (lbc_apply(ng)%east(j)) THEN
884 cff=dt2d*grid(ng)%pm(iend,j)
885 cff1=sqrt(g*(grid(ng)%h(iend,j)+ &
886 & zeta(iend,j,know)))
887 cx=cff*cff1
888 cff2=1.0_r8/(1.0_r8+cx)
889# ifdef MASKING
890!^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend+1,j,kout)* &
891!^ & GRID(ng)%rmask(Iend+1,j)
892!^
893 ad_zeta(iend+1,j,kout)=ad_zeta(iend+1,j,kout)* &
894 & grid(ng)%rmask(iend+1,j)
895# endif
896!^ tl_zeta(Iend+1,j,kout)=tl_cff2*(zeta(Iend+1,j,know)+ &
897!^ & Cx*zeta(Iend,j,kout))+ &
898!^ & cff2*(tl_zeta(Iend+1,j,know)+ &
899!^ & tl_Cx*zeta(Iend,j,kout)+ &
900!^ & Cx*tl_zeta(Iend,j,kout))
901!^
902 adfac=cff2*ad_zeta(iend+1,j,kout)
903 ad_zeta(iend ,j,kout)=ad_zeta(iend ,j,kout)+cx*adfac
904 ad_zeta(iend+1,j,know)=ad_zeta(iend+1,j,know)+adfac
905 ad_cx=ad_cx+zeta(iend,j,kout)*adfac
906 ad_cff2=ad_cff2+ &
907 & (zeta(iend+1,j,know)+ &
908 & cx*zeta(iend,j,kout))*ad_zeta(iend+1,j,kout)
909 ad_zeta(iend+1,j,kout)=0.0_r8
910!^ tl_cff2=-cff2*cff2*tl_Cx
911!^
912 ad_cx=ad_cx-cff2*cff2*ad_cff2
913 ad_cff2=0.0_r8
914!^ tl_Cx=cff*tl_cff1
915!^
916 ad_cff1=ad_cff1+cff*ad_cx
917 ad_cx=0.0_r8
918!^ tl_cff1=0.5_r8*g*(GRID(ng)%tl_h(Iend,j)+ &
919!^ & tl_zeta(Iend,j,know))/cff1
920!^
921 adfac=0.5_r8*g*ad_cff1/cff1
922 grid(ng)%ad_h(iend,j)=grid(ng)%ad_h(iend,j)+adfac
923 ad_zeta(iend,j,know)=ad_zeta(iend,j,know)+adfac
924 ad_cff1=0.0_r8
925 END IF
926 END DO
927!
928! Eastern edge, clamped boundary condition.
929!
930 ELSE IF (ad_lbc(ieast,isfsur,ng)%clamped) THEN
931 DO j=jstr,jend
932 IF (lbc_apply(ng)%east(j)) THEN
933# ifdef MASKING
934!^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend+1,j,kout)* &
935!^ & GRID(ng)%rmask(Iend+1,j)
936!^
937 ad_zeta(iend+1,j,kout)=ad_zeta(iend+1,j,kout)* &
938 & grid(ng)%rmask(iend+1,j)
939# endif
940# ifdef ADJUST_BOUNDARY
941 IF (lobc(ieast,isfsur,ng)) THEN
942!^ tl_zeta(Iend+1,j,kout)=BOUNDARY(ng)%tl_zeta_east(j)
943!^
944 boundary(ng)%ad_zeta_east(j)=boundary(ng)% &
945 & ad_zeta_east(j)+ &
946 & ad_zeta(iend+1,j,kout)
947 ad_zeta(iend+1,j,kout)=0.0_r8
948 ELSE
949!^ tl_zeta(Iend+1,j,kout)=0.0_r8
950!^
951 ad_zeta(iend+1,j,kout)=0.0_r8
952 END IF
953# else
954!^ tl_zeta(Iend+1,j,kout)=0.0_r8
955!^
956 ad_zeta(iend+1,j,kout)=0.0_r8
957# endif
958 END IF
959 END DO
960!
961! Eastern edge, gradient boundary condition.
962!
963 ELSE IF (ad_lbc(ieast,isfsur,ng)%gradient) THEN
964 DO j=jstr,jend
965 IF (lbc_apply(ng)%east(j)) THEN
966# ifdef MASKING
967!^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend+1,j,kout)* &
968!^ & GRID(ng)%rmask(Iend+1,j)
969!^
970 ad_zeta(iend+1,j,kout)=ad_zeta(iend+1,j,kout)* &
971 & grid(ng)%rmask(iend+1,j)
972# endif
973!^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend,j,kout)
974!^
975 ad_zeta(iend ,j,kout)=ad_zeta(iend ,j,kout)+ &
976 & ad_zeta(iend+1,j,kout)
977 ad_zeta(iend+1,j,kout)=0.0_r8
978 END IF
979 END DO
980!
981! Eastern edge, closed boundary condition.
982!
983 ELSE IF (ad_lbc(ieast,isfsur,ng)%closed) THEN
984 DO j=jstr,jend
985 IF (lbc_apply(ng)%east(j)) THEN
986# ifdef MASKING
987!^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend+1,j,kout)* &
988!^ & GRID(ng)%rmask(Iend+1,j)
989!^
990 ad_zeta(iend+1,j,kout)=ad_zeta(iend+1,j,kout)* &
991 & grid(ng)%rmask(iend+1,j)
992# endif
993!^ tl_zeta(Iend+1,j,kout)=tl_zeta(Iend,j,kout)
994!^
995 ad_zeta(iend ,j,kout)=ad_zeta(iend ,j,kout)+ &
996 & ad_zeta(iend+1,j,kout)
997 ad_zeta(iend+1,j,kout)=0.0_r8
998 END IF
999 END DO
1000 END IF
1001 END IF
1002!
1003!-----------------------------------------------------------------------
1004! Lateral boundary conditions at the western edge.
1005!-----------------------------------------------------------------------
1006!
1007 IF (domain(ng)%Western_Edge(tile)) THEN
1008!
1009! Western edge, implicit upstream radiation condition.
1010!
1011 IF (ad_lbc(iwest,isfsur,ng)%radiation) THEN
1012 IF (iic(ng).ne.0) THEN
1013 DO j=jstr,jend
1014 IF (lbc_apply(ng)%west(j)) THEN
1015# if defined CELERITY_READ && defined FORWARD_READ
1016 IF (ad_lbc(iwest,isfsur,ng)%nudging) THEN
1017 IF (boundary(ng)%zeta_west_Cx(j).eq.0.0_r8) THEN
1018 tau=fsobc_in(ng,iwest)
1019 ELSE
1020 tau=fsobc_out(ng,iwest)
1021 END IF
1022 tau=tau*dt2d
1023 END IF
1024 cx=boundary(ng)%zeta_west_Cx(j)
1025# ifdef RADIATION_2D
1026 ce=boundary(ng)%zeta_west_Ce(j)
1027# else
1028 ce=0.0_r8
1029# endif
1030 cff=boundary(ng)%zeta_west_C2(j)
1031# endif
1032# ifdef MASKING
1033!^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr-1,j,kout)* &
1034!^ & GRID(ng)%rmask(Istr-1,j)
1035!^
1036 ad_zeta(istr-1,j,kout)=ad_zeta(istr-1,j,kout)* &
1037 & grid(ng)%rmask(istr-1,j)
1038# endif
1039 IF (ad_lbc(iwest,isfsur,ng)%nudging) THEN
1040!^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr-1,j,kout)- &
1041!^ & tau*tl_zeta(Istr-1,j,know)
1042!^
1043 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)- &
1044 & tau*ad_zeta(istr-1,j,kout)
1045
1046 END IF
1047!^ tl_zeta(Istr-1,j,kout)=(cff*tl_zeta(Istr-1,j,know)+ &
1048!^ & Cx *tl_zeta(Istr ,j,kout)- &
1049!^ & MAX(Ce,0.0_r8)* &
1050!^ & tl_grad(Istr-1,j )- &
1051!^ & MIN(Ce,0.0_r8)* &
1052!^ & tl_grad(Istr-1,j+1))/ &
1053!^ & (cff+Cx)
1054!^
1055 adfac=ad_zeta(istr-1,j,kout)/(cff+cx)
1056 ad_grad(istr-1,j )=ad_grad(istr-1,j )- &
1057 & max(ce,0.0_r8)*adfac
1058 ad_grad(istr-1,j+1)=ad_grad(istr-1,j+1)- &
1059 & min(ce,0.0_r8)*adfac
1060 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)+cff*adfac
1061 ad_zeta(istr ,j,kout)=ad_zeta(istr ,j,kout)+cx *adfac
1062 ad_zeta(istr-1,j,kout)=0.0_r8
1063 END IF
1064 END DO
1065 END IF
1066!
1067! Western edge, explicit Chapman boundary condition.
1068!
1069 ELSE IF (ad_lbc(iwest,isfsur,ng)%Chapman_explicit) THEN
1070 DO j=jstr,jend
1071 IF (lbc_apply(ng)%west(j)) THEN
1072 cff=dt2d*grid(ng)%pm(istr,j)
1073 cff1=sqrt(g*(grid(ng)%h(istr,j)+ &
1074 & zeta(istr,j,know)))
1075 cx=cff*cff1
1076# ifdef MASKING
1077!^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr-1,j,kout)* &
1078!^ & GRID(ng)%rmask(Istr-1,j)
1079!^
1080 ad_zeta(istr-1,j,kout)=ad_zeta(istr-1,j,kout)* &
1081 & grid(ng)%rmask(istr-1,j)
1082# endif
1083!^ tl_zeta(Istr-1,j,kout)=(1.0_r8-Cx)*tl_zeta(Istr-1,j,know)+&
1084!^ & tl_Cx*(zeta(Istr-1,j,know)+ &
1085!^ & zeta(Istr ,j,know))+ &
1086!^ & Cx*tl_zeta(Istr,j,know)
1087!^
1088 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)+ &
1089 & (1.0_r8-cx)*ad_zeta(istr-1,j,kout)
1090 ad_cx=ad_cx+(zeta(istr-1,j,know)+ &
1091 & zeta(istr ,j,know))*ad_zeta(istr-1,j,kout)
1092 ad_zeta(istr,j,know)=ad_zeta(istr,j,know)+ &
1093 & cx*ad_zeta(istr-1,j,kout)
1094 ad_zeta(istr-1,j,kout)=0.0_r8
1095!^ tl_Cx=cff*tl_cff1
1096!^
1097 ad_cff1=ad_cff1+cff*ad_cx
1098 ad_cx=0.0_r8
1099!^ tl_cff1=0.5_r8*g*(GRID(ng)%tl_h(Istr,j)+ &
1100!^ & tl_zeta(Istr,j,know))/cff1
1101!^
1102 adfac=0.5_r8*g*ad_cff1/cff1
1103 grid(ng)%ad_h(istr,j)=grid(ng)%ad_h(istr,j)+adfac
1104 ad_zeta(istr,j,know)=ad_zeta(istr,j,know)+adfac
1105 ad_cff1=0.0_r8
1106 END IF
1107 END DO
1108!
1109! Western edge, implicit Chapman boundary condition.
1110!
1111 ELSE IF (ad_lbc(iwest,isfsur,ng)%Chapman_implicit) THEN
1112 DO j=jstr,jend
1113 IF (lbc_apply(ng)%west(j)) THEN
1114 cff=dt2d*grid(ng)%pm(istr,j)
1115 cff1=sqrt(g*(grid(ng)%h(istr,j)+ &
1116 & zeta(istr,j,know)))
1117 cx=cff*cff1
1118 cff2=1.0_r8/(1.0_r8+cx)
1119# ifdef MASKING
1120!^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr-1,j,kout)* &
1121!^ & GRID(ng)%rmask(Istr-1,j)
1122!^
1123 ad_zeta(istr-1,j,kout)=ad_zeta(istr-1,j,kout)* &
1124 & grid(ng)%rmask(istr-1,j)
1125# endif
1126!^ tl_zeta(Istr-1,j,kout)=tl_cff2*(zeta(Istr-1,j,know)+ &
1127!^ & Cx*zeta(Istr,j,kout))+ &
1128!^ & cff2*(tl_zeta(Istr-1,j,know)+ &
1129!^ & tl_Cx*zeta(Istr,j,kout)+ &
1130!^ & Cx*tl_zeta(Istr,j,kout))
1131!^
1132 adfac=cff2*ad_zeta(istr-1,j,kout)
1133 ad_zeta(istr-1,j,know)=ad_zeta(istr-1,j,know)+adfac
1134 ad_zeta(istr ,j,kout)=ad_zeta(istr ,j,kout)+cx*adfac
1135 ad_cx=ad_cx+zeta(istr,j,kout)*adfac
1136 ad_cff2=ad_cff2+ &
1137 & (zeta(istr-1,j,know)+ &
1138 & cx*zeta(istr,j,kout))*ad_zeta(istr-1,j,kout)
1139 ad_zeta(istr-1,j,kout)=0.0_r8
1140!^ tl_cff2=-cff2*cff2*tl_Cx
1141!^
1142 ad_cx=ad_cx-cff2*cff2*ad_cff2
1143 ad_cff2=0.0_r8
1144!^ tl_Cx=cff*tl_cff1
1145!^
1146 ad_cff1=ad_cff1+cff*ad_cx
1147 ad_cx=0.0_r8
1148!^ tl_cff1=0.5_r8*g*(GRID(ng)%tl_h(Istr,j)+ &
1149!^ & tl_zeta(Istr,j,know))/cff1
1150!^
1151 adfac=0.5_r8*g*ad_cff1/cff1
1152 grid(ng)%ad_h(istr,j)=grid(ng)%ad_h(istr,j)+adfac
1153 ad_zeta(istr,j,know)=ad_zeta(istr,j,know)+adfac
1154 ad_cff1=0.0_r8
1155 END IF
1156 END DO
1157!
1158! Western edge, clamped boundary condition.
1159!
1160 ELSE IF (ad_lbc(iwest,isfsur,ng)%clamped) THEN
1161 DO j=jstr,jend
1162 IF (lbc_apply(ng)%west(j)) THEN
1163# ifdef MASKING
1164!^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr-1,j,kout)* &
1165!^ & GRID(ng)%rmask(Istr-1,j)
1166!^
1167 ad_zeta(istr-1,j,kout)=ad_zeta(istr-1,j,kout)* &
1168 & grid(ng)%rmask(istr-1,j)
1169# endif
1170# ifdef ADJUST_BOUNDARY
1171 IF (lobc(iwest,isfsur,ng)) THEN
1172!^ tl_zeta(Istr-1,j,kout)=BOUNDARY(ng)%tl_zeta_west(j)
1173!^
1174 boundary(ng)%ad_zeta_west(j)=boundary(ng)% &
1175 & ad_zeta_west(j)+ &
1176 & ad_zeta(istr-1,j,kout)
1177 ad_zeta(istr-1,j,kout)=0.0_r8
1178 ELSE
1179!^ tl_zeta(Istr-1,j,kout)=0.0_r8
1180!^
1181 ad_zeta(istr-1,j,kout)=0.0_r8
1182 END IF
1183# else
1184!^ tl_zeta(Istr-1,j,kout)=0.0_r8
1185!^
1186 ad_zeta(istr-1,j,kout)=0.0_r8
1187# endif
1188 END IF
1189 END DO
1190!
1191! Western edge, gradient boundary condition.
1192!
1193 ELSE IF (ad_lbc(iwest,isfsur,ng)%gradient) THEN
1194 DO j=jstr,jend
1195 IF (lbc_apply(ng)%west(j)) THEN
1196# ifdef MASKING
1197!^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr-1,j,kout)* &
1198!^ & GRID(ng)%rmask(Istr-1,j)
1199!^
1200 ad_zeta(istr-1,j,kout)=ad_zeta(istr-1,j,kout)* &
1201 & grid(ng)%rmask(istr-1,j)
1202# endif
1203!^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr,j,kout)
1204!^
1205 ad_zeta(istr ,j,kout)=ad_zeta(istr ,j,kout)+ &
1206 & ad_zeta(istr-1,j,kout)
1207 ad_zeta(istr-1,j,kout)=0.0_r8
1208 END IF
1209 END DO
1210!
1211! Western edge, closed boundary condition.
1212!
1213 ELSE IF (ad_lbc(iwest,isfsur,ng)%closed) THEN
1214 DO j=jstr,jend
1215 IF (lbc_apply(ng)%west(j)) THEN
1216# ifdef MASKING
1217!^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr-1,j,kout)* &
1218!^ & GRID(ng)%rmask(Istr-1,j)
1219!^
1220 ad_zeta(istr-1,j,kout)=ad_zeta(istr-1,j,kout)* &
1221 & grid(ng)%rmask(istr-1,j)
1222# endif
1223!^ tl_zeta(Istr-1,j,kout)=tl_zeta(Istr,j,kout)
1224!^
1225 ad_zeta(istr ,j,kout)=ad_zeta(istr ,j,kout)+ &
1226 & ad_zeta(istr-1,j,kout)
1227 ad_zeta(istr-1,j,kout)=0.0_r8
1228 END IF
1229 END DO
1230 END IF
1231 END IF
1232
1233 RETURN
type(t_boundary), dimension(:), allocatable boundary
type(t_apply), dimension(:), allocatable lbc_apply
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
integer isfsur
type(t_lbc), dimension(:,:,:), allocatable ad_lbc
Definition mod_param.F:378
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, dimension(:), allocatable iic
real(r8), dimension(:), allocatable dcrit
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 fsobc_out
integer, parameter isouth
real(dp), dimension(:), allocatable dtfast
integer, parameter ieast
real(dp) g
integer, parameter inorth
real(dp), dimension(:,:), allocatable fsobc_in

References mod_param::ad_lbc, mod_boundary::boundary, mod_scalars::dcrit, mod_param::domain, mod_scalars::dtfast, mod_scalars::ewperiodic, mod_scalars::fsobc_in, mod_scalars::fsobc_out, mod_scalars::g, mod_grid::grid, mod_scalars::ieast, mod_scalars::iic, mod_scalars::inorth, mod_ncparam::isfsur, mod_scalars::isouth, mod_scalars::iwest, mod_boundary::lbc_apply, mod_scalars::lobc, mod_scalars::nsperiodic, and mod_scalars::predictor_2d_step.

Referenced by ad_ini_fields_mod::ad_ini_zeta_tile(), ad_ini_fields_mod::ad_out_zeta_tile(), ad_step2d_mod::ad_step2d_tile(), and ad_zetabc().

Here is the caller graph for this function: