[MITgcm-support] Problem with tides in the OBCS

Jody Klymak jklymak at uvic.ca
Sat Jun 10 22:31:38 EDT 2017


On 10 Jun 2017, at 18:10, Dimitris Menemenlis wrote:

> Mea culpa.  Lazy programming.
> For time being I will add warning as I have no time to fix and check 
> right now.
>
> We have used useOBCStide for regional studies, and to first order it 
> works OK,
> i.e., it generates tidal ellipses as expected and it preserves volume 
> in regional
> domain, since flow perpendicular to boundaries is correctly specified. 
>  I think that
> the flow parallel to boundaries will only impact the flow inside the 
> regional domain
> via horizontal viscosity term, so probably a second order impact in 
> most set-ups.

I certainly agree that you’ll get a tidal response.  I’m just unsure 
if it will be accurate.  I’d have though fV would matter (at an E/W 
boundary) and that dV/dy would have an effect.  But you know the 
numerics better than I do, and maybe neither of those terms enter into 
the equations just inside the boundary.

I looked at the only source I know of this, [Carter and Merrifield 
(2008)](http://doi.org/10.1016/j.ocemod.2007.04.003) and they recommend 
using a Flather BC that uses both U (normal flow) and \eta. They do test 
just setting U, but they had f=0, and they set dV/dx = 0, so they 
didn’t do nothing for V.

So, I’m not at all sure what you have done is incorrect, but it’d be 
fun if someone tested it.  If I get time the next few weeks I may have a 
quick try at it, and I’ll report back.

Thanks,   Jody



>
> Dimitris Menemenlis
>
>> On Jun 2, 2017, at 12:43 PM, Jody Klymak <jklymak at uvic.ca> wrote:
>>
>> On 2 Jun 2017, at 4:52, Bahman Alavi wrote:
>>
>>
>> I'm trying to input Tide on my open boundaries; the first question 
>> is: Is
>> this a good idea to enabling "useOBCStide"?
>>
>> OK, no offense to whoever wrote obcs_add_tides.F, but I don’t 
>> really see how it is supposed to work except if the Coriolis 
>> parameter f=0. Tidal boundary conditions at (for instance) the 
>> western boundary are only specified in the east-west velocity, OBWu. 
>> Unfortunately, it doesn’t also specify OBWv, and therefore assumes 
>> the north-south velocity is zero. So while a tidal signal will almost 
>> certainly result, I don’t see how it will produce the correct tidal 
>> signal as it will constantly drive a rectilinear flow through the 
>> western boundary, when really it should have a tidal ellipse in both 
>> u and v. Maybe I haven’t fully though through the physics of this 
>> BC, and forcing v=0 doesn’t matter, but intuitively it seems wrong.
>>
>> Given this, I’d not use this code, or I’d modify it to also 
>> accept along-boundary velocities.
>>
>> Conversely, you can force tides using the following in data.obcs:
>>
>>  useOBCSprescribe=.TRUE.,
>>  OBEuFile = 'Ue.bin',
>>  OBWuFile = 'Uw.bin',
>>  OBEtFile = 'Te.bin',
>>  OBWtFile = 'Tw.bin',
>> (where I have not specified OBEvFile, but you should unless f=0) and 
>> in data
>>
>> # 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.,
>> You can see examples of how to generate the files at 
>> https://github.com/jklymak/MITgcmExampleSteadyGauss 
>> <https://github.com/jklymak/MITgcmExampleSteadyGauss>
>>
>> The second question; I extracted "amplitude"s & "phase"s of my 
>> boundary
>> points. At first I was getting errors, i thought maybe my units for 
>> am/ph
>> are wrong. I've searched and finally found the units in file
>> "obcs_add_tide.F" which is somehow weird. For amplitude: "m/s" and 
>> for
>> phase "s". Is this correct? and how should i change to these units?
>>
>> They want the “phase” specified as a dt = phase*T/2/pi, where 
>> phase is the tidal phase in radians, and T is the tidal period.
>>
>> Cheers, Jody
>>
>>
>> Thanks in advance
>> _______________________________________________
>> MITgcm-support mailing list
>> MITgcm-support at mitgcm.org
>> http://mitgcm.org/mailman/listinfo/mitgcm-support 
>> <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/20170610/b63b9e53/attachment.htm>


More information about the MITgcm-support mailing list