namespace mom_regridding

Overview

Generates vertical grids as part of the ALE algorithm. More…

namespace mom_regridding {

// global variables

character(len=*), parameter, public regriddingcoordinatemodedoc =                  " LAYER - Isopycnal or stacked shallow water layers\n"//                 " ZSTAR, Z* - stretched geopotential z*\n"//                 " SIGMA_SHELF_ZSTAR - stretched geopotential z* ignoring shelf\n"//                 " SIGMA - terrain following coordinates\n"//                 " RHO   - continuous isopycnal\n"//                 " HYCOM1 - HyCOM-like hybrid coordinate\n"//                 " SLIGHT - stretched coordinates above continuous isopycnal\n"//                 " ADAPTIVE - optimize for smooth neutral density surfaces";
character(len=*), parameter, public regriddinginterpschemedoc =                  " P1M_H2(2nd-order accurate)\n"//                 " P1M_H4(2nd-order accurate)\n"//                 " P1M_IH4(2nd-order accurate)\n"//                 " PLM(2nd-order accurate)\n"//                 " PPM_H4(3rd-order accurate)\n"//                 " PPM_IH4(3rd-order accurate)\n"//                 " P3M_IH4IH3(4th-order accurate)\n"//                 " P3M_IH6IH5(4th-order accurate)\n"//                 " PQM_IH4IH3(4th-order accurate)\n"//                 " PQM_IH6IH5(5th-order accurate)";
character(len=*), parameter, public regriddingdefaultinterpscheme = "P1M_H2";
logical, parameter, public regriddingdefaultboundaryextrapolation = .false.;
real, parameter, public regriddingdefaultminthickness = 1.e-3;

// global functions

subroutine, public initialize_regridding(
    CS CS,
    GV GV,
    US US,
    max_depth max_depth,
    param_file param_file,
    mdl mdl,
    coord_mode coord_mode,
    param_prefix param_prefix,
    param_suffix param_suffix
    );

subroutine, public end_regridding(CS CS);

subroutine, public regridding_main(
    remapCS remapCS,
    CS CS,
    G G,
    GV GV,
    h h,
    tv tv,
    h_new h_new,
    dzInterface dzInterface,
    frac_shelf_h frac_shelf_h,
    conv_adjust conv_adjust
    );

subroutine, public check_remapping_grid(G G, GV GV, h h, dzInterface dzInterface, msg msg);
subroutine, public check_grid_column(nk nk, depth depth, h h, dzInterface dzInterface, msg msg);
subroutine, public inflate_vanished_layers_old(CS CS, G G, GV GV, h h);

real function, dimension(nk), public uniformresolution(
    nk nk,
    coordMode coordMode,
    maxDepth maxDepth,
    rhoLight rhoLight,
    rhoHeavy rhoHeavy
    );

subroutine, public setcoordinateresolution(dz dz, CS CS, scale scale);
subroutine, public set_target_densities_from_gv(GV GV, CS CS);
subroutine, public set_target_densities(CS CS, rho_int rho_int);
subroutine, public set_regrid_max_depths(CS CS, max_depths max_depths, units_to_H units_to_H);
subroutine, public set_regrid_max_thickness(CS CS, max_h max_h, units_to_H units_to_H);
real function, dimension(cs%nk), public getcoordinateresolution(CS CS, undo_scaling undo_scaling);
real function, dimension(cs%nk+1), public getcoordinateinterfaces(CS CS, undo_scaling undo_scaling);
character(len=20) function, public getcoordinateunits(CS CS);
character(len=20) function, public getcoordinateshortname(CS CS);

subroutine, public set_regrid_params(
    CS CS,
    boundary_extrapolation boundary_extrapolation,
    min_thickness min_thickness,
    old_grid_weight old_grid_weight,
    interp_scheme interp_scheme,
    depth_of_time_filter_shallow depth_of_time_filter_shallow,
    depth_of_time_filter_deep depth_of_time_filter_deep,
    compress_fraction compress_fraction,
    dz_min_surface dz_min_surface,
    nz_fixed_surface nz_fixed_surface,
    Rho_ML_avg_depth Rho_ML_avg_depth,
    nlay_ML_to_interior nlay_ML_to_interior,
    fix_haloclines fix_haloclines,
    halocline_filt_len halocline_filt_len,
    halocline_strat_tol halocline_strat_tol,
    integrate_downward_for_e integrate_downward_for_e,
    adaptTimeRatio adaptTimeRatio,
    adaptZoom adaptZoom,
    adaptZoomCoeff adaptZoomCoeff,
    adaptBuoyCoeff adaptBuoyCoeff,
    adaptAlpha adaptAlpha,
    adaptDoMin adaptDoMin
    );

integer function, public get_regrid_size(CS CS);
type(zlike_cs) function, public get_zlike_cs(CS CS);
type(sigma_cs) function, public get_sigma_cs(CS CS);
type(rho_cs) function, public get_rho_cs(CS CS);
real function, dimension(cs%nk), public getstaticthickness(CS CS, SSH SSH, depth depth);

} // namespace mom_regridding

Detailed Documentation

Generates vertical grids as part of the ALE algorithm.

A vertical grid is defined solely by the cell thicknesses, \(h\). Most calculations in this module start with the coordinate at the bottom of the column set to -depth, and use a increasing value of coordinate with decreasing k. This is consistent with the rest of MOM6 that uses position, \(z\) which is a negative quantity for most of the ocean.

A change in grid is define through a change in position of the interfaces:

\[z^n_{k+1/2} = z^{n-1}_{k+1/2} + \Delta z_{k+1/2}\]

with the positive upward coordinate convention

\[z_{k-1/2} = z_{k+1/2} + h_k\]

so that

\[h^n_k = h^{n-1}_k + ( \Delta z_{k-1/2} - \Delta z_{k+1/2} )\]

Original date of creation: 2008.06.09 by L. White

Global Variables

character(len=*), parameter, public regriddingcoordinatemodedoc =                  " LAYER - Isopycnal or stacked shallow water layers\n"//                 " ZSTAR, Z* - stretched geopotential z*\n"//                 " SIGMA_SHELF_ZSTAR - stretched geopotential z* ignoring shelf\n"//                 " SIGMA - terrain following coordinates\n"//                 " RHO   - continuous isopycnal\n"//                 " HYCOM1 - HyCOM-like hybrid coordinate\n"//                 " SLIGHT - stretched coordinates above continuous isopycnal\n"//                 " ADAPTIVE - optimize for smooth neutral density surfaces"

Documentation for coordinate options.

character(len=*), parameter, public regriddinginterpschemedoc =                  " P1M_H2(2nd-order accurate)\n"//                 " P1M_H4(2nd-order accurate)\n"//                 " P1M_IH4(2nd-order accurate)\n"//                 " PLM(2nd-order accurate)\n"//                 " PPM_H4(3rd-order accurate)\n"//                 " PPM_IH4(3rd-order accurate)\n"//                 " P3M_IH4IH3(4th-order accurate)\n"//                 " P3M_IH6IH5(4th-order accurate)\n"//                 " PQM_IH4IH3(4th-order accurate)\n"//                 " PQM_IH6IH5(5th-order accurate)"

Documentation for regridding interpolation schemes.

character(len=*), parameter, public regriddingdefaultinterpscheme = "P1M_H2"

Default interpolation scheme.

logical, parameter, public regriddingdefaultboundaryextrapolation = .false.

Default mode for boundary extrapolation.

real, parameter, public regriddingdefaultminthickness = 1.e-3

Default minimum thickness for some coordinate generation modes.

Global Functions

subroutine, public initialize_regridding(
    CS CS,
    GV GV,
    US US,
    max_depth max_depth,
    param_file param_file,
    mdl mdl,
    coord_mode coord_mode,
    param_prefix param_prefix,
    param_suffix param_suffix
    )

Initialization and configures a regridding control structure based on customizable run-time parameters.

Parameters:

cs

Regridding control structure

gv

Ocean vertical grid structure

us

A dimensional unit scaling type

max_depth

The maximum depth of the ocean [Z ~> m].

param_file

Parameter file

mdl

Name of calling module.

coord_mode

Coordinate mode

param_prefix

String to prefix to parameter names. If empty, causes main model parameters to be used.

param_suffix

String to append to parameter names.

subroutine, public end_regridding(CS CS)

Deallocation of regridding memory.

Parameters:

cs

Regridding control structure

subroutine, public regridding_main(
    remapCS remapCS,
    CS CS,
    G G,
    GV GV,
    h h,
    tv tv,
    h_new h_new,
    dzInterface dzInterface,
    frac_shelf_h frac_shelf_h,
    conv_adjust conv_adjust
    )

Dispatching regridding routine for orchestrating regridding & remapping.

Parameters:

remapcs

Remapping parameters and options

cs

Regridding control structure

g

Ocean grid structure

gv

Ocean vertical grid structure

h

Current 3D grid obtained after the last time step

tv

Thermodynamical variables (T, S, …)

h_new

New 3D grid consistent with target coordinate

dzinterface

The change in position of each interface

frac_shelf_h

Fractional ice shelf coverage

conv_adjust

If true, do convective adjustment

subroutine, public check_remapping_grid(
    G G,
    GV GV,
    h h,
    dzInterface dzInterface,
    msg msg
    )

Check that the total thickness of two grids match.

Parameters:

g

Grid structure

gv

Ocean vertical grid structure

h

Layer thicknesses [H ~> m or kg m-2]

dzinterface

Change in interface positions [H ~> m or kg m-2]

msg

Message to append to errors

subroutine, public check_grid_column(
    nk nk,
    depth depth,
    h h,
    dzInterface dzInterface,
    msg msg
    )

Check that the total thickness of new and old grids are consistent.

Parameters:

nk

Number of cells

depth

Depth of bottom [Z ~> m] or arbitrary units

h

Cell thicknesses [Z ~> m] or arbitrary units

dzinterface

Change in interface positions (same units as h)

msg

Message to append to errors

subroutine, public inflate_vanished_layers_old(CS CS, G G, GV GV, h h)

Parameters:

cs

Regridding control structure

g

The ocean’s grid structure

gv

The ocean’s vertical grid structure

h

Layer thicknesses [H ~> m or kg m-2]

real function, dimension(nk), public uniformresolution(
    nk nk,
    coordMode coordMode,
    maxDepth maxDepth,
    rhoLight rhoLight,
    rhoHeavy rhoHeavy
    )

Return a uniform resolution vector in the units of the coordinate.

Parameters:

nk

Number of cells in source grid

coordmode

A string indicating the coordinate mode. See the documenttion for regrid_consts for the recognized values.

maxdepth

The range of the grid values in some modes

rholight

The minimum value of the grid in RHO mode

rhoheavy

The maximum value of the grid in RHO mode

Returns:

The returned uniform resolution grid.

subroutine, public setcoordinateresolution(dz dz, CS CS, scale scale)

Set the fixed resolution data.

Parameters:

dz

A vector of vertical grid spacings

cs

Regridding control structure

scale

A scaling factor converting dz to coordRes

subroutine, public set_target_densities_from_gv(GV GV, CS CS)

Set target densities based on the old Rlay variable.

Parameters:

gv

Ocean vertical grid structure

cs

Regridding control structure

subroutine, public set_target_densities(CS CS, rho_int rho_int)

Set target densities based on vector of interface values.

Parameters:

cs

Regridding control structure

rho_int

Interface densities

subroutine, public set_regrid_max_depths(
    CS CS,
    max_depths max_depths,
    units_to_H units_to_H
    )

Set maximum interface depths based on a vector of input values.

Parameters:

cs

Regridding control structure

max_depths

Maximum interface depths, in arbitrary units

units_to_h

A conversion factor for max_depths into H units

subroutine, public set_regrid_max_thickness(
    CS CS,
    max_h max_h,
    units_to_H units_to_H
    )

Set maximum layer thicknesses based on a vector of input values.

Parameters:

cs

Regridding control structure

max_h

Maximum interface depths, in arbitrary units

units_to_h

A conversion factor for max_h into H units

real function, dimension(cs%nk), public getcoordinateresolution(
    CS CS,
    undo_scaling undo_scaling
    )

Query the fixed resolution data.

Parameters:

cs

Regridding control structure

undo_scaling

If present and true, undo any internal rescaling of the resolution data.

real function, dimension(cs%nk+1), public getcoordinateinterfaces(
    CS CS,
    undo_scaling undo_scaling
    )

Query the target coordinate interface positions.

Parameters:

cs

Regridding control structure

undo_scaling

If present and true, undo any internal rescaling of the resolution data.

Returns:

Interface positions in target coordinate

character(len=20) function, public getcoordinateunits(CS CS)

Query the target coordinate units.

Parameters:

cs

Regridding control structure

character(len=20) function, public getcoordinateshortname(CS CS)

Query the short name of the coordinate.

Parameters:

cs

Regridding control structure

subroutine, public set_regrid_params(
    CS CS,
    boundary_extrapolation boundary_extrapolation,
    min_thickness min_thickness,
    old_grid_weight old_grid_weight,
    interp_scheme interp_scheme,
    depth_of_time_filter_shallow depth_of_time_filter_shallow,
    depth_of_time_filter_deep depth_of_time_filter_deep,
    compress_fraction compress_fraction,
    dz_min_surface dz_min_surface,
    nz_fixed_surface nz_fixed_surface,
    Rho_ML_avg_depth Rho_ML_avg_depth,
    nlay_ML_to_interior nlay_ML_to_interior,
    fix_haloclines fix_haloclines,
    halocline_filt_len halocline_filt_len,
    halocline_strat_tol halocline_strat_tol,
    integrate_downward_for_e integrate_downward_for_e,
    adaptTimeRatio adaptTimeRatio,
    adaptZoom adaptZoom,
    adaptZoomCoeff adaptZoomCoeff,
    adaptBuoyCoeff adaptBuoyCoeff,
    adaptAlpha adaptAlpha,
    adaptDoMin adaptDoMin
    )

Can be used to set any of the parameters for MOM_regridding.

Parameters:

cs

Regridding control structure

boundary_extrapolation

Extrapolate in boundary cells

min_thickness

Minimum thickness allowed when building the new grid [H ~> m or kg m-2]

old_grid_weight

Weight given to old coordinate when time-filtering grid

interp_scheme

Interpolation method for state-dependent coordinates

depth_of_time_filter_shallow

Depth to start cubic [H ~> m or kg m-2]

depth_of_time_filter_deep

Depth to end cubic [H ~> m or kg m-2]

compress_fraction

Fraction of compressibility to add to potential density

dz_min_surface

The fixed resolution in the topmost SLight_nkml_min layers [H ~> m or kg m-2]

nz_fixed_surface

The number of fixed-thickness layers at the top of the model

rho_ml_avg_depth

Averaging depth over which to determine mixed layer potential density [H ~> m or kg m-2]

nlay_ml_to_interior

Number of layers to offset the mixed layer density to find resolved stratification [nondim]

fix_haloclines

Detect regions with much weaker stratification in the coordinate

halocline_filt_len

Length scale over which to filter T & S when looking for spuriously unstable water mass profiles [m]

halocline_strat_tol

Value of the stratification ratio that defines a problematic halocline region.

integrate_downward_for_e

If true, integrate for interface positions downward from the top.

adapttimeratio

Ratio of the ALE timestep to the grid timescale [nondim].

adaptzoom

Depth of near-surface zooming region [H ~> m or kg m-2].

adaptzoomcoeff

Coefficient of near-surface zooming diffusivity [nondim].

adaptbuoycoeff

Coefficient of buoyancy diffusivity [nondim].

adaptalpha

Scaling factor on optimization tendency [nondim].

adaptdomin

If true, make a HyCOM-like mixed layer by preventing interfaces from being shallower than the depths specified by the regridding coordinate.

integer function, public get_regrid_size(CS CS)

Returns the number of levels/layers in the regridding control structure.

Parameters:

cs

Regridding control structure

type(zlike_cs) function, public get_zlike_cs(CS CS)

This returns a copy of the zlike_CS stored in the regridding control structure.

Parameters:

cs

Regridding control structure

type(sigma_cs) function, public get_sigma_cs(CS CS)

This returns a copy of the sigma_CS stored in the regridding control structure.

Parameters:

cs

Regridding control structure

type(rho_cs) function, public get_rho_cs(CS CS)

This returns a copy of the rho_CS stored in the regridding control structure.

Parameters:

cs

Regridding control structure

real function, dimension(cs%nk), public getstaticthickness(
    CS CS,
    SSH SSH,
    depth depth
    )

Return coordinate-derived thicknesses for fixed coordinate systems.

Parameters:

cs

Regridding control structure

ssh

The sea surface height, in the same units as depth

depth

The maximum depth of the grid, often [Z ~> m]

Returns:

The returned thicknesses in the units of depth