/
pca.py
287 lines (250 loc) · 12.1 KB
/
pca.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
"""
PCA and related signal decomposition methods for tedana
"""
import logging
import os.path as op
from numbers import Number
import numpy as np
import pandas as pd
from scipy import stats
from sklearn.decomposition import PCA
from tedana import metrics, utils, io
from tedana.decomposition import ma_pca
from tedana.stats import computefeats2
from tedana.selection import kundu_tedpca
LGR = logging.getLogger(__name__)
RepLGR = logging.getLogger('REPORT')
RefLGR = logging.getLogger('REFERENCES')
def low_mem_pca(data):
"""
Run Singular Value Decomposition (SVD) on input data.
Parameters
----------
data : (S [*E] x T) array_like
Optimally combined (S x T) or full multi-echo (S*E x T) data.
Returns
-------
u : (S [*E] x C) array_like
Component weight map for each component.
s : (C,) array_like
Variance explained for each component.
v : (C x T) array_like
Component timeseries.
"""
from sklearn.decomposition import IncrementalPCA
ppca = IncrementalPCA(n_components=(data.shape[-1] - 1))
ppca.fit(data)
v = ppca.components_.T
s = ppca.explained_variance_
u = np.dot(np.dot(data, v), np.diag(1. / s))
return u, s, v
def tedpca(data_cat, data_oc, combmode, mask, adaptive_mask, t2sG,
ref_img, tes, algorithm='mdl', kdaw=10., rdaw=1.,
out_dir='.', verbose=False, low_mem=False):
"""
Use principal components analysis (PCA) to identify and remove thermal
noise from multi-echo data.
Parameters
----------
data_cat : (S x E x T) array_like
Input functional data
data_oc : (S x T) array_like
Optimally combined time series data
combmode : {'t2s', 'paid'} str
How optimal combination of echos should be made, where 't2s' indicates
using the method of Posse 1999 and 'paid' indicates using the method of
Poser 2006
mask : (S,) array_like
Boolean mask array
adaptive_mask : (S,) array_like
Array where each value indicates the number of echoes with good signal
for that voxel. This mask may be thresholded; for example, with values
less than 3 set to 0.
For more information on thresholding, see `make_adaptive_mask`.
t2sG : (S,) array_like
Map of voxel-wise T2* estimates.
ref_img : :obj:`str` or img_like
Reference image to dictate how outputs are saved to disk
tes : :obj:`list`
List of echo times associated with `data_cat`, in milliseconds
algorithm : {'kundu', 'kundu-stabilize', 'mdl', 'aic', 'kic', float}, optional
Method with which to select components in TEDPCA. PCA
decomposition with the mdl, kic and aic options are based on a Moving Average
(stationary Gaussian) process and are ordered from most to least aggressive
(see Li et al., 2007).
If a float is provided, then it is assumed to represent percentage of variance
explained (0-1) to retain from PCA.
Default is 'mdl'.
kdaw : :obj:`float`, optional
Dimensionality augmentation weight for Kappa calculations. Must be a
non-negative float, or -1 (a special value). Default is 10.
rdaw : :obj:`float`, optional
Dimensionality augmentation weight for Rho calculations. Must be a
non-negative float, or -1 (a special value). Default is 1.
out_dir : :obj:`str`, optional
Output directory.
verbose : :obj:`bool`, optional
Whether to output files from fitmodels_direct or not. Default: False
low_mem : :obj:`bool`, optional
Whether to use incremental PCA (for low-memory systems) or not.
This is only compatible with the "kundu" or "kundu-stabilize" algorithms.
Default: False
Returns
-------
kept_data : (S x T) :obj:`numpy.ndarray`
Dimensionally reduced optimally combined functional data
n_components : :obj:`int`
Number of components retained from PCA decomposition
Notes
-----
====================== =================================================
Notation Meaning
====================== =================================================
:math:`\\kappa` Component pseudo-F statistic for TE-dependent
(BOLD) model.
:math:`\\rho` Component pseudo-F statistic for TE-independent
(artifact) model.
:math:`v` Voxel
:math:`V` Total number of voxels in mask
:math:`\\zeta` Something
:math:`c` Component
:math:`p` Something else
====================== =================================================
Steps:
1. Variance normalize either multi-echo or optimally combined data,
depending on settings.
2. Decompose normalized data using PCA or SVD.
3. Compute :math:`{\\kappa}` and :math:`{\\rho}`:
.. math::
{\\kappa}_c = \\frac{\\sum_{v}^V {\\zeta}_{c,v}^p * \
F_{c,v,R_2^*}}{\\sum {\\zeta}_{c,v}^p}
{\\rho}_c = \\frac{\\sum_{v}^V {\\zeta}_{c,v}^p * \
F_{c,v,S_0}}{\\sum {\\zeta}_{c,v}^p}
4. Some other stuff. Something about elbows.
5. Classify components as thermal noise if they meet both of the
following criteria:
- Nonsignificant :math:`{\\kappa}` and :math:`{\\rho}`.
- Nonsignificant variance explained.
Outputs:
This function writes out several files:
====================== =================================================
Filename Content
====================== =================================================
pca_decomposition.json PCA component table.
pca_mixing.tsv PCA mixing matrix.
pca_components.nii.gz Component weight maps.
====================== =================================================
See Also
--------
:func:`tedana.utils.make_adaptive_mask` : The function used to create the ``adaptive_mask``
parameter.
"""
if algorithm == 'kundu':
alg_str = ("followed by the Kundu component selection decision "
"tree (Kundu et al., 2013)")
RefLGR.info("Kundu, P., Brenowitz, N. D., Voon, V., Worbe, Y., "
"Vértes, P. E., Inati, S. J., ... & Bullmore, E. T. "
"(2013). Integrated strategy for improving functional "
"connectivity mapping using multiecho fMRI. Proceedings "
"of the National Academy of Sciences, 110(40), "
"16187-16192.")
elif algorithm == 'kundu-stabilize':
alg_str = ("followed by the 'stabilized' Kundu component "
"selection decision tree (Kundu et al., 2013)")
RefLGR.info("Kundu, P., Brenowitz, N. D., Voon, V., Worbe, Y., "
"Vértes, P. E., Inati, S. J., ... & Bullmore, E. T. "
"(2013). Integrated strategy for improving functional "
"connectivity mapping using multiecho fMRI. Proceedings "
"of the National Academy of Sciences, 110(40), "
"16187-16192.")
elif isinstance(algorithm, Number):
alg_str = (
"in which the number of components was determined based on a "
"variance explained threshold")
else:
alg_str = ("based on the PCA component estimation with a Moving Average"
"(stationary Gaussian) process (Li et al., 2007)")
RefLGR.info("Li, Y.O., Adalı, T. and Calhoun, V.D., (2007). "
"Estimating the number of independent components for "
"functional magnetic resonance imaging data. "
"Human brain mapping, 28(11), pp.1251-1266.")
RepLGR.info("Principal component analysis {0} was applied to "
"the optimally combined data for dimensionality "
"reduction.".format(alg_str))
n_samp, n_echos, n_vols = data_cat.shape
LGR.info('Computing PCA of optimally combined multi-echo data')
data = data_oc[mask, :]
data_z = ((data.T - data.T.mean(axis=0)) / data.T.std(axis=0)).T # var normalize ts
data_z = (data_z - data_z.mean()) / data_z.std() # var normalize everything
if algorithm in ['mdl', 'aic', 'kic']:
data_img = io.new_nii_like(ref_img, utils.unmask(data, mask))
mask_img = io.new_nii_like(ref_img, mask.astype(int))
voxel_comp_weights, varex, varex_norm, comp_ts = ma_pca.ma_pca(
data_img, mask_img, algorithm)
elif isinstance(algorithm, Number):
ppca = PCA(copy=False, n_components=algorithm, svd_solver="full")
ppca.fit(data_z)
comp_ts = ppca.components_.T
varex = ppca.explained_variance_
voxel_comp_weights = np.dot(np.dot(data_z, comp_ts),
np.diag(1. / varex))
varex_norm = varex / varex.sum()
elif low_mem:
voxel_comp_weights, varex, comp_ts = low_mem_pca(data_z)
varex_norm = varex / varex.sum()
else:
ppca = PCA(copy=False, n_components=(n_vols - 1))
ppca.fit(data_z)
comp_ts = ppca.components_.T
varex = ppca.explained_variance_
voxel_comp_weights = np.dot(np.dot(data_z, comp_ts),
np.diag(1. / varex))
varex_norm = varex / varex.sum()
# Compute Kappa and Rho for PCA comps
# Normalize each component's time series
vTmixN = stats.zscore(comp_ts, axis=0)
comptable, _, _, _ = metrics.dependence_metrics(
data_cat, data_oc, comp_ts, adaptive_mask, tes, ref_img,
reindex=False, mmixN=vTmixN, algorithm=None,
label='mepca_', out_dir=out_dir, verbose=verbose)
# varex_norm from PCA overrides varex_norm from dependence_metrics,
# but we retain the original
comptable['estimated normalized variance explained'] = \
comptable['normalized variance explained']
comptable['normalized variance explained'] = varex_norm
# write component maps to 4D image
comp_ts_z = stats.zscore(comp_ts, axis=0)
comp_maps = utils.unmask(computefeats2(data_oc, comp_ts_z, mask), mask)
io.filewrite(comp_maps, op.join(out_dir, 'pca_components.nii.gz'), ref_img)
# Select components using decision tree
if algorithm == 'kundu':
comptable = kundu_tedpca(comptable, n_echos, kdaw, rdaw, stabilize=False)
elif algorithm == 'kundu-stabilize':
comptable = kundu_tedpca(comptable, n_echos, kdaw, rdaw, stabilize=True)
else:
alg_str = "variance explained-based" if isinstance(algorithm, Number) else algorithm
LGR.info('Selected {0} components with {1} dimensionality '
'detection'.format(comptable.shape[0], alg_str))
comptable['classification'] = 'accepted'
comptable['rationale'] = ''
# Save decomposition
comp_names = [io.add_decomp_prefix(comp, prefix='pca', max_value=comptable.index.max())
for comp in comptable.index.values]
mixing_df = pd.DataFrame(data=comp_ts, columns=comp_names)
mixing_df.to_csv(op.join(out_dir, 'pca_mixing.tsv'), sep='\t', index=False)
comptable['Description'] = 'PCA fit to optimally combined data.'
mmix_dict = {}
mmix_dict['Method'] = ('Principal components analysis implemented by '
'sklearn. Components are sorted by variance '
'explained in descending order. '
'Component signs are flipped to best match the '
'data.')
io.save_comptable(comptable, op.join(out_dir, 'pca_decomposition.json'),
label='pca', metadata=mmix_dict)
acc = comptable[comptable.classification == 'accepted'].index.values
n_components = acc.size
voxel_kept_comp_weighted = (voxel_comp_weights[:, acc] * varex[None, acc])
kept_data = np.dot(voxel_kept_comp_weighted, comp_ts[:, acc].T)
kept_data = stats.zscore(kept_data, axis=1) # variance normalize time series
kept_data = stats.zscore(kept_data, axis=None) # variance normalize everything
return kept_data, n_components