関数interpolate2D_3p マニュアル

(The documentation of function interpolate2D_3p)

Last Update: 2023/4/27


◆機能・用途(Purpose)

最寄りの3地点での2次元データを用いて 指定した任意の地点での値を計算する。
Compute the value of a 2-D data at an arbitrary specified point from the data values at nearest three points.


◆形式(Format)

#include <interpolate.h>
inline double interpolate2D_3p
(const int N,
 const double ∗data_x,const double ∗data_y, const double ∗data_values,
 const double target_x,const double target_y)


◆引数(Arguments)

N 2次元データのサンプル数。
The number of samples of the 2-D data.
data_x 2次元データの\(x\)座標のリスト。 \(N\)個の要素から成る配列として与える。
A list of the \(x\)-coordinates of the 2-D data, given as an array of \(N\) components.
data_y 2次元データの\(y\)座標のリスト。 \(N\)個の要素から成る配列として与える。
A list of the \(y\)-coordinates of the 2-D data, given as an array of \(N\) components.
data_values 引数data_x, data_yで与えた各地点\((x,y)\)における 2次元データの値のリスト。 \(N\)個の要素から成る配列として与える。
A list of the values of the 2-D data at each \((x,y)\) specified by the arguments data_x and data_y; an array of \(N\) components.
target_x 補間値を計算したい地点の\(x\)座標。
The \(x\)-coordinate of the location where the interpolated value is needed.
target_y 補間値を計算したい地点の\(y\)座標。
The \(y\)-coordinate of the location where the interpolated value is needed.


◆戻り値(Return value)

引数target_x, target_yで指定した地点\((x,y)\)における 2次元データの補間値。
The interpolated value of the 2-D data at the location \((x,y)\) specified by arguments target_x and target_y.


◆使用例(Example)

double data_x[]={0.0,2.0,1.0,0.0,2.0};
double data_y[]={0.0,0.0,1.0,3.0,4.0};
double data_values[]={1.2,3.4,5.6,7.8,9.0};
double value=interpolate2D_3p(5,data_x,data_y,data_values,1.5,2.0);

この例では
という5地点で定義された2次元データを元に 位置(1.5,2)での補間値を計算する。
In this example, the interpolated value at a location (1.5, 2) is computed from a 2-D data defined at the following five points:


◆計算式(Formula)

補間値を求めたい地点を\(P_0(x_0,y_0)\)とし、 2次元データが定義されている\(P_0\)の最寄りの3地点を \(P_1(x_1,y_1)\), \(P_2(x_2,y_2)\), \(P_3(x_3,y_3)\)、 それらの地点での2次元データの値を\(h_1\), \(h_2\), \(h_3\)とする。 この関数では2次元データの値が位置\((x,y)\)の線形関数として \[\begin{equation} h(x,y)=h_1+c_x(x-x_1)+c_y(y-y_1) \label{eq.h} \end{equation}\] のように変化すると仮定する。 ここで\(c_x\), \(c_y\)は未知定数である。 点\(P_2\)において\(h(x,y)=h_2\)となる条件から \[\begin{equation} h_2=h_1+c_x(x_2-x_1)+c_y(y_2-y_1) \label{eq.h2} \end{equation}\] 点\(P_3\)において\(h(x,y)=h_3\)となる条件から \[\begin{equation} h_3=h_1+c_x(x_3-x_1)+c_y(y_3-y_1) \label{eq.h3} \end{equation}\] であり、これらは未知数\(c_x\), \(c_y\)に対する連立方程式として行列形式で \[\begin{equation} \begin{pmatrix} x_2-x_1 & y_2-y_1 \\ x_3-x_1 & y_3-y_1 \end{pmatrix} \begin{pmatrix} c_x \\ c_y \end{pmatrix} = \begin{pmatrix} h_2-h_1 \\ h_3-h_1 \end{pmatrix} \label{eq.cx_cy.matrix} \end{equation}\] と書ける。その解はクラメルの公式を用いて \[\begin{eqnarray} c_x &=& \frac{\begin{vmatrix} h_2-h_1 & y_2-y_1 \\ h_3-h_1 & y_3-y_1 \end{vmatrix}} {\begin{vmatrix} x_2-x_1 & y_2-y_1 \\ x_3-x_1 & y_3-y_1 \end{vmatrix}} \nonumber \\ &=& \frac{(h_2-h_1)(y_3-y_1)-(h_3-h_1)(y_2-y_1)} {(x_2-x_1)(y_3-y_1)-(x_3-x_1)(y_2-y_1)} \nonumber \\ &=& \frac{(h_2y_3-h_2y_1-h_1y_3+h_1y_1)-(h_3y_2-h_3y_1-h_1y_2+h_1y_1)} {(x_2y_3-x_2y_1-x_1y_3+x_1y_1)-(x_3y_2-x_3y_1-x_1y_2+x_1y_1)} \nonumber \\ &=& \frac{(h_1y_2-h_2y_1)+(h_2y_3-h_3y_2)+(h_3y_1-h_1y_3)} {(x_1y_2-x_2y_1)+(x_2y_3-x_3y_2)+(x_3y_1-x_1y_3)} \label{eq.cx} \end{eqnarray}\] \[\begin{eqnarray} c_y &=& \frac{\begin{vmatrix} x_2-x_1 & h_2-h_1 \\ x_3-x_1 & h_3-h_1 \end{vmatrix}} {\begin{vmatrix} x_2-x_1 & y_2-y_1 \\ x_3-x_1 & y_3-y_1 \end{vmatrix}} \nonumber \\ &=& \frac{(x_2-x_1)(h_3-h_1)-(x_3-x_1)(h_2-h_1)} {(x_2-x_1)(y_3-y_1)-(x_3-x_1)(y_2-y_1)} \nonumber \\ &=& \frac{(x_2h_3-x_2h_1-x_1h_3+x_1h_1)-(x_3h_2-x_3h_1-x_1h_2+x_1h_1)} {(x_2y_3-x_2y_1-x_1y_3+x_1y_1)-(x_3y_2-x_3y_1-x_1y_2+x_1y_1)} \nonumber \\ &=& \frac{(x_1h_2-x_2h_1)+(x_2h_3-x_3h_2)+(x_3h_1-x_1h_3)} {(x_1y_2-x_2y_1)+(x_2y_3-x_3y_2)+(x_3y_1-x_1y_3)} \label{eq.cy} \end{eqnarray}\] と求められる。 この\(c_x\), \(c_y\)を用いて(\ref{eq.h})式から計算される \[\begin{equation} h_0=h(x_0,y_0)=h_1+c_x(x_0-x_1)+c_y(y_0-y_1) \label{eq.h0} \end{equation}\] が求める\(P_0\)での補間値である。
Let \(P_0(x_0,y_0)\) be the location where the interpolated value is needed, and \(P_1(x_1,y_1)\), \(P_2(x_2,y_2)\), and \(P_3(x_3,y_3)\) be the three closest locations from \(P_0\) where the 2-D data are defined. Let the 2-D data values at the three points be \(h_1\), \(h_2\), and \(h_3\). This function assumes that the 2-D data value varies linearly with the location \((x,y)\) as Eq. (\ref{eq.h}), where \(c_x\) and \(c_y\) are unknown constants. The requirements of \(h(x,y)=h_2\) at \(P_2\) and \(h(x,y)=h_3\) at \(P_3\) can be written as Eqs. (\ref{eq.h2}) and (\ref{eq.h3}), respectively. These equations can be regarded as a simultaneous equation for the unknowns \(c_x\) and \(c_y\), and can be writtein in a matrix form as (\ref{eq.cx_cy.matrix}). The solution can be computed using the Cramer's formula as Eqs. (\ref{eq.cx}) and (\ref{eq.cy}). Using these \(c_x\) and \(c_y\), the interpolated value at \(P_0\) is computed from Eq. (\ref{eq.h}) as (\ref{eq.h0}).