DMP_BBO library
Public Member Functions | Static Public Member Functions | Protected Member Functions | Friends | List of all members

Model parameters for the GMR function approximator. More...

#include <ModelParametersGMR.hpp>

Inheritance diagram for ModelParametersGMR:
Inheritance graph
[legend]
Collaboration diagram for ModelParametersGMR:
Collaboration graph
[legend]

Public Member Functions

 ModelParametersGMR (std::vector< double > priors, std::vector< Eigen::VectorXd > means, std::vector< Eigen::MatrixXd > covars, int n_output_dims=1)
 Constructor for the model parameters of the GMR function approximator.
 
 ModelParametersGMR (std::vector< double > priors, std::vector< Eigen::VectorXd > mu_xs, std::vector< Eigen::VectorXd > mu_ys, std::vector< Eigen::MatrixXd > sigma_xs, std::vector< Eigen::MatrixXd > sigma_ys, std::vector< Eigen::MatrixXd > sigma_x_ys)
 Constructor for the model parameters of the GMR function approximator.
 
 ModelParametersGMR (int n_observations, std::vector< double > priors, std::vector< Eigen::VectorXd > means, std::vector< Eigen::MatrixXd > covars, int n_output_dims=1)
 Constructor for the model parameters of the GMR function approximator (Used by the incremental training).
 
 ModelParametersGMR (int n_observations, std::vector< double > priors, std::vector< Eigen::VectorXd > mu_xs, std::vector< Eigen::VectorXd > mu_ys, std::vector< Eigen::MatrixXd > sigma_xs, std::vector< Eigen::MatrixXd > sigma_ys, std::vector< Eigen::MatrixXd > sigma_x_ys)
 Constructor for the model parameters of the GMR function approximator (Used by the incremental training).
 
unsigned int getNumberOfGaussians (void) const
 Get the number of Gaussians in the GMM. More...
 
virtual int getExpectedOutputDim (void) const
 The expected dimensionality of the output data. More...
 
virtual int getExpectedInputDim (void) const
 The expected dimensionality of the input data. More...
 
std::string toString (void) const
 Returns a string representation of the object. More...
 
ModelParametersclone (void) const
 Return a pointer to a deep copy of the ModelParameters object. More...
 
void getSelectableParameters (std::set< std::string > &selected_values_labels) const
 Return all the names of the parameter types that can be selected. More...
 
void getParameterVectorMask (const std::set< std::string > selected_values_labels, Eigen::VectorXi &selected_mask) const
 Get a mask for selecting parameters. More...
 
void getParameterVectorAll (Eigen::VectorXd &all_values) const
 Return a vector that returns all available parameter values. More...
 
int getParameterVectorAllSize (void) const
 Get the size of the parameter values vector when it contains all available parameter values. More...
 
bool saveGMM (std::string directory, bool overwrite=false) const
 Save a Gaussian mixture model to a directory; useful for debugging. More...
 
void toMatrix (Eigen::MatrixXd &gmm_as_matrix) const
 This function represents the Gaussian Mixture Model as one big matrix. More...
 
bool saveGMMToMatrix (std::string filename, bool overwrite=false) const
 Save the GMM as a matrix in an ASCII file. More...
 
UnifiedModeltoUnifiedModel (void) const
 Convert these model parameters to unified model parameters. More...
 
- Public Member Functions inherited from ModelParameters
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 Serialize class data members to boost archive. More...
 
- Public Member Functions inherited from Parameterizable
virtual ~Parameterizable (void)
 Destructor.
 
virtual int getParameterVectorSelectedSize (void) const
 Get the size of the vector of selected parameters, as returned by getParameterVectorSelected(. More...
 
virtual void getParameterVectorSelected (Eigen::VectorXd &values, bool normalized=false) const
 Get the values of the selected parameters in one vector. More...
 
virtual void getParameterVectorSelectedNormalized (Eigen::VectorXd &values) const
 Get the normalized values of the selected parameters in one vector. More...
 
void getParameterVectorSelectedMinMax (Eigen::VectorXd &min, Eigen::VectorXd &max) const
 Get the minimum and maximum of the selected parameters in one vector. More...
 
void getParameterVectorAllMinMax (Eigen::VectorXd &min, Eigen::VectorXd &max) const
 Get the minimum and maximum values of the current parameter vector. More...
 
void getParameterVectorSelectedRanges (Eigen::VectorXd &ranges) const
 Get the ranges of the selected parameters, i.e. More...
 
virtual void setParameterVectorSelected (const Eigen::VectorXd &values, bool normalized=false)
 Set all the values of the selected parameters with one vector. More...
 
virtual void setParameterVectorSelectedNormalized (const Eigen::VectorXd &values)
 Set all the values of the selected parameters with one vector of normalized values. More...
 
virtual void setSelectedParameters (const std::set< std::string > &selected_values_labels)
 Determine which subset of parameters is represented in the vector returned by Parameterizable::getParameterVectorSelected. More...
 
void setSelectedParametersOne (std::string selected)
 Set the parameters that are currently selected. More...
 
void setParameterVectorModifier (std::string modifier, bool new_value)
 Turn certain modifiers on or off. More...
 
void setVectorLengthsPerDimension (const Eigen::VectorXi &lengths_per_dimension)
 The vector (VectorXd) with parameter values can be split into different parts (as vector<VectorXd>; this function specifices the length of each sub-vector. More...
 
Eigen::VectorXi getVectorLengthsPerDimension (void) const
 Get the specified length of each vector in each dimension. More...
 
void getParameterVectorSelected (std::vector< Eigen::VectorXd > &values, bool normalized=false) const
 Get the values of the selected parameters in one vector. More...
 
void setParameterVectorSelected (const std::vector< Eigen::VectorXd > &values, bool normalized=false)
 Set all the values of the selected parameters with a vector of vectors. More...
 

Static Public Member Functions

static bool saveGMM (std::string directory, const std::vector< Eigen::VectorXd > &centers, const std::vector< Eigen::MatrixXd > &covars, bool overwrite=false, int iter=-1)
 Save a Gaussian mixture model to a directory; useful for debugging. More...
 
static ModelParametersGMRfromMatrix (const Eigen::MatrixXd &gmm_matrix)
 Initialize a GMM from a matrix. More...
 
static ModelParametersGMRloadGMMFromMatrix (std::string filename)
 Load the GMM from a matrix in an ASCII file. More...
 

Protected Member Functions

void setParameterVectorAll (const Eigen::VectorXd &values)
 Set all available parameter values with one vector. More...
 

Friends

class FunctionApproximatorGMR
 
class boost::serialization::access
 Give boost serialization access to private members. More...
 

Detailed Description

Model parameters for the GMR function approximator.

Definition at line 38 of file ModelParametersGMR.hpp.

Member Function Documentation

unsigned int getNumberOfGaussians ( void  ) const
inline

Get the number of Gaussians in the GMM.

Returns
The number of Gaussians in the GMM.

Definition at line 73 of file ModelParametersGMR.hpp.

74  {
75  return priors_.size();
76  }
virtual int getExpectedOutputDim ( void  ) const
inlinevirtual

The expected dimensionality of the output data.

For now, we only consider 1-dimensional output by default.

Returns
Expected dimensionality of the output data

Reimplemented from ModelParameters.

Definition at line 78 of file ModelParametersGMR.hpp.

79  {
80  assert(means_y_.size()>0); // This is also checked in the constructor
81  return means_y_[0].size();
82  }
virtual int getExpectedInputDim ( void  ) const
inlinevirtual

The expected dimensionality of the input data.

Returns
Expected dimensionality of the input data

Implements ModelParameters.

Definition at line 84 of file ModelParametersGMR.hpp.

85  {
86  assert(means_x_.size()>0); // This is also checked in the constructor
87  return means_x_[0].size();
88  };

Here is the call graph for this function:

string toString ( void  ) const
virtual

Returns a string representation of the object.

Returns
A string representation of the object.

Implements ModelParameters.

Definition at line 236 of file ModelParametersGMR.cpp.

237 {
238  RETURN_STRING_FROM_BOOST_SERIALIZATION_XML("ModelParametersGMR");
239 };
#define RETURN_STRING_FROM_BOOST_SERIALIZATION_XML(name)
Macro to convert the boost XML serialization of an object into a string.
ModelParameters * clone ( void  ) const
virtual

Return a pointer to a deep copy of the ModelParameters object.

Returns
Pointer to a deep copy

Implements ModelParameters.

Definition at line 196 of file ModelParametersGMR.cpp.

197 {
198  std::vector<double> priors;
199  std::vector<VectorXd> means_x;
200  std::vector<VectorXd> means_y;
201  std::vector<MatrixXd> covars_x;
202  std::vector<MatrixXd> covars_y;
203  std::vector<MatrixXd> covars_y_x;
204  int n_observations = n_observations_;
205 
206  for (size_t i = 0; i < priors_.size(); i++)
207  {
208  priors.push_back(priors_[i]);
209  means_x.push_back(VectorXd(means_x_[i]));
210  means_y.push_back(VectorXd(means_y_[i]));
211  covars_x.push_back(MatrixXd(covars_x_[i]));
212  covars_y.push_back(MatrixXd(covars_y_[i]));
213  covars_y_x.push_back(MatrixXd(covars_y_x_[i]));
214  }
215 
216  return new ModelParametersGMR(n_observations, priors, means_x, means_y, covars_x, covars_y, covars_y_x);
217 }
ModelParametersGMR(std::vector< double > priors, std::vector< Eigen::VectorXd > means, std::vector< Eigen::MatrixXd > covars, int n_output_dims=1)
Constructor for the model parameters of the GMR function approximator.
void getSelectableParameters ( std::set< std::string > &  selected_values_labels) const
virtual

Return all the names of the parameter types that can be selected.

Parameters
[out]selected_values_labelsNames of the parameter types that can be selected
Todo:
Determine which parameters should be modifiable in GMR.

Implements Parameterizable.

Definition at line 403 of file ModelParametersGMR.cpp.

404 {
405  selected_values_labels = set<string>();
406  // selected_values_labels.insert("centers");
407  // selected_values_labels.insert("priors");
408  // selected_values_labels.insert("slopes");
409  // selected_values_labels.insert("biases");
410  // selected_values_labels.insert("inverse_covars_l");
411 }
void getParameterVectorMask ( const std::set< std::string >  selected_values_labels,
Eigen::VectorXi &  selected_mask 
) const
virtual

Get a mask for selecting parameters.

Parameters
[in]selected_values_labelsLabels of the selected parameter values
[out]selected_maskA mask indicating indices of selected parameters. 0 indicates not selected, >0 indicates selected.

For instance, if the parameters consists of centers, widths and slopes the parameter values vector will be something like

    centers     widths    slopes
[ 100 110 120 10 10 10 0.4 0.7 0.4 ]

In this case, if selected_values_labels contains "centers" and "slopes", the mask will be:

    centers     widths    slopes
[   1   1   1  0  0  0   3   3   3 ]

The '0' indicates that these parameters are not selected. The other ones have different numbers so that they may be discerned from one another (as required in Parameterizable::getParameterVectorSelectedMinMax for instance.

Implements Parameterizable.

Definition at line 413 of file ModelParametersGMR.cpp.

414 {
415  // selected_mask.resize(getParameterVectorAllSize());
416  // selected_mask.fill(0);
417 
418  // int offset = 0;
419  // int size;
420 
421  // size = centers_.size() * centers_[0].size();
422  // if (selected_values_labels.find("centers")!=selected_values_labels.end())
423  // selected_mask.segment(offset,size).fill(1);
424  // offset += size;
425 
426  // size = priors_.size();
427  // if (selected_values_labels.find("priors")!=selected_values_labels.end())
428  // selected_mask.segment(offset,size).fill(2);
429  // offset += size;
430 
431  // size = slopes_.size() * slopes_[0].rows() * slopes_[0].cols();
432  // if (selected_values_labels.find("slopes")!=selected_values_labels.end())
433  // selected_mask.segment(offset,size).fill(3);
434  // offset += size;
435 
436 
437  // size = biases_.size() * biases_[0].size();
438  // if (selected_values_labels.find("biases")!=selected_values_labels.end())
439  // selected_mask.segment(offset,size).fill(4);
440  // offset += size;
441 
442  // size = inverseCovarsL_.size() * (inverseCovarsL_[0].rows() * (inverseCovarsL_[0].cols() + 1))/2;
443  // if (selected_values_labels.find("inverse_covars_l")!=selected_values_labels.end())
444  // selected_mask.segment(offset,size).fill(5);
445  // offset += size;
446 
447  // assert(offset == getParameterVectorAllSize());
448 }
void getParameterVectorAll ( Eigen::VectorXd &  values) const
virtual

Return a vector that returns all available parameter values.

Parameters
[out]valuesAll available parameter values in one vector.
Remarks
Contrast this with Parameterizable::getParameterVectorSelected, which return only the SELECTED parameter values. Selecting parameters is done with Parameterizable::setSelectedParameters

Implements Parameterizable.

Definition at line 451 of file ModelParametersGMR.cpp.

452 {
453  // values.resize(getParameterVectorAllSize());
454  // int offset = 0;
455 
456  // for (size_t i = 0; i < centers_.size(); i++)
457  // {
458  // values.segment(offset, centers_[i].size()) = centers_[i];
459  // offset += centers_[i].size();
460  // }
461 
462  // for (size_t i = 0; i < centers_.size(); i++)
463  // {
464  // values[offset] = priors_[i];
465  // offset += 1;
466  // }
467 
468  // for (size_t i = 0; i < slopes_.size(); i++)
469  // {
470  // for (int col = 0; col < slopes_[i].cols(); col++)
471  // {
472  // values.segment(offset, slopes_[i].rows()) = slopes_[i].col(col);
473  // offset += slopes_[i].rows();
474  // }
475  // }
476 
477  // for (size_t i = 0; i < centers_.size(); i++)
478  // {
479  // values.segment(offset, biases_[i].size()) = biases_[i];
480  // offset += biases_[i].size();
481  // }
482 
483  // for (size_t i = 0; i < inverseCovarsL_.size(); i++)
484  // {
485  // for (int row = 0; row < inverseCovarsL_[i].rows(); row++)
486  // for (int col = 0; col < row + 1; col++)
487  // {
488  // values[offset] = inverseCovarsL_[i](row, col);
489  // offset += 1;
490  // }
491  // }
492 
493  // assert(offset == getParameterVectorAllSize());
494 
495 };
int getParameterVectorAllSize ( void  ) const
inlinevirtual

Get the size of the parameter values vector when it contains all available parameter values.

Returns
The size of the parameter vector

For instance, if the parameters consists of centers, widths and slopes the parameter values vector will be something like

    centers     widths    slopes
[ 100 110 120 10 10 10 0.4 0.7 0.4 ]

then getParameterVectorAllSize() will return 9

Implements Parameterizable.

Definition at line 102 of file ModelParametersGMR.hpp.

103  {
104  return all_values_vector_size_;
105  }

Here is the call graph for this function:

bool saveGMM ( std::string  directory,
const std::vector< Eigen::VectorXd > &  centers,
const std::vector< Eigen::MatrixXd > &  covars,
bool  overwrite = false,
int  iter = -1 
)
static

Save a Gaussian mixture model to a directory; useful for debugging.

Parameters
[in]directoryDirectory to save to
[in]centersCenters of the Gaussians
[in]covarsCovariance matrices of the Gaussians
[in]overwriteWhether to overwrite existing files or not.
[in]iterIteration number when running Expectation-Maximization. Allows the GMM to be stored with a different filename (in the same directory) at each iteration.
Returns
true if successful, false otherwise

Definition at line 241 of file ModelParametersGMR.cpp.

242 {
243  for (size_t i_gau = 0; i_gau < centers.size(); i_gau++)
244  {
245  stringstream stream;
246  stream << "gmm";
247  if (iter>=0)
248  stream << "_iter" << setw(2) << setfill('0') << iter;
249  stream << "_mu" << setw(3) << setfill('0') << i_gau << ".txt";
250  string filename = stream.str();
251  if (!saveMatrix(directory, filename, centers[i_gau], overwrite))
252  return false;
253  //cout << " filename=" << filename << endl;
254 
255  stringstream stream2;
256  stream2 << "gmm";
257  if (iter>=0)
258  stream2 << "_iter" << setw(2) << setfill('0') << iter;
259  stream2 << "_covar" << setw(3) << setfill('0') << i_gau << ".txt";
260  filename = stream2.str();
261  //cout << " filename=" << filename << endl;
262  if (!saveMatrix(directory, filename, covars[i_gau], overwrite))
263  return false;
264  }
265  return true;
266 }
bool saveMatrix(std::string filename, Eigen::Matrix< Scalar, RowsAtCompileTime, ColsAtCompileTime > matrix, bool overwrite=false)
Save an Eigen matrix to an ASCII file.

Here is the call graph for this function:

bool saveGMM ( std::string  directory,
bool  overwrite = false 
) const

Save a Gaussian mixture model to a directory; useful for debugging.

Parameters
[in]directoryDirectory to save to
[in]overwriteWhether to overwrite existing files or not.
Returns
true if successful, false otherwise

Definition at line 268 of file ModelParametersGMR.cpp.

269 {
270  if (save_directory.empty())
271  return true;
272 
273  //MatrixXd inputs;
274  //FunctionApproximator::generateInputsGrid(min, max, n_samples_per_dim, inputs);
275  //saveMatrix(save_directory,"n_samples_per_dim.txt",n_samples_per_dim,overwrite);
276 
277  int n_gaussians = means_x_.size();
278  int n_dims_in = means_x_[0].size();
279  int n_dims_out = means_y_[0].size();
280  int n_dims_gmm = n_dims_in + n_dims_out;
281 
282 
283  std::vector<VectorXd> means(n_gaussians);
284  std::vector<MatrixXd> covars(n_gaussians);
285  for (int i_gau = 0; i_gau < n_gaussians; i_gau++)
286  {
287  means[i_gau] = VectorXd(n_dims_gmm);
288  means[i_gau].segment(0, n_dims_in) = means_x_[i_gau];
289  means[i_gau].segment(n_dims_in, n_dims_out) = means_y_[i_gau];
290 
291  covars[i_gau] = MatrixXd(n_dims_gmm,n_dims_gmm);
292  covars[i_gau].fill(0);
293  covars[i_gau].block(0, 0, n_dims_in, n_dims_in) = covars_x_[i_gau];
294  covars[i_gau].block(n_dims_in, n_dims_in, n_dims_out, n_dims_out) = covars_y_[i_gau];
295  covars[i_gau].block(n_dims_in, 0, n_dims_out, n_dims_in) = covars_y_x_[i_gau];
296  covars[i_gau].block(0, n_dims_in, n_dims_in, n_dims_out) = covars_y_x_[i_gau].transpose();
297  }
298 
299  saveGMM(save_directory,means,covars,overwrite);
300 
301 
302  return true;
303 }
static bool saveGMM(std::string directory, const std::vector< Eigen::VectorXd > &centers, const std::vector< Eigen::MatrixXd > &covars, bool overwrite=false, int iter=-1)
Save a Gaussian mixture model to a directory; useful for debugging.

Here is the call graph for this function:

void toMatrix ( Eigen::MatrixXd &  gmm_as_matrix) const

This function represents the Gaussian Mixture Model as one big matrix.

Parameters
[out]gmm_as_matrixThe Gaussian Mixture Model as one big matrix

In combination with ModelParametersGMR::saveGMMToMatrix() and ModelParametersGMR::loadGMMFromMatrix(), this hacky function allowed for easier exchange with some Matlab code we wrote.

This matrix contains the following rows (example for two Gaussians, with dimensionality 3)

 0. meta_data  (n_gaussians and n_output_dims)
 1. meta_data  (n_observations)
__________________________
 2. prior1 and E1
 3. mean1      (size: 3)
 4. covar1     (row1 of covar matrix)
 5. covar1     (row2 of covar matrix)
 6. covar1     (row3 of covar matrix)
 __________________________
 7. prior2 and E2
 8. mean2      (size: 3)
 9. covar2     (row1 of covar matrix)
 10. covar2     (row2 of covar matrix)
 11. covar2     (row3 of covar matrix)
 Thus, the number of rows in the Matrix is 2 (for meta-data) + n_gaussians * (1+1+n_dims)

Definition at line 305 of file ModelParametersGMR.cpp.

306 {
307  int n_gaussians = means_x_.size();
308  assert(n_gaussians>0);
309  int n_dims_in = means_x_[0].size();
310  int n_dims_out = means_y_[0].size();
311  int n_dims_gmm = n_dims_in + n_dims_out;
312 
313  int n_rows = 2; // First row contains n_gaussians, n_output_dims, second row contains n_observations and responsability
314  for (int i_gau = 0; i_gau < n_gaussians; i_gau++)
315  {
316  n_rows += 1; // Add one row for the prior
317  n_rows += 1; // Add one row for the mean
318  n_rows += n_dims_gmm; // For the covariance matrix
319  }
320 
321  gmm_as_matrix = MatrixXd::Zero(n_rows,n_dims_gmm);
322 
323  gmm_as_matrix(0,0) = n_gaussians;
324  gmm_as_matrix(0,1) = n_dims_out;
325  gmm_as_matrix(1,0) = n_observations_;
326 
327  VectorXd mean = VectorXd(n_dims_gmm);
328  MatrixXd covar = MatrixXd(n_dims_gmm,n_dims_gmm);
329  int cur_row = 2;
330  for (int i_gau = 0; i_gau < n_gaussians; i_gau++)
331  {
332  mean.segment(0, n_dims_in) = means_x_[i_gau];
333  mean.segment(n_dims_in, n_dims_out) = means_y_[i_gau];
334 
335  covar.block(0, 0, n_dims_in, n_dims_in) = covars_x_[i_gau];
336  covar.block(n_dims_in, n_dims_in, n_dims_out, n_dims_out) = covars_y_[i_gau];
337  covar.block(n_dims_in, 0, n_dims_out, n_dims_in) = covars_y_x_[i_gau];
338  covar.block(0, n_dims_in, n_dims_in, n_dims_out) = covars_y_x_[i_gau].transpose();
339 
340  gmm_as_matrix(cur_row,0) = priors_[i_gau];
341  gmm_as_matrix.row(cur_row+1) = mean;
342  gmm_as_matrix.block(cur_row+2,0,n_dims_gmm,n_dims_gmm) = covar;
343 
344  cur_row += 1 + 1 + n_dims_gmm;
345  }
346 }
ModelParametersGMR * fromMatrix ( const Eigen::MatrixXd &  gmm_matrix)
static

Initialize a GMM from a matrix.

See also
toMatrix() for the format of the file
Parameters
[in]gmm_matrixThe GMM, represented as a matrix.
Returns
The model parameters of the GMM.

Definition at line 348 of file ModelParametersGMR.cpp.

349 {
350  int n_dims_gmm = gmm_matrix.cols();
351  int n_rows = gmm_matrix.rows();
352  assert(n_dims_gmm>1);
353  assert(n_rows>0);
354 
355  int n_gaussians = gmm_matrix(0,0);
356  int n_dims_out = gmm_matrix(0,1);
357  int n_observations = gmm_matrix(1,0);
358 
359  assert(n_rows == (2+ (n_gaussians*(1+1+n_dims_gmm))));
360 
361  vector<double> priors(n_gaussians);
362  vector<VectorXd> means(n_gaussians);
363  vector<MatrixXd> covars(n_gaussians);
364 
365  int cur_row = 2;
366  for (int i_gau = 0; i_gau < n_gaussians; i_gau++)
367  {
368  priors[i_gau] = gmm_matrix(cur_row,0);
369 
370  means[i_gau] = gmm_matrix.row(cur_row+1);
371 
372  covars[i_gau] = gmm_matrix.block(cur_row+2,0,n_dims_gmm,n_dims_gmm);
373 
374  cur_row += 1 + 1 + n_dims_gmm;
375  }
376 
377  return new ModelParametersGMR(n_observations,priors,means,covars,n_dims_out);
378 }
ModelParametersGMR(std::vector< double > priors, std::vector< Eigen::VectorXd > means, std::vector< Eigen::MatrixXd > covars, int n_output_dims=1)
Constructor for the model parameters of the GMR function approximator.
bool saveGMMToMatrix ( std::string  filename,
bool  overwrite = false 
) const

Save the GMM as a matrix in an ASCII file.

See also
toMatrix() for the format of the file
Parameters
[in]filenameThe name of the file to save the GMM to.
[in]overwriteWhether to overwrite the file if it already exists.
Returns
true if saving was successful, false otherwise.

Definition at line 380 of file ModelParametersGMR.cpp.

381 {
382  if (filename.empty())
383  return true;
384 
385  MatrixXd gmm_as_matrix;
386  toMatrix(gmm_as_matrix);
387 
388  if (!saveMatrix(filename, gmm_as_matrix, overwrite))
389  return false;
390 
391  return true;
392 }
bool saveMatrix(std::string filename, Eigen::Matrix< Scalar, RowsAtCompileTime, ColsAtCompileTime > matrix, bool overwrite=false)
Save an Eigen matrix to an ASCII file.
void toMatrix(Eigen::MatrixXd &gmm_as_matrix) const
This function represents the Gaussian Mixture Model as one big matrix.

Here is the call graph for this function:

ModelParametersGMR * loadGMMFromMatrix ( std::string  filename)
static

Load the GMM from a matrix in an ASCII file.

See also
toMatrix() for the format of the file
Parameters
[in]filenameThe name of the file to load the GMM from.
Returns
The model parameters of the GMM.

Definition at line 394 of file ModelParametersGMR.cpp.

395 {
396  MatrixXd gmm_matrix;
397  if (!loadMatrix(filename, gmm_matrix))
398  return NULL;
399 
400  return ModelParametersGMR::fromMatrix(gmm_matrix);
401 }
static ModelParametersGMR * fromMatrix(const Eigen::MatrixXd &gmm_matrix)
Initialize a GMM from a matrix.
bool loadMatrix(std::string filename, Eigen::Matrix< Scalar, RowsAtCompileTime, ColsAtCompileTime > &m)
Load an Eigen matrix from an ASCII file.
Definition: EigenFileIO.tpp:33

Here is the call graph for this function:

UnifiedModel * toUnifiedModel ( void  ) const
virtual

Convert these model parameters to unified model parameters.

Returns
Unified model parameter representation (NULL if not implemented for a particular subclass)

Implements ModelParameters.

Definition at line 548 of file ModelParametersGMR.cpp.

549 {
550  int n_gaussians = means_x_.size();
551 
552  // This copying is not necessary. It is just done to show which variable in GMR relates to which
553  // variable in the unified model.
554  vector<VectorXd> centers = means_x_;
555  vector<MatrixXd> covars = covars_x_;
556 
557  vector<VectorXd> slopes(n_gaussians);
558  vector<double> offsets(n_gaussians);
559  VectorXd rest;
560  for (int i_gau=0; i_gau<n_gaussians; i_gau++)
561  {
562  slopes[i_gau] = (covars_y_x_[i_gau] * covars_x_inv_[i_gau]).transpose();
563 
564  assert(means_y_[i_gau].size()==1); // Only works for 1D y output for now
565  offsets[i_gau] = means_y_[i_gau][0] - slopes[i_gau].dot(means_x_[i_gau]);
566  }
567 
568 
569  bool normalized_basis_functions = true;
570  bool lines_pivot_at_max_activation = false;
571 
572  return new UnifiedModel(centers, covars, slopes, offsets, priors_, normalized_basis_functions,lines_pivot_at_max_activation);
573 
574 }
void setParameterVectorAll ( const Eigen::VectorXd &  values)
protectedvirtual

Set all available parameter values with one vector.

Parameters
[in]valuesAll available parameter values in one vector.
Remarks
Contrast this with Parameterizable::setParameterVectorSelected, which sets only the SELECTED parameter values. Selecting parameters is done with Parameterizable::setSelectedParameters

Implements Parameterizable.

Definition at line 497 of file ModelParametersGMR.cpp.

498 {
499  // if (all_values_vector_size_ != values.size())
500  // {
501  // cerr << __FILE__ << ":" << __LINE__ << ": values is of wrong size." << endl;
502  // return;
503  // }
504 
505  // int offset = 0;
506 
507  // for (size_t i = 0; i < centers_.size(); i++)
508  // {
509  // centers_[i] = values.segment(offset, centers_[i].size());
510  // offset += centers_[i].size();
511  // }
512  // for (size_t i = 0; i < centers_.size(); i++)
513  // {
514  // priors_[i] = values[offset];
515  // offset += 1;
516  // }
517 
518  // for (size_t i = 0; i < slopes_.size(); i++)
519  // {
520  // for (int col = 0; col < slopes_[i].cols(); col++)
521  // {
522  // slopes_[i].col(col) = values.segment(offset, slopes_[i].rows());
523  // offset += slopes_[i].rows();
524  // }
525  // }
526 
527  // for (size_t i = 0; i < centers_.size(); i++)
528  // {
529  // biases_[i] = values.segment(offset, biases_[i].size());
530  // offset += biases_[i].size();
531  // }
532 
533  // for (size_t i = 0; i < inverseCovarsL_.size(); i++)
534  // {
535  // for (int row = 0; row < inverseCovarsL_[i].rows(); row++)
536  // for (int col = 0; col < row + 1; col++)
537  // {
538  // inverseCovarsL_[i](row, col) = values[offset];
539  // offset += 1;
540  // }
541  // }
542 
543  // assert(offset == getParameterVectorAllSize());
544 
545 };

Friends And Related Function Documentation

friend class boost::serialization::access
friend

Give boost serialization access to private members.

Definition at line 214 of file ModelParametersGMR.hpp.


The documentation for this class was generated from the following files: