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"
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) |
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.
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().
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:
where = datapoints in the new pose, = average position vector for the points in new pose, = datapoints in the old pose, = average position vector for the points in old pose.
[out] | G | The matrix relating point orientation in old pose to new pose |
[out] | P | The matrix describing point distribution, in new pose |
[in] | datapoints_new | Vector 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_base | Vector 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().
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().
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.
[in] | datapoints_base | The co-ordinates of the points in the base frame |
[in] | datapoints_new | The co-ordinates of the points in the new frame |
[out] | rotmat | The 2-dimensional C++ Vector representing the Rotation matrix |
[out] | p_vec | The 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().
void getavgvec | ( | vector< double > & | datavec, |
vector< double > & | avgvec, | ||
int | np | ||
) |
Get mean position vector.
Referenced by get_trans_matrix().
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().
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:
Referenced by get_betas().
int sol_betas_f | ( | gsl_vector * | betas, |
void * | g_vals, | ||
gsl_vector * | err | ||
) |
Equations to be solved using least squares.
Equations are :
Referenced by get_betas().
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().