[MITgcm-support] dissipation and diffusion rates
Martin Losch
Martin.Losch at awi.de
Wed Oct 2 03:58:56 EDT 2019
Hi Yangxin,
1) Um_Diss, etc is “guDissip” in the code, which contains only the explicit contributions. With implicit viscosity, you need the extra terms that are derived from VISrI_Um and Vm (the vertical diffusive fluxes of momentum, that are computed with an implicit method).
2) Um_Diss is the tendency due to dissipation, but VISrE/I_Um is the vertical flux (du/dz*area). To get to a tendency, you need to compute (VISrE_Um[k]-VISrE_Um[k+1])/volumeOfGridCell.
3) the same is true for 'DFrE_TH' and ‘DFrI_TH’. They are vertical fluxes that already include the area. The units are correct to, as far as I can see, because the diffusion term (d/dz kappa dtheta/dz) should have degC/s, so kappa dtheta/dz has the unis degC/s * m, and kappa dtheta/dz * area then degC/s m^3, you need to do the same computation as for VISrE/I_Um to get the tendencies.
Martin
> On 27. Sep 2019, at 19:46, Yangxin He <y67he at uwaterloo.ca> wrote:
>
> Hi there,
>
> I am currently running a simple oceanic internal wave simulation. I have implemented Smagorinsky for horizontal viscosity, pp81 for the vertical viscosity and diffusion. I have the horizontal diffusion set to 0. Part of my data file is attached in the end.
> I would like to output dissipation/diffusion rates using the diagnostic package. I know I should have Um_Diss, Vm_Diss and Wm_Diss. This way my dissipation rate for the kinetic energy will be
> Kt=U*(Um_Diss)+V*(Vm_Diss)+W*(Wm_Diss).
> pp81 can only be used if implicitViscosity and implicitDiffusion are on. Few questions here
> 1) From reading some of the previous posts, I got the impression that if implicitViscosity is on, I need to add more terms to my expression of Kt. Why?
>
> 2) If 1) is true, I will need 'VISrE_Um', 'VISrI_Um', 'VISrE_Vm' and ''VISrI_Vm'', but their units are 'm^4/s^2' instead of 'm/s^2'. How do I understand this? and how can I add them to Kt?
>
> 3) For the diffusion rate of the density equation, I will need both 'DFrE_TH' and 'DFrI_TH'? Their units are degC.m^3/s, but I am expecting degC/s. Does this need some scaling? or my understanding is not quite right?
>
> viscAz=1.E-5,
> viscAh=0.,
> no_slip_sides=.FALSE.,
> no_slip_bottom=.FALSE.,
> diffKhT=0.,
> diffKzT=1.E-5,
> #new parameters for smag
> viscC2Smag=2.0,
> # viscC4Smag=1.0,
> f0=6.5E-5,
> beta=0.,
> eosType='LINEAR',
> tAlpha=2.E-4,
> sBeta =0.,
> rhonil=1028.0,
> rhoConst=1028.0,
> gravity=9.81,
> rigidLid=.TRUE.,
> implicitFreeSurface=.FALSE.,
> # exactConserv=.TRUE.
> nonHydrostatic=.FALSE.,
> # minimum cell fraction. This reduces steppiness from Jody
> hFacMin=0.1,
> # implicSurfPress=0.5,
> # implicDiv2DFlow=0.5,
> # nonlinFreeSurf=3,
> # hFacInf=0.2,
> # hFacSup=1.8,
> #- not safe to use globalFiles in multi-processors runs
> globalFiles=.TRUE.,
> readBinaryPrec=64,
> writeBinaryPrec=64,
> writeStatePrec=64,
> saltStepping=.FALSE.,
> #parameters for pp81
> implicitDiffusion=.TRUE.,
> implicitViscosity=.TRUE.,
> #change advection shceme
> tempAdvScheme=33,
> staggerTimeStep=.TRUE.,
> &
>
>
> Thanks
>
> Yangxin
>
>
> _______________________________________________
> MITgcm-support mailing list
> MITgcm-support at mitgcm.org
> http://mailman.mitgcm.org/mailman/listinfo/mitgcm-support
More information about the MITgcm-support
mailing list