MNE-CPP 0.1.9
A Framework for Electrophysiology
Loading...
Searching...
No Matches
detecttrigger.cpp
Go to the documentation of this file.
1//=============================================================================================================
35//=============================================================================================================
36// INCLUDES
37//=============================================================================================================
38
39#include "detecttrigger.h"
40
41//=============================================================================================================
42// QT INCLUDES
43//=============================================================================================================
44
45#include <QMapIterator>
46
47//=============================================================================================================
48// USED NAMESPACES
49//=============================================================================================================
50
51using namespace Eigen;
52using namespace RTPROCESSINGLIB;
53
54//=============================================================================================================
55// DEFINE RTPROCESSINGLIB GLOBAL METHODS
56//=============================================================================================================
57
58QList<MatrixXi> RTPROCESSINGLIB::toEventMatrix(QMap<int,QList<QPair<int,double> > > mapTriggers)
59{
60 QList<MatrixXi> lMatDetectedTrigger;
61
62 QMapIterator<int,QList<QPair<int,double> > > idx(mapTriggers);
63
64 while (idx.hasNext()) {
65 idx.next();
66 MatrixXi matDetectedTrigger(idx.value().size(),3);
67
68 for(int i = 0; i < idx.value().size(); ++i) {
69 matDetectedTrigger(i,0) = idx.value().at(i).first;
70 matDetectedTrigger(i,1) = 0;
71 matDetectedTrigger(i,2) = idx.value().at(i).second;
72 }
73
74 lMatDetectedTrigger << matDetectedTrigger;
75 }
76
77 return lMatDetectedTrigger;
78}
79
80//=============================================================================================================
81
82QMap<int,QList<QPair<int,double> > > RTPROCESSINGLIB::detectTriggerFlanksMax(const MatrixXd &data,
83 const QList<int>& lTriggerChannels,
84 int iOffsetIndex,
85 double dThreshold,
86 bool bRemoveOffset,
87 int iBurstLengthSamp)
88{
89 QMap<int,QList<QPair<int,double> > > qMapDetectedTrigger;
90
91 //Find all triggers above threshold in the data block
92 for(int i = 0; i < lTriggerChannels.size(); ++i) {
93 int iChIdx = lTriggerChannels.at(i);
94
95 //Add empty list to map
96 QList<QPair<int,double> > temp;
97 qMapDetectedTrigger.insert(iChIdx, temp);
98
99 //detect the actual triggers in the current data matrix
100 if(iChIdx > data.rows() || iChIdx < 0) {
101 return qMapDetectedTrigger;
102 }
103
104 //Find positive maximum in data vector.
105 for(int j = 0; j < data.cols(); ++j) {
106 double dMatVal = bRemoveOffset ? data(iChIdx,j) - data(iChIdx,0) : data(iChIdx,j);
107
108 if(dMatVal >= dThreshold) {
109 QPair<int,double> pair;
110 pair.first = iOffsetIndex+j;
111 pair.second = data(iChIdx,j);
112
113 qMapDetectedTrigger[iChIdx].append(pair);
114
115 j += iBurstLengthSamp;
116 }
117 }
118 }
119
120 return qMapDetectedTrigger;
121}
122
123//=============================================================================================================
124
125QList<QPair<int,double> > RTPROCESSINGLIB::detectTriggerFlanksMax(const MatrixXd &data,
126 int iTriggerChannelIdx,
127 int iOffsetIndex,
128 double dThreshold,
129 bool bRemoveOffset,
130 int iBurstLengthSamp)
131{
132 QList<QPair<int,double> > lDetectedTriggers;
133
134 //Find all triggers above threshold in the data block
135 //detect the actual triggers in the current data matrix
136 if(iTriggerChannelIdx > data.rows() || iTriggerChannelIdx < 0) {
137 return lDetectedTriggers;
138 }
139
140 //Find positive maximum in data vector.
141 for(int j = 0; j < data.cols(); ++j) {
142 double dMatVal = bRemoveOffset ? data(iTriggerChannelIdx,j) - data(iTriggerChannelIdx,0) : data(iTriggerChannelIdx,j);
143
144 if(dMatVal >= dThreshold) {
145 QPair<int,double> pair;
146 pair.first = iOffsetIndex+j;
147 pair.second = data(iTriggerChannelIdx,j);
148
149 lDetectedTriggers.append(pair);
150
151 j += iBurstLengthSamp;
152 }
153 }
154
155 return lDetectedTriggers;
156}
157
158//=============================================================================================================
159
160QMap<int,QList<QPair<int,double> > > RTPROCESSINGLIB::detectTriggerFlanksGrad(const MatrixXd& data,
161 const QList<int>& lTriggerChannels,
162 int iOffsetIndex,
163 double dThreshold,
164 bool bRemoveOffset,
165 const QString& type,
166 int iBurstLengthSamp)
167{
168 QMap<int,QList<QPair<int,double> > > qMapDetectedTrigger;
169 RowVectorXd tGradient = RowVectorXd::Zero(data.cols());
170
171 //Find all triggers above threshold in the data block
172 for(int i = 0; i < lTriggerChannels.size(); ++i) {
173 int iChIdx = lTriggerChannels.at(i);
174
175 //Add empty list to map
176 QList<QPair<int,double> > temp;
177 qMapDetectedTrigger.insert(iChIdx, temp);
178
179 //detect the actual triggers in the current data matrix
180 if(iChIdx > data.rows() || iChIdx < 0) {
181 return qMapDetectedTrigger;
182 }
183
184 //Compute gradient
185 for(int t = 1; t<tGradient.cols(); t++) {
186 tGradient(t) = data(iChIdx,t)-data(iChIdx,t-1);
187 }
188
189 // If falling flanks are to be detected flip the gradient's sign
190 if(type == "Falling") {
191 tGradient = tGradient * -1;
192 }
193
194 //Find positive maximum in gradient vector. This position is equal to the rising trigger flank.
195 for(int j = 0; j < tGradient.cols(); ++j) {
196 double dMatVal = bRemoveOffset ? tGradient(j) - data(iChIdx,0) : tGradient(j);
197
198 if(dMatVal >= dThreshold) {
199 QPair<int,double> pair;
200 pair.first = iOffsetIndex+j;
201 pair.second = tGradient(j);
202
203 qMapDetectedTrigger[iChIdx].append(pair);
204
205 j += iBurstLengthSamp;
206 }
207 }
208 }
209
210 return qMapDetectedTrigger;
211}
212
213//=============================================================================================================
214
215QList<QPair<int,double> > RTPROCESSINGLIB::detectTriggerFlanksGrad(const MatrixXd &data,
216 int iTriggerChannelIdx,
217 int iOffsetIndex,
218 double dThreshold,
219 bool bRemoveOffset,
220 const QString& type,
221 int iBurstLengthSamp)
222{
223 QList<QPair<int,double> > lDetectedTriggers;
224
225 RowVectorXd tGradient = RowVectorXd::Zero(data.cols());
226
227 //detect the actual triggers in the current data matrix
228 if(iTriggerChannelIdx > data.rows() || iTriggerChannelIdx < 0) {
229 return lDetectedTriggers;
230 }
231
232 //Compute gradient
233 for(int t = 1; t < tGradient.cols(); ++t) {
234 tGradient(t) = data(iTriggerChannelIdx,t) - data(iTriggerChannelIdx,t-1);
235 }
236
237 //If falling flanks are to be detected flip the gradient's sign
238 if(type == "Falling") {
239 tGradient = tGradient * -1;
240 }
241
242 //Find all triggers above threshold in the data block
243 for(int j = 0; j < tGradient.cols(); ++j) {
244 double dMatVal = bRemoveOffset ? tGradient(j) - data(iTriggerChannelIdx,0) : tGradient(j);
245
246 if(dMatVal >= dThreshold) {
247 QPair<int,double> pair;
248 pair.first = iOffsetIndex+j;
249 pair.second = tGradient(j);
250
251 lDetectedTriggers.append(pair);
252
253 j += iBurstLengthSamp;
254 }
255 }
256
257 return lDetectedTriggers;
258}
DetectTrigger declarations.
RTPROCESINGSHARED_EXPORT QMap< int, QList< QPair< int, double > > > detectTriggerFlanksMax(const Eigen::MatrixXd &data, const QList< int > &lTriggerChannels, int iOffsetIndex, double dThreshold, bool bRemoveOffset, int iBurstLengthSamp=100)
RTPROCESINGSHARED_EXPORT QList< Eigen::MatrixXi > toEventMatrix(QMap< int, QList< QPair< int, double > > > mapTriggers)
RTPROCESINGSHARED_EXPORT QMap< int, QList< QPair< int, double > > > detectTriggerFlanksGrad(const Eigen::MatrixXd &data, const QList< int > &lTriggerChannels, int iOffsetIndex, double dThreshold, bool bRemoveOffset, const QString &type, int iBurstLengthSamp=100)