v2.0.0
Loading...
Searching...
No Matches
dipole_fit.cpp
Go to the documentation of this file.
1
2
3#include "dipole_fit.h"
5#include "guess_data.h"
6
7#include <string.h>
8#include <memory>
9
10using namespace INVERSELIB;
11using namespace MNELIB;
12using namespace FWDLIB;
13
14#if defined(_WIN32) || defined(_WIN64)
15#define snprintf _snprintf
16#define vsnprintf _vsnprintf
17#define strcasecmp _stricmp
18#define strncasecmp _strnicmp
19#endif
20
21//=============================================================================================================
22// STATIC DEFINITIONS
23//=============================================================================================================
24
25#ifndef PROGRAM_VERSION
26#define PROGRAM_VERSION "1.00"
27#endif
28
29#ifndef FAIL
30#define FAIL -1
31#endif
32
33#ifndef OK
34#define OK 0
35#endif
36
37#ifndef TRUE
38#define TRUE 1
39#endif
40
41#ifndef FALSE
42#define FALSE 0
43#endif
44
45#define BIG_TIME 1e6
46
47#define SEG_LEN 10.0
48
49#define EPS_VALUES 0.05
50
51//=============================================================================================================
52// STATIC DEFINITIONS ToDo make members
53//=============================================================================================================
54
55#define MALLOC(x,t) (t *)malloc((x)*sizeof(t))
56#define REALLOC(x,y,t) (t *)((x == NULL) ? malloc((y)*sizeof(t)) : realloc((x),(y)*sizeof(t)))
57#define FREE(x) if ((char *)(x) != NULL) free((char *)(x))
58
59#define ALLOC_CMATRIX(x,y) mne_cmatrix((x),(y))
60
61#define FREE_CMATRIX(m) mne_free_cmatrix((m))
62
63static void matrix_error(int kind, int nr, int nc)
64
65{
66 if (kind == 1)
67 printf("Failed to allocate memory pointers for a %d x %d matrix\n",nr,nc);
68 else if (kind == 2)
69 printf("Failed to allocate memory for a %d x %d matrix\n",nr,nc);
70 else
71 printf("Allocation error for a %d x %d matrix\n",nr,nc);
72 if (sizeof(void *) == 4) {
73 printf("This is probably because you seem to be using a computer with 32-bit architecture.\n");
74 printf("Please consider moving to a 64-bit platform.");
75 }
76 printf("Cannot continue. Sorry.\n");
77 exit(1);
78}
79
80float **mne_cmatrix (int nr,int nc)
81
82{
83 int i;
84 float **m;
85 float *whole;
86
87 m = MALLOC(nr,float *);
88 if (!m) matrix_error(1,nr,nc);
89 whole = MALLOC(nr*nc,float);
90 if (!whole) matrix_error(2,nr,nc);
91
92 for(i=0;i<nr;i++)
93 m[i] = whole + i*nc;
94 return m;
95}
96
97void mne_free_cmatrix (float **m)
98{
99 if (m) {
100 FREE(*m);
101 FREE(m);
102 }
103}
104
105//char *mne_compose_mne_name(const char *path, const char *filename)
107// * Compose a filename under the "$MNE_ROOT" directory
108// */
109//{
110// char *res;
111// char *mne_root;
112
113// if (filename == NULL) {
114// qCritical("No file name specified to mne_compose_mne_name");
115// return NULL;
116// }
117// mne_root = getenv(MNE_ENV_ROOT);
118// if (mne_root == NULL || strlen(mne_root) == 0) {
119// qCritical("Environment variable MNE_ROOT not set");
120// return NULL;
121// }
122// if (path == NULL || strlen(path) == 0) {
123// res = MALLOC(strlen(mne_root)+strlen(filename)+2,char);
124// strcpy(res,mne_root);
125// strcat(res,"/");
126// strcat(res,filename);
127// }
128// else {
129// res = MALLOC(strlen(mne_root)+strlen(filename)+strlen(path)+3,char);
130// strcpy(res,mne_root);
131// strcat(res,"/");
132// strcat(res,path);
133// strcat(res,"/");
134// strcat(res,filename);
135// }
136// return res;
137//}
138
139//============================= misc_util.c =============================
140
141//============================= mne_named_matrix.c =============================
142
143void mne_free_name_list(char **list, int nlist)
144/*
145 * Free a name list array
146 */
147{
148 int k;
149 if (list == NULL || nlist == 0)
150 return;
151 for (k = 0; k < nlist; k++) {
152#ifdef FOO
153 printf("%d %s\n",k,list[k]);
154#endif
155 FREE(list[k]);
156 }
157 FREE(list);
158 return;
159}
160
161//============================= mne_ch_selections.c =============================
162
163/*
164 * Mandatory allocation functions
165 */
166static mneChSelection new_ch_selection()
167
168{
169 mneChSelection newsel = new MNEChSelection();
170
171 newsel->name.clear();
172 newsel->chdef.clear();
173 newsel->chspick.clear();
174 newsel->chspick_nospace.clear();
175 newsel->ndef = 0;
176 newsel->nchan = 0;
178 return newsel;
179}
180
181mneChSelection mne_ch_selection_these(const QString& selname, const QStringList& names, int nch)
182/*
183 * Give an explicit list of interesting channels
184 */
185{
186 int c;
187 mneChSelection sel;
188
189 sel = new_ch_selection();
190 sel->name = selname;
191// sel->chdef;
192 sel->ndef = nch;
194
195 for (c = 0; c < nch; c++)
196 sel->chdef.append(names[c]);
197
198 return sel;
199}
200
201static void omit_spaces(QStringList names, int nnames)
202
203{
204 char *c,*cc;
205 for (int k = 0; k < names.size(); k++) {
206 while( names[k].startsWith(" ") ) {
207 names[k].remove(0,1);
208 }
209 }
210 return;
211}
212
214 MNERawData* data)
215/*
216 * Make the channel picking real easy
217 */
218{
219 int c,rc,d;
220 MNERawInfo* info;
221 int nch;
222 QString dash;
223
224 if (!sel || !data)
225 return 0;
226
227 info = data->info.get();
228 sel->chspick.clear();
229 sel->chspick_nospace.clear();
230 /*
231 * Expansion of possible regular expressions must be added eventually
232 */
233 sel->chspick = sel->chdef;
234 sel->chspick_nospace = sel->chdef;
235 omit_spaces(sel->chspick_nospace,sel->ndef);
236 sel->nchan = sel->ndef;
237
238 sel->pick.setConstant(sel->nchan, -1);
239 sel->pick_deriv.setConstant(sel->nchan, -1);
240 sel->ch_kind.setConstant(sel->nchan, -1);
241
242 for (c = 0; c < sel->nchan; c++) {
243 for (rc = 0; rc < info->nchan; rc++) {
244 if (QString::compare(sel->chspick[c],info->chInfo[rc].ch_name,Qt::CaseInsensitive) == 0 ||
245 QString::compare(sel->chspick_nospace[c],info->chInfo[rc].ch_name,Qt::CaseInsensitive) == 0) {
246 sel->pick[c] = rc;
247 sel->ch_kind[c] = info->chInfo[rc].kind;
248 break;
249 }
250 }
251 }
252 /*
253 * Maybe the derivations will help
254 */
255 sel->nderiv = 0;
256 if (data->deriv_matched) {
257 QStringList deriv_names = data->deriv_matched->deriv_data->rowlist;
258 int nderiv = data->deriv_matched->deriv_data->nrow;
259
260 for (c = 0; c < sel->nchan; c++) {
261 if (sel->pick[c] == -1) {
262 for (d = 0; d < nderiv; d++) {
263 if (QString::compare(sel->chspick[c],deriv_names[d],Qt::CaseInsensitive) == 0 &&
264 data->deriv_matched->valid.size() > 0 && data->deriv_matched->valid[d]) {
265 sel->pick_deriv[c] = d;
266 sel->ch_kind[c] = data->deriv_matched->chs[d].kind;
267 sel->nderiv++;
268 break;
269 }
270 }
271 }
272 }
273 }
274 /*
275 * Try simple channels again without the part after dashes
276 */
277 for (c = 0; c < sel->nchan; c++) {
278 if (sel->pick[c] == -1 && sel->pick_deriv[c] == -1) {
279 for (rc = 0; rc < info->nchan; rc++) {
280 dash = QString(info->chInfo[rc].ch_name).mid(QString(info->chInfo[rc].ch_name).indexOf("-")+1);// strchr(info->chInfo[rc].ch_name,'-');
281 if (!dash.isNull()) {
282// *dash = '\0';
283 if (QString::compare(sel->chspick[c],info->chInfo[rc].ch_name,Qt::CaseInsensitive) == 0 ||
284 QString::compare(sel->chspick_nospace[c],info->chInfo[rc].ch_name,Qt::CaseInsensitive) == 0) {
285// *dash = '-';
286 sel->pick[c] = rc;
287 sel->ch_kind[c] = info->chInfo[rc].kind;
288 break;
289 }
290// *dash = '-';
291 }
292 }
293 }
294 }
295 for (c = 0, nch = 0; c < sel->nchan; c++) {
296 if (sel->pick[c] >= 0)
297 nch++;
298 }
299 if (sel->nderiv > 0)
300 printf("Selection %c%s%c has %d matched derived channels.\n",'"',sel->name.toUtf8().constData(),'"',sel->nderiv);
301 return nch;
302}
303
304//============================= mne_get_values.c =============================
305
306int mne_get_values_from_data (float time, /* Interesting time point */
307 float integ, /* Time integration */
308 float **data, /* The data values (time by time) */
309 int nsamp, /* How many time points? */
310 int nch, /* How many channels */
311 float tmin, /* Time of first sample */
312 float sfreq, /* Sampling frequency */
313 int use_abs, /* Use absolute values */
314 float *value) /* The picked values */
315/*
316 * Pick a signal value using linear interpolation
317 */
318{
319 int n1,n2,k;
320 float s1,s2;
321 float f1,f2;
322 float sum;
323 float width;
324 int ch;
325
326 for (ch = 0; ch < nch; ch++) {
327 /*
328 * Find out the correct samples
329 */
330 if (std::fabs(sfreq*integ) < EPS_VALUES) { /* This is the single-sample case */
331 s1 = sfreq*(time - tmin);
332 n1 = floor(s1);
333 f1 = 1.0 + n1 - s1;
334 if (n1 < 0 || n1 > nsamp-1) {
335 printf("Sample value out of range %d (0..%d)",n1,nsamp-1);
336 return(-1);
337 }
338 /*
339 * Avoid rounding error
340 */
341 if (n1 == nsamp-1) {
342 if (std::fabs(f1-1.0) < 1e-3)
343 f1 = 1.0;
344 }
345 if (f1 < 1.0 && n1 > nsamp-2) {
346 printf("Sample value out of range %d (0..%d) %.4f",n1,nsamp-1,f1);
347 return(-1);
348 }
349 if (f1 < 1.0) {
350 if (use_abs)
351 sum = f1*std::fabs(data[n1][ch]) + (1.0-f1)*std::fabs(data[n1+1][ch]);
352 else
353 sum = f1*data[n1][ch] + (1.0-f1)*data[n1+1][ch];
354 }
355 else {
356 if (use_abs)
357 sum = std::fabs(data[n1][ch]);
358 else
359 sum = data[n1][ch];
360 }
361 }
362 else { /* Multiple samples */
363 s1 = sfreq*(time - 0.5*integ - tmin);
364 s2 = sfreq*(time + 0.5*integ - tmin);
365 n1 = ceil(s1); n2 = floor(s2);
366 if (n2 < n1) { /* We are within one sample interval */
367 n1 = floor(s1);
368 if (n1 < 0 || n1 > nsamp-2)
369 return (-1);
370 f1 = s1 - n1;
371 f2 = s2 - n1;
372 if (use_abs)
373 sum = 0.5*((f1+f2)*std::fabs(data[n1+1][ch]) + (2.0-f1-f2)*std::fabs(data[n1][ch]));
374 else
375 sum = 0.5*((f1+f2)*data[n1+1][ch] + (2.0-f1-f2)*data[n1][ch]);
376 }
377 else {
378 f1 = n1 - s1;
379 f2 = s2 - n2;
380 if (n1 < 0 || n1 > nsamp-1) {
381 printf("Sample value out of range %d (0..%d)",n1,nsamp-1);
382 return(-1);
383 }
384 if (n2 < 0 || n2 > nsamp-1) {
385 printf("Sample value out of range %d (0..%d)",n2,nsamp-1);
386 return(-1);
387 }
388 if (f1 != 0.0 && n1 < 1)
389 return(-1);
390 if (f2 != 0.0 && n2 > nsamp-2)
391 return(-1);
392 sum = 0.0;
393 width = 0.0;
394 if (n2 > n1) { /* Do the whole intervals */
395 if (use_abs) {
396 sum = 0.5 * std::fabs(data[n1][ch]);
397 for (k = n1+1; k < n2; k++)
398 sum = sum + std::fabs(data[k][ch]);
399 sum = sum + 0.5 * std::fabs(data[n2][ch]);
400 }
401 else {
402 sum = 0.5*data[n1][ch];
403 for (k = n1+1; k < n2; k++)
404 sum = sum + data[k][ch];
405 sum = sum + 0.5*data[n2][ch];
406 }
407 width = n2 - n1;
408 }
409 /*
410 * Add tails
411 */
412 if (use_abs) {
413 if (f1 != 0.0)
414 sum = sum + 0.5 * f1 * (f1 * std::fabs(data[n1-1][ch]) + (2.0 - f1) * std::fabs(data[n1][ch]));
415 if (f2 != 0.0)
416 sum = sum + 0.5 * f2 * (f2 * std::fabs(data[n2+1][ch]) + (2.0 - f2) * std::fabs(data[n2][ch]));
417 }
418 else {
419 if (f1 != 0.0)
420 sum = sum + 0.5*f1*(f1*data[n1-1][ch] + (2.0-f1)*data[n1][ch]);
421 if (f2 != 0.0)
422 sum = sum + 0.5*f2*(f2*data[n2+1][ch] + (2.0-f2)*data[n2][ch]);
423 }
424 width = width + f1 + f2;
425 sum = sum/width;
426 }
427 }
428 value[ch] = sum;
429 }
430 return (0);
431}
432
433int mne_get_values_from_data_ch (float time, /* Interesting time point */
434 float integ, /* Time integration */
435 float **data, /* The data values (channel by channel) */
436 int nsamp, /* How many time points? */
437 int nch, /* How many channels */
438 float tmin, /* Time of first sample */
439 float sfreq, /* Sampling frequency */
440 int use_abs, /* Use absolute values */
441 float *value) /* The picked values */
442/*
443 * Pick a signal value using linear interpolation
444 */
445{
446 int n1,n2,k;
447 float s1,s2;
448 float f1,f2;
449 float sum;
450 float width;
451 int ch;
452
453 for (ch = 0; ch < nch; ch++) {
454 /*
455 * Find out the correct samples
456 */
457 if (std::fabs(sfreq * integ) < EPS_VALUES) { /* This is the single-sample case */
458 s1 = sfreq*(time - tmin);
459 n1 = floor(s1);
460 f1 = 1.0 + n1 - s1;
461 if (n1 < 0 || n1 > nsamp-1)
462 return(-1);
463 if (f1 < 1.0 && n1 > nsamp-2)
464 return(-1);
465 if (f1 < 1.0) {
466 if (use_abs)
467 sum = f1 * std::fabs(data[ch][n1]) + (1.0 - f1) * std::fabs(data[ch][n1+1]);
468 else
469 sum = f1*data[ch][n1] + (1.0-f1)*data[ch][n1+1];
470 }
471 else {
472 if (use_abs)
473 sum = std::fabs(data[ch][n1]);
474 else
475 sum = data[ch][n1];
476 }
477 }
478 else { /* Multiple samples */
479 s1 = sfreq*(time - 0.5*integ - tmin);
480 s2 = sfreq*(time + 0.5*integ - tmin);
481 n1 = ceil(s1); n2 = floor(s2);
482 if (n2 < n1) { /* We are within one sample interval */
483 n1 = floor(s1);
484 if (n1 < 0 || n1 > nsamp-2)
485 return (-1);
486 f1 = s1 - n1;
487 f2 = s2 - n1;
488 if (use_abs)
489 sum = 0.5*((f1+f2)*std::fabs(data[ch][n1+1]) + (2.0-f1-f2)*std::fabs(data[ch][n1]));
490 else
491 sum = 0.5*((f1+f2)*data[ch][n1+1] + (2.0-f1-f2)*data[ch][n1]);
492 }
493 else {
494 f1 = n1 - s1;
495 f2 = s2 - n2;
496 if (n1 < 0 || n1 > nsamp-1 || n2 < 0 || n2 > nsamp-1)
497 return(-1);
498 if (f1 != 0.0 && n1 < 1)
499 return(-1);
500 if (f2 != 0.0 && n2 > nsamp-2)
501 return(-1);
502 sum = 0.0;
503 width = 0.0;
504 if (n2 > n1) { /* Do the whole intervals */
505 if (use_abs) {
506 sum = 0.5 * std::fabs(data[ch][n1]);
507 for (k = n1+1; k < n2; k++)
508 sum = sum + std::fabs(data[ch][k]);
509 sum = sum + 0.5 * std::fabs(data[ch][n2]);
510 }
511 else {
512 sum = 0.5*data[ch][n1];
513 for (k = n1+1; k < n2; k++)
514 sum = sum + data[ch][k];
515 sum = sum + 0.5*data[ch][n2];
516 }
517 width = n2 - n1;
518 }
519 /*
520 * Add tails
521 */
522 if (use_abs) {
523 if (f1 != 0.0)
524 sum = sum + 0.5 * f1 * (f1 * std::fabs(data[ch][n1-1]) + (2.0 - f1) * std::fabs(data[ch][n1]));
525 if (f2 != 0.0)
526 sum = sum + 0.5 * f2 * (f2 * std::fabs(data[ch][n2+1]) + (2.0 - f2) * std::fabs(data[ch][n2]));
527 }
528 else {
529 if (f1 != 0.0)
530 sum = sum + 0.5*f1*(f1*data[ch][n1-1]+ (2.0-f1)*data[ch][n1]);
531 if (f2 != 0.0)
532 sum = sum + 0.5*f2*(f2*data[ch][n2+1] + (2.0-f2)*data[ch][n2]);
533 }
534 width = width + f1 + f2;
535 sum = sum/width;
536 }
537 }
538 value[ch] = sum;
539 }
540 return (0);
541}
542
543//=============================================================================================================
544// DEFINE MEMBER METHODS
545//=============================================================================================================
546
548: settings(p_settings)
549{
550}
551
552//=============================================================================================================
553//todo split in initFit where the settings are handed over and the actual fit
555{
556 std::unique_ptr<GuessData> guess;
557 ECDSet set;
558 FwdEegSphereModel* eeg_model = NULL;
559 DipoleFitData* fit_data = NULL;
560 MNEMeasData* data = NULL;
561 MNERawData* raw = NULL;
562 mneChSelection sel = NULL;
563
564 printf("---- Setting up...\n\n");
565 if (settings->include_eeg) {
566 if ((eeg_model = FwdEegSphereModel::setup_eeg_sphere_model(settings->eeg_model_file,settings->eeg_model_name,settings->eeg_sphere_rad)) == NULL)
567 goto out;
568 }
569
570 if ((fit_data = DipoleFitData::setup_dipole_fit_data( settings->mriname,
571 settings->measname,
572 settings->bemname.isEmpty() ? NULL : settings->bemname.toUtf8().data(),
573 &settings->r0,
574 eeg_model,
575 settings->accurate,
576 settings->badname,
577 settings->noisename,
578 settings->grad_std,
579 settings->mag_std,
580 settings->eeg_std,
581 settings->mag_reg,
582 settings->grad_reg,
583 settings->eeg_reg,
584 settings->diagnoise,
585 settings->projnames,
586 settings->include_meg,
587 settings->include_eeg)) == NULL)
588 goto out;
589
590 fit_data->fit_mag_dipoles = settings->fit_mag_dipoles;
591 if (settings->is_raw) {
592 int c;
593 float t1,t2;
594
595 printf("\n---- Opening a raw data file...\n\n");
596 if ((raw = MNERawData::open_file(settings->measname.isEmpty() ? NULL : settings->measname.toUtf8().data(),TRUE,FALSE,settings->filter)) == NULL)
597 goto out;
598 /*
599 * A channel selection is needed to access the data
600 */
601 sel = mne_ch_selection_these("fit",fit_data->ch_names,fit_data->nmeg+fit_data->neeg);
603 for (c = 0; c < sel->nchan; c++)
604 if (sel->pick[c] < 0) {
605 qCritical ("All desired channels were not available");
606 goto out;
607 }
608 printf("\tChannel selection created.\n");
609 /*
610 * Let's be a little generous here
611 */
612 t1 = raw->first_samp/raw->info->sfreq;
613 t2 = (raw->first_samp+raw->nsamp-1)/raw->info->sfreq;
614 if (settings->tmin < t1 + settings->integ)
615 settings->tmin = t1 + settings->integ;
616 if (settings->tmax > t2 - settings->integ)
617 settings->tmax = t2 - settings->integ;
618 if (settings->tstep < 0)
619 settings->tstep = 1.0/raw->info->sfreq;
620
621 printf("\tOpened raw data file %s : %d MEG and %d EEG \n",
622 settings->measname.toUtf8().data(),fit_data->nmeg,fit_data->neeg);
623 }
624 else {
625 printf("\n---- Reading data...\n\n");
626 if ((data = MNEMeasData::mne_read_meas_data(settings->measname,
627 settings->setno,
628 NULL,
629 NULL,
630 fit_data->ch_names,
631 fit_data->nmeg+fit_data->neeg)) == NULL)
632 goto out;
633 if (settings->do_baseline)
634 data->adjust_baselines(settings->bmin,settings->bmax);
635 else
636 printf("\tNo baseline setting in effect.\n");
637 if (settings->tmin < data->current->tmin + settings->integ/2.0)
638 settings->tmin = data->current->tmin + settings->integ/2.0;
639 if (settings->tmax > data->current->tmin + (data->current->np-1)*data->current->tstep - settings->integ/2.0)
640 settings->tmax = data->current->tmin + (data->current->np-1)*data->current->tstep - settings->integ/2.0;
641 if (settings->tstep < 0)
642 settings->tstep = data->current->tstep;
643
644 printf("\tRead data set %d from %s : %d MEG and %d EEG \n",
645 settings->setno,settings->measname.toUtf8().data(),fit_data->nmeg,fit_data->neeg);
646 if (!settings->noisename.isEmpty()) {
647 printf("\nScaling the noise covariance...\n");
648 if (DipoleFitData::scale_noise_cov(fit_data,data->current->nave) == FAIL)
649 goto out;
650 }
651 }
652
653 /*
654 * Proceed to computing the fits
655 */
656 printf("\n---- Computing the forward solution for the guesses...\n\n");
657 guess.reset(new GuessData( settings->guessname,
658 settings->guess_surfname,
659 settings->guess_mindist, settings->guess_exclude, settings->guess_grid, fit_data));
660 if (!guess)
661 goto out;
662
663 fprintf (stderr,"\n---- Fitting : %7.1f ... %7.1f ms (step: %6.1f ms integ: %6.1f ms)\n\n",
664 1000*settings->tmin,1000*settings->tmax,1000*settings->tstep,1000*settings->integ);
665
666 if (raw) {
667 if (fit_dipoles_raw(settings->measname,raw,sel,fit_data,guess.release(),settings->tmin,settings->tmax,settings->tstep,settings->integ,settings->verbose) == FAIL)
668 goto out;
669 }
670 else {
671 if (fit_dipoles(settings->measname,data,fit_data,guess.release(),settings->tmin,settings->tmax,settings->tstep,settings->integ,settings->verbose,set) == FAIL)
672 goto out;
673 }
674 printf("%d dipoles fitted\n",set.size());
675
676out : {
677 return set;
678 }
679}
680
681//=============================================================================================================
682
683int DipoleFit::fit_dipoles( const QString& dataname, MNEMeasData* data, DipoleFitData* fit, GuessData* guess, float tmin, float tmax, float tstep, float integ, int verbose, ECDSet& p_set)
684{
685 float *one = MALLOC(data->nchan,float);
686 float time;
687 ECDSet set;
688 ECD dip;
689 int s;
690 int report_interval = 10;
691
692 set.dataname = dataname;
693
694 printf("Fitting...%c",verbose ? '\n' : '\0');
695 for (s = 0, time = tmin; time < tmax; s++, time = tmin + s*tstep) {
696 /*
697 * Pick the data point
698 */
699 if (mne_get_values_from_data(time,integ,data->current->data,data->current->np,data->nchan,data->current->tmin,
700 1.0/data->current->tstep,FALSE,one) == FAIL) {
701 printf("Cannot pick time: %7.1f ms\n",1000*time);
702 continue;
703 }
704
705 if (!DipoleFitData::fit_one(fit,guess,time,one,verbose,dip))
706 printf("t = %7.1f ms : %s\n",1000*time,"error (tbd: catch)");
707 else {
708 set.addEcd(dip);
709 if (verbose)
710 dip.print(stdout);
711 else {
712 if (set.size() % report_interval == 0)
713 printf("%d..",set.size());
714 }
715 }
716 }
717 if (!verbose)
718 printf("[done]\n");
719 FREE(one);
720 p_set = set;
721 return OK;
722}
723
724//=============================================================================================================
725
726int DipoleFit::fit_dipoles_raw(const QString& dataname, MNERawData* raw, mneChSelection sel, DipoleFitData* fit, GuessData* guess, float tmin, float tmax, float tstep, float integ, int verbose, ECDSet& p_set)
727{
728 float *one = MALLOC(sel->nchan,float);
729 float sfreq = raw->info->sfreq;
730 float myinteg = integ > 0.0 ? 2*integ : 0.1;
731 int overlap = ceil(myinteg*sfreq);
732 int length = SEG_LEN*sfreq;
733 int step = length - overlap;
734 int stepo = step + overlap/2;
735 int start = raw->first_samp;
736 int s,picks;
737 float time,stime;
738 float **data = ALLOC_CMATRIX(sel->nchan,length);
739 ECD dip;
740 ECDSet set;
741 int report_interval = 10;
742
743 set.dataname = dataname;
744
745 /*
746 * Load the initial data segment
747 */
748 stime = start/sfreq;
749 if (raw->pick_data_filt(sel,start,length,data) == FAIL)
750 goto bad;
751 printf("Fitting...%c",verbose ? '\n' : '\0');
752 for (s = 0, time = tmin; time < tmax; s++, time = tmin + s*tstep) {
753 picks = time*sfreq - start;
754 if (picks > stepo) { /* Need a new data segment? */
755 start = start + step;
756 if (raw->pick_data_filt(sel,start,length,data) == FAIL)
757 goto bad;
758 picks = time*sfreq - start;
759 stime = start/sfreq;
760 }
761 /*
762 * Get the values
763 */
764 if (mne_get_values_from_data_ch (time,integ,data,length,sel->nchan,stime,sfreq,FALSE,one) == FAIL) {
765 printf("Cannot pick time: %8.3f s\n",time);
766 continue;
767 }
768 /*
769 * Fit
770 */
771 if (!DipoleFitData::fit_one(fit,guess,time,one,verbose,dip))
772 qWarning() << "Error";
773 else {
774 set.addEcd(dip);
775 if (verbose)
776 dip.print(stdout);
777 else {
778 if (set.size() % report_interval == 0)
779 printf("%d..",set.size());
780 }
781 }
782 }
783 if (!verbose)
784 printf("[done]\n");
785 FREE_CMATRIX(data);
786 FREE(one);
787 p_set = set;
788 return OK;
789
790bad : {
791 FREE_CMATRIX(data);
792 FREE(one);
793 return FAIL;
794 }
795}
796
797//=============================================================================================================
798
799int DipoleFit::fit_dipoles_raw(const QString& dataname, MNERawData* raw, mneChSelection sel, DipoleFitData* fit, GuessData* guess, float tmin, float tmax, float tstep, float integ, int verbose)
800{
801 ECDSet set;
802 return fit_dipoles_raw(dataname, raw, sel, fit, guess, tmin, tmax, tstep, integ, verbose, set);
803}
#define MNE_CH_SELECTION_USER
Definition mne_types.h:119
#define MNE_CH_SELECTION_UNKNOWN
Definition mne_types.h:117
#define FREE(x)
#define TRUE
#define FALSE
#define OK
#define FAIL
MNE Meas Data Set (MNEMeasDataSet) class declaration.
GuessData class declaration.
Dipole Fit class declaration.
#define FREE(x)
#define EPS_VALUES
int mne_get_values_from_data(float time, float integ, float **data, int nsamp, int nch, float tmin, float sfreq, int use_abs, float *value)
int mne_ch_selection_assign_chs(mneChSelection sel, MNERawData *data)
mneChSelection mne_ch_selection_these(const QString &selname, const QStringList &names, int nch)
int mne_get_values_from_data_ch(float time, float integ, float **data, int nsamp, int nch, float tmin, float sfreq, int use_abs, float *value)
float ** mne_cmatrix(int nr, int nc)
void mne_free_cmatrix(float **m)
void mne_free_name_list(char **list, int nlist)
#define ALLOC_CMATRIX(x, y)
#define SEG_LEN
#define FREE_CMATRIX(m)
#define MALLOC(x, t)
Core MNE data structures (source spaces, source estimates, hemispheres).
MNEChSelection * mneChSelection
Inverse source estimation (MNE, dSPM, sLORETA, dipole fitting).
Forward modelling (BEM, MEG/EEG lead fields).
Definition compute_fwd.h:95
Multi-layer spherical head model for EEG forward computation.
static FwdEegSphereModel * setup_eeg_sphere_model(const QString &eeg_model_file, QString eeg_model_name, float eeg_sphere_rad)
Measurement data container holding multiple data sets for MNE inverse computations.
void adjust_baselines(float bmin, float bmax)
static MNEMeasData * mne_read_meas_data(const QString &name, int set, MNEInverseOperator *op, MNELIB::MNENamedMatrix *fwd, const QStringList &namesp, int nnamesp)
MNEMeasDataSet * current
static int fit_dipoles_raw(const QString &dataname, MNELIB::MNERawData *raw, MNELIB::mneChSelection sel, DipoleFitData *fit, GuessData *guess, float tmin, float tmax, float tstep, float integ, int verbose, ECDSet &p_set)
DipoleFit(DipoleFitSettings *p_settings)
ECDSet calculateFit() const
static int fit_dipoles(const QString &dataname, MNEMeasData *data, DipoleFitData *fit, GuessData *guess, float tmin, float tmax, float tstep, float integ, int verbose, ECDSet &p_set)
Dipole Fit Data implementation.
static DipoleFitData * setup_dipole_fit_data(const QString &mriname, const QString &measname, const QString &bemname, Eigen::Vector3f *r0, FWDLIB::FwdEegSphereModel *eeg_model, int accurate_coils, const QString &badname, const QString &noisename, float grad_std, float mag_std, float eeg_std, float mag_reg, float grad_reg, float eeg_reg, int diagnoise, const QList< QString > &projnames, int include_meg, int include_eeg)
static bool fit_one(DipoleFitData *fit, GuessData *guess, float time, float *B, int verbose, ECD &res)
static int scale_noise_cov(DipoleFitData *f, int nave)
Dipole Fit setting implementation.
Single equivalent current dipole with position, orientation, amplitude, and goodness-of-fit.
Definition ecd.h:73
void print(FILE *f) const
Definition ecd.cpp:88
Holds a set of Electric Current Dipoles.
Definition ecd_set.h:81
void addEcd(const ECD &p_ecd)
Definition ecd_set.cpp:162
qint32 size() const
Definition ecd_set.h:193
Precomputed guess point grid with forward fields for initial dipole position candidates.
Definition guess_data.h:79
Eigen::VectorXi pick_deriv
A comprehensive raw data structure.
std::unique_ptr< MNELIB::MNERawInfo > info
std::unique_ptr< MNELIB::MNEDeriv > deriv_matched
static MNERawData * open_file(const QString &name, int omit_skip, int allow_maxshield, const MNEFilterDef &filter)
int pick_data_filt(mneChSelection sel, int firsts, int ns, float **picked)
Information about raw data in a FIFF file.
QList< FIFFLIB::FiffChInfo > chInfo