[MITgcm-support] Balance of momentum

Ryan Abernathey ryan.abernathey at gmail.com
Fri Dec 6 18:28:15 EST 2013


JMC,
What you say about the etan diagnostics makes perfect sense. I misunderstood the earlier email and thought there was some missing term that could only be output through hacking the code. I never had to do this in my momentum budgets, so I was kind of confused / worried. 
-Ryan

Sent from my iPhone

On Dec 6, 2013, at 5:48 PM, Jean-Michel Campin <jmc at ocean.mit.edu> wrote:

> Hi Marco, Dave, Martin and Ryan,
> 
> I would like to clarify few things regarding diagnostics of the
> horizontal momentum equation. I did not find very easy to follow
> the discussion started ~ a week ago (the subject has changed, not 
> always easy to distinguish quickly the answer from the quoting of 
> the original message and long lines are not easy to read from the 
> mitgcm-support archive web-site). But I hope this will not add more
> confusion.
> 
> The different horizontal momentum terms are split into different
> diagnostics, mainly to reflect the time-stepping options (forward
> in time = explicit / backward in time = implicit, and going
> through the Adams-Bashforth or not). A finer decomposition can
> be obtained with a 2nd (larger) set of diagnostics that are specific
> to the formulation used (flux-form or vector-invariant).
> 
> 1) the surface pressure horizontal gradient needs to be added.
> this was mentioned in the Emily's message that Martin (on Nov 28) pointed to:
> http://mitgcm.org/pipermail/mitgcm-support/2010-December/006921.html
> and it's needed whether or not the free-surface or the rigid-lid is used
> (in this later case, this is the pressure from the lid).
> 
> To Ryan: Having many diagnostics slows down the code (it's not 
> clear how much, but it does) and when some term can be easily and
> precisely recomputed from other diags, it's a matter of choice
> if we decide to add a diagnostics or not. 
> In this case, from a single 2-D diagnostics output (ETAN), one can 
> compute the horizontal gradient in both direction (DXC & DYC are 
> standard grid output); but we would need two 3-D diagnostics (to get it 
> properly masked) for the contribution at each level of the du/dt 
> and dv/dt equations. So, it's not very clear to me that these two
> 3-D diagnostics should be added (may be we can continue this discussion
> elsewhere). 
> 
> 2) Diagnostics Um_Advec and Um_Advec already contain the Coriolis
> term, unless one uses the CD-Scheme (and in this case, 
> Um_Cori and Vm_Cori needs to be added). 
> The main reason for this is that when using the vector-invariant form, 
> Coriolis term can be incorporated to the absolute vorticity (and then 
> treated with the advection), and in order to make the diagnostics more 
> similar between the 2 momentum formulations, the same convention was 
> retained also for the flux-form.
> Note that, for the same reasons, the metrics-terms are also part
> of Um_Advec and Um_Advec.
> 
> Now, here is where I get confused: In Emily's message, the CD-Scheme is 
> not used, but momentum budget is closed by adding Um_Cori to to Um_Advec ?
> 
> 3) Diagnostics Um_Diss and Vm_Diss contain all the explicit (forward
> in time) dissipation terms (horizontal viscosity, explicit vertical
> viscosity, side drag, bottom friction).
> If one uses implicit vertical viscosity (implicitViscosity=.TRUE.,),
> the vertical viscosity tendency needs to be added (as explained in
> Emily's message).
> 
> 4) forcing (diagnostics Um_Ext, Vm_Ext) and baroclinic pressure gradient
> (diagnostics Um_dPHdx and Vm_dPHdy) needs to be added separately.
> If filters are used (Shapiro, FFT), their contribution needs to be
> diagnosed separately (e.g., using SHAP_dU and SHAP_dV ).
> 
> 5) The Adams-Bashforth will alter the momentum budget (although over
> a long period, most cancel in time-averaging) and, depending on the 
> particular set of parameter (e.g., staggerTimeStep, forcing_In_AB, 
> momDissip_In_AB ) it's not always negligible. This AB effect can
> be diagnosed using the 2 diagnostics AB_gU and AB_gV (added 2 years ago).
> 
> So that the momentum budget for the U componentic can be diagnosed with:
> 
> [du/dt] = - gravity * ( ETAN(i) - ETAN(i-1) ) / DXC
>        + Um_dPHdx
>        + Um_Advec   (+ Um_Cori : but only if using CD-Scheme )
>        + Um_Diss    (+ Implicit vertical viscosity tendency)
>        + Um_Ext
>        + AB_gU
> 
> Some remarks:
> a) When using the flux-form formulation, the non-linear free-surface can introduce 
>   a small difference (rescaling) which is not currently diagnosed.
> b) the 3-D Coriolis ( in 2 Omega cos(Phi) between U and W) is not currently
>  diagnosed separately (and is not part of Um_Cori) but is part of
>  Um_Advec. Also using the quasiHydrostatic option might alter some of the 
>  diagnostics (e.g., what is described as the hydrostatic pressure is no 
>  longer strictly hydrostatic).
> c) when trying to do a budget over just 1 time-step, one needs to be aware 
>  of the precise time-discretisation of the the first term which is evaluated 
>  in the model as 
>   - gravity * {  implicSurfPress * (d.EtaN/dx)^n+1 
>                +(1-implicSurfPress)*(d.EtaN/dx)^n }
>  (using the default implicSurfPress=1 , it might be be easier to use 
>   snap-shot output Eta.[iter+1] to get EtaN^n+1 )
> 
> Cheers,
> Jean-Michel
> 
> On Thu, Dec 05, 2013 at 10:23:04AM +0000, David Munday wrote:
>> Hi Marco,
>> 
>> Um_HydroP isn't in my available_diagnostics, I can't find it anywhere in the documentation, and grepping the code turns up nothing. Is it the same as Um_dPHdx?
>> 
>> In general, I'd expect Um_Advec and Um_dPHdx to be similar in magnitude and of opposite sign, due to the dominance of geostrophic balance. Assuming Um_HydroP is Um_dPHdx by another name, then its magnitude is about right, but its sign is wrong (or Um_Advec's is). Taking the difference and adding Um_Diss gives me -4.3197E-9, so you're still missing a term.
>> 
>> Have you output every single Um_* diagnostic available, including the ones you think are zero, the implicit/explicit viscosity, and the gradient in SSH?
>> 
>> Best wishes,
>> 
>> Dave
>> 
>> 
>> On 4 Dec 2013, at 21:44, Marco Reale wrote:
>> 
>> Hi David,
>> 
>> and thanks a lot for your response: I followed your suggestion.
>> 
>> So finally  the equation will be :
>> 
>> TOTUTEND/86400=Um_Advec+Um_diss+Um_HydroP;
>> 
>> each of these members of the equation are 4D arrays : following a point to check the balance at t = 30, level=12, i=9, j=30 and I got the following results:
>> 
>> TOTUTEND/86400= -8.0561e-10
>> Um_Advec= 2.5398e-07
>> Um_diss= 7.7003e-09
>> Um_HydroP=2.4196e-07
>> 
>> the balance of this term is 5.0365e-07
>> 
>> 
>> So I have done something wrong: I’m working with rough data without taking into account volume of cells, etc.
>> 
>> Any suggestions?
>> 
>> Marco
> 
>> _______________________________________________
>> 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



More information about the MITgcm-support mailing list