mom_regularize_layers module reference

Provides regularization of layers in isopycnal mode.

More…

Data Types

regularize_layers_cs

This control structure holds parameters used by the MOM_regularize_layers module.

Functions/Subroutines

regularize_layers()

This subroutine partially steps the bulk mixed layer model.

regularize_surface()

This subroutine ensures that there is a degree of horizontal smoothness in the depths of the near-surface interfaces.

find_deficit_ratios()

This subroutine determines the amount by which the harmonic mean thickness at velocity points differ from the arithmetic means, relative to the arithmetic means, after eliminating thickness variations that are solely due to topography and aggregating all interior layers into one.

regularize_layers_init()

Initializes the regularize_layers control structure.

Detailed Description

Provides regularization of layers in isopycnal mode.

Type Documentation

type mom_regularize_layers/regularize_layers_cs

This control structure holds parameters used by the MOM_regularize_layers module.

Type fields:
  • % initialized [logical] :: True if this control structure has been initialized.

  • % regularize_surface_layers [logical] :: If true, vertically restructure the near-surface layers when they have too much lateral variations to allow for sensible lateral barotropic transports.

  • % reg_sfc_detrain [logical] :: If true, allow the buffer layers to detrain into the interior as a part of the restructuring when regularize_surface_layers is true.

  • % density_match_tol [real] :: A relative tolerance for how well the densities must match with the target densities during detrainment when regularizing the near-surface layers [nondim].

  • % sufficient_adjustment [real] :: The fraction of the target entrainment of mass to the mixed and buffer layers that is enough for one timestep when regularizing the near-surface layers [nondim]. No more mass will be sought from deeper layers in the interior after this fraction is exceeded.

  • % h_def_tol1 [real] :: The value of the relative thickness deficit at which to start modifying the structure, 0.5 by default (or a thickness ratio of 5.83) [nondim].

  • % h_def_tol2 [real] :: The value of the relative thickness deficit at which to the structure modification is in full force, now 20% of the way from h_def_tol1 to 1 [nondim].

  • % h_def_tol3 [real] :: The value of the relative thickness deficit at which to start detrainment from the buffer layers to the interior, now 30% of the way from h_def_tol1 to 1 [nondim].

  • % h_def_tol4 [real] :: The value of the relative thickness deficit at which to do detrainment from the buffer layers to the interior at full force, now 50% of the way from h_def_tol1 to 1 [nondim].

  • % hmix_min [real] :: The minimum mixed layer thickness [H ~> m or kg m-2].

  • % time [type(time_type),pointer] :: A pointer to the ocean model’s clock.

  • % diag [type( diag_ctrl ),pointer] :: A structure that is used to regulate the timing of diagnostic output.

  • % answer_date [integer] :: The vintage of the order of arithmetic and expressions in this module’s calculations. Values below 20190101 recover the answers from the end of 2018, while higher values use updated and more robust forms of the same expressions.

  • % debug [logical] :: If true, do more thorough checks for debugging purposes.

  • % id_def_rat [integer] :: A diagnostic ID.

Function/Subroutine Documentation

subroutine mom_regularize_layers/regularize_layers(h, tv, dt, ea, eb, G, GV, US, CS)

This subroutine partially steps the bulk mixed layer model. The following processes are executed, in the order listed.

Parameters:
  • g :: [inout] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure.

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

  • tv :: [inout] A structure containing pointers to any available thermodynamic fields. Absent fields have NULL pointers.

  • dt :: [in] Time increment [T ~> s].

  • ea :: [inout] The amount of fluid moved downward into a

  • eb :: [inout] The amount of fluid moved upward into a layer

  • us :: [in] A dimensional unit scaling type

  • cs :: [in] Regularize layer control structure

Call to:

id_clock_pass mom_error_handler::mom_error regularize_surface

Called from:

mom_diabatic_driver::layered_diabatic

subroutine mom_regularize_layers/regularize_surface(h, tv, dt, ea, eb, G, GV, US, CS)

This subroutine ensures that there is a degree of horizontal smoothness in the depths of the near-surface interfaces.

Parameters:
  • g :: [inout] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure.

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

  • tv :: [inout] A structure containing pointers to any available thermodynamic fields. Absent fields have NULL pointers.

  • dt :: [in] Time increment [T ~> s].

  • ea :: [inout] The amount of fluid moved downward into a

  • eb :: [inout] The amount of fluid moved upward into a layer

  • us :: [in] A dimensional unit scaling type

  • cs :: [in] Regularize layer control structure

Call to:

mom_eos::eos_domain find_deficit_ratios mom_error_handler::mom_error

Called from:

regularize_layers

subroutine mom_regularize_layers/find_deficit_ratios(e, def_rat_u, def_rat_v, G, GV, CS, h)

This subroutine determines the amount by which the harmonic mean thickness at velocity points differ from the arithmetic means, relative to the arithmetic means, after eliminating thickness variations that are solely due to topography and aggregating all interior layers into one.

Parameters:
  • g :: [in] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure.

  • e :: [in] Interface depths [H ~> m or kg m-2]

  • def_rat_u :: [out] The thickness deficit ratio at u points,

  • def_rat_v :: [out] The thickness deficit ratio at v points,

  • cs :: [in] Regularize layer control structure

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

Called from:

regularize_surface

subroutine mom_regularize_layers/regularize_layers_init(Time, G, GV, param_file, diag, CS)

Initializes the regularize_layers control structure.

Parameters:
  • time :: [in] The current model time.

  • g :: [in] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure.

  • param_file :: [in] A structure to parse for run-time parameters.

  • diag :: [inout] A structure that is used to regulate diagnostic output.

  • cs :: [inout] Regularize layer control structure

Call to:

id_clock_pass

Called from:

mom_diabatic_driver::diabatic_driver_init