Gridpp 0.6.0
A post-processing library for gridded weather forecasts
gridpp.h
Go to the documentation of this file.
1 #ifndef GRIDPP_API_H
2 #define GRIDPP_API_H
3 #include <vector>
4 #include <string>
5 #include <boost/geometry.hpp>
6 #include <boost/geometry/geometries/point.hpp>
7 #include <boost/geometry/geometries/box.hpp>
8 #include <boost/geometry/index/rtree.hpp>
9 #include <boost/math/distributions/gamma.hpp>
10 #include <boost/math/distributions/normal.hpp>
11 #ifdef _OPENMP
12  #include <omp.h>
13 #endif
14 #include <exception>
15 
16 #define GRIDPP_VERSION "0.6.0"
17 #define __version__ GRIDPP_VERSION
18 
19 namespace gridpp {
23  // Preferred vector types
24  typedef std::vector<float> vec;
25  typedef std::vector<vec> vec2;
26  typedef std::vector<vec2> vec3;
27  typedef std::vector<int> ivec;
28  typedef std::vector<ivec> ivec2;
29  typedef std::vector<ivec2> ivec3;
30 
31  // Only use when double is required
32  typedef std::vector<double> dvec;
33  typedef std::vector<dvec> dvec2;
41  static const float MV = NAN;
43  static const float MV_CML = -999;
45  static const float pi = 3.14159265;
47  static const double radius_earth = 6.378137e6;
49  static const float lapse_rate=0.0065;
51  static const float standard_surface_temperature = 288.15;
53  static const float gravit = 9.80665;
55  static const float molar_mass = 0.0289644;
57  static const float gas_constant_mol = 8.31447;
59  static const float gas_constant_si = 287.05;
62  class KDTree;
63  class Points;
64  class Grid;
65  class Point;
66  class Nearest;
67  class StructureFunction;
68  class Transform;
69 
72  OneToOne = 0,
73  MeanSlope = 10,
74  NearestSlope = 20,
75  Zero = 30,
76  Unchanged = 40,
77  };
78 
80  enum Statistic {
81  Mean = 0,
82  Min = 10,
83  Median = 20,
84  Max = 30,
85  Quantile = 40,
86  Std = 50,
87  Variance = 60,
88  Sum = 70,
89  Count = 80,
90  Unknown = -1
91  };
92 
94  enum Metric {
95  Ets = 0,
96  Ts = 1,
97  Kss = 20,
98  Pc = 30,
99  Bias = 40,
100  Hss = 50,
101  };
102 
105  Qq = 0,
107  Additive = 20,
108  };
109 
112  Geodetic = 0,
113  Cartesian = 1,
114  };
115 
118  MinMax = 0,
120  };
121 
123  enum Downscaler {
124  Nearest = 0,
125  Bilinear = 1,
126  };
127 
130  Lt = 0,
131  Leq = 10,
132  Gt = 20,
133  Geq = 30,
134  };
135 
152  vec2 optimal_interpolation(const Grid& bgrid,
153  const vec2& background,
154  const Points& points,
155  const vec& pobs,
156  const vec& pratios,
157  const vec& pbackground,
158  const StructureFunction& structure,
159  int max_points,
160  bool allow_extrapolation=true);
161 
173  vec optimal_interpolation(const Points& bpoints,
174  const vec& background,
175  const Points& points,
176  const vec& pobs,
177  const vec& pratios,
178  const vec& pbackground,
179  const StructureFunction& structure,
180  int max_points,
181  bool allow_extrapolation=true);
182 
196  vec2 optimal_interpolation_full(const Grid& bgrid,
197  const vec2& background,
198  const vec2& bvariance,
199  const Points& points,
200  const vec& obs,
201  const vec& obs_variance,
202  const vec& background_at_points,
203  const vec& bvariance_at_points,
204  const StructureFunction& structure,
205  int max_points,
206  vec2& analysis_variance,
207  bool allow_extrapolation=true);
208 
222  vec optimal_interpolation_full(const Points& bpoints,
223  const vec& background,
224  const vec& bvariance,
225  const Points& points,
226  const vec& obs,
227  const vec& obs_variance,
228  const vec& background_at_points,
229  const vec& bvariance_at_points,
230  const StructureFunction& structure,
231  int max_points,
232  vec& analysis_sigmas,
233  bool allow_extrapolation=true);
234 
243  vec3 optimal_interpolation_ensi(const Grid& bgrid,
244  const vec3& background,
245  const Points& points,
246  const vec& pobs,
247  const vec& psigmas,
248  const vec2& pbackground,
249  const StructureFunction& structure,
250  int max_points,
251  bool allow_extrapolation=true);
252 
253  vec2 optimal_interpolation_ensi(const Points& bpoints,
254  const vec2& background,
255  const Points& points,
256  const vec& pobs,
257  const vec& psigmas,
258  const vec2& pbackground,
259  const StructureFunction& structure,
260  int max_points,
261  bool allow_extrapolation=true);
262 
275  vec2 local_distribution_correction(const Grid& bgrid,
276  const vec2& background,
277  const Points& points,
278  const vec& pobs,
279  const vec& pbackground,
280  const StructureFunction& structure,
281  float min_quantile,
282  float max_quantile,
283  int min_points=0);
284 
297  vec2 local_distribution_correction(const Grid& bgrid,
298  const vec2& background,
299  const Points& points,
300  const vec2& pobs,
301  const vec2& pbackground,
302  const StructureFunction& structure,
303  float min_quantile,
304  float max_quantile,
305  int min_points=0);
306 
313  vec2 fill(const Grid& igrid, const vec2& input, const Points& points, const vec& radii, float value, bool outside);
314 
315 
324  vec2 doping_square(const Grid& igrid, const vec2& background, const Points& points, const vec& observations, const ivec& halfwidths, float max_elev_diff=gridpp::MV);
325 
334  vec2 doping_circle(const Grid& igrid, const vec2& background, const Points& points, const vec& observations, const vec& radii, float max_elev_diff=gridpp::MV);
335 
346  vec gamma_inv(const vec& levels, const vec& shape, const vec& scale);
347 
360  vec2 neighbourhood(const vec2& input, int halfwidth, Statistic statistic);
361 
367  vec2 neighbourhood(const vec3& input, int halfwidth, Statistic statistic);
368 
374  vec2 neighbourhood_quantile(const vec2& input, float quantile, int halfwidth);
375 
381  vec2 neighbourhood_quantile(const vec3& input, float quantile, int halfwidth);
382 
389  vec2 neighbourhood_quantile_fast(const vec2& input, float quantile, int halfwidth, const vec& thresholds);
390 
397  vec2 neighbourhood_quantile_fast(const vec3& input, float quantile, int halfwidth, const vec& thresholds);
398 
405  vec2 neighbourhood_quantile_fast(const vec2& input, const vec2& quantile, int halfwidth, const vec& thresholds);
406 
413  vec2 neighbourhood_quantile_fast(const vec3& input, const vec2& quantile, int halfwidth, const vec& thresholds);
414 
420  vec2 neighbourhood_brute_force(const vec2& input, int halfwidth, Statistic statistic);
421 
427  vec2 neighbourhood_brute_force(const vec3& input, int halfwidth, Statistic statistic);
428 
433  vec get_neighbourhood_thresholds(const vec2& input, int num_thresholds);
434 
439  vec get_neighbourhood_thresholds(const vec3& input, int num_thresholds);
440 
443  vec2 neighbourhood_ens(const vec3& input, int halfwidth, Statistic statistic);
446  vec2 neighbourhood_quantile_ens(const vec3& input, float quantile, int halfwidth);
449  vec2 neighbourhood_quantile_ens_fast(const vec3& input, float quantile, int radius, const vec& thresholds);
464  vec quantile_mapping_curve(const vec& ref, const vec& fcst, vec& output_fcst, vec quantiles=vec());
465 
474  vec metric_optimizer_curve(const vec& ref, const vec& fcst, const vec& thresholds, Metric metric, vec& output_fcst);
475 
484  float apply_curve(float fcst, const vec& curve_ref, const vec& curve_fcst, Extrapolation policy_below, Extrapolation policy_above);
485 
494  vec apply_curve(const vec& fcst, const vec& curve_ref, const vec& curve_fcst, Extrapolation policy_below, Extrapolation policy_above);
495 
504  vec2 apply_curve(const vec2& fcst, const vec& curve_ref, const vec& curve_fcst, Extrapolation policy_below, Extrapolation policy_above);
505 
514  vec2 apply_curve(const vec2& fcst, const vec3& curve_ref, const vec3& curve_fcst, Extrapolation policy_below, Extrapolation policy_above);
515 
522  vec monotonize_curve(vec curve_ref, vec curve_fcst, vec& output_fcst);
523 
524  float get_optimal_threshold(const vec& ref, const vec& fcst, float threshold, Metric metric);
525 
526  float calc_score(float a, float b, float c, float d, Metric metric);
527  float calc_score(const vec& ref, const vec& fcst, float threshold, Metric metric);
528  float calc_score(const vec& ref, const vec& fcst, float threshold, float fthreshold, Metric metric);
529 
544  vec2 downscaling(const Grid& igrid, const Grid& ogrid, const vec2& ivalues, Downscaler downscaler);
545  vec3 downscaling(const Grid& igrid, const Grid& ogrid, const vec3& ivalues, Downscaler downscaler);
546  vec downscaling(const Grid& igrid, const Points& opoints, const vec2& ivalues, Downscaler downscaler);
547  vec2 downscaling(const Grid& igrid, const Points& opoints, const vec3& ivalues, Downscaler downscaler);
548  // vec downscaling(const Points& ipoints, const Points& opoints, const vec& ivalues, Downscaler downscaler);
549  // vec2 downscaling(const Points& ipoints, const Points& opoints, const vec2& ivalues, Downscaler downscaler);
550  // vec2 downscaling(const Points& ipoints, const Grid& ogrid, const vec& ivalues, Downscaler downscaler);
551  // vec3 downscaling(const Points& ipoints, const Grid& ogrid, const vec2& ivalues, Downscaler downscaler);
552 
559  vec2 nearest(const Grid& igrid, const Grid& ogrid, const vec2& ivalues);
560  vec3 nearest(const Grid& igrid, const Grid& ogrid, const vec3& ivalues);
561 
568  vec nearest(const Grid& igrid, const Points& opoints, const vec2& ivalues);
569  vec2 nearest(const Grid& igrid, const Points& opoints, const vec3& ivalues);
570 
577  vec nearest(const Points& ipoints, const Points& opoints, const vec& ivalues);
578  vec2 nearest(const Points& ipoints, const Points& opoints, const vec2& ivalues);
579 
586  vec2 nearest(const Points& ipoints, const Grid& ogrid, const vec& ivalues);
587  vec3 nearest(const Points& ipoints, const Grid& ogrid, const vec2& ivalues);
588 
597  vec2 downscale_probability(const Grid& igrid, const Grid& ogrid, const vec3& ivalues, const vec2& threshold, const ComparisonOperator& comparison_operator);
598 
610  vec2 mask_threshold_downscale_consensus(const Grid& igrid, const Grid& ogrid, const vec3& ivalues_true, const vec3& ivalues_false, const vec3& theshold_values, const vec2& threshold, const ComparisonOperator& comparison_operator, const Statistic& statistic);
611 
618  vec2 bilinear(const Grid& igrid, const Grid& ogrid, const vec2& ivalues);
619  vec3 bilinear(const Grid& igrid, const Grid& ogrid, const vec3& ivalues);
620 
627  vec bilinear(const Grid& igrid, const Points& opoints, const vec2& ivalues);
628  vec2 bilinear(const Grid& igrid, const Points& opoints, const vec3& ivalues);
629 
630  vec2 simple_gradient(const Grid& igrid, const Grid& ogrid, const vec2& ivalues, float elev_gradient, Downscaler downscaler=Nearest);
631  vec3 simple_gradient(const Grid& igrid, const Grid& ogrid, const vec3& ivalues, float elev_gradient, Downscaler downscaler=Nearest);
632 
633  vec simple_gradient(const Grid& igrid, const Points& opoints, const vec2& ivalues, float elev_gradient, Downscaler downscaler=Nearest);
634  vec2 simple_gradient(const Grid& igrid, const Points& opoints, const vec3& ivalues, float elev_gradient, Downscaler downscaler=Nearest);
635 
643  vec2 full_gradient(const Grid& igrid, const Grid& ogrid, const vec2& ivalues, const vec2& elev_gradient, const vec2& laf_gradient=vec2(), Downscaler downscaler=Nearest);
644  vec3 full_gradient(const Grid& igrid, const Grid& ogrid, const vec3& ivalues, const vec3& elev_gradient, const vec3& laf_gradient, Downscaler downscaler=Nearest);
645  vec full_gradient(const Grid& igrid, const Points& opoints, const vec2& ivalues, const vec2& elev_gradient, const vec2& laf_gradient, Downscaler downscaler=Nearest);
646  vec2 full_gradient(const Grid& igrid, const Points& opoints, const vec3& ivalues, const vec3& elev_gradient, const vec3& laf_gradient, Downscaler downscaler=Nearest);
647 
648  /* Elevation and land area fraction downscaling with debug output fields
649  */
650  vec3 full_gradient_debug(const Grid& igrid, const Grid& ogrid, const vec2& ivalues, const vec2& elev_gradient, const vec2& laf_gradient=vec2(), Downscaler downscaler=Nearest);
651 
660  vec2 smart(const Grid& igrid, const Grid& ogrid, const vec2& ivalues, int num, const StructureFunction& structure);
674  vec count(const Grid& grid, const Points& points, float radius);
675 
682  vec2 count(const Grid& igrid, const Grid& ogrid, float radius);
683 
690  vec2 count(const Points& points, const Grid& grid, float radius);
691 
698  vec count(const Points& ipoints, const Points& opoints, float radius);
699 
706  vec distance(const Grid& grid, const Points& points, int num=1);
707 
714  vec2 distance(const Grid& igrid, const Grid& ogrid, int num=1);
715 
722  vec2 distance(const Points& points, const Grid& grid, int num=1);
723 
730  vec distance(const Points& ipoints, const Points& opoint, int num=1);
731 
736  vec2 fill_missing(const vec2& values);
737 
746  vec2 gridding(const Grid& grid, const Points& points, const vec& values, float radius, int min_num, Statistic statistic);
747 
755  vec2 gridding_nearest(const Grid& grid, const Points& points, const vec& values, int min_num, gridpp::Statistic statistic);
756 
769  float dewpoint(float temperature, float relative_humidity);
770 
776  vec dewpoint(const vec& temperature, const vec& relative_humidity);
777 
785  float pressure(float ielev, float oelev, float ipressure, float itemperature=288.15);
786 
794  vec pressure(const vec& ielev, const vec& oelev, const vec& ipressure, const vec& itemperature);
795 
804  float sea_level_pressure(float ps, float altitude, float temperature, float rh=gridpp::MV, float dewpoint=gridpp::MV);
805 
814  vec sea_level_pressure(const vec& ps, const vec& altitude, const vec& temperature, const vec& rh, const vec& dewpoint);
815 
821  float qnh(float pressure, float altitude);
822 
828  vec qnh(const vec& pressure, const vec& altitude);
829 
835  float relative_humidity(float temperature, float dewpoint);
836 
842  vec relative_humidity(const vec& temperature, const vec& dewpoint);
843 
850  float wetbulb(float temperature, float pressure, float relative_humidity);
851 
858  vec wetbulb(const vec& temperature, const vec& pressure, const vec& relative_humidity);
859 
865  float wind_speed(float xwind, float ywind);
866 
872  vec wind_speed(const vec& xwind, const vec& ywind);
873 
880  float wind_direction(float xwind, float ywind);
881 
887  vec wind_direction(const vec& xwind, const vec& ywind);
888 
896  void set_omp_threads(int num);
897 
899  int get_omp_threads();
900 
902  void initialize_omp();
903 
908  // ivec regression(const vec& x, const vec& y);
909 
911  void set_debug_level(int level);
912 
913  static int _debug_level = 0;
914 
916  int get_debug_level();
917 
919  Statistic get_statistic(std::string name);
920 
924  std::string version();
925 
926  double clock();
927  void debug(std::string string);
928  void warning(std::string string);
929  void error(std::string string);
930  void future_deprecation_warning(std::string function, std::string other="");
931  bool is_valid(float value);
932  float calc_statistic(const vec& array, Statistic statistic);
933  float calc_quantile(const vec& array, float quantile);
934  vec calc_statistic(const vec2& array, Statistic statistic);
935  vec calc_quantile(const vec2& array, float quantile=MV);
936 
942  vec2 calc_quantile(const vec3& array, const vec2& quantile);
943  int num_missing_values(const vec2& iArray);
944 
950  int get_lower_index(float iX, const vec& iValues);
951 
957  int get_upper_index(float iX, const vec& iValues);
958 
970  float interpolate(float x, const vec& iX, const vec& iY);
971 
973  ivec2 init_ivec2(int Y, int X, int value);
974  vec2 init_vec2(int Y, int X, float value=MV);
975 
977  ivec3 init_ivec3(int Y, int X, int E, int value);
978  vec3 init_vec3(int Y, int X, int E, float value=MV);
979 
986  vec calc_even_quantiles(const vec& values, int num);
987 
997  vec2 calc_gradient(const vec2& base, const vec2& values, GradientType gradient_type, int halfwidth, int min_num=2, float min_range=gridpp::MV, float default_gradient=0);
1010  vec2 neighbourhood_search(const vec2& array, const vec2& search_array, int halfwidth, float search_target_min, float search_target_max, float search_delta, const ivec2& apply_array=ivec2());
1011 
1021  vec2 window(const vec2& array, int length, gridpp::Statistic statistic, bool before=false, bool keep_missing=false, bool missing_edges=true);
1022 
1025  bool compatible_size(const Grid& grid, const vec2& v);
1026  bool compatible_size(const Grid& grid, const vec3& v);
1027  bool compatible_size(const Points& points, const vec& v);
1028  bool compatible_size(const Points& points, const vec2& v);
1029  bool compatible_size(const vec2& a, const vec2& b);
1030  bool compatible_size(const vec2& a, const vec3& b);
1031  bool compatible_size(const vec3& a, const vec3& b);
1032 
1042  bool point_in_rectangle(const Point& A, const Point& B, const Point& C, const Point& D, const Point& m );
1043 
1049  float* test_array(float* v, int n);
1051  float test_vec_input(const vec& input);
1053  int test_ivec_input(const ivec& input);
1055  float test_vec2_input(const vec2& input);
1057  float test_vec3_input(const vec3& input);
1059  vec test_vec_output();
1060  ivec test_ivec_output();
1062  vec2 test_vec2_output();
1063  ivec2 test_ivec2_output();
1065  vec3 test_vec3_output();
1066  ivec3 test_ivec3_output();
1067 
1069  float test_vec_argout(vec& distances);
1071  float test_vec2_argout(vec2& distances);
1072 
1074 
1076  static const float swig_default_value = -1;
1077 
1081  class Point {
1082  public:
1090  Point(float lat, float lon, float elev=MV, float laf=MV, CoordinateType type=Geodetic);
1091  float lat;
1092  float lon;
1093  float elev;
1094  float laf;
1096  };
1097 
1099  class KDTree {
1100  public:
1101  KDTree(vec lats, vec lons, CoordinateType type=Geodetic);
1102  KDTree& operator=(KDTree other);
1103  KDTree(const KDTree& other);
1105 
1110  int get_nearest_neighbour(float lat, float lon, bool include_match=true) const;
1111 
1117  ivec get_neighbours(float lat, float lon, float radius, bool include_match=true) const;
1118 
1125  ivec get_neighbours_with_distance(float lat, float lon, float radius, vec& distances, bool include_match=true) const;
1126 
1132  int get_num_neighbours(float lat, float lon, float radius, bool include_match=true) const;
1133 
1139  ivec get_closest_neighbours(float lat, float lon, int num, bool include_match=true) const;
1140 
1141 
1149  bool convert_coordinates(const vec& lats, const vec& lons, vec& x_coords, vec& y_coords, vec& z_coords) const;
1150 
1158  bool convert_coordinates(float lat, float lon, float& x_coord, float& y_coord, float& z_coord) const;
1159  static float deg2rad(float deg);
1160  static float rad2deg(float deg);
1161  static float calc_distance(float lat1, float lon1, float lat2, float lon2, CoordinateType type=Geodetic);
1162  static float calc_distance(float x0, float y0, float z0, float x1, float y1, float z1);
1163  static float calc_distance_fast(float lat1, float lon1, float lat2, float lon2, CoordinateType type=Geodetic);
1164  static float calc_distance_fast(const Point& p1, const Point& p2);
1165  vec get_lats() const;
1166  vec get_lons() const;
1167  int size() const;
1168  CoordinateType get_coordinate_type() const;
1169  protected:
1170  typedef boost::geometry::model::point<float, 3, boost::geometry::cs::cartesian> point;
1171  typedef std::pair<point, unsigned> value;
1172  typedef boost::geometry::model::box<point> box;
1173  boost::geometry::index::rtree< value, boost::geometry::index::quadratic<16> > mTree;
1174  vec mLats;
1175  vec mLons;
1177 
1178  struct within_radius {
1179  public:
1180  within_radius(point p, float radius, bool include_match);
1181  bool operator()(value const& v) const;
1182  private:
1183  float radius;
1184  point p;
1185  bool include_match;
1186  };
1187  struct is_not_equal {
1188  public:
1189  is_not_equal(point p);
1190  bool operator()(value const& v) const;
1191  private:
1192  point p;
1193  };
1194  };
1195 
1197  class Points {
1198  public:
1199  Points();
1207  Points(vec lats, vec lons, vec elevs=vec(), vec lafs=vec(), CoordinateType type=Geodetic);
1208  Points(KDTree tree, vec elevs=vec(), vec lafs=vec());
1209  Points& operator=(Points other);
1210  Points(const Points& other);
1211  // Returns -1 if there are no neighbours
1212  int get_nearest_neighbour(float lat, float lon, bool include_match=true) const;
1213  ivec get_neighbours(float lat, float lon, float radius, bool include_match=true) const;
1214  ivec get_neighbours_with_distance(float lat, float lon, float radius, vec& distances, bool include_match=true) const;
1215  int get_num_neighbours(float lat, float lon, float radius, bool include_match=true) const;
1216  ivec get_closest_neighbours(float lat, float lon, int num, bool include_match=true) const;
1217 
1218  vec get_lats() const;
1219  vec get_lons() const;
1220  vec get_elevs() const;
1221  vec get_lafs() const;
1222  int size() const;
1223  ivec get_in_domain_indices(const Grid& grid) const;
1224  Points get_in_domain(const Grid& grid) const;
1225  CoordinateType get_coordinate_type() const;
1226  Point get_point(int index) const;
1227  Points subset(const ivec& indices) const;
1228  private:
1229  KDTree mTree;
1230  vec mLats;
1231  vec mLons;
1232  vec mElevs;
1233  vec mLafs;
1234  };
1235 
1237  class Grid {
1238  public:
1239  Grid();
1240 
1248  Grid(vec2 lats, vec2 lons, vec2 elevs=vec2(), vec2 lafs=vec2(), CoordinateType type=Geodetic);
1249  ivec get_nearest_neighbour(float lat, float lon, bool include_match=true) const;
1250  ivec2 get_neighbours(float lat, float lon, float radius, bool include_match=true) const;
1251  ivec2 get_neighbours_with_distance(float lat, float lon, float radius, vec& distances, bool include_match=true) const;
1252  int get_num_neighbours(float lat, float lon, float radius, bool include_match=true) const;
1253  ivec2 get_closest_neighbours(float lat, float lon, int num, bool include_match=true) const;
1254 
1255  bool get_box(float lat, float lon, int& Y1_out, int& X1_out, int& Y2_out, int& X2_out) const;
1256 
1258  Points to_points() const;
1259 
1260  vec2 get_lats() const;
1261  vec2 get_lons() const;
1262  vec2 get_elevs() const;
1263  vec2 get_lafs() const;
1264  ivec size() const;
1265  CoordinateType get_coordinate_type() const;
1266  Point get_point(int y_index, int x_index) const;
1267  private:
1268  KDTree mTree;
1269  int mX;
1270  vec2 get_2d(vec input) const;
1271  ivec get_indices(int index) const;
1272  ivec2 get_indices(ivec indices) const;
1273  vec2 mLats;
1274  vec2 mLons;
1275  vec2 mElevs;
1276  vec2 mLafs;
1277  };
1278  class not_implemented_exception: public std::logic_error
1279  {
1280  public:
1281  not_implemented_exception() : std::logic_error("Function not yet implemented") { };
1282  };
1283 
1286  public:
1287  StructureFunction(float localization_distance=0);
1289  virtual float corr(const Point& p1, const Point& p2) const = 0;
1290  virtual float corr_background(const Point& p1, const Point& p2) const;
1294  virtual float localization_distance(const Point& p) const;
1295  virtual StructureFunction* clone() const = 0;
1296  static const float default_min_rho;
1297  protected:
1302  float barnes_rho(float dist, float length) const;
1303 
1308  float cressman_rho(float dist, float length) const;
1310  };
1312  public:
1319  MultipleStructure(const StructureFunction& structure_h, const StructureFunction& structure_v, const StructureFunction& structure_w);
1320  float corr(const Point& p1, const Point& p2) const;
1321  StructureFunction* clone() const;
1322  float localization_distance(const Point& p) const;
1323  private:
1324  StructureFunction* m_structure_h;
1325  StructureFunction* m_structure_v;
1326  StructureFunction* m_structure_w;
1327  };
1330  public:
1337  BarnesStructure(float h, float v=0, float w=0, float hmax=MV);
1345  BarnesStructure(Grid grid, vec2 h, vec2 v, vec2 w, float min_rho=StructureFunction::default_min_rho);
1346  float corr(const Point& p1, const Point& p2) const;
1347  StructureFunction* clone() const;
1348  float localization_distance(const Point& p) const;
1349  private:
1350  Grid m_grid;
1351  vec2 mH;
1352  vec2 mV;
1353  vec2 mW;
1354  float m_min_rho;
1355  bool m_is_spatial;
1356  };
1357 
1360  public:
1361  CressmanStructure(float h, float v=0, float w=0);
1362  float corr(const Point& p1, const Point& p2) const;
1363  StructureFunction* clone() const;
1364  private:
1365  float mH;
1366  float mV;
1367  float mW;
1368  };
1369 
1371  public:
1376  CrossValidation(StructureFunction& structure, float dist=MV);
1377  float corr(const Point& p1, const Point& p2) const;
1378  float corr_background(const Point& p1, const Point& p2) const;
1379  StructureFunction* clone() const;
1380  float localization_distance(const Point& p) const;
1381  private:
1382  StructureFunction* m_structure;
1383  float m_dist;
1384  };
1385 
1386  class Transform {
1387  public:
1388  // Note these cannot be pure virtual, otherwise SWIG does not expose
1389  // the vector functions (with the same name) in python. Therefore, make sure these
1390  // functions are overloaded in the subclass implementation
1391  virtual float forward(float value) const;
1392  virtual float backward(float value) const;
1393 
1394  vec forward(const vec& input) const;
1395  vec backward(const vec& input) const;
1396  vec2 forward(const vec2& input) const;
1397  vec2 backward(const vec2& input) const;
1398  vec3 forward(const vec3& input) const;
1399  vec3 backward(const vec3& input) const;
1400  };
1401  class Identity : public Transform {
1402  public:
1403  // SWIG requires these "using" statements to enable the vectorized versions in the
1404  // subclasses
1405  using Transform::forward;
1406  using Transform::backward;
1407  float forward(float value) const;
1408  float backward(float value) const;
1409  };
1410  class Log : public Transform {
1411  public:
1412  using Transform::forward;
1413  using Transform::backward;
1414  float forward(float value) const;
1415  float backward(float value) const;
1416  };
1417  class BoxCox : public Transform {
1418  public:
1419  BoxCox(float threshold);
1420  using Transform::forward;
1421  using Transform::backward;
1422  float forward(float value) const;
1423  float backward(float value) const;
1424  private:
1425  float mThreshold;
1426  };
1427  class Gamma : public Transform {
1431  public:
1432  Gamma(float shape, float scale, float tolerance=0.01);
1433  using Transform::forward;
1434  using Transform::backward;
1435  float forward(float value) const;
1436  float backward(float value) const;
1437  private:
1438  float m_tolerance;
1439  boost::math::gamma_distribution<> m_gamma_dist;
1440  boost::math::normal m_norm_dist;
1441  };
1442 };
1443 #endif
ComparisonOperator
Types of comparison operators.
Definition: gridpp.h:129
Definition: gridpp.h:1278
float test_vec2_input(const vec2 &input)
Testing function for 2D input vector.
Definition: swig.cpp:25
vec2 nearest(const Grid &igrid, const Grid &ogrid, const vec2 &ivalues)
Nearest neighbour dowscaling grid to grid.
Definition: nearest.cpp:7
std::string version()
The gridpp version.
Definition: gridpp.cpp:8
vec2 calc_gradient(const vec2 &base, const vec2 &values, GradientType gradient_type, int halfwidth, int min_num=2, float min_range=gridpp::MV, float default_gradient=0)
Computes gradients based on values in neighbourhood.
Definition: calc_gradient.cpp:6
Definition: gridpp.h:1178
Bias.
Definition: gridpp.h:99
vec2 neighbourhood_quantile_ens_fast(const vec3 &input, float quantile, int radius, const vec &thresholds)
Deprecated: Compute neighbourhood quantiles fast on ensemble field.
Definition: neighbourhood.cpp:549
float m_localization_distance
Definition: gridpp.h:1309
float laf
Definition: gridpp.h:1094
Greater or equal than, >=.
Definition: gridpp.h:133
vec2 optimal_interpolation(const Grid &bgrid, const vec2 &background, const Points &points, const vec &pobs, const vec &pratios, const vec &pbackground, const StructureFunction &structure, int max_points, bool allow_extrapolation=true)
Optimal interpolation for a deterministic gridded field.
Definition: oi.cpp:26
static const float MV_CML
Missing value indicator in gridpp command-line tool.
Definition: gridpp.h:43
Helper class for Grid and Points.
Definition: gridpp.h:1099
void test_not_implemented_exception()
Definition: swig.cpp:98
vec2 bilinear(const Grid &igrid, const Grid &ogrid, const vec2 &ivalues)
Bilinear downscaling grid to grid.
Definition: bilinear.cpp:27
static const float standard_surface_temperature
Temperature at surface in standard atmosphere [K].
Definition: gridpp.h:51
Minimum of values.
Definition: gridpp.h:82
vec2 window(const vec2 &array, int length, gridpp::Statistic statistic, bool before=false, bool keep_missing=false, bool missing_edges=true)
Compute window statistics.
Definition: window.cpp:6
float dewpoint(float temperature, float relative_humidity)
Calculate dewpoint temperature from temperature and relative humidity.
Definition: humidity.cpp:5
Maximum of values.
Definition: gridpp.h:84
float wind_speed(float xwind, float ywind)
Diagnose wind speed from its components.
Definition: wind.cpp:6
Definition: gridpp.h:1401
CorrectionType
Method for statistical correction.
Definition: gridpp.h:104
std::vector< int > ivec
Definition: gridpp.h:27
ivec2 test_ivec2_output()
Definition: swig.cpp:69
ivec test_ivec_output()
Definition: swig.cpp:65
vec2 test_vec2_output()
Testing function for 2D output vector.
Definition: swig.cpp:49
Unknown statistic.
Definition: gridpp.h:90
Metric
Binary verification metrics.
Definition: gridpp.h:94
Simple structure function based on distance, elevation, and land area fraction.
Definition: gridpp.h:1329
virtual float forward(float value) const
Definition: transform.cpp:6
GradientType
Types of methods to calculate the gradient.
Definition: gridpp.h:117
Definition: gridpp.h:1370
vec2 full_gradient(const Grid &igrid, const Grid &ogrid, const vec2 &ivalues, const vec2 &elev_gradient, const vec2 &laf_gradient=vec2(), Downscaler downscaler=Nearest)
Compute Downscale.
Definition: gradient.cpp:5
vec distance(const Grid &grid, const Points &points, int num=1)
For each point, calculates the distance to nearest gridpoint.
Definition: distance.cpp:6
double clock()
Definition: util.cpp:234
static const float gas_constant_mol
Universal Gas Constant [kg*m^2*s^-2/(K*mol)].
Definition: gridpp.h:57
vec2 downscale_probability(const Grid &igrid, const Grid &ogrid, const vec3 &ivalues, const vec2 &threshold, const ComparisonOperator &comparison_operator)
Nearest neighbour downscaling grid to grid and probability in one.
Definition: downscale_probability.cpp:7
Extrapolation
Methods for extrapolating outside a curve.
Definition: gridpp.h:71
float test_vec_argout(vec &distances)
Testing function for 1D vector treated as output.
Definition: swig.cpp:86
void set_debug_level(int level)
Set the verbosity of debug messages.
Definition: gridpp.cpp:67
CoordinateType type
Definition: gridpp.h:1095
not_implemented_exception()
Definition: gridpp.h:1281
static const float gas_constant_si
Universal Gas Constant [J/(kg*K)].
Definition: gridpp.h:59
Proportion correct.
Definition: gridpp.h:98
float sea_level_pressure(float ps, float altitude, float temperature, float rh=gridpp::MV, float dewpoint=gridpp::MV)
Convert Surface Pressure to Sea Level Pressure.
Definition: pressure.cpp:28
Population variance of values.
Definition: gridpp.h:87
Point(float lat, float lon, float elev=MV, float laf=MV, CoordinateType type=Geodetic)
Constructor.
Definition: point.cpp:5
vec2 smart(const Grid &igrid, const Grid &ogrid, const vec2 &ivalues, int num, const StructureFunction &structure)
Smart neighbour downscaling grid to grid.
Definition: smart.cpp:12
void set_omp_threads(int num)
Set the number of OpenMP threads to use.
Definition: gridpp.cpp:54
float relative_humidity(float temperature, float dewpoint)
Calculate relative humidity from temperature and dewpoint temperature.
Definition: humidity.cpp:33
ivec3 test_ivec3_output()
Definition: swig.cpp:76
vec2 local_distribution_correction(const Grid &bgrid, const vec2 &background, const Points &points, const vec &pobs, const vec &pbackground, const StructureFunction &structure, float min_quantile, float max_quantile, int min_points=0)
Correction of a gridded field ensuring the distribution of values nearby match that of observations...
Definition: local_distribution_correction.cpp:18
Statistic get_statistic(std::string name)
Convert name of a statistic enum.
Definition: gridpp.cpp:11
int get_lower_index(float iX, const vec &iValues)
Find the index in a vector that is equal or just below a value.
Definition: util.cpp:291
vec2 neighbourhood_quantile(const vec2 &input, float quantile, int halfwidth)
Computes a quantile in a sliding square neighbourhood.
Definition: neighbourhood.cpp:534
boost::geometry::model::point< float, 3, boost::geometry::cs::cartesian > point
Definition: gridpp.h:1170
Threat score.
Definition: gridpp.h:96
vec2 fill(const Grid &igrid, const vec2 &input, const Points &points, const vec &radii, float value, bool outside)
Fill in values inside or outside a set of circles (useful for masking)
Definition: fill.cpp:6
Represents a vector of locations and their metadata.
Definition: gridpp.h:1197
vec2 doping_circle(const Grid &igrid, const vec2 &background, const Points &points, const vec &observations, const vec &radii, float max_elev_diff=gridpp::MV)
Insert observations into gridded field using a circle.
Definition: doping.cpp:50
vec2 neighbourhood_quantile_fast(const vec2 &input, float quantile, int halfwidth, const vec &thresholds)
Fast and approximate neighbourhood quantile.
Definition: neighbourhood.cpp:296
std::pair< point, unsigned > value
Definition: gridpp.h:1171
std::vector< ivec2 > ivec3
Definition: gridpp.h:29
static const float molar_mass
Molar Mass of Dry Air [kg/mol].
Definition: gridpp.h:55
Definition: gridpp.h:1386
Keep values the way they were.
Definition: gridpp.h:76
Greater than, >
Definition: gridpp.h:132
Definition: gridpp.h:1187
vec2 neighbourhood_brute_force(const vec2 &input, int halfwidth, Statistic statistic)
Spatial neighbourhood filter without any shortcuts.
Definition: neighbourhood.cpp:528
ivec3 init_ivec3(int Y, int X, int E, int value)
Initialize a vector of size Y, X, E, with a given value.
Definition: util.cpp:436
std::vector< vec > vec2
Definition: gridpp.h:25
Latitude and longitude.
Definition: gridpp.h:112
Represents a 2D grid of locations and their metadata.
Definition: gridpp.h:1237
static int _debug_level
Definition: gridpp.h:913
Simple structure function based on distance, elevation, and land area fraction.
Definition: gridpp.h:1359
vec2 neighbourhood(const vec2 &input, int halfwidth, Statistic statistic)
Spatial neighbourhood filter, computing a statistic for a sliding square window.
Definition: neighbourhood.cpp:28
Covariance structure function.
Definition: gridpp.h:1285
int get_upper_index(float iX, const vec &iValues)
Find the index in a vector that is equal or just above a value.
Definition: util.cpp:310
Heidke skill score.
Definition: gridpp.h:100
Continue past the end-points using the mean slope of the curve.
Definition: gridpp.h:73
void debug(std::string string)
Definition: util.cpp:205
void initialize_omp()
Sets the number of OpenMP threads to 1 if OMP_NUM_THREADS undefined.
Definition: gridpp.cpp:42
float get_optimal_threshold(const vec &ref, const vec &fcst, float threshold, Metric metric)
Definition: metric_optimizer.cpp:131
static const float swig_default_value
Default value used to fill array in SWIG testing functions.
Definition: gridpp.h:1076
vec count(const Grid &grid, const Points &points, float radius)
For each point, counts the number of gridpoints within the radius.
Definition: count.cpp:6
Lower than, <.
Definition: gridpp.h:130
Continue past the end-points using the slope of the two lowermost or uppermost points in the curve...
Definition: gridpp.h:74
vec2 optimal_interpolation_full(const Grid &bgrid, const vec2 &background, const vec2 &bvariance, const Points &points, const vec &obs, const vec &obs_variance, const vec &background_at_points, const vec &bvariance_at_points, const StructureFunction &structure, int max_points, vec2 &analysis_variance, bool allow_extrapolation=true)
Optimal interpolation for a deterministic gridded field including analysis variance.
Definition: oi.cpp:326
boost::geometry::model::box< point > box
Definition: gridpp.h:1172
float test_vec2_argout(vec2 &distances)
Testing function for 2D vector treated as output.
Definition: swig.cpp:91
vec2 gridding_nearest(const Grid &grid, const Points &points, const vec &values, int min_num, gridpp::Statistic statistic)
Assign each point to nearest neighbour in grid and aggregate values.
Definition: gridding.cpp:31
float test_vec_input(const vec &input)
Testing function for 1D input vector.
Definition: swig.cpp:12
CoordinateType
Types of coordinates for position of points.
Definition: gridpp.h:111
vec test_vec_output()
Testing function for 1D output vector.
Definition: swig.cpp:45
vec quantile_mapping_curve(const vec &ref, const vec &fcst, vec &output_fcst, vec quantiles=vec())
Create quantile mapping calibration curve.
Definition: quantile_mapping.cpp:5
vec3 init_vec3(int Y, int X, int E, float value=MV)
Definition: util.cpp:427
vec mLats
Definition: gridpp.h:1174
Definition: gridpp.h:1311
bool point_in_rectangle(const Point &A, const Point &B, const Point &C, const Point &D, const Point &m)
Checks if a point is located inside a rectangle formed by 4 points.
Definition: util.cpp:511
float wetbulb(float temperature, float pressure, float relative_humidity)
Calculate wetbulb temperature from temperature, pressure, and relative humidity.
Definition: humidity.cpp:91
vec2 mask_threshold_downscale_consensus(const Grid &igrid, const Grid &ogrid, const vec3 &ivalues_true, const vec3 &ivalues_false, const vec3 &theshold_values, const vec2 &threshold, const ComparisonOperator &comparison_operator, const Statistic &statistic)
Nearest neighbour downscaling grid to grid and threshold and consensus in one.
Definition: mask_threshold_downscale_consensus.cpp:7
std::vector< double > dvec
Definition: gridpp.h:32
int num_missing_values(const vec2 &iArray)
Definition: util.cpp:196
float lon
Definition: gridpp.h:1092
std::vector< dvec > dvec2
Definition: gridpp.h:33
vec2 gridding(const Grid &grid, const Points &points, const vec &values, float radius, int min_num, Statistic statistic)
Aggregate points onto a grid.
Definition: gridding.cpp:6
Mean of values.
Definition: gridpp.h:83
Statistic
Statistical operations to reduce a vector to a scalar.
Definition: gridpp.h:80
int test_ivec_input(const ivec &input)
Testing function for 1D input vector.
Definition: swig.cpp:18
vec mLons
Definition: gridpp.h:1175
bool is_valid(float value)
Definition: util.cpp:20
CoordinateType mType
Definition: gridpp.h:1176
Definition: gridpp.h:19
vec2 init_vec2(int Y, int X, float value=MV)
Definition: util.cpp:415
vec calc_even_quantiles(const vec &values, int num)
Get reasonably spaced quantiles from a vector of values, ignoring duplicate values but including the ...
Definition: util.cpp:241
vec2 neighbourhood_ens(const vec3 &input, int halfwidth, Statistic statistic)
Deprecated: Compute neighbourhood statistic on ensemble field.
Definition: neighbourhood.cpp:541
A quantile from values.
Definition: gridpp.h:85
float elev
Definition: gridpp.h:1093
vec3 optimal_interpolation_ensi(const Grid &bgrid, const vec3 &background, const Points &points, const vec &pobs, const vec &psigmas, const vec2 &pbackground, const StructureFunction &structure, int max_points, bool allow_extrapolation=true)
Optimal interpolation using a structure function based on an ensemble See Lussana et al 2019 (DOI: 10...
Definition: oi_ensi.cpp:33
static const float pi
Mathematical constant pi.
Definition: gridpp.h:45
float apply_curve(float fcst, const vec &curve_ref, const vec &curve_fcst, Extrapolation policy_below, Extrapolation policy_above)
Apply arbitrary calibration curve to a single value.
Definition: curve.cpp:6
vec metric_optimizer_curve(const vec &ref, const vec &fcst, const vec &thresholds, Metric metric, vec &output_fcst)
Create calibration curve that optimizes a metric.
Definition: metric_optimizer.cpp:105
static const float default_min_rho
Definition: gridpp.h:1296
virtual float backward(float value) const
Definition: transform.cpp:9
float calc_quantile(const vec &array, float quantile)
Definition: util.cpp:93
Definition: gridpp.h:124
float interpolate(float x, const vec &iX, const vec &iY)
Piecewise linear interpolation If x is outside the range of iX, then the min/max value of iY is used...
Definition: util.cpp:329
Multiplicative.
Definition: gridpp.h:106
ivec2 init_ivec2(int Y, int X, int value)
Initialize a vector of size Y, X, with a given value.
Definition: util.cpp:421
float qnh(float pressure, float altitude)
Diagnose QNH from pressure and altitude.
Definition: qnh.cpp:6
Standard deviation of values.
Definition: gridpp.h:86
Definition: gridpp.h:118
static const float MV
Missing value indicator.
Definition: gridpp.h:41
Definition: gridpp.h:119
Downscaler
Types of simple downscaling methods.
Definition: gridpp.h:123
vec get_neighbourhood_thresholds(const vec2 &input, int num_thresholds)
Calculate appropriate approximation thresholds for neighbourhood quantile.
Definition: neighbourhood.cpp:243
float pressure(float ielev, float oelev, float ipressure, float itemperature=288.15)
Calculate pressure at a new elevation.
Definition: pressure.cpp:5
float calc_score(float a, float b, float c, float d, Metric metric)
Definition: metric_optimizer.cpp:209
float wind_direction(float xwind, float ywind)
Diagnose wind direction from its components.
Definition: wind.cpp:20
X and Y.
Definition: gridpp.h:113
vec monotonize_curve(vec curve_ref, vec curve_fcst, vec &output_fcst)
Ensure calibration curve is monotonic, by removing points.
Definition: curve.cpp:126
vec gamma_inv(const vec &levels, const vec &shape, const vec &scale)
Definition: distribution.cpp:5
vec3 test_vec3_output()
Testing function for 3D output vector.
Definition: swig.cpp:55
float calc_statistic(const vec &array, Statistic statistic)
Definition: util.cpp:23
float test_vec3_input(const vec3 &input)
Testing function for 3D input vector.
Definition: swig.cpp:34
void future_deprecation_warning(std::string function, std::string other="")
Definition: util.cpp:226
boost::geometry::index::rtree< value, boost::geometry::index::quadratic< 16 > > mTree
Definition: gridpp.h:1173
vec2 neighbourhood_quantile_ens(const vec3 &input, float quantile, int halfwidth)
Deprecated: Compute neighbourhood quantiles on ensemble field.
Definition: neighbourhood.cpp:545
Represents a single point in some coordinate system.
Definition: gridpp.h:1081
Mean of values.
Definition: gridpp.h:81
void warning(std::string string)
Definition: util.cpp:209
int get_omp_threads()
Get the number of OpenMP threads currently set.
Definition: gridpp.cpp:60
Lower or equal than, <=.
Definition: gridpp.h:131
Additive.
Definition: gridpp.h:107
vec2 fill_missing(const vec2 &values)
Fill in missing values based on nearby values.
Definition: fill.cpp:43
static const float lapse_rate
Constant Lapse Rate moist air standard atmosphere [K/m].
Definition: gridpp.h:49
static const float gravit
Gravitational acceleration [m/s^2].
Definition: gridpp.h:53
float * test_array(float *v, int n)
Special function whose presense is needed for SWIG.
Definition: swig.cpp:6
vec2 neighbourhood_search(const vec2 &array, const vec2 &search_array, int halfwidth, float search_target_min, float search_target_max, float search_delta, const ivec2 &apply_array=ivec2())
Find suitable value in neighbourhood based on a search criteria.
Definition: neighbourhood_search.cpp:7
vec2 downscaling(const Grid &igrid, const Grid &ogrid, const vec2 &ivalues, Downscaler downscaler)
Generic downscaler for simple methods that don&#39;t use information other than the grids.
Definition: downscaling.cpp:21
Continue past the end-points using a slope of 0.
Definition: gridpp.h:75
Continue past the end-points using a slope of 1.
Definition: gridpp.h:72
Definition: gridpp.h:125
std::vector< float > vec
Definition: gridpp.h:24
int get_debug_level()
Get the currently set level of debug messages.
Definition: gridpp.cpp:71
bool compatible_size(const Grid &grid, const vec2 &v)
Check if the grid is the same size as the 2D vector.
Definition: util.cpp:367
vec3 full_gradient_debug(const Grid &igrid, const Grid &ogrid, const vec2 &ivalues, const vec2 &elev_gradient, const vec2 &laf_gradient=vec2(), Downscaler downscaler=Nearest)
Definition: gradient.cpp:276
vec2 doping_square(const Grid &igrid, const vec2 &background, const Points &points, const vec &observations, const ivec &halfwidths, float max_elev_diff=gridpp::MV)
Insert observations into gridded field using a square box.
Definition: doping.cpp:5
Hannsen-Kuiper skill score.
Definition: gridpp.h:97
Count of values.
Definition: gridpp.h:89
float lat
Definition: gridpp.h:1091
Equitable threat score.
Definition: gridpp.h:95
Sum of values.
Definition: gridpp.h:88
static const double radius_earth
Radius of the earth [m].
Definition: gridpp.h:47
Definition: gridpp.h:1410
vec2 simple_gradient(const Grid &igrid, const Grid &ogrid, const vec2 &ivalues, float elev_gradient, Downscaler downscaler=Nearest)
Definition: simple_gradient.cpp:23
Definition: gridpp.h:1417
Definition: gridpp.h:1427
void error(std::string string)
Definition: util.cpp:213
Quantile mapping.
Definition: gridpp.h:105
KDTree(CoordinateType type=Geodetic)
Definition: gridpp.h:1104
std::vector< ivec > ivec2
Definition: gridpp.h:28
std::vector< vec2 > vec3
Definition: gridpp.h:26