reconciliation_lib 0.1
Library for data reconciliation algorithms
Loading...
Searching...
No Matches
dvrlib Namespace Reference

Classes

struct  assertion_error
 
class  func
 Abstract function object mapping argtype to restype. More...
 
struct  gsl_exception
 Exception type thrown by GSL when gsl_enable_exceptions() is active. More...
 
class  matrix
 RAII wrapper around a heap-allocated gsl_matrix. More...
 
class  matrix_view
 Non-owning view into a rectangular sub-region of a matrix. More...
 
class  ostream_format_guard
 
struct  permutation
 
class  recon_system
 Encapsulates a system to be reconciled, including variables, confidence intervals, covariance coefficient, etc. More...
 
class  vector
 RAII wrapper around a heap-allocated gsl_vector. More...
 
class  vector_view
 Non-owning view into a contiguous slice of a vector or matrix row/column. More...
 

Functions

void gsl_err_handler (const char *reason, const char *file, int line, int gsl_errno)
 
void gsl_enable_exceptions ()
 Install a GSL error handler that throws gsl_exception instead of aborting.
 
void set_vdi_print_format (bool enable)
 Switch between standard and VDI 2048 pretty-print format for vectors and matrices written to std::ostream.
 
vector operator* (double d, const vector &src)
 Scalar multiplication d*src (commutative with vector::operator*).
 
std::ostream & print_std_format (std::ostream &out, const vector &vec)
 
std::ostream & print_vdi_format (std::ostream &out, const vector &vec)
 
std::ostream & operator<< (std::ostream &out, const vector &vec)
 Write a vector to an output stream.
 
std::ostream & print_std_format (std::ostream &out, const matrix &mat)
 
std::ostream & print_vdi_format (std::ostream &out, const matrix &mat)
 
std::ostream & operator<< (std::ostream &out, const matrix &mat)
 Write a matrix to an output stream.
 
matrix operator* (double d, const matrix &src)
 Scalar multiplication d*src (commutative with matrix::operator*).
 
std::ostream & operator<< (std::ostream &out, const gsl_vector &v)
 
std::ostream & operator<< (std::ostream &out, const gsl_matrix &m)
 
std::ostream & operator<< (std::ostream &out, const gsl_matrix_view &mv)
 
void lin_cov_update_Streit (const matrix &S_x, const matrix &F, matrix &S_v)
 
void lin_cov_update_Zander (const matrix &S_x, const matrix &F, matrix &S_v)
 Compute the update of the covariance matrix.
 
void lin_cov_update (const matrix &S_x, const matrix &F, matrix &S_v)
 Propagate measurement covariance through a linear constraint matrix.
 
void lin_recon (const vector &r, const matrix &S_x, const matrix &F, vector &v)
 Solve the linear reconciliation problem.
 
void lin_recon_update (const vector &r, const matrix &S_x_inv, const matrix &F, const vector &v, vector &dv)
 Compute the reconciliation correction step for non-linear iteration.
 
int recon (const vector &x, const matrix &S_x, const func< vector, vector > &f, const func< vector, matrix > &J, vector &v, matrix &S_v, double eps=1e-6, int maxiter=50)
 Solve the non-linear reconciliation problem by iterative linearisation.
 
double confint2var (double confint)
 Convert a 95% confidence interval into a variance.
 
double var2confint (double var)
 Convert a variance into a 95% confidence interval.
 
void extract_confidence (const matrix &S_xnew, vector &conf_results)
 Extract 95% confidence intervals from a covariance matrix diagonal.
 

Variables

static bool g_vdi_format = false
 
static constexpr auto z_alpha_95 = 1.959963984540054
 

Function Documentation

◆ confint2var()

double dvrlib::confint2var ( double  confint)

Convert a 95% confidence interval into a variance.

Uses the relation confint = z_{0.975} * sqrt(var) where z_{0.975} = Φ⁻¹(0.975) ≈ 1.96.

◆ extract_confidence()

void dvrlib::extract_confidence ( const matrix S_xnew,
vector conf_results 
)

Extract 95% confidence intervals from a covariance matrix diagonal.

Parameters
[in]S_xnewCovariance matrix of the reconciled values.
[out]conf_resultsConfidence intervals (one per variable).

◆ gsl_enable_exceptions()

void dvrlib::gsl_enable_exceptions ( )

Install a GSL error handler that throws gsl_exception instead of aborting.

◆ gsl_err_handler()

void dvrlib::gsl_err_handler ( const char reason,
const char file,
int  line,
int  gsl_errno 
)

◆ lin_cov_update()

void dvrlib::lin_cov_update ( const matrix S_x,
const matrix F,
matrix S_v 
)

Propagate measurement covariance through a linear constraint matrix.

Computes S_v = F * S_x * F^T, the covariance of the reconciled values given the measurement covariance S_x and the constraint matrix F.

Parameters
[in]S_xMeasurement covariance matrix.
[in]FLinear constraint (Jacobian) matrix.
[out]S_vResulting covariance of the reconciled values.

◆ lin_cov_update_Streit()

void dvrlib::lin_cov_update_Streit ( const matrix S_x,
const matrix F,
matrix S_v 
)

◆ lin_cov_update_Zander()

void dvrlib::lin_cov_update_Zander ( const matrix S_x,
const matrix F,
matrix S_v 
)

Compute the update of the covariance matrix.

The update is computed by the following formula

\[ S_v = G S_F G^T \]

where

\[ G = P_m * Z^{-1} * P_l \]

and

\[ S_F = F_m * S_x * F_m. \]

$P_m$ and $P_l$ indicate projection matrices on the space of measured variables and on the space of Lagrange multipliers, respectively.

Parameters
[in]S_xThe input covariance matrix
[in]FThe matrix of constraints
[out]S_vThe covariance of the updates

◆ lin_recon()

void dvrlib::lin_recon ( const vector r,
const matrix S_x,
const matrix F,
vector v 
)

Solve the linear reconciliation problem.

Computes the reconciled values v that satisfy the linear constraints encoded in F, given raw measurements r and measurement covariance S_x.

Parameters
[in]rRaw measurement vector.
[in]S_xMeasurement covariance matrix.
[in]FLinear constraint matrix.
[out]vReconciled values.

◆ lin_recon_update()

void dvrlib::lin_recon_update ( const vector r,
const matrix S_x_inv,
const matrix F,
const vector v,
vector dv 
)

Compute the reconciliation correction step for non-linear iteration.

Given the current iterate v, computes the Newton-like correction dv to reduce the residual of the constraint equations.

Parameters
[in]rRaw measurement vector.
[in]S_x_invInverse of the measurement covariance matrix.
[in]FConstraint Jacobian evaluated at the current iterate.
[in]vCurrent reconciled values.
[out]dvCorrection to apply: next iterate is v + dv.

◆ operator*() [1/2]

matrix dvrlib::operator* ( double  d,
const matrix src 
)

Scalar multiplication d*src (commutative with matrix::operator*).

◆ operator*() [2/2]

vector dvrlib::operator* ( double  d,
const vector src 
)

Scalar multiplication d*src (commutative with vector::operator*).

◆ operator<<() [1/5]

std::ostream & dvrlib::operator<< ( std::ostream &  out,
const gsl_matrix m 
)

◆ operator<<() [2/5]

std::ostream & dvrlib::operator<< ( std::ostream &  out,
const gsl_matrix_view mv 
)

◆ operator<<() [3/5]

std::ostream & dvrlib::operator<< ( std::ostream &  out,
const gsl_vector v 
)

◆ operator<<() [4/5]

std::ostream & dvrlib::operator<< ( std::ostream &  out,
const matrix mat 
)

Write a matrix to an output stream.

Format depends on set_vdi_print_format().

◆ operator<<() [5/5]

std::ostream & dvrlib::operator<< ( std::ostream &  out,
const vector vec 
)

Write a vector to an output stream.

Format depends on set_vdi_print_format().

◆ print_std_format() [1/2]

std::ostream & dvrlib::print_std_format ( std::ostream &  out,
const matrix mat 
)

◆ print_std_format() [2/2]

std::ostream & dvrlib::print_std_format ( std::ostream &  out,
const vector vec 
)

◆ print_vdi_format() [1/2]

std::ostream & dvrlib::print_vdi_format ( std::ostream &  out,
const matrix mat 
)

◆ print_vdi_format() [2/2]

std::ostream & dvrlib::print_vdi_format ( std::ostream &  out,
const vector vec 
)

◆ recon()

int dvrlib::recon ( const vector x,
const matrix S_x,
const func< vector, vector > &  f,
const func< vector, matrix > &  J,
vector v,
matrix S_v,
double  eps = 1e-6,
int  maxiter = 50 
)

Solve the non-linear reconciliation problem by iterative linearisation.

Iterates a Newton-type update until the correction dv is smaller than eps (in the L2 norm) or maxiter iterations are reached.

Parameters
[in]xRaw measurement vector.
[in]S_xMeasurement covariance matrix.
[in]fConstraint residual function f(v) = 0.
[in]JJacobian of f with respect to v.
[out]vReconciled values.
[out]S_vCovariance of the reconciled values.
[in]epsConvergence tolerance on the correction norm (default 1e-6).
[in]maxiterMaximum number of iterations (default 50).
Returns
Number of iterations performed, or -1 if the method did not converge.

◆ set_vdi_print_format()

void dvrlib::set_vdi_print_format ( bool  enable)

Switch between standard and VDI 2048 pretty-print format for vectors and matrices written to std::ostream.

Parameters
[in]enabletrue selects VDI format, false selects standard format.

◆ var2confint()

double dvrlib::var2confint ( double  var)

Convert a variance into a 95% confidence interval.

Uses the relation confint = z_{0.975} * sqrt(var) where z_{0.975} = Φ⁻¹(0.975) ≈ 1.96.

Variable Documentation

◆ g_vdi_format

bool dvrlib::g_vdi_format = false
static

◆ z_alpha_95

constexpr auto dvrlib::z_alpha_95 = 1.959963984540054
staticconstexpr