ATSuite C++  v1.0
Scientific C++ routines originally developed by Alexis Tantet
Classes | Typedefs | Functions
transferOperator.hpp File Reference

Calculate discretized approximation of transfer operators from time series. More...

#include <iostream>
#include <vector>
#include <gsl/gsl_vector.h>
#include <gsl/gsl_vector_uint.h>
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_matrix_uint.h>
#include <Eigen/Dense>
#include <Eigen/Sparse>
#include <ATSuite/atmatrix.hpp>
#include <ATSuite/atmarkov.hpp>
#include <ATSuite/atio.hpp>

Go to the source code of this file.

Classes

class  Grid
 Grid class. More...
 
class  transferOperator
 Transfer operator class. More...
 

Typedefs

typedef Eigen::Triplet< double > triplet
 Eigen triplet of double.
 
typedef std::vector< triplettripletVector
 STD vector of Eigen triplet of double.
 
typedef Eigen::Triplet< size_t > tripletUInt
 Eigen triplet of integer.
 
typedef std::vector< tripletUInttripletUIntVector
 STD vector of Eigen triplet of integer.
 
typedef Eigen::SparseMatrix< double, Eigen::ColMajor > SpMatCSC
 Eigen CSC matrix of double type.
 
typedef Eigen::SparseMatrix< double, Eigen::RowMajor > SpMatCSR
 Eigen CSR matrix of double type.
 

Functions

gsl_matrix_uint * getGridMemMatrix (const gsl_matrix *, const gsl_matrix *, const Grid *)
 Get membership matrix from initial and final states for a grid. More...
 
gsl_matrix_uint * getGridMemMatrix (const gsl_matrix *, const Grid *, const size_t)
 Get the grid membership vector from a single long trajectory. More...
 
gsl_vector_uint * getGridMemVector (const gsl_matrix *, const Grid *)
 Get the grid membership matrix from a single long trajectory. More...
 
gsl_matrix_uint * memVector2memMatrix (const gsl_vector_uint *, size_t)
 Get the grid membership matrix from the membership vector for a given lag.
 
gsl_matrix_uint * memVectorList2memMatrix (const std::vector< gsl_vector_uint * > *, size_t)
 Concatenate a list of membership vectors into one membership matrix. More...
 
int getBoxMembership (const gsl_vector *, const Grid *)
 Get membership to a grid box of a single realization. More...
 
tripletUIntVectorgetTransitionCountTriplet (const gsl_matrix_uint *, size_t)
 Get triplet vector from membership matrix. More...
 
gsl_matrix_uint * memVector2memMatrix (gsl_vector_uint *gridMemVect, size_t tauStep)
 

Detailed Description

Calculate discretized approximation of transfer operators from time series. The result is given as forward and backward Markov transition matrices and their distributions.

Definition in file transferOperator.hpp.

Function Documentation

int getBoxMembership ( const gsl_vector *  state,
const Grid grid 
)

Get membership to a grid box of a single realization.

Parameters
[in]stateGSL vector of a single state.
[in]gridPointer to Grid object.
Returns
Box index to which the state belongs.

Definition at line 716 of file transferOperator.hpp.

References Grid::dim, Grid::gridBounds, and Grid::N.

Referenced by getGridMemMatrix(), getGridMemVector(), and transferOperator::transferOperator().

gsl_matrix_uint * getGridMemMatrix ( const gsl_matrix *  initStates,
const gsl_matrix *  finalStates,
const Grid grid 
)

Get membership matrix from initial and final states for a grid.

Parameters
[in]initStatesGSL matrix of initial states.
[in]finalStatesGSL matrix of final states.
[in]gridPointer to Grid object.
Returns
GSL grid membership matrix.

Definition at line 563 of file transferOperator.hpp.

References Grid::dim, getBoxMembership(), and Grid::N.

Referenced by transferOperator::transferOperator().

gsl_matrix_uint * getGridMemMatrix ( const gsl_matrix *  states,
const Grid grid,
const size_t  tauStep 
)

Get the grid membership matrix from a single long trajectory.

Parameters
[in]statesGSL matrix of states for each time step.
[in]gridPointer to Grid object.
[in]tauStepLag used to calculate the transitions.
Returns
GSL grid membership matrix.

Definition at line 658 of file transferOperator.hpp.

References getGridMemVector(), and memVector2memMatrix().

gsl_vector_uint * getGridMemVector ( const gsl_matrix *  states,
const Grid grid 
)

Get the grid membership vector from a single long trajectory.

Parameters
[in]statesGSL matrix of states for each time step.
[in]gridPointer to Grid object.
Returns
GSL grid membership vector.

Definition at line 601 of file transferOperator.hpp.

References Grid::dim, and getBoxMembership().

Referenced by getGridMemMatrix(), and transferOperator::transferOperator().

tripletUIntVector * getTransitionCountTriplet ( const gsl_matrix_uint *  gridMem,
size_t  N 
)

Get the triplet vector counting the transitions from pairs of grid boxes from the grid membership matrix.

Parameters
[in]gridMemGrid membership matrix.
[in]NSize of the grid.
Returns
Triplet vector counting the transitions.

Definition at line 757 of file transferOperator.hpp.

Referenced by transferOperator::transferOperator(), and Grid::~Grid().

gsl_matrix_uint* memVector2memMatrix ( gsl_vector_uint *  gridMemVect,
size_t  tauStep 
)

Get the grid membership matrix from the membership vector for a given lag.

Parameters
[in]gridMemVectGrid membership vector of a long trajectory for a grid.
[in]tauStepLag used to calculate the transitions.
Returns
GSL grid membership matrix.

Definition at line 634 of file transferOperator.hpp.

gsl_matrix_uint * memVectorList2memMatrix ( const std::vector< gsl_vector_uint * > *  memList,
size_t  tauStep 
)

Concatenate a list of membership vectors into one membership matrix.

Parameters
[in]memListSTD vector of membership GSL vectors each of them associated with a single long trajectory.
[in]tauStepLag used to calculate the transitions.
Returns
GSL grid membership matrix.

Definition at line 680 of file transferOperator.hpp.

References memVector2memMatrix().

Referenced by transferOperator::transferOperator().