MNE-CPP  0.1.9
A Framework for Electrophysiology
interpolation.cpp
Go to the documentation of this file.
1 //=============================================================================================================
35 //=============================================================================================================
36 // INCLUDES
37 //=============================================================================================================
38 
39 #include "interpolation.h"
40 
41 //=============================================================================================================
42 // QT INCLUDES
43 //=============================================================================================================
44 
45 #include <QSet>
46 #include <QDebug>
47 
48 //=============================================================================================================
49 // EIGEN INCLUDES
50 //=============================================================================================================
51 
52 //=============================================================================================================
53 // USED NAMESPACES
54 //=============================================================================================================
55 
56 using namespace DISP3DLIB;
57 using namespace Eigen;
58 
59 //=============================================================================================================
60 // DEFINE GLOBAL METHODS
61 //=============================================================================================================
62 
63 //=============================================================================================================
64 // INITIALIZE STATIC MEMBER
65 //=============================================================================================================
66 
67 //=============================================================================================================
68 // DEFINE MEMBER METHODS
69 //=============================================================================================================
70 
71 QSharedPointer<SparseMatrix<float> > Interpolation::createInterpolationMat(const QVector<int> &vecProjectedSensors,
72  const QSharedPointer<MatrixXd> matDistanceTable,
73  double (*interpolationFunction) (double),
74  const double dCancelDist,
75  const QVector<int> &vecExcludeIndex)
76 {
77 
78  if(matDistanceTable->rows() == 0 && matDistanceTable->cols() == 0) {
79  qDebug() << "[WARNING] Interpolation::createInterpolationMat - received an empty distance table.";
80  return QSharedPointer<SparseMatrix<float> >::create();
81  }
82 
83  // initialization
84  QSharedPointer<Eigen::SparseMatrix<float> > matInterpolationMatrix = QSharedPointer<SparseMatrix<float> >::create(matDistanceTable->rows(), vecProjectedSensors.size());
85 
86  // temporary helper structure for filling sparse matrix
87  QVector<Triplet<float> > vecNonZeroEntries;
88  const qint32 iRows = matInterpolationMatrix->rows();
89  const qint32 iCols = matInterpolationMatrix->cols();
90 
91  // insert all sensor nodes into set for faster lookup during later computation. Also consider bad channels here.
92  QSet<qint32> sensorLookup;
93  int idx = 0;
94 
95  for(const qint32& s : vecProjectedSensors){
96  if(!vecExcludeIndex.contains(idx)){
97  sensorLookup.insert(s);
98  }
99  idx++;
100  }
101 
102  // main loop: go through all rows of distance table and calculate weights
103  for (qint32 r = 0; r < iRows; ++r) {
104  if (sensorLookup.contains(r) == false) {
105  // "normal" node, i.e. one which was not assigned a sensor
106  // bLoThreshold: stores the indizes that point to distances which are below the passed distance threshold (dCancelDist)
107  QVector<QPair<qint32, float> > vecBelowThresh;
108  float dWeightsSum = 0.0;
109  const RowVectorXd& rowVec = matDistanceTable->row(r);
110 
111  for (qint32 c = 0; c < iCols; ++c) {
112  const float dDist = rowVec[c];
113 
114  if (dDist < dCancelDist) {
115  const float dValueWeight = std::fabs(1.0 / interpolationFunction(dDist));
116  dWeightsSum += dValueWeight;
117  vecBelowThresh.push_back(qMakePair(c, dValueWeight));
118  }
119  }
120 
121  for (const QPair<qint32, float> &qp : vecBelowThresh) {
122  vecNonZeroEntries.push_back(Eigen::Triplet<float> (r, qp.first, qp.second / dWeightsSum));
123  }
124  } else {
125  // a sensor has been assigned to this node, we do not need to interpolate anything
126  //(final vertex signal is equal to sensor input signal, thus factor 1)
127  const int iIndexInSubset = vecProjectedSensors.indexOf(r);
128 
129  vecNonZeroEntries.push_back(Eigen::Triplet<float> (r, iIndexInSubset, 1));
130  }
131  }
132 
133  matInterpolationMatrix->setFromTriplets(vecNonZeroEntries.begin(), vecNonZeroEntries.end());
134 
135  return matInterpolationMatrix;
136 }
137 
138 //=============================================================================================================
139 
140 VectorXf Interpolation::interpolateSignal(const QSharedPointer<SparseMatrix<float> > matInterpolationMatrix,
141  const QSharedPointer<VectorXf> &vecMeasurementData)
142 {
143  if (matInterpolationMatrix->cols() != vecMeasurementData->rows()) {
144  qDebug() << "[WARNING] Interpolation::interpolateSignal - Dimension mismatch. Return null pointer...";
145  return VectorXf();
146  }
147 
148  VectorXf pOutVec = *matInterpolationMatrix * (*vecMeasurementData);
149 
150  return pOutVec;
151 }
152 
153 //=============================================================================================================
154 
155 VectorXf Interpolation::interpolateSignal(const SparseMatrix<float> &matInterpolationMatrix,
156  const VectorXf &vecMeasurementData)
157 {
158  if (matInterpolationMatrix.cols() != vecMeasurementData.rows()) {
159  qDebug() << "[WARNING] Interpolation::interpolateSignal - Dimension mismatch. Return null pointer...";
160  return VectorXf();
161  }
162 
163  VectorXf pOutVec = matInterpolationMatrix * vecMeasurementData;
164 
165  return pOutVec;
166 }
167 
168 //=============================================================================================================
169 
170 double Interpolation::linear(const double dIn)
171 {
172  return dIn;
173 }
174 
175 //=============================================================================================================
176 
177 double Interpolation::gaussian(const double dIn)
178 {
179  return exp(-((dIn * dIn) / 2.0));
180 }
181 
182 //=============================================================================================================
183 
184 double Interpolation::square(const double dIn)
185 {
186  return std::max((-(1.0f / 9.0f) * (dIn * dIn) + 1), 0.0);
187 }
188 
189 //=============================================================================================================
190 
191 double Interpolation::cubic(const double dIn)
192 {
193  return dIn * dIn * dIn;
194 }
DISP3DLIB::Interpolation::cubic
static double cubic(const double dIn)
Definition: interpolation.cpp:191
DISP3DLIB::Interpolation::interpolateSignal
static Eigen::VectorXf interpolateSignal(const QSharedPointer< Eigen::SparseMatrix< float > > matInterpolationMatrix, const QSharedPointer< Eigen::VectorXf > &vecMeasurementData)
DISP3DLIB::Interpolation::linear
static double linear(const double dIn)
Definition: interpolation.cpp:170
DISP3DLIB::Interpolation::gaussian
static double gaussian(const double dIn)
Definition: interpolation.cpp:177
DISP3DLIB::Interpolation::square
static double square(const double dIn)
Definition: interpolation.cpp:184
DISP3DLIB::Interpolation::createInterpolationMat
static QSharedPointer< Eigen::SparseMatrix< float > > createInterpolationMat(const QVector< int > &vecProjectedSensors, const QSharedPointer< Eigen::MatrixXd > matDistanceTable, double(*interpolationFunction)(double), const double dCancelDist=FLOAT_INFINITY, const QVector< int > &vecExcludeIndex=QVector< int >())
Definition: interpolation.cpp:71
interpolation.h
Interpolation class declaration.