Titanlib 0.3.0-dev3
Library for quality control algorithms
Enumerations | Classes
titanlib Namespace Reference

Typedefs

Short-hand notation for vectors of different dimensions sizes

typedef std::vector< int > ivec
 
typedef std::vector< float > vec
 
typedef std::vector< double > dvec
 
typedef std::vector< vecvec2
 

Functions

Spatial checks

Checks that consider the spatial properties of observations

ivec sct (const Points &points, const vec &values, int num_min, int num_max, float inner_radius, float outer_radius, int num_iterations, int num_min_prof, float min_elev_diff, float min_horizontal_scale, float vertical_scale, const vec &pos, const vec &neg, const vec &eps2, vec &prob_gross_error, vec &rep)
 Spatial Consistency Test. More...
 
ivec sct_resistant (const Points &points, const vec &values, const ivec &obs_to_check, const vec &background_values, BackgroundType background_elab_type, int num_min_outer, int num_max_outer, float inner_radius, float outer_radius, int num_iterations, int num_min_prof, float min_elev_diff, float min_horizontal_scale, float max_horizontal_scale, int kth_closest_obs_horizontal_scale, float vertical_scale, const vec &value_mina, const vec &value_maxa, const vec &value_minv, const vec &value_maxv, const vec &eps2, const vec &tpos, const vec &tneg, bool debug, vec &scores)
 Spatial Consistency Test (SCT) - resistant to outliers. More...
 
ivec fgt (const Points &points, const vec &values, const ivec &obs_to_check, const vec &background_values, const vec &background_uncertainties, BackgroundType background_elab_type, int num_min_outer, int num_max_outer, float inner_radius, float outer_radius, int num_iterations, int num_min_prof, float min_elev_diff, const vec &value_mina, const vec &value_maxa, const vec &value_minv, const vec &value_maxv, const vec &tpos, const vec &tneg, bool debug, vec &scores)
 First Guess Test (FGT) - simplified (without OI) SCT. More...
 
ivec range_check (const vec &values, const vec &min, const vec &max)
 Range check. More...
 
ivec range_check_climatology (const Points &points, const vec &values, int unixtime, const vec &pos, const vec &neg)
 
ivec buddy_check (const Points &points, const vec &values, const vec &radius, const ivec &num_min, float threshold, float max_elev_diff, float elev_gradient, float min_std, int num_iterations, const ivec &obs_to_check=ivec())
 Buddy check. More...
 
ivec buddy_event_check (const Points &points, const vec &values, const vec &radius, const ivec &num_min, float event_threshold, float threshold, float max_elev_diff, float elev_gradient, int num_iterations, const ivec &obs_to_check=ivec())
 
ivec isolation_check (const Points &points, int num_min, float radius, float vertical_radius=MV)
 Isolation check. More...
 
ivec duplicate_check (const Points &points, float radius, float vertical_range=titanlib::MV)
 Duplicate check. More...
 
ivec metadata_check (const Points &points, bool check_lat=true, bool check_lon=true, bool check_elev=true, bool check_laf=true)
 
Timeserie methods

Functions that operate on timeseries of observations

vec lag_reduction_filter (const vec &times, const vec &values, float a=0.5, float b=0.5, float k1=0.25, float k2=0.25, int n=10)
 Method by McCarthy 1973 https://doi.org/10.1175/1520-0450(1973)012%3C0211:AMFCAT%3E2.0.CO;2. 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

std::string version ()
 
double clock ()
 
bool is_valid (float value)
 
bool convert_coordinates (const vec &lats, const vec &lons, vec &x_coords, vec &y_coords, vec &z_coords)
 Convert lat/lons to 3D cartesian coordinates with the centre of the earth as the origin. More...
 
bool convert_coordinates (float lat, float lon, float &x_coord, float &y_coord, float &z_coord)
 Same as above, but convert a single lat/lon to 3D cartesian coordinates. More...
 
vec interpolate_to_points (const vec2 &input_lats, const vec2 &input_lons, const vec2 &input_values, const vec &output_lats, const vec &output_lons)
 
float deg2rad (float deg)
 
float calc_distance (float lat1, float lon1, float lat2, float lon2)
 
float calc_distance (float x0, float y0, float z0, float x1, float y1, float z1)
 
float compute_quantile (double quantile, const vec &array)
 
float find_k_closest (const vec &array, int k)
 
vec subset (const vec &array, const ivec &indices)
 
Points subset (const Points &input, const ivec &indices)
 
vec background (const vec &elevs, const vec &values, int num_min_prof, float min_elev_diff, float value_minp, float value_maxp, BackgroundType background_type, const vec &external_background_values, const ivec &indices_global_outer, bool debug)
 Background (first guess) calculations at observation locations. More...
 
bool invert_matrix (const boost::numeric::ublas::matrix< float > &input, boost::numeric::ublas::matrix< float > &inverse)
 
bool set_indices (const ivec &indices_global_outer_guess, const ivec &obs_test, const ivec &dqcflags, const vec &dist_outer_guess, float inner_radius, int test_just_this, ivec &indices_global_outer, ivec &indices_global_test, ivec &indices_outer_inner, ivec &indices_outer_test, ivec &indices_inner_test)
 
SWIG testing functions

Functions for testing the SWIG interface.

Not useful for any other purpose.

float * test_array (float *v, int n)
 Required for SWIG only. More...
 
void test_not_implemented_exception ()
 

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  BackgroundType {
  VerticalProfile = 0, VerticalProfileTheilSen = 1, MeanOuterCircle = 2, MedianOuterCircle = 3,
  External = 4
}
 
enum  CoordinateType { Geodetic = 0, Cartesian = 1 }
 Types of coordinates for position of points. More...
 

Classes

class  Dataset
 Represents point and their observed values. More...
 
class  KDTree
 
class  not_implemented_exception
 
class  Point
 
class  Points
 
struct  sort_pair_first
 

Typedef Documentation

◆ dvec

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

◆ ivec

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

◆ vec

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

◆ vec2

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

Enumeration Type Documentation

◆ BackgroundType

Enumerator
VerticalProfile 
VerticalProfileTheilSen 
MeanOuterCircle 
MedianOuterCircle 
External 

◆ CoordinateType

Types of coordinates for position of points.

Enumerator
Geodetic 

Latitude and longitude.

Cartesian 

X and Y.

Function Documentation

◆ background()

vec titanlib::background ( const vec elevs,
const vec values,
int  num_min_prof,
float  min_elev_diff,
float  value_minp,
float  value_maxp,
BackgroundType  background_type,
const vec external_background_values,
const ivec indices_global_outer,
bool  debug 
)

Background (first guess) calculations at observation locations.

Parameters
background_typebackground elaboration type
elevselevations (m above mean sea level)
valuesobserved values
num_min_profminimum number of observations needed to compute a vertical profile different from the pseudo-adiabatic lapse rate (-0.0065 degC/m)
min_elev_diffwhen computing vertical profiles, the elevation difference between the 95th percentile and 5th percentile must be greater than this parameter. Otherwise use the pseudo-adiabatic lapse rate
value_minpminimum plausible value
value_maxpmaximum plausible value
external_background_valuesexternal background values (used when background_elab_type=external)
indices_global_outervector of positions of matches of (vector with the observations belonging to the outer circle) in the (global vector) (dimension is the number of observations in the outer circle)
debugverbose mode (true / false)

◆ buddy_check()

ivec titanlib::buddy_check ( const Points points,
const vec values,
const vec radius,
const ivec num_min,
float  threshold,
float  max_elev_diff,
float  elev_gradient,
float  min_std,
int  num_iterations,
const ivec obs_to_check = ivec() 
)

Buddy check.

Compares a station to all its neighbours within a certain distance

Parameters
valuesvector of observation values
radiussearch radius [m]
num_minthe minimum number of buddies a station can have
thresholdthe threshold for flagging a station
max_elev_diffthe maximum difference in elevation for a buddy (if negative will not check for heigh difference)
elev_gradientlinear elevation gradient with height. For temperature, use something like -0.0065.
min_std
num_iterations
obs_to_checkthe observations that will be checked (since can pass in observations that will not be checked)
flagsvector of return flags

◆ buddy_event_check()

ivec titanlib::buddy_event_check ( const Points points,
const vec values,
const vec radius,
const ivec num_min,
float  event_threshold,
float  threshold,
float  max_elev_diff,
float  elev_gradient,
int  num_iterations,
const ivec obs_to_check = ivec() 
)

◆ calc_distance() [1/2]

float titanlib::calc_distance ( float  lat1,
float  lon1,
float  lat2,
float  lon2 
)

◆ calc_distance() [2/2]

float titanlib::calc_distance ( float  x0,
float  y0,
float  z0,
float  x1,
float  y1,
float  z1 
)

◆ clock()

double titanlib::clock ( )
Returns
The current UTC time (in seconds since 1970-01-01 00:00:00Z)

◆ compute_quantile()

float titanlib::compute_quantile ( double  quantile,
const vec array 
)

◆ convert_coordinates() [1/2]

bool titanlib::convert_coordinates ( const vec lats,
const vec lons,
vec x_coords,
vec y_coords,
vec z_coords 
)

Convert lat/lons to 3D cartesian coordinates with the centre of the earth as the origin.

Parameters
latsvector of latitudes [deg]
lonsvector of longitudes [deg]
x_coordsvector of x-coordinates [m]
y_coordsvector of y-coordinates [m]
z_coordsvector of z-coordinates [m]

◆ convert_coordinates() [2/2]

bool titanlib::convert_coordinates ( float  lat,
float  lon,
float &  x_coord,
float &  y_coord,
float &  z_coord 
)

Same as above, but convert a single lat/lon to 3D cartesian coordinates.

Parameters
latlatitude [deg]
lonlongitude [deg]
x_coordx-coordinate [m]
y_coordy-coordinate [m]
z_coordz-coordinate [m]

◆ deg2rad()

float titanlib::deg2rad ( float  deg)

◆ duplicate_check()

ivec titanlib::duplicate_check ( const Points points,
float  radius,
float  vertical_range = titanlib::MV 
)

Duplicate check.

Checks Flags duplicate stations. Keeps the first one when duplicates.

Parameters
radiusStations within this radius are considered duplicates [m]
vertical_rangeStations must also be within this vertical range to be considered duplicates [m]. Disable if MV.

◆ fgt()

ivec titanlib::fgt ( const Points points,
const vec values,
const ivec obs_to_check,
const vec background_values,
const vec background_uncertainties,
BackgroundType  background_elab_type,
int  num_min_outer,
int  num_max_outer,
float  inner_radius,
float  outer_radius,
int  num_iterations,
int  num_min_prof,
float  min_elev_diff,
const vec value_mina,
const vec value_maxa,
const vec value_minv,
const vec value_maxv,
const vec tpos,
const vec tneg,
bool  debug,
vec scores 
)

First Guess Test (FGT) - simplified (without OI) SCT.

Parameters
pointsInput points
valuesobserved values to check (and/or to use)
obs_to_checkObservations that will be checked (since can pass in observations that will not be checked). 1=check the corresponding observation
background_valuesexternal background value (not used if background_elab_type!=external)
background_uncertaintiesuncertainty of the external background value (not used if background_elab_type!=external, optional when background_elab_type=external)
background_elab_typeone of: vertical_profile, vertical_profile_Theil_Sen, mean_outer_circle, external
num_min_outerMinimum number of observations inside the outer circle to compute FGT
num_max_outerMaximum number of observations inside the outer circle used
num_min_profMinimum number of observations to compute vertical profile
inner_radiusRadius for flagging [m]
outer_radiusRadius for computing OI and background [m]
num_iterationsNumber of FGT iterations
min_elev_diffMinimum elevation difference to compute vertical profile [m]
value_minaMinimum admissible value
value_maxaMaximum admissible value
value_minvMinimum valid value
value_maxvMaximum valid value
tposFGT-score threshold. Positive deviation allowed
tnegFGT-score threshold. Negative deviation allowed
debugVerbose output
scoresFGT-score. The higher the score, the more likely is the presence of a gross measurement error
Returns
flags

◆ find_k_closest()

float titanlib::find_k_closest ( const vec array,
int  k 
)

◆ initialize_omp()

void titanlib::initialize_omp ( )

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

◆ interpolate_to_points()

vec titanlib::interpolate_to_points ( const vec2 input_lats,
const vec2 input_lons,
const vec2 input_values,
const vec output_lats,
const vec output_lons 
)

◆ invert_matrix()

bool titanlib::invert_matrix ( const boost::numeric::ublas::matrix< float > &  input,
boost::numeric::ublas::matrix< float > &  inverse 
)

◆ is_valid()

bool titanlib::is_valid ( float  value)

◆ isolation_check()

ivec titanlib::isolation_check ( const Points points,
int  num_min,
float  radius,
float  vertical_radius = MV 
)

Isolation check.

Checks that a station is not located alone

Parameters
num_minrequired number of observations
radiussearch radius [m]
vertical_radiusVertical search radius [m]
flagsvector of return flags

◆ lag_reduction_filter()

vec titanlib::lag_reduction_filter ( const vec times,
const vec values,
float  a = 0.5,
float  b = 0.5,
float  k1 = 0.25,
float  k2 = 0.25,
int  n = 10 
)

◆ metadata_check()

ivec titanlib::metadata_check ( const Points points,
bool  check_lat = true,
bool  check_lon = true,
bool  check_elev = true,
bool  check_laf = true 
)

◆ range_check()

ivec titanlib::range_check ( const vec values,
const vec min,
const vec max 
)

Range check.

Checks observation is within the ranges given

Parameters
valuesvector of observations
minmin allowed value
maxmax allowed value
flagsvector of return flags

◆ range_check_climatology()

ivec titanlib::range_check_climatology ( const Points points,
const vec values,
int  unixtime,
const vec pos,
const vec neg 
)

◆ sct()

ivec titanlib::sct ( const Points points,
const vec values,
int  num_min,
int  num_max,
float  inner_radius,
float  outer_radius,
int  num_iterations,
int  num_min_prof,
float  min_elev_diff,
float  min_horizontal_scale,
float  vertical_scale,
const vec pos,
const vec neg,
const vec eps2,
vec prob_gross_error,
vec rep 
)

Spatial Consistency Test.

Parameters
num_min_profMinimum number of observations to compute vertical profile
inner_radiusRadius for flagging [m]
outer_radiusRadius for computing OI and background [m]
min_elev_diffMinimum elevation difference to compute vertical profile [m]
min_horizontal_scaleMinimum horizontal decorrelation length [m]
vertical_scaleVertical decorrelation length [m]
posPositive deviation allowed
negNegative deviation allowed
eps2
prob_gross_errorProbability of gross error for each observation
repCoefficient of representativity
Returns
flags

◆ sct_resistant()

ivec titanlib::sct_resistant ( const Points points,
const vec values,
const ivec obs_to_check,
const vec background_values,
BackgroundType  background_elab_type,
int  num_min_outer,
int  num_max_outer,
float  inner_radius,
float  outer_radius,
int  num_iterations,
int  num_min_prof,
float  min_elev_diff,
float  min_horizontal_scale,
float  max_horizontal_scale,
int  kth_closest_obs_horizontal_scale,
float  vertical_scale,
const vec value_mina,
const vec value_maxa,
const vec value_minv,
const vec value_maxv,
const vec eps2,
const vec tpos,
const vec tneg,
bool  debug,
vec scores 
)

Spatial Consistency Test (SCT) - resistant to outliers.

Parameters
pointsInput points
valuesobserved values to check (and/or to use)
obs_to_checkObservations that will be checked (since can pass in observations that will not be checked). 1=check the corresponding observation
background_valuesexternal background value (not used if background_elab_type!=external)
background_elab_typeone of: vertical_profile, vertical_profile_Theil_Sen, mean_outer_circle, external
num_min_outerMinimum number of observations inside the outer circle to compute SCT
num_max_outerMaximum number of observations inside the outer circle used
num_min_profMinimum number of observations to compute vertical profile
inner_radiusRadius for flagging [m]
outer_radiusRadius for computing OI and background [m]
num_iterationsNumber of SCT iterations
min_elev_diffMinimum elevation difference to compute vertical profile [m]
min_horizontal_scaleMinimum horizontal decorrelation length [m]
max_horizontal_scaleMaximum horizontal decorrelation length [m]
kth_closest_obs_horizontal_scaleNumber of closest observations to consider in the adaptive estimation of the horizontal decorrelation length
vertical_scaleVertical decorrelation length [m]
value_minaMinimum admissible value
value_maxaMaximum admissible value
value_minvMinimum valid value
value_maxvMaximum valid value
eps2ratio between observation and background error variances
tposSCT-score threshold. Positive deviation allowed
tnegSCT-score threshold. Negative deviation allowed
debugVerbose output
scoresSCT-score. The higher the score, the more likely is the presence of a gross measurement error
Returns
flags

◆ set_indices()

bool titanlib::set_indices ( const ivec indices_global_outer_guess,
const ivec obs_test,
const ivec dqcflags,
const vec dist_outer_guess,
float  inner_radius,
int  test_just_this,
ivec indices_global_outer,
ivec indices_global_test,
ivec indices_outer_inner,
ivec indices_outer_test,
ivec indices_inner_test 
)

◆ set_omp_threads()

void titanlib::set_omp_threads ( int  num)

Set the number of OpenMP threads to use.

Overwrides OMP_NUM_THREAD env variable.

◆ subset() [1/2]

vec titanlib::subset ( const vec array,
const ivec indices 
)

◆ subset() [2/2]

Points titanlib::subset ( const Points input,
const ivec indices 
)

◆ test_array()

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

Required for SWIG only.

◆ test_not_implemented_exception()

void titanlib::test_not_implemented_exception ( )

◆ version()

std::string titanlib::version ( )
Returns
Titanlib version

Variable Documentation

◆ MV

const float titanlib::MV = NAN
static

Missing value indicator.

◆ MV_CML

const float titanlib::MV_CML = -999
static

Missing value indicator in gridpp command-line tool.

◆ pi

const float titanlib::pi = 3.14159265
static

Mathematical constant pi.

◆ radius_earth

const double titanlib::radius_earth = 6.378137e6
static

Radius of the earth [m].