Module inter_module with interpolation routines. More...
Functions/Subroutines | |
subroutine, private | lininter (n, ref_array, ref_value, array, res, flin) |
Linear interpolation in lin-log space. More... | |
subroutine, public | interp1d (n, xp, xb, yp, res, flin, itype) |
Interface for 1D interpolation routines. More... | |
subroutine, private | cubinter (n, xp, xb, yp, res, flin) |
Cubic interpolation in log-log space. Here, it is assumed that. More... | |
subroutine, private | pchipinter (n, xp, xb, yp, res, flin) |
PCHIP interpolation in lin-log space. More... | |
subroutine, private | akimainter (n, xp, xb, yp, res, flin) |
Interpolation using prescription of Akima. More... | |
subroutine, private | makimainter (n, xp, xb, yp, res, flin) |
Interpolation using modified prescription of Akima. More... | |
subroutine, public | get_indice_2d (xval, yval, x, y, x_dim, y_dim, indices) |
Get indices of the grid for given x- and y- values. More... | |
subroutine, public | calc_derivative_2d (f, x, y, dfx, dfy, dfxy, x_dim, y_dim) |
Calculate derivatives at all points in a given grid. More... | |
real(r_kind) function, public | lininter_2d (xin, yin, x, y, f, x_dim, y_dim, indice, extrapolation) |
Interpolate linearly on 2 Dimensional grid. More... | |
real(r_kind) function, public | cubinter_2d (xin, yin, x, y, f, dfx, dfy, dfxy, x_dim, y_dim, indice, extrapolation) |
Cubic interpolate on 2-Dimensional grid. More... | |
subroutine, public | inverse_interp1d (n, xp, xb, yp, res, flin, itype) |
The inverse of the 1D interpolation function. More... | |
Module inter_module with interpolation routines.
This routines are used to interpolate the hydrodynamic quantaties (e.g., temperature, density,..)
|
private |
Interpolation using prescription of Akima.
The interpolation is a spline interpolation that tries to avoid overshoots.
Here, it is assumed that the x and y arrays are given in the form of {x_i, log(y_i)}, but before interpolating the x coordinate, it is transformed to a log space: {log(x_i), log(y_i)}. An optional argument flin allows to control the return value:
[in] | n | array sizes |
[in] | xp | array of x-values |
[in] | xb | value of x |
[in] | yp | array of log(y)-values |
[out] | res | interpolated value of y |
[in] | flin | if present, return log(y) instead of y |
Definition at line 363 of file inter_module.f90.
subroutine, public inter_module::calc_derivative_2d | ( | real(r_kind), dimension(x_dim, y_dim), intent(in) | f, |
real(r_kind), dimension(x_dim), intent(in) | x, | ||
real(r_kind), dimension(y_dim), intent(in) | y, | ||
real(r_kind), dimension(x_dim,y_dim), intent(inout) | dfx, | ||
real(r_kind), dimension(x_dim,y_dim), intent(inout) | dfy, | ||
real(r_kind), dimension(x_dim,y_dim), intent(inout) | dfxy, | ||
integer, intent(in) | x_dim, | ||
integer, intent(in) | y_dim | ||
) |
Calculate derivatives at all points in a given grid.
These derivatives are necessary for the cubic interpolation. The subroutine returns the derivative df/dx, df/dy, and ddf/dxdy. They are second order accurate, i.e.,
\[ df/dx = \frac{h_1^2 f(x+h_2) + (h_2^2-h_1^2)f(x)-h_2^2 f(x-h_1) }{h_1 h_2 (h_1+h_2)} \]
[in] | y_dim | Dimension of x- and y- |
[in] | f | Function values |
[in] | x | x values |
[in] | y | y values |
[in,out] | dfx | Derivative with respect to x |
[in,out] | dfy | Derivative with respect to y |
[in,out] | dfxy | Derivative with respect to x and y |
Definition at line 622 of file inter_module.f90.
|
private |
Cubic interpolation in log-log space. Here, it is assumed that.
the x and y arrays are given in the form of {x_i, log(y_i)}, but before interpolating the x coordinate, it is transformed to a log space: {log(x_i), log(y_i)}. An optional argument flin allows to control the return value:
Edited:
[in] | n | array sizes |
[in] | xp | array of x-values |
[in] | xb | value of x |
[in] | yp | array of log(y)-values |
[out] | res | interpolated value of y |
[in] | flin | if present, return log(y) instead of y |
Definition at line 160 of file inter_module.f90.
real(r_kind) function, public inter_module::cubinter_2d | ( | real(r_kind), intent(in) | xin, |
real(r_kind), intent(in) | yin, | ||
real(r_kind), dimension(x_dim), intent(in) | x, | ||
real(r_kind), dimension(y_dim), intent(in) | y, | ||
real(r_kind), dimension(x_dim,y_dim), intent(in) | f, | ||
real(r_kind), dimension(x_dim,y_dim), intent(in) | dfx, | ||
real(r_kind), dimension(x_dim,y_dim), intent(in) | dfy, | ||
real(r_kind), dimension(x_dim,y_dim), intent(in) | dfxy, | ||
integer, intent(in) | x_dim, | ||
integer, intent(in) | y_dim, | ||
integer, dimension(2,2), intent(in), optional | indice, | ||
integer, intent(in), optional | extrapolation | ||
) |
Cubic interpolate on 2-Dimensional grid.
This interpolation is a bicubic interpolation. For values outside the grid, either a constant extrapolation is performed or an error is raised.
[in] | y_dim | Dimension of x- and y- |
[in] | yin | Function input values |
[in] | x | X values |
[in] | y | Y values |
[in] | dfxy | Function and derivatives |
[in] | extrapolation | 1: None; 2: constant |
[in] | indice | Index where to interpolate |
Definition at line 794 of file inter_module.f90.
subroutine, public inter_module::get_indice_2d | ( | real(r_kind), intent(in) | xval, |
real(r_kind), intent(in) | yval, | ||
real(r_kind), dimension(x_dim), intent(in) | x, | ||
real(r_kind), dimension(y_dim), intent(in) | y, | ||
integer, intent(in) | x_dim, | ||
integer, intent(in) | y_dim, | ||
integer, dimension(2,2), intent(out) | indices | ||
) |
Get indices of the grid for given x- and y- values.
[in] | y_dim | Dimension of x- and y- |
[in] | yval | Function values |
[in] | x | X-grid values |
[in] | y | Y-grid values |
[out] | indices | Index in the grid |
Definition at line 586 of file inter_module.f90.
subroutine, public inter_module::interp1d | ( | integer, intent(in) | n, |
real(r_kind), dimension(n), intent(in) | xp, | ||
real(r_kind), intent(in) | xb, | ||
real(r_kind), dimension(n), intent(in) | yp, | ||
real(r_kind), intent(out) | res, | ||
logical, intent(in), optional | flin, | ||
integer, intent(in), optional | itype | ||
) |
Interface for 1D interpolation routines.
This routine is a wrapper for the 1D interpolation routines lininter, cubinter, akimainter, makimainter, and pchipinter. The interpolation type is controlled by the optional argument itype. If itype is absent, interp_mode is used. The optional argument flin allows to control the return value:
[in] | n | array sizes |
[in] | xp | array of x-values |
[in] | xb | value of x |
[in] | yp | array of log(y)-values |
[out] | res | interpolated value of y |
[in] | flin | if present, return log(y) instead of y |
[in] | itype | Type of interpolation (1: linear, 2: cubic) |
Definition at line 98 of file inter_module.f90.
subroutine, public inter_module::inverse_interp1d | ( | integer, intent(in) | n, |
real(r_kind), dimension(n), intent(in) | xp, | ||
real(r_kind), intent(in) | xb, | ||
real(r_kind), dimension(n), intent(in) | yp, | ||
real(r_kind), intent(inout) | res, | ||
logical, intent(in), optional | flin, | ||
integer, intent(in), optional | itype | ||
) |
The inverse of the 1D interpolation function.
Finds the inverse of cubinter and returns the value of x. This is done by using a Newton-Raphson method.
Edited:
[in] | n | array sizes |
[in] | xp | array of log(x)-values |
[in] | xb | value of x |
[in] | yp | array of y-values |
[in,out] | res | interpolated value of y |
[in] | flin | flag for linear interpolation |
[in] | itype | interpolation type |
Definition at line 930 of file inter_module.f90.
|
private |
Linear interpolation in lin-log space.
This routine assumes that the x and y arrays are given in the form of {x_i, log(y_i)}. An optional argument flin allows to control the return value:
Edited:
[in] | n | array sizes |
[in] | ref_array | array of x-values {x_i} |
[in] | array | array of log(y)-values {log(y_i)} |
[in] | ref_value | value of log(x) |
[out] | res | interpolated value of y |
[in] | flin | if present, return log(y) instead of y |
Definition at line 39 of file inter_module.f90.
real(r_kind) function, public inter_module::lininter_2d | ( | real(r_kind), intent(in) | xin, |
real(r_kind), intent(in) | yin, | ||
real(r_kind), dimension(x_dim), intent(in) | x, | ||
real(r_kind), dimension(y_dim), intent(in) | y, | ||
real(r_kind), dimension(x_dim,y_dim), intent(in) | f, | ||
integer, intent(in) | x_dim, | ||
integer, intent(in) | y_dim, | ||
integer, dimension(2,2), intent(in), optional | indice, | ||
integer, intent(in), optional | extrapolation | ||
) |
Interpolate linearly on 2 Dimensional grid.
Bilinear interpolation
[in] | y_dim | Dimension of x- and y- |
[in] | yin | Function input values |
[in] | x | X values |
[in] | y | Y values |
[in] | f | Function and derivatives |
[in] | extrapolation | 1: None; 2: constant |
[in] | indice | Index where to interpolate |
Definition at line 698 of file inter_module.f90.
|
private |
Interpolation using modified prescription of Akima.
The interpolation is a spline interpolation that tries to avoid overshoots. It is doing this slightly more radically than the Akima interpolation, but not as strong as the pchip interpolation.
Here, it is assumed that the x and y arrays are given in the form of {x_i, log(y_i)}, but before interpolating the x coordinate, it is transformed to a log space: {log(x_i), log(y_i)}. An optional argument flin allows to control the return value:
[in] | n | array sizes |
[in] | xp | array of x-values |
[in] | xb | value of x |
[in] | yp | array of log(y)-values |
[out] | res | interpolated value of y |
[in] | flin | if present, return log(y) instead of y |
Definition at line 480 of file inter_module.f90.
|
private |
PCHIP interpolation in lin-log space.
The Piecewise Cubic Hermite Interpolating Polynomial (PCHIP) interpolation is a method for interpolating data points where every section between points is monoton. The interpolation avoids overshoots and undershoots and is close to the linear interpolation, but having the first derivative continuous.
[in] | n | array sizes |
[in] | xp | array of x-values |
[in] | xb | value of x |
[in] | yp | array of log(y)-values |
[out] | res | interpolated value of y |
[in] | flin | if present, return log(y) instead of y |
Definition at line 242 of file inter_module.f90.