v2.0.0
Loading...
Searching...
No Matches
sts_cluster.h
Go to the documentation of this file.
1//=============================================================================================================
34
35#ifndef STS_CLUSTER_H
36#define STS_CLUSTER_H
37
38//=============================================================================================================
39// INCLUDES
40//=============================================================================================================
41
42#include "sts_global.h"
43#include "sts_types.h"
44
45//=============================================================================================================
46// QT INCLUDES
47//=============================================================================================================
48
49#include <QVector>
50#include <QPair>
51
52//=============================================================================================================
53// EIGEN INCLUDES
54//=============================================================================================================
55
56#include <Eigen/Core>
57#include <Eigen/SparseCore>
58
59//=============================================================================================================
60// DEFINE NAMESPACE STSLIB
61//=============================================================================================================
62
63namespace STSLIB
64{
65
66//=============================================================================================================
71 Eigen::MatrixXd matTObs;
72 QVector<double> vecClusterStats;
73 QVector<double> vecClusterPvals;
74 Eigen::MatrixXi matClusterIds;
76};
77
78//=============================================================================================================
85{
86public:
87 //=========================================================================================================
102 const QVector<Eigen::MatrixXd>& dataA,
103 const QVector<Eigen::MatrixXd>& dataB,
104 const Eigen::SparseMatrix<int>& adjacency,
105 int nPermutations = 1024,
106 double clusterAlpha = 0.05,
107 double pThreshold = 0.05,
109
110 //=========================================================================================================
128 const QVector<Eigen::MatrixXd>& data,
129 const Eigen::SparseMatrix<int>& adjacency,
130 double threshold,
131 int nPermutations,
132 StatsTailType tail);
133
134 //=========================================================================================================
151 const QVector<QVector<Eigen::MatrixXd>>& conditions,
152 const Eigen::SparseMatrix<int>& adjacency,
153 double threshold,
154 int nPermutations);
155
156 //=========================================================================================================
173 static Eigen::MatrixXd tfce(
174 const Eigen::MatrixXd& statMap,
175 const Eigen::SparseMatrix<int>& adjacency,
176 double E = 0.5,
177 double H = 2.0,
178 int nSteps = 100);
179
180private:
181 //=========================================================================================================
185 static Eigen::MatrixXd computeTMap(
186 const QVector<Eigen::MatrixXd>& dataA,
187 const QVector<Eigen::MatrixXd>& dataB);
188
189 //=========================================================================================================
195 static QPair<Eigen::MatrixXi, QVector<double>> findClusters(
196 const Eigen::MatrixXd& tMap,
197 double threshold,
198 const Eigen::SparseMatrix<int>& adjacency,
199 StatsTailType tail);
200
201 //=========================================================================================================
205 static double permuteOnce(
206 const QVector<Eigen::MatrixXd>& allData,
207 int nA,
208 const Eigen::SparseMatrix<int>& adjacency,
209 double threshold,
210 StatsTailType tail);
211
212 //=========================================================================================================
216 static double inverseTCdf(double p, int df);
217
218 //=========================================================================================================
222 static Eigen::MatrixXd computeOneSampleTMap(
223 const QVector<Eigen::MatrixXd>& data);
224
225 //=========================================================================================================
229 static Eigen::MatrixXd computeFMap(
230 const QVector<QVector<Eigen::MatrixXd>>& conditions);
231
232 //=========================================================================================================
238 static QPair<Eigen::MatrixXi, QVector<double>> findClustersFlat(
239 const Eigen::MatrixXd& statMap,
240 double threshold,
241 const Eigen::SparseMatrix<int>& adjacency,
242 bool positiveOnly);
243
244 //=========================================================================================================
248 static double permuteOnceOneSample(
249 const QVector<Eigen::MatrixXd>& data,
250 const Eigen::SparseMatrix<int>& adjacency,
251 double threshold,
252 StatsTailType tail);
253
254 //=========================================================================================================
258 static double permuteOnceFTest(
259 const QVector<Eigen::MatrixXd>& allData,
260 const QVector<int>& groupSizes,
261 const Eigen::SparseMatrix<int>& adjacency,
262 double threshold);
263};
264
265} // namespace STSLIB
266
267#endif // STS_CLUSTER_H
Stats library type enumerations.
stats library export/import macros.
#define STSSHARED_EXPORT
Definition sts_global.h:50
Statistical testing (t-tests, F-tests, cluster permutation, multiple comparison correction).
StatsTailType
Definition sts_types.h:49
Eigen::MatrixXd matTObs
Definition sts_cluster.h:71
Eigen::MatrixXi matClusterIds
Definition sts_cluster.h:74
QVector< double > vecClusterStats
Definition sts_cluster.h:72
QVector< double > vecClusterPvals
Definition sts_cluster.h:73
Cluster permutation testing.
Definition sts_cluster.h:85
static StatsClusterResult permutationTest(const QVector< Eigen::MatrixXd > &dataA, const QVector< Eigen::MatrixXd > &dataB, const Eigen::SparseMatrix< int > &adjacency, int nPermutations=1024, double clusterAlpha=0.05, double pThreshold=0.05, StatsTailType tail=StatsTailType::Both)
static Eigen::MatrixXd tfce(const Eigen::MatrixXd &statMap, const Eigen::SparseMatrix< int > &adjacency, double E=0.5, double H=2.0, int nSteps=100)
static StatsClusterResult oneSamplePermutationTest(const QVector< Eigen::MatrixXd > &data, const Eigen::SparseMatrix< int > &adjacency, double threshold, int nPermutations, StatsTailType tail)
static StatsClusterResult fTestPermutationTest(const QVector< QVector< Eigen::MatrixXd > > &conditions, const Eigen::SparseMatrix< int > &adjacency, double threshold, int nPermutations)