namespace regrid_interp¶
Overview¶
Vertical interpolation for regridding. More…
namespace regrid_interp { // global variables integer, parameter, public degree_max = 5; real, parameter, public nr_offset = 1e-6; integer, parameter, public nr_iterations = 8; real, parameter, public nr_tolerance = 1e-12; integer, parameter interpolation_p1m_h2 = 0; // global functions subroutine, public regridding_set_ppolys( CS CS, densities densities, n0 n0, h0 h0, ppoly0_E ppoly0_E, ppoly0_S ppoly0_S, ppoly0_coefs ppoly0_coefs, degree degree, h_neglect h_neglect, h_neglect_edge h_neglect_edge ); subroutine, public interpolate_grid( n0 n0, h0 h0, x0 x0, ppoly0_E ppoly0_E, ppoly0_coefs ppoly0_coefs, target_values target_values, degree degree, n1 n1, h1 h1, x1 x1 ); subroutine, public build_and_interpolate_grid( CS CS, densities densities, n0 n0, h0 h0, x0 x0, target_values target_values, n1 n1, h1 h1, x1 x1, h_neglect h_neglect, h_neglect_edge h_neglect_edge ); subroutine, public set_interp_scheme(CS CS, interp_scheme interp_scheme); subroutine, public set_interp_extrap(CS CS, extrap extrap); } // namespace regrid_interp
Detailed Documentation¶
Vertical interpolation for regridding.
Global Variables¶
integer, parameter, public degree_max = 5
Interpolant degrees.
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.
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)
integer, parameter interpolation_p1m_h2 = 0
O(h^2)
Global Functions¶
subroutine, public regridding_set_ppolys( CS CS, densities densities, n0 n0, h0 h0, ppoly0_E ppoly0_E, ppoly0_S ppoly0_S, ppoly0_coefs ppoly0_coefs, degree degree, h_neglect h_neglect, h_neglect_edge 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 |
Interpolation control structure |
densities |
Actual cell densities |
n0 |
Number of cells on source grid |
h0 |
cell widths on source grid |
ppoly0_e |
Edge value of polynomial |
ppoly0_s |
Edge slope of polynomial |
ppoly0_coefs |
Coefficients of polynomial |
degree |
The degree of the polynomials |
h_neglect |
A negligibly small width for the purpose of cell reconstructions in the same units as h0. |
h_neglect_edge |
A negligibly small width for the purpose of edge value calculations in the same units as h0. |
subroutine, public interpolate_grid( n0 n0, h0 h0, x0 x0, ppoly0_E ppoly0_E, ppoly0_coefs ppoly0_coefs, target_values target_values, degree degree, n1 n1, h1 h1, x1 x1 )
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 |
Number of points on source grid |
h0 |
Thicknesses of source grid cells |
x0 |
Source interface positions |
ppoly0_e |
Edge values of interpolating polynomials |
ppoly0_coefs |
Coefficients of interpolating polynomials |
target_values |
Target values of interfaces |
degree |
Degree of interpolating polynomials |
n1 |
Number of points on target grid |
h1 |
Thicknesses of target grid cells |
x1 |
Target interface positions |
subroutine, public build_and_interpolate_grid( CS CS, densities densities, n0 n0, h0 h0, x0 x0, target_values target_values, n1 n1, h1 h1, x1 x1, h_neglect h_neglect, h_neglect_edge h_neglect_edge )
Build a grid by interpolating for target values.
Parameters:
cs |
A control structure for regrid_interp |
densities |
Input cell densities [kg m-3] |
target_values |
Target values of interfaces |
n0 |
The number of points on the input grid |
h0 |
Initial cell widths |
x0 |
Source interface positions |
n1 |
The number of points on the output grid |
h1 |
Output cell widths |
x1 |
Target interface positions |
h_neglect |
A negligibly small width for the purpose of cell reconstructions in the same units as h0. |
h_neglect_edge |
A negligibly small width for the purpose of edge value calculations in the same units as h0. |
subroutine, public set_interp_scheme(CS CS, interp_scheme interp_scheme)
Store the interpolation_scheme value in the interp_CS based on the input string.
Parameters:
cs |
A control structure for regrid_interp |
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” |
subroutine, public set_interp_extrap(CS CS, extrap extrap)
Store the boundary_extrapolation value in the interp_CS.
Parameters:
cs |
A control structure for regrid_interp |
extrap |
Indicate whether high-order boundary extrapolation should be used in boundary cells |