MeteoIODoc 20240426.aefd3c94
libinterpol2D.h
Go to the documentation of this file.
1// SPDX-License-Identifier: LGPL-3.0-or-later
2/***********************************************************************************/
3/* Copyright 2009 WSL Institute for Snow and Avalanche Research SLF-DAVOS */
4/***********************************************************************************/
5/* This file is part of MeteoIO.
6 MeteoIO is free software: you can redistribute it and/or modify
7 it under the terms of the GNU Lesser General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 MeteoIO is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU Lesser General Public License for more details.
15
16 You should have received a copy of the GNU Lesser General Public License
17 along with MeteoIO. If not, see <http://www.gnu.org/licenses/>.
18*/
23#ifndef INTERPOL2D_H
24#define INTERPOL2D_H
25
30#include <vector>
31
32namespace mio {
33
46 public:
48 typedef enum REG_TYPES {
50 R_LIN
52
53 static void stdPressure(const DEMObject& dem, Grid2DObject& grid);
54 static void constant(const double& value, const DEMObject& dem, Grid2DObject& grid);
55 static void IDW(const std::vector<double>& vecData_in, const std::vector<StationData>& vecStations_in,
56 const DEMObject& dem, Grid2DObject& grid, const double& scale, const double& alpha=1.);
57 static void LocalLapseIDW(const std::vector<double>& vecData_in,
58 const std::vector<StationData>& vecStations_in,
59 const DEMObject& dem, const size_t& nrOfNeighbors, const double& MaxDistance,
60 Grid2DObject& grid, const double& scale, const double& alpha=1.);
61 static void ListonWind(const DEMObject& i_dem, Grid2DObject& VW, Grid2DObject& DW, const double& eta);
62 static void CurvatureCorrection(DEMObject& dem, const Grid2DObject& ta, Grid2DObject& grid);
63 static void SteepSlopeRedistribution(const DEMObject& dem, const Grid2DObject& ta, Grid2DObject& grid);
64 static void PrecipSnow(const DEMObject& dem, const Grid2DObject& ta, Grid2DObject& grid);
65 static void ODKriging(const std::vector<double>& vecData,
66 const std::vector<StationData>& vecStations,
67 const DEMObject& dem, const Fit1D& variogram, Grid2DObject& grid);
68
69 static void RyanWind(const DEMObject& dem, Grid2DObject& VW, Grid2DObject& DW);
70 static void Winstral(const DEMObject& dem, const Grid2DObject& TA, const double& dmax, const double& in_bearing, Grid2DObject& grid);
71 static void Winstral(const DEMObject& dem, const Grid2DObject& TA, const Grid2DObject& DW, const Grid2DObject& VW, const double& dmax, Grid2DObject& grid);
72
73 static bool allZeroes(const std::vector<double>& vecData);
74
75 static double getTanMaxSlope(const Grid2DObject& dem, const double& dmin, const double& dmax, const double& bearing, const size_t& ii, const size_t& jj);
76 private:
77 //generic functions
78 static double InvHorizontalDistance(const double& X1, const double& Y1, const double& X2, const double& Y2);
79 static double HorizontalDistance(const double& X1, const double& Y1, const double& X2, const double& Y2);
80 static double HorizontalDistance(const DEMObject& dem, const int& i, const int& j,
81 const double& X2, const double& Y2);
82 static void getNeighbors(const double& x, const double& y,
83 const std::vector<StationData>& vecStations,
84 std::vector< std::pair<double, size_t> >& list);
85 static void buildPositionsVectors(const std::vector<StationData>& vecStations,
86 std::vector<double>& vecEastings, std::vector<double>& vecNorthings);
87
88 //core methods
89 static double IDWCore(const double& x, const double& y,
90 const std::vector<double>& vecData_in,
91 const std::vector<double>& vecEastings, const std::vector<double>& vecNorthings, const double& scale, const double& alpha=1.);
92 static double IDWCore(const std::vector<double>& vecData_in, const std::vector<double>& vecDistance_sq, const double& scale, const double& alpha=1.);
93 static double LLIDW_pixel(const size_t& i, const size_t& j,
94 const std::vector<double>& vecData_in,
95 const std::vector<StationData>& vecStations_in,
96 const DEMObject& dem, const size_t& nrOfNeighbors, const double& MaxDistance, const double& scale, const double& alpha=1.);
97
98 static void steepestDescentDisplacement(const DEMObject& dem, const Grid2DObject& grid, const size_t& ii, const size_t& jj, char &d_i_dest, char &d_j_dest);
99 static double depositAroundCell(const DEMObject& dem, const size_t& ii, const size_t& jj, const double& precip, Grid2DObject &grid);
100
101 static void WinstralSX(const DEMObject& dem, const double& dmax, const double& in_bearing, Grid2DObject& grid);
102 static void WinstralSX(const DEMObject& dem, const double& dmax, const Grid2DObject& DW, Grid2DObject& grid);
103
104 //weighting methods
105 static double weightInvDist(const double& d2);
106 static double weightInvDistSqrt(const double& d2);
107 static double weightInvDist2(const double& d2);
108 double weightInvDistN(const double& d2);
109 double dist_pow; //power for the weighting method weightInvDistN
110};
111
112} //end namespace
113
114#endif
A class to represent DEMs and automatically compute some properties. This class stores elevation grid...
Definition: DEMObject.h:40
A class to perform 1D regressions.
Definition: libfit1D.h:159
A class to represent 2D Grids. Typical application as DEM or Landuse Model.
Definition: Grid2DObject.h:42
A class to perform 2D spatial interpolations. Each parameter to be interpolated declares which interp...
Definition: libinterpol2D.h:45
reg_types
Keywords for selecting the regression algorithm to use.
Definition: libinterpol2D.h:48
@ R_CST
no elevation dependence (ie: constant)
Definition: libinterpol2D.h:49
@ R_LIN
linear elevation dependence
Definition: libinterpol2D.h:50
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,...
Definition: libinterpol2D.cc:617
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 eac...
Definition: libinterpol2D.cc:235
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...
Definition: libinterpol2D.cc:39
static void constant(const double &value, const DEMObject &dem, Grid2DObject &grid)
Grid filling function: This implementation fills the grid with a constant value.
Definition: libinterpol2D.cc:179
static void PrecipSnow(const DEMObject &dem, const Grid2DObject &ta, Grid2DObject &grid)
Distribute precipitation in a way that reflects snow redistribution on the ground,...
Definition: libinterpol2D.cc:576
static void stdPressure(const DEMObject &dem, Grid2DObject &grid)
Grid filling function: This implementation builds a standard air pressure as a function of the elevat...
Definition: libinterpol2D.cc:160
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,...
Definition: libinterpol2D.cc:348
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 l...
Definition: libinterpol2D.cc:668
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 exposu...
Definition: libinterpol2D.cc:799
static void CurvatureCorrection(DEMObject &dem, const Grid2DObject &ta, Grid2DObject &grid)
Distribute precipitation in a way that reflects snow redistribution on the ground,...
Definition: libinterpol2D.cc:426
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,...
Definition: libinterpol2D.cc:953
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 ...
Definition: libinterpol2D.cc:515
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....
Definition: libinterpol2D.cc:309
double bearing(std::string bearing_str)
Converts a string bearing to a compass bearing.
Definition: IOUtils.cc:76
Definition: Config.cc:31