mom_remapping module reference

Provides column-wise vertical remapping functions.

More…

Data Types

remapping_cs

Container for remapping parameters.

Functions/Subroutines

remapping_set_param()

Set parameters within remapping object.

extract_member_remapping_cs()

buildgridfromh()

Calculate edge coordinate x from cell width h.

remapping_core_h()

Remaps column of values u0 on grid h0 to grid h1 assuming the top edge is aligned.

remapping_core_w()

Remaps column of values u0 on grid h0 to implied grid h1 where the interfaces of h1 differ from those of h0 by dx.

build_reconstructions_1d()

Creates polynomial reconstructions of u0 on the source grid h0.

check_reconstructions_1d()

Checks that edge values and reconstructions satisfy bounds.

remap_via_sub_cells()

Remaps column of n0 values u0 on grid h0 to grid h1 with n1 cells by calculating the n0+n1+1 sub-integrals of the intersection of h0 and h1, and the summing the appropriate integrals into the h1*u1 values.

interpolate_column()

Linearly interpolate interface data, u_src, from grid h_src to a grid h_dest.

reintegrate_column()

Conservatively calculate integrated data, uh_dest, on grid h_dest, from layer-integrated data, uh_src, on grid h_src.

average_value_ppoly()

Returns the average value of a reconstruction within a single source cell, i0, between the non-dimensional positions xa and xb (xa<=xb) with dimensional separation dh.

check_remapped_values()

This subroutine checks for sufficient consistence in the extrema and total amounts on the old and new grids.

measure_input_bounds()

Measure totals and bounds on source grid.

measure_output_bounds()

Measure totals and bounds on destination grid.

dzfromh1h2()

Calculates the change in interface positions based on h1 and h2.

initialize_remapping()

Constructor for remapping control structure.

setreconstructiontype()

Changes the method of reconstruction Use this routine to parse a string parameter specifying the reconstruction and re-allocates work arrays appropriately.

end_remapping()

Destrcutor for remapping control structure.

remapping_unit_tests()

Runs unit tests on remapping functions.

test_answer()

Returns true if any cell of u and u_true are not identical.

test_interp()

Returns true if a test of interpolate_column() produces the wrong answer.

test_reintegrate()

Returns true if a test of reintegrate_column() produces the wrong answer.

dumpgrid()

Convenience function for printing grid to screen.

Detailed Description

Provides column-wise vertical remapping functions.

Type Documentation

type mom_remapping/remapping_cs

Container for remapping parameters.

Type fields:
  • % remapping_scheme [integer] :: Determines which reconstruction to use.

  • % degree [integer] :: Degree of polynomial reconstruction.

  • % boundary_extrapolation [logical] :: If true, extrapolate boundaries.

  • % check_reconstruction [logical] :: If true, reconstructions are checked for consistency.

  • % check_remapping [logical] :: If true, the result of remapping are checked for conservation and bounds.

  • % force_bounds_in_subcell [logical] :: If true, the intermediate values used in remapping are forced to be bounded.

  • % answer_date [integer] :: The vintage of the expressions to use for remapping. Values below 20190101 result in the use of older, less accurate expressions.

Function/Subroutine Documentation

subroutine mom_remapping/remapping_set_param(CS, remapping_scheme, boundary_extrapolation, check_reconstruction, check_remapping, force_bounds_in_subcell, answers_2018, answer_date)

Set parameters within remapping object.

Parameters:
  • cs :: [inout] Remapping control structure

  • remapping_scheme :: [in] Remapping scheme to use

  • boundary_extrapolation :: [in] Indicate to extrapolate in boundary cells

  • check_reconstruction :: [in] Indicate to check reconstructions

  • check_remapping :: [in] Indicate to check results of remapping

  • force_bounds_in_subcell :: [in] Force subcells values to be bounded

  • answers_2018 :: [in] If true use older, less accurate expressions.

  • answer_date :: [in] The vintage of the expressions to use

Call to:

setreconstructiontype

Called from:

initialize_remapping

subroutine mom_remapping/extract_member_remapping_cs(CS, remapping_scheme, degree, boundary_extrapolation, check_reconstruction, check_remapping, force_bounds_in_subcell)
Parameters:
  • cs :: [in] Control structure for remapping module

  • remapping_scheme :: [out] Determines which reconstruction scheme to use

  • degree :: [out] Degree of polynomial reconstruction

  • boundary_extrapolation :: [out] If true, extrapolate boundaries

  • check_reconstruction :: [out] If true, reconstructions are checked for consistency.

  • check_remapping :: [out] If true, the result of remapping are checked for conservation and bounds.

  • force_bounds_in_subcell :: [out] If true, the intermediate values used in remapping are forced to be bounded.

subroutine mom_remapping/buildgridfromh(nz, h, x)

Calculate edge coordinate x from cell width h.

Parameters:
  • nz :: [in] Number of cells

  • h :: [in] Cell widths [H]

  • x :: [inout] Edge coordinates starting at x(1)=0 [H]

Called from:

remapping_unit_tests

subroutine mom_remapping/remapping_core_h(CS, n0, h0, u0, n1, h1, u1, h_neglect, h_neglect_edge, PCM_cell)

Remaps column of values u0 on grid h0 to grid h1 assuming the top edge is aligned.

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

  • n0 :: [in] Number of cells on source grid

  • h0 :: [in] Cell widths on source grid [H]

  • u0 :: [in] Cell averages on source grid [A]

  • n1 :: [in] Number of cells on target grid

  • h1 :: [in] Cell widths on target grid [H]

  • u1 :: [out] Cell averages on target grid [A]

  • h_neglect :: [in] A negligibly small width for the purpose of cell reconstructions in the same units as h0 [H]

  • h_neglect_edge :: [in] A negligibly small width for the purpose of edge value calculations in the same units as h0 [H]

  • pcm_cell :: [in] If present, use PCM remapping for cells in the source grid where this is true.

Call to:

build_reconstructions_1d check_reconstructions_1d check_remapped_values remap_via_sub_cells

Called from:

mom_ale_sponge::apply_ale_sponge mom_oda_incupd::apply_oda_incupd mom_oda_driver_mod::apply_oda_tracer_increments coord_hycom::build_hycom1_column coord_rho::build_rho_column_iteratively mom_oda_incupd::calc_oda_increments mom_tidal_mixing::calculate_cvmix_tidal mom_state_initialization::cut_off_column_top mom_oda_driver_mod::set_prior_tracer mom_wave_speed::wave_speed mom_wave_speed::wave_speeds

subroutine mom_remapping/remapping_core_w(CS, n0, h0, u0, n1, dx, u1, h_neglect, h_neglect_edge)

Remaps column of values u0 on grid h0 to implied grid h1 where the interfaces of h1 differ from those of h0 by dx.

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

  • n0 :: [in] Number of cells on source grid

  • h0 :: [in] Cell widths on source grid [H]

  • u0 :: [in] Cell averages on source grid [A]

  • n1 :: [in] Number of cells on target grid

  • dx :: [in] Cell widths on target grid [H]

  • u1 :: [out] Cell averages on target grid [A]

  • h_neglect :: [in] A negligibly small width for the purpose of cell reconstructions in the same units as h0 [H].

  • h_neglect_edge :: [in] A negligibly small width for the purpose of edge value calculations in the same units as h0 [H].

Call to:

build_reconstructions_1d check_reconstructions_1d check_remapped_values remap_via_sub_cells

Called from:

remapping_unit_tests

subroutine mom_remapping/build_reconstructions_1d(CS, n0, h0, u0, ppoly_r_coefs, ppoly_r_E, ppoly_r_S, iMethod, h_neglect, h_neglect_edge, PCM_cell)

Creates polynomial reconstructions of u0 on the source grid h0.

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

  • n0 :: [in] Number of cells on source grid

  • h0 :: [in] Cell widths on source grid [H]

  • u0 :: [in] Cell averages on source grid [A]

  • ppoly_r_coefs :: [out] Coefficients of polynomial [A]

  • ppoly_r_e :: [out] Edge value of polynomial [A]

  • ppoly_r_s :: [out] Edge slope of polynomial [A H-1]

  • imethod :: [out] Integration method

  • h_neglect :: [in] A negligibly small width for the purpose of cell reconstructions in the same units as h0 [H]

  • h_neglect_edge :: [in] A negligibly small width for the purpose of edge value calculations in the same units as h0 [H]

  • pcm_cell :: [in] If present, use PCM remapping for cells from the source grid where this is true.

Call to:

regrid_edge_values::edge_slopes_implicit_h3 regrid_edge_values::edge_slopes_implicit_h5 mom_hybgen_remap::hybgen_plm_coefs mom_hybgen_remap::hybgen_ppm_coefs mom_hybgen_remap::hybgen_weno_coefs integration_pcm integration_plm integration_ppm integration_pqm mom_error_handler::mom_error pcm_functions::pcm_reconstruction plm_functions::plm_boundary_extrapolation plm_functions::plm_reconstruction ppm_functions::ppm_monotonicity pqm_functions::pqm_boundary_extrapolation_v1 pqm_functions::pqm_reconstruction remapping_pcm remapping_plm remapping_plm_hybgen remapping_ppm_cw remapping_ppm_h4 remapping_ppm_hybgen remapping_ppm_ih4 remapping_pqm_ih4ih3 remapping_pqm_ih6ih5 remapping_weno_hybgen

Called from:

remapping_core_h remapping_core_w

subroutine mom_remapping/check_reconstructions_1d(n0, h0, u0, deg, boundary_extrapolation, ppoly_r_coefs, ppoly_r_E, ppoly_r_S)

Checks that edge values and reconstructions satisfy bounds.

Parameters:
  • n0 :: [in] Number of cells on source grid

  • h0 :: [in] Cell widths on source grid [H]

  • u0 :: [in] Cell averages on source grid [A]

  • deg :: [in] Degree of polynomial reconstruction

  • boundary_extrapolation :: [in] Extrapolate at boundaries if true

  • ppoly_r_coefs :: [out] Coefficients of polynomial [A]

  • ppoly_r_e :: [out] Edge value of polynomial [A]

  • ppoly_r_s :: [out] Edge slope of polynomial [A H-1]

Call to:

mom_error_handler::mom_error

Called from:

remapping_core_h remapping_core_w

subroutine mom_remapping/remap_via_sub_cells(n0, h0, u0, ppoly0_E, ppoly0_coefs, n1, h1, method, force_bounds_in_subcell, u1, uh_err, ah_sub, aisub_src, aiss, aise)

Remaps column of n0 values u0 on grid h0 to grid h1 with n1 cells by calculating the n0+n1+1 sub-integrals of the intersection of h0 and h1, and the summing the appropriate integrals into the h1*u1 values. h0 and h1 must have the same units.

Parameters:
  • n0 :: [in] Number of cells in source grid

  • h0 :: [in] Source grid widths (size n0) [H]

  • u0 :: [in] Source cell averages (size n0) [A]

  • ppoly0_e :: [in] Edge value of polynomial [A]

  • ppoly0_coefs :: [in] Coefficients of polynomial [A]

  • n1 :: [in] Number of cells in target grid

  • h1 :: [in] Target grid widths (size n1) [H]

  • method :: [in] Remapping scheme to use

  • force_bounds_in_subcell :: [in] Force sub-cell values to be bounded

  • u1 :: [out] Target cell averages (size n1) [A]

  • uh_err :: [out] Estimate of bound on error in sum of u*h [A H]

  • ah_sub :: [out] Overlapping sub-cell thicknesses, h_sub [H]

  • aisub_src :: [out] i_sub_src

  • aiss :: [out] isrc_start

  • aise :: [out] isrc_ens

Call to:

average_value_ppoly measure_input_bounds measure_output_bounds mom_error_handler::mom_error

Called from:

remapping_core_h remapping_core_w remapping_unit_tests

subroutine mom_remapping/interpolate_column(nsrc, h_src, u_src, ndest, h_dest, u_dest, mask_edges)

Linearly interpolate interface data, u_src, from grid h_src to a grid h_dest.

Parameters:
  • nsrc :: [in] Number of source cells

  • h_src :: [in] Thickness of source cells [H]

  • u_src :: [in] Values at source cell interfaces [A]

  • ndest :: [in] Number of destination cells

  • h_dest :: [in] Thickness of destination cells [H]

  • u_dest :: [inout] Interpolated value at destination cell interfaces [A]

  • mask_edges :: [in] If true, mask the values outside of massless layers at the top and bottom of the column.

Called from:

test_interp mom_diag_remap::vertically_interpolate_field mom_wave_speed::wave_speeds

subroutine mom_remapping/reintegrate_column(nsrc, h_src, uh_src, ndest, h_dest, uh_dest)

Conservatively calculate integrated data, uh_dest, on grid h_dest, from layer-integrated data, uh_src, on grid h_src.

Parameters:
  • nsrc :: [in] Number of source cells

  • h_src :: [in] Thickness of source cells [H]

  • uh_src :: [in] Values at source cell interfaces [A H]

  • ndest :: [in] Number of destination cells

  • h_dest :: [in] Thickness of destination cells [H]

  • uh_dest :: [inout] Interpolated value at destination cell interfaces [A H]

Called from:

test_reintegrate mom_diag_remap::vertically_reintegrate_field

function mom_remapping/average_value_ppoly(n0, u0, ppoly0_E, ppoly0_coefs, method, i0, xa, xb) [real]

Returns the average value of a reconstruction within a single source cell, i0, between the non-dimensional positions xa and xb (xa<=xb) with dimensional separation dh.

Parameters:
  • n0 :: [in] Number of cells in source grid

  • u0 :: [in] Cell means [A]

  • ppoly0_e :: [in] Edge value of polynomial [A]

  • ppoly0_coefs :: [in] Coefficients of polynomial [A]

  • method :: [in] Remapping scheme to use

  • i0 :: [in] Source cell index

  • xa :: [in] Non-dimensional start position within source cell [nondim]

  • xb :: [in] Non-dimensional end position within source cell [nondim]

Call to:

integration_pcm integration_plm integration_ppm integration_pqm mom_error_handler::mom_error

Called from:

mom_neutral_diffusion::neutral_surface_t_eval remap_via_sub_cells

subroutine mom_remapping/check_remapped_values(n0, h0, u0, ppoly_r_E, deg, ppoly_r_coefs, n1, h1, u1, iMethod, uh_err, caller)

This subroutine checks for sufficient consistence in the extrema and total amounts on the old and new grids.

Parameters:
  • n0 :: [in] Number of cells on source grid

  • h0 :: [in] Cell widths on source grid [H]

  • u0 :: [in] Cell averages on source grid [A]

  • ppoly_r_e :: [in] Edge values of polynomial fits [A]

  • deg :: [in] Degree of the piecewise polynomial reconstrution

  • ppoly_r_coefs :: [in] Coefficients of the piecewise polynomial reconstructions [A]

  • n1 :: [in] Number of cells on target grid

  • h1 :: [in] Cell widths on target grid [H]

  • u1 :: [in] Cell averages on target grid [A]

  • imethod :: [in] An integer indicating the integration method used

  • uh_err :: [in] A bound on the error in the sum of u*h as estimated by the remapping code [H A]

  • caller :: [in] The name of the calling routine.

Call to:

measure_input_bounds measure_output_bounds mom_error_handler::mom_error

Called from:

remapping_core_h remapping_core_w

subroutine mom_remapping/measure_input_bounds(n0, h0, u0, edge_values, h0tot, h0err, u0tot, u0err, u0min, u0max)

Measure totals and bounds on source grid.

Parameters:
  • n0 :: [in] Number of cells on source grid

  • h0 :: [in] Cell widths on source grid [H]

  • u0 :: [in] Cell averages on source grid [A]

  • edge_values :: [in] Cell edge values on source grid [A]

  • h0tot :: [out] Sum of cell widths [H]

  • h0err :: [out] Magnitude of round-off error in h0tot [H]

  • u0tot :: [out] Sum of cell widths times values [H A]

  • u0err :: [out] Magnitude of round-off error in u0tot [H A]

  • u0min :: [out] Minimum value in reconstructions of u0 [A]

  • u0max :: [out] Maximum value in reconstructions of u0 [A]

Called from:

check_remapped_values remap_via_sub_cells

subroutine mom_remapping/measure_output_bounds(n1, h1, u1, h1tot, h1err, u1tot, u1err, u1min, u1max)

Measure totals and bounds on destination grid.

Parameters:
  • n1 :: [in] Number of cells on destination grid

  • h1 :: [in] Cell widths on destination grid [H]

  • u1 :: [in] Cell averages on destination grid [A]

  • h1tot :: [out] Sum of cell widths [H]

  • h1err :: [out] Magnitude of round-off error in h1tot [H]

  • u1tot :: [out] Sum of cell widths times values [H A]

  • u1err :: [out] Magnitude of round-off error in u1tot [H A]

  • u1min :: [out] Minimum value in reconstructions of u1 [A]

  • u1max :: [out] Maximum value in reconstructions of u1 [A]

Called from:

check_remapped_values remap_via_sub_cells

subroutine mom_remapping/dzfromh1h2(n1, h1, n2, h2, dx)

Calculates the change in interface positions based on h1 and h2.

Parameters:
  • n1 :: [in] Number of cells on source grid

  • h1 :: [in] Cell widths of source grid (size n1) [H]

  • n2 :: [in] Number of cells on target grid

  • h2 :: [in] Cell widths of target grid (size n2) [H]

  • dx :: [out] Change in interface position (size n2+1) [H]

Called from:

mom_ale::ale_remap_scalar remapping_unit_tests

subroutine mom_remapping/initialize_remapping(CS, remapping_scheme, boundary_extrapolation, check_reconstruction, check_remapping, force_bounds_in_subcell, answers_2018, answer_date)

Constructor for remapping control structure.

Parameters:
  • cs :: [inout] Remapping control structure

  • remapping_scheme :: [in] Remapping scheme to use

  • boundary_extrapolation :: [in] Indicate to extrapolate in boundary cells

  • check_reconstruction :: [in] Indicate to check reconstructions

  • check_remapping :: [in] Indicate to check results of remapping

  • force_bounds_in_subcell :: [in] Force subcells values to be bounded

  • answers_2018 :: [in] If true use older, less accurate expressions.

  • answer_date :: [in] The vintage of the expressions to use

Call to:

remapping_set_param

Called from:

mom_oda_driver_mod::init_oda mom_ale_sponge::initialize_ale_sponge::initialize_ale_sponge_fixed mom_ale_sponge::initialize_ale_sponge::initialize_ale_sponge_varying mom_oda_incupd::initialize_oda_incupd mom_tracer_initialization_from_z::mom_initialize_tracer_from_z mom_open_boundary::open_boundary_config mom_tidal_mixing::read_tidal_constituents remapping_unit_tests mom_state_initialization::trim_for_ice mom_wave_speed::wave_speed_init

subroutine mom_remapping/setreconstructiontype(string, CS)

Changes the method of reconstruction Use this routine to parse a string parameter specifying the reconstruction and re-allocates work arrays appropriately. It is called from initialize_remapping but can be called from an external module too.

Parameters:
  • string :: [in] String to parse for method

  • cs :: [inout] Remapping control structure

Call to:

mom_error_handler::mom_error remapping_pcm remapping_plm remapping_plm_hybgen remapping_ppm_cw remapping_ppm_h4 remapping_ppm_hybgen remapping_ppm_ih4 remapping_pqm_ih4ih3 remapping_pqm_ih6ih5 remapping_weno_hybgen mom_string_functions::uppercase

Called from:

remapping_set_param

subroutine mom_remapping/end_remapping(CS)

Destrcutor for remapping control structure.

Parameters:

cs :: [inout] Remapping control structure

function mom_remapping/remapping_unit_tests(verbose) [logical]

Runs unit tests on remapping functions. Should only be called from a single/root thread Returns True if a test fails, otherwise False.

Parameters:

verbose :: [in] If true, write results to stdout

Call to:

buildgridfromh dumpgrid dzfromh1h2 initialize_remapping integration_plm integration_ppm pcm_functions::pcm_reconstruction plm_functions::plm_reconstruction remap_via_sub_cells remapping_attic::remapping_attic_unit_tests remapping_core_w mom_io::stdout test_answer test_interp test_reintegrate

Called from:

mom_unit_tests::unit_tests

function mom_remapping/test_answer(verbose, n, u, u_true, label, tol) [logical]

Returns true if any cell of u and u_true are not identical. Returns false otherwise.

Parameters:
  • verbose :: [in] If true, write results to stdout

  • n :: [in] Number of cells in u

  • u :: [in] Values to test [A]

  • u_true :: [in] Values to test against (correct answer) [A]

  • label :: [in] Message

  • tol :: [in] The tolerance for differences between u and u_true [A]

Call to:

mom_io::stderr mom_io::stdout

Called from:

remapping_unit_tests

function mom_remapping/test_interp(verbose, msg, nsrc, h_src, u_src, ndest, h_dest, u_true) [logical]

Returns true if a test of interpolate_column() produces the wrong answer. produces the wrong answer.

Parameters:
  • verbose :: [in] If true, write results to stdout

  • msg :: [in] Message to label test

  • nsrc :: [in] Number of source cells

  • h_src :: [in] Thickness of source cells [H]

  • u_src :: [in] Values at source cell interfaces [A]

  • ndest :: [in] Number of destination cells

  • h_dest :: [in] Thickness of destination cells [H]

  • u_true :: [in] Correct value at destination cell interfaces [A]

Call to:

interpolate_column mom_io::stderr mom_io::stdout

Called from:

remapping_unit_tests

function mom_remapping/test_reintegrate(verbose, msg, nsrc, h_src, uh_src, ndest, h_dest, uh_true) [logical]

Returns true if a test of reintegrate_column() produces the wrong answer. produces the wrong answer.

Parameters:
  • verbose :: [in] If true, write results to stdout

  • msg :: [in] Message to label test

  • nsrc :: [in] Number of source cells

  • h_src :: [in] Thickness of source cells [H]

  • uh_src :: [in] Values of source cell stuff [A H]

  • ndest :: [in] Number of destination cells

  • h_dest :: [in] Thickness of destination cells [H]

  • uh_true :: [in] Correct value of destination cell stuff [A H]

Call to:

reintegrate_column mom_io::stderr mom_io::stdout

Called from:

remapping_unit_tests

subroutine mom_remapping/dumpgrid(n, h, x, u)

Convenience function for printing grid to screen.

Parameters:
  • n :: [in] Number of cells

  • h :: [in] Cell thickness [H]

  • x :: [in] Interface delta [H]

  • u :: [in] Cell average values [A]

Call to:

mom_io::stdout

Called from:

remapping_unit_tests