PUBLIC INTERFACE ~ PUBLIC ROUTINES ~ NAMELIST

Module ocean_bihgen_friction_mod

Contact:  S. M. Griffies
Reviewers: 
Change History: WebCVS Log


OVERVIEW

This module computes thickness weighted and density weighted time tendency for horizontal velocity arising from biharmonic friction.

This module computes thickness weighted and density weighted time tendency for horizontal velocity arising from biharmonic friction. The viscosity used to determine the strength of the tendency can be a general function of space and time as specified by the Smagorinsky approach as well as a grid-scale dependent background viscosity. The form of the friction operator can be isotropic or anisotropic in the lateral plane.


OTHER MODULES USED

       constants_mod
diag_manager_mod
fms_mod
mpp_domains_mod
mpp_mod
ocean_domains_mod
ocean_obc_mod
ocean_operators_mod
ocean_parameters_mod
ocean_types_mod
ocean_util_mod
ocean_workspace_mod

PUBLIC INTERFACE

ocean_bihgen_friction_init:
bihgen_friction:
ncar_boundary_scale_read:
ncar_boundary_scale_create:
BDX_EU_smag:
BDY_NU_smag:
bihgen_viscosity_check:
bihgen_reynolds_check:
compute_neptune_velocity:


PUBLIC ROUTINES

  1. ocean_bihgen_friction_init

    DESCRIPTION
    Initialize the horizontal biharmonic friction module by registering fields for diagnostic output and performing some numerical checks to see that viscosity is set appropriately.


  2. bihgen_friction

    DESCRIPTION
    This subroutine computes the time tendency for horizontal velocity (i.e., the acceleration) from horizontal biharmonic friction. The algorithm is derived from a functional approach that ensures kinetic energy is consistenty dissipated for all flow configurations. The triad do-loops are expanded in order to enhance the ability of cache-based machines to keep most of the variables on-cache. Fundamental to the scheme are the rates of horizontal deformation
    horizontal tension = DT = (dy)(u/dy)_x - (dx)(v/dx)_y
    horizontal strain = DS = (dx)(u/dx)_y + (dy)(v/dy)_x
    Units of the tension and strain are sec^-1. Four tensions and four strains are computed for each velocity point,
    corresponding to the four triads surrounding the point.
    The following notation is used to distinguish the triads:
    (0,1)=northwest triad (1,1)=northeast triad,
    (0,0)=southwest triad, (1,0)=southeast triad A triad contributes when at least one of its velocities is not a land point. In order to obtain the correct tension and strain next to boundaries, tension and strain should not be masked with umask. As shown in Griffies and Hallberg (2000), a biharmonic operator with a nonconstant viscosity is guaranteed to dissipate kinetic energy *only* when using the sqrt of the biharmonic viscosity at each of the two stages of the algorithm. The sqrt approach is employed here.


  3. ncar_boundary_scale_read

    DESCRIPTION
    Read in the 3d ncar boundary scaling field and use this to rescale the background viscosities. To use this routine, we need to already have generated the field ncar_rescale using the routine ncar_boundary_scale_create. The advantage of reading ncar_rescale is that we do not need to introduce any global 2d arrays required for ncar_boundary_scale_create. So the idea is to pay the price once by running ncar_boundary_scale_create, save ncar_rescale, then read that field in during subsequent runs through ncar_boundary_scale_read. Here are the steps: 1/ run one time with ncar_boundary_scaling_read=.false. and ncar_boundary_scaling=.true. Be sure that the field ncar_rescale is saved in diagnostic table. To ensure answers agree whether reading ncar_rescale or creating it during initialization, it is necessary to save ncar_rescale using the double precision option in the diagnostic table (packing=1). 2/ extract field ncar_rescale from the diagnostics output and place into its own file INPUT/ncar_rescale.nc example extraction using ncks: ncks -v ncar_rescale 19900101.ocean_month.nc ncar_rescale.nc 3/ set ncar_boundary_scaling_read=.true. and ncar_boundary_scaling=.true., and now run the model reading in ncar_rescale rather than regenerating it during each initialization (which can be a bottleneck for large models on huge processor counts). 4/ As a check that all is fine, save ncar_rescale as a diagnostic for both the create and the read stage and make sure they agree. Also, all checksums should agree whether reading in ncar_rescale or creating it each initialization, so long as the ncar_rescale.nc was saved with double precision (see step 1/ above).


  4. ncar_boundary_scale_create

    DESCRIPTION
    Recale the background viscosities to be larger in the western boundary regions. The algorithm is taken directly from the anisotropic_ncar routine in ocean_lapgen_friction.F90. NOTE: The nearest western boundary computations are done along the model i-grid lines. Therefore, viscosity based on these are only approximate in the high Northern Hemisphere when using generalized coordinates with coordinate pole(s) shifted onto land.


  5. BDX_EU_smag

    DESCRIPTION
    Compute backwards Derivative in X of a quantity defined on the east face of a U-cell. Slightly modified version of BDX_EU used in ocean_operators.F90. If input is a(i,j) then output is defined at (i-1/2,j). BDX_EU_smag changes dimensions by m^-3


    INPUT
    a    field defined on the east face of a U-cell
       [real, dimension(isd:ied,jsd:jed)]

  6. BDY_NU_smag

    DESCRIPTION
    Compute backwards Derivative in Y of a quantity defined on the north face of a U-cell. Slightly modified version of BDY_EU used in ocean_operators.F90. If input is a(i,j) then output is defined at (i,j-1/2). BDY_EU_smag changes dimensions by m^-3


    INPUT
    a    field defined on the north face of a U-cell
       [real, dimension(isd:ied,jsd:jed)]

  7. bihgen_viscosity_check

    DESCRIPTION
    Subroutine to perform linear stability check for the biharmonic operator given a value for the horizontal biharmonic viscosity.


  8. bihgen_reynolds_check

    DESCRIPTION
    Subroutine to compute the biharmonic grid Reynolds number. Large Reynolds numbers indicate regions where solution may experience some grid noise due to lack of enough horizontal friction.


  9. compute_neptune_velocity

    DESCRIPTION
    Compute Neptune velocity. Method follows that of Maltrud and Holloway, 2008: Implementing biharmonic neptune in a global eddying ocean model, Ocean Modelling, vol. 21, pages 22-34. March 2012 Stephen.Griffies



NAMELIST

&ocean_bihgen_friction_nml

use_this_module
Must be true to use this module. Default is false.
[logical]
debug_this_module
For debugging by printing checksums.
[logical]
k_smag_iso
This is the dimensionless Smagorinsky coefficient used to set the scale of the Smagorinsky isotropic viscosity.
[real, units: dimensionless]
k_smag_aniso
This is the dimensionless Smagorinsky coefficient used to set the scale of the Smagorinsky anisotropic viscosity.
[real, units: dimensionless]
vel_micom_iso
Velocity scale that is used for computing the MICOM isotropic viscosity.
[real, units: m/sec]
vel_micom_aniso
Velocity scale that is used for computing the MICOM anisotropic viscosity.
[real, units: m/sec]
visc_crit_scale
Scaling factor used to determine the critical viscosity, above which the viscosity is not allowed to reach. Use visc_crit_scale < 1.0 for cases where the visc_crit from linear stability allows for still too large of a viscosity. Use visc_crit_scale>1.0 when wish to allow for larger viscosity. Default is visc_crit_scale=1.0.
[real, units: dimensionless]
equatorial_zonal
Orient the anisotropic friction within a latitudinal band according to zonal direction.
[real]
equatorial_zonal_lat
Latitudinal band to use the zonal friction orientation.
[real]
eq_vel_micom_iso
Velocity scale that is used for computing the MICOM isotropic viscosity within a user specified equatorial band.
[real]
eq_vel_micom_aniso
Velocity scale that is used for computing the MICOM anisotropic viscosity within a user specified equatorial band.
[real]
eq_lat_micom
Equatorial latitude band (degrees) within which the MICOM viscosity is set according to eq_vel_micom_iso and eq_vel_micom_aniso.
[real]
bottom_5point
To alleviate problems with small partial cells, it is often necessary to reduce the operator to the traditional 5-point Laplacian at the ocean bottom. This logical implements this mixing. Default bottom_5point=.false.
[logical]
vel_micom_bottom
Velocity scale that is used for computing the MICOM viscosity for 5point Laplacian at the bottom.
[real, units: m/sec]
ncar_boundary_scaling
To enhance the velocity scale used in western boundaries for the isotropic and anisotropic background viscosities, we compute a scaling using the algorithm from the laplacian NCAR anisotropic scheme. Default ncar_boundary_scaling=.false.
[logical]
ncar_boundary_scaling_read
To read in the ncar boundary scaling field rather than generating it during initialization. Generating during initialization can be a bottle-neck on fine resolution models since there are some global 2d fields needed. So if the rescaling is produced once and then saved, it can be read in during subsequent runs without incurring the slowdown of re-generating the scalings. Default ncar_boundary_scaling_read=.false.
[logical]
ncar_rescale_power
For determining rescaling of the viscosity so to enhance the friction near the western boundaries. Default ncar_rescale_power=1.
[integer, units: dimensionless]
ncar_vconst_4
Inverse damping length for exponential falloff of the velocity scale as move eastward away from western boundary. Default ncar_vconst_4=2.e-8.
[real, units: 1/cm]
ncar_vconst_5
For determining number of grid points in boundary calculation. Default ncar_vconst_5=3.
[integer, units: dimensionless]
neptune
Set to true for computing friction relative to Neptune barotropic velocity. Default neptune=.false.
[logical]
neptune_length_eq
Length scale used to compute Neptune velocity at equator. Default neptune_length_eq= 4.2e3 from Maltrud and Holloway.
[real, units: m]
neptune_length_pole
Length scale used to compute Neptune velocity at pole. Default neptune_length_pole= 17.0e3 from Maltrud and Holloway.
[real, units: m]
neptune_depth_min
Minimum depth scale used for computing Neptune velocity. Default neptune_depth_min=100.0 from Maltrud and Holloway.
[real, units: m]
neptune_smooth
For doing a horizontal 1-2-1 smoothing on the diagnosed neptune velocity scale. Default neptune_smooth=.true.
[logical]
neptune_smooth_num
Number of smoothing passes for neptune velocity. Default neptune_smooth_num=1.
[integer]
neptune_scaling
Overall scaling parameter to help tune neptune. Default neptune_scaling=1.0
[real]
visc_diverge_scaling
Dimensionless scaling used for divergence based viscosity. Default visc_diverge_scaling=0.0 turns off the scheme. visc_diverge_scaling=10.0 produces sensible viscosities for 1-degree model. May need tuning for different resolutions.
[real]
use_side_drag_friction
For converting friction at U-cells next to walls into a drag law, as per Deremble et al. Use cdbot_array from ocean_core/ocean_bbc.F90 to compute drag force. Default use_side_drag_friction=.false.
[logical]
side_drag_friction_scaling
Dimensionless scaling used for cdbot_array when setting side drag friction. So the effective side dragy coefficient is side_drag_friction_scaling*cdbot_array. Default side_drag_friction_scaling=1.0.
[real]
side_drag_friction_uvmag_max
Maximum magnitude of horizontal velocity used to compute the side drag friction. This parameter can be useful especially for pressure models where the bottom cells can be quite thin and subject to sporadic large magnitudes. We do the same thing with bottom drag calculations. Default side_drag_friction_uvmag_max=10.0.
[real, units: m/s]
side_drag_friction_max
Maximum magnitude of the side drag induced friction. This parameter can be useful especially for pressure models where the bottom cells can be quite thin and subject to sporadic large magnitudes. We do the same thing with bottom drag calculations. Default side_drag_friction_max=1.0.
[real, units: N/m^2]


REFERENCES

  1. S.M. Griffies and R.W. Hallberg, 2000: Biharmonic friction with a Smagorinsky viscosity for use in large-scale eddy-permitting ocean models Monthly Weather Review, vol. 128, pages 2935-2946.
  2. R.D. Smith and J.C. McWilliams, 2003: Anisotropic horizontal viscosity for ocean models, Ocean Modelling, vol. 5, pages 129-156.
  3. Maltrud and Holloway, 2008: Implementing biharmonic neptune in a global eddying ocean model, Ocean Modelling, vol. 21, pages 22-34.
  4. Deremble, Hogg, Berloff, and Dewar, 2011: On the application of no-slip lateral boundary conditions to coarsely resolved ocean models, Ocean Modelling.
  5. S.M. Griffies, 2004: Fundamentals of Ocean Climate Models Princeton University Press
  6. Griffies, Elements of MOM (2012)


NOTES

The ocean model can generally run with both Laplacian and biharmonic friction enabled at the same time. Such has been found useful for some eddying ocean simulations.


top