[MITgcm-support] Regional high-res model configuration

Stanislav Martyanov martyanov.sd at gmail.com
Wed Apr 22 12:28:45 EDT 2020


Thank you very much for your answers, Martin! You have greatly helped!

Now it's time to re-configure my previous setup and try it again..

Best regards,
Stanislav Martyanov

ср, 22 апр. 2020 г. в 11:16, Martin Losch <Martin.Losch at awi.de>:

> Hi
>
> this is the post about tuning viscosities, that I couldn’t find earlier:
> http://mailman.mitgcm.org/pipermail/mitgcm-support/2006-June/004086.html
>
> M.
>
> > On 22. Apr 2020, at 09:38, Martin Losch <Martin.Losch at awi.de> wrote:
> >
> > Hi again,
> >
> >> On 21. Apr 2020, at 19:15, Stanislav Martyanov <martyanov.sd at gmail.com>
> wrote:
> >>
> >> Hi Martin!
> >>
> >> Thanks for your reply!
> >>
> >> Right to the business:
> >>
> >>>> I would use viscAhGrid or viscAhGridMin for this (something small)
> >> To be honest, I hardly understand the nature of viscAhGrid parameter...
> While viscAh can be seen as some kind of  constant physical background
> viscosity not resolved by the main horizontal turbulence parameterization
> used (Smagorinsky, Leith, etc),  viscAhGrid seems like some numerical
> trick. From its formula in MOM_CALC_VISC it can be concluded that it has
> something to do with grid size (just like Smagorinsky).
> > Please have a look at the corresponding section the online documentation
> <
> https://mitgcm.readthedocs.io/en/latest/algorithm/algorithm.html#nonlinear-viscosities-for-large-eddy-simulation>,
> and also at the header of mom_calc_visc.F, where it should become clear how
> all of this is handled.
> >>
> >>>> this is only 2D Smagorinsky. I would use 2D Leith instead.
> >> The only reason I chose Smagorinsky parameterization was its dependence
> on grid resolution and dynamical regime (velocity gradients are involved).
> But your alternative (viscAhGrid + viscC2leith), perhaps, does the same in
> a combined manner.
> > As far as I know (and I am by no means a specialist on turbulence),
> Smagorinsky parameterises 3D turbulence (which is kind of funny to then do
> it only in 2D by default in the code), and there’s even code to do that.
> The Leith scheme parameterises 2D turbulence and as long as the aspect
> ratio of your grid is so large (dx/dz) that makes more sense to me.
> >>
> >>>> setting viscAhGridMax (something a little smaller than 1)
> >>>> may be necessary
> >> In MOM_CALC_VISC  the recommended value for viscAhGridMax is 1. What
> value would you use?
> > that’s the maximum allowed for numberical “cfl-like” stability. It’s
> usually OK.
> >>
> >>>> depending on your advection scheme, horizontal diffusivities
> >>>> can be zero (I would use 77, 33 or 7 for advection)
> >> I use the advection scheme 33 (3-order DST with Sweby flux limiter) for
> T, S and seaice.
> >>
> >> Does the next set of parameters seem better than the previous one?
> >> viscAhGrid = 0.5, viscAhgridmax = 1.0, useFullLeith = .TRUE.,
> viscC2leith = 1.0, viscC2LeithD = 1.0, viscAr = 1.0e-5.
> >> Diffusivities: diffKhT = 1.0, diffKhS= 1.0, diffKrT = 1.0e-5, diffKrS =
> 1.0e-5.
> >> Here I additionally set viscC2LeithD as it is recommented in the manual
> (a damping specifically targeting purely divergent instabilities near the
> gridscale).
> > As these these are all additive, the gridMax =1 is important.
> viscAhGrid=0.5 is too large for a small background. See header of
> mom_calc_visc.F to estimate how large this will be in dimensional units
> (0.25*L**2*viscAhGrid/deltaT = 0.25 * 1e6 * 0.5 / deltaT in your case).
> Horizontal diffusivities for T and S are not necessary with your advection
> scheme.
> > In the end you’ll have to experiment with all of that. See this post:
> >>
> >>>> Don’t use the CD scheme, if necessary use small biharmonic viscosity
> >>>> instead. Implicit vertical advection is probably not required either
> >> OK, the CD_scheme is omitted then, thank you! But can Implicit vertical
> advection significally alter the results? I just hoped to use it in order
> to somewhat increase the time step.
> > You’ll have to test that. There’s some discussion on this list. I cannot
> find the exact post I was looking for, though it was a good set of
> instructions by Baylor how to tune the viscosity parameters.
> >>
> >>>> Why do you want to use the MY82 package?
> >> The only reason is that I previously used it and familiar with it.
> Moreover, the MY is a well-known mixing scheme, and both its good and bad
> sides are also documented.
> > Make sure that it is really is what you think it is. This is not the
> 2.5-level closure MY scheme, but something simpler.
> >
> >>>> I’d prefer GGL90 or KPP (although at your resolution, KPP might be
> noisy (see this thread <
> http://mailman.mitgcm.org/pipermail/mitgcm-support/2020-April/012427.html
> >).
> >> Previously our team used the KPP for another configuration, but for now
> I decided to omit it. I read the thread you mentioned, and it was the last
> straw to find another scheme :-)
> >>
> >> As for GGL90, I did not previously used it at all. From the original
> paper by Gaspar et al, 1990, I know that it has only one prognostic
> equation (for TKE) and that the turbulent length scale is diagnosed rather
> than prognosed. From your experience, can you  recommend it over the MY82?
> > GGL90 is definitely more sophisticated taht MY82. You’ll have to try it
> out.
> >>
> >>>> but it should (like GGL90) remove the instabilities,
> >>>> because it’s also Richardson number based
> >> Yes, I know, and that is why I was surprized by the comment in the code
> "convective adjustment might be needed even with ggl90". Do you use GGL90
> with or without convAdj or  ivdc_kappa?
> > I never paid attention to that, but I guess it should work without cAdj
> or ivdc_kappa.
> >>
> >> I also know that you are the author of the SEAICE package for MITgcm,
> so could you tell, is it OK for now to use the newly added feature
> SEAICE_ALLOW_MOM_ADVECTION and corresponding run-time parameters
> SEAICEmomAdvection=TRUE, SEAICEuseMetricTerms = TRUE, and
> SEAICEuseFluxForm=TRUE if my configuration uses curvilinear grid (and
> that's why I must use vectorInvariant for momentum)? The vectorInvariant
> does not use the metric terms for momentum, should I use them for SEAICE on
> my grid? And as I understood, there is no alternative for
> SEAICEuseFluxForm?
> > I have never tested the option SEAICE_ALLOW_MOM_ADVECTION properly, but
> just implemented it upon request. Please go ahead and try it.
> > As far as I understand, the momentum advection is only correct on
> curvilinear grids with the vector invariant option, but that only applies
> to the momentum advection. The entire stress formulation uses the metric
> terms by default (you should not have to set it explicitly), but recently I
> found that I may have missed a bunch of terms because of insufficient
> tensor analysis background at my end. With grids that are
> quasi-rectangular, this is probably not much of an issue.
> > SEAICEuseFluxForm only applies to the advection scheme 2 for tracers
> (HEFF/AREA), which I wouldn’t use anyway (use 77 or 33).
> >
> > Martin
> >
> >>
> >> Stanislav
> >>
> >>
> >>
> >>
> >>
> >> вт, 21 апр. 2020 г. в 14:31, Martin Losch <Martin.Losch at awi.de>:
> >> Hi Stanislav,
> >>
> >> generally your parameters look good, so comments below
> >>
> >>> On 21. Apr 2020, at 11:23, Stanislav Martyanov <martyanov.sd at gmail.com>
> wrote:
> >>>
> >>> Hello, dear colleagues!
> >>>
> >>> I am configuring the MITgcm for the Kara Sea region. With the
> curvilinear grid, the model's resolution is: dx=500-1200 m, dz=2 m in upper
> layers, and up to 20 m in the deepest. The minimum model depth is set equal
> to 5 m. I use the pure r vertical coordinate (not r*) for now. Because of
> the curvilinear grid, I use vectorinvariant package. Other packages used
> are: seaice, my82, exf, cal, obcs.
> >>>
> >>> Having carefully read (once again) the newest version of the
> online-manual and browsed the code, I have configured almost all necessary
> model parameters, but still some questions remain. They are qiute
> specialized, so I could not find the needed information in puplished
> papers. But for those experienced with MITgcm, they might be easy to
> answer. I will appreciate any advise!
> >>>
> >>> 1) The newly added feature with partial cells (#undef
> EXCLUDE_PCELL_MIX_CODE in CPP_OPTIONS). The parameter pCellMix_select = 0
> by default, which means that the enhanced mixing is OFF in bottom and
> surface layers. The parameters intended for the enhanced mixing in the
> inner layers are also OFF by default: interViscAr_pCell and
> interDiffKr_pCell = FALSE. Are this option and the corresponding run-time
> parameters really important for such configuration (Kara Sea)? Any
> experience with them?
> >> I have no experience with this, but I think you can safely keep this
> turned off.
> >>>
> >>> 2) I use the following viscocity and diffusivity values, intending to
> use the Smagorinsky viscocity for lateral turbulent exchange:
> >>> viscAh = 1.0,
> >> I would use viscAhGrid or viscAhGridMin for this (something small)
> >> this is only 2D Smagorinsky. I would use 2D Leith instead.
> >>> viscC2smag = 3.0,
> >> setting viscAhGridMax (something a little smaller than 1) may be
> necessary
> >>> viscAr = 1.0e-5,
> >> depending on your advection scheme, horizontal diffusivities can be
> zero (I would use 77, 33 or 7 for advection)
> >>> diffKhT = 1.0,
> >>> diffKhS= 1.0,
> >>> diffKrT = 1.0e-5,
> >>> diffKrS = 1.0e-5,
> >>> But for diffusivity there is no Smagorinsky analogue. What approach
> can be the most suitable in this case?
> >>>
> >>> 3) In the paper "ECCO version 4 - an integrated framework... " by
> Forget et al., 2015, the authors decided to switch OFF the C-D sheme and to
> switch ON the Implicit vertical advection option. But they investigated the
> global scales, not high-res. From your experience, can their advise be
> applyed to the high-res ocean simulations?
> >> Don’t use the CD scheme, if necessary use small biharmonic viscosity
> instead. Implicit vertical advection is probably not required either
> >>>
> >>> 4) The MY82 package does not allow using the convective adjustment
> (cAdjFreq) and ivdc_kappa. Does it mean that the MY in MITgcm works fine
> with vertical instability cases? Actually, it should, as follows from J.
> Mellor's papers, but this question arose when I was investigating the GGL90
> code, where I found the line "convective adjustment might be needed even
> with ggl90"…
> >> Why do you want to use the MY82 package? I’d prefer GGL90 or KPP
> (although at your resolution, KPP might be noisy (see this thread <
> http://mailman.mitgcm.org/pipermail/mitgcm-support/2020-April/012427.html
> >).
> >> You can just comment out the stop-statement in the code and use MY82
> with other convective adjustment, but it should (like GGL90) remove the
> instabilities, because it’s also Richardson number based. MY82 is not
> really used that much (at all?).
> >>
> >> Martin
> >>>
> >>> Once again, I will appreciate any advise!
> >>>
> >>> PS: Sorry for large wall-of-text here, but I tried to be concrete. :-)
> >>> _______________________________________________
> >>> MITgcm-support mailing list
> >>> MITgcm-support at mitgcm.org
> >>> http://mailman.mitgcm.org/mailman/listinfo/mitgcm-support
> >>
> >> _______________________________________________
> >> MITgcm-support mailing list
> >> MITgcm-support at mitgcm.org
> >> http://mailman.mitgcm.org/mailman/listinfo/mitgcm-support
> >> _______________________________________________
> >> MITgcm-support mailing list
> >> MITgcm-support at mitgcm.org
> >> http://mailman.mitgcm.org/mailman/listinfo/mitgcm-support
> >
> > _______________________________________________
> > MITgcm-support mailing list
> > MITgcm-support at mitgcm.org
> > http://mailman.mitgcm.org/mailman/listinfo/mitgcm-support
>
> _______________________________________________
> MITgcm-support mailing list
> MITgcm-support at mitgcm.org
> http://mailman.mitgcm.org/mailman/listinfo/mitgcm-support
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://mailman.mitgcm.org/pipermail/mitgcm-support/attachments/20200422/66497231/attachment-0001.html>


More information about the MITgcm-support mailing list