Opened 16 years ago

Closed 16 years ago

Last modified 16 years ago

#149 closed upgrade (Done)

Time dependent horizontal diffusion/viscosity coefficients

Reported by: arango Owned by: arango
Priority: major Milestone: Release ROMS/TOMS 3.2
Component: Nonlinear Version: 3.2
Keywords: Cc:

Description

The horizontal diffusion (t3dmix2_s.h, t3dmix4_s.h, t3dmix2_geo.h, t3dmix4_geo.h, t3dmix2_iso.h, t3dmix4_iso.h) and viscosity (uv3dmix2_s.h, uv3dmix4_s.h, uv3dmix2_geo.h, uv3dmix4_geo.h) routines were updated to include time-dependent mixing coefficients.

Few new CPP options are introduced:

   TS_U3ADV_SPLIT      use if 3rd-order upstream split tracer advection  
   TS_SMAGORINSKY      use to turn ON or OFF Smagorinsky-like diffusion

   UV_U3ADV_SPLIT      use if 3rd-order upstream split momentum advection  
   UV_SMAGORINSKY      use to turn ON or OFF Smagorinsky-like diffusion

The primary goal of this update is to correct for the spurious numerical diapycnal diffusion of the tracer advection in coarse basin scale applications. Since this requires a velocity-dependent 3D horizontal diffusion coefficient, a Smagorinsky-like mixing coefficient is also coded. However, its not recommended to use this option on tracers since it is not clear, physically, how to compute the Smagorinsky tracer diffusion from momentum viscosity, as discussed by Griffies and Hallberg (2000). It is added for pedagogical reasons. The Smagorisky mixing on momentum makes more sense.

The 3rd-order upstream split advection (TS_U3ADV_SPLIT) can be used to correct for the spurious diapycnal diffusion of the 3rd-order upstream advection operator in terrain-following coordinates. The advection operator is split between advective and diffusive components, as shown by Holland et al. (1998), Webb et al. (1998), and Marchesiello et al. (2008). They show that the split operator has a 4th-order centered differences term plus a 3rd-order hyper-diffusion to correct for spurious diapycnal mixing.

The total (time dependent plus time independent) horizontal mixing coefficient are computed in new routine hmixing.F. Many thanks to Patrick Marchesiello for providing his codes from which hmixing.F is based.

Notice that several internal cpp options are activated in globaldefs.h to facilitate the split of the advection operator:

#ifdef TS_U3ADV_SPLIT
# define DIFF_3DCOEF
# ifdef TS_U3HADVECTION
#  undef TS_U3HADVECTION
# endif
# ifndef TS_C4HADVECTION
#  define TS_C4HADVECTION
# endif
# ifndef TS_C4VADVECTION
#  define TS_C4VADVECTION
# endif
# ifndef TS_DIF4
#  define TS_DIF4
# endif
# ifdef TS_DIF2
#  undef TS_DIF2
# endif
# ifdef TS_SMAGORINSKY
#  undef TS_SMAGORINSKY
# endif
#endif

It is also recommended to use the rotated mixing tensor along geopotentials (MIX_GEO_TS) for the biharmonic operator.

The modifications to the tangent linear, representer, and adjoint models will be updated soon.

Change History (2)

comment:1 by arango, 16 years ago

Resolution: Done
Status: newclosed

Also corrected a small problem roms_export.F and esmf_roms.F. Additional logic was needed for array Outfield, as reported in the forum.

comment:2 by Barbara, 16 years ago

The "SMAGORINSKI" in cppdefs.h should be changed into "SMAGORINSKY" as those in other codes.

Note: See TracTickets for help on using tickets.