mom_spatial_means module reference

Functions and routines to take area, volume, mass-weighted, layerwise, zonal or meridional means.

More…

Functions/Subroutines

global_area_mean()

Return the global area mean of a variable.

global_area_integral()

Return the global area integral of a variable, by default using the masked area from the grid, but an alternate could be used instead.

global_layer_mean()

Return the layerwise global thickness-weighted mean of a variable.

global_volume_mean()

Find the global thickness-weighted mean of a variable.

global_mass_integral()

Find the global mass-weighted integral of a variable.

global_i_mean()

Determine the global mean of a field along rows of constant i, returning it in a 1-d array using the local indexing.

global_j_mean()

Determine the global mean of a field along rows of constant j, returning it in a 1-d array using the local indexing.

adjust_area_mean_to_zero()

Adjust 2d array such that area mean is zero without moving the zero contour.

Detailed Description

Functions and routines to take area, volume, mass-weighted, layerwise, zonal or meridional means.

Function/Subroutine Documentation

function mom_spatial_means/global_area_mean(var, G, scale) [real]

Return the global area mean of a variable. This uses reproducing sums.

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

  • var :: [in] The variable to average

  • scale :: [in] A rescaling factor for the variable

Called from

mom_forcing_type::forcing_diagnostics mom_generic_tracer::mom_generic_tracer_column_physics mom_generic_tracer::mom_generic_tracer_surface_state mom_internal_tides::sum_en

function mom_spatial_means/global_area_integral(var, G, scale, area) [real]

Return the global area integral of a variable, by default using the masked area from the grid, but an alternate could be used instead. This uses reproducing sums.

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

  • var :: [in] The variable to integrate

  • scale :: [in] A rescaling factor for the variable

  • area :: [in] The alternate area to use, including any required masking [L2 ~> m2].

Return

undefined :: The returned area integral, usually in the units of var times [m2].

Called from

mom_ice_shelf::add_shelf_flux mom_forcing_type::forcing_diagnostics mom_diagnostics::post_surface_thermo_diags

function mom_spatial_means/global_layer_mean(var, h, G, GV, scale) [real]

Return the layerwise global thickness-weighted mean of a variable. This uses reproducing sums.

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

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

  • var :: [in] The variable to average

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

  • scale :: [in] A rescaling factor for the variable

function mom_spatial_means/global_volume_mean(var, h, G, GV, scale) [real]

Find the global thickness-weighted mean of a variable. This uses reproducing sums.

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

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

  • var :: [in] The variable being averaged

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

  • scale :: [in] A rescaling factor for the variable

Return

undefined :: The thickness-weighted average of var

Called from

mom_diagnostics::calculate_diagnostic_fields

function mom_spatial_means/global_mass_integral(h, G, GV, var, on_PE_only, scale) [real]

Find the global mass-weighted integral of a variable. This uses reproducing sums.

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

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

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

  • var :: [in] The variable being integrated

  • on_pe_only :: [in] If present and true, the sum is only done on the local PE, and it is not order invariant.

  • scale :: [in] A rescaling factor for the variable

Return

undefined :: The mass-weighted integral of var (or 1) in kg times the units of var

subroutine mom_spatial_means/global_i_mean(array, i_mean, G, mask, scale, tmp_scale)

Determine the global mean of a field along rows of constant i, returning it in a 1-d array using the local indexing. This uses reproducing sums.

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

  • array :: [in] The variable being averaged

  • i_mean :: [out] Global mean of array along its i-axis

  • mask :: [in] An array used for weighting the i-mean

  • scale :: [in] A rescaling factor for the output variable

  • tmp_scale :: [in] A rescaling factor for the internal calculations that is removed from the output

Call to

mom_error_handler::mom_error mom_coms::query_efp_overflow_error mom_coms::reset_efp_overflow_error

Called from

mom_sponge::apply_sponge

subroutine mom_spatial_means/global_j_mean(array, j_mean, G, mask, scale, tmp_scale)

Determine the global mean of a field along rows of constant j, returning it in a 1-d array using the local indexing. This uses reproducing sums.

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

  • array :: [in] The variable being averaged

  • j_mean :: [out] Global mean of array along its j-axis

  • mask :: [in] An array used for weighting the j-mean

  • scale :: [in] A rescaling factor for the output variable

  • tmp_scale :: [in] A rescaling factor for the internal calculations that is removed from the output

Call to

mom_error_handler::mom_error mom_coms::query_efp_overflow_error mom_coms::reset_efp_overflow_error

subroutine mom_spatial_means/adjust_area_mean_to_zero(array, G, scaling, unit_scale)

Adjust 2d array such that area mean is zero without moving the zero contour.

Parameters
  • g :: [in] Grid structure

  • array :: [inout] 2D array to be adjusted

  • scaling :: [out] The scaling factor used

  • unit_scale :: [in] A rescaling factor for the variable

Called from

mom_surface_forcing_gfdl::convert_iob_to_fluxes