[MITgcm-support] Orlanski and mass conservation

samar khatiwala spk at ldeo.columbia.edu
Mon Aug 2 10:41:16 EDT 2004


Hi Stefano

I may be wrong here, and its been ages since I worked with the
obcs package, but mass conservation was one of problems I
encountered when I first wrote this package. (Both Sonya Legg
and Jake Gebbie are more current with this package and may have
different experiences to report.) In my simulations,
there was a mean flow in the problem and this had to be accounted
for in computing OBEu etc. I don't know the details of your experiment,
but in my problem the OBCS were applied to the "wave" part of the
flow, i.e., deviations from a vertical average. To implement this,
in say orlanski_east, I have code like:

...
               OBEu(J,K,bi,bj) = ubarnew + uWave
C--------------Average uWave vertically---------
               xArea =
     &           _dyG(I_obc,J,bi,bj)*drF(K)*_hFacW(I_obc,J,K,bi,bj)
               uWaveAvg(J) = uWaveAvg(J) + uWave*xArea
               xAsum(J) = xAsum(J) + xArea
C-----------------------------------------------
...
C     Subtract vertical average of uWave
      DO J=1-Oly,sNy+Oly
         I_obc=OB_Ie(J,bi,bj)
         IF (I_obc.ne.0) THEN
            uWaveAvg(J)=uWaveAvg(J)/xAsum(J)
         ENDIF
      ENDDO
      DO K=1,Nr
         DO J=1-Oly,sNy+Oly
            I_obc=OB_Ie(J,bi,bj)
            IF (I_obc.ne.0) THEN
               OBEu(J,K,bi,bj) = OBEu(J,K,bi,bj) - uWaveAvg(J)
            ENDIF
         ENDDO
      ENDDO

Perhaps you may have to do something similar.

Samar

On Mon, 2 Aug 2004, Stefano Querin wrote:

> Hi,
> I am making some tests using the Orlanski OBCs. I recently observed that when a non-negligible velocity field reaches the OB, the mass conservation in the domain is no longer satisfied.
> For example, in the Plume on slope verification experiment, when the plume reaches the boundary, there is a volume increase in the basin (eta reaches 180 m in few days!). In this case, I simply made a longer run to see the evolution of the plume for some more days.
> I observed similar problems even when I used the flag: exactConserv=.TRUE.
> To be more clear I attach three pictures showing the problem:
> the first two show the evolution of eta and u-velocity (varying with depth) in time, the third is a sketch of the temperature profile when the plume "reaches" the OB (after 2.5 days; the horizontal lenghtscale is of tens of meters).
> My questions are:
> - (first of all!) did I make some mistakes?
> - is there a way to control the mass fluxes in the domain (in the orlanski_east(north, south, west).F routines there is no mention to eta)?
> - if I impose the Orlanski condition on the eastern boundary (as prescribed in data.obcs), do I (implicitly) impose the same condition on the western OB (u-velocity is non zero on that OB)?
> Thank you very much for any suggestion.
>
> Sincerely,
>
> Stefano Querin
>
> P.S.: (for Patrick) any news about the active OBCs?



More information about the MITgcm-support mailing list