[MITgcm-support] open boundaries!

Matthew Mazloff mmazloff at MIT.EDU
Tue Oct 11 14:51:47 EDT 2005


Hi Van Thinh,

I have no idea why you have bad results.  If you are following what was 
done in exp4 and only changing the domain and obcs forcing function 
things should work.  Did you make sure useobcsbalance=.FALSE. and that 
the frequency is given in radians/second?
You will have to give me more info as to what is wrong with the solution 
if you want more advice.  Perhaps you should print out the prescribed 
forcing for the first few time steps to see if it is what you want it to 
be, other things may be going wrong.

-Matt


Van Thinh Nguyen wrote:

>Hi Matt,
>
>Thanks so much for helping me about OBs.
>I just would like to have more questions about this problem:
>
> I exactly follow your hints, however I just edited 'obcs_calc.F', in
>which I wrote my function at OBWu and OBEu
>(e.g. OBWu(J,K,bi,bj)=Uinflow*cos(futureTime*freq)
> where freq=M2=1.4053e-4, Uinflow=0.2
>I am not sure how the program pickup futureTime here!)
> So I've got very bad results, it was something confused!
>
>I've also done the same way as exp4, of course with different domain.
>
>Could you or someone experienced with OB give me a hint!
>
>Thanks a lot,
>
>Van Thinh
>--------------
>
>On Thu, 29 Sep 2005, Matthew Robert Mazloff wrote:
>
>  
>
>>Hi Van Thinh,
>>
>>As you know, you tell the model you have an open western boundary by writing in
>>data.obcs:
>> OB_Iwest=   55*1
>>
>>You can then prescribe a u velocity component by adding to this file:
>>
>>useOBCSprescribe=.TRUE.,
>>OBNuFile='UNBC.bin',
>>
>>This can have multiple records in it...and then in data.exf you add
>>obcsWstartdate1     = 19920101,
>>obcsWstartdate2     = 000000,
>>obcsWperiod         = 2635200.0,
>>
>>which tells the model to switch to the next velocity record every 2635200.0
>>seconds.
>>
>>A better way to impliment your function would be to write a new routine.  The
>>model sets the U velocity on the open western boundary to
>>OBWu.  This is read in when the routine OBCS_PRESCRIBE_READ executes
>>
>> call EXF_SET_OBCS_YZ(  OBWu, OBWu0, OBWu1, OBWufile, 'u'
>>     I                     , fac, first, changed, count0, count1
>>     I                     , mycurrenttime, mycurrentiter, mythid )
>>
>>Here you could change this call to
>>
>> call VAN_THINHS_FUNCTION(  OBWu, OBWu0, OBWu1, OBWufile, 'u'
>>     I                     , fac, first, changed, count0, count1
>>     I                     , mycurrenttime, mycurrentiter, mythid )
>>
>>and in the routine VAN_THINHS_FUNCTION, OBWu is set from your function.  The
>>same can be done for the eastern boundary.
>>
>>Good luck,
>>Matt
>>
>>
>>
>>Van Thinh Nguyen wrote:
>>
>>    
>>
>>>Hi Stefano and Patrick,
>>>
>>>I took a look on MITgcm Support Archives to find how to set up OBCs,
>>>however I still confuse in setting up some parameters.
>>>
>>>My test case is almost the same as exp4:
>>>
>>>Flow in a channel over a bump with
>>>- north and south boundaries are walls
>>>- west & east are OBs
>>>I would like to set the tidal forcing at west OB with velocity component
>>>U=Uo*cos(omega*t). (Uo=0.025; omega=2*pi/43200)
>>>
>>>Could you please give a hint how can I set these inputs?
>>>
>>>Thanks so much for your helps,
>>>
>>>Van Thinh
>>>----------------
>>>
>>>
>>>On Mon, 26 Sep 2005, Stefano Querin wrote:
>>>
>>>      
>>>
>>>>Hi Van-Thinh,
>>>>maybe the best thing you can do is to take a look to the files in the CVS
>>>>repository (or directly in your /pkg/obcs folder if you have a recent
>>>>version of the code). Some routines have short and useful comments that help
>>>>the user in the preliminar obcs tests.
>>>>Basically (and in very few words), you can impose:
>>>>- constant OBCs (obcs_calc.F sets u, v, w and eta to zero, and T, S to tRef
>>>>and sRef);
>>>>- radiating (passive) OBCs for momentum and tracers following Orlanski (see
>>>>Orlanski, 1976);
>>>>- time varying (active) OBCs using the obcs package (set
>>>>useOBCSprescribe=.TRUE.,) together with exf and cal.
>>>>
>>>>You can also define a sponge (or better, nudging?) layer on the OB (setting
>>>>useOBCSsponge=.TRUE., sponge thickness and relaxation parameters) and/or
>>>>prescribe the net flow across the OB using the useOBCSbalance option.
>>>>
>>>>A useful example is in:
>>>>
>>>>http://mitgcm.org/cgi-bin/viewcvs.cgi/MITgcm_contrib/heimbach/obcs_prescribe/
>>>>
>>>>As the old obcs/README said, happy open-boundarying!
>>>>
>>>>Cheers,
>>>>
>>>>Stefano
>>>>
>>>>
>>>>----- Original Message -----
>>>>From: "Van Thinh Nguyen" <vtnguyen at moisie.math.uwaterloo.ca>
>>>>To: <mitgcm-support at mitgcm.org>
>>>>Sent: Monday, September 26, 2005 6:26 PM
>>>>Subject: [MITgcm-support] open boundaries!
>>>>
>>>>
>>>>        
>>>>
>>>>>Dear MITgcm-supporters & users,
>>>>>
>>>>>Could someone please explain for me more detail about how to set up an
>>>>>open boundary conditions. For instance, in exp4:
>>>>>
>>>>>when we set in data.obcs:
>>>>>
>>>>>"OB_Jnorth=80*-1" to define a face where OB is appied, however I still
>>>>>don't know which boundary condition is appied here (Neumann?).
>>>>>
>>>>>And what is the "useOrlanski" in OB option?
>>>>>
>>>>>There are probably some docu around in MITgcm users' manuals (I actually
>>>>>took a look but haven't found), if someone know please give me a hint?
>>>>>
>>>>>Thanks a lot!
>>>>>
>>>>>Van Thinh
>>>>>_______________________________________________
>>>>>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
>>>
>>>      
>>>
>>_______________________________________________
>>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/20051011/f1ba14b3/attachment.htm>


More information about the MITgcm-support mailing list