|
MOM6
|
Vertical interpolation for regridding.
Data Types | |
| type | interp_cs_type |
| Control structure for regrid_interp module. More... | |
Functions/Subroutines | |
| subroutine, public | 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. More... | |
| subroutine, public | interpolate_grid (n0, h0, x0, ppoly0_E, ppoly0_coefs, target_values, degree, n1, h1, x1, answers_2018) |
| Given target values (e.g., density), build new grid based on polynomial. More... | |
| subroutine, public | 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. More... | |
| real function | get_polynomial_coordinate (N, h, x_g, edge_values, ppoly_coefs, target_value, degree, answers_2018) |
| Given a target value, find corresponding coordinate for given polynomial. More... | |
| integer function | interpolation_scheme (interp_scheme) |
| Numeric value of interpolation_scheme corresponding to scheme name. More... | |
| subroutine, public | set_interp_scheme (CS, interp_scheme) |
| Store the interpolation_scheme value in the interp_CS based on the input string. More... | |
| subroutine, public | set_interp_extrap (CS, extrap) |
| Store the boundary_extrapolation value in the interp_CS. More... | |
Variables | |
| integer, parameter | interpolation_p1m_h2 = 0 |
| O(h^2) | |
| integer, parameter | interpolation_p1m_h4 = 1 |
| O(h^2) | |
| integer, parameter | interpolation_p1m_ih4 = 2 |
| O(h^2) | |
| integer, parameter | interpolation_plm = 3 |
| O(h^2) | |
| integer, parameter | interpolation_ppm_h4 = 4 |
| O(h^3) | |
| integer, parameter | interpolation_ppm_ih4 = 5 |
| O(h^3) | |
| integer, parameter | interpolation_p3m_ih4ih3 = 6 |
| O(h^4) | |
| integer, parameter | interpolation_p3m_ih6ih5 = 7 |
| O(h^4) | |
| integer, parameter | interpolation_pqm_ih4ih3 = 8 |
| O(h^4) | |
| integer, parameter | interpolation_pqm_ih6ih5 = 9 |
| O(h^5) | |
| real, parameter, public | nr_offset = 1e-6 |
| When the N-R algorithm produces an estimate that lies outside [0,1], the estimate is set to be equal to the boundary location, 0 or 1, plus or minus an offset, respectively, when the derivative is zero at the boundary [nondim]. | |
| integer, parameter, public | nr_iterations = 8 |
| Maximum number of Newton-Raphson iterations. Newton-Raphson iterations are used to build the new grid by finding the coordinates associated with target densities and interpolations of degree larger than 1. | |
| real, parameter, public | nr_tolerance = 1e-12 |
| Tolerance for Newton-Raphson iterations (stop when increment falls below this) [nondim]. | |
| integer, parameter | degree_1 = 1 |
| Interpolant degrees. | |
| integer, parameter | degree_2 = 2 |
| Interpolant degrees. | |
| integer, parameter | degree_3 = 3 |
| Interpolant degrees. | |
| integer, parameter | degree_4 = 4 |
| Interpolant degrees. | |
| integer, parameter, public | degree_max = 5 |
| Interpolant degrees. | |
| subroutine, public regrid_interp::build_and_interpolate_grid | ( | type(interp_cs_type), intent(in) | CS, |
| real, dimension(n0), intent(in) | densities, | ||
| integer, intent(in) | n0, | ||
| real, dimension(n0), intent(in) | h0, | ||
| real, dimension(n0+1), intent(in) | x0, | ||
| real, dimension(n1+1), intent(in) | target_values, | ||
| integer, intent(in) | n1, | ||
| real, dimension(n1), intent(inout) | h1, | ||
| real, dimension(n1+1), intent(inout) | x1, | ||
| real, intent(in), optional | h_neglect, | ||
| real, intent(in), optional | h_neglect_edge | ||
| ) |
Build a grid by interpolating for target values.
| [in] | cs | A control structure for regrid_interp |
| [in] | n0 | The number of points on the input grid |
| [in] | n1 | The number of points on the output grid |
| [in] | densities | Input cell densities [R ~> kg m-3] |
| [in] | target_values | Target values of interfaces [R ~> kg m-3] |
| [in] | h0 | Initial cell widths [H] |
| [in] | x0 | Source interface positions [H] |
| [in,out] | h1 | Output cell widths [H] |
| [in,out] | x1 | Target interface positions [H] |
| [in] | h_neglect | A negligibly small width for the purpose of cell reconstructions [H] in the same units as h0. |
| [in] | h_neglect_edge | A negligibly small width for the purpose of edge value calculations [H] in the same units as h0. |
Definition at line 307 of file regrid_interp.F90.
|
private |
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.
| [in] | n | Number of grid cells |
| [in] | h | Grid cell thicknesses [H] |
| [in] | x_g | Grid interface locations [H] |
| [in] | edge_values | Edge values of interpolating polynomials [A] |
| [in] | ppoly_coefs | Coefficients of interpolating polynomials [A] |
| [in] | target_value | Target value to find position for [A] |
| [in] | degree | Degree of the interpolating polynomials |
| [in] | answers_2018 | If true use older, less acccurate expressions. |
Definition at line 352 of file regrid_interp.F90.
| subroutine, public regrid_interp::interpolate_grid | ( | integer, intent(in) | n0, |
| real, dimension(n0), intent(in) | h0, | ||
| real, dimension(n0+1), intent(in) | x0, | ||
| real, dimension(n0,2), intent(in) | ppoly0_E, | ||
| real, dimension(n0,degree_max+1), intent(in) | ppoly0_coefs, | ||
| real, dimension(n1+1), intent(in) | target_values, | ||
| integer, intent(in) | degree, | ||
| integer, intent(in) | n1, | ||
| real, dimension(n1), intent(inout) | h1, | ||
| real, dimension(n1+1), intent(inout) | x1, | ||
| logical, intent(in), optional | answers_2018 | ||
| ) |
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.
| [in] | n0 | Number of points on source grid |
| [in] | n1 | Number of points on target grid |
| [in] | h0 | Thicknesses of source grid cells [H] |
| [in] | x0 | Source interface positions [H] |
| [in] | ppoly0_e | Edge values of interpolating polynomials [A] |
| [in] | ppoly0_coefs | Coefficients of interpolating polynomials [A] |
| [in] | target_values | Target values of interfaces [A] |
| [in] | degree | Degree of interpolating polynomials |
| [in,out] | h1 | Thicknesses of target grid cells [H] |
| [in,out] | x1 | Target interface positions [H] |
| [in] | answers_2018 | If true use older, less acccurate expressions. |
Definition at line 270 of file regrid_interp.F90.
|
private |
Numeric value of interpolation_scheme corresponding to scheme name.
| [in] | interp_scheme | Name of the interpolation scheme Valid values include "P1M_H2", "P1M_H4", "P1M_IH2", "PLM", "PPM_H4", "PPM_IH4", "P3M_IH4IH3", "P3M_IH6IH5", "PQM_IH4IH3", and "PQM_IH6IH5" |
Definition at line 487 of file regrid_interp.F90.
| subroutine, public regrid_interp::regridding_set_ppolys | ( | type(interp_cs_type), intent(in) | CS, |
| real, dimension(n0), intent(in) | densities, | ||
| integer, intent(in) | n0, | ||
| real, dimension(n0), intent(in) | h0, | ||
| real, dimension(n0,2), intent(inout) | ppoly0_E, | ||
| real, dimension(n0,2), intent(inout) | ppoly0_S, | ||
| real, dimension(n0,degree_max+1), intent(inout) | ppoly0_coefs, | ||
| integer, intent(inout) | degree, | ||
| real, intent(in), optional | h_neglect, | ||
| real, intent(in), optional | 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).
| [in] | cs | Interpolation control structure |
| [in] | n0 | Number of cells on source grid |
| [in] | densities | Actual cell densities [A] |
| [in] | h0 | cell widths on source grid [H] |
| [in,out] | ppoly0_e | Edge value of polynomial [A] |
| [in,out] | ppoly0_s | Edge slope of polynomial [A H-1] |
| [in,out] | ppoly0_coefs | Coefficients of polynomial [A] |
| [in,out] | degree | The degree of the polynomials |
| [in] | h_neglect | A negligibly small width for the purpose of cell reconstructions [H] in the same units as h0. |
| [in] | h_neglect_edge | A negligibly small width for the purpose of edge value calculations [H] in the same units as h0. |
Definition at line 78 of file regrid_interp.F90.
| subroutine, public regrid_interp::set_interp_extrap | ( | type(interp_cs_type), intent(inout) | CS, |
| logical, intent(in) | extrap | ||
| ) |
Store the boundary_extrapolation value in the interp_CS.
| [in,out] | cs | A control structure for regrid_interp |
| [in] | extrap | Indicate whether high-order boundary extrapolation should be used in boundary cells |
Definition at line 519 of file regrid_interp.F90.
| subroutine, public regrid_interp::set_interp_scheme | ( | type(interp_cs_type), intent(inout) | CS, |
| character(len=*), intent(in) | interp_scheme | ||
| ) |
Store the interpolation_scheme value in the interp_CS based on the input string.
| [in,out] | cs | A control structure for regrid_interp |
| [in] | interp_scheme | Name of the interpolation scheme Valid values include "P1M_H2", "P1M_H4", "P1M_IH2", "PLM", "PPM_H4", "PPM_IH4", "P3M_IH4IH3", "P3M_IH6IH5", "PQM_IH4IH3", and "PQM_IH6IH5" |
Definition at line 509 of file regrid_interp.F90.