MOM6
|
Data Types | |
type | interp_cs_type |
Functions/Subroutines | |
subroutine, public | regridding_set_ppolys (CS, densities, n0, h0, ppoly0_E, ppoly0_S, ppoly0_coefficients, degree) |
Given the set of target values and cell densities, this routine 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). More... | |
subroutine, public | interpolate_grid (n0, h0, x0, ppoly0_E, ppoly0_coefficients, target_values, degree, n1, h1, x1) |
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) |
real function | get_polynomial_coordinate (N, h, x_g, ppoly_E, ppoly_coefficients, target_value, degree) |
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) |
subroutine, public | set_interp_extrap (CS, extrapolation) |
Variables | |
integer, parameter | interpolation_p1m_h2 = 0 |
O(h^2) More... | |
integer, parameter | interpolation_p1m_h4 = 1 |
O(h^2) More... | |
integer, parameter | interpolation_p1m_ih4 = 2 |
O(h^2) More... | |
integer, parameter | interpolation_plm = 3 |
O(h^2) More... | |
integer, parameter | interpolation_ppm_h4 = 4 |
O(h^3) More... | |
integer, parameter | interpolation_ppm_ih4 = 5 |
O(h^3) More... | |
integer, parameter | interpolation_p3m_ih4ih3 = 6 |
O(h^4) More... | |
integer, parameter | interpolation_p3m_ih6ih5 = 7 |
O(h^4) More... | |
integer, parameter | interpolation_pqm_ih4ih3 = 8 |
O(h^4) More... | |
integer, parameter | interpolation_pqm_ih6ih5 = 9 |
O(h^5) More... | |
integer, parameter | degree_1 = 1 |
List of interpolant degrees. More... | |
integer, parameter | degree_2 = 2 |
integer, parameter | degree_3 = 3 |
integer, parameter | degree_4 = 4 |
integer, parameter, public | degree_max = 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. More... | |
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. More... | |
real, parameter, public | nr_tolerance = 1e-12 |
Tolerance for Newton-Raphson iterations (stop when increment falls below this) More... | |
subroutine, public regrid_interp::build_and_interpolate_grid | ( | type(interp_cs_type), intent(in) | CS, |
real, dimension(:), intent(in) | densities, | ||
integer, intent(in) | n0, | ||
real, dimension(:), intent(in) | h0, | ||
real, dimension(:), intent(in) | x0, | ||
real, dimension(:), intent(in) | target_values, | ||
integer, intent(in) | n1, | ||
real, dimension(:), intent(inout) | h1, | ||
real, dimension(:), intent(inout) | x1 | ||
) |
Definition at line 279 of file regrid_interp.F90.
References interpolate_grid(), and regridding_set_ppolys().
Referenced by coord_hycom::build_hycom1_column(), and coord_rho::build_rho_column().
|
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 |
[in] | x_g | Grid interface locations |
[in] | ppoly_e | Edge values of interpolating polynomials |
[in] | ppoly_coefficients | Coefficients of interpolating polynomials |
[in] | target_value | Target value to find position for |
[in] | degree | Degree of the interpolating polynomials |
Definition at line 313 of file regrid_interp.F90.
References mom_error_handler::mom_error(), nr_iterations, nr_offset, and nr_tolerance.
Referenced by interpolate_grid().
subroutine, public regrid_interp::interpolate_grid | ( | integer, intent(in) | n0, |
real, dimension(:), intent(in) | h0, | ||
real, dimension(:), intent(in) | x0, | ||
real, dimension(:,:), intent(in) | ppoly0_E, | ||
real, dimension(:,:), intent(in) | ppoly0_coefficients, | ||
real, dimension(:), intent(in) | target_values, | ||
integer, intent(in) | degree, | ||
integer, intent(in) | n1, | ||
real, dimension(:), intent(inout) | h1, | ||
real, dimension(:), intent(inout) | 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.
[in] | n0 | Number of points on source grid |
[in] | h0 | Thicknesses of source grid cells |
[in] | x0 | Source interface positions |
[in] | ppoly0_e | Edge values of interpolating polynomials |
[in] | ppoly0_coefficients | Coefficients of interpolating polynomials |
[in] | target_values | Target values of interfaces |
[in] | degree | Degree of interpolating polynomials |
[in] | n1 | Number of points on target grid |
[in,out] | h1 | Thicknesses of target grid cells |
[in,out] | x1 | Target interface positions |
Definition at line 247 of file regrid_interp.F90.
References get_polynomial_coordinate().
Referenced by build_and_interpolate_grid().
|
private |
Numeric value of interpolation_scheme corresponding to scheme name.
[in] | interp_scheme | Name of interpolation scheme |
Definition at line 449 of file regrid_interp.F90.
References interpolation_p1m_h2, interpolation_p1m_h4, interpolation_p1m_ih4, interpolation_p3m_ih4ih3, interpolation_p3m_ih6ih5, interpolation_plm, interpolation_ppm_h4, interpolation_ppm_ih4, interpolation_pqm_ih4ih3, interpolation_pqm_ih6ih5, mom_error_handler::mom_error(), and mom_string_functions::uppercase().
Referenced by set_interp_scheme().
subroutine, public regrid_interp::regridding_set_ppolys | ( | type(interp_cs_type), intent(in) | CS, |
real, dimension(:), intent(in) | densities, | ||
integer, intent(in) | n0, | ||
real, dimension(:), intent(in) | h0, | ||
real, dimension(:,:), intent(inout) | ppoly0_E, | ||
real, dimension(:,:), intent(inout) | ppoly0_S, | ||
real, dimension(:,:), intent(inout) | ppoly0_coefficients, | ||
integer, intent(inout) | degree | ||
) |
Given the set of target values and cell densities, this routine 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] | densities | Actual cell densities |
[in] | n0 | Number of cells on source grid |
[in] | h0 | cell widths on source grid |
[in,out] | ppoly0_e | Edge value of polynomial |
[in,out] | ppoly0_s | Edge slope of polynomial |
[in,out] | ppoly0_coefficients | Coefficients of polynomial |
[in,out] | degree | The degree of the polynomials |
Definition at line 73 of file regrid_interp.F90.
References degree_1, degree_2, degree_3, degree_4, regrid_edge_slopes::edge_slopes_implicit_h3(), regrid_edge_slopes::edge_slopes_implicit_h5(), regrid_edge_values::edge_values_implicit_h4(), regrid_edge_values::edge_values_implicit_h6(), interpolation_p1m_h2, interpolation_p1m_h4, interpolation_p1m_ih4, interpolation_p3m_ih4ih3, interpolation_p3m_ih6ih5, interpolation_plm, 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_boundary_extrapolation(), ppm_functions::ppm_reconstruction(), pqm_functions::pqm_boundary_extrapolation_v1(), and pqm_functions::pqm_reconstruction().
Referenced by build_and_interpolate_grid().
subroutine, public regrid_interp::set_interp_extrap | ( | type(interp_cs_type), intent(inout) | CS, |
logical, intent(in) | extrapolation | ||
) |
Definition at line 475 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 | ||
) |
Definition at line 468 of file regrid_interp.F90.
References interpolation_scheme().
|
private |
List of interpolant degrees.
Definition at line 49 of file regrid_interp.F90.
Referenced by regridding_set_ppolys().
|
private |
Definition at line 49 of file regrid_interp.F90.
Referenced by regridding_set_ppolys().
|
private |
Definition at line 49 of file regrid_interp.F90.
Referenced by regridding_set_ppolys().
|
private |
Definition at line 49 of file regrid_interp.F90.
Referenced by regridding_set_ppolys().
integer, parameter, public regrid_interp::degree_max = 5 |
Definition at line 50 of file regrid_interp.F90.
integer, parameter regrid_interp::interpolation_p1m_h2 = 0 |
O(h^2)
Definition at line 37 of file regrid_interp.F90.
Referenced by interpolation_scheme(), and regridding_set_ppolys().
|
private |
O(h^2)
Definition at line 38 of file regrid_interp.F90.
Referenced by interpolation_scheme(), and regridding_set_ppolys().
|
private |
O(h^2)
Definition at line 39 of file regrid_interp.F90.
Referenced by interpolation_scheme(), and regridding_set_ppolys().
|
private |
O(h^4)
Definition at line 43 of file regrid_interp.F90.
Referenced by interpolation_scheme(), and regridding_set_ppolys().
|
private |
O(h^4)
Definition at line 44 of file regrid_interp.F90.
Referenced by interpolation_scheme(), and regridding_set_ppolys().
|
private |
O(h^2)
Definition at line 40 of file regrid_interp.F90.
Referenced by interpolation_scheme(), and regridding_set_ppolys().
|
private |
O(h^3)
Definition at line 41 of file regrid_interp.F90.
Referenced by interpolation_scheme(), and regridding_set_ppolys().
|
private |
O(h^3)
Definition at line 42 of file regrid_interp.F90.
Referenced by interpolation_scheme(), and regridding_set_ppolys().
|
private |
O(h^4)
Definition at line 45 of file regrid_interp.F90.
Referenced by interpolation_scheme(), and regridding_set_ppolys().
|
private |
O(h^5)
Definition at line 46 of file regrid_interp.F90.
Referenced by interpolation_scheme(), and regridding_set_ppolys().
integer, parameter, public regrid_interp::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.
Definition at line 59 of file regrid_interp.F90.
Referenced by get_polynomial_coordinate(), and coord_slight::rho_interfaces_col().
real, parameter, public regrid_interp::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.
Definition at line 55 of file regrid_interp.F90.
Referenced by get_polynomial_coordinate().
real, parameter, public regrid_interp::nr_tolerance = 1e-12 |
Tolerance for Newton-Raphson iterations (stop when increment falls below this)
Definition at line 61 of file regrid_interp.F90.
Referenced by get_polynomial_coordinate(), and coord_slight::rho_interfaces_col().