v2.0.0
Loading...
Searching...
No Matches
interpolation.cpp
Go to the documentation of this file.
1//=============================================================================================================
34
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// STL INCLUDES
50//=============================================================================================================
51
52#include <unordered_set>
53
54//=============================================================================================================
55// USED NAMESPACES
56//=============================================================================================================
57
58using namespace DISP3DRHILIB;
59using namespace Eigen;
60
61//=============================================================================================================
62// DEFINE MEMBER METHODS
63//=============================================================================================================
64
65QSharedPointer<SparseMatrix<float> > Interpolation::createInterpolationMat(const VectorXi &vecProjectedSensors,
66 const QSharedPointer<MatrixXd> matDistanceTable,
67 double (*interpolationFunction) (double),
68 const double dCancelDist,
69 const VectorXi &vecExcludeIndex)
70{
71 if(matDistanceTable->rows() == 0 && matDistanceTable->cols() == 0) {
72 qDebug() << "[WARNING] Interpolation::createInterpolationMat - received an empty distance table.";
73 return QSharedPointer<SparseMatrix<float> >::create();
74 }
75
76 QSharedPointer<Eigen::SparseMatrix<float> > matInterpolationMatrix = QSharedPointer<SparseMatrix<float> >::create(matDistanceTable->rows(), static_cast<int>(vecProjectedSensors.size()));
77
78 QVector<Triplet<float> > vecNonZeroEntries;
79 const qint32 iRows = matInterpolationMatrix->rows();
80 const qint32 iCols = matInterpolationMatrix->cols();
81
82 // Build exclude set for O(1) lookup
83 std::unordered_set<int> excludeSet;
84 for(Eigen::Index i = 0; i < vecExcludeIndex.size(); ++i) {
85 excludeSet.insert(vecExcludeIndex[i]);
86 }
87
88 QSet<qint32> sensorLookup;
89 for(Eigen::Index idx = 0; idx < vecProjectedSensors.size(); ++idx){
90 if(excludeSet.count(static_cast<int>(idx)) == 0){
91 sensorLookup.insert(vecProjectedSensors[idx]);
92 }
93 }
94
95 for (qint32 r = 0; r < iRows; ++r) {
96 if (sensorLookup.contains(r) == false) {
97 QVector<QPair<qint32, float> > vecBelowThresh;
98 float dWeightsSum = 0.0;
99 const RowVectorXd& rowVec = matDistanceTable->row(r);
100
101 for (qint32 c = 0; c < iCols; ++c) {
102 const float dDist = rowVec[c];
103
104 if (dDist < dCancelDist) {
105 const float dValueWeight = std::fabs(1.0 / interpolationFunction(dDist));
106 dWeightsSum += dValueWeight;
107 vecBelowThresh.push_back(qMakePair(c, dValueWeight));
108 }
109 }
110
111 for (const QPair<qint32, float> &qp : vecBelowThresh) {
112 vecNonZeroEntries.push_back(Eigen::Triplet<float> (r, qp.first, qp.second / dWeightsSum));
113 }
114 } else {
115 // Find index of r in vecProjectedSensors
116 int iIndexInSubset = 0;
117 for(Eigen::Index k = 0; k < vecProjectedSensors.size(); ++k) {
118 if(vecProjectedSensors[k] == r) {
119 iIndexInSubset = static_cast<int>(k);
120 break;
121 }
122 }
123 vecNonZeroEntries.push_back(Eigen::Triplet<float> (r, iIndexInSubset, 1));
124 }
125 }
126
127 matInterpolationMatrix->setFromTriplets(vecNonZeroEntries.begin(), vecNonZeroEntries.end());
128
129 return matInterpolationMatrix;
130}
131
132//=============================================================================================================
133
134VectorXf Interpolation::interpolateSignal(const QSharedPointer<SparseMatrix<float> > matInterpolationMatrix,
135 const QSharedPointer<VectorXf> &vecMeasurementData)
136{
137 if (matInterpolationMatrix->cols() != vecMeasurementData->rows()) {
138 qDebug() << "[WARNING] Interpolation::interpolateSignal - Dimension mismatch. Return null pointer...";
139 return VectorXf();
140 }
141
142 VectorXf pOutVec = *matInterpolationMatrix * (*vecMeasurementData);
143 return pOutVec;
144}
145
146//=============================================================================================================
147
148VectorXf Interpolation::interpolateSignal(const SparseMatrix<float> &matInterpolationMatrix,
149 const VectorXf &vecMeasurementData)
150{
151 if (matInterpolationMatrix.cols() != vecMeasurementData.rows()) {
152 qDebug() << "[WARNING] Interpolation::interpolateSignal - Dimension mismatch. Return null pointer...";
153 return VectorXf();
154 }
155
156 VectorXf pOutVec = matInterpolationMatrix * vecMeasurementData;
157 return pOutVec;
158}
159
160//=============================================================================================================
161
162double Interpolation::linear(const double dIn)
163{
164 return dIn;
165}
166
167//=============================================================================================================
168
169double Interpolation::gaussian(const double dIn)
170{
171 return exp(-((dIn * dIn) / 2.0));
172}
173
174//=============================================================================================================
175
176double Interpolation::square(const double dIn)
177{
178 return std::max((-(1.0f / 9.0f) * (dIn * dIn) + 1), 0.0);
179}
180
181//=============================================================================================================
182
183double Interpolation::cubic(const double dIn)
184{
185 return dIn * dIn * dIn;
186}
Interpolation class declaration.
3-D brain visualisation using the Qt RHI rendering backend.
static double gaussian(const double dIn)
gaussian Gaussian interpolation function (sigma=1).
static Eigen::VectorXf interpolateSignal(const QSharedPointer< Eigen::SparseMatrix< float > > matInterpolationMatrix, const QSharedPointer< Eigen::VectorXf > &vecMeasurementData)
interpolateSignal Interpolates sensor data using the weight matrix (shared pointer version).
static double linear(const double dIn)
linear Identity interpolation function.
static double cubic(const double dIn)
cubic Cubic hyperbola interpolation function.
static QSharedPointer< Eigen::SparseMatrix< float > > createInterpolationMat(const Eigen::VectorXi &vecProjectedSensors, const QSharedPointer< Eigen::MatrixXd > matDistanceTable, double(*interpolationFunction)(double), const double dCancelDist=FLOAT_INFINITY, const Eigen::VectorXi &vecExcludeIndex=Eigen::VectorXi())
createInterpolationMat Calculates the weight matrix for interpolation.
static double square(const double dIn)
square Negative parabola interpolation function with y-offset of 1.