Gridpp 0.5.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 #ifdef _OPENMP
10  #include <omp.h>
11 #endif
12 #include <exception>
13 
14 #define GRIDPP_VERSION "0.5.0"
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  };
63 
65  enum Statistic {
66  Mean = 0,
67  Min = 10,
68  Median = 20,
69  Max = 30,
70  Quantile = 40,
71  Std = 50,
72  Variance = 60,
73  Sum = 70,
74  Unknown = -1
75  };
76 
78  enum Metric {
79  Ets = 0,
80  Ts = 1,
81  Kss = 20,
82  Pc = 30,
83  Bias = 40,
84  Hss = 50,
85  };
86 
89  Qq = 0,
91  Additive = 20,
92  };
93 
96  Geodetic = 0,
97  Cartesian = 1,
98  };
99 
115  vec2 optimal_interpolation(const Grid& bgrid,
116  const vec2& background,
117  const Points& points,
118  const vec& pobs,
119  const vec& pratios,
120  const vec& pbackground,
121  const StructureFunction& structure,
122  int max_points);
123 
134  vec optimal_interpolation(const Points& bpoints,
135  const vec& background,
136  const Points& points,
137  const vec& pobs,
138  const vec& pratios,
139  const vec& pbackground,
140  const StructureFunction& structure,
141  int max_points);
142 
145  vec2 optimal_interpolation_transform(const Grid& bgrid,
146  const vec2& background,
147  float bsigma,
148  const Points& points,
149  const vec& pobs,
150  const vec& psigmas,
151  const vec& pbackground,
152  const StructureFunction& structure,
153  int max_points,
154  const Transform& transform);
155 
156  vec optimal_interpolation_transform(const Points& bpoints,
157  const vec& background,
158  float bsigma,
159  const Points& points,
160  const vec& pobs,
161  const vec& psigmas,
162  const vec& pbackground,
163  const StructureFunction& structure,
164  int max_points,
165  const Transform& transform);
166 
175  vec3 optimal_interpolation_ensi(const Grid& bgrid,
176  const vec3& background,
177  const Points& points,
178  const vec& pobs,
179  const vec& psigmas,
180  const vec2& pbackground,
181  const StructureFunction& structure,
182  int max_points);
183 
184  vec2 optimal_interpolation_ensi(const Points& bpoints,
185  const vec2& background,
186  const Points& points,
187  const vec& pobs,
188  const vec& psigmas,
189  const vec2& pbackground,
190  const StructureFunction& structure,
191  int max_points);
192 
199  vec2 fill(const Grid& igrid, const vec2& input, const Points& points, const vec& radii, float value, bool outside);
212  vec2 neighbourhood(const vec2& input, int halfwidth, Statistic statistic);
213 
219  vec2 neighbourhood(const vec3& input, int halfwidth, Statistic statistic);
220 
226  vec2 neighbourhood_quantile(const vec2& input, float quantile, int halfwidth);
227 
233  vec2 neighbourhood_quantile(const vec3& input, float quantile, int halfwidth);
234 
241  vec2 neighbourhood_quantile_fast(const vec2& input, float quantile, int halfwidth, const vec& thresholds);
242 
249  vec2 neighbourhood_quantile_fast(const vec3& input, float quantile, int halfwidth, const vec& thresholds);
250 
257  vec2 neighbourhood_quantile_fast(const vec2& input, const vec2& quantile, int halfwidth, const vec& thresholds);
258 
265  vec2 neighbourhood_quantile_fast(const vec3& input, const vec2& quantile, int halfwidth, const vec& thresholds);
266 
272  vec2 neighbourhood_brute_force(const vec2& input, int halfwidth, Statistic statistic);
273 
279  vec2 neighbourhood_brute_force(const vec3& input, int halfwidth, Statistic statistic);
280 
285  vec get_neighbourhood_thresholds(const vec2& input, int num_thresholds);
286 
291  vec get_neighbourhood_thresholds(const vec3& input, int num_thresholds);
292 
295  vec2 neighbourhood_ens(const vec3& input, int halfwidth, Statistic statistic);
298  vec2 neighbourhood_quantile_ens(const vec3& input, float quantile, int halfwidth);
301  vec2 neighbourhood_quantile_ens_fast(const vec3& input, float quantile, int radius, const vec& thresholds);
315  vec2 quantile_mapping_curve(const vec& ref, const vec& fcst, vec quantiles=vec());
316 
325  vec2 metric_optimizer_curve(const vec& ref, const vec& fcst, const vec& thresholds, Metric metric);
334  vec apply_curve(const vec& fcst, const vec2& curve, Extrapolation policy_below, Extrapolation policy_above);
335 
343  vec2 apply_curve(const vec2& fcst, const vec2& curve, Extrapolation policy_below, Extrapolation policy_above);
344 
349  vec2 monotonize_curve(vec2 curve);
350 
351  float get_optimal_threshold(const vec& ref, const vec& fcst, float threshold, Metric metric);
352 
353  float calc_score(float a, float b, float c, float d, Metric metric);
354  float calc_score(const vec& ref, const vec& fcst, float threshold, Metric metric);
355  float calc_score(const vec& ref, const vec& fcst, float threshold, float fthreshold, Metric metric);
356 
357  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);
358  // Apply correction based on multiple timesteps
359  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);
360 
374  vec2 nearest(const Grid& igrid, const Grid& ogrid, const vec2 ivalues);
375  vec3 nearest(const Grid& igrid, const Grid& ogrid, const vec3 ivalues);
376 
383  vec nearest(const Grid& igrid, const Points& opoints, const vec2 ivalues);
384  vec2 nearest(const Grid& igrid, const Points& opoints, const vec3 ivalues);
385 
392  vec2 bilinear(const Grid& igrid, const Grid& ogrid, const vec2 ivalues);
393 
400  vec bilinear(const Grid& igrid, const Points& opoints, const vec2 ivalues);
401 
402  vec2 simple_gradient(const Grid& igrid, const Grid& ogrid, const vec2 ivalues, float elev_gradient);
403  vec simple_gradient(const Grid& igrid, const Points& opoints, const vec2 ivalues, float elev_gradient);
404 
413  vec2 smart(const Grid& igrid, const Grid& ogrid, const vec2& ivalues, int num, const StructureFunction& structure);
427  vec count(const Grid& grid, const Points& points, float radius);
428 
435  vec2 count(const Points& points, const Grid& grid, float radius);
436 
443  vec distance(const Grid& grid, const Points& points, int num=1);
444 
451  vec2 distance(const Grid& igrid, const Grid& ogrid, int num=1);
452 
459  vec2 distance(const Points& points, const Grid& grid, int num=1);
460 
465  vec2 fill_missing(const vec2& values);
466 
475  vec2 gridding(const Grid& grid, const Points& points, const vec& values, float radius, int min_num, Statistic statistic);
476 
489  float dewpoint(float temperature, float relative_humidity);
490 
496  vec dewpoint(const vec& temperature, const vec& relative_humidity);
497 
505  float pressure(float ielev, float oelev, float ipressure, float itemperature=288.15);
506 
514  vec pressure(const vec& ielev, const vec& oelev, const vec& ipressure, const vec& itemperature);
515 
521  float qnh(float pressure, float altitude);
522 
528  vec qnh(const vec& pressure, const vec& altitude);
529 
535  float relative_humidity(float temperature, float dewpoint);
536 
542  vec relative_humidity(const vec& temperature, const vec& dewpoint);
543 
550  float wetbulb(float temperature, float pressure, float relative_humidity);
551 
558  vec wetbulb(const vec& temperature, const vec& pressure, const vec& relative_humidity);
559 
565  float wind_speed(float xwind, float ywind);
566 
572  vec wind_speed(const vec& xwind, const vec& ywind);
573 
580  float wind_direction(float xwind, float ywind);
581 
587  vec wind_direction(const vec& xwind, const vec& ywind);
588 
596  void set_omp_threads(int num);
597 
599  void initialize_omp();
600 
605  // vec2 calc_gradient(const vec2& values, const vec2& aux, int radius);
606  // ivec regression(const vec& x, const vec& y);
607 
609  Statistic get_statistic(std::string name);
610 
614  std::string version();
615 
616  double clock();
617  void debug(std::string string);
618  void warning(std::string string);
619  void error(std::string string);
620  void future_deprecation_warning(std::string function, std::string other="");
621  bool is_valid(float value);
622  float calc_statistic(const vec& array, Statistic statistic);
623  float calc_quantile(const vec& array, float quantile);
624  vec calc_statistic(const vec2& array, Statistic statistic);
625  vec calc_quantile(const vec2& array, float quantile=MV);
626  int num_missing_values(const vec2& iArray);
627 
633  int get_lower_index(float iX, const std::vector<float>& iValues);
634 
640  int get_upper_index(float iX, const std::vector<float>& iValues);
641 
649  float interpolate(float x, const std::vector<float>& iX, const std::vector<float>& iY);
650 
652  ivec2 init_ivec2(int Y, int X, int value);
653  vec2 init_vec2(int Y, int X, float value=MV);
654 
656  ivec3 init_ivec3(int Y, int X, int E, int value);
657  vec3 init_vec3(int Y, int X, int E, float value=MV);
658 
665  vec calc_even_quantiles(const vec& values, int num);
666 
676  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);
677 
679  bool compatible_size(const Grid& grid, const vec2& v);
680  bool compatible_size(const Grid& grid, const vec3& v);
681 
691  bool point_in_rectangle(const Point& A, const Point& B, const Point& C, const Point& D, const Point& m );
692 
698  float* test_array(float* v, int n);
700  float test_vec_input(const vec& input);
702  int test_ivec_input(const ivec& input);
704  float test_vec2_input(const vec2& input);
706  float test_vec3_input(const vec3& input);
708  vec test_vec_output();
709  ivec test_ivec_output();
711  vec2 test_vec2_output();
712  ivec2 test_ivec2_output();
714  vec3 test_vec3_output();
715  ivec3 test_ivec3_output();
716 
718  float test_vec_argout(vec& distances);
720  float test_vec2_argout(vec2& distances);
721 
723 
725  static const float swig_default_value = -1;
726 
730  class Point {
731  public:
739  Point(float lat, float lon, float elev=MV, float laf=MV, CoordinateType type=Geodetic);
740  float lat;
741  float lon;
742  float elev;
743  float laf;
745  };
748  public:
749  StructureFunction(float localization_distance);
751  virtual float corr(const Point& p1, const Point& p2) const = 0;
752  virtual float corr_background(const Point& p1, const Point& p2) const;
756  float localization_distance() const;
757  virtual StructureFunction* clone() const = 0;
758  protected:
763  float barnes_rho(float dist, float length) const;
764 
769  float cressman_rho(float dist, float length) const;
771  };
774  public:
781  BarnesStructure(float h, float v=0, float w=0, float hmax=MV);
782  float corr(const Point& p1, const Point& p2) const;
783  StructureFunction* clone() const;
784  private:
785  float mH;
786  float mV;
787  float mW;
788  };
789 
792  public:
793  CressmanStructure(float h, float v=0, float w=0);
794  float corr(const Point& p1, const Point& p2) const;
795  StructureFunction* clone() const;
796  private:
797  float mH;
798  float mV;
799  float mW;
800  };
802  public:
807  CrossValidation(StructureFunction& structure, float dist=MV);
808  float corr(const Point& p1, const Point& p2) const;
809  float corr_background(const Point& p1, const Point& p2) const;
810  StructureFunction* clone() const;
811  private:
812  StructureFunction* m_structure;
813  float m_dist;
814  };
815 
816  class Transform {
817  public:
818  // Note these cannot be pure virtual, otherwise SWIG does not expose
819  // the vector functions (with the same name) in python. Therefore, make sure these
820  // functions are overloaded in the subclass implementation
821  virtual float forward(float value) const;
822  virtual float backward(float value) const;
823 
824  vec forward(const vec& input) const;
825  vec backward(const vec& input) const;
826  vec2 forward(const vec2& input) const;
827  vec2 backward(const vec2& input) const;
828  vec3 forward(const vec3& input) const;
829  vec3 backward(const vec3& input) const;
830  };
831  class Identity : public Transform {
832  public:
833  // SWIG requires these "using" statements to enable the vectorized versions in the
834  // subclasses
835  using Transform::forward;
836  using Transform::backward;
837  float forward(float value) const;
838  float backward(float value) const;
839  };
840  class Log : public Transform {
841  public:
842  using Transform::forward;
843  using Transform::backward;
844  float forward(float value) const;
845  float backward(float value) const;
846  };
847  class BoxCox : public Transform {
848  public:
849  BoxCox(float threshold);
850  using Transform::forward;
851  using Transform::backward;
852  float forward(float value) const;
853  float backward(float value) const;
854  private:
855  float mThreshold;
856  };
857 
859  class KDTree {
860  public:
861  KDTree(vec lats, vec lons, CoordinateType type=Geodetic);
862  KDTree& operator=(KDTree other);
863  KDTree(const KDTree& other);
864  KDTree() {};
865 
870  int get_nearest_neighbour(float lat, float lon) const;
871 
877  ivec get_neighbours(float lat, float lon, float radius) const;
878 
885  ivec get_neighbours_with_distance(float lat, float lon, float radius, vec& distances) const;
886 
892  int get_num_neighbours(float lat, float lon, float radius) const;
893 
899  ivec get_closest_neighbours(float lat, float lon, int num) const;
900 
901 
909  bool convert_coordinates(const vec& lats, const vec& lons, vec& x_coords, vec& y_coords, vec& z_coords) const;
910 
918  bool convert_coordinates(float lat, float lon, float& x_coord, float& y_coord, float& z_coord) const;
919  static float deg2rad(float deg);
920  static float rad2deg(float deg);
921  static float calc_distance(float lat1, float lon1, float lat2, float lon2, CoordinateType type=Geodetic);
922  static float calc_distance(float x0, float y0, float z0, float x1, float y1, float z1);
923  static float calc_distance_fast(float lat1, float lon1, float lat2, float lon2, CoordinateType type=Geodetic);
924  static float calc_distance_fast(const Point& p1, const Point& p2);
925  vec get_lats() const;
926  vec get_lons() const;
927  int size() const;
928  CoordinateType get_coordinate_type() const;
929  protected:
930  typedef boost::geometry::model::point<float, 3, boost::geometry::cs::cartesian> point;
931  typedef std::pair<point, unsigned> value;
932  typedef boost::geometry::model::box<point> box;
933  boost::geometry::index::rtree< value, boost::geometry::index::quadratic<16> > mTree;
934  vec mLats;
935  vec mLons;
937  };
938 
940  class Points {
941  public:
942  Points();
950  Points(vec lats, vec lons, vec elevs=vec(), vec lafs=vec(), CoordinateType type=Geodetic);
951  Points(KDTree tree, vec elevs=vec(), vec lafs=vec());
952  Points& operator=(Points other);
953  Points(const Points& other);
954  // Returns -1 if there are no neighbours
955  int get_nearest_neighbour(float lat, float lon) const;
956  ivec get_neighbours(float lat, float lon, float radius) const;
957  ivec get_neighbours_with_distance(float lat, float lon, float radius, vec& distances) const;
958  int get_num_neighbours(float lat, float lon, float radius) const;
959  ivec get_closest_neighbours(float lat, float lon, int num) const;
960 
961  vec get_lats() const;
962  vec get_lons() const;
963  vec get_elevs() const;
964  vec get_lafs() const;
965  int size() const;
966  ivec get_in_domain_indices(const Grid& grid) const;
967  Points get_in_domain(const Grid& grid) const;
968  CoordinateType get_coordinate_type() const;
969  private:
970  KDTree mTree;
971  vec mLats;
972  vec mLons;
973  vec mElevs;
974  vec mLafs;
975  };
976 
978  class Grid {
979  public:
980  Grid();
981 
989  Grid(vec2 lats, vec2 lons, vec2 elevs=vec2(), vec2 lafs=vec2(), CoordinateType type=Geodetic);
990  ivec get_nearest_neighbour(float lat, float lon) const;
991  ivec2 get_neighbours(float lat, float lon, float radius) const;
992  ivec2 get_neighbours_with_distance(float lat, float lon, float radius, vec& distances) const;
993  int get_num_neighbours(float lat, float lon, float radius) const;
994  ivec2 get_closest_neighbours(float lat, float lon, int num) const;
995 
996  bool get_box(float lat, float lon, int& Y1_out, int& X1_out, int& Y2_out, int& X2_out) const;
997 
999  Points to_points() const;
1000 
1001  vec2 get_lats() const;
1002  vec2 get_lons() const;
1003  vec2 get_elevs() const;
1004  vec2 get_lafs() const;
1005  ivec size() const;
1006  CoordinateType get_coordinate_type() const;
1007  private:
1008  KDTree mTree;
1009  int mX;
1010  vec2 get_2d(vec input) const;
1011  ivec get_indices(int index) const;
1012  ivec2 get_indices(ivec indices) const;
1013  vec2 mLats;
1014  vec2 mLons;
1015  vec2 mElevs;
1016  vec2 mLafs;
1017  };
1018  class not_implemented_exception: public std::logic_error
1019  {
1020  public:
1021  not_implemented_exception() : std::logic_error("Function not yet implemented") { };
1022  };
1023 };
1024 #endif
Definition: gridpp.h:1018
float test_vec2_input(const vec2 &input)
Testing function for 2D input vector.
Definition: swig.cpp:25
std::string version()
The gridpp version.
Definition: gridpp.cpp:8
Bias.
Definition: gridpp.h:83
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:530
float laf
Definition: gridpp.h:743
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:859
void test_not_implemented_exception()
Definition: swig.cpp:98
Minimum of values.
Definition: gridpp.h:67
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:69
float wind_speed(float xwind, float ywind)
Diagnose wind speed from its components.
Definition: wind.cpp:6
Definition: gridpp.h:831
CorrectionType
Method for statistical correction.
Definition: gridpp.h:88
std::vector< int > ivec
Definition: gridpp.h:25
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:74
Metric
Binary verification metrics.
Definition: gridpp.h:78
Simple structure function based on distance, elevation, and land area fraction.
Definition: gridpp.h:773
vec2 metric_optimizer_curve(const vec &ref, const vec &fcst, const vec &thresholds, Metric metric)
Create calibration curve that optimizes a metric.
Definition: metric_optimizer.cpp:105
virtual float forward(float value) const
Definition: transform.cpp:5
Definition: gridpp.h:801
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:178
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.
Definition: swig.cpp:86
CoordinateType type
Definition: gridpp.h:744
not_implemented_exception()
Definition: gridpp.h:1021
Proportion correct.
Definition: gridpp.h:82
Population variance of values.
Definition: gridpp.h:72
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:51
float relative_humidity(float temperature, float dewpoint)
Calculate relative humidity from temperature and dewpoint temperature.
Definition: humidity.cpp:32
ivec3 test_ivec3_output()
Definition: swig.cpp:76
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)
Definition: correction.cpp:251
vec2 nearest(const Grid &igrid, const Grid &ogrid, const vec2 ivalues)
Nearest neighbour dowscaling grid to grid.
Definition: nearest.cpp:7
Statistic get_statistic(std::string name)
Convert name of a statistic enum.
Definition: gridpp.cpp:11
vec2 monotonize_curve(vec2 curve)
Ensure calibration curve is monotonic, by removing points.
Definition: curve.cpp:82
vec2 neighbourhood_quantile(const vec2 &input, float quantile, int halfwidth)
Computes a quantile in a sliding square neighbourhood.
Definition: neighbourhood.cpp:515
boost::geometry::model::point< float, 3, boost::geometry::cs::cartesian > point
Definition: gridpp.h:930
Threat score.
Definition: gridpp.h:80
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.
Definition: fill.cpp:6
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)
Optimal interpolation using a structure function based on an ensemble See Lussana et al 2019 (DOI: 10...
Definition: oi_ensi.cpp:33
Represents a vector of locations and their metadata.
Definition: gridpp.h:940
vec2 neighbourhood_quantile_fast(const vec2 &input, float quantile, int halfwidth, const vec &thresholds)
Fast and approximate neighbourhood quantile.
Definition: neighbourhood.cpp:293
std::pair< point, unsigned > value
Definition: gridpp.h:931
std::vector< ivec2 > ivec3
Definition: gridpp.h:27
Definition: gridpp.h:816
vec2 simple_gradient(const Grid &igrid, const Grid &ogrid, const vec2 ivalues, float elev_gradient)
Definition: simple_gradient.cpp:27
vec2 bilinear(const Grid &igrid, const Grid &ogrid, const vec2 ivalues)
Bilinear downscaling grid to grid.
Definition: bilinear.cpp:18
vec2 neighbourhood_brute_force(const vec2 &input, int halfwidth, Statistic statistic)
Spatial neighbourhood filter without any shortcuts.
Definition: neighbourhood.cpp:509
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:328
std::vector< vec > vec2
Definition: gridpp.h:23
Latitude and longitude.
Definition: gridpp.h:96
Represents a 2D grid of locations and their metadata.
Definition: gridpp.h:978
Simple structure function based on distance, elevation, and land area fraction.
Definition: gridpp.h:791
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:747
Heidke skill score.
Definition: gridpp.h:84
Continue past the end-points using the mean slope of the curve.
Definition: gridpp.h:59
void debug(std::string string)
Definition: util.cpp:149
void initialize_omp()
Sets the number of OpenMP threads to 1 if OMP_NUM_THREADS undefined.
Definition: gridpp.cpp:39
float get_optimal_threshold(const vec &ref, const vec &fcst, float threshold, Metric metric)
Definition: metric_optimizer.cpp:129
static const float swig_default_value
Default value used to fill array in SWIG testing functions.
Definition: gridpp.h:725
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
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.
Definition: util.cpp:254
Continue past the end-points using the slope of the two lowermost or uppermost points in the curve...
Definition: gridpp.h:60
boost::geometry::model::box< point > box
Definition: gridpp.h:932
float test_vec2_argout(vec2 &distances)
Testing function for 2D vector treated as output.
Definition: swig.cpp:91
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:95
vec test_vec_output()
Testing function for 1D output vector.
Definition: swig.cpp:45
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...
Definition: util.cpp:273
vec3 init_vec3(int Y, int X, int E, float value=MV)
Definition: util.cpp:319
vec mLats
Definition: gridpp.h:934
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:403
float wetbulb(float temperature, float pressure, float relative_humidity)
Calculate wetbulb temperature from temperature, pressure, and relative humidity.
Definition: humidity.cpp:89
std::vector< double > dvec
Definition: gridpp.h:30
int num_missing_values(const vec2 &iArray)
Definition: util.cpp:140
float lon
Definition: gridpp.h:741
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.
Definition: gridding.cpp:6
Mean of values.
Definition: gridpp.h:68
Statistic
Statistical operations to reduce a vector to a scalar.
Definition: gridpp.h:65
vec apply_curve(const vec &fcst, const vec2 &curve, Extrapolation policy_below, Extrapolation policy_above)
Apply arbitrary calibration curve to 1D forecasts.
Definition: curve.cpp:6
vec2 quantile_mapping_curve(const vec &ref, const vec &fcst, vec quantiles=vec())
Create quantile mapping calibration curve.
Definition: quantile_mapping.cpp:5
int test_ivec_input(const ivec &input)
Testing function for 1D input vector.
Definition: swig.cpp:18
vec mLons
Definition: gridpp.h:935
bool is_valid(float value)
Definition: util.cpp:20
CoordinateType mType
Definition: gridpp.h:936
Definition: gridpp.h:17
vec2 init_vec2(int Y, int X, float value=MV)
Definition: util.cpp:307
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:185
vec2 neighbourhood_ens(const vec3 &input, int halfwidth, Statistic statistic)
Deprecated: Compute neighbourhood statistic on ensemble field.
Definition: neighbourhood.cpp:522
A quantile from values.
Definition: gridpp.h:70
float elev
Definition: gridpp.h:742
vec2 optimal_interpolation_transform(const Grid &bgrid, const vec2 &background, float bsigma, const Points &points, const vec &pobs, const vec &psigmas, const vec &pbackground, const StructureFunction &structure, int max_points, const Transform &transform)
This is its own function because the variance parameterization is different than in the non-transform...
Definition: oi.cpp:305
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.
Definition: calc_gradient.cpp:5
static const float pi
Mathematical constant pi.
Definition: gridpp.h:43
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.
Definition: util.cpp:235
KDTree()
Definition: gridpp.h:864
virtual float backward(float value) const
Definition: transform.cpp:8
float calc_quantile(const vec &array, float quantile)
Definition: util.cpp:91
Multiplicative.
Definition: gridpp.h:90
ivec2 init_ivec2(int Y, int X, int value)
Initialize a vector of size Y, X, with a given value.
Definition: util.cpp:313
float qnh(float pressure, float altitude)
Diagnose QNH from pressure and altitude.
Definition: qnh.cpp:6
Standard deviation of values.
Definition: gridpp.h:71
static const float MV
Missing value indicator.
Definition: gridpp.h:39
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)
Optimal interpolation for a deterministic gridded field.
Definition: oi.cpp:25
vec get_neighbourhood_thresholds(const vec2 &input, int num_thresholds)
Calculate appropriate approximation thresholds for neighbourhood quantile.
Definition: neighbourhood.cpp:240
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:200
float wind_direction(float xwind, float ywind)
Diagnose wind direction from its components.
Definition: wind.cpp:19
X and Y.
Definition: gridpp.h:97
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:170
boost::geometry::index::rtree< value, boost::geometry::index::quadratic< 16 > > mTree
Definition: gridpp.h:933
vec2 neighbourhood_quantile_ens(const vec3 &input, float quantile, int halfwidth)
Deprecated: Compute neighbourhood quantiles on ensemble field.
Definition: neighbourhood.cpp:526
Represents a single point in some coordinate system.
Definition: gridpp.h:730
Mean of values.
Definition: gridpp.h:66
void warning(std::string string)
Definition: util.cpp:153
float mLocalizationDistance
Definition: gridpp.h:770
Additive.
Definition: gridpp.h:91
vec2 fill_missing(const vec2 &values)
Fill in missing values based on nearby values.
Definition: fill.cpp:43
float * test_array(float *v, int n)
Special function whose presense is needed for SWIG.
Definition: swig.cpp:6
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
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.
Definition: util.cpp:301
Hannsen-Kuiper skill score.
Definition: gridpp.h:81
float lat
Definition: gridpp.h:740
Equitable threat score.
Definition: gridpp.h:79
Sum of values.
Definition: gridpp.h:73
static const double radius_earth
Radius of the earth [m].
Definition: gridpp.h:45
Definition: gridpp.h:840
Definition: gridpp.h:847
void error(std::string string)
Definition: util.cpp:157
Quantile mapping.
Definition: gridpp.h:89
std::vector< ivec > ivec2
Definition: gridpp.h:26
std::vector< vec2 > vec3
Definition: gridpp.h:24