[MITgcm-support] SSH drift and obcs again

Martin Losch mlosch at awi-bremerhaven.de
Tue Oct 25 04:18:35 EDT 2005


Hi Tom,

unfortunately I am not really the authority on the obcs-pkg (and as far 
as I remember from a few discussions at MIT, there isn't really anyone 
who feels particularily responsible for all the details), so my 
comments on this topic may or may not be valuable. A few things I think 
I can explain, though:

obcs_apply_uv and solve_for pressure (the hydrostatic part only):
the normal flux across a boundary is specified at velocity points 
(v-points at the southern boundary), but cg2d_b/x is defined at 
pressure/tracer points in solve_for_pressure. So at a southern open 
boundary at Ny=1, cg2d_b/x(iy=1) = 0, and v(iy=1)=v(iy=2)=obsv, so that 
the flux across the boundary does not have a gradient along the 
direction of the flux (which would also show up in cg2d_b). I makes 
sense to use the velocity mask hFacS on inside of the boundary to make 
sure that topography across the boundary does not destroy this 
condition (I guess you shouldn't have a topography gradient across the 
boundary for other reasons, but here the problem is explicitly 
avoided).

I don't understand, where the drift in Eta comes from. I have a box 
with four open boundaries (no topography at all), and I observe drift 
in mean-ssh myself, but I don't worry too much about it, because I only 
want to integrate for 40days and I can actually attribute the drift to 
a flux divergence over my domain of my horizontal boundary velocities. 
So, I am not much help there. Maybe I should look into that problem a 
little more myself because the drift is significant. However, I can 
alway substract a mean value  and then SSH is OK. I may be wrong about 
this, but the the matrix system that cg2d solve is subject to Neumann 
conditions, so that is determined up to a constant, which is the mean 
SSH. This mean is usually zero, but small (numerical) inconsistencies 
(e.g., overall flux divergence) may make this mean non-zero. What 
happens if you substract a mean from the SSH after each timestep, will 
the model still explode?

You will still have corner points, even if you move the open boundary 
one grid-point inwards. Now your corner points will be at (2,2), 
(2,Ny-1), (Nx-1,2), (Nx-1,Ny-1). That shouldn't change anything about 
the numerics (cleverly, the MITgcm has overlap regions beyond Nx,Ny 
(o:).

Now, I don't know how and if all this works with a nonlinear free 
surface or non-hydrostatic code ...

Martin

On Oct 24, 2005, at 8:54 PM, Thomas Haine wrote:

> Hi Martin, Chris,
>
> OK, I understand that bcs revert to periodic when they're unspecified
> and there's no wall in the bathy file. But then why does my case i)
> work?  E/W obcs are unspecified here and there is an open periodic
> boundary in the east/west direction too.
>
> Also, when I set the corner obc speeds as you suggest:
>
> a) Case i) before still works the same. No surprise.
>
> b) With all 4 obcs open (the case I'm really interested in): the source
> term is O(1) again and eta is drifting but it integrates for many steps
> OK. Why is there a cg2d source?
>
> If I remove the corner obc points altogether in data.obcs (i.e. obcs
> from 2 -> NX-1, not 1 -> NX etc) I get similar results to b) except it
> dies after 12 steps (w overflow).
>
> More generally, I don't understand how obcs_apply_uv.F is supposed to
> work. Why is the mask shifted relative to the grid in the 2nd obc
> assignment, e.g.:
>
>         vFld(I,OB_Js(I,bi,bj)+1,K,bi,bj)=OBSv(I,K,bi,bj)
>      &  *_maskS(I,OB_Js(I,bi,bj)+1,K,bi,bj)
>
>         vFld(I,OB_Js(I,bi,bj),K,bi,bj)=OBSv(I,K,bi,bj)
>      &  *_maskS(I,OB_Js(I,bi,bj)+1,K,bi,bj)
>
> Is this consistent with solve_for_pressure.F where cg2d_b and cg2d_x 
> are
> set at the obcs points themselves, not 1 step in on the S/W bdy, or at
> multiple points as in obcs_apply_uv?  I can't see how.
>
> Tom.
>
> _______________________________________________
> MITgcm-support mailing list
> MITgcm-support at mitgcm.org
> http://mitgcm.org/mailman/listinfo/mitgcm-support




More information about the MITgcm-support mailing list