24 SUBROUTINE u2dbc (ng, tile, kout)
33 integer,
intent(in) :: ng, tile, kout
40 & lbi, ubi, lbj, ubj, &
41 & imins, imaxs, jmins, jmaxs, &
52 & LBi, UBi, LBj, UBj, &
53 & IminS, ImaxS, JminS, JmaxS, &
68 integer,
intent(in) :: ng, tile
69 integer,
intent(in) :: LBi, UBi, LBj, UBj
70 integer,
intent(in) :: IminS, ImaxS, JminS, JmaxS
71 integer,
intent(in) :: krhs, kstp, kout
74 real(r8),
intent(in) :: vbar(LBi:,LBj:,:)
75 real(r8),
intent(in) :: zeta(LBi:,LBj:,:)
77 real(r8),
intent(inout) :: ubar(LBi:,LBj:,:)
79 real(r8),
intent(in) :: vbar(LBi:UBi,LBj:UBj,:)
80 real(r8),
intent(in) :: zeta(LBi:UBi,LBj:UBj,:)
82 real(r8),
intent(inout) :: ubar(LBi:UBi,LBj:UBj,:)
90 real(r8),
parameter :: eps = 1.0e-20_r8
93 real(r8) :: bry_pgr, bry_cor, bry_str, bry_val
94 real(r8) :: cff, cff1, cff2, dt2d, dUde, dUdt, dUdx
95 real(r8) :: obc_in, obc_out, tau
97 real(r8),
dimension(IminS:ImaxS,JminS:JmaxS) :: grad
99#include "set_bounds.h"
105 IF (first_2d_step)
THEN
120 IF (
domain(ng)%Western_Edge(tile))
THEN
126 grad(istr ,j)=ubar(istr ,j ,know)- &
127 & ubar(istr ,j-1,know)
128 grad(istr+1,j)=ubar(istr+1,j ,know)- &
129 & ubar(istr+1,j-1,know)
133 dudt=ubar(istr+1,j,know)-ubar(istr+1,j,kout)
134 dudx=ubar(istr+1,j,know)-ubar(istr+2,j,know)
139 & (
clima(ng)%M2nudgcof(istr-1,j)+ &
140 &
clima(ng)%M2nudgcof(istr ,j))
141 obc_in =
obcfac(ng)*obc_out
146 IF ((dudt*dudx).lt.0.0_r8)
THEN
154 IF ((dudt*dudx).lt.0.0_r8) dudt=0.0_r8
155 IF ((dudt*(grad(istr+1,j )+ &
156 & grad(istr+1,j+1))).gt.0.0_r8)
THEN
159 dude=grad(istr+1,j+1)
161 cff=dudt/max(dudx*dudx+dude*dude,eps)
162 cx=min(1.0_r8,cff*dudx)
164 ce=min(1.0_r8,max(-1.0_r8,cff*dude))
168#if defined CELERITY_WRITE && defined FORWARD_WRITE
173 ubar(istr,j,kout)=(1.0_r8-cx)*ubar(istr,j,know)+ &
174 & cx*ubar(istr+1,j,know)- &
175 & max(ce,0.0_r8)*grad(istr,j )- &
176 & min(ce,0.0_r8)*grad(istr,j+1)
179 ubar(istr,j,kout)=ubar(istr,j,kout)+ &
184 ubar(istr,j,kout)=ubar(istr,j,kout)* &
185 &
grid(ng)%umask(istr,j)
195#if defined SSH_TIDES && !defined UV_TIDES
197 bry_pgr=-
g*(zeta(istr,j,know)- &
199 & 0.5_r8*
grid(ng)%pm(istr,j)
201 bry_pgr=-
g*(zeta(istr ,j,know)- &
202 & zeta(istr-1,j,know))* &
203 & 0.5_r8*(
grid(ng)%pm(istr-1,j)+ &
204 &
grid(ng)%pm(istr ,j))
207 bry_cor=0.125_r8*(vbar(istr-1,j ,know)+ &
208 & vbar(istr-1,j+1,know)+ &
209 & vbar(istr ,j ,know)+ &
210 & vbar(istr ,j+1,know))* &
211 & (
grid(ng)%f(istr-1,j)+ &
212 &
grid(ng)%f(istr ,j))
216 cff1=1.0_r8/(0.5_r8*(
grid(ng)%h(istr-1,j)+ &
217 & zeta(istr-1,j,know)+ &
218 &
grid(ng)%h(istr ,j)+ &
219 & zeta(istr ,j,know)))
220 bry_str=cff1*(
forces(ng)%sustr(istr,j)- &
221 &
forces(ng)%bustr(istr,j))
222 cx=1.0_r8/sqrt(
g*0.5_r8*(
grid(ng)%h(istr-1,j)+ &
223 & zeta(istr-1,j,know)+ &
224 &
grid(ng)%h(istr ,j)+ &
225 & zeta(istr ,j,know)))
226 cff2=
grid(ng)%om_u(istr,j)*cx
228 bry_val=ubar(istr+1,j,know)+ &
235 cff=1.0_r8/(0.5_r8*(
grid(ng)%h(istr-1,j)+ &
236 & zeta(istr-1,j,know)+ &
237 &
grid(ng)%h(istr ,j)+ &
238 & zeta(istr ,j,know)))
240 ubar(istr,j,kout)=bry_val- &
241 & cx*(0.5_r8*(zeta(istr-1,j,know)+ &
242 & zeta(istr ,j,know))- &
245 ubar(istr,j,kout)=ubar(istr,j,kout)* &
246 &
grid(ng)%umask(istr,j)
256 ubar(istr,j,kout)=
boundary(ng)%ubar_west(j)
258 ubar(istr,j,kout)=ubar(istr,j,kout)* &
259 &
grid(ng)%umask(istr,j)
269 ubar(istr,j,kout)=ubar(istr+1,j,kout)
271 ubar(istr,j,kout)=ubar(istr,j,kout)* &
272 &
grid(ng)%umask(istr,j)
283 bry_pgr=-
g*(zeta(istr,j,know)- &
285 & 0.5_r8*
grid(ng)%pm(istr,j)
287 bry_pgr=-
g*(zeta(istr ,j,know)- &
288 & zeta(istr-1,j,know))* &
289 & 0.5_r8*(
grid(ng)%pm(istr-1,j)+ &
290 &
grid(ng)%pm(istr ,j))
293 bry_cor=0.125_r8*(vbar(istr-1,j ,know)+ &
294 & vbar(istr-1,j+1,know)+ &
295 & vbar(istr ,j ,know)+ &
296 & vbar(istr ,j+1,know))* &
297 & (
grid(ng)%f(istr-1,j)+ &
298 &
grid(ng)%f(istr ,j))
302 cff=1.0_r8/(0.5_r8*(
grid(ng)%h(istr-1,j)+ &
303 & zeta(istr-1,j,know)+ &
304 &
grid(ng)%h(istr ,j)+ &
305 & zeta(istr ,j,know)))
306 bry_str=cff*(
forces(ng)%sustr(istr,j)- &
307 &
forces(ng)%bustr(istr,j))
308 ubar(istr,j,kout)=ubar(istr,j,know)+ &
313 ubar(istr,j,kout)=ubar(istr,j,kout)* &
314 &
grid(ng)%umask(istr,j)
324 ubar(istr,j,kout)=0.0_r8
334 IF (
domain(ng)%Eastern_Edge(tile))
THEN
340 grad(iend ,j)=ubar(iend ,j ,know)- &
341 & ubar(iend ,j-1,know)
342 grad(iend+1,j)=ubar(iend+1,j ,know)- &
343 & ubar(iend+1,j-1,know)
347 dudt=ubar(iend,j,know)-ubar(iend ,j,kout)
348 dudx=ubar(iend,j,know)-ubar(iend-1,j,know)
353 & (
clima(ng)%M2nudgcof(iend ,j)+ &
354 &
clima(ng)%M2nudgcof(iend+1,j))
355 obc_in =
obcfac(ng)*obc_out
360 IF ((dudt*dudx).lt.0.0_r8)
THEN
368 IF ((dudt*dudx).lt.0.0_r8) dudt=0.0_r8
369 IF ((dudt*(grad(iend,j )+ &
370 & grad(iend,j+1))).gt.0.0_r8)
THEN
375 cff=dudt/max(dudx*dudx+dude*dude,eps)
376 cx=min(1.0_r8,cff*dudx)
378 ce=min(1.0_r8,max(-1.0_r8,cff*dude))
382#if defined CELERITY_WRITE && defined FORWARD_WRITE
387 ubar(iend+1,j,kout)=(1.0_r8-cx)*ubar(iend+1,j,know)+ &
388 & cx*ubar(iend,j,know)- &
389 & max(ce,0.0_r8)*grad(iend+1,j )- &
390 & min(ce,0.0_r8)*grad(iend+1,j+1)
393 ubar(iend+1,j,kout)=ubar(iend+1,j,kout)+ &
395 & ubar(iend+1,j,know))
398 ubar(iend+1,j,kout)=ubar(iend+1,j,kout)* &
399 &
grid(ng)%umask(iend+1,j)
409#if defined SSH_TIDES && !defined UV_TIDES
412 & zeta(iend,j,know))* &
413 & 0.5_r8*
grid(ng)%pm(iend,j)
415 bry_pgr=-
g*(zeta(iend+1,j,know)- &
416 & zeta(iend ,j,know))* &
417 & 0.5_r8*(
grid(ng)%pm(iend ,j)+ &
418 &
grid(ng)%pm(iend+1,j))
421 bry_cor=0.125_r8*(vbar(iend ,j ,know)+ &
422 & vbar(iend ,j+1,know)+ &
423 & vbar(iend+1,j ,know)+ &
424 & vbar(iend+1,j+1,know))* &
425 & (
grid(ng)%f(iend ,j)+ &
426 &
grid(ng)%f(iend+1,j))
430 cff1=1.0_r8/(0.5_r8*(
grid(ng)%h(iend ,j)+ &
431 & zeta(iend ,j,know)+ &
432 &
grid(ng)%h(iend+1,j)+ &
433 & zeta(iend+1,j,know)))
434 bry_str=cff1*(
forces(ng)%sustr(iend+1,j)- &
435 &
forces(ng)%bustr(iend+1,j))
436 cx=1.0_r8/sqrt(
g*0.5_r8*(
grid(ng)%h(iend+1,j)+ &
437 & zeta(iend+1,j,know)+ &
438 &
grid(ng)%h(iend ,j)+ &
439 & zeta(iend ,j,know)))
440 cff2=
grid(ng)%om_u(iend+1,j)*cx
442 bry_val=ubar(iend,j,know)+ &
449 cff=1.0_r8/(0.5_r8*(
grid(ng)%h(iend ,j)+ &
450 & zeta(iend ,j,know)+ &
451 &
grid(ng)%h(iend+1,j)+ &
452 & zeta(iend+1,j,know)))
454 ubar(iend+1,j,kout)=bry_val+ &
455 & cx*(0.5_r8*(zeta(iend ,j,know)+ &
456 & zeta(iend+1,j,know))- &
459 ubar(iend+1,j,kout)=ubar(iend+1,j,kout)* &
460 &
grid(ng)%umask(iend+1,j)
470 ubar(iend+1,j,kout)=
boundary(ng)%ubar_east(j)
472 ubar(iend+1,j,kout)=ubar(iend+1,j,kout)* &
473 &
grid(ng)%umask(iend+1,j)
483 ubar(iend+1,j,kout)=ubar(iend,j,kout)
485 ubar(iend+1,j,kout)=ubar(iend+1,j,kout)* &
486 &
grid(ng)%umask(iend+1,j)
498 & zeta(iend,j,know))* &
499 & 0.5_r8*
grid(ng)%pm(iend,j)
501 bry_pgr=-
g*(zeta(iend+1,j,know)- &
502 & zeta(iend ,j,know))* &
503 & 0.5_r8*(
grid(ng)%pm(iend ,j)+ &
504 &
grid(ng)%pm(iend+1,j))
507 bry_cor=0.125_r8*(vbar(iend ,j ,know)+ &
508 & vbar(iend ,j+1,know)+ &
509 & vbar(iend+1,j ,know)+ &
510 & vbar(iend+1,j+1,know))* &
511 & (
grid(ng)%f(iend ,j)+ &
512 &
grid(ng)%f(iend+1,j))
516 cff=1.0_r8/(0.5_r8*(
grid(ng)%h(iend ,j)+ &
517 & zeta(iend ,j,know)+ &
518 &
grid(ng)%h(iend+1,j)+ &
519 & zeta(iend+1,j,know)))
520 bry_str=cff*(
forces(ng)%sustr(iend+1,j)- &
521 &
forces(ng)%bustr(iend+1,j))
522 ubar(iend+1,j,kout)=ubar(iend+1,j,know)+ &
527 ubar(iend+1,j,kout)=ubar(iend+1,j,kout)* &
528 &
grid(ng)%umask(iend+1,j)
538 ubar(iend+1,j,kout)=0.0_r8
548 IF (
domain(ng)%Southern_Edge(tile))
THEN
554 grad(i,jstr-1)=ubar(i+1,jstr-1,know)- &
555 & ubar(i ,jstr-1,know)
556 grad(i,jstr )=ubar(i+1,jstr ,know)- &
557 & ubar(i ,jstr ,know)
561 dudt=ubar(i,jstr,know)-ubar(i,jstr ,kout)
562 dude=ubar(i,jstr,know)-ubar(i,jstr+1,know)
567 & (
clima(ng)%M2nudgcof(i-1,jstr-1)+ &
568 &
clima(ng)%M2nudgcof(i ,jstr-1))
569 obc_in =
obcfac(ng)*obc_out
574 IF ((dudt*dude).lt.0.0_r8)
THEN
582 IF ((dudt*dude).lt.0.0_r8) dudt=0.0_r8
583 IF ((dudt*(grad(i-1,jstr)+ &
584 & grad(i ,jstr))).gt.0.0_r8)
THEN
589 cff=dudt/max(dudx*dudx+dude*dude,eps)
591 cx=min(1.0_r8,max(-1.0_r8,cff*dudx))
595 ce=min(1.0_r8,cff*dude)
596#if defined CELERITY_WRITE && defined FORWARD_WRITE
601 ubar(i,jstr-1,kout)=(1.0_r8-ce)*ubar(i,jstr-1,know)+ &
602 & ce*ubar(i,jstr,know)- &
603 & max(cx,0.0_r8)*grad(i-1,jstr-1)- &
604 & min(cx,0.0_r8)*grad(i ,jstr-1)
607 ubar(i,jstr-1,kout)=ubar(i,jstr-1,kout)+ &
608 & tau*(
boundary(ng)%ubar_south(i)- &
609 & ubar(i,jstr-1,know))
612 ubar(i,jstr-1,kout)=ubar(i,jstr-1,kout)* &
613 &
grid(ng)%umask(i,jstr-1)
624 cff=dt2d*0.5_r8*(
grid(ng)%pn(i-1,jstr)+ &
625 &
grid(ng)%pn(i ,jstr))
626 cff1=sqrt(
g*0.5_r8*(
grid(ng)%h(i-1,jstr)+ &
627 & zeta(i-1,jstr,know)+ &
628 &
grid(ng)%h(i ,jstr)+ &
629 & zeta(i ,jstr,know)))
631 cff2=1.0_r8/(1.0_r8+ce)
632 ubar(i,jstr-1,kout)=cff2*(ubar(i,jstr-1,know)+ &
633 & ce*ubar(i,jstr,kout))
635 ubar(i,jstr-1,kout)=ubar(i,jstr-1,kout)* &
636 &
grid(ng)%umask(i,jstr-1)
646 ubar(i,jstr-1,kout)=
boundary(ng)%ubar_south(i)
648 ubar(i,jstr-1,kout)=ubar(i,jstr-1,kout)* &
649 &
grid(ng)%umask(i,jstr-1)
659 ubar(i,jstr-1,kout)=ubar(i,jstr,kout)
661 ubar(i,jstr-1,kout)=ubar(i,jstr-1,kout)* &
662 &
grid(ng)%umask(i,jstr-1)
680 ubar(i,jstr-1,kout)=
gamma2(ng)*ubar(i,jstr,kout)
682 ubar(i,jstr-1,kout)=ubar(i,jstr-1,kout)* &
683 &
grid(ng)%umask(i,jstr-1)
694 IF (
domain(ng)%Northern_Edge(tile))
THEN
700 grad(i,jend )=ubar(i+1,jend ,know)- &
701 & ubar(i ,jend ,know)
702 grad(i,jend+1)=ubar(i+1,jend+1,know)- &
703 & ubar(i ,jend+1,know)
707 dudt=ubar(i,jend,know)-ubar(i,jend ,kout)
708 dude=ubar(i,jend,know)-ubar(i,jend-1,know)
713 & (
clima(ng)%M2nudgcof(i-1,jend+1)+ &
714 &
clima(ng)%M2nudgcof(i ,jend+1))
715 obc_in =
obcfac(ng)*obc_out
720 IF ((dudt*dude).lt.0.0_r8)
THEN
728 IF ((dudt*dude).lt.0.0_r8) dudt=0.0_r8
729 IF ((dudt*(grad(i-1,jend)+ &
730 & grad(i ,jend))).gt.0.0_r8)
THEN
735 cff=dudt/max(dudx*dudx+dude*dude,eps)
737 cx=min(1.0_r8,max(-1.0_r8,cff*dudx))
741 ce=min(1.0_r8,cff*dude)
742#if defined CELERITY_WRITE && defined FORWARD_WRITE
747 ubar(i,jend+1,kout)=(1.0_r8-ce)*ubar(i,jend+1,know)+ &
748 & ce*ubar(i,jend,know)- &
749 & max(cx,0.0_r8)*grad(i-1,jend+1)- &
750 & min(cx,0.0_r8)*grad(i ,jend+1)
753 ubar(i,jend+1,kout)=ubar(i,jend+1,kout)+ &
754 & tau*(
boundary(ng)%ubar_north(i)- &
755 & ubar(i,jend+1,know))
758 ubar(i,jend+1,kout)=ubar(i,jend+1,kout)* &
759 &
grid(ng)%umask(i,jend+1)
770 cff=dt2d*0.5_r8*(
grid(ng)%pn(i-1,jend)+ &
771 &
grid(ng)%pn(i ,jend))
772 cff1=sqrt(
g*0.5_r8*(
grid(ng)%h(i-1,jend)+ &
773 & zeta(i-1,jend,know)+ &
774 &
grid(ng)%h(i ,jend)+ &
775 & zeta(i ,jend,know)))
777 cff2=1.0_r8/(1.0_r8+ce)
778 ubar(i,jend+1,kout)=cff2*(ubar(i,jend+1,know)+ &
779 & ce*ubar(i,jend,kout))
781 ubar(i,jend+1,kout)=ubar(i,jend+1,kout)* &
782 &
grid(ng)%umask(i,jend+1)
792 ubar(i,jend+1,kout)=
boundary(ng)%ubar_north(i)
794 ubar(i,jend+1,kout)=ubar(i,jend+1,kout)* &
795 &
grid(ng)%umask(i,jend+1)
805 ubar(i,jend+1,kout)=ubar(i,jend,kout)
807 ubar(i,jend+1,kout)=ubar(i,jend+1,kout)* &
808 &
grid(ng)%umask(i,jend+1)
826 ubar(i,jend+1,kout)=
gamma2(ng)*ubar(i,jend,kout)
828 ubar(i,jend+1,kout)=ubar(i,jend+1,kout)* &
829 &
grid(ng)%umask(i,jend+1)
841 IF (
domain(ng)%SouthWest_Corner(tile))
THEN
844 ubar(istr,jstr-1,kout)=0.5_r8*(ubar(istr+1,jstr-1,kout)+ &
845 & ubar(istr ,jstr ,kout))
848 IF (
domain(ng)%SouthEast_Corner(tile))
THEN
851 ubar(iend+1,jstr-1,kout)=0.5_r8*(ubar(iend ,jstr-1,kout)+ &
852 & ubar(iend+1,jstr ,kout))
855 IF (
domain(ng)%NorthWest_Corner(tile))
THEN
858 ubar(istr,jend+1,kout)=0.5_r8*(ubar(istr ,jend ,kout)+ &
859 & ubar(istr+1,jend+1,kout))
862 IF (
domain(ng)%NorthEast_Corner(tile))
THEN
865 ubar(iend+1,jend+1,kout)=0.5_r8*(ubar(iend+1,jend ,kout)+ &
866 & ubar(iend ,jend+1,kout))
878 IF (
domain(ng)%Western_Edge(tile))
THEN
881 cff1=abs(abs(
grid(ng)%umask_wet(istr,j))-1.0_r8)
882 cff2=0.5_r8+dsign(0.5_r8,ubar(istr,j,kout))* &
883 &
grid(ng)%umask_wet(istr,j)
884 cff=0.5_r8*
grid(ng)%umask_wet(istr,j)*cff1+ &
886 ubar(istr,j,kout)=ubar(istr,j,kout)*cff
890 IF (
domain(ng)%Eastern_Edge(tile))
THEN
893 cff1=abs(abs(
grid(ng)%umask_wet(iend+1,j))-1.0_r8)
894 cff2=0.5_r8+dsign(0.5_r8,ubar(iend+1,j,kout))* &
895 &
grid(ng)%umask_wet(iend+1,j)
896 cff=0.5_r8*
grid(ng)%umask_wet(iend+1,j)*cff1+ &
898 ubar(iend+1,j,kout)=ubar(iend+1,j,kout)*cff
905 IF (
domain(ng)%Southern_Edge(tile))
THEN
908 cff1=abs(abs(
grid(ng)%umask_wet(i,jstr-1))-1.0_r8)
909 cff2=0.5_r8+dsign(0.5_r8,ubar(i,jstr-1,kout))* &
910 &
grid(ng)%umask_wet(i,jstr-1)
911 cff=0.5_r8*
grid(ng)%umask_wet(i,jstr-1)*cff1+ &
913 ubar(i,jstr-1,kout)=ubar(i,jstr-1,kout)*cff
917 IF (
domain(ng)%Northern_Edge(tile))
THEN
920 cff1=abs(abs(
grid(ng)%umask_wet(i,jend+1))-1.0_r8)
921 cff2=0.5_r8+dsign(0.5_r8,ubar(i,jend+1,kout))* &
922 &
grid(ng)%umask_wet(i,jend+1)
923 cff=0.5_r8*
grid(ng)%umask_wet(i,jend+1)*cff1+ &
925 ubar(i,jend+1,kout)=ubar(i,jend+1,kout)*cff
932 IF (
domain(ng)%SouthWest_Corner(tile))
THEN
935 cff1=abs(abs(
grid(ng)%umask_wet(istr,jstr-1))-1.0_r8)
936 cff2=0.5_r8+dsign(0.5_r8,ubar(istr,jstr-1,kout))* &
937 &
grid(ng)%umask_wet(istr,jstr-1)
938 cff=0.5_r8*
grid(ng)%umask_wet(istr,jstr-1)*cff1+ &
940 ubar(istr,jstr-1,kout)=ubar(istr,jstr-1,kout)*cff
943 IF (
domain(ng)%SouthEast_Corner(tile))
THEN
946 cff1=abs(abs(
grid(ng)%umask_wet(iend+1,jstr-1))-1.0_r8)
947 cff2=0.5_r8+dsign(0.5_r8,ubar(iend+1,jstr-1,kout))* &
948 &
grid(ng)%umask_wet(iend+1,jstr-1)
949 cff=0.5_r8*
grid(ng)%umask_wet(iend+1,jstr-1)*cff1+ &
951 ubar(iend+1,jstr-1,kout)=ubar(iend+1,jstr-1,kout)*cff
954 IF (
domain(ng)%NorthWest_Corner(tile))
THEN
957 cff1=abs(abs(
grid(ng)%umask_wet(istr,jend+1))-1.0_r8)
958 cff2=0.5_r8+dsign(0.5_r8,ubar(istr,jend+1,kout))* &
959 &
grid(ng)%umask_wet(istr,jend+1)
960 cff=0.5_r8*
grid(ng)%umask_wet(istr,jend+1)*cff1+ &
962 ubar(istr,jend+1,kout)=ubar(istr,jend+1,kout)*cff
965 IF (
domain(ng)%NorthEast_Corner(tile))
THEN
968 cff1=abs(abs(
grid(ng)%umask_wet(iend+1,jend+1))-1.0_r8)
969 cff2=0.5_r8+dsign(0.5_r8,ubar(iend+1,jend+1,kout))* &
970 &
grid(ng)%umask_wet(iend+1,jend+1)
971 cff=0.5_r8*
grid(ng)%umask_wet(iend+1,jend+1)+cff1+ &
973 ubar(iend+1,jend+1,kout)=ubar(iend+1,jend+1,kout)*cff
type(t_boundary), dimension(:), allocatable boundary
type(t_apply), dimension(:), allocatable lbc_apply
type(t_clima), dimension(:), allocatable clima
type(t_forces), dimension(:), allocatable forces
type(t_grid), dimension(:), allocatable grid
type(t_ocean), dimension(:), allocatable ocean
type(t_lbc), dimension(:,:,:), allocatable lbc
type(t_domain), dimension(:), allocatable domain
logical, dimension(:), allocatable lnudgem2clm
logical, dimension(:), allocatable ewperiodic
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 inorth
real(dp), dimension(:,:), allocatable m2obc_in
integer, dimension(:), allocatable kstp
integer, dimension(:), allocatable krhs
subroutine, public u2dbc_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, krhs, kstp, kout, ubar, vbar, zeta)
subroutine, public u2dbc(ng, tile, kout)