[MITgcm-support] diagnosing gmredi fluxes for tracers
Christoph Voelker
christoph.voelker at awi.de
Wed Oct 30 03:35:48 EDT 2013
Hi Jean-Michel,
I am indeed using implicitDiffusion=.TRUE. and my data.gmredi contains
that I am using the Large, Damabasoglu, Doney tapering scheme, and the
Visback variable kappa. I have attached data.gmredi and
GMREDI_OPTIONS.h, so if you could have a quick look I'd be grateful!
Cheers, Christoph
On 10/29/13 11:05 PM, Jean-Michel Campin wrote:
> Hi Christoph,
>
> There are ways to get diagnostics for some of the different pieces of GM-Redi.
> But it would help if we know what parameter you are using (from data.gmredi)
> and if you are using implicitDiffusion=.TRUE., (from main paramter file "data").
>
> The fact that you have:
>> #define GM_BOLUS_ADVEC
>> in GMREDI_OPTIONS.h.
> does not mean that you are using the advective form (it's just that the
> code is compiled, and you can use it or not).
>
> Cheers,
> Jean-Michel
>
> On Mon, Oct 28, 2013 at 11:34:05PM +0100, Christoph Voelker wrote:
>> Dear MITgcm'ers,
>>
>> I have a run of a biogeochemical model (i.e. with several additional
>> passive tracers) and I am trying to diagnose (and hopefully understand)
>> what the GMREDI eddy parameterization contributes to my tracer fluxes.
>> In doing so I got somewhat confused my the many options and parameter
>> settings that are possible in that package; I found it hard to follow
>> how the advective and diffusive fluxes are really calculated.
>>
>> My settings include
>> #define GM_VISBECK_VARIABLE_K
>> #define GM_BOLUS_ADVEC
>> in GMREDI_OPTIONS.h.
>> so, I am using the bolus velocity form, not the skew flux form, with
>> variable kappa_GM, a la Visbeck et al. '97. In addition I am using the
>> Large, Danabasoglu, Doney tapering scheme.
>>
>> My questions concern the diagnostics associated with the tracer fluxes.
>>
>> - Firstly, I assume that the explicit diffusive terms like DFrETr01
>> contain the Redi isopycnal diffusion, perhaps plus the additional
>> diffusion induced e.g. by the KPP sheme that I am using. Is this right?
>>
>> - Secondly, is the bolus advective flux part of the advective flux
>> diagnostics, like ADVrTr01? Or is this only the non-bolus advection?
>>
>> - And finally, is there a way to diagnose the advective tracer fluxes
>> from the bolus velocity alone, i.e. separate it from the background
>> velocity advection, in a diagnostic variable, or does one have to write
>> out the bolus stream function terms GM_PsiX and GM_PsiY, and then
>> compute the induced additional advection offline?
>>
>> Sorry if my questions are daft, but I really first tried to understand
>> all this from the code,
>>
>> Thanks, Christoph
>>
>> --
>> Christoph Voelker
>> Alfred Wegener Institute for Polar and Marine Research
>> Am Handelshafen 12
>> 27570 Bremerhaven, Germany
>> e: Christoph.Voelker at awi.de
>> t: +49 471 4831 1848
>>
>>
>> _______________________________________________
>> 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
--
Christoph Voelker
Alfred Wegener Institute for Polar and Marine Research
Am Handelshafen 12
27570 Bremerhaven, Germany
e: Christoph.Voelker at awi.de
t: +49 471 4831 1848
-------------- next part --------------
# GM+Redi package parameters:
#-from MOM :
# GM_background_K: G & Mc.W diffusion coefficient
# GM_maxSlope : max slope of isopycnals
# GM_Scrit : transition for scaling diffusion coefficient
# GM_Sd : half width scaling for diffusion coefficient
# GM_taper_scheme: slope clipping or one of the tapering schemes
# GM_Kmin_horiz : horizontal diffusion minimum value
#-Option parameters (needs to "define" options in GMREDI_OPTIONS.h")
# GM_isopycK : isopycnal diffusion coefficient (default=GM_background_K)
# GM_AdvForm : turn on GM Advective form (default=Skew flux form)
&GM_PARM01
GM_background_K = 10.,
GM_taper_scheme = 'ldd97',
GM_maxSlope = 1.e-2,
GM_Kmin_horiz = 50.,
GM_Scrit = 4.e-3,
GM_Sd = 1.e-3,
# parameters taken from Gnanadesikan et al 2006, JCL
GM_Visbeck_alpha = 0.07,
GM_Visbeck_length = 50.e+3,
GM_Visbeck_depth = 2000.,
GM_Visbeck_maxval_K= 600.,
&end
-------------- next part --------------
C $Header: /u/gcmpack/MITgcm/pkg/gmredi/GMREDI_OPTIONS.h,v 1.9 2005/07/29 18:24:35 edhill Exp $
C $Name: $
C CPP options file for GM/Redi package
C
C Use this file for selecting options within the GM/Redi package
#ifndef GMREDI_OPTIONS_H
#define GMREDI_OPTIONS_H
#include "PACKAGES_CONFIG.h"
#ifdef ALLOW_GMREDI
#include "CPP_OPTIONS.h"
#undef GMREDI_COS_SCALING
C Designed to simplify the Ajoint code:
C exclude the clipping/tapering part of the code that is not used
C #define GM_EXCLUDE_CLIPPING
C #define GM_EXCLUDE_AC02_TAP
C #define GM_EXCLUDE_TAPERING
C This allows to use Visbeck et al formulation to compute K_GM+Redi
#define GM_VISBECK_VARIABLE_K
C This allows the leading diagonal (top two rows) to be non-unity
C (a feature required when tapering adiabatically).
#define GM_NON_UNITY_DIAGONAL
C Allows to use different values of K_GM and K_Redi ; also to
C be used with the advective form (Bolus velocity) of GM
#define GM_EXTRA_DIAGONAL
C Allows to use the advective form (Bolus velocity) of GM
C instead of the Skew-Flux form (=default)
#define GM_BOLUS_ADVEC
#endif /* ALLOW_GMREDI */
#endif /* GMREDI_OPTIONS_H */
CEH3 ;;; Local Variables: ***
CEH3 ;;; mode:fortran ***
CEH3 ;;; End: ***
More information about the MITgcm-support
mailing list