Titanlib 0.3.1
Library for quality control algorithms
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 #ifdef _OPENMP
10  #include <omp.h>
11 #endif
12 #include <exception>
13 
14 #define GRIDPP_VERSION "0.6.0.dev4"
15 #define __version__ GRIDPP_VERSION
16 
17 namespace gridpp {
21  // Preferred vector types
22  typedef std::vector<float> vec;
23  typedef std::vector<vec> vec2;
24  typedef std::vector<vec2> vec3;
25  typedef std::vector<int> ivec;
26  typedef std::vector<ivec> ivec2;
27  typedef std::vector<ivec2> ivec3;
28 
29  // Only use when double is required
30  typedef std::vector<double> dvec;
31  typedef std::vector<dvec> dvec2;
39  static const float MV = NAN;
41  static const float MV_CML = -999;
43  static const float pi = 3.14159265;
45  static const double radius_earth = 6.378137e6;
48  class KDTree;
49  class Points;
50  class Grid;
51  class Point;
52  class Nearest;
53  class StructureFunction;
54  class Transform;
55 
58  OneToOne = 0,
59  MeanSlope = 10,
60  NearestSlope = 20,
61  Zero = 30,
62  Unchanged = 40,
63  };
64 
66  enum Statistic {
67  Mean = 0,
68  Min = 10,
69  Median = 20,
70  Max = 30,
71  Quantile = 40,
72  Std = 50,
73  Variance = 60,
74  Sum = 70,
75  Unknown = -1
76  };
77 
79  enum Metric {
80  Ets = 0,
81  Ts = 1,
82  Kss = 20,
83  Pc = 30,
84  Bias = 40,
85  Hss = 50,
86  };
87 
90  Qq = 0,
92  Additive = 20,
93  };
94 
97  Geodetic = 0,
98  Cartesian = 1,
99  };
100 
117  vec2 optimal_interpolation(const Grid& bgrid,
118  const vec2& background,
119  const Points& points,
120  const vec& pobs,
121  const vec& pratios,
122  const vec& pbackground,
123  const StructureFunction& structure,
124  int max_points,
125  bool allow_extrapolation=true);
126 
138  vec optimal_interpolation(const Points& bpoints,
139  const vec& background,
140  const Points& points,
141  const vec& pobs,
142  const vec& pratios,
143  const vec& pbackground,
144  const StructureFunction& structure,
145  int max_points,
146  bool allow_extrapolation=true);
147 
161  vec2 optimal_interpolation_full(const Grid& bgrid,
162  const vec2& background,
163  const vec2& bvariance,
164  const Points& points,
165  const vec& obs,
166  const vec& obs_variance,
167  const vec& background_at_points,
168  const vec& bvariance_at_points,
169  const StructureFunction& structure,
170  int max_points,
171  vec2& analysis_variance,
172  bool allow_extrapolation=true);
173 
187  vec optimal_interpolation_full(const Points& bpoints,
188  const vec& background,
189  const vec& bvariance,
190  const Points& points,
191  const vec& obs,
192  const vec& obs_variance,
193  const vec& background_at_points,
194  const vec& bvariance_at_points,
195  const StructureFunction& structure,
196  int max_points,
197  vec& analysis_sigmas,
198  bool allow_extrapolation=true);
199 
208  vec3 optimal_interpolation_ensi(const Grid& bgrid,
209  const vec3& background,
210  const Points& points,
211  const vec& pobs,
212  const vec& psigmas,
213  const vec2& pbackground,
214  const StructureFunction& structure,
215  int max_points,
216  bool allow_extrapolation=true);
217 
218  vec2 optimal_interpolation_ensi(const Points& bpoints,
219  const vec2& background,
220  const Points& points,
221  const vec& pobs,
222  const vec& psigmas,
223  const vec2& pbackground,
224  const StructureFunction& structure,
225  int max_points,
226  bool allow_extrapolation=true);
227 
240  vec2 local_distribution_correction(const Grid& bgrid,
241  const vec2& background,
242  const Points& points,
243  const vec& pobs,
244  const vec& pbackground,
245  const StructureFunction& structure,
246  float min_quantile,
247  float max_quantile,
248  int min_points=0);
249 
250  vec2 local_distribution_correction(const Grid& bgrid,
251  const vec2& background,
252  const Points& points,
253  const vec2& pobs,
254  const vec2& pbackground,
255  const StructureFunction& structure,
256  float min_quantile,
257  float max_quantile,
258  int min_points=0);
259 
266  vec2 fill(const Grid& igrid, const vec2& input, const Points& points, const vec& radii, float value, bool outside);
267  vec2 doping(const Grid& igrid, const vec2& input, const Points& points, const vec& values, int half_width, float max_elev_diff=gridpp::MV);
268 
281  vec2 neighbourhood(const vec2& input, int halfwidth, Statistic statistic);
282 
288  vec2 neighbourhood(const vec3& input, int halfwidth, Statistic statistic);
289 
295  vec2 neighbourhood_quantile(const vec2& input, float quantile, int halfwidth);
296 
302  vec2 neighbourhood_quantile(const vec3& input, float quantile, int halfwidth);
303 
310  vec2 neighbourhood_quantile_fast(const vec2& input, float quantile, int halfwidth, const vec& thresholds);
311 
318  vec2 neighbourhood_quantile_fast(const vec3& input, float quantile, int halfwidth, const vec& thresholds);
319 
326  vec2 neighbourhood_quantile_fast(const vec2& input, const vec2& quantile, int halfwidth, const vec& thresholds);
327 
334  vec2 neighbourhood_quantile_fast(const vec3& input, const vec2& quantile, int halfwidth, const vec& thresholds);
335 
341  vec2 neighbourhood_brute_force(const vec2& input, int halfwidth, Statistic statistic);
342 
348  vec2 neighbourhood_brute_force(const vec3& input, int halfwidth, Statistic statistic);
349 
354  vec get_neighbourhood_thresholds(const vec2& input, int num_thresholds);
355 
360  vec get_neighbourhood_thresholds(const vec3& input, int num_thresholds);
361 
364  vec2 neighbourhood_ens(const vec3& input, int halfwidth, Statistic statistic);
367  vec2 neighbourhood_quantile_ens(const vec3& input, float quantile, int halfwidth);
370  vec2 neighbourhood_quantile_ens_fast(const vec3& input, float quantile, int radius, const vec& thresholds);
385  vec quantile_mapping_curve(const vec& ref, const vec& fcst, vec& output_fcst, vec quantiles=vec());
386 
395  vec2 metric_optimizer_curve(const vec& ref, const vec& fcst, const vec& thresholds, Metric metric);
396 
405  float apply_curve(float fcst, const vec& curve_ref, const vec& curve_fcst, Extrapolation policy_below, Extrapolation policy_above);
406 
415  vec apply_curve(const vec& fcst, const vec& curve_ref, const vec& curve_fcst, Extrapolation policy_below, Extrapolation policy_above);
416 
425  vec2 apply_curve(const vec2& fcst, const vec& curve_ref, const vec& curve_fcst, Extrapolation policy_below, Extrapolation policy_above);
426 
435  vec2 apply_curve(const vec2& fcst, const vec3& curve_ref, const vec3& curve_fcst, Extrapolation policy_below, Extrapolation policy_above);
436 
443  vec monotonize_curve(vec curve_ref, vec curve_fcst, vec& output_fcst);
444 
445  float get_optimal_threshold(const vec& ref, const vec& fcst, float threshold, Metric metric);
446 
447  float calc_score(float a, float b, float c, float d, Metric metric);
448  float calc_score(const vec& ref, const vec& fcst, float threshold, Metric metric);
449  float calc_score(const vec& ref, const vec& fcst, float threshold, float fthreshold, Metric metric);
450 
451  vec2 correction(const Grid& rgrid, const vec2& rvalues, const Points& npoints, const vec& nvalues, float mean_radius, float outer_radius, float inner_radius, int min_num, int max_num, CorrectionType type, ivec2& count);
452  // Apply correction based on multiple timesteps
453  vec2 correction(const Grid& rgrid, const vec3& rvalues, const Points& npoints, const vec2& nvalues, const vec2& apply_values, float mean_radius, float outer_radius, float inner_radius, int min_num, int max_num, CorrectionType type, ivec2& count);
454 
468  vec2 nearest(const Grid& igrid, const Grid& ogrid, const vec2& ivalues);
469  vec3 nearest(const Grid& igrid, const Grid& ogrid, const vec3& ivalues);
470 
477  vec nearest(const Grid& igrid, const Points& opoints, const vec2& ivalues);
478  vec2 nearest(const Grid& igrid, const Points& opoints, const vec3& ivalues);
479 
486  vec nearest(const Points& ipoints, const Points& opoints, const vec& ivalues);
487  vec2 nearest(const Points& ipoints, const Points& opoints, const vec2& ivalues);
488 
495  vec2 nearest(const Points& ipoints, const Grid& ogrid, const vec& ivalues);
496  vec3 nearest(const Points& ipoints, const Grid& ogrid, const vec2& ivalues);
497 
504  vec2 bilinear(const Grid& igrid, const Grid& ogrid, const vec2& ivalues);
505  vec3 bilinear(const Grid& igrid, const Grid& ogrid, const vec3& ivalues);
506 
513  vec bilinear(const Grid& igrid, const Points& opoints, const vec2& ivalues);
514  vec2 bilinear(const Grid& igrid, const Points& opoints, const vec3& ivalues);
515 
516  vec2 simple_gradient(const Grid& igrid, const Grid& ogrid, const vec2& ivalues, float elev_gradient);
517  vec3 simple_gradient(const Grid& igrid, const Grid& ogrid, const vec3& ivalues, float elev_gradient);
518 
519  vec simple_gradient(const Grid& igrid, const Points& opoints, const vec2& ivalues, float elev_gradient);
520  vec2 simple_gradient(const Grid& igrid, const Points& opoints, const vec3& ivalues, float elev_gradient);
521 
530  vec2 smart(const Grid& igrid, const Grid& ogrid, const vec2& ivalues, int num, const StructureFunction& structure);
544  vec count(const Grid& grid, const Points& points, float radius);
545 
552  vec2 count(const Points& points, const Grid& grid, float radius);
553 
560  vec distance(const Grid& grid, const Points& points, int num=1);
561 
568  vec2 distance(const Grid& igrid, const Grid& ogrid, int num=1);
569 
576  vec2 distance(const Points& points, const Grid& grid, int num=1);
577 
582  vec2 fill_missing(const vec2& values);
583 
592  vec2 gridding(const Grid& grid, const Points& points, const vec& values, float radius, int min_num, Statistic statistic);
593 
606  float dewpoint(float temperature, float relative_humidity);
607 
613  vec dewpoint(const vec& temperature, const vec& relative_humidity);
614 
622  float pressure(float ielev, float oelev, float ipressure, float itemperature=288.15);
623 
631  vec pressure(const vec& ielev, const vec& oelev, const vec& ipressure, const vec& itemperature);
632 
638  float qnh(float pressure, float altitude);
639 
645  vec qnh(const vec& pressure, const vec& altitude);
646 
652  float relative_humidity(float temperature, float dewpoint);
653 
659  vec relative_humidity(const vec& temperature, const vec& dewpoint);
660 
667  float wetbulb(float temperature, float pressure, float relative_humidity);
668 
675  vec wetbulb(const vec& temperature, const vec& pressure, const vec& relative_humidity);
676 
682  float wind_speed(float xwind, float ywind);
683 
689  vec wind_speed(const vec& xwind, const vec& ywind);
690 
697  float wind_direction(float xwind, float ywind);
698 
704  vec wind_direction(const vec& xwind, const vec& ywind);
705 
713  void set_omp_threads(int num);
714 
716  void initialize_omp();
717 
722  // vec2 calc_gradient(const vec2& values, const vec2& aux, int radius);
723  // ivec regression(const vec& x, const vec& y);
724 
726  Statistic get_statistic(std::string name);
727 
731  std::string version();
732 
733  double clock();
734  void debug(std::string string);
735  void warning(std::string string);
736  void error(std::string string);
737  void future_deprecation_warning(std::string function, std::string other="");
738  bool is_valid(float value);
739  float calc_statistic(const vec& array, Statistic statistic);
740  float calc_quantile(const vec& array, float quantile);
741  vec calc_statistic(const vec2& array, Statistic statistic);
742  vec calc_quantile(const vec2& array, float quantile=MV);
743 
749  vec2 calc_quantile(const vec3& array, const vec2& quantile);
750  int num_missing_values(const vec2& iArray);
751 
757  int get_lower_index(float iX, const std::vector<float>& iValues);
758 
764  int get_upper_index(float iX, const std::vector<float>& iValues);
765 
777  float interpolate(float x, const std::vector<float>& iX, const std::vector<float>& iY);
778 
780  ivec2 init_ivec2(int Y, int X, int value);
781  vec2 init_vec2(int Y, int X, float value=MV);
782 
784  ivec3 init_ivec3(int Y, int X, int E, int value);
785  vec3 init_vec3(int Y, int X, int E, float value=MV);
786 
793  vec calc_even_quantiles(const vec& values, int num);
794 
804  vec2 calc_gradient(const Grid& grid, const vec2& base, const vec2& values, int radius, int min_num=2, float min_range=0, float default_gradient=0);
805 
807  bool compatible_size(const Grid& grid, const vec2& v);
808  bool compatible_size(const Grid& grid, const vec3& v);
809  bool compatible_size(const Points& points, const vec& v);
810  bool compatible_size(const Points& points, const vec2& v);
811  bool compatible_size(const vec2& a, const vec3& b);
812  bool compatible_size(const vec3& a, const vec3& b);
813 
823  bool point_in_rectangle(const Point& A, const Point& B, const Point& C, const Point& D, const Point& m );
824 
830  float* test_array(float* v, int n);
832  float test_vec_input(const vec& input);
834  int test_ivec_input(const ivec& input);
836  float test_vec2_input(const vec2& input);
838  float test_vec3_input(const vec3& input);
840  vec test_vec_output();
841  ivec test_ivec_output();
843  vec2 test_vec2_output();
844  ivec2 test_ivec2_output();
846  vec3 test_vec3_output();
847  ivec3 test_ivec3_output();
848 
850  float test_vec_argout(vec& distances);
852  float test_vec2_argout(vec2& distances);
853 
855 
857  static const float swig_default_value = -1;
858 
862  class Point {
863  public:
871  Point(float lat, float lon, float elev=MV, float laf=MV, CoordinateType type=Geodetic);
872  float lat;
873  float lon;
874  float elev;
875  float laf;
877  };
880  public:
881  StructureFunction(float localization_distance);
883  virtual float corr(const Point& p1, const Point& p2) const = 0;
884  virtual float corr_background(const Point& p1, const Point& p2) const;
888  float localization_distance() const;
889  virtual StructureFunction* clone() const = 0;
890  protected:
895  float barnes_rho(float dist, float length) const;
896 
901  float cressman_rho(float dist, float length) const;
903  };
905  public:
912  MultipleStructure(const StructureFunction& structure_h, const StructureFunction& structure_v, const StructureFunction& structure_w);
913  float corr(const Point& p1, const Point& p2) const;
914  StructureFunction* clone() const;
915  private:
916  StructureFunction* m_structure_h;
917  StructureFunction* m_structure_v;
918  StructureFunction* m_structure_w;
919  };
922  public:
929  BarnesStructure(float h, float v=0, float w=0, float hmax=MV);
930  float corr(const Point& p1, const Point& p2) const;
931  StructureFunction* clone() const;
932  private:
933  float mH;
934  float mV;
935  float mW;
936  };
937 
940  public:
941  CressmanStructure(float h, float v=0, float w=0);
942  float corr(const Point& p1, const Point& p2) const;
943  StructureFunction* clone() const;
944  private:
945  float mH;
946  float mV;
947  float mW;
948  };
950  public:
955  CrossValidation(StructureFunction& structure, float dist=MV);
956  float corr(const Point& p1, const Point& p2) const;
957  float corr_background(const Point& p1, const Point& p2) const;
958  StructureFunction* clone() const;
959  private:
960  StructureFunction* m_structure;
961  float m_dist;
962  };
963 
964  class Transform {
965  public:
966  // Note these cannot be pure virtual, otherwise SWIG does not expose
967  // the vector functions (with the same name) in python. Therefore, make sure these
968  // functions are overloaded in the subclass implementation
969  virtual float forward(float value) const;
970  virtual float backward(float value) const;
971 
972  vec forward(const vec& input) const;
973  vec backward(const vec& input) const;
974  vec2 forward(const vec2& input) const;
975  vec2 backward(const vec2& input) const;
976  vec3 forward(const vec3& input) const;
977  vec3 backward(const vec3& input) const;
978  };
979  class Identity : public Transform {
980  public:
981  // SWIG requires these "using" statements to enable the vectorized versions in the
982  // subclasses
983  using Transform::forward;
984  using Transform::backward;
985  float forward(float value) const;
986  float backward(float value) const;
987  };
988  class Log : public Transform {
989  public:
990  using Transform::forward;
991  using Transform::backward;
992  float forward(float value) const;
993  float backward(float value) const;
994  };
995  class BoxCox : public Transform {
996  public:
997  BoxCox(float threshold);
998  using Transform::forward;
999  using Transform::backward;
1000  float forward(float value) const;
1001  float backward(float value) const;
1002  private:
1003  float mThreshold;
1004  };
1005 
1007  class KDTree {
1008  public:
1009  KDTree(vec lats, vec lons, CoordinateType type=Geodetic);
1010  KDTree& operator=(KDTree other);
1011  KDTree(const KDTree& other);
1012  KDTree() {};
1013 
1018  int get_nearest_neighbour(float lat, float lon, bool include_match=true) const;
1019 
1025  ivec get_neighbours(float lat, float lon, float radius, bool include_match=true) const;
1026 
1033  ivec get_neighbours_with_distance(float lat, float lon, float radius, vec& distances, bool include_match=true) const;
1034 
1040  int get_num_neighbours(float lat, float lon, float radius, bool include_match=true) const;
1041 
1047  ivec get_closest_neighbours(float lat, float lon, int num, bool include_match=true) const;
1048 
1049 
1057  bool convert_coordinates(const vec& lats, const vec& lons, vec& x_coords, vec& y_coords, vec& z_coords) const;
1058 
1066  bool convert_coordinates(float lat, float lon, float& x_coord, float& y_coord, float& z_coord) const;
1067  static float deg2rad(float deg);
1068  static float rad2deg(float deg);
1069  static float calc_distance(float lat1, float lon1, float lat2, float lon2, CoordinateType type=Geodetic);
1070  static float calc_distance(float x0, float y0, float z0, float x1, float y1, float z1);
1071  static float calc_distance_fast(float lat1, float lon1, float lat2, float lon2, CoordinateType type=Geodetic);
1072  static float calc_distance_fast(const Point& p1, const Point& p2);
1073  vec get_lats() const;
1074  vec get_lons() const;
1075  int size() const;
1076  CoordinateType get_coordinate_type() const;
1077  protected:
1078  typedef boost::geometry::model::point<float, 3, boost::geometry::cs::cartesian> point;
1079  typedef std::pair<point, unsigned> value;
1080  typedef boost::geometry::model::box<point> box;
1081  boost::geometry::index::rtree< value, boost::geometry::index::quadratic<16> > mTree;
1082  vec mLats;
1083  vec mLons;
1085 
1086  struct within_radius {
1087  public:
1088  within_radius(point p, float radius, bool include_match);
1089  bool operator()(value const& v) const;
1090  private:
1091  float radius;
1092  point p;
1093  bool include_match;
1094  };
1095  struct is_not_equal {
1096  public:
1097  is_not_equal(point p);
1098  bool operator()(value const& v) const;
1099  private:
1100  point p;
1101  };
1102  };
1103 
1105  class Points {
1106  public:
1107  Points();
1115  Points(vec lats, vec lons, vec elevs=vec(), vec lafs=vec(), CoordinateType type=Geodetic);
1116  Points(KDTree tree, vec elevs=vec(), vec lafs=vec());
1117  Points& operator=(Points other);
1118  Points(const Points& other);
1119  // Returns -1 if there are no neighbours
1120  int get_nearest_neighbour(float lat, float lon, bool include_match=true) const;
1121  ivec get_neighbours(float lat, float lon, float radius, bool include_match=true) const;
1122  ivec get_neighbours_with_distance(float lat, float lon, float radius, vec& distances, bool include_match=true) const;
1123  int get_num_neighbours(float lat, float lon, float radius, bool include_match=true) const;
1124  ivec get_closest_neighbours(float lat, float lon, int num, bool include_match=true) const;
1125 
1126  vec get_lats() const;
1127  vec get_lons() const;
1128  vec get_elevs() const;
1129  vec get_lafs() const;
1130  int size() const;
1131  ivec get_in_domain_indices(const Grid& grid) const;
1132  Points get_in_domain(const Grid& grid) const;
1133  CoordinateType get_coordinate_type() const;
1134  Point get_point(int index) const;
1135  Points subset(const ivec& indices) const;
1136  private:
1137  KDTree mTree;
1138  vec mLats;
1139  vec mLons;
1140  vec mElevs;
1141  vec mLafs;
1142  };
1143 
1145  class Grid {
1146  public:
1147  Grid();
1148 
1156  Grid(vec2 lats, vec2 lons, vec2 elevs=vec2(), vec2 lafs=vec2(), CoordinateType type=Geodetic);
1157  ivec get_nearest_neighbour(float lat, float lon, bool include_match=true) const;
1158  ivec2 get_neighbours(float lat, float lon, float radius, bool include_match=true) const;
1159  ivec2 get_neighbours_with_distance(float lat, float lon, float radius, vec& distances, bool include_match=true) const;
1160  int get_num_neighbours(float lat, float lon, float radius, bool include_match=true) const;
1161  ivec2 get_closest_neighbours(float lat, float lon, int num, bool include_match=true) const;
1162 
1163  bool get_box(float lat, float lon, int& Y1_out, int& X1_out, int& Y2_out, int& X2_out) const;
1164 
1166  Points to_points() const;
1167 
1168  vec2 get_lats() const;
1169  vec2 get_lons() const;
1170  vec2 get_elevs() const;
1171  vec2 get_lafs() const;
1172  ivec size() const;
1173  CoordinateType get_coordinate_type() const;
1174  Point get_point(int y_index, int x_index) const;
1175  private:
1176  KDTree mTree;
1177  int mX;
1178  vec2 get_2d(vec input) const;
1179  ivec get_indices(int index) const;
1180  ivec2 get_indices(ivec indices) const;
1181  vec2 mLats;
1182  vec2 mLons;
1183  vec2 mElevs;
1184  vec2 mLafs;
1185  };
1186  class not_implemented_exception: public std::logic_error
1187  {
1188  public:
1189  not_implemented_exception() : std::logic_error("Function not yet implemented") { };
1190  };
1191 };
1192 #endif
Definition: gridpp.h:1186
float test_vec2_input(const vec2 &input)
Testing function for 2D input vector.
vec subset(const vec &array, const ivec &indices)
Definition: util.cpp:172
vec2 nearest(const Grid &igrid, const Grid &ogrid, const vec2 &ivalues)
Nearest neighbour dowscaling grid to grid.
std::string version()
The gridpp version.
Definition: gridpp.h:1086
Bias.
Definition: gridpp.h:84
vec2 neighbourhood_quantile_ens_fast(const vec3 &input, float quantile, int radius, const vec &thresholds)
Deprecated: Compute neighbourhood quantiles fast on ensemble field.
float laf
Definition: gridpp.h:875
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.
static const float MV_CML
Missing value indicator in gridpp command-line tool.
Definition: gridpp.h:41
Helper class for Grid and Points.
Definition: gridpp.h:1007
void test_not_implemented_exception()
vec2 bilinear(const Grid &igrid, const Grid &ogrid, const vec2 &ivalues)
Bilinear downscaling grid to grid.
Minimum of values.
Definition: gridpp.h:68
float dewpoint(float temperature, float relative_humidity)
Calculate dewpoint temperature from temperature and relative humidity.
Statistic get_statistic(std::string name)
Convert name of a statistic enum.
Maximum of values.
Definition: gridpp.h:70
float wind_speed(float xwind, float ywind)
Diagnose wind speed from its components.
Definition: gridpp.h:979
CorrectionType
Method for statistical correction.
Definition: gridpp.h:89
std::vector< int > ivec
Definition: gridpp.h:25
ivec2 test_ivec2_output()
ivec test_ivec_output()
vec2 test_vec2_output()
Testing function for 2D output vector.
virtual float forward(float value) const
Unknown statistic.
Definition: gridpp.h:75
Metric
Binary verification metrics.
Definition: gridpp.h:79
Simple structure function based on distance, elevation, and land area fraction.
Definition: gridpp.h:921
vec2 metric_optimizer_curve(const vec &ref, const vec &fcst, const vec &thresholds, Metric metric)
Create calibration curve that optimizes a metric.
Definition: gridpp.h:949
vec background(const vec &elevs, const vec &values, int num_min_prof, float min_elev_diff, float value_minp, float value_maxp, BackgroundType background_type, const vec &external_background_values, const ivec &indices_global_outer, bool debug)
Background (first guess) calculations at observation locations.
Definition: background.cpp:35
vec distance(const Grid &grid, const Points &points, int num=1)
For each point, calculates the distance to nearest gridpoint.
double clock()
Extrapolation
Methods for extrapolating outside a curve.
Definition: gridpp.h:57
float test_vec_argout(vec &distances)
Testing function for 1D vector treated as output.
CoordinateType type
Definition: gridpp.h:876
not_implemented_exception()
Definition: gridpp.h:1189
Proportion correct.
Definition: gridpp.h:83
bool convert_coordinates(const vec &lats, const vec &lons, vec &x_coords, vec &y_coords, vec &z_coords)
Convert lat/lons to 3D cartesian coordinates with the centre of the earth as the origin.
Definition: util.cpp:59
Population variance of values.
Definition: gridpp.h:73
Point(float lat, float lon, float elev=MV, float laf=MV, CoordinateType type=Geodetic)
Constructor.
vec2 smart(const Grid &igrid, const Grid &ogrid, const vec2 &ivalues, int num, const StructureFunction &structure)
Smart neighbour downscaling grid to grid.
void set_omp_threads(int num)
Set the number of OpenMP threads to use.
float relative_humidity(float temperature, float dewpoint)
Calculate relative humidity from temperature and dewpoint temperature.
ivec3 test_ivec3_output()
vec2 correction(const Grid &rgrid, const vec2 &rvalues, const Points &npoints, const vec &nvalues, float mean_radius, float outer_radius, float inner_radius, int min_num, int max_num, CorrectionType type, ivec2 &count)
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...
vec2 neighbourhood_quantile(const vec2 &input, float quantile, int halfwidth)
Computes a quantile in a sliding square neighbourhood.
boost::geometry::model::point< float, 3, boost::geometry::cs::cartesian > point
Definition: gridpp.h:1078
Threat score.
Definition: gridpp.h:81
vec2 doping(const Grid &igrid, const vec2 &input, const Points &points, const vec &values, int half_width, float max_elev_diff=gridpp::MV)
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.
Represents a vector of locations and their metadata.
Definition: gridpp.h:1105
vec2 neighbourhood_quantile_fast(const vec2 &input, float quantile, int halfwidth, const vec &thresholds)
Fast and approximate neighbourhood quantile.
std::pair< point, unsigned > value
Definition: gridpp.h:1079
std::vector< ivec2 > ivec3
Definition: gridpp.h:27
Definition: gridpp.h:964
Keep values the way they were.
Definition: gridpp.h:62
Definition: gridpp.h:1095
vec2 neighbourhood_brute_force(const vec2 &input, int halfwidth, Statistic statistic)
Spatial neighbourhood filter without any shortcuts.
ivec3 init_ivec3(int Y, int X, int E, int value)
Initialize a vector of size Y, X, E, with a given value.
std::vector< vec > vec2
Definition: gridpp.h:23
Latitude and longitude.
Definition: gridpp.h:97
Represents a 2D grid of locations and their metadata.
Definition: gridpp.h:1145
Simple structure function based on distance, elevation, and land area fraction.
Definition: gridpp.h:939
vec2 neighbourhood(const vec2 &input, int halfwidth, Statistic statistic)
Spatial neighbourhood filter, computing a statistic for a sliding square window.
Covariance structure function.
Definition: gridpp.h:879
Heidke skill score.
Definition: gridpp.h:85
Continue past the end-points using the mean slope of the curve.
Definition: gridpp.h:59
void debug(std::string string)
void initialize_omp()
Sets the number of OpenMP threads to 1 if OMP_NUM_THREADS undefined.
float get_optimal_threshold(const vec &ref, const vec &fcst, float threshold, Metric metric)
virtual float backward(float value) const
static const float swig_default_value
Default value used to fill array in SWIG testing functions.
Definition: gridpp.h:857
vec count(const Grid &grid, const Points &points, float radius)
For each point, counts the number of gridpoints within the radius.
int get_upper_index(float iX, const std::vector< float > &iValues)
Find the index in a vector that is equal or just above a value.
Continue past the end-points using the slope of the two lowermost or uppermost points in the curve...
Definition: gridpp.h:60
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.
boost::geometry::model::box< point > box
Definition: gridpp.h:1080
float test_vec2_argout(vec2 &distances)
Testing function for 2D vector treated as output.
float test_vec_input(const vec &input)
Testing function for 1D input vector.
CoordinateType
Types of coordinates for position of points.
Definition: gridpp.h:96
vec test_vec_output()
Testing function for 1D output vector.
float interpolate(float x, const std::vector< float > &iX, const std::vector< float > &iY)
Piecewise linear interpolation If x is outside the range of iX, then the min/max value of iY is used...
vec quantile_mapping_curve(const vec &ref, const vec &fcst, vec &output_fcst, vec quantiles=vec())
Create quantile mapping calibration curve.
vec3 init_vec3(int Y, int X, int E, float value=MV)
vec mLats
Definition: gridpp.h:1082
Definition: gridpp.h:904
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.
float wetbulb(float temperature, float pressure, float relative_humidity)
Calculate wetbulb temperature from temperature, pressure, and relative humidity.
std::vector< double > dvec
Definition: gridpp.h:30
int num_missing_values(const vec2 &iArray)
float lon
Definition: gridpp.h:873
std::vector< dvec > dvec2
Definition: gridpp.h:31
vec2 gridding(const Grid &grid, const Points &points, const vec &values, float radius, int min_num, Statistic statistic)
Aggregate points onto a grid.
Mean of values.
Definition: gridpp.h:69
Statistic
Statistical operations to reduce a vector to a scalar.
Definition: gridpp.h:66
int test_ivec_input(const ivec &input)
Testing function for 1D input vector.
vec mLons
Definition: gridpp.h:1083
bool is_valid(float value)
CoordinateType mType
Definition: gridpp.h:1084
Definition: gridpp.h:17
vec2 init_vec2(int Y, int X, float value=MV)
vec calc_even_quantiles(const vec &values, int num)
Get reasonably spaced quantiles from a vector of values, ignoring duplicate values but including the ...
vec2 neighbourhood_ens(const vec3 &input, int halfwidth, Statistic statistic)
Deprecated: Compute neighbourhood statistic on ensemble field.
A quantile from values.
Definition: gridpp.h:71
float elev
Definition: gridpp.h:874
vec2 calc_gradient(const Grid &grid, const vec2 &base, const vec2 &values, int radius, int min_num=2, float min_range=0, float default_gradient=0)
Computes gradients based on values in neighbourhood.
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...
static const float pi
Mathematical constant pi.
Definition: gridpp.h:43
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.
int get_lower_index(float iX, const std::vector< float > &iValues)
Find the index in a vector that is equal or just below a value.
KDTree()
Definition: gridpp.h:1012
float calc_quantile(const vec &array, float quantile)
Multiplicative.
Definition: gridpp.h:91
ivec2 init_ivec2(int Y, int X, int value)
Initialize a vector of size Y, X, with a given value.
float qnh(float pressure, float altitude)
Diagnose QNH from pressure and altitude.
Standard deviation of values.
Definition: gridpp.h:72
static const float MV
Missing value indicator.
Definition: gridpp.h:39
vec get_neighbourhood_thresholds(const vec2 &input, int num_thresholds)
Calculate appropriate approximation thresholds for neighbourhood quantile.
float pressure(float ielev, float oelev, float ipressure, float itemperature=288.15)
Calculate pressure at a new elevation.
float calc_score(float a, float b, float c, float d, Metric metric)
float wind_direction(float xwind, float ywind)
Diagnose wind direction from its components.
X and Y.
Definition: gridpp.h:98
vec monotonize_curve(vec curve_ref, vec curve_fcst, vec &output_fcst)
Ensure calibration curve is monotonic, by removing points.
vec3 test_vec3_output()
Testing function for 3D output vector.
float calc_statistic(const vec &array, Statistic statistic)
float test_vec3_input(const vec3 &input)
Testing function for 3D input vector.
void future_deprecation_warning(std::string function, std::string other="")
boost::geometry::index::rtree< value, boost::geometry::index::quadratic< 16 > > mTree
Definition: gridpp.h:1081
vec2 neighbourhood_quantile_ens(const vec3 &input, float quantile, int halfwidth)
Deprecated: Compute neighbourhood quantiles on ensemble field.
Represents a single point in some coordinate system.
Definition: gridpp.h:862
Mean of values.
Definition: gridpp.h:67
void warning(std::string string)
float mLocalizationDistance
Definition: gridpp.h:902
Additive.
Definition: gridpp.h:92
vec2 fill_missing(const vec2 &values)
Fill in missing values based on nearby values.
float * test_array(float *v, int n)
Special function whose presense is needed for SWIG.
float deg2rad(float deg)
Definition: util.cpp:108
Continue past the end-points using a slope of 0.
Definition: gridpp.h:61
Continue past the end-points using a slope of 1.
Definition: gridpp.h:58
vec2 simple_gradient(const Grid &igrid, const Grid &ogrid, const vec2 &ivalues, float elev_gradient)
std::vector< float > vec
Definition: gridpp.h:22
bool compatible_size(const Grid &grid, const vec2 &v)
Check if the grid is the same size as the 2D vector.
Hannsen-Kuiper skill score.
Definition: gridpp.h:82
float lat
Definition: gridpp.h:872
Equitable threat score.
Definition: gridpp.h:80
Sum of values.
Definition: gridpp.h:74
static const double radius_earth
Radius of the earth [m].
Definition: gridpp.h:45
Definition: gridpp.h:988
Definition: gridpp.h:995
void error(std::string string)
Quantile mapping.
Definition: gridpp.h:90
float calc_distance(float lat1, float lon1, float lat2, float lon2)
Definition: util.cpp:83
std::vector< ivec > ivec2
Definition: gridpp.h:26
std::vector< vec2 > vec3
Definition: gridpp.h:24