Gridpp 0.5.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
 
typedef std::vector< vecvec2
 
typedef std::vector< vec2vec3
 
typedef std::vector< int > ivec
 
typedef std::vector< ivecivec2
 
typedef std::vector< ivec2ivec3
 
typedef std::vector< double > dvec
 
typedef std::vector< dvecdvec2
 

Functions

Data assimilation methods

Functions that merge observations with a background field

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. More...
 
vec optimal_interpolation (const Points &bpoints, const vec &background, const Points &points, const vec &pobs, const vec &pratios, const vec &pbackground, const StructureFunction &structure, int max_points)
 Optimal interpolation for a deterministic vector of points. More...
 
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-transformed case. More...
 
vec optimal_interpolation_transform (const Points &bpoints, const vec &background, float bsigma, const Points &points, const vec &pobs, const vec &psigmas, const vec &pbackground, const StructureFunction &structure, int max_points, const Transform &transform)
 
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.1002/qj.3646) More...
 
vec2 optimal_interpolation_ensi (const Points &bpoints, const vec2 &background, const Points &points, const vec &pobs, const vec &psigmas, const vec2 &pbackground, const StructureFunction &structure, int max_points)
 
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. 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 for an ensemble of fields. 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. 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 neighbourhood filter without any shortcuts. More...
 
vec get_neighbourhood_thresholds (const vec2 &input, int num_thresholds)
 Calculate appropriate approximation thresholds for 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 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

vec2 quantile_mapping_curve (const vec &ref, const vec &fcst, vec quantiles=vec())
 Create quantile mapping calibration curve. More...
 
vec2 metric_optimizer_curve (const vec &ref, const vec &fcst, const vec &thresholds, Metric metric)
 Create calibration curve that optimizes a metric. More...
 
vec apply_curve (const vec &fcst, const vec2 &curve, Extrapolation policy_below, Extrapolation policy_above)
 Apply arbitrary calibration curve to 1D forecasts. More...
 
vec2 apply_curve (const vec2 &fcst, const vec2 &curve, Extrapolation policy_below, Extrapolation policy_above)
 Apply arbitrary calibration curve to 2D forecasts. More...
 
vec2 monotonize_curve (vec2 curve)
 Ensure calibration curve is monotonic, by removing points. More...
 
float get_optimal_threshold (const vec &ref, const vec &fcst, float threshold, Metric metric)
 
float calc_score (float a, float b, float c, float d, Metric metric)
 
float calc_score (const vec &ref, const vec &fcst, float threshold, Metric metric)
 
float calc_score (const vec &ref, const vec &fcst, float threshold, float fthreshold, Metric metric)
 
vec2 correction (const Grid &rgrid, const vec2 &rvalues, const Points &npoints, const vec &nvalues, float mean_radius, float outer_radius, float inner_radius, int min_num, int max_num, CorrectionType type, ivec2 &count)
 
vec2 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)
 
Downscaling methods

Functions that interpolate data from one grid to another

vec2 nearest (const Grid &igrid, const Grid &ogrid, const vec2 ivalues)
 Nearest neighbour dowscaling grid to grid. More...
 
vec3 nearest (const Grid &igrid, const Grid &ogrid, const vec3 ivalues)
 
vec nearest (const Grid &igrid, const Points &opoints, const vec2 ivalues)
 Nearest neighbour dowscaling grid to point. More...
 
vec2 nearest (const Grid &igrid, const Points &opoints, const vec3 ivalues)
 
vec2 bilinear (const Grid &igrid, const Grid &ogrid, const vec2 ivalues)
 Bilinear downscaling grid to grid. More...
 
vec bilinear (const Grid &igrid, const Points &opoints, const vec2 ivalues)
 Bilinear downscaling grid to points. More...
 
vec2 simple_gradient (const Grid &igrid, const Grid &ogrid, const vec2 ivalues, float elev_gradient)
 
vec simple_gradient (const Grid &igrid, const Points &opoints, const vec2 ivalues, float elev_gradient)
 
vec2 smart (const Grid &igrid, const Grid &ogrid, const vec2 &ivalues, int num, const StructureFunction &structure)
 Smart neighbour downscaling grid to grid. 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 Points &points, const Grid &grid, float radius)
 For each gridpoint, 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...
 
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...
 
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 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...
 
void initialize_omp ()
 Sets the number of OpenMP threads to 1 if OMP_NUM_THREADS undefined. More...
 
Utilities

Helper functions

Statistic get_statistic (std::string name)
 Convert name of a statistic enum. More...
 
std::string version ()
 The gridpp version. More...
 
double clock ()
 
void debug (std::string string)
 
void warning (std::string string)
 
void error (std::string string)
 
void future_deprecation_warning (std::string function, std::string other="")
 
bool is_valid (float value)
 
float calc_statistic (const vec &array, Statistic statistic)
 
float calc_quantile (const vec &array, float quantile)
 
vec calc_statistic (const vec2 &array, Statistic statistic)
 
vec calc_quantile (const vec2 &array, float quantile=MV)
 
int num_missing_values (const vec2 &iArray)
 
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. More...
 
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. More...
 
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. More...
 
ivec2 init_ivec2 (int Y, int X, int value)
 Initialize a vector of size Y, X, with a given value. More...
 
vec2 init_vec2 (int Y, int X, float value=MV)
 
ivec3 init_ivec3 (int Y, int X, int E, int value)
 Initialize a vector of size Y, X, E, with a given value. More...
 
vec3 init_vec3 (int Y, int X, int E, float value=MV)
 
vec calc_even_quantiles (const vec &values, int num)
 Get reasonably spaced quantiles from a vector of values, ignoring duplicate values but including the first number after duplicated values. More...
 
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. 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)
 
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...
 

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...
 

Enumerations

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  Extrapolation { OneToOne = 0, MeanSlope = 10, NearestSlope = 20, Zero = 30 }
 Methods for extrapolating outside a curve. 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,
  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
 
class  CressmanStructure
 Simple structure function based on distance, elevation, and land area fraction. More...
 
class  CrossValidation
 
class  Grid
 Represents a 2D grid of locations and their metadata. More...
 
class  Identity
 
class  KDTree
 Helper class for Grid and Points. More...
 
class  Log
 
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
 

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

◆ dvec2

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

◆ ivec

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

◆ ivec2

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

◆ ivec3

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

◆ vec

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

◆ vec2

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

◆ vec3

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

Enumeration Type Documentation

◆ 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.

◆ 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.

◆ 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.

Unknown 

Unknown statistic.

Function Documentation

◆ apply_curve() [1/2]

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

Apply arbitrary calibration curve to 1D forecasts.

Parameters
fcst1D vector of forecast values
curveCalibration curve
policy_belowExtrapolation policy below curve
policy_aboveExtrapolation policy above curve
Returns
Calibrated forecasts

◆ apply_curve() [2/2]

vec2 gridpp::apply_curve ( const vec2 fcst,
const vec2 curve,
gridpp::Extrapolation  policy_below,
gridpp::Extrapolation  policy_above 
)

Apply arbitrary calibration curve to 2D forecasts.

Parameters
fcst2D grid of forecast values
curveCalibration curve
policy_belowExtrapolation policy below curve
policy_aboveExtrapolation policy above curve
Returns
Calibrated forecasts

◆ bilinear() [1/2]

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

Bilinear downscaling grid to grid.

Parameters
igridInput grid
ogridOutput grid to downscale to
ivalues2D vector of values on the input grid
Returns
Values on the output grid

◆ bilinear() [2/2]

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

Bilinear downscaling grid to points.

Parameters
igridInput grid
ogridOutput points to downscale to
ivalues2D vector of values on the input grid
Returns
Values for the output points

◆ calc_even_quantiles()

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

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

Include the lowest and highest values.

Parameters
valuesvector of values (unsorted, and no invalid values)
numnumber of thresholds to get

◆ calc_gradient()

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

Parameters
gridGrid
baseDependent variable. Missing values are not used.
valuesIndependent variable. Missing values are not used.
radiusNeighbourhood radius 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 is not met

◆ calc_quantile() [1/2]

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

◆ calc_quantile() [2/2]

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

◆ calc_score() [1/3]

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

◆ calc_score() [2/3]

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

◆ calc_score() [3/3]

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

◆ calc_statistic() [1/2]

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

◆ calc_statistic() [2/2]

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

◆ clock()

double gridpp::clock ( )

◆ compatible_size() [1/2]

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

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

◆ compatible_size() [2/2]

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

◆ correction() [1/2]

vec2 gridpp::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,
gridpp::CorrectionType  type,
ivec2 count 
)

◆ correction() [2/2]

vec2 gridpp::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,
gridpp::CorrectionType  type,
ivec2 count 
)

◆ count() [1/2]

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

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

Parameters
gridGrid
pointsPoints
radiusRadius [m]
Returns
Number of gridpoints

◆ count() [2/2]

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

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

Parameters
gridGrid
pointsPoints
radiusRadius [m]
Returns
Number of points

◆ debug()

void gridpp::debug ( std::string  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
temperatureTemperatures [K]
relative_humidityRelative humidities [1]
Returns
Dewpoint temperatures [K]

◆ distance() [1/3]

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
Distance [m] to nearest gridpoint for each point

◆ distance() [2/3]

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
ogridOutput grid
numNumber of points
Returns
Distance [m] to nearest gridpoint for each gridpoint

◆ distance() [3/3]

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
Distance [m] to nearest gridpoint for each point

◆ error()

void gridpp::error ( std::string  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.

Parameters
inputDeterministic values with dimensions Y, X
radiiCircle radii for each point
valueFill in this value
outsideif True, fill outside circles, if False, fill inside circles

◆ fill_missing()

vec2 gridpp::fill_missing ( const vec2 values)

Fill in missing values based on nearby values.

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

◆ future_deprecation_warning()

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

◆ get_lower_index()

int gridpp::get_lower_index ( float  iX,
const std::vector< float > &  iValues 
)

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

Parameters
iXLookup value
iValuesLookup vector. Must be sorted.
Returns
The index into iValues that is equal 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 neighbourhood quantile.

Parameters
input2D (Y, X) array of values
num_thresholdsNumber of 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 (Y, X, T) array of values
num_thresholdsNumber of thresholds

◆ get_optimal_threshold()

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

◆ get_statistic()

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

Convert name of a statistic enum.

◆ get_upper_index()

int gridpp::get_upper_index ( float  iX,
const std::vector< float > &  iValues 
)

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

Parameters
iXLookup value
iValuesLookup vector. Must be sorted.
Returns
The index into iValues that is equal or just above iX

◆ gridding()

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 MV where there are not enough observations

Parameters
gridGrid to aggregate to
pointsPoints with values
valuesValues at points
radiusCircle radius for aggregate points [m]
min_numMinimum number of points in radius to create a value
statisticStatistic to compute on points within radius

◆ init_ivec2()

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

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

◆ init_ivec3()

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

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

◆ init_vec2()

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

◆ init_vec3()

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

◆ initialize_omp()

void gridpp::initialize_omp ( )

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

◆ interpolate()

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

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

◆ is_valid()

bool gridpp::is_valid ( float  value)

◆ metric_optimizer_curve()

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

Create calibration curve that optimizes a metric.

Parameters
refReference values (observations)
fcstForecast values
thresholdsThresholds for computing optimal values for
metricA Metric to optimize for
Returns
Calibration curve

◆ monotonize_curve()

vec2 gridpp::monotonize_curve ( vec2  curve)

Ensure calibration curve is monotonic, by removing points.

Parameters
curveCalibration curve
Returns
Monotonic calibration curve

◆ nearest() [1/4]

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

Nearest neighbour dowscaling grid to grid.

Parameters
igridInput grid
ogridOutput grid to downscale to
ivalues2D vector of values on the input grid
Returns
Values on the output grid

◆ nearest() [2/4]

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

◆ nearest() [3/4]

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

Nearest neighbour dowscaling grid to point.

Parameters
igridInput grid
ogridOutput points to downscale to
ivalues2D vector of values on the input grid
Returns
Values for the output points

◆ nearest() [4/4]

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

◆ 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 grid of values
halfwidthFilter halfwidth in number of gridpoints
statisticStatistic to compute

◆ neighbourhood() [2/2]

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

Spatial neighbourhood filter for an ensemble of fields.

Parameters
input3D grid of values with dimensions (Y, X, E)
halfwidthFilter halfwidth in number of gridpoints
statisticStatistic to compute

◆ 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 grid of values
halfwidthFilter halfwidth in number of gridpoints
operationone of min, mean, median, max

◆ neighbourhood_brute_force() [2/2]

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

Spatial neighbourhood filter without any shortcuts.

This is quite slow and is only useful for testing.

Parameters
input3D grid of values with dimensions (Y, X, E)
halfwidthFilter halfwidth in number of gridpoints
operationone of min, mean, median, max

◆ 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 grid of values
quantileQuantile to compute (between 0 and 1)
halfwidthFilter halfwidth in number of gridpoints

◆ neighbourhood_quantile() [2/2]

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

Computes a quantile in a sliding square neighbourhood for an ensemble of fields.

Parameters
input3D grid of values with dimensions (Y, X, E)
quantileQuantile to compute (between 0 and 1)
halfwidthFilter halfwidth in number of gridpoints

◆ 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 grid of values
quantileQuantile to compute (between 0 and 1)
halfwidthFilter halfwidth in number of gridpoints
thresholdsVector of thresholds to use to approximate value

◆ 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 grid of values with dimensions (Y, X, E)
quantileQuantile to compute (between 0 and 1)
halfwidthFilter halfwidth in number of gridpoints
thresholdsVector of thresholds to use to approximate value

◆ 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.

Parameters
input2D grid of values
quantile2D grid quantiles to compute (between 0 and 1)
halfwidthFilter halfwidth in number of gridpoints
thresholdsVector of thresholds to use to approximate value

◆ 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 grid of values with dimensions (Y, X, E)
quantile2D grid quantiles to compute (between 0 and 1)
halfwidthFilter halfwidth in number of gridpoints
thresholdsVector of thresholds to use to approximate value

◆ num_missing_values()

int gridpp::num_missing_values ( const vec2 iArray)

◆ optimal_interpolation() [1/2]

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

Parameters
bgridGrid of background field
background2D field of background values
pointsPoints of observations
pobsVector of observations
pratiosVector of ratio of observation error variance to background variance
pbackgroundBackground with observation operator
structureStructure function
max_pointsMaximum number of observations to use inside localization zone; Use 0 to disable

◆ optimal_interpolation() [2/2]

vec gridpp::optimal_interpolation ( const Points bpoints,
const vec background,
const Points points,
const vec pobs,
const vec pratios,
const vec pbackground,
const StructureFunction structure,
int  max_points 
)

Optimal interpolation for a deterministic vector of points.

Parameters
bpointsPoints of background field
background1D field of background values
pointsPoints of observations
pobsVector of observations
pratiosVector of ratio of observation error variance to background variance
pbackgroundBackground with observation operator
structureStructure function
max_pointsMaximum number of observations to use inside localization zone; Use 0 to disable

◆ optimal_interpolation_ensi() [1/2]

vec3 gridpp::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.1002/qj.3646)

Parameters
input3D field of background values (Y, X, E)
bgridgrid corresponding to input
pobsvector of observations
pcivector of ci values
pointsobservation points

◆ optimal_interpolation_ensi() [2/2]

vec2 gridpp::optimal_interpolation_ensi ( const Points bpoints,
const vec2 background,
const Points points,
const vec pobs,
const vec psigmas,
const vec2 pbackground,
const StructureFunction structure,
int  max_points 
)

◆ optimal_interpolation_transform() [1/2]

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

◆ optimal_interpolation_transform() [2/2]

vec gridpp::optimal_interpolation_transform ( const Points bpoints,
const vec background,
float  bsigma,
const Points points,
const vec pobs,
const vec psigmas,
const vec pbackground,
const StructureFunction structure,
int  max_points,
const Transform transform 
)

◆ 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
oelevElevation at new point
ipressurePressure at start point
itemperatureTemperature at start point
Returns
Pressure at new point

◆ 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
ielevElevations at start point
oelevElevations at new point
ipressurePressures at start point
itemperatureTemperatures at start point
Returns
Pressures at new points

◆ 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
pressurePressures at points [pa]
altitudeAltitudes of points [m]
Returns
QNH [pa]

◆ quantile_mapping_curve()

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

Create quantile mapping calibration curve.

Parameters
refReference values (observations)
fcstForecast values
quantilesVector of quantiles to extract. If empty, use all values.
Returns
Calibration curve

◆ 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
temperatureTemperatures [K]
dewpointDewpoint temperatures [K]
Returns
Relative humidities [1]

◆ set_omp_threads()

void gridpp::set_omp_threads ( int  num)

Set the number of OpenMP threads to use.

Overwrides OMP_NUM_THREAD env variable.

◆ simple_gradient() [1/2]

vec2 gridpp::simple_gradient ( const Grid igrid,
const Grid ogrid,
const vec2  ivalues,
float  elev_gradient 
)

◆ simple_gradient() [2/2]

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

◆ smart()

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

Smart neighbour downscaling grid to grid.

Parameters
igridInput grid
ogridOutput points to downscale to
ivalues2D vector of values on the input grid
numNumber of neighbours to average
structureStructure function for determining similarity
Returns
Values for the output points

◆ 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)

◆ 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]
Relativehumidity [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
temperatureTemperatures [K]
pressureAir pressures [pa]
Relativehumidities [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
xwindX-components of wind [any unit]
ywindY-components of wind [any unit]
Returns
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
xwindX-components of wind [any unit]
ywindY-components of wind [any unit]
Returns
Wind speeds [any unit]

Variable Documentation

◆ 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].

◆ 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.