[MITgcm-devel] Hack to increase mixing near bottom/surface
Daniel Goldberg
dngoldberg at gmail.com
Tue Jan 27 11:58:08 EST 2015
Hi J-M,
On Tue, Jan 27, 2015 at 4:19 PM, Jean-Michel Campin <jmc at ocean.mit.edu>
wrote:
> Hi Dan,
>
> On Mon, Jan 26, 2015 at 09:43:42AM +0000, Daniel Goldberg wrote:
> > Hi J-M and Martin
> >
> > I have looked at the s/r's and the tests and this looks very compelling.
> > Regarding interViscAr_pCell and interDiffKr_pCell, these changes do not
> > seem related to increasing viscosity and diffusivity in thin cells, but
> > rather to improve upon discretisation that already exists in cells just
> > above the bottom or below the surface, is this correct?
> You are right.
>
> > I see that with pCellMix_select, you override the cell visc and diff with
> > imposed values (increased by a factor of how much smaller hfacC is than
> > some threshold value).
> using pCellMix_select is not able to decrease visc or diff, can only result
> in an increase of visc and/or diff (use a MAX(,) in modified code)
>
No, i did not think this was the case
>
> > 2 questions on this -- (1) should there not be an
> > upper limit to prevent infinite diffusivity and viscosity, and (2) would
> > these replace any values calculated by a physical package like KPP or
> GGL90?
> initial hFac cannot be smaller than hFacMin ; and with Non-Lin Free-Surf
> (NLFS),
> hFac cannot become smaller than hFacInf. So will not produce infinite
> value.
> And regarding (2): in file "pCellMix.notes" it is written:
> > limitations of the present hack:
> > - should not be used with NH
> > - should not be used with KPP
> and with other than KPP mixing scheme, not really "replacing" GGL90 value
> (see comment above); Below ice-shelf, GGL90 is probably not going to
> get enough mixing since the friction stress is not account for (as the
> surface wind stress is over open ocean) as a source of TKE.
>
sorry, i did not notice the MAX(,) -- so i withdraw what i said about
"replacing"
> And would need more test regarding the use of pCellMix_select with GGL90
> over open-ocean (when using NLFS and thin surface level thickness).
>
> > I also am unsure where vertical viscosities and diffusivities are
> > physically located. Are they located at vertical (upward-facing) cell
> > faces? In the center of these faces or at the boundaries?
> diffusivity kappaRTr is at same location as vertical velocity ;
> viscosity kappaRU, kappaRV are above uVel,vVel, at interface between 2
> levels.
>
> > Further questions below --
> >
> > On Fri, Jan 16, 2015 at 7:05 PM, Jean-Michel Campin <jmc at ocean.mit.edu>
> > wrote:
> >
> > > Hi Dan and Martin,
> > >
> > > I made some simple test with a hack to increase mixing near above
> bottom
> > > and/or
> > > below surface when the bottom/surface grid cell is thin (small hFac).
> > >
> > > I put some files here: http://mitgcm.org/~jmc/pCellMix/
> > > with some documentation in file "pCellMix.notes"
> > > and some plots from my simple test.
> > >
> > > the hack: it should move to the main code at some point, but
> > > 1) need more tests (and feedback)
> > > 2) it involves changing more S/R (> 7 but only 2 right now) since it
> > > needs to be pushed where vertical diffusion and viscosity are applied:
> > > ( mom_u/v_implicit_r.F mom_u/v_rviscflux.F
> > > impldiff.F gad_diff_r.F gad_implicit_r.F )
> > > and passing more arguments to these S/R (-> also changing the calling
> S/R)
> > >
> > > Some remarks:
> > > 1) we should probably always use selectBotDragQuadr=1 or 2 instead of
> the
> > > original
> > > discretisation (selectBotDragQuadr=0), see
> > > http://mitgcm.org/~jmc/pCellMix/fig_wt4.ps
> > > 2) accounting for partial cell in the interior (interViscAr_pCell=T,
> > > interDiffKr_pCell=T) is safe when using implicit visc/diff ;
> > > but also accounting for hFac in bottom friction (bottomVisc_pCell=T)
> when
> > > using no-slip-bottom=T can be unstable; and we don't have code for
> > > implicit
> > > bottom friction + it would not work with present solve_for_pressure
> code.
> > >
> >
> > I am confused by the 3rd set of tests -- you seem to be comparing the
> > effect of drag formulation, not the partial cell treatment of viscosity.
> It depends on your point of view: If I compare fig_wt3.ps next to
> fig_wt4.ps,
> I have both effects.
>
> > Also -- I was unaware that vertical viscosity and diffusivity could be
> made
> > implicit, could you tell me which parameters/CPP options enable this?
> this is documented in the manual (section 2.6) but the run-time params are
> not listed there; you can find them in the KPP section (6.4.2.3.2);
> it's also used in many experiments (lab_sea, vermix, ...) including in few
> tutorials (e.g., tutorial_global_oce_biogeo)
>
> > > Regarding how this could be use with pkg/shelfice (to avoid using
> > > SHELFICEboundaryLayer):
> > >
> >
> > This change seems ready to test with the Shelfice package, in that the
> > surface implementation is coded for ksurfC>1..
> >
> >
> > > 1) in isomip test experiment, no_slip_shelfice=F and this is good,
> > > otherwise
> > > the velocity below ice-shelf is strongly reduced when hFac is small
> > > (and it's
> > > physical) and make the melt-rate very dependent on hFac.
> > > 2) using a type of quadratic drag: the equivalent of
> selectBotDragQuadr >
> > > 0 is not
> > > coded in shelfice_u/v_drag.F but should not be hard (+ I added the
> new
> > > argument)
> > >
>
> > There is a parameter call SHELFICEDragQuadratic and it does seem to make
> > drag quadratic (but it looks as though it can only be implicit, and
> > probably not implemented the same way as selectBotDragQuadr > 0)
> SHELFICEDragQuadratic is documented in manual (sec 6.6.3.3.2) and
> default value is bottomDragQuadratic (unless SHELFICEuseGammaFrict=T);
> the 3 isomip tests all use SHELFICEDragQuadratic, as reported
> in STDOUT, with no linear drag nor no-slip BC.
> I don't know exactly what you mean by "it can only be implicit".
> My point (2) is exactly this: selectBotDragQuadr > 0 is not
> (yet) coded in shelfice_u/v_drag.F
>
i meant to say *explicit* not implicit -- apologies. and this was just
based on my inference of what the S/R does. I guess my point was that
SHELFICEDragQuadratic .ne. 0 implements a quadratic dependence of
uDragTerms on velocity, no? Is the point that lines 84-87 of
shelfice_u_drag do not make use of kappaRU, and therefore pCellMix will
have no effect?
>
> > > 3) might need to test different averaging (e.g. wet-point) of the
> velocity
> > > to the grid-cell center (in shelfice_thermodynamics.F)
> > > I don't know if there is a diagnostic for this (but could be useful)
> > >
>
> > Should e.g. 0.5 * (uFld(i,j) + uFld(i+1,j)) replace uFld(i,j)? I believe
> > these currently are the velocities in ksurfC(i,j). Should any thought be
> > given to circumstances when not all are available? What if, for instance,
> > maskW(i+1,j,ksurfC(i,j)) = 0 due to variable ice shelf topography?
> I was refering to shelfice_thermodynamics.F, lines 252-257 where
> uLoc,vLoc are computed. A wet-point averaging (like with
> selectBotDragQuadr=2
> in mom_u_bottomdrag.F) would be different, in particular when there
> is a "step" with one of the 2 maskW(i,ksurfC(i)) or maskW(i+1,ksurfC(i))
> is zero, as you mentionned.
>
> Cheers,
> Jean-Michel
>
> > Thanks for beginning this work,
> > Dan
> >
> >
> >
> > > 4) would be interesting to try SHELFICEboundaryLayer=F and
> > > pCellMix_select=10 ;
> > >
> > > Cheers,
> > > Jean-Michel
> > >
> > > _______________________________________________
> > > MITgcm-devel mailing list
> > > MITgcm-devel at mitgcm.org
> > > http://mitgcm.org/mailman/listinfo/mitgcm-devel
> > >
> >
> >
> >
> > --
> >
> > Daniel Goldberg, PhD
> > Lecturer in Glaciology
> > School of Geosciences, University of Edinburgh
> > Geography Building, Drummond Street, Edinburgh EH8 9XP
> >
> >
> > em: D <dgoldber at mit.edu>an.Goldberg at ed.ac.uk
> > web: http://ocean.mit.edu/~dgoldberg
>
> > _______________________________________________
> > MITgcm-devel mailing list
> > MITgcm-devel at mitgcm.org
> > http://mitgcm.org/mailman/listinfo/mitgcm-devel
>
>
> _______________________________________________
> MITgcm-devel mailing list
> MITgcm-devel at mitgcm.org
> http://mitgcm.org/mailman/listinfo/mitgcm-devel
>
--
Daniel Goldberg, PhD
Lecturer in Glaciology
School of Geosciences, University of Edinburgh
Geography Building, Drummond Street, Edinburgh EH8 9XP
em: D <dgoldber at mit.edu>an.Goldberg at ed.ac.uk
web: http://ocean.mit.edu/~dgoldberg
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://mitgcm.org/pipermail/mitgcm-devel/attachments/20150127/c12a439a/attachment-0001.htm>
More information about the MITgcm-devel
mailing list