v2.0.0
Loading...
Searching...
No Matches
artifact_detect.cpp
Go to the documentation of this file.
1//=============================================================================================================
33
34//=============================================================================================================
35// INCLUDES
36//=============================================================================================================
37
38#include "artifact_detect.h"
39#include "iirfilter.h"
40
41//=============================================================================================================
42// FIFF INCLUDES
43//=============================================================================================================
44
45#include <fiff/fiff_constants.h>
46#include <fiff/fiff_ch_info.h>
47
48//=============================================================================================================
49// EIGEN INCLUDES
50//=============================================================================================================
51
52#include <Eigen/Core>
53
54//=============================================================================================================
55// QT INCLUDES
56//=============================================================================================================
57
58#include <QDebug>
59
60//=============================================================================================================
61// C++ INCLUDES
62//=============================================================================================================
63
64#include <cmath>
65#include <algorithm>
66
67//=============================================================================================================
68// USED NAMESPACES
69//=============================================================================================================
70
71using namespace UTILSLIB;
72using namespace FIFFLIB;
73using namespace Eigen;
74
75//=============================================================================================================
76// PRIVATE HELPERS
77//=============================================================================================================
78
79RowVectorXd ArtifactDetect::bandpassFilter(const RowVectorXd& vecSignal,
80 double dSFreq,
81 double dLow,
82 double dHigh,
83 int iOrder)
84{
85 QVector<IirBiquad> sos;
86
87 if (dLow < 1e-3) {
88 // Low-pass only
89 sos = IirFilter::designButterworth(iOrder, IirFilter::LowPass, dHigh, 0.0, dSFreq);
90 } else {
91 sos = IirFilter::designButterworth(iOrder, IirFilter::BandPass, dLow, dHigh, dSFreq);
92 }
93
94 return IirFilter::applyZeroPhase(vecSignal, sos);
95}
96
97//=============================================================================================================
98
99QVector<int> ArtifactDetect::findPeaks(const RowVectorXd& vecSignal,
100 double dThreshold,
101 int iMinDist)
102{
103 QVector<int> peaks;
104 const int N = static_cast<int>(vecSignal.size());
105 if (N < 3) return peaks;
106
107 int lastPeak = -iMinDist - 1;
108
109 for (int i = 1; i < N - 1; ++i) {
110 double v = vecSignal(i);
111 if (v < dThreshold) continue;
112
113 // Local maximum check
114 if (v > vecSignal(i - 1) && v >= vecSignal(i + 1)) {
115 // Enforce minimum distance
116 if (i - lastPeak >= iMinDist) {
117 peaks.append(i);
118 lastPeak = i;
119 } else if (!peaks.isEmpty() && vecSignal(peaks.last()) < v) {
120 // Replace last peak if this one is higher and within the distance window
121 peaks.last() = i;
122 lastPeak = i;
123 }
124 }
125 }
126
127 return peaks;
128}
129
130//=============================================================================================================
131// PUBLIC DEFINITIONS
132//=============================================================================================================
133
134QVector<int> ArtifactDetect::detectEcg(const MatrixXd& matData,
135 const FiffInfo& fiffInfo,
136 double dSFreq,
137 const EcgParams& params)
138{
139 // ---- Find ECG channel ----
140 int ecgIdx = -1;
141 for (int i = 0; i < fiffInfo.nchan; ++i) {
142 if (fiffInfo.chs[i].kind == FIFFV_ECG_CH) {
143 ecgIdx = i;
144 break;
145 }
146 }
147
148 RowVectorXd ecgSignal;
149
150 if (ecgIdx >= 0) {
151 // Use dedicated ECG channel (already calibrated)
152 ecgSignal = matData.row(ecgIdx).cast<double>();
153 } else {
154 // ---- Synthetic ECG from MEG channels ----
155 // The cardiac artifact is visible on most MEG channels. A robust synthetic ECG is
156 // obtained by summing the absolute values of all magnetometer channels (gradiometer
157 // baseline artefact reduces them; magnetometers show the global field pattern).
158 qWarning() << "ArtifactDetect::detectEcg: No ECG channel found. "
159 "Synthesising ECG proxy from MEG magnetometers.";
160
161 QVector<int> magIdx;
162 for (int i = 0; i < fiffInfo.nchan; ++i) {
163 if (fiffInfo.chs[i].kind == FIFFV_MEG_CH) {
164 // Distinguish magnetometers from gradiometers by coil type:
165 // magnetometers have a single integration point — heuristically identified
166 // as FIFFV_COIL_MAG type; use channel unit as a proxy (T vs T/m).
167 // Use unit: FIFF_UNIT_T = 112, FIFF_UNIT_T_M = 201
168 if (fiffInfo.chs[i].unit == 112) { // Tesla — magnetometer
169 magIdx.append(i);
170 }
171 }
172 }
173
174 if (magIdx.isEmpty()) {
175 // Fall back: any MEG channel
176 for (int i = 0; i < fiffInfo.nchan; ++i) {
177 if (fiffInfo.chs[i].kind == FIFFV_MEG_CH) {
178 magIdx.append(i);
179 }
180 }
181 }
182
183 if (magIdx.isEmpty()) {
184 qWarning() << "ArtifactDetect::detectEcg: No MEG channels found. Cannot detect ECG.";
185 return {};
186 }
187
188 const int nSamp = static_cast<int>(matData.cols());
189 ecgSignal = RowVectorXd::Zero(nSamp);
190 for (int idx : magIdx) {
191 ecgSignal += matData.row(idx).cwiseAbs().cast<double>();
192 }
193 ecgSignal /= static_cast<double>(magIdx.size());
194 }
195
196 // ---- Band-pass filter to isolate QRS complex ----
197 RowVectorXd filtered = bandpassFilter(ecgSignal, dSFreq,
198 params.dFilterLow, params.dFilterHigh,
199 params.iFilterOrder);
200
201 // ---- Adaptive threshold: fraction of the peak-to-peak amplitude ----
202 double sigMin = filtered.minCoeff();
203 double sigMax = filtered.maxCoeff();
204 double threshold = params.dThreshFactor * (sigMax - sigMin) + sigMin;
205
206 // Minimum inter-peak distance in samples
207 int iMinDist = static_cast<int>(std::round(params.dMinRRSec * dSFreq));
208 iMinDist = std::max(iMinDist, 1);
209
210 return findPeaks(filtered, threshold, iMinDist);
211}
212
213//=============================================================================================================
214
215QVector<int> ArtifactDetect::detectEog(const MatrixXd& matData,
216 const FiffInfo& fiffInfo,
217 double dSFreq,
218 const EogParams& params)
219{
220 // ---- Find EOG channel(s) ----
221 QVector<int> eogIdx;
222 for (int i = 0; i < fiffInfo.nchan; ++i) {
223 if (fiffInfo.chs[i].kind == FIFFV_EOG_CH) {
224 eogIdx.append(i);
225 }
226 }
227
228 if (eogIdx.isEmpty()) {
229 qWarning() << "ArtifactDetect::detectEog: No EOG channel found.";
230 return {};
231 }
232
233 // Select the EOG channel with the largest peak-to-peak amplitude
234 int bestIdx = eogIdx[0];
235 double bestPtp = 0.0;
236 for (int idx : eogIdx) {
237 RowVectorXd row = matData.row(idx).cast<double>();
238 double ptp = row.maxCoeff() - row.minCoeff();
239 if (ptp > bestPtp) {
240 bestPtp = ptp;
241 bestIdx = idx;
242 }
243 }
244
245 RowVectorXd eogSignal = matData.row(bestIdx).cast<double>();
246
247 // ---- Low-pass filter ----
248 RowVectorXd filtered = bandpassFilter(eogSignal, dSFreq,
249 0.0, params.dFilterHigh,
250 params.iFilterOrder);
251
252 // ---- Detect excursions above ±threshold ----
253 // Work on the absolute value so both positive (downward blinks, depending on EOG polarity)
254 // and negative deflections are detected.
255 RowVectorXd absFiltered = filtered.cwiseAbs();
256
257 int iMinDist = static_cast<int>(std::round(params.dMinGapSec * dSFreq));
258 iMinDist = std::max(iMinDist, 1);
259
260 return findPeaks(absFiltered, params.dThresholdV, iMinDist);
261}
Fiff constants.
#define FIFFV_EOG_CH
#define FIFFV_MEG_CH
#define FIFFV_ECG_CH
FiffChInfo class declaration.
Declaration of IirFilter — Butterworth IIR filter design and application using numerically stable sec...
Declaration of ArtifactDetect — ECG and EOG physiological artifact event detection.
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
static QVector< int > detectEcg(const Eigen::MatrixXd &matData, const FIFFLIB::FiffInfo &fiffInfo, double dSFreq, const EcgParams &params=EcgParams())
ArtifactDetectEogParams EogParams
ArtifactDetectEcgParams EcgParams
static QVector< int > detectEog(const Eigen::MatrixXd &matData, const FIFFLIB::FiffInfo &fiffInfo, double dSFreq, const EogParams &params=EogParams())
static Eigen::RowVectorXd applyZeroPhase(const Eigen::RowVectorXd &vecData, const QVector< IirBiquad > &sos)
static QVector< IirBiquad > designButterworth(int iOrder, FilterType type, double dCutoffLow, double dCutoffHigh, double dSFreq)
FIFF measurement file information.
Definition fiff_info.h:86
QList< FiffChInfo > chs