inter_module Module Reference

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, reverse_type)
 The inverse of the 1D interpolation function. More...
 
subroutine, private inverse_nr (n, xp, xb, yp, res, flin, itype)
 The inverse of the 1D interpolation function. More...
 
subroutine, private inverse_brent (n, xp, xb, yp, res, flin, itype)
 The inverse of the 1D interpolation function. More...
 

Detailed Description

Module inter_module with interpolation routines.

This routines are used to interpolate the hydrodynamic quantaties (e.g., temperature, density,..)

Function/Subroutine Documentation

◆ akimainter()

subroutine, private inter_module::akimainter ( 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 
)
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:

  • if flin is absent: return exp(log(y)) = y;
  • if flin is present: return log(y).
See also
Akima, Journal of the ACM, 17: 589–602, 1970, Wikipedia: Akima Spline
Author
M. Reichert
Date
30.04.2023
Parameters
[in]narray sizes
[in]xparray of x-values
[in]xbvalue of x
[in]yparray of log(y)-values
[out]resinterpolated value of y
[in]flinif present, return log(y) instead of y

Definition at line 364 of file inter_module.f90.

◆ calc_derivative_2d()

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)} \]

Author
: M. Reichert
Date
06.07.22
Parameters
[in]y_dimDimension of x- and y-
[in]fFunction values
[in]xx values
[in]yy values
[in,out]dfxDerivative with respect to x
[in,out]dfyDerivative with respect to y
[in,out]dfxyDerivative with respect to x and y

Definition at line 623 of file inter_module.f90.

◆ cubinter()

subroutine, private inter_module::cubinter ( 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 
)
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:

  • if flin is absent: return exp(log(y)) = y;
  • if flin is present: return log(y).

Edited:

  • 14.01.14
Parameters
[in]narray sizes
[in]xparray of x-values
[in]xbvalue of x
[in]yparray of log(y)-values
[out]resinterpolated value of y
[in]flinif present, return log(y) instead of y

Definition at line 161 of file inter_module.f90.

◆ cubinter_2d()

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.

See also
weak_inter_bilinear, weak_inter, readweak_logft
Author
: Moritz Reichert
Date
: 14.03.22
See also
Bicubic interpolation
Parameters
[in]y_dimDimension of x- and y-
[in]yinFunction input values
[in]xX values
[in]yY values
[in]dfxyFunction and derivatives
[in]extrapolation1: None; 2: constant
[in]indiceIndex where to interpolate
Returns
Interpolated value (output)

Definition at line 795 of file inter_module.f90.

Here is the call graph for this function:

◆ get_indice_2d()

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.

Author
: M. Reichert
Date
: 06.07.22
Parameters
[in]y_dimDimension of x- and y-
[in]yvalFunction values
[in]xX-grid values
[in]yY-grid values
[out]indicesIndex in the grid

Definition at line 587 of file inter_module.f90.

◆ interp1d()

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:

  • if flin is absent: return exp(log(y)) = y;
  • if flin is present: return log(y).
Author
M. Reichert
Date
01.05.2023
Parameters
[in]narray sizes
[in]xparray of x-values
[in]xbvalue of x
[in]yparray of log(y)-values
[out]resinterpolated value of y
[in]flinif present, return log(y) instead of y
[in]itypeType of interpolation (1: linear, 2: cubic)

Definition at line 99 of file inter_module.f90.

Here is the call graph for this function:

◆ inverse_brent()

subroutine, private inter_module::inverse_brent ( 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)  flin,
integer, intent(in)  itype 
)
private

The inverse of the 1D interpolation function.

Finds the inverse of interp1d and returns the value of x using Brent's method, a hybrid root-finding algorithm that combines the robustness of the bisection method with the speed of inverse quadratic interpolation and the secant method. Brent's method ensures that the root is bracketed at every iteration, providing a reliable and efficient way to converge on the solution. For more details, see the wikipedia article.

Note
The result depends on xb, because x_i does not necessarily have to be monotonic. The result should be the closest value to xb.
Author
M. Reichert
Date
16.10.2024
Parameters
[in]nsize of the array
[in]xparray of x-values
[in]yparray of y-values
[in]xbtarget x-value for root finding
[in,out]resinterpolated root result
[in]flinflag for linear interpolation
[in]itypeinterpolation type

Definition at line 1085 of file inter_module.f90.

Here is the call graph for this function:

◆ inverse_interp1d()

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,
integer, intent(in), optional  reverse_type 
)

The inverse of the 1D interpolation function.

Interface to the inverse interpolation functions using either Newton-Raphson or Brent's method.

Author
M. Reichert
Date
16.10.2024
See also
timestep_module::restrict_timestep
Parameters
[in]narray sizes
[in]xparray of log(x)-values
[in]xbvalue of x
[in]yparray of y-values
[in,out]resinterpolated value of y
[in]flinflag for linear interpolation
[in]itypeinterpolation type
[in]reverse_typeflag for reverse interpolation (1: NR, 2: Brent)

Definition at line 925 of file inter_module.f90.

Here is the call graph for this function:

◆ inverse_nr()

subroutine, private inter_module::inverse_nr ( 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)  flin,
integer, intent(in)  itype 
)
private

The inverse of the 1D interpolation function.

Finds the inverse of interp1d and returns the value of x. This is done by using a Newton-Raphson method.

Note
The result depends on xb, because x_i does not necessarily have to be monotonic. The result should be the closest value to xb. This is however not always the case in which brents method might be more suitable.
Author
M. Reichert

Edited:

  • 30.01.23
Parameters
[in]narray sizes
[in]xparray of log(x)-values
[in]xbvalue of x
[in]yparray of y-values
[in,out]resinterpolated value of y
[in]flinflag for linear interpolation
[in]itypeinterpolation type

Definition at line 997 of file inter_module.f90.

Here is the call graph for this function:

◆ lininter()

subroutine, private inter_module::lininter ( integer, intent(in)  n,
real(r_kind), dimension(n), intent(in)  ref_array,
real(r_kind), intent(in)  ref_value,
real(r_kind), dimension(n), intent(in)  array,
real(r_kind), intent(out)  res,
logical, intent(in), optional  flin 
)
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:

  • if flin is absent: return exp(log(y)) = y;
  • if flin is present: return log(y).

Edited:

  • 11.01.14
Parameters
[in]narray sizes
[in]ref_arrayarray of x-values {x_i}
[in]arrayarray of log(y)-values {log(y_i)}
[in]ref_valuevalue of log(x)
[out]resinterpolated value of y
[in]flinif present, return log(y) instead of y

Definition at line 40 of file inter_module.f90.

◆ lininter_2d()

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

See also
weak_inter
Author
: M. Reichert
Date
: 06.07.22
Parameters
[in]y_dimDimension of x- and y-
[in]yinFunction input values
[in]xX values
[in]yY values
[in]fFunction and derivatives
[in]extrapolation1: None; 2: constant
[in]indiceIndex where to interpolate
Returns
Interpolated value (output)

Definition at line 699 of file inter_module.f90.

Here is the call graph for this function:

◆ makimainter()

subroutine, private inter_module::makimainter ( 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 
)
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:

  • if flin is absent: return exp(log(y)) = y;
  • if flin is present: return log(y).
See also
https://blogs.mathworks.com/cleve/2019/04/29/makima-piecewise-cubic-interpolation/
Author
M. Reichert
Date
30.04.2023
Parameters
[in]narray sizes
[in]xparray of x-values
[in]xbvalue of x
[in]yparray of log(y)-values
[out]resinterpolated value of y
[in]flinif present, return log(y) instead of y

Definition at line 481 of file inter_module.f90.

◆ pchipinter()

subroutine, private inter_module::pchipinter ( 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 
)
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.

See also
F. N. Fritsch and R. E. Carlson, Monotone Piecewise Cubic Interpolation, SIAM Journal on Numerical Analysis, 17(2), 1980, pp. 238-246. Steffen 1990
Author
: M. Reichert
Date
: 01.05.2023
Parameters
[in]narray sizes
[in]xparray of x-values
[in]xbvalue of x
[in]yparray of log(y)-values
[out]resinterpolated value of y
[in]flinif present, return log(y) instead of y

Definition at line 243 of file inter_module.f90.

Here is the call graph for this function: