forked from davidackerman/nnmf
-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathestimatenoise.m
executable file
·404 lines (393 loc) · 12.9 KB
/
estimatenoise.m
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
function noisevar = estimatenoise(X,varargin)
% estimatenoise: additive noise estimation from a time series
% usage: noisevar = estimatenoise(X)
% usage: noisevar = estimatenoise(X,dim)
% usage: noisevar = estimatenoise(X,t)
% usage: noisevar = estimatenoise(X,t,dim)
%
%
% Assume the time series has additive, Gaussian noise on top
% of a smoothly varying signal. The noise component variance
% is extracted. I'd see this code as useful for Kalman filter
% estimation, where one needs an estimate of the noise variance
% in data. It can also be applied to spline fitting problems,
% where one can control the regularization parameter to match
% a known noise variance.
%
%
% arguments: (input)
% X - (REQUIRED) (numeric vector of length n, or any n-dimensional
% numeric array)
%
% X must be at least of length 5 in the specified dimension.
%
% WARNING: Estimation of any variance will be suspect with
% few data points in your sample.
%
% If X is complex, then I treat the problem as two parallel
% sequences in parallel. The variance returned will be the
% sum of the variances for the real and complex parts.
%
% t - (OPTIONAL) - list of "time" steps if X was sampled at an
% non-uniform spacing. If supplied, then t must be a vector
% of the same length as X is in dimension dim.
%
% dim - specifies the dimension of the array to act on
% If dim is not provided, then the first non-singleton
% dimension will be used.
%
%
% arguments: (output)
% noisevar - scalar (or vector or array as appropriate)
% estimated variance of additive noise in the data
% along dimension dim.
%
%
% Example:
% Simple linear data, with purely additive N(0,1) gaussian noise
%
% t = 0:10000;
% x = t + randn(size(t));
% mv = estimatenoise(x)
% mv =
% 1.0166
%
% Example:
% Gaussian noise added to a sine wave (Nominal variance = 0.01)
% t = linspace(0,1,1000)';
% x = sin(t*50) + randn(size(t))/10;
%
% mv = estimatenoise(x)
%
% mv =
% 0.0096887
%
% Example:
% Pure gaussian noise, with a nominal variance of 9. (Note that var
% would have been a better estimator for this particular case...)
%
% mv = estimatenoise(3*randn(2,3,1000),3)
% mv =
% 9.6584 8.2696 8.632
% 9.2404 8.5346 9.7725
%
% Example:
% A piecewise constant function with multiple discontinuities.
% The true noise variance should be 0.01.
%
% t = linspace(0,1,1000);
% X = round(cos(t*6*pi)) + randn(size(t))/10;
% plot(t,X)
%
% var(X)
% ans =
% 0.68256
%
% estimatenoise(X)
% ans =
% 0.010882
%
% Example:
% Test if estimatenoise is able to recover the variance of a
% normally distributed random sample with unit variance. (Yes,
% it will be much slower than var.)
%
% mean(estimatenoise(randn(1000,1000)))
% ans =
% 1.0002
%
% Example:
% Use of extimatenoise on a problem with non-uniform spacing
% The actual noise variance was 1.0 here. Perform the
% operation 1000 times, then look at the median estimate of
% the variance. How well did we do?
%
% t = sort([randn(1,100) , randn(1,100)+5]);
% X = repmat(sin(t*5)*100,1000,1) + randn(1000,length(t));
%
% Estimatenoise is clearly wrong when the spacing is ignored
% median(estimatenoise(X,2))
% ans =
% 16.438
%
% Supplying the sampling "times", we get quite a reasonable result.
% median(estimatenoise(X,2,t))
% ans =
% 1.1307
%
% Example:
% Use of extimatenoise on a problem with replicates. In this
% example, each point will be replicated up to 3 times.
% The actual noise variance was again 1.0 here. Look at the
% increase in times required for estimatenoise when t is
% supplied.
%
% t = sort([0:1:100, 1:2:100, 1:4:100]);
% X = repmat(sin(t/10)*100,1000,1) + randn(1000,length(t));
%
% Estimatenoise is clearly wrong when the non-uniform spacing
% is ignored.
%
% tic,median(estimatenoise(X,2)),toc
% ans =
% 4.2056
% Elapsed time is 2.690135 seconds.
%
% Supplying the sampling "times", we get quite a reasonable
% result. The cost in time is not quite 2x.
%
% tic,median(estimatenoise(X,2,t)),toc
% ans =
% 1.0116
% Elapsed time is 4.864486 seconds.
%
%
% Methodology behind estimatenoise:
% Assume that the time series X is composed of a smoothly varying
% function, plus additive gaussian noise with mean zero and
% variance sigma^2. Then at any point it can be represented
% by a low order polynomial (think of it as a truncated local
% Taylor series approximation). The finite differences applied
% using filter was designed to kill off low order polynomial
% terms, while leaving behind a correlated series of points,
% but the variance of that series will be unchanged at sigma^2.
%
% Estimatenoise is designed to be robust to a few arbitrary
% spikes or discontinuities in the function or its derivatives.
% This was achieved by trimming off the tails of the distributions,
% then using percentiles to back out the desired variance.
%
% This robustness has a penalty however, as it takes a bit more
% time to run, and there is also an observed bias in the results.
% This bias is predictable as a function of the series length, so
% one of the last steps in estimatenoise backs out that predicted
% bias from the estimated variance. As a test of this,
%
% X = randn(100000,20);
% nv = estimatenoise(X,2);
% mean(nv)
% ans =
% 0.98355
%
% Yes, var would have been a better estimator in this particular
% case, but estimatenoise was designed to solve a more difficult
% problem.
%
% When a non-uniform time series spacing is supplied, estimatenoise
% works the same way, but filter cannot be used to do the
% computations efficiently, since the filter coefficients vary
% along that series. Here null is used to compute the time
% varying filter coefficients.
%
%
% See also: std, var
%
%
% Author: John D'Errico
% e-mail: [email protected]
% Release: 2.0
% Release date: 8/25/07
% if no input arguments provided, just dump out the help
if nargin<1
help estimatenoise
return
end
% complex sequences are split into real and imaginary parts.
if ~isreal(X)
noisevar = estimatenoise(real(X)) + estimatenoise(imag(X));
return
end
% get the size of X, so we can pick the value for
% dim if it was not supplied. We will also need to
% make sure the length in that dimension is large
% enough.
sx = size(X);
% were t and/or dim supplied? In which order?
if (nargin>3)
error('Too many inputs supplied. Maximum of 3 inputs are allowed.')
elseif nargin == 1
% use defaults
dim = [];
t = [];
elseif (nargin == 2) && isempty(varargin{1})
% use defaults
dim = [];
t = [];
elseif (nargin==2) && (length(varargin{1}) == 1)
% dim was supplied
dim = varargin{1};
t = [];
elseif (nargin==2)
% t must have been specified
t = varargin{1};
t = t(:);
dim = [];
elseif (nargin==3) && (length(varargin{1}) == 1)
% dim and t were supplied, in that order
dim = varargin{1};
t = varargin{2};
t = t(:);
elseif (nargin==3) && (length(varargin{2}) == 1)
% t and dim were supplied, in that order
dim = varargin{2};
t = varargin{1};
t = t(:);
else
% there is a problem with the arguments
error('A and t arguments are of inconsistent size')
end
% get dim if it was not supplied
if isempty(dim)
dim = find(sx~=1,1,'first');
if isempty(dim)
error('X did not have at least one dimension of length >1')
end
end
% check the length in the dim dimension
nX = sx(dim);
if nX<5
error('The length of X in the specified dimension was less than 5')
end
% permute the dimensions so that the chosen dimension
% of X is moved to the end of the line.
ndim = length(sx);
nx = 1:ndim;
nx(dim) = [];
nx = [nx,dim];
Xp = permute(X,nx);
% then just reshape it to be a 2-d array
Xp = reshape(Xp,[],sx(dim));
% was t actually equally spaced anyway? If it was, then we
% want to use the equal spaced code.
equispaced = true;
if ~isempty(t)
% first sort t, just in case
[t,tags] = sort(t);
% and shuffle Xp in the first dimension
Xp = Xp(:,tags);
% check for equal spacing in t
dt = diff(t);
% average spacing
avespace = (t(end) - t(1))/(nX-1);
if (avespace == 0)
error('Invalid t vector: t was identically zero.')
end
tol = 10*eps(avespace);
if tol < (max(dt) - min(dt))
% unequal spacing
equispaced = false;
end
end
% estimate the measurement variability in the (now) second dimension
% later we will undo this reshape.
% The idea here is to form a linear combination of successive elements
% of the series. If the underlying form is locally nearly linear, then
% a [1 -2 1] combination (for equally spaced data) will leave only
% the noise remaining. Next, if we assume the measurement noise was
% iid, N(0,s^2), then we can try to back out the noise variance.
fda{1} = [1 -1];
fda{2} = [1 -2 1];
fda{3} = [1 -3 3 -1];
fda{4} = [1 -4 6 -4 1];
fda{5} = [1 -5 10 -10 5 -1];
fda{6} = [1 -6 15 -20 15 -6 1];
nfda = length(fda);
for i = 1:nfda
% normalize to unit norm
fda{i} = fda{i}/norm(fda{i});
end
% compute an interquantile range, like the distance between the 25%
% and 75% points. This trims off the trash at each end, potentially
% corrupted if there are discontinuities in the curve. It also deals
% simply with a non-zero mean in this data. Actually do this for
% several different interquantile ranges, then take a median.
% NOTE: While I could have used other methods for the final variance
% estimation, this method was chosen to avoid outlier issues when
% the curve may have isolated discontinuities in function value or
% a derivative.
% The following points correspond to the central 90, 80, 75, 70, 65,
% 60, 55, 50, 45, 40, 35, 30, 25, and 20 percent ranges.
perc = [0.05, 0.1:0.025:0.40];
z = erfinv((1-perc)*2-1)*sqrt(2);
sigmaest = nan(size(Xp,1),nfda);
for i = 1:nfda
% Apply each difference to the data series with convn
% the 'valid' option will trim off the junk at the ends
% from the convolution.
% did we have equal spacing?
if equispaced
% If so, then conv will do the trick.
% These convolutions will yield noise of the desired variance,
% although it will be colored noise.
noisedata = conv2(Xp,fda{i},'valid');
ntrim = size(noisedata,2);
else
% unequal spacing. do it the hard way
F = 1./gamma(1:i);
ntrim = nX - i;
noisedata = zeros(size(Xp,1),ntrim);
for j = 1:ntrim
tj = t(j+(0:i));
tj = tj - mean(tj);
if all(abs(tj) <= tol)
% be careful if all points were reps!
coef = rem((0:i)',2)*2-1;
coef = coef - mean(coef);
% and normalize
coef = coef/norm(coef);
else
% use null to choose the linear combination
% of the elements of Xp.
tj = tj/norm(tj);
A = repmat(tj,1,i).^repmat(0:(i-1),i+1,1);
A = A.*repmat(F,i+1,1);
nullvecs = null(A');
% already normalized to have unit norm by null
coef = nullvecs(:,1);
end
% form the appropriate linear combination of Xp
noisedata(:,j) = Xp(:,j+(0:i))*coef;
end
end % if equispaced
% were there enough points to even try anything?
if ntrim >= 2
% sorting will provide the necessary percentiles after
% interpolation.
noisedata = sort(noisedata,2);
p = 0.5 + (1:ntrim)';
p = p/(ntrim + 0.5);
Q = zeros(size(noisedata,1),length(perc));
for j = 1:length(perc)
Qj = interp1q(p,noisedata',1-perc(j)) - interp1q(p,noisedata',perc(j));
Qj = Qj/(2*z(j));
Q(:,j) = Qj';
end
% Trim off any nans first, since if the series was short enough,
% some of those percentiles were not present.
wasnan = isnan(Q(1,:));
Q(:,wasnan) = [];
% Our noise std estimate is given by the median of the interquantile
% range(s). This is an ad hoc, but hopefully effective, way of
% estimating the measurement noise present in the signal.
sigmaest(:,i) = median(Q,2);
end
end % for i = 1:nfda
% drop those estimates which failed for lack of enough data
sigmaest(:,isnan(sigmaest(1,:))) = [];
% use median of these estimates to get a noise estimate.
noisevar = median(sigmaest,2).^2;
% Use an adhoc correction to remove the bias in the noise estimate.
% This correction was determined by examination of a large number of
% random samples.
noisevar = noisevar/(1+15*(sx(dim)+1.225)^-1.245);
% The first order difference might be used to guesstimate the
% process noise.
% if nargout>1
% processvar = max(0,sigmaest(:,1).^2 - noisevar);
% end
% finally, reshape the result to be consistent with the input shape
sx(dim) = 1;
noisevar = reshape(noisevar,sx);
% if nargout>1
% processvar = reshape(processvar,sx);
% end