coord_rho module reference

Regrid columns for the continuous isopycnal (rho) coordinate.

More…

Data Types

rho_cs

Control structure containing required parameters for the rho coordinate.

Functions/Subroutines

init_coord_rho()

Initialise a rho_CS with pointers to parameters.

end_coord_rho()

This subroutine deallocates memory in the control structure for the coord_rho() module.

set_rho_params()

This subroutine can be used to set the parameters for the coord_rho() module.

build_rho_column()

Build a rho coordinate column.

build_rho_column_iteratively()

Iteratively build a rho coordinate column.

copy_finite_thicknesses()

Copy column thicknesses with vanished layers removed.

old_inflate_layers_1d()

Inflate vanished layers to finite (nonzero) width.

Detailed Description

Regrid columns for the continuous isopycnal (rho) coordinate.

Type Documentation

type coord_rho/rho_cs

Control structure containing required parameters for the rho coordinate.

Type fields:
  • % nk [integer] :: Number of layers.

  • % min_thickness [real] :: Minimum thickness allowed for layers, often in [H ~> m or kg m-2].

  • % ref_pressure [real] :: Reference pressure for density calculations [R L2 T-2 ~> Pa].

  • % integrate_downward_for_e [logical] :: If true, integrate for interface positions from the top downward. If false, integrate from the bottom upward, as does the rest of the model.

  • % target_density [real(:),allocatable] :: Nominal density of interfaces [R ~> kg m-3].

  • % interp_cs [type( interp_cs_type )] :: Interpolation control structure.

Function/Subroutine Documentation

subroutine coord_rho/init_coord_rho(CS, nk, ref_pressure, target_density, interp_CS)

Initialise a rho_CS with pointers to parameters.

Parameters:
  • cs :: Unassociated pointer to hold the control structure

  • nk :: [in] Number of layers in the grid

  • ref_pressure :: [in] Coordinate reference pressure [R L2 T-2 ~> Pa]

  • target_density :: [in] Nominal density of interfaces [R ~> kg m-3]

  • interp_cs :: [in] Controls for interpolation

Call to:

mom_error_handler::mom_error

subroutine coord_rho/end_coord_rho(CS)

This subroutine deallocates memory in the control structure for the coord_rho() module. module.

Parameters:

cs :: Coordinate control structure

subroutine coord_rho/set_rho_params(CS, min_thickness, integrate_downward_for_e, interp_CS, ref_pressure)

This subroutine can be used to set the parameters for the coord_rho() module. module.

Parameters:
  • cs :: Coordinate control structure

  • min_thickness :: [in] Minimum allowed thickness [H ~> m or kg m-2]

  • integrate_downward_for_e :: [in] If true, integrate for interface positions from the top downward. If false, integrate from the bottom upward, as does the rest of the model.

  • ref_pressure :: [in] The reference pressure for density-dependent coordinates [R L2 T-2 ~> Pa]

  • interp_cs :: [in] Controls for interpolation

Call to:

mom_error_handler::mom_error

subroutine coord_rho/build_rho_column(CS, nz, depth, h, T, S, eqn_of_state, z_interface, z_rigid_top, eta_orig, h_neglect, h_neglect_edge)

Build a rho coordinate column.

  1. Density profiles are calculated on the source grid.

  2. Positions of target densities (for interfaces) are found by interpolation.

Parameters:
  • cs :: [in] coord_rho() control structure control structure

  • nz :: [in] Number of levels on source grid (i.e. length of h, T, S)

  • depth :: [in] Depth of ocean bottom (positive downward) [H ~> m or kg m-2]

  • h :: [in] Layer thicknesses [H ~> m or kg m-2]

  • t :: [in] Temperature for source column [C ~> degC]

  • s :: [in] Salinity for source column [S ~> ppt]

  • eqn_of_state :: [in] Equation of state structure

  • z_interface :: [inout] Absolute positions of interfaces [H ~> m or kg m-2]

  • z_rigid_top :: [in] The height of a rigid top (positive upward in the same units as depth) [H ~> m or kg m-2]

  • eta_orig :: [in] The actual original height of the top in the same units as depth) [H ~> m or kg m-2]

  • h_neglect :: [in] A negligibly small width for the purpose of cell reconstructions [H ~> m or kg m-2]

  • h_neglect_edge :: [in] A negligibly small width for the purpose of edge value calculations [H ~> m or kg m-2]

Call to:

regrid_interp::build_and_interpolate_grid copy_finite_thicknesses old_inflate_layers_1d

Called from:

mom_diag_remap::diag_remap_update

subroutine coord_rho/build_rho_column_iteratively(CS, remapCS, nz, depth, h, T, S, eqn_of_state, zInterface, h_neglect, h_neglect_edge, dev_tol)

Iteratively build a rho coordinate column.

The algorithm operates as follows within each column:

  1. Given T & S within each layer, the layer densities are computed.

  2. Based on these layer densities, a global density profile is reconstructed (this profile is monotonically increasing and may be discontinuous)

  3. The new grid interfaces are determined based on the target interface densities.

  4. T & S are remapped onto the new grid.

  5. Return to step 1 until convergence or until the maximum number of iterations is reached, whichever comes first.

Parameters:
  • cs :: [in] Regridding control structure

  • remapcs :: [in] Remapping parameters and options

  • nz :: [in] Number of levels

  • depth :: [in] Depth of ocean bottom [Z ~> m]

  • h :: [in] Layer thicknesses in Z coordinates [Z ~> m]

  • t :: [in] T for column [C ~> degC]

  • s :: [in] S for column [S ~> ppt]

  • eqn_of_state :: [in] Equation of state structure

  • zinterface :: [inout] Absolute positions of interfaces [Z ~> m]

  • h_neglect :: [in] A negligibly small width for the purpose of cell reconstructions in the same units as h [Z ~> m]

  • h_neglect_edge :: [in] A negligibly small width for the purpose of edge value calculations in the same units as h [Z ~> m]

  • dev_tol :: [in] The tolerance for the deviation between successive grids for determining when the iterative solver has converged [Z ~> m]

Call to:

regrid_interp::build_and_interpolate_grid copy_finite_thicknesses old_inflate_layers_1d mom_remapping::remapping_core_h

subroutine coord_rho/copy_finite_thicknesses(nk, h_in, thresh, nout, h_out, mapping)

Copy column thicknesses with vanished layers removed.

Parameters:
  • nk :: [in] Number of layer for h_in, T_in, S_in

  • h_in :: [in] Thickness of input column [H ~> m or kg m-2] or [Z ~> m]

  • thresh :: [in] Thickness threshold defining vanished layers [H ~> m or kg m-2] or [Z ~> m]

  • nout :: [out] Number of non-vanished layers

  • h_out :: [out] Thickness of output column [H ~> m or kg m-2] or [Z ~> m]

  • mapping :: [out] Index of k-out corresponding to k-in

Called from:

build_rho_column build_rho_column_iteratively

subroutine coord_rho/old_inflate_layers_1d(min_thickness, nk, h)

Inflate vanished layers to finite (nonzero) width.

Parameters:
  • min_thickness :: [in] Minimum allowed thickness [H ~> m or kg m-2] or other units

  • nk :: [in] Number of layers in the grid

  • h :: [inout] Layer thicknesses [H ~> m or kg m-2] or other units

Called from:

build_rho_column build_rho_column_iteratively mom_regridding::inflate_vanished_layers_old