Typedefs | |
Short-hand notation for vectors of different dimensions sizes | |
typedef std::vector< float > | vec |
typedef std::vector< vec > | vec2 |
typedef std::vector< vec2 > | vec3 |
typedef std::vector< int > | ivec |
typedef std::vector< ivec > | ivec2 |
typedef std::vector< ivec2 > | ivec3 |
typedef std::vector< double > | dvec |
typedef std::vector< dvec > | dvec2 |
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, bool allow_extrapolation=true) |
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, 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 &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 &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_sigmas, 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 &points, const vec &pobs, const vec &psigmas, const vec2 &pbackground, const StructureFunction &structure, int max_points, bool allow_extrapolation=true) |
Optimal interpolation using a structure function based on an ensemble See Lussana et al 2019 (DOI: 10.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, bool allow_extrapolation=true) |
vec2 | local_distribution_correction (const Grid &bgrid, const vec2 &background, const Points &points, const vec &pobs, const vec &pbackground, const StructureFunction &structure, float min_quantile, float max_quantile, int min_points=0) |
Correction of a gridded field ensuring the distribution of values nearby match that of observations. More... | |
vec2 | local_distribution_correction (const Grid &bgrid, const vec2 &background, const Points &points, const vec2 &pobs, const vec2 &pbackground, const StructureFunction &structure, float min_quantile, float max_quantile, int min_points=0) |
Version with multiple number of timesteps. 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) |
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 | |
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. 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) |
Downscaling methods | |
Functions that interpolate data from one grid to another | |
vec2 | downscaling (const Grid &igrid, const Grid &ogrid, const vec2 &ivalues, Downscaler downscaler) |
Generic downscaler for simple methods that don't use information other than the grids. More... | |
vec3 | downscaling (const Grid &igrid, const Grid &ogrid, const vec3 &ivalues, Downscaler downscaler) |
vec | downscaling (const Grid &igrid, const Points &opoints, const vec2 &ivalues, Downscaler downscaler) |
vec2 | downscaling (const Grid &igrid, const Points &opoints, const vec3 &ivalues, Downscaler downscaler) |
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) |
vec | nearest (const Points &ipoints, const Points &opoints, const vec &ivalues) |
Nearest neighbour dowscaling point to point. More... | |
vec2 | nearest (const Points &ipoints, const Points &opoints, const vec2 &ivalues) |
vec2 | nearest (const Points &ipoints, const Grid &ogrid, const vec &ivalues) |
Nearest neighbour dowscaling point to grid. More... | |
vec3 | nearest (const Points &ipoints, const Grid &ogrid, const vec2 &ivalues) |
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) |
Nearest neighbour downscaling grid to grid and threshold and consensus in one. More... | |
vec2 | bilinear (const Grid &igrid, const Grid &ogrid, const vec2 &ivalues) |
Bilinear downscaling grid to grid. More... | |
vec3 | bilinear (const Grid &igrid, const Grid &ogrid, const vec3 &ivalues) |
vec | bilinear (const Grid &igrid, const Points &opoints, const vec2 &ivalues) |
Bilinear downscaling grid to points. More... | |
vec2 | bilinear (const Grid &igrid, const Points &opoints, const vec3 &ivalues) |
vec2 | simple_gradient (const Grid &igrid, const Grid &ogrid, const vec2 &ivalues, float elev_gradient, Downscaler downscaler=Nearest) |
vec3 | simple_gradient (const Grid &igrid, const Grid &ogrid, const vec3 &ivalues, float elev_gradient, Downscaler downscaler=Nearest) |
vec | simple_gradient (const Grid &igrid, const Points &opoints, const vec2 &ivalues, float elev_gradient, Downscaler downscaler=Nearest) |
vec2 | simple_gradient (const Grid &igrid, const Points &opoints, const vec3 &ivalues, float elev_gradient, Downscaler downscaler=Nearest) |
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) |
vec | full_gradient (const Grid &igrid, const Points &opoints, const vec2 &ivalues, const vec2 &elev_gradient, const vec2 &laf_gradient, Downscaler downscaler=Nearest) |
vec2 | full_gradient (const Grid &igrid, const Points &opoints, const vec3 &ivalues, const vec3 &elev_gradient, const vec3 &laf_gradient, Downscaler downscaler=Nearest) |
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. 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... | |
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... | |
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, 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 | Gamma |
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 | 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 messages. More... | |
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) |
vec2 | calc_quantile (const vec3 &array, const vec2 &quantile) |
Compute quantile with 2D varying quantile. More... | |
int | num_missing_values (const vec2 &iArray) |
int | get_lower_index (float iX, const vec &iValues) |
Find the index in a vector that is equal or just below a value. More... | |
int | get_upper_index (float iX, const vec &iValues) |
Find the index in a vector that is equal 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... | |
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 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 based on a search criteria. 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. 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 | compatible_size (const Points &points, const vec &v) |
bool | compatible_size (const Points &points, const vec2 &v) |
bool | compatible_size (const vec2 &a, const vec2 &b) |
bool | compatible_size (const vec2 &a, const vec3 &b) |
bool | compatible_size (const vec3 &a, const vec3 &b) |
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 std::vector<double> gridpp::dvec |
typedef std::vector<dvec> gridpp::dvec2 |
typedef std::vector<int> gridpp::ivec |
typedef std::vector<ivec> gridpp::ivec2 |
typedef std::vector<ivec2> gridpp::ivec3 |
typedef std::vector<float> gridpp::vec |
typedef std::vector<vec> gridpp::vec2 |
typedef std::vector<vec2> gridpp::vec3 |
enum gridpp::Downscaler |
Methods for extrapolating outside a curve.
enum gridpp::GradientType |
enum gridpp::Metric |
enum gridpp::Statistic |
Statistical operations to reduce a vector to a scalar.
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.
fcst | input forecast |
curve_ref | Reference quantiles |
curve_fcst | Forecast quantiles |
policy_below | Extrapolation policy below curve |
policy_above | Extrapolation policy above curve |
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.
fcst | 1D vector of forecast values |
curve_ref | Reference quantiles |
curve_fcst | Forecast quantiles |
policy_below | Extrapolation policy below curve |
policy_above | Extrapolation policy above curve |
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.
fcst | 2D grid of forecast values |
curve_ref | Reference quantiles |
curve_fcst | Forecast quantiles |
policy_below | Extrapolation policy below curve |
policy_above | Extrapolation policy above curve |
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.
fcst | 2D grid of forecast values |
curve_ref | Reference quantiles (Y, X, Q) |
curve_fcst | Forecast quantiles (Y, X, Q) |
policy_below | Extrapolation policy below curve |
policy_above | Extrapolation policy above curve |
Bilinear downscaling grid to grid.
igrid | Input grid |
ogrid | Output grid to downscale to |
ivalues | 2D vector of values on the input grid |
Bilinear downscaling grid to points.
igrid | Input grid |
ogrid | Output points to downscale to |
ivalues | 2D vector of values on the input grid |
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.
values | vector of values (unsorted, and no invalid values) |
num | number of thresholds to get |
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.
grid | Grid |
base | Dependent variable. Missing values are not used. |
values | Independent variable. Missing values are not used. |
radius | Neighbourhood radius in number of gridpoints |
min_nim | Minimum number of points required to compute gradient |
min_range | Minimum range of base to compute gradient |
default_gradient | Use this gradient if minimum number is not met |
float gridpp::calc_quantile | ( | const vec & | array, |
float | quantile | ||
) |
Compute quantile with 2D varying quantile.
array | Input array with dimensions (T, Y, X) |
quantile | Quantile array with dimensions (Y, X) |
float gridpp::calc_score | ( | float | a, |
float | b, | ||
float | c, | ||
float | d, | ||
gridpp::Metric | metric | ||
) |
float gridpp::calc_score | ( | const vec & | ref, |
const vec & | fcst, | ||
float | threshold, | ||
gridpp::Metric | metric | ||
) |
float gridpp::calc_score | ( | const vec & | ref, |
const vec & | fcst, | ||
float | threshold, | ||
float | fthreshold, | ||
gridpp::Metric | metric | ||
) |
float gridpp::calc_statistic | ( | const vec & | array, |
gridpp::Statistic | statistic | ||
) |
vec gridpp::calc_statistic | ( | const vec2 & | array, |
gridpp::Statistic | statistic | ||
) |
double gridpp::clock | ( | ) |
Check if the grid is the same size as the 2D vector.
If True, they are compatible, if false they are incompatible
For each gridpoint, counts the number of gridpoints within the radius.
igrid | Input grid |
ogrid | Output grid |
radius | Radius [m] |
For each point, counts the number of points within the radius.
ipoints | Input points |
opoints | Output points |
radius | Radius [m] |
void gridpp::debug | ( | std::string | string | ) |
float gridpp::dewpoint | ( | float | temperature, |
float | relative_humidity | ||
) |
Calculate dewpoint temperature from temperature and relative humidity.
temperature | Temperature [K] |
relative_humidity | Relative humidity [1] |
Vector version of dewpoint calculation.
temperature | Temperatures [K] |
relative_humidity | Relative humidities [1] |
For each output gridpoint, calculate the distance to nearest input gridpoint.
grid | Grid |
ogrid | Output grid |
num | Number of points |
For each output point, calculates the distance to nearest input point.
ipoints | Input points |
opoints | Output points |
num | Number of points |
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.
grid | Grid |
background | Deterministic values with dimensions Y, X |
points | Points representing observations |
observations | Vector of observations |
radii | Radius of circle where observations are inserted for each point [m] |
max_elev_diff | Only insert where elevation difference between grid and point is less than this value |
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.
grid | Grid |
background | Deterministic values with dimensions Y, X |
points | Points representing observations |
observations | Vector of observations |
halfwidths | Half width of square (in number of grid points) where observations are inserted for each point |
max_elev_diff | Only insert where elevation difference between grid and point is less than this value |
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.
igrid | Input grid |
ogrid | Output grid to downscale to |
ivalues | 3D vector of values on the input grid (Y, X, E) |
threshold | 2D vector of threshold values |
comparison_operator | lower than, lower or equal than, greater than, great or equal than |
vec2 gridpp::downscaling | ( | const Grid & | igrid, |
const Grid & | ogrid, | ||
const vec2 & | ivalues, | ||
Downscaler | downscaler | ||
) |
Generic downscaler for simple methods that don't use information other than the grids.
igrid | Input grid |
ogrid | Output grid to downscale to |
ivalues | 2D vector of values on the input grid |
downscaler | Downscaling method |
vec3 gridpp::downscaling | ( | const Grid & | igrid, |
const Grid & | ogrid, | ||
const vec3 & | ivalues, | ||
Downscaler | downscaler | ||
) |
vec gridpp::downscaling | ( | const Grid & | igrid, |
const Points & | opoints, | ||
const vec2 & | ivalues, | ||
Downscaler | downscaler | ||
) |
vec2 gridpp::downscaling | ( | const Grid & | igrid, |
const Points & | opoints, | ||
const vec3 & | ivalues, | ||
Downscaler | downscaler | ||
) |
void gridpp::error | ( | std::string | string | ) |
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)
input | Deterministic values with dimensions Y, X |
radii | Circle radii for each point |
value | Fill in this value |
outside | if True, fill outside circles, if False, fill inside circles |
Fill in missing values based on nearby values.
values | 2D array of values |
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.
igrid | input grid |
ogrid | output grid |
ivalues | values from igrid |
elev_gradient | elevation gradient |
laf_gradient | land area fraction gradient |
vec3 gridpp::full_gradient | ( | const Grid & | igrid, |
const Grid & | ogrid, | ||
const vec3 & | ivalues, | ||
const vec3 & | elev_gradient, | ||
const vec3 & | laf_gradient, | ||
Downscaler | downscaler = Nearest |
||
) |
vec gridpp::full_gradient | ( | const Grid & | igrid, |
const Points & | opoints, | ||
const vec2 & | ivalues, | ||
const vec2 & | elev_gradient, | ||
const vec2 & | laf_gradient, | ||
Downscaler | downscaler = Nearest |
||
) |
vec2 gridpp::full_gradient | ( | const Grid & | igrid, |
const Points & | opoints, | ||
const vec3 & | ivalues, | ||
const vec3 & | elev_gradient, | ||
const vec3 & | laf_gradient, | ||
Downscaler | downscaler = Nearest |
||
) |
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 |
||
) |
void gridpp::future_deprecation_warning | ( | std::string | function, |
std::string | other = "" |
||
) |
shape | Shape parameter of gamma distribution |
scale | Scale parameter of gamma distribution |
levels | Quantile levels to retrieve |
int gridpp::get_debug_level | ( | ) |
Get the currently set level of debug messages.
int gridpp::get_lower_index | ( | float | iX, |
const vec & | iValues | ||
) |
Find the index in a vector that is equal or just below a value.
iX | Lookup value |
iValues | Lookup vector. Must be sorted. |
Calculate appropriate approximation thresholds for neighbourhood quantile.
input | 2D (Y, X) array of values |
num_thresholds | Number of thresholds |
Calculate appropriate approximation thresholds for neighbourhood quantile based on an * ensemble.
input | 3D (Y, X, T) array of values |
num_thresholds | Number of thresholds |
int gridpp::get_omp_threads | ( | ) |
Get the number of OpenMP threads currently set.
float gridpp::get_optimal_threshold | ( | const vec & | ref, |
const vec & | fcst, | ||
float | threshold, | ||
gridpp::Metric | metric | ||
) |
gridpp::Statistic gridpp::get_statistic | ( | std::string | name | ) |
Convert name of a statistic enum.
int gridpp::get_upper_index | ( | float | iX, |
const vec & | iValues | ||
) |
Find the index in a vector that is equal or just above a value.
iX | Lookup value |
iValues | Lookup vector. Must be sorted. |
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.
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 MV where there are not enough observations.
ivec2 gridpp::init_ivec2 | ( | int | Y, |
int | X, | ||
int | value | ||
) |
Initialize a vector of size Y, X, with a given value.
ivec3 gridpp::init_ivec3 | ( | int | Y, |
int | X, | ||
int | E, | ||
int | value | ||
) |
Initialize a vector of size Y, X, E, with a given value.
void gridpp::initialize_omp | ( | ) |
Sets the number of OpenMP threads to 1 if OMP_NUM_THREADS undefined.
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.
x | Interpolation to this point |
iX | Vector of x-axis values. Vector must be sorted. |
iY | Vector of y-axis values corresponding to iX. |
bool gridpp::is_valid | ( | float | value | ) |
vec2 gridpp::local_distribution_correction | ( | const Grid & | bgrid, |
const vec2 & | background, | ||
const Points & | points, | ||
const vec & | pobs, | ||
const vec & | pbackground, | ||
const StructureFunction & | structure, | ||
float | min_quantile, | ||
float | max_quantile, | ||
int | min_points = 0 |
||
) |
Correction of a gridded field ensuring the distribution of values nearby match that of observations.
This is an experimental method.
bgrid | grid corresponding to input |
background | 2D field of background values (Y, X) |
points | observation points |
pobs | vector of observations |
pbackground | vector of background values at points |
structure | structure function specifying correlation between points |
min_quantile | truncate quantile map below this quantile |
max_quantile | truncate quantile map above this quantile |
max_points | maximum number of points used within localization radius (not used at moment) |
vec2 gridpp::local_distribution_correction | ( | const Grid & | bgrid, |
const vec2 & | background, | ||
const Points & | points, | ||
const vec2 & | pobs, | ||
const vec2 & | pbackground, | ||
const StructureFunction & | structure, | ||
float | min_quantile, | ||
float | max_quantile, | ||
int | min_points = 0 |
||
) |
Version with multiple number of timesteps.
Pool all observations across time in when computing the calibration curve.
bgrid | grid corresponding to input |
background | 2D field of background values (Y, X) |
points | observation points |
pobs | 2D vector of observations with dimensions (T, N) |
pbackground | vector of background values at points with dimensions (T, N) |
structure | structure function specifying correlation between points |
min_quantile | truncate quantile map below this quantile |
max_quantile | truncate quantile map above this quantile |
max_points | maximum number of points used within localization radius (not used at moment) |
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 | ||
) |
Nearest neighbour downscaling grid to grid and threshold and consensus in one.
igrid | Input grid |
ogrid | Output grid to downscale to |
ivalues_true | 3D vector of values on the input grid (Y, X, E) |
ivalues_false | 3D vector of values on the input grid (Y, X, E) |
threshold_values | 3D vector of values (Y, X, E), which defines the mask array for ivalues after applying the theshold |
threshold | 2D vector of threshold values |
comparison_operator | lower than, lower or equal than, greater than, great or equal than |
statistic | statistic to compute over the ensemble dimension |
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.
ref | Reference values (observations) |
fcst | Forecast values |
thresholds | Thresholds for computing optimal values for |
metric | A Metric to optimize for |
Ensure calibration curve is monotonic, by removing points.
curve_ref | Reference quantiles |
curve_fcst | Forecast quantiles |
output_fcst | New forecast quantiles |
Nearest neighbour dowscaling grid to grid.
igrid | Input grid |
ogrid | Output grid to downscale to |
ivalues | 2D vector of values on the input grid |
Nearest neighbour dowscaling grid to point.
igrid | Input grid |
ogrid | Output points to downscale to |
ivalues | 2D vector of values on the input grid |
Nearest neighbour dowscaling point to point.
ipoints | Input points |
opoints | Output points to downscale to |
ivalues | 2D vector of values on the input grid |
Nearest neighbour dowscaling point to grid.
ipoints | Input points |
ogrid | Output points to downscale to |
ivalues | 2D vector of values on the input grid |
vec2 gridpp::neighbourhood | ( | const vec2 & | input, |
int | halfwidth, | ||
gridpp::Statistic | statistic | ||
) |
Spatial neighbourhood filter, computing a statistic for a sliding square window.
input | 2D grid of values |
halfwidth | Filter halfwidth in number of gridpoints |
statistic | Statistic to compute |
vec2 gridpp::neighbourhood | ( | const vec3 & | input, |
int | halfwidth, | ||
gridpp::Statistic | statistic | ||
) |
Spatial neighbourhood filter for an ensemble of fields.
input | 3D grid of values with dimensions (Y, X, E) |
halfwidth | Filter halfwidth in number of gridpoints |
statistic | Statistic to compute |
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.
input | 2D grid of values |
halfwidth | Filter halfwidth in number of gridpoints |
operation | one of min, mean, median, max |
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.
input | 3D grid of values with dimensions (Y, X, E) |
halfwidth | Filter halfwidth in number of gridpoints |
operation | one of min, mean, median, max |
vec2 gridpp::neighbourhood_ens | ( | const vec3 & | input, |
int | halfwidth, | ||
gridpp::Statistic | statistic | ||
) |
Deprecated: Compute neighbourhood statistic on ensemble field.
Computes a quantile in a sliding square neighbourhood.
input | 2D grid of values |
quantile | Quantile to compute (between 0 and 1) |
halfwidth | Filter halfwidth in number of gridpoints |
Computes a quantile in a sliding square neighbourhood for an ensemble of fields.
input | 3D grid of values with dimensions (Y, X, E) |
quantile | Quantile to compute (between 0 and 1) |
halfwidth | Filter halfwidth in number of gridpoints |
Deprecated: Compute neighbourhood quantiles on ensemble field.
vec2 gridpp::neighbourhood_quantile_ens_fast | ( | const vec3 & | input, |
float | quantile, | ||
int | radius, | ||
const vec & | thresholds | ||
) |
Deprecated: Compute neighbourhood quantiles fast on ensemble field.
vec2 gridpp::neighbourhood_quantile_fast | ( | const vec2 & | input, |
float | quantile, | ||
int | halfwidth, | ||
const vec & | thresholds | ||
) |
Fast and approximate neighbourhood quantile.
input | 2D grid of values |
quantile | Quantile to compute (between 0 and 1) |
halfwidth | Filter halfwidth in number of gridpoints |
thresholds | Vector of thresholds to use to approximate value |
vec2 gridpp::neighbourhood_quantile_fast | ( | const vec3 & | input, |
float | quantile, | ||
int | halfwidth, | ||
const vec & | thresholds | ||
) |
Fast and approximate neighbourhood quantile for ensemble of fields.
input | 3D grid of values with dimensions (Y, X, E) |
quantile | Quantile to compute (between 0 and 1) |
halfwidth | Filter halfwidth in number of gridpoints |
thresholds | Vector of thresholds to use to approximate value |
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.
input | 2D grid of values |
quantile | 2D grid quantiles to compute (between 0 and 1) |
halfwidth | Filter halfwidth in number of gridpoints |
thresholds | Vector of thresholds to use to approximate value |
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.
input | 3D grid of values with dimensions (Y, X, E) |
quantile | 2D grid quantiles to compute (between 0 and 1) |
halfwidth | Filter halfwidth in number of gridpoints |
thresholds | Vector of thresholds to use to approximate value |
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 based on a search criteria.
If search value is within a criteria range, then the most suitable point is used. This is the nearest value of any point within the search_target range; or if no point fulfills this, the point with the highest search value.
base | base values (e.g elevation) |
values | values to compute gradients for (e.g. temperatures) |
gradient_type | what gradient type to compute |
halfwidth | neighbourhood halfwidth to compute gradient for |
min_num | minimum number of valid points needed to compute gradient |
min_range | minimum range of base values to compute gradient (for LinearRegression, this is the standard deviation of values |
default_gradient | The gradient to use when a gradient cannot be computed |
int gridpp::num_missing_values | ( | const vec2 & | iArray | ) |
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, | ||
bool | allow_extrapolation = true |
||
) |
Optimal interpolation for a deterministic gridded field.
bgrid | Grid of background field |
background | 2D field of background values |
points | Points of observations |
pobs | Vector of observations |
pratios | Vector of ratio of observation error variance to background variance |
pbackground | Background with observation operator |
structure | Structure function |
max_points | Maximum number of observations to use inside localization zone; Use 0 to disable |
allow_extrapolation | Allow OI to extrapolate increments outside increments at observations |
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, | ||
bool | allow_extrapolation = true |
||
) |
Optimal interpolation for a deterministic vector of points.
bpoints | Points of background field |
background | 1D field of background values |
points | Points of observations |
pobs | Vector of observations |
pratios | Vector of ratio of observation error variance to background variance |
pbackground | Background with observation operator |
structure | Structure function |
max_points | Maximum number of observations to use inside localization zone; Use 0 to disable |
allow_extrapolation | Allow OI to extrapolate increments outside increments at observations |
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, | ||
bool | allow_extrapolation = true |
||
) |
Optimal interpolation using a structure function based on an ensemble See Lussana et al 2019 (DOI: 10.1002/qj.3646)
input | 3D field of background values (Y, X, E) |
bgrid | grid corresponding to input |
pobs | vector of observations |
pci | vector of ci values |
points | observation points |
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, | ||
bool | allow_extrapolation = true |
||
) |
vec2 gridpp::optimal_interpolation_full | ( | const Grid & | bgrid, |
const vec2 & | background, | ||
const vec2 & | bvariance, | ||
const Points & | points, | ||
const vec & | obs, | ||
const vec & | obs_variance, | ||
const vec & | background_at_points, | ||
const vec & | bvariance_at_points, | ||
const StructureFunction & | structure, | ||
int | max_points, | ||
vec2 & | analysis_variance, | ||
bool | allow_extrapolation = true |
||
) |
Optimal interpolation for a deterministic gridded field including analysis variance.
bgrid | Grid of background field |
background | 2D field of background values |
bvariance | Variance of background field |
points | Points of observations |
obs | Vector of observations |
obs_variance | Variance of observations |
background_at_points | Background interpolated to observation points |
bvariance_at_points | Variance of background interpolated to observation points |
structure | Structure function |
max_points | Maximum number of observations to use inside localization zone; Use 0 to disable |
allow_extrapolation | Allow OI to extrapolate increments outside increments at observations |
vec gridpp::optimal_interpolation_full | ( | const Points & | bpoints, |
const vec & | background, | ||
const vec & | bvariance, | ||
const Points & | points, | ||
const vec & | obs, | ||
const vec & | obs_variance, | ||
const vec & | background_at_points, | ||
const vec & | bvariance_at_points, | ||
const StructureFunction & | structure, | ||
int | max_points, | ||
vec & | analysis_sigmas, | ||
bool | allow_extrapolation = true |
||
) |
Optimal interpolation for a deterministic vector of points including analysis variance.
bpoints | Points of background field |
background | 1D field of background values |
bvariance | Variance of background field |
points | Points of observations |
obs | Vector of observations |
obs_variance | Variance of observations |
background_at_points | Background interpolated to observation points |
bvariance_at_points | Variance of background interpolated to observation points |
structure | Structure function |
max_points | Maximum number of observations to use inside localization zone; Use 0 to disable |
allow_extrapolation | Allow OI to extrapolate increments outside increments at observations |
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)
A | A point in the rectangle |
B | A point in the rectangle |
C | A point in the rectangle |
D | A point in the rectangle |
m | The point to test if it is inside |
float gridpp::pressure | ( | float | ielev, |
float | oelev, | ||
float | ipressure, | ||
float | itemperature = 288.15 |
||
) |
Calculate pressure at a new elevation.
ielev | Elevation at start point |
oelev | Elevation at new point |
ipressure | Pressure at start point |
itemperature | Temperature at start point |
vec gridpp::pressure | ( | const vec & | ielev, |
const vec & | oelev, | ||
const vec & | ipressure, | ||
const vec & | itemperature | ||
) |
Calculate Vector version of pressure calculation.
ielev | Elevations at start point |
oelev | Elevations at new point |
ipressure | Pressures at start point |
itemperature | Temperatures at start point |
float gridpp::qnh | ( | float | pressure, |
float | altitude | ||
) |
Diagnose QNH from pressure and altitude.
pressure | Pressure at point [pa] |
altitude | Altitude of point [m] |
Vector version of QNH calculation.
pressure | Pressures at points [pa] |
altitude | Altitudes of points [m] |
vec gridpp::quantile_mapping_curve | ( | const vec & | ref, |
const vec & | fcst, | ||
vec & | output_fcst, | ||
vec | quantiles = vec() |
||
) |
Create quantile mapping calibration curve.
ref | Reference values (observations) |
fcst | Forecast values |
output_fcst | Output forecast quantiles |
quantiles | Vector of quantiles to extract. If empty, use all values. |
float gridpp::relative_humidity | ( | float | temperature, |
float | dewpoint | ||
) |
Calculate relative humidity from temperature and dewpoint temperature.
temperature | Temperature [K] |
dewpoint | Dewpoint temperature [K] |
Vector version of relative humidity calculation.
temperature | Temperatures [K] |
dewpoint | Dewpoint temperatures [K] |
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.
ps | Surface pressure [pa] |
altitude | Station altitude above sea level [m] |
temperature | 2m temperature [K] |
rh | 2m Relative humidity [1] |
dewpoint | 2m Dewpoint Temperature at station [K] |
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.
ps | Surface pressures [pa] |
altitude | Station altitudes above sea level [m] |
temperature | 2m temperatures [K] |
rh | 2m Relative humidities [1] |
dewpoint | 2m Dewpoint Temperatures at stations [K] |
void gridpp::set_debug_level | ( | int | level | ) |
Set the verbosity of debug messages.
Use 0 for no messages.
void gridpp::set_omp_threads | ( | int | num | ) |
Set the number of OpenMP threads to use.
Overrides OMP_NUM_THREAD env variable.
vec2 gridpp::simple_gradient | ( | const Grid & | igrid, |
const Grid & | ogrid, | ||
const vec2 & | ivalues, | ||
float | elev_gradient, | ||
Downscaler | downscaler = Nearest |
||
) |
vec3 gridpp::simple_gradient | ( | const Grid & | igrid, |
const Grid & | ogrid, | ||
const vec3 & | ivalues, | ||
float | elev_gradient, | ||
Downscaler | downscaler = Nearest |
||
) |
vec gridpp::simple_gradient | ( | const Grid & | igrid, |
const Points & | opoints, | ||
const vec2 & | ivalues, | ||
float | elev_gradient, | ||
Downscaler | downscaler = Nearest |
||
) |
vec2 gridpp::simple_gradient | ( | const Grid & | igrid, |
const Points & | opoints, | ||
const vec3 & | ivalues, | ||
float | elev_gradient, | ||
Downscaler | downscaler = Nearest |
||
) |
vec2 gridpp::smart | ( | const Grid & | igrid, |
const Grid & | ogrid, | ||
const vec2 & | ivalues, | ||
int | num, | ||
const StructureFunction & | structure | ||
) |
Smart neighbour downscaling grid to grid.
igrid | Input grid |
ogrid | Output points to downscale to |
ivalues | 2D vector of values on the input grid |
num | Number of neighbours to average |
structure | Structure function for determining similarity |
float * gridpp::test_array | ( | float * | v, |
int | n | ||
) |
Special function whose presense is needed for SWIG.
ivec2 gridpp::test_ivec2_output | ( | ) |
ivec3 gridpp::test_ivec3_output | ( | ) |
int gridpp::test_ivec_input | ( | const ivec & | input | ) |
Testing function for 1D input vector.
ivec gridpp::test_ivec_output | ( | ) |
void gridpp::test_not_implemented_exception | ( | ) |
float gridpp::test_vec2_argout | ( | vec2 & | distances | ) |
Testing function for 2D vector treated as output.
float gridpp::test_vec2_input | ( | const vec2 & | input | ) |
Testing function for 2D input vector.
vec2 gridpp::test_vec2_output | ( | ) |
Testing function for 2D output vector.
float gridpp::test_vec3_input | ( | const vec3 & | input | ) |
Testing function for 3D input vector.
vec3 gridpp::test_vec3_output | ( | ) |
Testing function for 3D output vector.
float gridpp::test_vec_argout | ( | vec & | distances | ) |
Testing function for 1D vector treated as output.
float gridpp::test_vec_input | ( | const vec & | input | ) |
Testing function for 1D input vector.
vec gridpp::test_vec_output | ( | ) |
Testing function for 1D output vector.
std::string gridpp::version | ( | ) |
The gridpp version.
void gridpp::warning | ( | std::string | string | ) |
float gridpp::wetbulb | ( | float | temperature, |
float | pressure, | ||
float | relative_humidity | ||
) |
Calculate wetbulb temperature from temperature, pressure, and relative humidity.
temperature | Temperature [K] |
pressure | Air pressure [pa] |
Relative | humidity [1] |
vec gridpp::wetbulb | ( | const vec & | temperature, |
const vec & | pressure, | ||
const vec & | relative_humidity | ||
) |
Vector version of wetbulb calculation.
temperature | Temperatures [K] |
pressure | Air pressures [pa] |
Relative | humidities [1] |
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
xwind | X-component of wind [any unit] |
ywind | Y-component of wind [any unit] |
Vector version of wind direction calculation.
xwind | X-components of wind [any unit] |
ywind | Y-components of wind [any unit] |
float gridpp::wind_speed | ( | float | xwind, |
float | ywind | ||
) |
Diagnose wind speed from its components.
xwind | X-component of wind [any unit] |
ywind | Y-component of wind [any unit] |
Vector version of wind speed calculation.
xwind | X-components of wind [any unit] |
ywind | Y-components of wind [any unit] |
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.
array | input array with dimensions (case, time) |
length | window length in number of timesteps |
statistic | statistic to apply to window |
before | if true, make the window end at the particular time. If false, centre it. |
keep_missing | if true, window value will be missing if one or more values in window are missing |
missing_edges | if true put missing values at the edges, where window overshoots the edge |
|
static |
|
static |
Universal Gas Constant [kg*m^2*s^-2/(K*mol)].
|
static |
Universal Gas Constant [J/(kg*K)].
|
static |
Gravitational acceleration [m/s^2].
|
static |
Constant Lapse Rate moist air standard atmosphere [K/m].
|
static |
Molar Mass of Dry Air [kg/mol].
|
static |
Missing value indicator.
|
static |
Missing value indicator in gridpp command-line tool.
|
static |
Mathematical constant pi.
|
static |
Radius of the earth [m].
|
static |
Temperature at surface in standard atmosphere [K].
|
static |
Default value used to fill array in SWIG testing functions.
Not useful for any other purpose.