-
Notifications
You must be signed in to change notification settings - Fork 1.3k
/
maxwell.py
1797 lines (1609 loc) · 79.2 KB
/
maxwell.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
# -*- coding: utf-8 -*-
# Authors: Mark Wronkiewicz <wronk.mark@gmail.com>
# Eric Larson <larson.eric.d@gmail.com>
# Jussi Nurminen <jnu@iki.fi>
# License: BSD (3-clause)
from functools import partial
from math import factorial
from os import path as op
import numpy as np
from scipy import linalg
from .. import __version__
from ..bem import _check_origin
from ..chpi import quat_to_rot, rot_to_quat
from ..transforms import (_str_to_frame, _get_trans, Transform, apply_trans,
_find_vector_rotation, _cart_to_sph, _get_n_moments,
_sph_to_cart_partials, _deg_ord_idx,
_sh_complex_to_real, _sh_real_to_complex, _sh_negate)
from ..forward import _concatenate_coils, _prep_meg_channels, _create_meg_coils
from ..surface import _normalize_vectors
from ..io.constants import FIFF
from ..io.proc_history import _read_ctc
from ..io.write import _generate_meas_id, _date_now
from ..io import _loc_to_coil_trans, BaseRaw
from ..io.pick import pick_types, pick_info
from ..utils import verbose, logger, _clean_names, warn, _time_mask, _pl
from ..fixes import _get_args, _safe_svd, _get_sph_harm
from ..externals.six import string_types
from ..channels.channels import _get_T1T2_mag_inds
# Note: Elekta uses single precision and some algorithms might use
# truncated versions of constants (e.g., μ0), which could lead to small
# differences between algorithms
@verbose
def maxwell_filter(raw, origin='auto', int_order=8, ext_order=3,
calibration=None, cross_talk=None, st_duration=None,
st_correlation=0.98, coord_frame='head', destination=None,
regularize='in', ignore_ref=False, bad_condition='error',
head_pos=None, st_fixed=True, st_only=False, mag_scale=100.,
verbose=None):
u"""Apply Maxwell filter to data using multipole moments.
.. warning:: Automatic bad channel detection is not currently implemented.
It is critical to mark bad channels before running Maxwell
filtering to prevent artifact spreading.
.. warning:: Maxwell filtering in MNE is not designed or certified
for clinical use.
Parameters
----------
raw : instance of mne.io.Raw
Data to be filtered
origin : array-like, shape (3,) | str
Origin of internal and external multipolar moment space in meters.
The default is ``'auto'``, which means a head-digitization-based
origin fit when ``coord_frame='head'``, and ``(0., 0., 0.)`` when
``coord_frame='meg'``.
int_order : int
Order of internal component of spherical expansion.
ext_order : int
Order of external component of spherical expansion.
calibration : str | None
Path to the ``'.dat'`` file with fine calibration coefficients.
File can have 1D or 3D gradiometer imbalance correction.
This file is machine/site-specific.
cross_talk : str | None
Path to the FIF file with cross-talk correction information.
st_duration : float | None
If not None, apply spatiotemporal SSS with specified buffer duration
(in seconds). Elekta's default is 10.0 seconds in MaxFilter™ v2.2.
Spatiotemporal SSS acts as implicitly as a high-pass filter where the
cut-off frequency is 1/st_dur Hz. For this (and other) reasons, longer
buffers are generally better as long as your system can handle the
higher memory usage. To ensure that each window is processed
identically, choose a buffer length that divides evenly into your data.
Any data at the trailing edge that doesn't fit evenly into a whole
buffer window will be lumped into the previous buffer.
st_correlation : float
Correlation limit between inner and outer subspaces used to reject
ovwrlapping intersecting inner/outer signals during spatiotemporal SSS.
coord_frame : str
The coordinate frame that the ``origin`` is specified in, either
``'meg'`` or ``'head'``. For empty-room recordings that do not have
a head<->meg transform ``info['dev_head_t']``, the MEG coordinate
frame should be used.
destination : str | array-like, shape (3,) | None
The destination location for the head. Can be ``None``, which
will not change the head position, or a string path to a FIF file
containing a MEG device<->head transformation, or a 3-element array
giving the coordinates to translate to (with no rotations).
For example, ``destination=(0, 0, 0.04)`` would translate the bases
as ``--trans default`` would in MaxFilter™ (i.e., to the default
head location).
regularize : str | None
Basis regularization type, must be "in" or None.
"in" is the same algorithm as the "-regularize in" option in
MaxFilter™.
ignore_ref : bool
If True, do not include reference channels in compensation. This
option should be True for KIT files, since Maxwell filtering
with reference channels is not currently supported.
bad_condition : str
How to deal with ill-conditioned SSS matrices. Can be "error"
(default), "warning", or "ignore".
head_pos : array | None
If array, movement compensation will be performed.
The array should be of shape (N, 10), holding the position
parameters as returned by e.g. `read_head_pos`.
.. versionadded:: 0.12
st_fixed : bool
If True (default), do tSSS using the median head position during the
``st_duration`` window. This is the default behavior of MaxFilter
and has been most extensively tested.
.. versionadded:: 0.12
st_only : bool
If True, only tSSS (temporal) projection of MEG data will be
performed on the output data. The non-tSSS parameters (e.g.,
``int_order``, ``calibration``, ``head_pos``, etc.) will still be
used to form the SSS bases used to calculate temporal projectors,
but the ouptut MEG data will *only* have temporal projections
performed. Noise reduction from SSS basis multiplication,
cross-talk cancellation, movement compensation, and so forth
will not be applied to the data. This is useful, for example, when
evoked movement compensation will be performed with
:func:`mne.epochs.average_movements`.
.. versionadded:: 0.12
mag_scale : float | str
The magenetometer scale-factor used to bring the magnetometers
to approximately the same order of magnitude as the gradiometers
(default 100.), as they have different units (T vs T/m).
Can be ``'auto'`` to use the reciprocal of the physical distance
between the gradiometer pickup loops (e.g., 0.0168 m yields
59.5 for VectorView).
.. versionadded:: 0.13
verbose : bool, str, int, or None
If not None, override default verbose level (see :func:`mne.verbose`
and :ref:`Logging documentation <tut_logging>` for more).
Returns
-------
raw_sss : instance of mne.io.Raw
The raw data with Maxwell filtering applied.
See Also
--------
mne.chpi.filter_chpi
mne.chpi.read_head_pos
mne.epochs.average_movements
Notes
-----
.. versionadded:: 0.11
Some of this code was adapted and relicensed (with BSD form) with
permission from Jussi Nurminen. These algorithms are based on work
from [1]_ and [2]_.
.. note:: This code may use multiple CPU cores, see the
:ref:`FAQ <faq_cpu>` for more information.
Compared to Elekta's MaxFilter™ software, the MNE Maxwell filtering
routines currently provide the following features:
+-----------------------------------------------------------------------------+-----+-----------+
| Feature | MNE | MaxFilter |
+=============================================================================+=====+===========+
| Maxwell filtering software shielding | X | X |
+-----------------------------------------------------------------------------+-----+-----------+
| Bad channel reconstruction | X | X |
+-----------------------------------------------------------------------------+-----+-----------+
| Cross-talk cancellation | X | X |
+-----------------------------------------------------------------------------+-----+-----------+
| Fine calibration correction (1D) | X | X |
+-----------------------------------------------------------------------------+-----+-----------+
| Fine calibration correction (3D) | X | |
+-----------------------------------------------------------------------------+-----+-----------+
| Spatio-temporal SSS (tSSS) | X | X |
+-----------------------------------------------------------------------------+-----+-----------+
| Coordinate frame translation | X | X |
+-----------------------------------------------------------------------------+-----+-----------+
| Regularization using information theory | X | X |
+-----------------------------------------------------------------------------+-----+-----------+
| Movement compensation (raw) | X | X |
+-----------------------------------------------------------------------------+-----+-----------+
| Movement compensation (:func:`epochs <mne.epochs.average_movements>`) | X | |
+-----------------------------------------------------------------------------+-----+-----------+
| :func:`cHPI subtraction <mne.chpi.filter_chpi>` | X | X |
+-----------------------------------------------------------------------------+-----+-----------+
| Double floating point precision | X | |
+-----------------------------------------------------------------------------+-----+-----------+
| Seamless processing of split (``-1.fif``) and concatenated files | X | |
+-----------------------------------------------------------------------------+-----+-----------+
| Certified for clinical use | | X |
+-----------------------------------------------------------------------------+-----+-----------+
| Automatic bad channel detection | | X |
+-----------------------------------------------------------------------------+-----+-----------+
| Head position estimation | | X |
+-----------------------------------------------------------------------------+-----+-----------+
Epoch-based movement compensation is described in [1]_.
Use of Maxwell filtering routines with non-Elekta systems is currently
**experimental**. Worse results for non-Elekta systems are expected due
to (at least):
* Missing fine-calibration and cross-talk cancellation data for
other systems.
* Processing with reference sensors has not been vetted.
* Regularization of components may not work well for all systems.
* Coil integration has not been optimized using Abramowitz/Stegun
definitions.
.. note:: Various Maxwell filtering algorithm components are covered by
patents owned by Elekta Oy, Helsinki, Finland.
These patents include, but may not be limited to:
- US2006031038 (Signal Space Separation)
- US6876196 (Head position determination)
- WO2005067789 (DC fields)
- WO2005078467 (MaxShield)
- WO2006114473 (Temporal Signal Space Separation)
These patents likely preclude the use of Maxwell filtering code
in commercial applications. Consult a lawyer if necessary.
Currently, in order to perform Maxwell filtering, the raw data must not
have any projectors applied. During Maxwell filtering, the spatial
structure of the data is modified, so projectors are discarded (unless
in ``st_only=True`` mode).
References
----------
.. [1] Taulu S. and Kajola M. "Presentation of electromagnetic
multichannel data: The signal space separation method,"
Journal of Applied Physics, vol. 97, pp. 124905 1-10, 2005.
http://lib.tkk.fi/Diss/2008/isbn9789512295654/article2.pdf
.. [2] Taulu S. and Simola J. "Spatiotemporal signal space separation
method for rejecting nearby interference in MEG measurements,"
Physics in Medicine and Biology, vol. 51, pp. 1759-1768, 2006.
http://lib.tkk.fi/Diss/2008/isbn9789512295654/article3.pdf
""" # noqa: E501
# There are an absurd number of different possible notations for spherical
# coordinates, which confounds the notation for spherical harmonics. Here,
# we purposefully stay away from shorthand notation in both and use
# explicit terms (like 'azimuth' and 'polar') to avoid confusion.
# See mathworld.wolfram.com/SphericalHarmonic.html for more discussion.
# Our code follows the same standard that ``scipy`` uses for ``sph_harm``.
# triage inputs ASAP to avoid late-thrown errors
if not isinstance(raw, BaseRaw):
raise TypeError('raw must be Raw, not %s' % type(raw))
_check_usable(raw)
_check_regularize(regularize)
st_correlation = float(st_correlation)
if st_correlation <= 0. or st_correlation > 1.:
raise ValueError('Need 0 < st_correlation <= 1., got %s'
% st_correlation)
if coord_frame not in ('head', 'meg'):
raise ValueError('coord_frame must be either "head" or "meg", not "%s"'
% coord_frame)
head_frame = True if coord_frame == 'head' else False
recon_trans = _check_destination(destination, raw.info, head_frame)
if st_duration is not None:
st_duration = float(st_duration)
if not 0. < st_duration <= raw.times[-1] + 1. / raw.info['sfreq']:
raise ValueError('st_duration (%0.1fs) must be between 0 and the '
'duration of the data (%0.1fs).'
% (st_duration, raw.times[-1]))
st_correlation = float(st_correlation)
st_duration = int(round(st_duration * raw.info['sfreq']))
if not 0. < st_correlation <= 1:
raise ValueError('st_correlation must be between 0. and 1.')
if not isinstance(bad_condition, string_types) or \
bad_condition not in ['error', 'warning', 'ignore']:
raise ValueError('bad_condition must be "error", "warning", or '
'"ignore", not %s' % bad_condition)
if raw.info['dev_head_t'] is None and coord_frame == 'head':
raise RuntimeError('coord_frame cannot be "head" because '
'info["dev_head_t"] is None; if this is an '
'empty room recording, consider using '
'coord_frame="meg"')
if st_only and st_duration is None:
raise ValueError('st_duration must not be None if st_only is True')
head_pos = _check_pos(head_pos, head_frame, raw, st_fixed,
raw.info['sfreq'])
_check_info(raw.info, sss=not st_only, tsss=st_duration is not None,
calibration=not st_only and calibration is not None,
ctc=not st_only and cross_talk is not None)
# Now we can actually get moving
logger.info('Maxwell filtering raw data')
add_channels = (head_pos[0] is not None) and not st_only
raw_sss, pos_picks = _copy_preload_add_channels(
raw, add_channels=add_channels)
del raw
if not st_only:
# remove MEG projectors, they won't apply now
_remove_meg_projs(raw_sss)
info = raw_sss.info
meg_picks, mag_picks, grad_picks, good_picks, mag_or_fine = \
_get_mf_picks(info, int_order, ext_order, ignore_ref)
# Magnetometers are scaled to improve numerical stability
coil_scale, mag_scale = _get_coil_scale(
meg_picks, mag_picks, grad_picks, mag_scale, info)
#
# Fine calibration processing (load fine cal and overwrite sensor geometry)
#
sss_cal = dict()
if calibration is not None:
calibration, sss_cal = _update_sensor_geometry(info, calibration,
ignore_ref)
mag_or_fine.fill(True) # all channels now have some mag-type data
# Determine/check the origin of the expansion
origin = _check_origin(origin, info, coord_frame, disp=True)
origin.setflags(write=False)
n_in, n_out = _get_n_moments([int_order, ext_order])
#
# Cross-talk processing
#
if cross_talk is not None:
sss_ctc = _read_ctc(cross_talk)
ctc_chs = sss_ctc['proj_items_chs']
meg_ch_names = [info['ch_names'][p] for p in meg_picks]
# checking for extra space ambiguity in channel names
# between old and new fif files
if meg_ch_names[0] not in ctc_chs:
ctc_chs = _clean_names(ctc_chs, remove_whitespace=True)
missing = sorted(list(set(meg_ch_names) - set(ctc_chs)))
if len(missing) != 0:
raise RuntimeError('Missing MEG channels in cross-talk matrix:\n%s'
% missing)
missing = sorted(list(set(ctc_chs) - set(meg_ch_names)))
if len(missing) > 0:
warn('Not all cross-talk channels in raw:\n%s' % missing)
ctc_picks = [ctc_chs.index(info['ch_names'][c])
for c in meg_picks[good_picks]]
assert len(ctc_picks) == len(good_picks) # otherwise we errored
ctc = sss_ctc['decoupler'][ctc_picks][:, ctc_picks]
# I have no idea why, but MF transposes this for storage..
sss_ctc['decoupler'] = sss_ctc['decoupler'].T.tocsc()
else:
sss_ctc = dict()
#
# Translate to destination frame (always use non-fine-cal bases)
#
exp = dict(origin=origin, int_order=int_order, ext_order=0)
all_coils = _prep_mf_coils(info, ignore_ref)
S_recon = _trans_sss_basis(exp, all_coils, recon_trans, coil_scale)
exp['ext_order'] = ext_order
# Reconstruct data from internal space only (Eq. 38), and rescale S_recon
S_recon /= coil_scale
if recon_trans is not None:
# warn if we have translated too far
diff = 1000 * (info['dev_head_t']['trans'][:3, 3] -
recon_trans['trans'][:3, 3])
dist = np.sqrt(np.sum(_sq(diff)))
if dist > 25.:
warn('Head position change is over 25 mm (%s) = %0.1f mm'
% (', '.join('%0.1f' % x for x in diff), dist))
# Reconstruct raw file object with spatiotemporal processed data
max_st = dict()
if st_duration is not None:
max_st.update(job=10, subspcorr=st_correlation,
buflen=st_duration / info['sfreq'])
logger.info(' Processing data using tSSS with st_duration=%s'
% max_st['buflen'])
st_when = 'before' if st_fixed else 'after' # relative to movecomp
else:
# st_duration from here on will act like the chunk size
st_duration = max(int(round(10. * info['sfreq'])), 1)
st_correlation = None
st_when = 'never'
st_duration = min(len(raw_sss.times), st_duration)
del st_fixed
# Generate time points to break up data into equal-length windows
read_lims = np.arange(0, len(raw_sss.times) + 1, st_duration)
if len(read_lims) == 1:
read_lims = np.concatenate([read_lims, [len(raw_sss.times)]])
if read_lims[-1] != len(raw_sss.times):
read_lims[-1] = len(raw_sss.times)
# len_last_buf < st_dur so fold it into the previous buffer
if st_correlation is not None and len(read_lims) > 2:
logger.info(' Spatiotemporal window did not fit evenly into '
'raw object. The final %0.2f seconds were lumped '
'onto the previous window.'
% ((read_lims[-1] - read_lims[-2]) / info['sfreq'],))
assert len(read_lims) >= 2
assert read_lims[0] == 0 and read_lims[-1] == len(raw_sss.times)
#
# Do the heavy lifting
#
# Figure out which transforms we need for each tSSS block
# (and transform pos[1] to times)
head_pos[1] = raw_sss.time_as_index(head_pos[1], use_rounding=True)
# Compute the first bit of pos_data for cHPI reporting
if info['dev_head_t'] is not None and head_pos[0] is not None:
this_pos_quat = np.concatenate([
rot_to_quat(info['dev_head_t']['trans'][:3, :3]),
info['dev_head_t']['trans'][:3, 3],
np.zeros(3)])
else:
this_pos_quat = None
_get_this_decomp_trans = partial(
_get_decomp, all_coils=all_coils,
cal=calibration, regularize=regularize,
exp=exp, ignore_ref=ignore_ref, coil_scale=coil_scale,
grad_picks=grad_picks, mag_picks=mag_picks, good_picks=good_picks,
mag_or_fine=mag_or_fine, bad_condition=bad_condition,
mag_scale=mag_scale)
S_decomp, pS_decomp, reg_moments, n_use_in = _get_this_decomp_trans(
info['dev_head_t'], t=0.)
reg_moments_0 = reg_moments.copy()
# Loop through buffer windows of data
n_sig = int(np.floor(np.log10(max(len(read_lims), 0)))) + 1
logger.info(' Processing %s data chunk%s of (at least) %0.1f sec'
% (len(read_lims) - 1, _pl(read_lims),
st_duration / info['sfreq']))
for ii, (start, stop) in enumerate(zip(read_lims[:-1], read_lims[1:])):
rel_times = raw_sss.times[start:stop]
t_str = '%8.3f - %8.3f sec' % tuple(rel_times[[0, -1]])
t_str += ('(#%d/%d)'
% (ii + 1, len(read_lims) - 1)).rjust(2 * n_sig + 5)
# Get original data
orig_data = raw_sss._data[meg_picks[good_picks], start:stop]
# This could just be np.empty if not st_only, but shouldn't be slow
# this way so might as well just always take the original data
out_meg_data = raw_sss._data[meg_picks, start:stop]
# Apply cross-talk correction
if cross_talk is not None:
orig_data = ctc.dot(orig_data)
out_pos_data = np.empty((len(pos_picks), stop - start))
# Figure out which positions to use
t_s_s_q_a = _trans_starts_stops_quats(head_pos, start, stop,
this_pos_quat)
n_positions = len(t_s_s_q_a[0])
# Set up post-tSSS or do pre-tSSS
if st_correlation is not None:
# If doing tSSS before movecomp...
resid = orig_data.copy() # to be safe let's operate on a copy
if st_when == 'after':
orig_in_data = np.empty((len(meg_picks), stop - start))
else: # 'before'
avg_trans = t_s_s_q_a[-1]
if avg_trans is not None:
# if doing movecomp
S_decomp_st, pS_decomp_st, _, n_use_in_st = \
_get_this_decomp_trans(avg_trans, t=rel_times[0])
else:
S_decomp_st, pS_decomp_st = S_decomp, pS_decomp
n_use_in_st = n_use_in
orig_in_data = np.dot(np.dot(S_decomp_st[:, :n_use_in_st],
pS_decomp_st[:n_use_in_st]),
resid)
resid -= np.dot(np.dot(S_decomp_st[:, n_use_in_st:],
pS_decomp_st[n_use_in_st:]), resid)
resid -= orig_in_data
# Here we operate on our actual data
proc = out_meg_data if st_only else orig_data
_do_tSSS(proc, orig_in_data, resid, st_correlation,
n_positions, t_str)
if not st_only or st_when == 'after':
# Do movement compensation on the data
for trans, rel_start, rel_stop, this_pos_quat in \
zip(*t_s_s_q_a[:4]):
# Recalculate bases if necessary (trans will be None iff the
# first position in this interval is the same as last of the
# previous interval)
if trans is not None:
S_decomp, pS_decomp, reg_moments, n_use_in = \
_get_this_decomp_trans(trans, t=rel_times[rel_start])
# Determine multipole moments for this interval
mm_in = np.dot(pS_decomp[:n_use_in],
orig_data[:, rel_start:rel_stop])
# Our output data
if not st_only:
out_meg_data[:, rel_start:rel_stop] = \
np.dot(S_recon.take(reg_moments[:n_use_in], axis=1),
mm_in)
if len(pos_picks) > 0:
out_pos_data[:, rel_start:rel_stop] = \
this_pos_quat[:, np.newaxis]
# Transform orig_data to store just the residual
if st_when == 'after':
# Reconstruct data using original location from external
# and internal spaces and compute residual
rel_resid_data = resid[:, rel_start:rel_stop]
orig_in_data[:, rel_start:rel_stop] = \
np.dot(S_decomp[:, :n_use_in], mm_in)
rel_resid_data -= np.dot(np.dot(S_decomp[:, n_use_in:],
pS_decomp[n_use_in:]),
rel_resid_data)
rel_resid_data -= orig_in_data[:, rel_start:rel_stop]
# If doing tSSS at the end
if st_when == 'after':
_do_tSSS(out_meg_data, orig_in_data, resid, st_correlation,
n_positions, t_str)
elif st_when == 'never' and head_pos[0] is not None:
logger.info(' Used % 2d head position%s for %s'
% (n_positions, _pl(n_positions), t_str))
raw_sss._data[meg_picks, start:stop] = out_meg_data
raw_sss._data[pos_picks, start:stop] = out_pos_data
# Update info
if not st_only:
info['dev_head_t'] = recon_trans # set the reconstruction transform
_update_sss_info(raw_sss, origin, int_order, ext_order, len(good_picks),
coord_frame, sss_ctc, sss_cal, max_st, reg_moments_0,
st_only)
logger.info('[done]')
return raw_sss
def _get_coil_scale(meg_picks, mag_picks, grad_picks, mag_scale, info):
"""Get the magnetometer scale factor."""
if isinstance(mag_scale, string_types):
if mag_scale != 'auto':
raise ValueError('mag_scale must be a float or "auto", got "%s"'
% mag_scale)
if len(mag_picks) in (0, len(meg_picks)):
mag_scale = 100. # only one coil type, doesn't matter
logger.info(' Setting mag_scale=%0.2f because only one '
'coil type is present' % mag_scale)
else:
# Find our physical distance between gradiometer pickup loops
# ("base line")
coils = _create_meg_coils(pick_info(info, meg_picks)['chs'],
'accurate')
grad_base = set(coils[pick]['base'] for pick in grad_picks)
if len(grad_base) != 1 or list(grad_base)[0] <= 0:
raise RuntimeError('Could not automatically determine '
'mag_scale, could not find one '
'proper gradiometer distance from: %s'
% list(grad_base))
grad_base = list(grad_base)[0]
mag_scale = 1. / grad_base
logger.info(' Setting mag_scale=%0.2f based on gradiometer '
'distance %0.2f mm' % (mag_scale, 1000 * grad_base))
mag_scale = float(mag_scale)
coil_scale = np.ones((len(meg_picks), 1))
coil_scale[mag_picks] = mag_scale
return coil_scale, mag_scale
def _remove_meg_projs(inst):
"""Remove inplace existing MEG projectors (assumes inactive)."""
meg_picks = pick_types(inst.info, meg=True, exclude=[])
meg_channels = [inst.ch_names[pi] for pi in meg_picks]
non_meg_proj = list()
for proj in inst.info['projs']:
if not any(c in meg_channels for c in proj['data']['col_names']):
non_meg_proj.append(proj)
inst.add_proj(non_meg_proj, remove_existing=True, verbose=False)
def _check_destination(destination, info, head_frame):
"""Triage our reconstruction trans."""
if destination is None:
return info['dev_head_t']
if not head_frame:
raise RuntimeError('destination can only be set if using the '
'head coordinate frame')
if isinstance(destination, string_types):
recon_trans = _get_trans(destination, 'meg', 'head')[0]
elif isinstance(destination, Transform):
recon_trans = destination
else:
destination = np.array(destination, float)
if destination.shape != (3,):
raise ValueError('destination must be a 3-element vector, '
'str, or None')
recon_trans = np.eye(4)
recon_trans[:3, 3] = destination
recon_trans = Transform('meg', 'head', recon_trans)
if recon_trans.to_str != 'head' or recon_trans.from_str != 'MEG device':
raise RuntimeError('Destination transform is not MEG device -> head, '
'got %s -> %s' % (recon_trans.from_str,
recon_trans.to_str))
return recon_trans
def _prep_mf_coils(info, ignore_ref=True):
"""Get all coil integration information loaded and sorted."""
coils, comp_coils = _prep_meg_channels(
info, accurate=True, elekta_defs=True, head_frame=False,
ignore_ref=ignore_ref, verbose=False)[:2]
mag_mask = _get_mag_mask(coils)
if len(comp_coils) > 0:
meg_picks = pick_types(info, meg=True, ref_meg=False, exclude=[])
ref_picks = pick_types(info, meg=False, ref_meg=True, exclude=[])
inserts = np.searchsorted(meg_picks, ref_picks)
# len(inserts) == len(comp_coils)
for idx, comp_coil in zip(inserts[::-1], comp_coils[::-1]):
coils.insert(idx, comp_coil)
# Now we have:
# [c['chname'] for c in coils] ==
# [info['ch_names'][ii]
# for ii in pick_types(info, meg=True, ref_meg=True)]
# Now coils is a sorted list of coils. Time to do some vectorization.
n_coils = len(coils)
rmags = np.concatenate([coil['rmag'] for coil in coils])
cosmags = np.concatenate([coil['cosmag'] for coil in coils])
ws = np.concatenate([coil['w'] for coil in coils])
cosmags *= ws[:, np.newaxis]
del ws
n_int = np.array([len(coil['rmag']) for coil in coils])
bins = np.repeat(np.arange(len(n_int)), n_int)
bd = np.concatenate(([0], np.cumsum(n_int)))
slice_map = dict((ii, slice(start, stop))
for ii, (start, stop) in enumerate(zip(bd[:-1], bd[1:])))
return rmags, cosmags, bins, n_coils, mag_mask, slice_map
def _trans_starts_stops_quats(pos, start, stop, this_pos_data):
"""Get all trans and limits we need."""
pos_idx = np.arange(*np.searchsorted(pos[1], [start, stop]))
used = np.zeros(stop - start, bool)
trans = list()
rel_starts = list()
rel_stops = list()
quats = list()
if this_pos_data is None:
avg_trans = None
else:
avg_trans = np.zeros(6)
for ti in range(-1, len(pos_idx)):
# first iteration for this block of data
if ti < 0:
rel_start = 0
rel_stop = pos[1][pos_idx[0]] if len(pos_idx) > 0 else stop
rel_stop = rel_stop - start
if rel_start == rel_stop:
continue # our first pos occurs on first time sample
# Don't calculate S_decomp here, use the last one
trans.append(None) # meaning: use previous
quats.append(this_pos_data)
else:
rel_start = pos[1][pos_idx[ti]] - start
if ti == len(pos_idx) - 1:
rel_stop = stop - start
else:
rel_stop = pos[1][pos_idx[ti + 1]] - start
trans.append(pos[0][pos_idx[ti]])
quats.append(pos[2][pos_idx[ti]])
assert 0 <= rel_start
assert rel_start < rel_stop
assert rel_stop <= stop - start
assert not used[rel_start:rel_stop].any()
used[rel_start:rel_stop] = True
rel_starts.append(rel_start)
rel_stops.append(rel_stop)
if this_pos_data is not None:
avg_trans += quats[-1][:6] * (rel_stop - rel_start)
assert used.all()
# Use weighted average for average trans over the window
if avg_trans is not None:
avg_trans /= (stop - start)
avg_trans = np.vstack([
np.hstack([quat_to_rot(avg_trans[:3]),
avg_trans[3:][:, np.newaxis]]),
[[0., 0., 0., 1.]]])
return trans, rel_starts, rel_stops, quats, avg_trans
def _do_tSSS(clean_data, orig_in_data, resid, st_correlation,
n_positions, t_str):
"""Compute and apply SSP-like projection vectors based on min corr."""
np.asarray_chkfinite(resid)
t_proj = _overlap_projector(orig_in_data, resid, st_correlation)
# Apply projector according to Eq. 12 in [2]_
msg = (' Projecting %2d intersecting tSSS components '
'for %s' % (t_proj.shape[1], t_str))
if n_positions > 1:
msg += ' (across %2d positions)' % n_positions
logger.info(msg)
clean_data -= np.dot(np.dot(clean_data, t_proj), t_proj.T)
def _copy_preload_add_channels(raw, add_channels):
"""Load data for processing and (maybe) add cHPI pos channels."""
raw = raw.copy()
if add_channels:
kinds = [FIFF.FIFFV_QUAT_1, FIFF.FIFFV_QUAT_2, FIFF.FIFFV_QUAT_3,
FIFF.FIFFV_QUAT_4, FIFF.FIFFV_QUAT_5, FIFF.FIFFV_QUAT_6,
FIFF.FIFFV_HPI_G, FIFF.FIFFV_HPI_ERR, FIFF.FIFFV_HPI_MOV]
out_shape = (len(raw.ch_names) + len(kinds), len(raw.times))
out_data = np.zeros(out_shape, np.float64)
msg = ' Appending head position result channels and '
if raw.preload:
logger.info(msg + 'copying original raw data')
out_data[:len(raw.ch_names)] = raw._data
raw._data = out_data
else:
logger.info(msg + 'loading raw data from disk')
raw._preload_data(out_data[:len(raw.ch_names)], verbose=False)
raw._data = out_data
assert raw.preload is True
off = len(raw.ch_names)
chpi_chs = [
dict(ch_name='CHPI%03d' % (ii + 1), logno=ii + 1,
scanno=off + ii + 1, unit_mul=-1, range=1., unit=-1,
kind=kinds[ii], coord_frame=FIFF.FIFFV_COORD_UNKNOWN,
cal=1e-4, coil_type=FIFF.FWD_COIL_UNKNOWN, loc=np.zeros(12))
for ii in range(len(kinds))]
raw.info['chs'].extend(chpi_chs)
raw.info._update_redundant()
raw.info._check_consistency()
assert raw._data.shape == (raw.info['nchan'], len(raw.times))
# Return the pos picks
pos_picks = np.arange(len(raw.ch_names) - len(chpi_chs),
len(raw.ch_names))
return raw, pos_picks
else:
if not raw.preload:
logger.info(' Loading raw data from disk')
raw.load_data(verbose=False)
else:
logger.info(' Using loaded raw data')
return raw, np.array([], int)
def _check_pos(pos, head_frame, raw, st_fixed, sfreq):
"""Check for a valid pos array and transform it to a more usable form."""
if pos is None:
return [None, np.array([-1])]
if not head_frame:
raise ValueError('positions can only be used if coord_frame="head"')
if not st_fixed:
warn('st_fixed=False is untested, use with caution!')
if not isinstance(pos, np.ndarray):
raise TypeError('pos must be an ndarray')
if pos.ndim != 2 or pos.shape[1] != 10:
raise ValueError('pos must be an array of shape (N, 10)')
t = pos[:, 0]
t_off = raw.first_samp / raw.info['sfreq']
if not np.array_equal(t, np.unique(t)):
raise ValueError('Time points must unique and in ascending order')
# We need an extra 1e-3 (1 ms) here because MaxFilter outputs values
# only out to 3 decimal places
if not _time_mask(t, tmin=t_off - 1e-3, tmax=None, sfreq=sfreq).all():
raise ValueError('Head position time points must be greater than '
'first sample offset, but found %0.4f < %0.4f'
% (t[0], t_off))
max_dist = np.sqrt(np.sum(pos[:, 4:7] ** 2, axis=1)).max()
if max_dist > 1.:
warn('Found a distance greater than 1 m (%0.3g m) from the device '
'origin, positions may be invalid and Maxwell filtering could '
'fail' % (max_dist,))
dev_head_ts = np.zeros((len(t), 4, 4))
dev_head_ts[:, 3, 3] = 1.
dev_head_ts[:, :3, 3] = pos[:, 4:7]
dev_head_ts[:, :3, :3] = quat_to_rot(pos[:, 1:4])
pos = [dev_head_ts, t - t_off, pos[:, 1:]]
return pos
def _get_decomp(trans, all_coils, cal, regularize, exp, ignore_ref,
coil_scale, grad_picks, mag_picks, good_picks, mag_or_fine,
bad_condition, t, mag_scale):
"""Get a decomposition matrix and pseudoinverse matrices."""
#
# Fine calibration processing (point-like magnetometers and calib. coeffs)
#
S_decomp = _get_s_decomp(exp, all_coils, trans, coil_scale, cal,
ignore_ref, grad_picks, mag_picks, good_picks,
mag_scale)
#
# Regularization
#
S_decomp, pS_decomp, sing, reg_moments, n_use_in = _regularize(
regularize, exp, S_decomp, mag_or_fine, t=t)
# Pseudo-inverse of total multipolar moment basis set (Part of Eq. 37)
cond = sing[0] / sing[-1]
logger.debug(' Decomposition matrix condition: %0.1f' % cond)
if bad_condition != 'ignore' and cond >= 1000.:
msg = 'Matrix is badly conditioned: %0.0f >= 1000' % cond
if bad_condition == 'error':
raise RuntimeError(msg)
else: # condition == 'warning':
warn(msg)
# Build in our data scaling here
pS_decomp *= coil_scale[good_picks].T
S_decomp /= coil_scale[good_picks]
return S_decomp, pS_decomp, reg_moments, n_use_in
def _get_s_decomp(exp, all_coils, trans, coil_scale, cal, ignore_ref,
grad_picks, mag_picks, good_picks, mag_scale):
"""Get S_decomp."""
S_decomp = _trans_sss_basis(exp, all_coils, trans, coil_scale)
if cal is not None:
# Compute point-like mags to incorporate gradiometer imbalance
grad_cals = _sss_basis_point(exp, trans, cal, ignore_ref, mag_scale)
# Add point like magnetometer data to bases.
S_decomp[grad_picks, :] += grad_cals
# Scale magnetometers by calibration coefficient
S_decomp[mag_picks, :] /= cal['mag_cals']
# We need to be careful about KIT gradiometers
S_decomp = S_decomp[good_picks]
return S_decomp
@verbose
def _regularize(regularize, exp, S_decomp, mag_or_fine, t, verbose=None):
"""Regularize a decomposition matrix."""
# ALWAYS regularize the out components according to norm, since
# gradiometer-only setups (e.g., KIT) can have zero first-order
# components
int_order, ext_order = exp['int_order'], exp['ext_order']
n_in, n_out = _get_n_moments([int_order, ext_order])
t_str = '%8.3f' % t
if regularize is not None: # regularize='in'
logger.info(' Computing regularization')
in_removes, out_removes = _regularize_in(
int_order, ext_order, S_decomp, mag_or_fine)
else:
in_removes = []
out_removes = _regularize_out(int_order, ext_order, mag_or_fine)
reg_in_moments = np.setdiff1d(np.arange(n_in), in_removes)
reg_out_moments = np.setdiff1d(np.arange(n_in, n_in + n_out),
out_removes)
n_use_in = len(reg_in_moments)
n_use_out = len(reg_out_moments)
reg_moments = np.concatenate((reg_in_moments, reg_out_moments))
S_decomp = S_decomp.take(reg_moments, axis=1)
pS_decomp, sing = _col_norm_pinv(S_decomp.copy())
if regularize is not None or n_use_out != n_out:
logger.info(' Using %s/%s harmonic components for %s '
'(%s/%s in, %s/%s out)'
% (n_use_in + n_use_out, n_in + n_out, t_str,
n_use_in, n_in, n_use_out, n_out))
return S_decomp, pS_decomp, sing, reg_moments, n_use_in
def _get_mf_picks(info, int_order, ext_order, ignore_ref=False):
"""Pick types for Maxwell filtering."""
# Check for T1/T2 mag types
mag_inds_T1T2 = _get_T1T2_mag_inds(info)
if len(mag_inds_T1T2) > 0:
warn('%d T1/T2 magnetometer channel types found. If using SSS, it is '
'advised to replace coil types using "fix_mag_coil_types".'
% len(mag_inds_T1T2))
# Get indices of channels to use in multipolar moment calculation
ref = not ignore_ref
meg_picks = pick_types(info, meg=True, ref_meg=ref, exclude=[])
meg_info = pick_info(info, meg_picks)
del info
good_picks = pick_types(meg_info, meg=True, ref_meg=ref, exclude='bads')
n_bases = _get_n_moments([int_order, ext_order]).sum()
if n_bases > len(good_picks):
raise ValueError('Number of requested bases (%s) exceeds number of '
'good sensors (%s)' % (str(n_bases), len(good_picks)))
recons = [ch for ch in meg_info['bads']]
if len(recons) > 0:
logger.info(' Bad MEG channels being reconstructed: %s' % recons)
else:
logger.info(' No bad MEG channels')
ref_meg = False if ignore_ref else 'mag'
mag_picks = pick_types(meg_info, meg='mag', ref_meg=ref_meg, exclude=[])
ref_meg = False if ignore_ref else 'grad'
grad_picks = pick_types(meg_info, meg='grad', ref_meg=ref_meg, exclude=[])
assert len(mag_picks) + len(grad_picks) == len(meg_info['ch_names'])
# Determine which are magnetometers for external basis purposes
mag_or_fine = np.zeros(len(meg_picks), bool)
mag_or_fine[mag_picks] = True
# KIT gradiometers are marked as having units T, not T/M (argh)
# We need a separate variable for this because KIT grads should be
# treated mostly like magnetometers (e.g., scaled by 100) for reg
mag_or_fine[np.array([ch['coil_type'] & 0xFFFF == FIFF.FIFFV_COIL_KIT_GRAD
for ch in meg_info['chs']], bool)] = False
msg = (' Processing %s gradiometers and %s magnetometers'
% (len(grad_picks), len(mag_picks)))
n_kit = len(mag_picks) - mag_or_fine.sum()
if n_kit > 0:
msg += ' (of which %s are actually KIT gradiometers)' % n_kit
logger.info(msg)
return meg_picks, mag_picks, grad_picks, good_picks, mag_or_fine
def _check_regularize(regularize):
"""Ensure regularize is valid."""
if not (regularize is None or (isinstance(regularize, string_types) and
regularize in ('in',))):
raise ValueError('regularize must be None or "in"')
def _check_usable(inst):
"""Ensure our data are clean."""
if inst.proj:
raise RuntimeError('Projectors cannot be applied to data during '
'Maxwell filtering.')
current_comp = inst.compensation_grade
if current_comp not in (0, None):
raise RuntimeError('Maxwell filter cannot be done on compensated '
'channels, but data have been compensated with '
'grade %s.' % current_comp)
def _col_norm_pinv(x):
"""Compute the pinv with column-normalization to stabilize calculation.
Note: will modify/overwrite x.
"""
norm = np.sqrt(np.sum(x * x, axis=0))
x /= norm
u, s, v = linalg.svd(x, full_matrices=False, overwrite_a=True,
**check_disable)
v /= norm
return np.dot(v.T * 1. / s, u.T), s
def _sq(x):
"""Square quickly."""
return x * x
def _check_finite(data):
"""Ensure data is finite."""
if not np.isfinite(data).all():
raise RuntimeError('data contains non-finite numbers')
def _sph_harm_norm(order, degree):
"""Compute normalization factor for spherical harmonics."""
# we could use scipy.special.poch(degree + order + 1, -2 * order)
# here, but it's slower for our fairly small degree
norm = np.sqrt((2 * degree + 1.) / (4 * np.pi))
if order != 0:
norm *= np.sqrt(factorial(degree - order) /
float(factorial(degree + order)))
return norm
def _concatenate_sph_coils(coils):
"""Concatenate MEG coil parameters for spherical harmoncs."""
rs = np.concatenate([coil['r0_exey'] for coil in coils])
wcoils = np.concatenate([coil['w'] for coil in coils])
ezs = np.concatenate([np.tile(coil['ez'][np.newaxis, :],
(len(coil['rmag']), 1))
for coil in coils])
bins = np.repeat(np.arange(len(coils)),
[len(coil['rmag']) for coil in coils])
return rs, wcoils, ezs, bins
_mu_0 = 4e-7 * np.pi # magnetic permeability
def _get_mag_mask(coils):
"""Get the coil_scale for Maxwell filtering."""
return np.array([coil['coil_class'] == FIFF.FWD_COILC_MAG
for coil in coils])
def _sss_basis_basic(exp, coils, mag_scale=100., method='standard'):
"""Compute SSS basis using non-optimized (but more readable) algorithms."""
int_order, ext_order = exp['int_order'], exp['ext_order']
origin = exp['origin']
# Compute vector between origin and coil, convert to spherical coords
if method == 'standard':