regrid_interp module reference

Vertical interpolation for regridding.

More…

Data Types

interp_cs_type

Control structure for regrid_interp() module.

Functions/Subroutines

regridding_set_ppolys()

Builds an interpolated profile for the densities within each grid cell.

interpolate_grid()

Given target values (e.g., density), build new grid based on polynomial.

build_and_interpolate_grid()

Build a grid by interpolating for target values.

get_polynomial_coordinate()

Given a target value, find corresponding coordinate for given polynomial.

interpolation_scheme()

Numeric value of interpolation_scheme corresponding to scheme name.

set_interp_scheme()

Store the interpolation_scheme value in the interp_CS based on the input string.

set_interp_extrap()

Store the boundary_extrapolation value in the interp_CS.

set_interp_answer_date()

Store the value of the answer_date in the interp_CS.

Detailed Description

Vertical interpolation for regridding.

Type Documentation

type regrid_interp/interp_cs_type

Control structure for regrid_interp() module. module.

Type fields:
  • % interpolation_scheme [integer] :: The following parameter is only relevant when used with the target interface densities regridding scheme. It indicates which interpolation to use to determine the grid.

  • % boundary_extrapolation [logical] :: Indicate whether high-order boundary extrapolation should be used within boundary cells.

  • % answer_date [integer] :: The vintage of the expressions to use for regridding.

Function/Subroutine Documentation

subroutine regrid_interp/regridding_set_ppolys(CS, densities, n0, h0, ppoly0_E, ppoly0_S, ppoly0_coefs, degree, h_neglect, h_neglect_edge)

Builds an interpolated profile for the densities within each grid cell.

It may happen that, given a high-order interpolator, the number of available layers is insufficient (e.g., there are two available layers for a third-order PPM ih4 scheme). In these cases, we resort to the simplest continuous linear scheme (P1M h2).

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

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

  • densities :: [in] Actual cell densities [A]

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

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

  • ppoly0_s :: [inout] Edge slope of polynomial [A H-1]

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

  • degree :: [inout] The degree of the polynomials

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

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

Call to:

degree_1 degree_2 degree_3 degree_4 regrid_edge_values::edge_slopes_implicit_h3 regrid_edge_values::edge_slopes_implicit_h5 interpolation_p1m_h2 interpolation_p1m_h4 interpolation_p1m_ih4 interpolation_p3m_ih4ih3 interpolation_p3m_ih6ih5 interpolation_plm interpolation_ppm_cw interpolation_ppm_h4 interpolation_ppm_ih4 interpolation_pqm_ih4ih3 interpolation_pqm_ih6ih5 p1m_functions::p1m_boundary_extrapolation p1m_functions::p1m_interpolation p3m_functions::p3m_boundary_extrapolation p3m_functions::p3m_interpolation plm_functions::plm_boundary_extrapolation plm_functions::plm_reconstruction ppm_functions::ppm_monotonicity pqm_functions::pqm_boundary_extrapolation_v1 pqm_functions::pqm_reconstruction

Called from:

build_and_interpolate_grid

subroutine regrid_interp/interpolate_grid(n0, h0, x0, ppoly0_E, ppoly0_coefs, target_values, degree, n1, h1, x1, answer_date)

Given target values (e.g., density), build new grid based on polynomial.

Given the grid ‘grid0’ and the piecewise polynomial interpolant ‘ppoly0’ (possibly discontinuous), the coordinates of the new grid ‘grid1’ are determined by finding the corresponding target interface densities.

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

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

  • h0 :: [in] Thicknesses of source grid cells [H]

  • x0 :: [in] Source interface positions [H]

  • ppoly0_e :: [in] Edge values of interpolating polynomials [A]

  • ppoly0_coefs :: [in] Coefficients of interpolating polynomials [A]

  • target_values :: [in] Target values of interfaces [A]

  • degree :: [in] Degree of interpolating polynomials

  • h1 :: [inout] Thicknesses of target grid cells [H]

  • x1 :: [inout] Target interface positions [H]

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

Call to:

get_polynomial_coordinate

Called from:

build_and_interpolate_grid

subroutine regrid_interp/build_and_interpolate_grid(CS, densities, n0, h0, x0, target_values, n1, h1, x1, h_neglect, h_neglect_edge)

Build a grid by interpolating for target values.

Parameters:
  • cs :: [in] A control structure for regrid_interp()

  • n0 :: [in] The number of points on the input grid

  • n1 :: [in] The number of points on the output grid

  • densities :: [in] Input cell densities [R ~> kg m-3]

  • target_values :: [in] Target values of interfaces [R ~> kg m-3]

  • h0 :: [in] Initial cell widths [H]

  • x0 :: [in] Source interface positions [H]

  • h1 :: [inout] Output cell widths [H]

  • x1 :: [inout] Target interface positions [H]

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

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

Call to:

interpolate_grid regridding_set_ppolys

Called from:

coord_rho::build_rho_column coord_rho::build_rho_column_iteratively

function regrid_interp/get_polynomial_coordinate(N, h, x_g, edge_values, ppoly_coefs, target_value, degree, answer_date) [real]

Given a target value, find corresponding coordinate for given polynomial.

Here, ‘ppoly’ is assumed to be a piecewise discontinuous polynomial of degree ‘degree’ throughout the domain defined by ‘grid’. A target value is given and we need to determine the corresponding grid coordinate to define the new grid.

If the target value is out of range, the grid coordinate is simply set to be equal to one of the boundary coordinates, which results in vanished layers near the boundaries.

IT IS ASSUMED THAT THE PIECEWISE POLYNOMIAL IS MONOTONICALLY INCREASING. IF THIS IS NOT THE CASE, THE NEW GRID MAY BE ILL-DEFINED.

It is assumed that the number of cells defining ‘grid’ and ‘ppoly’ are the same.

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

  • h :: [in] Grid cell thicknesses [H]

  • x_g :: [in] Grid interface locations [H]

  • edge_values :: [in] Edge values of interpolating polynomials [A]

  • ppoly_coefs :: [in] Coefficients of interpolating polynomials [A]

  • target_value :: [in] Target value to find position for [A]

  • degree :: [in] Degree of the interpolating polynomials

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

Return:

undefined :: The position of x_g at which target_value is found [H]

Call to:

mom_error_handler::mom_error nr_iterations nr_offset nr_tolerance

Called from:

interpolate_grid

function regrid_interp/interpolation_scheme(interp_scheme) [integer]

Numeric value of interpolation_scheme corresponding to scheme name.

Parameters:

interp_scheme :: [in] Name of the interpolation scheme Valid values include “P1M_H2”, “P1M_H4”, “P1M_IH2”, “PLM”, “PPM_CW”, “PPM_H4”, “PPM_IH4”, “P3M_IH4IH3”, “P3M_IH6IH5”, “PQM_IH4IH3”, and “PQM_IH6IH5”

Call to:

interpolation_p1m_h2 interpolation_p1m_h4 interpolation_p1m_ih4 interpolation_p3m_ih4ih3 interpolation_p3m_ih6ih5 interpolation_plm interpolation_ppm_cw interpolation_ppm_h4 interpolation_ppm_ih4 interpolation_pqm_ih4ih3 interpolation_pqm_ih6ih5 mom_error_handler::mom_error mom_string_functions::uppercase

Called from:

set_interp_scheme

subroutine regrid_interp/set_interp_scheme(CS, interp_scheme)

Store the interpolation_scheme value in the interp_CS based on the input string.

Parameters:
  • cs :: [inout] A control structure for regrid_interp()

  • interp_scheme :: [in] Name of the interpolation scheme Valid values include “P1M_H2”, “P1M_H4”, “P1M_IH2”, “PLM”, “PPM_CW”, “PPM_H4”, “PPM_IH4”, “P3M_IH4IH3”, “P3M_IH6IH5”, “PQM_IH4IH3”, and “PQM_IH6IH5”

Call to:

interpolation_scheme

Called from:

mom_regridding::set_regrid_params

subroutine regrid_interp/set_interp_extrap(CS, extrap)

Store the boundary_extrapolation value in the interp_CS.

Parameters:
  • cs :: [inout] A control structure for regrid_interp()

  • extrap :: [in] Indicate whether high-order boundary extrapolation should be used in boundary cells

Called from:

mom_regridding::set_regrid_params

subroutine regrid_interp/set_interp_answer_date(CS, answer_date)

Store the value of the answer_date in the interp_CS.

Parameters:
  • cs :: [inout] A control structure for regrid_interp()

  • answer_date :: [in] An integer encoding the vintage of the expressions to use for regridding

Called from:

mom_regridding::set_regrid_params