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

Functions/Subroutines

subroutine, public rp_zetabc (ng, tile, kout)
 
subroutine, public rp_zetabc_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, krhs, kstp, kout, zeta, tl_zeta)
 

Function/Subroutine Documentation

◆ rp_zetabc()

subroutine, public rp_zetabc_mod::rp_zetabc ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) kout )

Definition at line 28 of file rp_zetabc.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 rp_zetabc_tile (ng, tile, &
44 & lbi, ubi, lbj, ubj, &
45 & imins, imaxs, jmins, jmaxs, &
46 & krhs(ng), kstp(ng), kout, &
47 & ocean(ng) % zeta, &
48 & ocean(ng) % tl_zeta)
49
50 RETURN
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, dimension(:), allocatable kstp
integer, dimension(:), allocatable krhs

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

Here is the call graph for this function:

◆ rp_zetabc_tile()

subroutine, public rp_zetabc_mod::rp_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) tl_zeta )

Definition at line 54 of file rp_zetabc.F.

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

Referenced by rp_ini_fields_mod::rp_ini_zeta_tile(), rp_step2d_mod::rp_step2d_tile(), and rp_zetabc().

Here is the caller graph for this function: