WET_DRY and boundary conditions

Bug reports, work arounds and fixes

Moderators: arango, robertson

Post Reply
Message
Author
User avatar
kate
Posts: 4091
Joined: Wed Jul 02, 2003 5:29 pm
Location: CFOS/UAF, USA

WET_DRY and boundary conditions

#1 Unread post by kate »

I made this change to my code for negative values of h at the boundary. The sqrt(g*h) was not behaving well.

Code: Select all

diff --git a/ROMS/Nonlinear/u2dbc_im.F b/ROMS/Nonlinear/u2dbc_im.F
index ef55e36..eb3921d 100644
--- a/ROMS/Nonlinear/u2dbc_im.F
+++ b/ROMS/Nonlinear/u2dbc_im.F
@@ -311,7 +311,9 @@
               bry_val=BOUNDARY(ng)%ubar_west(j)
 #endif
               cff=1.0_r8/(0.5_r8*(GRID(ng)%h(Istr-1,j)+                 &
-     &                            GRID(ng)%h(Istr  ,j)))
+     &                            zeta(Istr-1,j,know)+                  &
+     &                            GRID(ng)%h(Istr  ,j)+                 &
+     &                            zeta(Istr  ,j,know)))
               cff1=SQRT(g*cff)
               Cx=dt2d*cff1*0.5_r8*(GRID(ng)%pm(Istr-1,j)+               &
      &                             GRID(ng)%pm(Istr  ,j))
@@ -624,7 +626,9 @@
               bry_val=BOUNDARY(ng)%ubar_east(j)
 #endif
               cff=1.0_r8/(0.5_r8*(GRID(ng)%h(Iend  ,j)+                 &
-     &                            GRID(ng)%h(Iend+1,j)))
+     &                            zeta(Iend  ,j,know)+                  &
+     &                            GRID(ng)%h(Iend+1,j)+                 &
+     &                            zeta(Iend+1,j,know)))
               cff1=SQRT(g*cff)
               Cx=dt2d*cff1*0.5_r8*(GRID(ng)%pm(Iend  ,j)+               &
      &                             GRID(ng)%pm(Iend+1,j))
diff --git a/ROMS/Nonlinear/v2dbc_im.F b/ROMS/Nonlinear/v2dbc_im.F
index 68b33a2..e4ee7a5 100644
--- a/ROMS/Nonlinear/v2dbc_im.F
+++ b/ROMS/Nonlinear/v2dbc_im.F
@@ -312,7 +312,9 @@
               bry_val=BOUNDARY(ng)%vbar_south(i)
 #endif
               cff=1.0_r8/(0.5_r8*(GRID(ng)%h(i,Jstr-1)+                 &
-     &                            GRID(ng)%h(i,Jstr  )))
+     &                            zeta(i,Jstr-1,know)+                  &
+     &                            GRID(ng)%h(i,Jstr  )+                 &
+     &                            zeta(i,Jstr  ,know)))
               cff1=SQRT(g*cff)
               Ce=dt2d*cff1*0.5_r8*(GRID(ng)%pn(i,Jstr-1)+               &
      &                             GRID(ng)%pn(i,Jstr  ))
@@ -625,7 +627,9 @@
               bry_val=BOUNDARY(ng)%vbar_north(i)
 #endif
               cff=1.0_r8/(0.5_r8*(GRID(ng)%h(i,Jend  )+                 &
-     &                            GRID(ng)%h(i,Jend+1)))
+     &                            zeta(i,Jend  ,know)+                  &
+     &                            GRID(ng)%h(i,Jend+1)+                 &
+     &                            zeta(i,Jend+1,know)))
               cff1=SQRT(g*cff)
               Ce=dt2d*cff1*0.5_r8*(GRID(ng)%pn(i,Jend  )+               &
      &                             GRID(ng)%pn(i,Jend+1))
Feel free to wrap it in #ifdef WET_DRY if you think it needs it.

jcwarner
Posts: 1200
Joined: Wed Dec 31, 2003 6:16 pm
Location: USGS, USA

Re: WET_DRY and boundary conditions

#2 Unread post by jcwarner »

there should be a call in
ROMS/Nonlinear/set_data
that compares the water levels to make sure they are above the h.

Oooh - is this for the new Shchepetkin BC?
-j

User avatar
kate
Posts: 4091
Joined: Wed Jul 02, 2003 5:29 pm
Location: CFOS/UAF, USA

Re: WET_DRY and boundary conditions

#3 Unread post by kate »

Yes, it's the new Shchepetkin boundary condition.

Post Reply