A class to perform 2D spatial interpolations. Each parameter to be interpolated declares which interpolation method to use. Then the class computes the interpolation for each 2D grid point, combining the inputs provided by the available data sources.
- Author
- Mathias Bavay
|
static void | stdPressure (const DEMObject &dem, Grid2DObject &grid) |
| Grid filling function: This implementation builds a standard air pressure as a function of the elevation. More...
|
|
static void | constant (const double &value, const DEMObject &dem, Grid2DObject &grid) |
| Grid filling function: This implementation fills the grid with a constant value. More...
|
|
static void | IDW (const std::vector< double > &vecData_in, const std::vector< StationData > &vecStations_in, const DEMObject &dem, Grid2DObject &grid, const double &scale, const double &alpha=1.) |
| Grid filling function: This implementation fills a grid using Inverse Distance Weighting. for example, the air temperatures measured at several stations would be given as values, the stations positions as positions and projected to a grid. No elevation detrending is performed, the DEM is only used for checking if a grid point is "nodata". More...
|
|
static void | LocalLapseIDW (const std::vector< double > &vecData_in, const std::vector< StationData > &vecStations_in, const DEMObject &dem, const size_t &nrOfNeighbors, const double &MaxDistance, Grid2DObject &grid, const double &scale, const double &alpha=1.) |
| Grid filling function: Similar to Interpol2D::LapseIDW but using a limited number of stations for each cell. More...
|
|
static void | ListonWind (const DEMObject &i_dem, Grid2DObject &VW, Grid2DObject &DW, const double &eta) |
| Grid filling function: This implementation fills a grid using a curvature and slope algorithm, as described in G. E. Liston and K. Elder, "A meteorological distribution system for high-resolution terrestrial modeling (MicroMet)", Journal of Hydrometeorology, 7.2, 2006. More...
|
|
static void | CurvatureCorrection (DEMObject &dem, const Grid2DObject &ta, Grid2DObject &grid) |
| Distribute precipitation in a way that reflects snow redistribution on the ground, according to (Huss, 2008) This method modifies the solid precipitation distribution according to the local slope and curvature. See "Quantitative evaluation of different hydrological modelling approaches in a partly glacierized Swiss watershed", Magnusson et All., Hydrological Processes, 2010, under review. and "Modelling runoff from highly glacierized alpine catchments in a changing climate", Huss et All., Hydrological Processes, 22, 3888-3902, 2008. More...
|
|
static void | SteepSlopeRedistribution (const DEMObject &dem, const Grid2DObject &ta, Grid2DObject &grid) |
| redistribute precip from steeper slopes to gentler slopes by following the steepest path from top to bottom and gradually depositing precip during descent More...
|
|
static void | PrecipSnow (const DEMObject &dem, const Grid2DObject &ta, Grid2DObject &grid) |
| Distribute precipitation in a way that reflects snow redistribution on the ground, according to (Huss, 2008) This method modifies the solid precipitation distribution according to the local slope and curvature: all pixels whose slope is greater than 60° will not receive any snow at all. All pixels whose slope is less than 40° will receive full snow and any pixel between 40° and 60° sees a linear correction between 100% and 0% snow. After this step, a curvature correction is applied: pixels having the minimu curvature see 50% snow more, pixels having the maximum curvature see 50% snow less and pixels ate the middle of the curvature range are unaffected. More...
|
|
static void | ODKriging (const std::vector< double > &vecData, const std::vector< StationData > &vecStations, const DEMObject &dem, const Fit1D &variogram, Grid2DObject &grid) |
| Ordinary Kriging matrix formulation This implements the matrix formulation of Ordinary Kriging, as shown (for example) in "Statistics for spatial data", Noel A. C. Cressie, John Wiley & Sons, revised edition, 1993, pp122. More...
|
|
static void | RyanWind (const DEMObject &dem, Grid2DObject &VW, Grid2DObject &DW) |
| compute the change of wind direction by the local terrain This is according to Ryan, "a mathematical model for diagnosis and prediction of surface
winds in mountainous terrain", 1977, journal of applied meteorology, 16, 6. More...
|
|
static void | Winstral (const DEMObject &dem, const Grid2DObject &TA, const double &dmax, const double &in_bearing, Grid2DObject &grid) |
| Alter a precipitation field with the Winstral Sx exposure coefficient This implements the wind exposure coefficient (Sx) for one bearing as in "Simulating wind fields and snow redistribution using terrain‐based parameters to model
snow accumulation and melt over a semi‐arid mountain catchment.", Winstral, Adam, and Danny Marks, Hydrological Processes 16.18 (2002), pp3585-3603. More...
|
|
static void | Winstral (const DEMObject &dem, const Grid2DObject &TA, const Grid2DObject &DW, const Grid2DObject &VW, const double &dmax, Grid2DObject &grid) |
|
static bool | allZeroes (const std::vector< double > &vecData) |
| check if the points measurements are all at zero This check can be performed to trigger optimizations: it is quicker to fill the grid directly with zeroes instead of running a complicated algorithm. More...
|
|
static double | getTanMaxSlope (const Grid2DObject &dem, const double &dmin, const double &dmax, const double &bearing, const size_t &ii, const size_t &jj) |
| compute the max slope angle looking toward the horizon in a given direction. The search distance is limited between dmin and dmax from the starting point (i,j) and the elvation difference must be at least 2 meters (otherwise it is considered flat). If a slope start to be positive/negative before turning negative/positive, the first one would be returned (so a pixel just behind a ridge would still see an uphill slope to the ridge even if the slope behind the ridge would be greater). More...
|
|
void mio::Interpol2D::ODKriging |
( |
const std::vector< double > & |
vecData, |
|
|
const std::vector< StationData > & |
vecStations, |
|
|
const DEMObject & |
dem, |
|
|
const Fit1D & |
variogram, |
|
|
Grid2DObject & |
grid |
|
) |
| |
|
static |
Ordinary Kriging matrix formulation This implements the matrix formulation of Ordinary Kriging, as shown (for example) in "Statistics for spatial data", Noel A. C. Cressie, John Wiley & Sons, revised edition, 1993, pp122.
First, Ordinary kriging assumes stationarity of the mean of all random variables. We start by solving the following system:
\begin{eqnarray*}
\mathbf{\lambda} & = & \mathbf{\Gamma_0^{-1}} \cdot \mathbf{\gamma^*} \\
\left[
\begin{array}{c}
\lambda_1 \\
\vdots \\
\lambda_i \\
\mu
\end{array}
\right]
&
=
&
{
\left[
\begin{array}{cccc}
\Gamma_{1,1} & \cdots & \Gamma_{1,i} & 1 \\
\vdots & \ddots & \vdots & \vdots \\
\Gamma_{i,1} & \cdots & \Gamma_{i,i} & 1 \\
1 & \cdots & 1 & 0
\end{array}
\right]
}^{-1}
\cdot
\left[
\begin{array}{c}
\gamma^*_1 \\
\vdots \\
\gamma^*_i \\
1
\end{array}
\right]
\end{eqnarray*}
where the \(\lambda_i\) are the interpolation weights (at each station i), \(\mu\) is the Lagrange multiplier (used to minimize the error), \(\Gamma_{i,j}\) is the covariance between the stations i and j and \(\gamma^*_i\) the covariances between the station i and the local position where the interpolation has to be computed. This covariance is computed based on distance, using the variogram that gives covariance = f(distance). The variogram is established by fitting a statistical model to all the (distance, covariance) points originating from the station measurements. The statistical model of the variogram enables computing the covariance for any distance, therefore it is possible to compute the \(\gamma^*_i\).
Once the \(\lambda_i\) have been computed, the locally interpolated value is computed as
\[
\mathbf{X^*} = \sum \mathbf{\lambda_i} * \mathbf{X_i}
\]
where \(X_i\) is the value measured at station i.
- Parameters
-
vecData | vector containing the values as measured at the stations |
vecStations | vector of stations |
dem | digital elevation model |
variogram | variogram regression model |
grid | 2D array of precipitation to fill |
- Author
- Mathias Bavay
Distribute precipitation in a way that reflects snow redistribution on the ground, according to (Huss, 2008) This method modifies the solid precipitation distribution according to the local slope and curvature: all pixels whose slope is greater than 60° will not receive any snow at all. All pixels whose slope is less than 40° will receive full snow and any pixel between 40° and 60° sees a linear correction between 100% and 0% snow. After this step, a curvature correction is applied: pixels having the minimu curvature see 50% snow more, pixels having the maximum curvature see 50% snow less and pixels ate the middle of the curvature range are unaffected.
For more, see "Quantitative evaluation of different hydrological modelling approaches in a partly glacierized Swiss watershed", Magnusson et All., Hydrological Processes, 2010, under review. and "Modelling runoff from highly glacierized alpine catchments in a changing climate", Huss et All., Hydrological Processes, 22, 3888-3902, 2008.
- Parameters
-
dem | array of elevations (dem). The slope must have been updated as it is required for the DEM analysis. |
ta | array of air temperatures used to determine if precipitation is rain or snow |
grid | 2D array of precipitation to fill |
- Author
- Florian Kobierska, Jan Magnusson and Mathias Bavay