dataprocess.h File Reference

Functions to process data. More...

#include <gsl/gsl_matrix.h>
#include <gsl/gsl_vector.h>
#include <gsl/gsl_blas.h>
#include <gsl/gsl_linalg.h>
#include <gsl/gsl_rng.h>
#include <gsl/gsl_randist.h>
#include <gsl/gsl_multifit_nlin.h>
#include <gsl/gsl_deriv.h>
#include "../functions/postprocess.h"
Include dependency graph for dataprocess.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

int get_trans_matrix (vector< double > &datapoints_base, vector< double > &datapoints_new, vector< vector< double > > &rotmat, vector< double > &p_vec)
 Takes in 3-d position data of a set of points in 2 poses and calculates the rotation matrix and translation vector.
void getavgvec (vector< double > &datavec, vector< double > &avgvec, int np)
 Get mean position vector.
void get_distri_matrices (gsl_matrix *G, gsl_matrix *P, vector< double > &datapoints_new, double *pav_ar, vector< double > &datapoints_base, double *a_ar)
 Get matrices that describe point distribution.
void get_matrix_adjoin (gsl_matrix *M, gsl_matrix *M_adj)
 Get matrix adjoint of given matrix (using gsl)
int sol_betas_fdf (gsl_vector *betas, void *g_vals, gsl_vector *err, gsl_matrix *J)
 Equation (and its Jacobian) to be solved using least squares.
int sol_betas_f (gsl_vector *betas, void *g_vals, gsl_vector *err)
 Equations to be solved using least squares.
int sol_betas_df (gsl_vector *betas, void *g_vals, gsl_matrix *J)
 Jacobian of set of equations to be solved using least squares.
int get_betas (double &b1, double &b2, double g1sq, double g2sq, double g3)
void getRot (gsl_matrix *R, gsl_matrix *G, gsl_matrix *GtG, double b1, double b2)

Detailed Description

Functions to process data.

This file includes function definitions for those functions that will be needed for post-processing the data, and combining it to obtain useful results regarding motion, position and orientation of collar and so on. Similar to postprocess.h, but with other functions that build on it.


Function Documentation

int get_betas ( double &  b1,
double &  b2,
double  g1sq,
double  g2sq,
double  g3 
)

References g_vals::g1sq, g_vals::g2sq, print_results_lm(), sol_betas_df(), sol_betas_f(), sol_betas_fdf(), and solve_gsl_lm().

Referenced by get_trans_matrix().

Here is the call graph for this function:

Here is the caller graph for this function:

void get_distri_matrices ( gsl_matrix *  G,
gsl_matrix *  P,
vector< double > &  datapoints_new,
double *  pav_ar,
vector< double > &  datapoints_base,
double *  a_ar 
)

Get matrices that describe point distribution.

The function returns two matrices G and P given by the equations:

\[ G = \frac{1}{m} \times \sum_{i=1}^m ((p_i-p)(a_i-a)') \]

\[ P = \frac{1}{m} \times \sum_{i=1}^m ((p_i-p)(p_i-p)') \]

where $p_i$ = datapoints in the new pose, $p$ = average position vector for the points in new pose, $ a_i$ = datapoints in the old pose, $ a$ = average position vector for the points in old pose.

Parameters:
[out]GThe matrix relating point orientation in old pose to new pose
[out]PThe matrix describing point distribution, in new pose
[in]datapoints_newVector containing position data of all points in new pose
[in]pav_ar= array containing average position vector for the points in new pose
[in]datapoints_baseVector containing position data of all points in base pose
[in]a_ar= array containing average position vector for the points in base pose

References vect2double().

Referenced by get_trans_matrix().

Here is the call graph for this function:

Here is the caller graph for this function:

void get_matrix_adjoin ( gsl_matrix *  M,
gsl_matrix *  M_adj 
)

Get matrix adjoint of given matrix (using gsl)

References colcat(), and gsl_vector_cross_product3D().

Referenced by get_trans_matrix(), and getRot().

Here is the call graph for this function:

Here is the caller graph for this function:

int get_trans_matrix ( vector< double > &  datapoints_base,
vector< double > &  datapoints_new,
vector< vector< double > > &  rotmat,
vector< double > &  p_vec 
)

Takes in 3-d position data of a set of points in 2 poses and calculates the rotation matrix and translation vector.

Parameters:
[in]datapoints_baseThe co-ordinates of the points in the base frame
[in]datapoints_newThe co-ordinates of the points in the new frame
[out]rotmatThe 2-dimensional C++ Vector representing the Rotation matrix
[out]p_vecThe 1-dimensional C++ Vector representing the translation vector

This function is recommended over calc_trans_mat() in case we have position data for more than 3 points, and errors in position data are significant. It uses method of Veldpaus, Woltring & Dortmans (1988) to calculate the best transformation matrix relating thigh marker set to the tibial marker set.

References get_betas(), get_distri_matrices(), get_matrix_adjoin(), getavgvec(), getRot(), gsl_det(), and vect2double().

Here is the call graph for this function:

void getavgvec ( vector< double > &  datavec,
vector< double > &  avgvec,
int  np 
)

Get mean position vector.

Referenced by get_trans_matrix().

Here is the caller graph for this function:

void getRot ( gsl_matrix *  R,
gsl_matrix *  G,
gsl_matrix *  GtG,
double  b1,
double  b2 
)

References get_matrix_adjoin(), and gsl_inverse().

Referenced by get_trans_matrix().

Here is the call graph for this function:

Here is the caller graph for this function:

int sol_betas_df ( gsl_vector *  betas,
void *  g_vals,
gsl_matrix *  J 
)

Jacobian of set of equations to be solved using least squares.

Jacobian is obtained as the matrix:

\[ J = \left[\begin{array}{cc} 2\beta_1 & -2 \\ -2g_3 & 2\beta_2 \end{array}\right] \]

Referenced by get_betas().

Here is the caller graph for this function:

int sol_betas_f ( gsl_vector *  betas,
void *  g_vals,
gsl_vector *  err 
)

Equations to be solved using least squares.

Equations are :

\[ \beta_1^2 - 2\beta_2 = g_1^2 \]

\[ \beta_2^2 -2\beta_1 g_3 = g_2^2 \]

Referenced by get_betas().

Here is the caller graph for this function:

int sol_betas_fdf ( gsl_vector *  betas,
void *  g_vals,
gsl_vector *  err,
gsl_matrix *  J 
)

Equation (and its Jacobian) to be solved using least squares.

Referenced by get_betas().

Here is the caller graph for this function:

 All Data Structures Files Functions Variables Typedefs Defines