[MITgcm-support] Tidal forcing with Orlanski boundary conditions.
Jody Klymak
jklymak at uvic.ca
Tue May 23 13:27:21 EDT 2017
Hi Hanut,
The boundary conditions need to be altered to force the flow through them as well. Orlanski will allow a mode-1 wave to pass, but without modification will set the barotropic velocity to zero. The body force should still move the water around in the interior, but your boundary conditions will be set to zero barotropic flow. I think ye olde internal wave examples modify `calc_obcs.F`, but you don’t need to do this anymore.
For a radiating internal wave problem I recommend using obcs or rbcs and force a sponged velocity at the boundaries, and I do not use orlanski etc. The sponge needs to be wide and gentle enough to absorb internal waves without reflection.
An example setup (with steady forcing, but how to setup oscilatory forcing should be clear from `gendata.py`) is at
https://github.com/jklymak/MITgcmExampleSteadyGauss <https://github.com/jklymak/MITgcmExampleSteadyGauss>
I think you can use the external forcing of the boundary with Orlanski BC instead of a sponge.
I think you can also do all this with a body force as well, but you would want to be sure that your body force and BC forcing are in sync.
If your domain is long compared to the surface tidal wave, you may also want to lag one end of the domain relative to the other to allow the tidal waves to be in sync.
Cheers, Jody
in data &PARAM03:
`
# Forcing for boundary condition files
periodicExternalForcing=.TRUE.,
# 1/12th M2 tidal period...
externForcingPeriod = 3720.,
# one tidal period (12 records need to be in OBEuFile, etc)
externForcingCycle = 44640.,
`
and in data.obcs:
`
# Open-boundaries
&OBCS_PARM01
OB_Iwest=1*1,
OB_Ieast=1*-1,
OB_Jnorth=80*0,
OB_Jsouth=80*0,
useOBCSsponge=.TRUE.,
useOBCSprescribe=.TRUE.,
OBEuFile = 'Ue.bin',
OBWuFile = 'Uw.bin',
OBEtFile = 'Te.bin',
OBWtFile = 'Tw.bin',
/
# Orlanski parameters
&OBCS_PARM02
/
# Sponge layer parameters
&OBCS_PARM03
Urelaxobcsinner=1000.,
Urelaxobcsbound=100.,
Vrelaxobcsinner=1000.0,
Vrelaxobcsbound=100.0,
spongeThickness=10,
/`
> On May 23, 2017, at 10:12 AM, Hanut Vemulapalli <hanut2 at gmail.com> wrote:
>
> Hi all,
>
> I am trying to simulate gravity waves produced by a gaussian topography in uniform stratificaiton. I have a 2D simulation with unit depth in the y-coordinate. I am giving tidal forcing to the model in code/external_forcing.F by adding the following lines .
>
> DO j=0,sNy+1
> DO i=1,sNx+1
> gU(i,j,kLev,bi,bj) = gU(i,j,kLev,bi,bj) + 0.02*0.00014*cos(0.00014*myTime)
> ENDDO
> ENDDO
>
>
> I am prescribing Orlanski boundary condition at the east and west boundaries.
> However, the velocity magnitude tends to very small values of the order 1E-5 times the tidal amplitude (i.e 1E-12) once I prescribe open boundary condition at the east and west boundaries.
> Can anyone tell me where I am going wrong.
>
> Thanks in advance.
>
> Best regards
> Hanut
>
>
>
>
> <https://mailtrack.io/>Sent with Mailtrack <https://mailtrack.io/install?source=signature&lang=en&referral=hanut2@gmail.com&idSignature=22>_______________________________________________
> 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/20170523/cfae1364/attachment.htm>
More information about the MITgcm-support
mailing list