Gridpp 0.7.0
A post-processing library for gridded weather forecasts
Enumerations | Classes
gridpp Namespace Reference

Typedefs

Short-hand notation for vectors of different dimensions sizes

typedef std::vector< float > vec
 1D float vector More...
 
typedef std::vector< vecvec2
 2D float vector More...
 
typedef std::vector< vec2vec3
 3D float vector More...
 
typedef std::vector< int > ivec
 1D integer vector More...
 
typedef std::vector< ivecivec2
 2D integer vector More...
 
typedef std::vector< ivec2ivec3
 3D integer vector More...
 
typedef std::vector< double > dvec
 1D double vector More...
 
typedef std::vector< dvecdvec2
 2D double vector More...
 

Functions

Data assimilation methods

Functions that merge observations with a background field

vec2 optimal_interpolation (const Grid &bgrid, const vec2 &background, const Points &obs_points, const vec &obs, const vec &variance_ratios, const vec &background_at_points, const StructureFunction &structure, int max_points, bool allow_extrapolation=true)
 Optimal interpolation for a deterministic gridded field. More...
 
vec optimal_interpolation (const Points &bpoints, const vec &background, const Points &obs_points, const vec &obs, const vec &variance_ratios, const vec &background_at_points, const StructureFunction &structure, int max_points, bool allow_extrapolation=true)
 Optimal interpolation for a deterministic vector of points. More...
 
vec2 optimal_interpolation_full (const Grid &bgrid, const vec2 &background, const vec2 &bvariance, const Points &obs_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. More...
 
vec optimal_interpolation_full (const Points &bpoints, const vec &background, const vec &bvariance, const Points &obs_points, const vec &obs, const vec &obs_variance, const vec &background_at_points, const vec &bvariance_at_points, const StructureFunction &structure, int max_points, vec &analysis_variance, bool allow_extrapolation=true)
 Optimal interpolation for a deterministic vector of points including analysis variance. More...
 
vec3 optimal_interpolation_ensi (const Grid &bgrid, const vec3 &background, const Points &obs_points, const vec &obs, const vec &obs_standard_deviations, const vec2 &background_at_points, const StructureFunction &structure, int max_points, bool allow_extrapolation=true)
 Optimal interpolation for an ensemble gridded field See Lussana et al 2019 (DOI: 10.1002/qj.3646) More...
 
vec2 optimal_interpolation_ensi (const Points &bpoints, const vec2 &background, const Points &obs_points, const vec &obs, const vec &obs_standard_deviations, const vec2 &background_at_points, const StructureFunction &structure, int max_points, bool allow_extrapolation=true)
 Optimal interpolation for an ensemble point field See Lussana et al 2019 (DOI: 10.1002/qj.3646) More...
 
vec2 local_distribution_correction (const Grid &bgrid, const vec2 &background, const Points &obs_points, const vec &obs, const vec &background_at_points, 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. More...
 
vec2 local_distribution_correction (const Grid &bgrid, const vec2 &background, const Points &obs_points, const vec2 &obs, const vec2 &background_at_points, 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. More...
 
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) More...
 
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. More...
 
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. More...
 
Distributions

Functions that extract values from probability distributions

vec gamma_inv (const vec &levels, const vec &shape, const vec &scale)
 Extract quantiles from a gamma distribution. More...
 
Spatial neighbourhood filters

Functions that apply neighbourhood filters on a gridded field

vec2 neighbourhood (const vec2 &input, int halfwidth, Statistic statistic)
 Spatial neighbourhood filter, computing a statistic for a sliding square window. More...
 
vec2 neighbourhood (const vec3 &input, int halfwidth, Statistic statistic)
 Spatial neighbourhood filter for an ensemble of fields. More...
 
vec2 neighbourhood_quantile (const vec2 &input, float quantile, int halfwidth)
 Computes a quantile in a sliding square neighbourhood. More...
 
vec2 neighbourhood_quantile (const vec3 &input, float quantile, int halfwidth)
 Computes a quantile in a sliding square neighbourhood across ensemble members. More...
 
vec2 neighbourhood_quantile_fast (const vec2 &input, float quantile, int halfwidth, const vec &thresholds)
 Fast and approximate neighbourhood quantile. More...
 
vec2 neighbourhood_quantile_fast (const vec3 &input, float quantile, int halfwidth, const vec &thresholds)
 Fast and approximate neighbourhood quantile for ensemble of fields. More...
 
vec2 neighbourhood_quantile_fast (const vec2 &input, const vec2 &quantile, int halfwidth, const vec &thresholds)
 Fast and approximate neighbourhood quantile, with spatially varying quantile levels. More...
 
vec2 neighbourhood_quantile_fast (const vec3 &input, const vec2 &quantile, int halfwidth, const vec &thresholds)
 Fast and approximate neighbourhood quantile for ensemble of fields, with spatially varying quantile. More...
 
vec2 neighbourhood_brute_force (const vec2 &input, int halfwidth, Statistic statistic)
 Spatial neighbourhood filter without any shortcuts. More...
 
vec2 neighbourhood_brute_force (const vec3 &input, int halfwidth, Statistic statistic)
 Spatial ensemble neighbourhood filter without any shortcuts. More...
 
vec get_neighbourhood_thresholds (const vec2 &input, int num_thresholds)
 Calculate appropriate approximation thresholds for fast neighbourhood quantile. More...
 
vec get_neighbourhood_thresholds (const vec3 &input, int num_thresholds)
 Calculate appropriate approximation thresholds for neighbourhood quantile based on an ensemble. More...
 
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. More...
 
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 that match a search criteria. More...
 
vec2 neighbourhood_ens (const vec3 &input, int halfwidth, Statistic statistic)
 Deprecated: Compute neighbourhood statistic on ensemble field. More...
 
vec2 neighbourhood_quantile_ens (const vec3 &input, float quantile, int halfwidth)
 Deprecated: Compute neighbourhood quantiles on ensemble field. More...
 
vec2 neighbourhood_quantile_ens_fast (const vec3 &input, float quantile, int radius, const vec &thresholds)
 Deprecated: Compute neighbourhood quantiles fast on ensemble field. More...
 
Calibration methods

Functions that apply statistical methods to data

vec quantile_mapping_curve (const vec &ref, const vec &fcst, vec &output_fcst, vec quantiles=vec())
 Create quantile mapping calibration curve. More...
 
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. More...
 
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. More...
 
vec apply_curve (const vec &fcst, const vec &curve_ref, const vec &curve_fcst, Extrapolation policy_below, Extrapolation policy_above)
 Apply arbitrary calibration curve to 1D forecasts. More...
 
vec2 apply_curve (const vec2 &fcst, const vec &curve_ref, const vec &curve_fcst, Extrapolation policy_below, Extrapolation policy_above)
 Apply arbitrary calibration curve to 2D forecasts. More...
 
vec2 apply_curve (const vec2 &fcst, const vec3 &curve_ref, const vec3 &curve_fcst, Extrapolation policy_below, Extrapolation policy_above)
 Apply arbitrary calibration curve to 2D forecasts with spatially varying QQ map. More...
 
vec monotonize_curve (vec curve_ref, vec curve_fcst, vec &output_fcst)
 Ensure calibration curve is monotonic, by removing points on the curve. More...
 
float get_optimal_threshold (const vec &curve_ref, const vec &curve_fcst, float threshold, Metric metric)
 Compute the optimal threshold to remap an input threshold to. More...
 
float calc_score (float a, float b, float c, float d, Metric metric)
 Compute the score for a 2x2 contingency table. More...
 
float calc_score (const vec &ref, const vec &fcst, float threshold, Metric metric)
 Compute the score for a 2x2 contingency table by thresholding. More...
 
float calc_score (const vec &ref, const vec &fcst, float threshold, float fthreshold, Metric metric)
 Compute the score for a 2x2 contingency table by thresholding and allowing a different threshold for the forecast. More...
 
Downscaling methods

Functions that interpolate data from one grid to another

vec2 downscaling (const Grid &igrid, const Grid &ogrid, const vec2 &ivalues, Downscaler downscaler)
 Downscale a gridded field. More...
 
vec3 downscaling (const Grid &igrid, const Grid &ogrid, const vec3 &ivalues, Downscaler downscaler)
 Downscale a gridded ensemble field. More...
 
vec downscaling (const Grid &igrid, const Points &opoints, const vec2 &ivalues, Downscaler downscaler)
 Downscale a gridded field to points. More...
 
vec2 downscaling (const Grid &igrid, const Points &opoints, const vec3 &ivalues, Downscaler downscaler)
 Downscale a gridded ensemble field to points. More...
 
vec2 nearest (const Grid &igrid, const Grid &ogrid, const vec2 &ivalues)
 Nearest neighbour dowscaling grid to grid for a deterministic field. More...
 
vec3 nearest (const Grid &igrid, const Grid &ogrid, const vec3 &ivalues)
 Nearest neighbour dowscaling grid to grid for an ensemble field. More...
 
vec nearest (const Grid &igrid, const Points &opoints, const vec2 &ivalues)
 Nearest neighbour dowscaling grid to point for a deterministic field. More...
 
vec2 nearest (const Grid &igrid, const Points &opoints, const vec3 &ivalues)
 Nearest neighbour dowscaling grid to point for an ensemble field. More...
 
vec nearest (const Points &ipoints, const Points &opoints, const vec &ivalues)
 Nearest neighbour dowscaling point to point for a deterministic field. More...
 
vec2 nearest (const Points &ipoints, const Points &opoints, const vec2 &ivalues)
 Nearest neighbour dowscaling point to point for an ensemble field. More...
 
vec2 nearest (const Points &ipoints, const Grid &ogrid, const vec &ivalues)
 Nearest neighbour dowscaling point to grid for a deterministic field. More...
 
vec3 nearest (const Points &ipoints, const Grid &ogrid, const vec2 &ivalues)
 Nearest neighbour dowscaling point to grid for an ensemble field. More...
 
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. More...
 
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)
 Masked ensemble downscaling and consensus. More...
 
vec2 mask_threshold_downscale_quantile (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 float quantile_level)
 Masked ensemble downscaling and quantile extraction. More...
 
vec2 bilinear (const Grid &igrid, const Grid &ogrid, const vec2 &ivalues)
 Bilinear downscaling grid to grid for a deterministic field. More...
 
vec3 bilinear (const Grid &igrid, const Grid &ogrid, const vec3 &ivalues)
 Bilinear dowscaling grid to grid for an ensemble field. More...
 
vec bilinear (const Grid &igrid, const Points &opoints, const vec2 &ivalues)
 Bilinear dowscaling grid to point for a deterministic field. More...
 
vec2 bilinear (const Grid &igrid, const Points &opoints, const vec3 &ivalues)
 Bilinear dowscaling grid to point for an ensemble field. More...
 
vec2 simple_gradient (const Grid &igrid, const Grid &ogrid, const vec2 &ivalues, float elev_gradient, Downscaler downscaler=Nearest)
 Elevation correction dowscaling grid to grid for a deterministic field. More...
 
vec3 simple_gradient (const Grid &igrid, const Grid &ogrid, const vec3 &ivalues, float elev_gradient, Downscaler downscaler=Nearest)
 Elevation correction dowscaling grid to grid for an ensemble field. More...
 
vec simple_gradient (const Grid &igrid, const Points &opoints, const vec2 &ivalues, float elev_gradient, Downscaler downscaler=Nearest)
 Elevation correction dowscaling grid to point for a deterministic field. More...
 
vec2 simple_gradient (const Grid &igrid, const Points &opoints, const vec3 &ivalues, float elev_gradient, Downscaler downscaler=Nearest)
 Elevation correction dowscaling grid to point for an ensemble field. More...
 
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. More...
 
vec3 full_gradient (const Grid &igrid, const Grid &ogrid, const vec3 &ivalues, const vec3 &elev_gradient, const vec3 &laf_gradient, Downscaler downscaler=Nearest)
 Spatially varying gradient dowscaling grid to grid for an ensemble field. More...
 
vec full_gradient (const Grid &igrid, const Points &opoints, const vec2 &ivalues, const vec2 &elev_gradient, const vec2 &laf_gradient, Downscaler downscaler=Nearest)
 Spatially varying gradient dowscaling grid to point for a deterministic field. More...
 
vec2 full_gradient (const Grid &igrid, const Points &opoints, const vec3 &ivalues, const vec3 &elev_gradient, const vec3 &laf_gradient, Downscaler downscaler=Nearest)
 Spatially varying gradient dowscaling grid to point for an ensemble field. More...
 
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)
 
vec2 smart (const Grid &igrid, const Grid &ogrid, const vec2 &ivalues, int num, const StructureFunction &structure)
 Smart neighbour downscaling grid to grid for a deterministic field. More...
 
Grid calculations

Functions that calculate statistics on a grid

vec count (const Grid &grid, const Points &points, float radius)
 For each point, counts the number of gridpoints within the radius. More...
 
vec2 count (const Grid &igrid, const Grid &ogrid, float radius)
 For each gridpoint, counts the number of gridpoints within the radius. More...
 
vec2 count (const Points &points, const Grid &grid, float radius)
 For each gridpoint, counts the number of points within the radius. More...
 
vec count (const Points &ipoints, const Points &opoints, float radius)
 For each point, counts the number of points within the radius. More...
 
vec distance (const Grid &grid, const Points &points, int num=1)
 For each point, calculates the distance to nearest gridpoint. More...
 
vec2 distance (const Grid &igrid, const Grid &ogrid, int num=1)
 For each output gridpoint, calculate the distance to nearest input gridpoint. More...
 
vec2 distance (const Points &points, const Grid &grid, int num=1)
 For each gridpoint, calculates the distance to nearest point. More...
 
vec distance (const Points &ipoints, const Points &opoint, int num=1)
 For each output point, calculates the distance to nearest input point. More...
 
vec2 fill_missing (const vec2 &values)
 Fill in missing values based on nearby values. More...
 
vec2 gridding (const Grid &grid, const Points &points, const vec &values, float radius, int min_num, Statistic statistic)
 Aggregate points onto a grid. More...
 
vec gridding (const Points &opoints, const Points &ipoints, const vec &values, float radius, int min_num, Statistic statistic)
 Aggregate points onto a points. More...
 
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. More...
 
vec gridding_nearest (const Points &opoints, const Points &ipoints, const vec &values, int min_num, gridpp::Statistic statistic)
 Assign each ipoint to nearest neighbour in opoints and aggregate values. More...
 
vec2 neighbourhood_score (const Grid &grid, const Points &points, const vec2 &fcst, const vec &ref, int half_width, gridpp::Metric metric, float threshold)
 Compute a score for a metric of all points within a radius. More...
 
Diagnosing meteorological variables

Functions that diagnose a meteorological variable based on other variables

float dewpoint (float temperature, float relative_humidity)
 Calculate dewpoint temperature from temperature and relative humidity. More...
 
vec dewpoint (const vec &temperature, const vec &relative_humidity)
 Vector version of dewpoint calculation. More...
 
float pressure (float ielev, float oelev, float ipressure, float itemperature=288.15)
 Calculate pressure at a new elevation. More...
 
vec pressure (const vec &ielev, const vec &oelev, const vec &ipressure, const vec &itemperature)
 Calculate Vector version of pressure calculation. More...
 
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. More...
 
vec sea_level_pressure (const vec &ps, const vec &altitude, const vec &temperature, const vec &rh, const vec &dewpoint)
 Vector version of Convert Surface Pressure to Sea Level Pressure. More...
 
float qnh (float pressure, float altitude)
 Diagnose QNH from pressure and altitude. More...
 
vec qnh (const vec &pressure, const vec &altitude)
 Vector version of QNH calculation. More...
 
float relative_humidity (float temperature, float dewpoint)
 Calculate relative humidity from temperature and dewpoint temperature. More...
 
vec relative_humidity (const vec &temperature, const vec &dewpoint)
 Vector version of relative humidity calculation. More...
 
float wetbulb (float temperature, float pressure, float relative_humidity)
 Calculate wetbulb temperature from temperature, pressure, and relative humidity. More...
 
vec wetbulb (const vec &temperature, const vec &pressure, const vec &relative_humidity)
 Vector version of wetbulb calculation. More...
 
float wind_speed (float xwind, float ywind)
 Diagnose wind speed from its components. More...
 
vec wind_speed (const vec &xwind, const vec &ywind)
 Vector version of wind speed calculation. More...
 
float wind_direction (float xwind, float ywind)
 Diagnose wind direction from its components. More...
 
vec wind_direction (const vec &xwind, const vec &ywind)
 Vector version of wind direction calculation. More...
 
OpenMP settings

Functions that configure OpenMP

void set_omp_threads (int num)
 Set the number of OpenMP threads to use. More...
 
int get_omp_threads ()
 Get the number of OpenMP threads currently set. More...
 
void initialize_omp ()
 Sets the number of OpenMP threads to 1 if OMP_NUM_THREADS undefined. More...
 

Variables

Constants

Functions that assimilate observations onto a gridded background

static const float MV = NAN
 Missing value indicator. More...
 
static const float MV_CML = -999
 Missing value indicator in gridpp command-line tool. More...
 
static const float pi = 3.14159265
 Mathematical constant pi. More...
 
static const double radius_earth = 6.378137e6
 Radius of the earth [m]. More...
 
static const float lapse_rate =0.0065
 Constant Lapse Rate moist air standard atmosphere [K/m]. More...
 
static const float standard_surface_temperature = 288.15
 Temperature at surface in standard atmosphere [K]. More...
 
static const float gravit = 9.80665
 Gravitational acceleration [m/s^2]. More...
 
static const float molar_mass = 0.0289644
 Molar Mass of Dry Air [kg/mol]. More...
 
static const float gas_constant_mol = 8.31447
 Universal Gas Constant [kg*m^2*s^-2/(K*mol)]. More...
 
static const float gas_constant_si = 287.05
 Universal Gas Constant [J/(kg*K)]. More...
 

Enumerations

enum  ComparisonOperator { Lt = 0, Leq = 10, Gt = 20, Geq = 30 }
 Types of comparison operators. More...
 
enum  CoordinateType { Geodetic = 0, Cartesian = 1 }
 Types of coordinates for position of points. More...
 
enum  CorrectionType { Qq = 0, Multiplicative = 10, Additive = 20 }
 Method for statistical correction. More...
 
enum  Downscaler { Nearest = 0, Bilinear = 1 }
 Types of simple downscaling methods. More...
 
enum  Extrapolation {
  OneToOne = 0, MeanSlope = 10, NearestSlope = 20, Zero = 30,
  Unchanged = 40
}
 Methods for extrapolating outside a curve. More...
 
enum  GradientType { MinMax = 0, LinearRegression = 10 }
 Types of methods to calculate the gradient. More...
 
enum  Metric {
  Ets = 0, Ts = 1, Kss = 20, Pc = 30,
  Bias = 40, Hss = 50
}
 Binary verification metrics. More...
 
enum  Statistic {
  Mean = 0, Min = 10, Median = 20, Max = 30,
  Quantile = 40, Std = 50, Variance = 60, Sum = 70,
  Count = 80, RandomChoice = 90, Unknown = -1
}
 Statistical operations to reduce a vector to a scalar. More...
 

Classes

class  BarnesStructure
 Simple structure function based on distance, elevation, and land area fraction. More...
 
class  BoxCox
 Box-Cox transformation. More...
 
class  CressmanStructure
 Simple structure function based on distance, elevation, and land area fraction. More...
 
class  CrossValidation
 
class  Gamma
 Gamma transformation. More...
 
class  Grid
 Represents a 2D grid of locations and their metadata. More...
 
class  Identity
 Identity transform, i.e. More...
 
class  KDTree
 Helper class for Grid and Points representing a tree of points. More...
 
class  Log
 Log transformation: output = log(input) More...
 
class  MultipleStructure
 
class  not_implemented_exception
 
class  Point
 Represents a single point in some coordinate system. More...
 
class  Points
 Represents a vector of locations and their metadata. More...
 
class  StructureFunction
 Covariance structure function. More...
 
class  Transform
 

Utilities


Helper functions

static int _debug_level = 0
 
void set_debug_level (int level)
 Set the verbosity of debug messages. More...
 
int get_debug_level ()
 Get the currently set level of debug messagess. More...
 
Statistic get_statistic (std::string name)
 Convert name of a statistic enums. More...
 
std::string version ()
 The gridpp version. More...
 
double clock ()
 The current time. More...
 
void debug (std::string string)
 Writes a debug message to standard out. More...
 
void warning (std::string string)
 Writes a warning message to standard out. More...
 
void error (std::string string)
 Writes an error message to standard out. More...
 
void future_deprecation_warning (std::string function, std::string other="")
 Writes an deprecation warning to standard out. More...
 
bool is_valid (float value)
 Checks if a value is valid. More...
 
float calc_statistic (const vec &array, Statistic statistic)
 Compute a statistic on a 1D vector. More...
 
float calc_quantile (const vec &array, float quantile_level)
 Compute a quantile from a 1D vector. More...
 
vec calc_statistic (const vec2 &array, Statistic statistic)
 Compute a statistic on a 2D vector. More...
 
vec calc_quantile (const vec2 &array, float quantile=MV)
 Compute a quantile from a 2D vector. More...
 
vec2 calc_quantile (const vec3 &array, const vec2 &quantile_levels)
 Compute quantile with 2D varying quantile level. More...
 
bool convert_coordinates (const vec &lats, const vec &lons, CoordinateType type, vec &x_coords, vec &y_coords, vec &z_coords)
 Convert lats/lons or 2D x/y coordinates to 3D cartesian coordinates with the centre of the earth as the origin. More...
 
bool convert_coordinates (float lat, float lon, CoordinateType type, float &x_coord, float &y_coord, float &z_coord)
 Convert lat/lon or 2D x/y coordinate to 3D cartesian coordinate with the centre of the earth as the origin. More...
 
bool is_valid_lat (float lat, CoordinateType type)
 Checks that a lat-coordinate is valid (based on the coordinate type) More...
 
bool is_valid_lon (float lon, CoordinateType type)
 Checks that a lon-coordinate is valid (based on the coordinate type) More...
 
int num_missing_values (const vec2 &iArray)
 Count the number of missing values in a vector. More...
 
int get_lower_index (float iX, const vec &iValues)
 Find the index in a vector that is equal to or just below a value. More...
 
int get_upper_index (float iX, const vec &iValues)
 Find the index in a vector that is equal to or just above a value. More...
 
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. More...
 
vec interpolate (const vec &x, const vec &iX, const vec &iY)
 Piecewise linear interpolation, vectorised version. More...
 
ivec2 init_ivec2 (int Y, int X, int value)
 Initialize a 2D integer vector of size Y, X, with a given value. More...
 
vec2 init_vec2 (int Y, int X, float value=MV)
 Initialize a 2D float vector of size Y, X, with a given value. More...
 
ivec3 init_ivec3 (int Y, int X, int E, int value)
 Initialize a 3D integer vector of size Y, X, E, with a given value. More...
 
vec3 init_vec3 (int Y, int X, int E, float value=MV)
 Initialize a 3D float vector of size Y, X, E, with a given value. More...
 
vec calc_even_quantiles (const vec &values, int num)
 Get reasonably spaced quantile levels from a vector of values, ignoring duplicate values but including the first number after duplicated values. More...
 
vec2 window (const vec2 &array, int length, gridpp::Statistic statistic, bool before=false, bool keep_missing=false, bool missing_edges=true)
 Compute window statistics across time, independently for each case. More...
 
bool compatible_size (const Grid &grid, const vec2 &v)
 Check if the grid is the same size as the 2D vector. More...
 
bool compatible_size (const Grid &grid, const vec3 &v)
 Check if the grid is compatible with 3D vector. More...
 
bool compatible_size (const Points &points, const vec &v)
 Check if the points are compatible with a 1D vector. More...
 
bool compatible_size (const Points &points, const vec2 &v)
 Check if the points are compatible with a 2D vector. More...
 
bool compatible_size (const vec2 &a, const vec2 &b)
 Check if two 2D vectors are compatible. More...
 
bool compatible_size (const vec2 &a, const vec3 &b)
 Check if a 2D and 3D vectors are compatible. More...
 
bool compatible_size (const vec3 &a, const vec3 &b)
 Check if two 3D vectors are compatible. More...
 
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. More...
 

SWIG testing functions


Functions for testing the SWIG interface.

Not useful for any other purpose.

static const float swig_default_value = -1
 Default value used to fill array in SWIG testing functions. More...
 
float * test_array (float *v, int n)
 Special function whose presense is needed for SWIG. More...
 
float test_vec_input (const vec &input)
 Testing function for 1D input vector. More...
 
int test_ivec_input (const ivec &input)
 Testing function for 1D input vector. More...
 
float test_vec2_input (const vec2 &input)
 Testing function for 2D input vector. More...
 
float test_vec3_input (const vec3 &input)
 Testing function for 3D input vector. More...
 
vec test_vec_output ()
 Testing function for 1D output vector. More...
 
ivec test_ivec_output ()
 
vec2 test_vec2_output ()
 Testing function for 2D output vector. More...
 
ivec2 test_ivec2_output ()
 
vec3 test_vec3_output ()
 Testing function for 3D output vector. More...
 
ivec3 test_ivec3_output ()
 
float test_vec_argout (vec &distances)
 Testing function for 1D vector treated as output. More...
 
float test_vec2_argout (vec2 &distances)
 Testing function for 2D vector treated as output. More...
 
void test_not_implemented_exception ()
 

Typedef Documentation

◆ dvec

typedef std::vector<double> gridpp::dvec

1D double vector

◆ dvec2

typedef std::vector<dvec> gridpp::dvec2

2D double vector

◆ ivec

typedef std::vector<int> gridpp::ivec

1D integer vector

◆ ivec2

typedef std::vector<ivec> gridpp::ivec2

2D integer vector

◆ ivec3

typedef std::vector<ivec2> gridpp::ivec3

3D integer vector

◆ vec

typedef std::vector<float> gridpp::vec

1D float vector

◆ vec2

typedef std::vector<vec> gridpp::vec2

2D float vector

◆ vec3

typedef std::vector<vec2> gridpp::vec3

3D float vector

Enumeration Type Documentation

◆ ComparisonOperator

Types of comparison operators.

Enumerator
Lt 

Lower than, <.

Leq 

Lower or equal than, <=.

Gt 

Greater than, >

Geq 

Greater or equal than, >=.

◆ CoordinateType

Types of coordinates for position of points.

Enumerator
Geodetic 

Latitude and longitude.

Cartesian 

X and Y.

◆ CorrectionType

Method for statistical correction.

Enumerator
Qq 

Quantile mapping.

Multiplicative 

Multiplicative.

Additive 

Additive.

◆ Downscaler

Types of simple downscaling methods.

Enumerator
Nearest 

Nearest neighour downscaler.

Bilinear 

Bilinear downscaler.

◆ Extrapolation

Methods for extrapolating outside a curve.

Enumerator
OneToOne 

Continue past the end-points using a slope of 1.

MeanSlope 

Continue past the end-points using the mean slope of the curve.

NearestSlope 

Continue past the end-points using the slope of the two lowermost or uppermost points in the curve.

Zero 

Continue past the end-points using a slope of 0.

Unchanged 

Keep values the way they were.

◆ GradientType

Types of methods to calculate the gradient.

Enumerator
MinMax 
LinearRegression 

◆ Metric

Binary verification metrics.

Enumerator
Ets 

Equitable threat score.

Ts 

Threat score.

Kss 

Hannsen-Kuiper skill score.

Pc 

Proportion correct.

Bias 

Bias.

Hss 

Heidke skill score.

◆ Statistic

Statistical operations to reduce a vector to a scalar.

Enumerator
Mean 

Mean of values.

Min 

Minimum of values.

Median 

Mean of values.

Max 

Maximum of values.

Quantile 

A quantile from values.

Std 

Standard deviation of values.

Variance 

Population variance of values.

Sum 

Sum of values.

Count 

Count of values.

RandomChoice 

Randomly pick a non-nan value.

Unknown 

Unknown statistic.

Function Documentation

◆ apply_curve() [1/4]

float gridpp::apply_curve ( float  fcst,
const vec curve_ref,
const vec curve_fcst,
gridpp::Extrapolation  policy_below,
gridpp::Extrapolation  policy_above 
)

Apply arbitrary calibration curve to a single value.

Parameters
fcst1D vector of forecast values to correct
curve_ref1D vector of reference curve values
curve_fcst1D vector of forecast curve values
policy_belowExtrapolation policy below curve
policy_aboveExtrapolation policy above curve
Returns
Calibrated forecasts

◆ apply_curve() [2/4]

vec gridpp::apply_curve ( const vec fcst,
const vec curve_ref,
const vec curve_fcst,
gridpp::Extrapolation  policy_below,
gridpp::Extrapolation  policy_above 
)

Apply arbitrary calibration curve to 1D forecasts.

Parameters
fcst1D vector of forecast values to correct
curve_ref1D vector of reference curve values
curve_fcst1D vector of forecast curve values
policy_belowExtrapolation policy below curve
policy_aboveExtrapolation policy above curve
Returns
1D vector of calibrated forecasts

◆ apply_curve() [3/4]

vec2 gridpp::apply_curve ( const vec2 fcst,
const vec curve_ref,
const vec curve_fcst,
gridpp::Extrapolation  policy_below,
gridpp::Extrapolation  policy_above 
)

Apply arbitrary calibration curve to 2D forecasts.

Parameters
fcst2D vector of forecast values to correct
curve_ref1D vector of reference curve values
curve_fcst1D vector of forecast curve values
policy_belowExtrapolation policy below curve
policy_aboveExtrapolation policy above curve
Returns
2D array of calibrated forecasts (Y, X)

◆ apply_curve() [4/4]

vec2 gridpp::apply_curve ( const vec2 fcst,
const vec3 curve_ref,
const vec3 curve_fcst,
gridpp::Extrapolation  policy_below,
gridpp::Extrapolation  policy_above 
)

Apply arbitrary calibration curve to 2D forecasts with spatially varying QQ map.

Parameters
fcst2D vector of forecast values to correct
curve_ref3D vector of reference curve values (Y, X, Curve)
curve_fcst3D vector of Forecast curve values (Y, X, Curve)
policy_belowExtrapolation policy below curve
policy_aboveExtrapolation policy above curve
Returns
2D vector of calibrated forecasts (Y, X)

◆ bilinear() [1/4]

vec2 gridpp::bilinear ( const Grid igrid,
const Grid ogrid,
const vec2 ivalues 
)

Bilinear downscaling grid to grid for a deterministic field.

Parameters
igridGrid of input (Yi, Xi)
ogridGrid of downscaled output (Yo, Xo)
ivalues2D vector of input values (Yi, Xi)
Returns
2D vector of output values (Yo, Xo)

◆ bilinear() [2/4]

vec3 gridpp::bilinear ( const Grid igrid,
const Grid ogrid,
const vec3 ivalues 
)

Bilinear dowscaling grid to grid for an ensemble field.

Parameters
igridGrid of input (Yi, Xi)
ogridGrid of downscaled output (Yo, Xo)
ivalues3D vector of input values (E, Yi, Xi)
Returns
3D vector of output values (E, Yo, Xo)

◆ bilinear() [3/4]

vec gridpp::bilinear ( const Grid igrid,
const Points opoints,
const vec2 ivalues 
)

Bilinear dowscaling grid to point for a deterministic field.

Parameters
igridGrid of input (Yi, Xi)
opointsPoints of downscaled output
ivalues2D vector of input values grid (Yi, Xi)
Returns
1D vector of output values

◆ bilinear() [4/4]

vec2 gridpp::bilinear ( const Grid igrid,
const Points opoints,
const vec3 ivalues 
)

Bilinear dowscaling grid to point for an ensemble field.

Parameters
igridGrid of input (Yi, Xi)
opointsPoints of downscaled output
ivalues3D vector of input values (E, Yi, Xi)
Returns
2D vector of output values (E, Points)

◆ calc_even_quantiles()

vec gridpp::calc_even_quantiles ( const vec values,
int  num 
)

Get reasonably spaced quantile levels from a vector of values, ignoring duplicate values but including the first number after duplicated values.

Include the lowest and highest values.

Parameters
values1D vector of values (unsorted, and no invalid values)
numnumber of thresholds to get
Returns
1D vector of sorted quantile levels

◆ calc_gradient()

vec2 gridpp::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.

Parameters
base2D vector of dependent values, missing values are ignored (Y, X)
values2D vector of independent values, missing values are ignored (Y, X)
gradient_typeWhat gradient type to compute
halfwidthNeighbourhood half width in number of gridpoints
min_nimMinimum number of points required to compute gradient
min_rangeMinimum range of base to compute gradient
default_gradientUse this gradient if minimum number of points not available
Returns
2D vector of gradients (Y, X)

◆ calc_quantile() [1/3]

float gridpp::calc_quantile ( const vec array,
float  quantile_level 
)

Compute a quantile from a 1D vector.

Parameters
array1D vector of values
quantile_levelQuantile level to compute
Returns
Computed quantile

◆ calc_quantile() [2/3]

vec gridpp::calc_quantile ( const vec2 array,
float  quantile = MV 
)

Compute a quantile from a 2D vector.

Parameters
array2D vector of values (Y, X)
quantile_levelQuantile level to compute
Returns
1D vector of computed quantiles (Y)

◆ calc_quantile() [3/3]

vec2 gridpp::calc_quantile ( const vec3 array,
const vec2 quantile_levels 
)

Compute quantile with 2D varying quantile level.

Parameters
array3D input vector (T, Y, X)
quantile2D vector of Quantile levels (Y, X)
Returns
Extracted quantiles (Y, X)

◆ calc_score() [1/3]

float gridpp::calc_score ( float  a,
float  b,
float  c,
float  d,
gridpp::Metric  metric 
)

Compute the score for a 2x2 contingency table.

Parameters
aFraction of hits
bFraction of false alarms
cFraction of misses
dFraction of correct rejections
metricMetric to conpute score for
Returns
Score

◆ calc_score() [2/3]

float gridpp::calc_score ( const vec ref,
const vec fcst,
float  threshold,
gridpp::Metric  metric 
)

Compute the score for a 2x2 contingency table by thresholding.

Parameters
ref1D vector of reference values
fcst1D vector forecast values
thresholdThreshold to compute score for
metricMetric to conpute score for
Returns
Score

◆ calc_score() [3/3]

float gridpp::calc_score ( const vec ref,
const vec fcst,
float  threshold,
float  fthreshold,
gridpp::Metric  metric 
)

Compute the score for a 2x2 contingency table by thresholding and allowing a different threshold for the forecast.

Parameters
ref1D vector of reference values
fcst1D vector forecast values
thresholdThreshold for the reference values
fthresholdThreshold for the forecast values
metricMetric to conpute score for
Returns
Score

◆ calc_statistic() [1/2]

float gridpp::calc_statistic ( const vec array,
gridpp::Statistic  statistic 
)

Compute a statistic on a 1D vector.

Parameters
array1D vector of values
statisticStatistic to compute
Returns
Computed value

◆ calc_statistic() [2/2]

vec gridpp::calc_statistic ( const vec2 array,
gridpp::Statistic  statistic 
)

Compute a statistic on a 2D vector.

Parameters
array2D vector of values (Y, X)
statisticStatistic to compute
Returns
1D vector of computed values (Y)

◆ clock()

double gridpp::clock ( )

The current time.

Returns
The current time in seconds since 1870-01-01 00:00:00Z [s]

◆ compatible_size() [1/7]

bool gridpp::compatible_size ( const Grid grid,
const vec2 v 
)

Check if the grid is the same size as the 2D vector.

Parameters
gridGrid (Y, X)
v2D vector (Y, X)
Returns
True if compatible

◆ compatible_size() [2/7]

bool gridpp::compatible_size ( const Grid grid,
const vec3 v 
)

Check if the grid is compatible with 3D vector.

Parameters
gridGrid (Y, X)
v3D vector (T, Y, X)
Returns
True if compatible

◆ compatible_size() [3/7]

bool gridpp::compatible_size ( const Points points,
const vec v 
)

Check if the points are compatible with a 1D vector.

Parameters
gridPoints
v1D vector
Returns
True if compatible

◆ compatible_size() [4/7]

bool gridpp::compatible_size ( const Points points,
const vec2 v 
)

Check if the points are compatible with a 2D vector.

Parameters
gridPoints (P)
v2D vector (T, P)
Returns
True if compatible

◆ compatible_size() [5/7]

bool gridpp::compatible_size ( const vec2 a,
const vec2 b 
)

Check if two 2D vectors are compatible.

Parameters
a2D vector (Y, X)
b2D vector (Y, X)
Returns
True if compatible

◆ compatible_size() [6/7]

bool gridpp::compatible_size ( const vec2 a,
const vec3 b 
)

Check if a 2D and 3D vectors are compatible.

Parameters
a2D vector (Y, X)
b3D vector (Y, X, E)
Returns
True if compatible

◆ compatible_size() [7/7]

bool gridpp::compatible_size ( const vec3 a,
const vec3 b 
)

Check if two 3D vectors are compatible.

Parameters
a3D vector (Y, X, E)
b3D vector (Y, X, E)
Returns
True if compatible

◆ convert_coordinates() [1/2]

bool gridpp::convert_coordinates ( const vec lats,
const vec lons,
CoordinateType  type,
vec x_coords,
vec y_coords,
vec z_coords 
)

Convert lats/lons or 2D x/y coordinates to 3D cartesian coordinates with the centre of the earth as the origin.

Parameters
latsvector of latitudes [deg] or 2D y-coordinates [m]
lonsvector of longitudes [deg] or 2D x-coordinates [m]
typeCoordinate type for lats and lons
x_coordsvector of x-coordinates [m]
y_coordsvector of y-coordinates [m]
z_coordsvector of z-coordinates [m]

◆ convert_coordinates() [2/2]

bool gridpp::convert_coordinates ( float  lat,
float  lon,
CoordinateType  type,
float &  x_coord,
float &  y_coord,
float &  z_coord 
)

Convert lat/lon or 2D x/y coordinate to 3D cartesian coordinate with the centre of the earth as the origin.

Parameters
latvector of latitudes [deg] or 2D y-coordinates [m]
lonvector of longitudes [deg] or 2D x-coordinates [m]
typeCoordinate type for lat and lon
x_coordx-coordinate [m]
y_coordy-coordinate [m]
z_coordz-coordinate [m]

◆ count() [1/4]

vec gridpp::count ( const Grid grid,
const Points points,
float  radius 
)

For each point, counts the number of gridpoints within the radius.

Parameters
gridInput grid
pointsOutput points
radiusRadius [m]
Returns
1D vector of number of gridpoints

◆ count() [2/4]

vec2 gridpp::count ( const Grid igrid,
const Grid ogrid,
float  radius 
)

For each gridpoint, counts the number of gridpoints within the radius.

Parameters
igridInput grid
ogridOutput grid
radiusRadius [m]
Returns
2D vector of number of gridpoints

◆ count() [3/4]

vec2 gridpp::count ( const Points points,
const Grid grid,
float  radius 
)

For each gridpoint, counts the number of points within the radius.

Parameters
pointsInput points
gridOutput grid
radiusRadius [m]
Returns
2D vector of number of points

◆ count() [4/4]

vec gridpp::count ( const Points ipoints,
const Points opoints,
float  radius 
)

For each point, counts the number of points within the radius.

Parameters
ipointsInput points (Pi)
opointsOutput points (Po)
radiusRadius [m]
Returns
1D vector of number of points (Po)

◆ debug()

void gridpp::debug ( std::string  string)

Writes a debug message to standard out.

Parameters
stringWrite this string

◆ dewpoint() [1/2]

float gridpp::dewpoint ( float  temperature,
float  relative_humidity 
)

Calculate dewpoint temperature from temperature and relative humidity.

Parameters
temperatureTemperature [K]
relative_humidityRelative humidity [1]
Returns
Dewpoint temperature [K]

◆ dewpoint() [2/2]

vec gridpp::dewpoint ( const vec temperature,
const vec relative_humidity 
)

Vector version of dewpoint calculation.

Parameters
temperature1D vector of temperature [K]
relative_humidity1D vector of relative humidity [1]
Returns
1D vector of dewpoint temperature [K]

◆ distance() [1/4]

vec gridpp::distance ( const Grid grid,
const Points points,
int  num = 1 
)

For each point, calculates the distance to nearest gridpoint.

Parameters
gridGrid
pointsPoints
numNumber of points
Returns
1D vector of distance [m] to nearest gridpoint for each point

◆ distance() [2/4]

vec2 gridpp::distance ( const Grid igrid,
const Grid ogrid,
int  num = 1 
)

For each output gridpoint, calculate the distance to nearest input gridpoint.

Parameters
gridGrid (Yi, Xi)
ogridOutput grid /Yo, Xo)
numNumber of points
Returns
2D vector of distance [m] to nearest gridpoint for each gridpoint (Yo, Xo)

◆ distance() [3/4]

vec2 gridpp::distance ( const Points points,
const Grid grid,
int  num = 1 
)

For each gridpoint, calculates the distance to nearest point.

Parameters
pointsPoints
gridGrid
numNumber of points
Returns
2D vector of distance [m] to nearest point for each gridpoint

◆ distance() [4/4]

vec gridpp::distance ( const Points ipoints,
const Points opoint,
int  num = 1 
)

For each output point, calculates the distance to nearest input point.

Parameters
ipointsInput points (Pi)
opointsOutput points (Po)
numNumber of points
Returns
1D vector of distance [m] to nearest point for each point (Po)

◆ doping_circle()

vec2 gridpp::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.

Parameters
gridGrid of background
background2D vector of background values with (Y, X)
pointsPoints of observations
observations1D vector of observations
radiiRadius of circle where observations are inserted for each point [m]
max_elev_diffOnly insert where elevation difference between grid and point is less than this value
Returns
2D vector of final field

◆ doping_square()

vec2 gridpp::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.

Parameters
gridGrid of background
background2D vector of background values (Y, X)
pointsPoints of observations
observations1D vector of observations
halfwidthsHalf width of square (in number of grid points) where observations are inserted for each point
max_elev_diffOnly insert where elevation difference between grid and point is less than this value [m]
Returns
2D vector of final field

◆ downscale_probability()

vec2 gridpp::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.

Parameters
igridGrid of input (Yi, Xi)
ogridGrid of downscaled output (Yo, Xo)
ivalues3D vector of input values (Yi, Xi, E)
threshold2D vector of threshold values (Yo, Xo)
comparison_operatorComparison operator to create probability
Returns
2D vector of output values (Yo, Xo)

◆ downscaling() [1/4]

vec2 gridpp::downscaling ( const Grid igrid,
const Grid ogrid,
const vec2 ivalues,
Downscaler  downscaler 
)

Downscale a gridded field.

Parameters
igridGrid of input (Yi, Xi)
ogridGrid of downscaled output (Yo, Xo)
ivalues2D vector of values on the input grid (Yi, Xi)
downscalerDownscaling method
Returns
2D vector of output values (Yo, Xo)

◆ downscaling() [2/4]

vec3 gridpp::downscaling ( const Grid igrid,
const Grid ogrid,
const vec3 ivalues,
Downscaler  downscaler 
)

Downscale a gridded ensemble field.

Parameters
igridGrid of input (Yi, Xi)
ogridGrid of downscaled output (Yo, Xo)
ivalues3D vector of values on the input grid (E, Yi, Xi)
downscalerDownscaling method
Returns
3D vector of output values (E, Yo, Xo)

◆ downscaling() [3/4]

vec gridpp::downscaling ( const Grid igrid,
const Points opoints,
const vec2 ivalues,
Downscaler  downscaler 
)

Downscale a gridded field to points.

Parameters
igridGrid of input (Yi, Xi)
opointsPoints of downscaled output
ivalues2D vector of values on the input grid (Yi, Xi)
downscalerDownscaling method
Returns
1D vector of output values

◆ downscaling() [4/4]

vec2 gridpp::downscaling ( const Grid igrid,
const Points opoints,
const vec3 ivalues,
Downscaler  downscaler 
)

Downscale a gridded ensemble field to points.

Parameters
igridGrid of input (Yi, Xi)
opointsPoints of downscaled output
ivalues3D vector of values on the input grid (E, Yi, Xi)
downscalerDownscaling method
Returns
2D vector of output values (E, Points)

◆ error()

void gridpp::error ( std::string  string)

Writes an error message to standard out.

Parameters
stringWrite this string

◆ fill()

vec2 gridpp::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)

Parameters
igridGrid of input
input2D vector of values (Y, X)
pointsPoints to fill in
radii1D vector of circle radii for each point [m]
valueValue to fill in
outsideIf True, fill outside circles, if False, fill inside circles
Returns
2D vector of final field

◆ fill_missing()

vec2 gridpp::fill_missing ( const vec2 values)

Fill in missing values based on nearby values.

Parameters
values2D vector of values
Returns
2D vector of values without any missing values

◆ full_gradient() [1/4]

vec2 gridpp::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.

Parameters
igridinput grid
ogridoutput grid
ivaluesvalues from igrid
elev_gradientelevation gradient
laf_gradientland area fraction gradientSpatially varying gradient dowscaling grid to grid for a deterministic field
igridGrid of input (Yi, Xi)
ogridGrid of downscaled output (Yo, Xo)
ivalues2D vector of input values (Yi, Xi)
elev_gradient2D vector of elevation gradients (Yi, Xi) [input unit/m]
laf_gradient2D vector of land-area-fraction gradients (Yi, Xi) [input unit/1]
downscalerDownscaling method
Returns
2D vector of output values (Yo, Xo)

◆ full_gradient() [2/4]

vec3 gridpp::full_gradient ( const Grid igrid,
const Grid ogrid,
const vec3 ivalues,
const vec3 elev_gradient,
const vec3 laf_gradient,
Downscaler  downscaler = Nearest 
)

Spatially varying gradient dowscaling grid to grid for an ensemble field.

Parameters
igridGrid of input (Yi, Xi)
ogridGrid of downscaled output (Yo, Xo)
ivalues3D vector of input values (E, Yi, Xi)
elev_gradient3D vector of elevation gradients (Yi, Xi, E) [input unit/m]
laf_gradient3D vector of land-area-fraction gradients (Yi, Xi, E) [input unit/1]
downscalerDownscaling method
Returns
3D vector of output values (E, Yo, Xo)

◆ full_gradient() [3/4]

vec gridpp::full_gradient ( const Grid igrid,
const Points opoints,
const vec2 ivalues,
const vec2 elev_gradient,
const vec2 laf_gradient,
Downscaler  downscaler = Nearest 
)

Spatially varying gradient dowscaling grid to point for a deterministic field.

Parameters
igridGrid of input (Yi, Xi)
opointsPoints of downscaled output
ivalues2D vector of input values grid (Yi, Xi)
elev_gradient2D vector of elevation gradients (Yi, Xi) [input unit/m]
laf_gradient2D vector of land-area-fraction gradients (Yi, Xi) [input unit/1]
downscalerDownscaling method
Returns
1D vector of output values

◆ full_gradient() [4/4]

vec2 gridpp::full_gradient ( const Grid igrid,
const Points opoints,
const vec3 ivalues,
const vec3 elev_gradient,
const vec3 laf_gradient,
Downscaler  downscaler = Nearest 
)

Spatially varying gradient dowscaling grid to point for an ensemble field.

Parameters
igridGrid of input (Yi, Xi)
opointsPoints of downscaled output
ivalues3D vector of input values (E, Yi, Xi)
elev_gradient3D vector of elevation gradients (E, Yi, Xi) [input unit/m]
laf_gradient3D vector of land-area-fraction gradients (E, Yi, Xi) [input unit/1]
downscalerDownscaling method
Returns
2D vector of output values (E, Points)

◆ full_gradient_debug()

vec3 gridpp::full_gradient_debug ( const Grid igrid,
const Grid ogrid,
const vec2 ivalues,
const vec2 elev_gradient,
const vec2 laf_gradient = vec2(),
Downscaler  downscaler = Nearest 
)

◆ future_deprecation_warning()

void gridpp::future_deprecation_warning ( std::string  function,
std::string  other = "" 
)

Writes an deprecation warning to standard out.

Parameters
stringName of the deprecated function
otherAdditional message to write

◆ gamma_inv()

vec gridpp::gamma_inv ( const vec levels,
const vec shape,
const vec scale 
)

Extract quantiles from a gamma distribution.

Parameters
levels1D vector of quantile levels to retrieve
shape1D vector of shape parameter of gamma distribution
scale1D vector of scale parameter of gamma distribution
Returns
1D vector of quantiles

◆ get_debug_level()

int gridpp::get_debug_level ( )

Get the currently set level of debug messagess.

Returns
Currently set debug level

◆ get_lower_index()

int gridpp::get_lower_index ( float  iX,
const vec iValues 
)

Find the index in a vector that is equal to or just below a value.

Parameters
iXLookup value
iValues1D vector of lookup values. Must be sorted.
Returns
The index into iValues that is equal to or just below iX

◆ get_neighbourhood_thresholds() [1/2]

vec gridpp::get_neighbourhood_thresholds ( const vec2 input,
int  num_thresholds 
)

Calculate appropriate approximation thresholds for fast neighbourhood quantile.

Parameters
input2D vector of values (Y, X)
num_thresholdsNumber of thresholds
Returns
1D vector of approximation thresholds

◆ get_neighbourhood_thresholds() [2/2]

vec gridpp::get_neighbourhood_thresholds ( const vec3 input,
int  num_thresholds 
)

Calculate appropriate approximation thresholds for neighbourhood quantile based on an ensemble.

Parameters
input3D vector of values (Y, X, T)
num_thresholdsNumber of thresholds
Returns
1D vector of approximation thresholds

◆ get_omp_threads()

int gridpp::get_omp_threads ( )

Get the number of OpenMP threads currently set.

Returns
Number of threads

◆ get_optimal_threshold()

float gridpp::get_optimal_threshold ( const vec curve_ref,
const vec curve_fcst,
float  threshold,
gridpp::Metric  metric 
)

Compute the optimal threshold to remap an input threshold to.

Parameters
curve_ref1D vector of reference curve values
curve_fcst1D vector of reference curve values
thresholdInput threshold
metricMetric to optimize
Returns
Optimal remapping threshold

◆ get_statistic()

gridpp::Statistic gridpp::get_statistic ( std::string  name)

Convert name of a statistic enums.

Parameters
nameName of a statistic
Returns
Statistic

◆ get_upper_index()

int gridpp::get_upper_index ( float  iX,
const vec iValues 
)

Find the index in a vector that is equal to or just above a value.

Parameters
iXLookup value
iValues1D vector of lookup values. Must be sorted.
Returns
The index into iValues that is equal to or just above iX

◆ gridding() [1/2]

vec2 gridpp::gridding ( const Grid grid,
const Points points,
const vec values,
float  radius,
int  min_num,
gridpp::Statistic  statistic 
)

Aggregate points onto a grid.

Writes gridpp::MV where there are not enough points.

Parameters
gridOutput grid to aggregate to
pointsInput points
values1D vector of input values
radiusCircle radius for aggregating points [m]
min_numMinimum number of points in radius to create a value
statisticStatistic to compute on points within radius
Returns
2D vector of output values

◆ gridding() [2/2]

vec gridpp::gridding ( const Points opoints,
const Points ipoints,
const vec values,
float  radius,
int  min_num,
gridpp::Statistic  statistic 
)

Aggregate points onto a points.

Writes gridpp::MV where there are not enough points.

Parameters
opointsOutput points to aggregate to
ipointsInput points
values1D vector of input Values
radiusCircle radius for aggregate points [m]
min_numMinimum number of points in radius to create a value
statisticStatistic to compute on points within radius
Returns
1D vector of output values

◆ gridding_nearest() [1/2]

vec2 gridpp::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.

Writes gridpp::MV where there are not enough observations.

Parameters
gridOutput grid to aggregate to
pointsInput points
valuesInput values
min_numMinimum number of points in a gridpoint to create a value
statisticStatistic to compute on points within gridpoint
Returns
2D vector of output values

◆ gridding_nearest() [2/2]

vec gridpp::gridding_nearest ( const Points opoints,
const Points ipoints,
const vec values,
int  min_num,
gridpp::Statistic  statistic 
)

Assign each ipoint to nearest neighbour in opoints and aggregate values.

Writes gridpp::MV where there are not enough observations.

Parameters
opointsOutput points to aggregate to
ipointsInput points
valuesInput values
min_numMinimum number of points in a gridpoint to create a value
statisticStatistic to compute on points within gridpoint
Returns
1D vector of output values

◆ init_ivec2()

ivec2 gridpp::init_ivec2 ( int  Y,
int  X,
int  value 
)

Initialize a 2D integer vector of size Y, X, with a given value.

Parameters
YSize of first dimension
XSize of second dimension value Initialize with this value
Returns
2D interger vector

◆ init_ivec3()

ivec3 gridpp::init_ivec3 ( int  Y,
int  X,
int  E,
int  value 
)

Initialize a 3D integer vector of size Y, X, E, with a given value.

Parameters
YSize of first dimension
XSize of second dimension
ESize of third dimension value Initialize with this value
Returns
3D integer vector

◆ init_vec2()

vec2 gridpp::init_vec2 ( int  Y,
int  X,
float  value = MV 
)

Initialize a 2D float vector of size Y, X, with a given value.

Parameters
YSize of first dimension
XSize of second dimension value Initialize with this value
Returns
2D float vector

◆ init_vec3()

vec3 gridpp::init_vec3 ( int  Y,
int  X,
int  E,
float  value = MV 
)

Initialize a 3D float vector of size Y, X, E, with a given value.

Parameters
YSize of first dimension
XSize of second dimension
ESize of third dimension value Initialize with this value
Returns
3D float vector

◆ initialize_omp()

void gridpp::initialize_omp ( )

Sets the number of OpenMP threads to 1 if OMP_NUM_THREADS undefined.

◆ interpolate() [1/2]

float gridpp::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.

If there are multiple identical x values, then the average of the y values at each end of the x-interval that is constant is used. The exception is if the constant interval is on one (only) of the edges of the interpolation curve. In that case, the y-value at the end of the interval further away from the boundary of the curve is used.

Parameters
xInterpolation to this point
iX1D vector of x-axis values. Vector must be sorted.
iY1D vector of y-axis values corresponding to iX.
Returns
Y value corresponding to x

◆ interpolate() [2/2]

vec gridpp::interpolate ( const vec x,
const vec iX,
const vec iY 
)

Piecewise linear interpolation, vectorised version.

Parameters
x1D vector of values to interpolate to
iX1D vector of x-axis values. Vector must be sorted.
iY1D vector of y-axis values corresponding to iX.
Returns
1D vector of y-axis values corresponding to x

◆ is_valid()

bool gridpp::is_valid ( float  value)

Checks if a value is valid.

Parameters
valueValue to check
Returns
True if the value is valid, False if it is NaN, Inf, or other invalid value

◆ is_valid_lat()

bool gridpp::is_valid_lat ( float  lat,
CoordinateType  type 
)

Checks that a lat-coordinate is valid (based on the coordinate type)

Parameters
latlatitude [deg] or y-coordinate [m]
typeCoordinate type for lat
Returns
True if lat is valid latitude or y-coordinate

◆ is_valid_lon()

bool gridpp::is_valid_lon ( float  lon,
CoordinateType  type 
)

Checks that a lon-coordinate is valid (based on the coordinate type)

Parameters
latlongitude [deg] or x-coordinate [m]
typeCoordinate type for lon
Returns
True if lon is valid latitude or x-coordinate

◆ local_distribution_correction() [1/2]

vec2 gridpp::local_distribution_correction ( const Grid bgrid,
const vec2 background,
const Points obs_points,
const vec obs,
const vec background_at_points,
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.

This is an experimental method.

Parameters
bgridGrid of background field
background2D vector of background values (Y, X)
obs_pointsPoints of observations
obs1D vector of observations
background_at_points1D vector of background values at observation points
structureStructure function
min_quantileTruncate quantile map below this quantile level
max_quantileTruncate quantile map above this quantile level
max_pointsMaximum number of points used within localization radius (not used at moment)
Returns
2D vector of corrected values

◆ local_distribution_correction() [2/2]

vec2 gridpp::local_distribution_correction ( const Grid bgrid,
const vec2 background,
const Points obs_points,
const vec2 obs,
const vec2 background_at_points,
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.

This is an experimental method. Version with multiple number of timesteps. Pool all observations across time in when computing the calibration curve.

Parameters
bgridGrid of background
background2D vector of background values (Y, X)
obs_pointsPoints of observations
obs2D vector of observations (Time, Points)
background_at_points2D vector of background values at observation points (Time, Points)
structurestructure function
min_quantileTruncate quantile map below this quantile
max_quantileTruncate quantile map above this quantile
max_pointsMaximum number of points used within localization radius (not used at moment)
Returns
2D vector of corrected values

◆ mask_threshold_downscale_consensus()

vec2 gridpp::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 
)

Masked ensemble downscaling and consensus.

This method computes a high-resolution consensus from a low-resolution ensemble, without downscaling all members. See the wiki for more details.

Parameters
igridGrid of input (Yi, Xi)
ogridGrid of downscaled output (Yo, Xo)
ivalues_true3D vector of input values to use when condition is true (Yi, Xi, E)
ivalues_false3D vector of input values to use when condition is false (Yi, Xi, E)
threshold_values3D vector of values (Yi, Xi, E) defining the mask for ivalues after applying the theshold
threshold2D vector of thresholds on output grid that is upscaled to input grid (Yo, Xo)
comparison_operatorComparison operator to create probability
statisticstatistic to compute over the ensemble dimension after downscaling
Returns
2D vector of output values (Yo, Xo)

◆ mask_threshold_downscale_quantile()

vec2 gridpp::mask_threshold_downscale_quantile ( 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 float  quantile_level 
)

Masked ensemble downscaling and quantile extraction.

This method extractes a quantile for a high-resolution grid from a low-resolution ensemble, without downscaling all members. See the wiki for more details.

Parameters
igridGrid of input (Yi, Xi)
ogridGrid of downscaled output (Yo, Xo)
ivalues_true3D vector of input values to use when condition is true (Yi, Xi, E)
ivalues_false3D vector of input values to use when condition is false (Yi, Xi, E)
threshold_values3D vector of values (Yi, Xi, E) defining the mask for ivalues after applying the theshold
threshold2D vector of thresholds on output grid that is upscaled to input grid (Yo, Xo)
comparison_operatorComparison operator to create probability
quantile_levelQuantile level (between 0 and 1) to extract
Returns
2D vector of output values (Yo, Xo)

◆ metric_optimizer_curve()

vec gridpp::metric_optimizer_curve ( const vec ref,
const vec fcst,
const vec thresholds,
gridpp::Metric  metric,
vec output_fcst 
)

Create calibration curve that optimizes a metric.

Parameters
ref1D vector of reference values (observations)
fcst1d vector of forecast values
thresholds1D vector of thresholds for computing optimal values for
metricA Metric to optimize for
output_fcst1D vector of output forecast quantiles
Returns
1D vector of output reference quantiles

◆ monotonize_curve()

vec gridpp::monotonize_curve ( vec  curve_ref,
vec  curve_fcst,
vec output_fcst 
)

Ensure calibration curve is monotonic, by removing points on the curve.

Parameters
curve_ref1D vector of reference curve values
curve_fcst1D vector of forecast curve values
output_fcst1D vector of output forecast curve values
Returns
1D vector of output reference curve values

◆ nearest() [1/8]

vec2 gridpp::nearest ( const Grid igrid,
const Grid ogrid,
const vec2 ivalues 
)

Nearest neighbour dowscaling grid to grid for a deterministic field.

Parameters
igridGrid of input (Yi, Xi)
ogridGrid of downscaled output (Yo, Xo)
ivalues2D vector of input values (Yi, Xi)
Returns
2D vector of output values (Yo, Xo)

◆ nearest() [2/8]

vec3 gridpp::nearest ( const Grid igrid,
const Grid ogrid,
const vec3 ivalues 
)

Nearest neighbour dowscaling grid to grid for an ensemble field.

Parameters
igridGrid of input (Yi, Xi)
ogridGrid of downscaled output (Yo, Xo)
ivalues3D vector of input values (E, Yi, Xi)
Returns
3D vector of output values (E, Yo, Xo)

◆ nearest() [3/8]

vec gridpp::nearest ( const Grid igrid,
const Points opoints,
const vec2 ivalues 
)

Nearest neighbour dowscaling grid to point for a deterministic field.

Parameters
igridGrid of input (Yi, Xi)
opointsPoints of downscaled output
ivalues2D vector of input values grid (Yi, Xi)
Returns
1D vector of output values

◆ nearest() [4/8]

vec2 gridpp::nearest ( const Grid igrid,
const Points opoints,
const vec3 ivalues 
)

Nearest neighbour dowscaling grid to point for an ensemble field.

Parameters
igridGrid of input (Yi, Xi)
opointsPoints of downscaled output
ivalues3D vector of input values (E, Yi, Xi)
Returns
2D vector of output values (E, Points)

◆ nearest() [5/8]

vec gridpp::nearest ( const Points ipoints,
const Points opoints,
const vec ivalues 
)

Nearest neighbour dowscaling point to point for a deterministic field.

Parameters
ipointsPoints of input
opointsPoints of downscaled output
ivalues1D vector of input values
Returns
Values 1D vector of output values

◆ nearest() [6/8]

vec2 gridpp::nearest ( const Points ipoints,
const Points opoints,
const vec2 ivalues 
)

Nearest neighbour dowscaling point to point for an ensemble field.

Parameters
ipointsPoints of input (Pi)
opointsPoints of downscaled output Po)
ivalues2D vector of input values (E, Pi)
Returns
2D vector of output values (E, Po)

◆ nearest() [7/8]

vec2 gridpp::nearest ( const Points ipoints,
const Grid ogrid,
const vec ivalues 
)

Nearest neighbour dowscaling point to grid for a deterministic field.

Parameters
ipointsPoints of input
ogridGrid of downscaled output
ivalues1D vector of values for the input points
Returns
2D vector of output values (Y, X)

◆ nearest() [8/8]

vec3 gridpp::nearest ( const Points ipoints,
const Grid ogrid,
const vec2 ivalues 
)

Nearest neighbour dowscaling point to grid for an ensemble field.

Parameters
ipointsPoints of input
ogridGrid of downscaled output
ivalues2D vector input values (E, Points)
Returns
3D vector of output values (E, Y, X)

◆ neighbourhood() [1/2]

vec2 gridpp::neighbourhood ( const vec2 input,
int  halfwidth,
gridpp::Statistic  statistic 
)

Spatial neighbourhood filter, computing a statistic for a sliding square window.

Parameters
input2D vector of input values (Y, X)
halfwidthFilter halfwidth in number of gridpoints
statisticStatistic to compute
Returns
2D vector of final field (Y, x)

◆ neighbourhood() [2/2]

vec2 gridpp::neighbourhood ( const vec3 input,
int  halfwidth,
gridpp::Statistic  statistic 
)

Spatial neighbourhood filter for an ensemble of fields.

Computes static across neighbourhood and all members

Parameters
input3D vector of input values (Y, X, E)
halfwidthFilter halfwidth in number of gridpoints
statisticStatistic to compute
Returns
2D vector of final field (Y, x)

◆ neighbourhood_brute_force() [1/2]

vec2 gridpp::neighbourhood_brute_force ( const vec2 input,
int  halfwidth,
gridpp::Statistic  statistic 
)

Spatial neighbourhood filter without any shortcuts.

This is quite slow and is only useful for testing.

Parameters
input2D vector of values (Y, X)
halfwidthFilter halfwidth in number of gridpoints
statisticStatistic to compute
Returns
2D vector of final field (Y, X)

◆ neighbourhood_brute_force() [2/2]

vec2 gridpp::neighbourhood_brute_force ( const vec3 input,
int  halfwidth,
gridpp::Statistic  statistic 
)

Spatial ensemble neighbourhood filter without any shortcuts.

This is quite slow and is only useful for testing.

Parameters
input3D vector of values (Y, X, E)
halfwidthFilter halfwidth in number of gridpoints
statisticStatistic to compute
Returns
2D vector of final field (Y, X)

◆ neighbourhood_ens()

vec2 gridpp::neighbourhood_ens ( const vec3 input,
int  halfwidth,
gridpp::Statistic  statistic 
)

Deprecated: Compute neighbourhood statistic on ensemble field.

Deprecated:
Use neighbourhood() function

◆ neighbourhood_quantile() [1/2]

vec2 gridpp::neighbourhood_quantile ( const vec2 input,
float  quantile,
int  halfwidth 
)

Computes a quantile in a sliding square neighbourhood.

Parameters
input2D vector of values
quantileQuantile to compute (between 0 and 1)
halfwidthFilter halfwidth in number of gridpoints
Returns
2D vector of final field (Y, x)

◆ neighbourhood_quantile() [2/2]

vec2 gridpp::neighbourhood_quantile ( const vec3 input,
float  quantile,
int  halfwidth 
)

Computes a quantile in a sliding square neighbourhood across ensemble members.

Parameters
input3D vector of values (Y, X, E)
quantileQuantile to compute (between 0 and 1)
halfwidthFilter halfwidth in number of gridpoints
Returns
2D vector of final field (Y, x)

◆ neighbourhood_quantile_ens()

vec2 gridpp::neighbourhood_quantile_ens ( const vec3 input,
float  quantile,
int  halfwidth 
)

Deprecated: Compute neighbourhood quantiles on ensemble field.

Deprecated:
Use neighbourhood_quantile() function

◆ neighbourhood_quantile_ens_fast()

vec2 gridpp::neighbourhood_quantile_ens_fast ( const vec3 input,
float  quantile,
int  radius,
const vec thresholds 
)

Deprecated: Compute neighbourhood quantiles fast on ensemble field.

Deprecated:
Use neighbourhood_quantile_fast() function

◆ neighbourhood_quantile_fast() [1/4]

vec2 gridpp::neighbourhood_quantile_fast ( const vec2 input,
float  quantile,
int  halfwidth,
const vec thresholds 
)

Fast and approximate neighbourhood quantile.

Parameters
input2D vector of values (Y, X)
quantileQuantile to compute (between 0 and 1)
halfwidthFilter halfwidth in number of gridpoints
thresholds1D vector of thresholds to use to approximate value
Returns
2D vector of final field (Y, X)

◆ neighbourhood_quantile_fast() [2/4]

vec2 gridpp::neighbourhood_quantile_fast ( const vec3 input,
float  quantile,
int  halfwidth,
const vec thresholds 
)

Fast and approximate neighbourhood quantile for ensemble of fields.

Parameters
input3D vector of values (Y, X, E)
quantileQuantile to compute (between 0 and 1)
halfwidthFilter halfwidth in number of gridpoints
thresholds1D vector of thresholds to use to approximate value
Returns
2D vector of final field (Y, X)

◆ neighbourhood_quantile_fast() [3/4]

vec2 gridpp::neighbourhood_quantile_fast ( const vec2 input,
const vec2 quantile,
int  halfwidth,
const vec thresholds 
)

Fast and approximate neighbourhood quantile, with spatially varying quantile levels.

Parameters
input2D vector of values (Y, X)
quantile2D vector quantile levels to compute (between 0 and 1) (Y, X)
halfwidthFilter halfwidth in number of gridpoints
thresholds1D vector of thresholds to use to approximate value
Returns
2D vector of final field (Y, X)

◆ neighbourhood_quantile_fast() [4/4]

vec2 gridpp::neighbourhood_quantile_fast ( const vec3 input,
const vec2 quantile,
int  halfwidth,
const vec thresholds 
)

Fast and approximate neighbourhood quantile for ensemble of fields, with spatially varying quantile.

Parameters
input3D vector of values with dimensions (Y, X, E)
quantile2D vector of quantile levels to compute (between 0 and 1) (Y, X)
halfwidthFilter halfwidth in number of gridpoints
thresholds1D vector of thresholds to use to approximate value
Returns
2D vector of final field (Y, X)

◆ neighbourhood_score()

vec2 gridpp::neighbourhood_score ( const Grid grid,
const Points points,
const vec2 fcst,
const vec ref,
int  half_width,
gridpp::Metric  metric,
float  threshold 
)

Compute a score for a metric of all points within a radius.

This is an experimental method.

◆ neighbourhood_search()

vec2 gridpp::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 that match a search criteria.

If values are found within the range, then the average of the values are used. If not, the closest point to the range is used, provided this is within search_delta of the criteria range. if no point fulfills this, the original value is not modified.

Parameters
input2D vector of input values, e.g. temperatures (Y, X)
search_array2D vector to search in, e.g. land area fraction (Y, X)
halfwidthNeighbourhood halfwidth to search in
search_target_minLower bound of criteria, e.g. a land area fraction value
search_target_MaxUpper bound of criteria, e.g. a land area fraction value
search_deltaMax distance outside target range where the most suitable point is still used
apply_array2D integer vector of 1 (correct this grid point) or 0 (don't correct)
Returns
2D vector of output values (Y, X)

◆ num_missing_values()

int gridpp::num_missing_values ( const vec2 iArray)

Count the number of missing values in a vector.

Parameters
array2D vector
Returns
Number of invalid values

◆ optimal_interpolation() [1/2]

vec2 gridpp::optimal_interpolation ( const Grid bgrid,
const vec2 background,
const Points obs_points,
const vec obs,
const vec variance_ratios,
const vec background_at_points,
const StructureFunction structure,
int  max_points,
bool  allow_extrapolation = true 
)

Optimal interpolation for a deterministic gridded field.

Parameters
bgridGrid of background field
background2D field of background values
obs_pointsPoints of observations
obs1D vector of observations
variance_ratios1D vector of ratio of observation to background error variance
background_at_points1D vector of background at observation points
structureStructure function
max_pointsMaximum number of observations to use inside localization zone; Use 0 to disable
allow_extrapolationAllow OI to extrapolate increments outside increments at observations points
Returns
2D vector of analised values

◆ optimal_interpolation() [2/2]

vec gridpp::optimal_interpolation ( const Points bpoints,
const vec background,
const Points obs_points,
const vec obs,
const vec variance_ratios,
const vec background_at_points,
const StructureFunction structure,
int  max_points,
bool  allow_extrapolation = true 
)

Optimal interpolation for a deterministic vector of points.

Parameters
bpointsPoints of background field
background1D field of background values
obs_pointsPoints of observations
obs1D vector of observations
variance_ratios1D vector of ratio of observation to background error variance
background_at_points1D vector of background at observation points
structureStructure function
max_pointsMaximum number of observations to use inside localization zone; Use 0 to disable
allow_extrapolationAllow OI to extrapolate increments outside increments at observations points
Returns
1D vector of analised values

◆ optimal_interpolation_ensi() [1/2]

vec3 gridpp::optimal_interpolation_ensi ( const Grid bgrid,
const vec3 background,
const Points obs_points,
const vec obs,
const vec obs_standard_deviations,
const vec2 background_at_points,
const StructureFunction structure,
int  max_points,
bool  allow_extrapolation = true 
)

Optimal interpolation for an ensemble gridded field See Lussana et al 2019 (DOI: 10.1002/qj.3646)

Parameters
bgridGrid of background field
background3D vector of background values (Y, X, E)
obs_pointsobservation points
obs1D vector of observations
obs_standard_deviations1D vector of observation standard deviations
background_at_points2D vector of background at observation points (L, E)
structureStructure function
max_pointsMaximum number of observations to use inside localization zone; Use 0 to disable
allow_extrapolationAllow OI to extrapolate increments outside increments at observations
Returns
3D vector of analised values (Y, X, E)

◆ optimal_interpolation_ensi() [2/2]

vec2 gridpp::optimal_interpolation_ensi ( const Points bpoints,
const vec2 background,
const Points obs_points,
const vec obs,
const vec obs_standard_deviations,
const vec2 background_at_points,
const StructureFunction structure,
int  max_points,
bool  allow_extrapolation = true 
)

Optimal interpolation for an ensemble point field See Lussana et al 2019 (DOI: 10.1002/qj.3646)

Parameters
bpointsPoints of background field
background2D vector of background values (L, E)
obs_pointsObservation points
obs1D vector of observations
obs_standard_deviations1D vector of observation standard deviations
background_at_points2D vector of background at observation points (L, E)
structureStructure function
max_pointsMaximum number of observations to use inside localization zone; Use 0 to disable
allow_extrapolationAllow OI to extrapolate increments outside increments at observations
Returns
2D vector of analised values (L, E)

◆ optimal_interpolation_full() [1/2]

vec2 gridpp::optimal_interpolation_full ( const Grid bgrid,
const vec2 background,
const vec2 bvariance,
const Points obs_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.

Parameters
bgridGrid of background field
background2D vector of background values
bvariance2D vector of background variances
obs_pointsPoints of observations
obs1D vector of observations
obs_variance1D vector of observation variances
background_at_points1D vector of background at observation points
bvariance_at_points1D vector of background variance at observation points
structureStructure function
max_pointsMaximum number of observations to use inside localization zone; Use 0 to disable
analysis_variance2D output vector of analysis variance
allow_extrapolationAllow OI to extrapolate increments outside increments at observations
Returns
2D vector of analysed values

◆ optimal_interpolation_full() [2/2]

vec gridpp::optimal_interpolation_full ( const Points bpoints,
const vec background,
const vec bvariance,
const Points obs_points,
const vec obs,
const vec obs_variance,
const vec background_at_points,
const vec bvariance_at_points,
const StructureFunction structure,
int  max_points,
vec analysis_variance,
bool  allow_extrapolation = true 
)

Optimal interpolation for a deterministic vector of points including analysis variance.

Parameters
bpointsPoints of background field
background1D vector of background values
bvariance1D vector of background variance
obs_pointsPoints of observations
obs1D vector of observations
obs_variance1D vector of observation variances
background_at_points1D vector of background at observation points
bvariance_at_points1D vector of background variance at observation points
structureStructure function
max_pointsMaximum number of observations to use inside localization zone; Use 0 to disable
analysis_variance1D output vector of analysis variance
allow_extrapolationAllow OI to extrapolate increments outside increments at observations
Returns
2D vector of analised values

◆ point_in_rectangle()

bool gridpp::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.

The 4 points must be provided in an order that draws out a rectangle (either clockwise or counter-clockwise).

Parameters
AA point in the rectangle
BA point in the rectangle
CA point in the rectangle
DA point in the rectangle
mThe point to test if it is inside
Returns
True if the point is inside, False otherwise

◆ pressure() [1/2]

float gridpp::pressure ( float  ielev,
float  oelev,
float  ipressure,
float  itemperature = 288.15 
)

Calculate pressure at a new elevation.

Parameters
ielevElevation at start point [m]
oelevElevation at new point [m]
ipressurePressure at start point [Pa]
itemperatureTemperature at start point [K]
Returns
Pressure at new point [Pa]

◆ pressure() [2/2]

vec gridpp::pressure ( const vec ielev,
const vec oelev,
const vec ipressure,
const vec itemperature 
)

Calculate Vector version of pressure calculation.

Parameters
ielev1D vector of elevation at start point [m]
oelev1D vector of elevation at new point [m]
ipressure1D vector of pressure at start point [Pa]
itemperature1D vector of temperature at start point [K]
Returns
1D vector of pressure at new points [Pa]

◆ qnh() [1/2]

float gridpp::qnh ( float  pressure,
float  altitude 
)

Diagnose QNH from pressure and altitude.

Parameters
pressurePressure at point [Pa]
altitudeAltitude of point [m]
Returns
QNH [Pa]

◆ qnh() [2/2]

vec gridpp::qnh ( const vec pressure,
const vec altitude 
)

Vector version of QNH calculation.

Parameters
pressure1D vector of Pressure [Pa]
altitude1D vector of altitude [m]
Returns
1D vector of QNH [Pa]

◆ quantile_mapping_curve()

vec gridpp::quantile_mapping_curve ( const vec ref,
const vec fcst,
vec output_fcst,
vec  quantiles = vec() 
)

Create quantile mapping calibration curve.

Parameters
ref1D vector of reference values (observations)
fcst1D vector of forecast values
output_fcst1D vector of output forecast quantiles
quantiles1D vector of quantile levels to extract. If empty, use all values.
Returns
1D vector of output reference quantiles

◆ relative_humidity() [1/2]

float gridpp::relative_humidity ( float  temperature,
float  dewpoint 
)

Calculate relative humidity from temperature and dewpoint temperature.

Parameters
temperatureTemperature [K]
dewpointDewpoint temperature [K]
Returns
Relative humidity [1]

◆ relative_humidity() [2/2]

vec gridpp::relative_humidity ( const vec temperature,
const vec dewpoint 
)

Vector version of relative humidity calculation.

Parameters
temperature1D vector of temperature [K]
dewpoint1D vector of dewpoint temperature [K]
Returns
1D vector of relative humidity [1]

◆ sea_level_pressure() [1/2]

float gridpp::sea_level_pressure ( float  ps,
float  altitude,
float  temperature,
float  rh = gridpp::MV,
float  dewpoint = gridpp::MV 
)

Convert Surface Pressure to Sea Level Pressure.

Parameters
psSurface pressure [Pa]
altitudeStation altitude above sea level [m]
temperature2m temperature [K]
rh2m Relative humidity [1]
dewpoint2m Dewpoint Temperature at station [K]
Returns
Sea Level Pressure [Pa]

◆ sea_level_pressure() [2/2]

vec gridpp::sea_level_pressure ( const vec ps,
const vec altitude,
const vec temperature,
const vec rh,
const vec dewpoint 
)

Vector version of Convert Surface Pressure to Sea Level Pressure.

Parameters
ps1D vector of Surface pressures [Pa]
altitude1D vector of station altitudes above sea level [m]
temperature1D vector of temperature [K]
rh1D vector of relative humidity [1]
dewpoint1D vector of dewpoint temperature [K]
Returns
1D vector of sea Level pressure [Pa]

◆ set_debug_level()

void gridpp::set_debug_level ( int  level)

Set the verbosity of debug messages.

Use 0 for no messages.s

Parameters
levelDebug level

◆ set_omp_threads()

void gridpp::set_omp_threads ( int  num)

Set the number of OpenMP threads to use.

Overrides OMP_NUM_THREAD env variable.

Parameters
numNumber of threads

◆ simple_gradient() [1/4]

vec2 gridpp::simple_gradient ( const Grid igrid,
const Grid ogrid,
const vec2 ivalues,
float  elev_gradient,
Downscaler  downscaler = Nearest 
)

Elevation correction dowscaling grid to grid for a deterministic field.

Parameters
igridGrid of input (Yi, Xi)
ogridGrid of downscaled output (Yo, Xo)
ivalues2D vector of input values (Yi, Xi)
elev_gradientElevation gradient [input unit/m]
downscalerDownscaling method
Returns
2D vector of output values (Yo, Xo)

◆ simple_gradient() [2/4]

vec3 gridpp::simple_gradient ( const Grid igrid,
const Grid ogrid,
const vec3 ivalues,
float  elev_gradient,
Downscaler  downscaler = Nearest 
)

Elevation correction dowscaling grid to grid for an ensemble field.

Parameters
igridGrid of input (Yi, Xi)
ogridGrid of downscaled output (Yo, Xo)
ivalues3D vector of input values (E, Yi, Xi)
elev_gradientElevation gradient [input unit/m]
downscalerDownscaling method
Returns
3D vector of output values (E, Yo, Xo)

◆ simple_gradient() [3/4]

vec gridpp::simple_gradient ( const Grid igrid,
const Points opoints,
const vec2 ivalues,
float  elev_gradient,
Downscaler  downscaler = Nearest 
)

Elevation correction dowscaling grid to point for a deterministic field.

Parameters
igridGrid of input (Yi, Xi)
opointsPoints of downscaled output
ivalues2D vector of input values grid (Yi, Xi)
elev_gradientElevation gradient [input unit/m]
downscalerDownscaling method
Returns
1D vector of output values

◆ simple_gradient() [4/4]

vec2 gridpp::simple_gradient ( const Grid igrid,
const Points opoints,
const vec3 ivalues,
float  elev_gradient,
Downscaler  downscaler = Nearest 
)

Elevation correction dowscaling grid to point for an ensemble field.

Parameters
igridGrid of input (Yi, Xi)
opointsPoints of downscaled output
ivalues3D vector of input values (E, Yi, Xi)
elev_gradientElevation gradient [input unit/m]
downscalerDownscaling method
Returns
2D vector of output values (E, Points)

◆ smart()

vec2 gridpp::smart ( const Grid igrid,
const Grid ogrid,
const vec2 ivalues,
int  num,
const StructureFunction structure 
)

Smart neighbour downscaling grid to grid for a deterministic field.

Parameters
igridGrid of input (Yi, Xi)
ogridGrid of downscaled output (Yo, Xo)
ivalues2D vector of input values (Yi, Xi)
numNumber of neighbours to average
structureStructure function for determining how similar neighbour gridpoints are
Returns
2D vector of output values (Yo, Xo)

◆ test_array()

float * gridpp::test_array ( float *  v,
int  n 
)

Special function whose presense is needed for SWIG.

◆ test_ivec2_output()

ivec2 gridpp::test_ivec2_output ( )

◆ test_ivec3_output()

ivec3 gridpp::test_ivec3_output ( )

◆ test_ivec_input()

int gridpp::test_ivec_input ( const ivec input)

Testing function for 1D input vector.

◆ test_ivec_output()

ivec gridpp::test_ivec_output ( )

◆ test_not_implemented_exception()

void gridpp::test_not_implemented_exception ( )

◆ test_vec2_argout()

float gridpp::test_vec2_argout ( vec2 distances)

Testing function for 2D vector treated as output.

◆ test_vec2_input()

float gridpp::test_vec2_input ( const vec2 input)

Testing function for 2D input vector.

◆ test_vec2_output()

vec2 gridpp::test_vec2_output ( )

Testing function for 2D output vector.

◆ test_vec3_input()

float gridpp::test_vec3_input ( const vec3 input)

Testing function for 3D input vector.

◆ test_vec3_output()

vec3 gridpp::test_vec3_output ( )

Testing function for 3D output vector.

◆ test_vec_argout()

float gridpp::test_vec_argout ( vec distances)

Testing function for 1D vector treated as output.

◆ test_vec_input()

float gridpp::test_vec_input ( const vec input)

Testing function for 1D input vector.

◆ test_vec_output()

vec gridpp::test_vec_output ( )

Testing function for 1D output vector.

◆ version()

std::string gridpp::version ( )

The gridpp version.

Returns
The gridpp version

◆ warning()

void gridpp::warning ( std::string  string)

Writes a warning message to standard out.

Parameters
stringWrite this string

◆ wetbulb() [1/2]

float gridpp::wetbulb ( float  temperature,
float  pressure,
float  relative_humidity 
)

Calculate wetbulb temperature from temperature, pressure, and relative humidity.

Parameters
temperatureTemperature [K]
pressureAir pressure [Pa]
relative_humidityRelative humidity [1]
Returns
Wetbulb temperature [K]

◆ wetbulb() [2/2]

vec gridpp::wetbulb ( const vec temperature,
const vec pressure,
const vec relative_humidity 
)

Vector version of wetbulb calculation.

Parameters
temperature1D vector of temperature [K]
pressure1D vector of air pressure [pa]
relative_humidity1D vector of Relative humidity [1]
Returns
Wetbulb temperatures [K]

◆ wind_direction() [1/2]

float gridpp::wind_direction ( float  xwind,
float  ywind 
)

Diagnose wind direction from its components.

If both xwind and ywind are 0, then direction is 180

Parameters
xwindX-component of wind [any unit]
ywindY-component of wind [any unit]
Returns
Wind direction [degrees]

◆ wind_direction() [2/2]

vec gridpp::wind_direction ( const vec xwind,
const vec ywind 
)

Vector version of wind direction calculation.

Parameters
xwind1D vector of X-component of wind [any unit]
ywind1D vector of Y-component of wind [any unit]
Returns
1D vector of wind direction [degrees]

◆ wind_speed() [1/2]

float gridpp::wind_speed ( float  xwind,
float  ywind 
)

Diagnose wind speed from its components.

Parameters
xwindX-component of wind [any unit]
ywindY-component of wind [any unit]
Returns
Wind speed [any unit]

◆ wind_speed() [2/2]

vec gridpp::wind_speed ( const vec xwind,
const vec ywind 
)

Vector version of wind speed calculation.

Parameters
xwind1D vector of X-component of wind [any unit]
ywind1D vector of Y-component of wind [any unit]
Returns
1D vector of wind speeds [any unit]

◆ window()

vec2 gridpp::window ( const vec2 array,
int  length,
gridpp::Statistic  statistic,
bool  before = false,
bool  keep_missing = false,
bool  missing_edges = true 
)

Compute window statistics across time, independently for each case.

Parameters
array2D vector of input values (Case, Time)
lengthWindow length in number of timesteps
statisticStatistic to apply to window
beforeIf true, make the window end at the particular time. If false, centre it.
keep_missingIf true, window value will be missing if one or more values in window are missing
missing_edgesIf true put missing values at the edges, where window overshoots the edge
Returns
2D vector of output values (Case, Time)

Variable Documentation

◆ _debug_level

int gridpp::_debug_level = 0
static

◆ gas_constant_mol

const float gridpp::gas_constant_mol = 8.31447
static

Universal Gas Constant [kg*m^2*s^-2/(K*mol)].

◆ gas_constant_si

const float gridpp::gas_constant_si = 287.05
static

Universal Gas Constant [J/(kg*K)].

◆ gravit

const float gridpp::gravit = 9.80665
static

Gravitational acceleration [m/s^2].

◆ lapse_rate

const float gridpp::lapse_rate =0.0065
static

Constant Lapse Rate moist air standard atmosphere [K/m].

◆ molar_mass

const float gridpp::molar_mass = 0.0289644
static

Molar Mass of Dry Air [kg/mol].

◆ MV

const float gridpp::MV = NAN
static

Missing value indicator.

◆ MV_CML

const float gridpp::MV_CML = -999
static

Missing value indicator in gridpp command-line tool.

◆ pi

const float gridpp::pi = 3.14159265
static

Mathematical constant pi.

◆ radius_earth

const double gridpp::radius_earth = 6.378137e6
static

Radius of the earth [m].

◆ standard_surface_temperature

const float gridpp::standard_surface_temperature = 288.15
static

Temperature at surface in standard atmosphere [K].

◆ swig_default_value

const float gridpp::swig_default_value = -1
static

Default value used to fill array in SWIG testing functions.

Not useful for any other purpose.