[MITgcm-support] diagnosing gmredi fluxes for tracers

Jean-Michel Campin jmc at ocean.mit.edu
Wed Nov 13 00:39:46 EST 2013


Hi Christoph,

Does not look like you got an answer on support to your question after 
this clarification (including data.gmredi & GMREDI_OPTIONS.h).
so if it's not too late, I will try to address your question.

In your data.gmredi, which contains this comment:
# GM_AdvForm     :      turn on GM Advective form       (default=Skew flux form)
GM_AdvForm is not set to TRUE, which means that you are using 
the Skew-Flux form of GM.
You can find some description in the manual about this part, 
in section "6.4.1.3 Griffies Skew Flux".

Now regarding pkg/diagnostics options: 
The isopycnal diffusion (Redi part) contribution is put into the
diffusive flux, splitted into the 3 directions and, 
if using implicitDiffusion=.TRUE., with a distinction
between explicit vertical flux (which contains the contribution
 Kwx d.T/dx + Kwy d/T/dy) and the implicit vertical flux 
(which contains the contribution Kwz d.T/dz).
For instance, for temperature (but it's available for any tracer):
   112 |DFrE_TH | 15 |       |WM      LR|degC.m^3/s      |Vertical Diffusive Flux of Pot.Temperature (Explicit part)  
   113 |DFxE_TH | 15 |   114 |UU      MR|degC.m^3/s      |Zonal      Diffusive Flux of Pot.Temperature
   114 |DFyE_TH | 15 |   113 |VV      MR|degC.m^3/s      |Meridional Diffusive Flux of Pot.Temperature
   115 |DFrI_TH | 15 |       |WM      LR|degC.m^3/s      |Vertical Diffusive Flux of Pot.Temperature (Implicit part)  
And if using the Skew-Flux form, the GM part (which double Kwx & Kwy
and cancel Kuz and Kvz if K_GM = K_iso ) will also fall into these diffusive
flux diagnostics.

With the advective form of GM ( GM_AdvForm=.TRUE.), and assuming the default 
GM_AdvSeparate=.FALSE. is used, the eulerian mean and the bolus transport
are combined together to form the residual mean which is used to advect all
tracers. And in this case, the GM contribution to the tracer tendency is 
part of the advective flux diagnostics, e.g. for temperature:
   109 |ADVr_TH | 15 |       |WM      LR|degC.m^3/s      |Vertical   Advective Flux of Pot.Temperature
   110 |ADVx_TH | 15 |   111 |UU      MR|degC.m^3/s      |Zonal      Advective Flux of Pot.Temperature 
   111 |ADVy_TH | 15 |   110 |VV      MR|degC.m^3/s      |Meridional Advective Flux of Pot.Temperature
Note that these diagnostics match exactly what the model computes
since they come directly out of the model advection subroutines.
But no distinction is made between the advective flux due to the eulerian mean 
and the one due to the bolus transport (hard to define anyway when 
the advection scheme is non-linear regarding the velocity).

For this reason, but only available for temperature (sorry), there is an
approximate estimation of the Bolus transport of temperature
(consistent with a centered in time and space advection scheme, like tempAdvScheme=2) 
which is diagnose (in pkg/gmredi/gmredi_xtransport.F & pkg/gmredi/gmredi_ytransport.F)
and available with the advective form:
   209 |GM_ubT  | 15 |   210 |UUr     MR|degC.m^3/s      |Zonal Mass-Weight Bolus Transp of Pot Temp
   210 |GM_vbT  | 15 |   209 |VVr     MR|degC.m^3/s      |Meridional Mass-Weight Bolus Transp of Pot Temp

An other alternative, since all Redi (+GM if skewFlux form) matrix coefficient are available
as diagnostics:
   197 |GM_Kux  | 15 |   198 |UU P    MR|m^2/s           |K_11 element (U.point, X.dir) of GM-Redi tensor
   198 |GM_Kvy  | 15 |   197 |VV P    MR|m^2/s           |K_22 element (V.point, Y.dir) of GM-Redi tensor
   199 |GM_Kuz  | 15 |   200 |UU      MR|m^2/s           |K_13 element (U.point, Z.dir) of GM-Redi tensor
   200 |GM_Kvz  | 15 |   199 |VV      MR|m^2/s           |K_23 element (V.point, Z.dir) of GM-Redi tensor
   201 |GM_Kwx  | 15 |   202 |UM      LR|m^2/s           |K_31 element (W.point, X.dir) of GM-Redi tensor
   202 |GM_Kwy  | 15 |   201 |VM      LR|m^2/s           |K_32 element (W.point, Y.dir) of GM-Redi tensor
   203 |GM_Kwz  | 15 |       |WM P    LR|m^2/s           |K_33 element (W.point, Z.dir) of GM-Redi tensor
would be to calculate (offline) an approximate evaluation (neglecting time correlation between tracer
gradient and Redi-GM matrix coeff) of gmredi flux, base on the time-mean tracer gradient and 
the time-mean GM-Redi coefficient.
 
Note that, but again for temperature only, some gmredi flux are available as diagnostics:
   206 |GM_KuzTz| 15 |   207 |UU      MR|degC.m^3/s      |Redi Off-diagonal Temperature flux: X component
   207 |GM_KvzTz| 15 |   206 |VV      MR|degC.m^3/s      |Redi Off-diagonal Temperature flux: Y component
   208 |GM_KwzTz| 15 |       |WM      LR|degC.m^3/s      |Redi main-diagonal vertical Temperature flux

Hope I did not add too much confusion here.
Cheers,
Jean-Michel

On Wed, Oct 30, 2013 at 08:35:48AM +0100, Christoph Voelker wrote:
> 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
> 



> _______________________________________________
> MITgcm-support mailing list
> MITgcm-support at mitgcm.org
> http://mitgcm.org/mailman/listinfo/mitgcm-support




More information about the MITgcm-support mailing list