[MITgcm-support] Anamalous Decrease in Sea Surface Level and SPONGE Layer

Himansu Pradhan oceancalling at gmail.com
Wed Sep 26 09:53:39 EDT 2012


@ Schaferkotter : thank you very much for the discussion regarding Open
Boundary Conditions and Orlanski scheme. Thanks for the Paper.

@ Martin : thank you .
 i have started using OBCS_balanceFacE/S/N/E.

my query here is : i am using "OBCS_balanceFacE/S/N/W = 1" along with
Orlanski BC ( useOrlanskiE/W/N/S =TRUE.) .... am i doing right here ????

 I have experimented this : The attached "fig1.jpg" has 3 experiments
(elevation at a location and elevation in the domain). Here i am not
loosing Sea Level as previously ...!!! relieved  !!! ....

Pls give a short explanation regarding the VALUES(1/-1/2/0) behind
OBCS_balanceFacE/S/N/W = 1/ -1/ 0/ 2 .

Thanking You Again...

On Mon, Sep 24, 2012 at 8:55 PM, michael schaferkotter <
schaferk at bellsouth.net> wrote:

> this phenomena is not surprising since the Orlanski scheme is based on the
> Sommerfield equation, which is known
> to not conserve volume.
>
> the proposal to balance each boundary separately will be problematic as
> you have no way of actually estimating
> the individual boundary coefficents.
>
> you might consider the method put forth:
>
> http://www.sciencedirect.com/science/article/pii/S1463500311001442
>
> On Sep 24, 2012, at 1:59 AM, Martin Losch wrote:
>
> > Himansu,
> >
> > the most likely cause for "loosing sea level" with obcs is an imbalance
> of the horizontal velocities through the boundaries. To convince yourself,
> compute the fluxes from averaged output of UVEL*DYG*DR*HFACW or
> /VVEL*DXG*DR*HFACES. The divergence of those should match your sea level
> change.
> > With Orlanski BCs you cannot prescribe velocities so that the fluxes are
> correct. I suggest that you use the balancing option described here:
> > <
> http://mitgcm.org/public/r2_manual/latest/online_documents/node236.html#SECTION00731560000000000000
> >
> >
> > M.
> >
> > On Sep 21, 2012, at 9:42 PM, Himansu Pradhan wrote:
> >
> >> Sorry for a belated answer.
> >>
> >> i too also confused why i am LOOSING sea level. Except the values for
> "Eta" i am getting a good circulation , temperature and salinity fields
> using Orlanski boundary condition.
> >>
> >> i am attaching my "input" and "code" directory for checking and
> suggesting me anything to move forward.
> >>
> >> I am trying the suggestation regarding SPONGE layer.
> >>
> >> Thank you
> >>
> >>
> >>
> >> On Thu, Sep 20, 2012 at 10:51 PM, Jody Klymak <jklymak at uvic.ca> wrote:
> >>
> >> Himansu,
> >>
> >> Not sure why you are losing sea level...
> >>
> >> For a sponge, you need the spatial extent to be at least a large
> fraction of a mode-1 tidal wavelength.  Unless you are in shallow water, 40
> km is not going to do it.  Second, your time constants are far too long to
> absorb internal waves.  They will bounce off that boundary no problem.  Use
> timescales that are closer to an hour at the outer boundary and 12 h at the
> inner side.
> >>
> >> I do something like
> >>
> >> # Sponge layer parameters
> >> &OBCS_PARM03
> >> Urelaxobcsinner=500.,
> >> Urelaxobcsbound=500.,
> >> Vrelaxobcsinner=500.0,
> >> Vrelaxobcsbound=500.0,
> >> spongeThickness=10,
> >> /
> >>
> >> but my sponge is closer to 150 km wide because I telescope the grid on
> the boundaries.
> >>
> >> Cheers,    Jody
> >>
> >> On Sep 20, 2012, at  9:44 AM, Himansu Pradhan wrote:
> >>
> >>> Hello MITgcm users,
> >>>
> >>> I have configured MITgcm for my region (174*200*17) to look after
> circulation and internal waves. dx = 2.5 km and dy = 2.5km approximatly .
> Tides are included  into the  region through momenteum equations to take
> note of the thermal/density oscillations. Orlanski OBC were applied at the
> boundary to take care for physics/dynamics.
> >>> ... here starts the PROBLEM ???? their is an ANAMALOUS DECREASE in the
> Sea Surface Level with time. ?????
> >>> .
> >>> so finding myself NOWHERE i wanted to include SPONGE layer at the
> BOUNDARY .------> Is this OK ??? .( need some suggestation/help HERE ).
> >>>
> >>> to proceed with I included :
> >>> #define ALLOW_OBCS_SPONGE   : ---> in ..../code/OBCS_OPTIONS.h
> >>>
> >>> and edited "data.obcs" :
> >>> # Open-boundaries
> >>> &OBCS_PARM01
> >>> OBCSfixTopo=.FALSE.,
> >>> OB_Jnorth=174*-1,
> >>> OB_Jsouth=174*1,
> >>> OB_Ieast=200*-1,
> >>> OB_Iwest=200*1,
> >>> # useOrlanskiN/S/E/S=.FALSE.,
> >>> useOBCSbalance=.FALSE.,
> >>> useOBCSprescribe = .FALSE.,
> >>> useOBCSsponge=.TRUE.,
> >>> &
> >>> # Orlanski parameters
> >>> &OBCS_PARM02
> >>> Cmax=0.45,
> >>> cVelTimeScale=1000.,
> >>> &
> >>> # Sponge-layer parameters
> >>> &OBCS_PARM03
> >>> Urelaxobcsinner=864000,
> >>> Urelaxobcsbound=43200,
> >>> Vrelaxobcsinner=864000,
> >>> Vrelaxobcsbound=43200,
> >>> spongeThickness=15,
> >>> &
> >>>
> >>> to start with. i had given : " # Sponge-layer parameters ;
>  &OBCS_PARM03 " .  values BLINDLY. it didn't HELP. The temperature and
> Circulation results were not satisfactory.
> >>>
> >>> I have gone trough the manual , . But need a bit more explanation
> regarding the physics of these paramaters.
> >>>
> >>> And supposingly if i want to run MITgcm with this current
> configuration for one month, ????? what should be the VALUES (of
> Sponge-layer parameters) to START with. ??????
> >>>
> >>> Any Suggestation/Help is highly appreciated.
> >>>
> >>> Thanking all in advance.
> >>>
> >>> Himansu.
> >>>
> >>>
> >>>
> >>> _______________________________________________
> >>> MITgcm-support mailing list
> >>> MITgcm-support at mitgcm.org
> >>> http://mitgcm.org/mailman/listinfo/mitgcm-support
> >>
> >> --
> >> Jody Klymak
> >> http://web.uvic.ca/~jklymak/
> >>
> >>
> >>
> >>
> >>
> >> _______________________________________________
> >> MITgcm-support mailing list
> >> MITgcm-support at mitgcm.org
> >> http://mitgcm.org/mailman/listinfo/mitgcm-support
> >>
> >> <FEBzI.zip>_______________________________________________
> >> MITgcm-support mailing list
> >> MITgcm-support at mitgcm.org
> >> http://mitgcm.org/mailman/listinfo/mitgcm-support
> >
> >
> > _______________________________________________
> > MITgcm-support mailing list
> > MITgcm-support at mitgcm.org
> > http://mitgcm.org/mailman/listinfo/mitgcm-support
>
>
> _______________________________________________
> MITgcm-support mailing list
> MITgcm-support at mitgcm.org
> http://mitgcm.org/mailman/listinfo/mitgcm-support
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://mitgcm.org/pipermail/mitgcm-support/attachments/20120926/78bee5d9/attachment.htm>


More information about the MITgcm-support mailing list