def demo_kalman_voltimeter(): """ Examples -------- >>> demo_kalman_voltimeter() >>> plt.close() """ V0 = 12 h = np.atleast_2d(1) # voltimeter measure the voltage itself q = 1e-9 # variance of process noise as the car operates r = 0.05**2 # variance of measurement error b = 0 # no system input u = 0 # no system input filt = Kalman(R=r, A=1, Q=q, H=h, B=b) # Generate random voltages and watch the filter operate. n = 50 truth = np.random.randn(n) * np.sqrt(q) + V0 z = truth + np.random.randn(n) * np.sqrt(r) # measurement x = np.zeros(n) for i, zi in enumerate(z): x[i] = filt(zi, u) # perform a Kalman filter iteration _hz = plt.plot(z, 'r.', label='observations') # a-posteriori state estimates: _hx = plt.plot(x, 'b-', label='Kalman output') _ht = plt.plot(truth, 'g-', label='true voltage') plt.legend() plt.title('Automobile Voltimeter Example')
def kde_demo4(N=50): """Demonstrate that the improved Sheather-Jones plug-in (hisj) is superior for 1D multimodal distributions KDEDEMO4 shows that the improved Sheather-Jones plug-in smoothing is a better compared to normal reference rules (in this case the hns) Examples -------- >>> kde_demo4() """ data = np.hstack((st.norm.rvs(loc=5, scale=1, size=(N,)), st.norm.rvs(loc=-5, scale=1, size=(N,)))) # x = np.linspace(1.5e-3, 5, 55) kde = KDE(data, kernel=Kernel('gauss', 'hns')) f = kde(output='plot', title='Ordinary KDE', plotflag=1) kde1 = KDE(data, kernel=Kernel('gauss', 'hisj')) f1 = kde1(output='plot', label='Ordinary KDE', plotflag=1) plt.figure(0) f.plot('r', label='hns={0}'.format(kde.hs)) # plt.figure(2) f1.plot('b', label='hisj={0}'.format(kde1.hs)) x = np.linspace(-9, 9) plt.plot(x, (st.norm.pdf(x, loc=-5, scale=1) + st.norm.pdf(x, loc=5, scale=1)) / 2, 'k:', label='True density') plt.legend()
def test_kalman(): V0 = 12 h = np.atleast_2d(1) # voltimeter measure the voltage itself q = 1e-9 # variance of process noise as the car operates r = 0.05 ** 2 # variance of measurement error b = 0 # no system input u = 0 # no system input filt = Kalman(R=r, A=1, Q=q, H=h, B=b) # Generate random voltages and watch the filter operate. n = 50 truth = np.random.randn(n) * np.sqrt(q) + V0 z = truth + np.random.randn(n) * np.sqrt(r) # measurement x = np.zeros(n) for i, zi in enumerate(z): x[i] = filt(zi, u) # perform a Kalman filter iteration _hz = plt.plot(z, 'r.', label='observations') # a-posteriori state estimates: _hx = plt.plot(x, 'b-', label='Kalman output') _ht = plt.plot(truth, 'g-', label='true voltage') plt.legend() plt.title('Automobile Voltimeter Example') plt.show('hold')
def kreg_demo1(hs=None, fast=True, fun='hisj'): """Compare KRegression to KernelReg from statsmodels.nonparametric Examples -------- >>> kreg_demo1() """ N = 100 # ei = np.random.normal(loc=0, scale=0.075, size=(N,)) ei = np.array([ -0.08508516, 0.10462496, 0.07694448, -0.03080661, 0.05777525, 0.06096313, -0.16572389, 0.01838912, -0.06251845, -0.09186784, -0.04304887, -0.13365788, -0.0185279, -0.07289167, 0.02319097, 0.06887854, -0.08938374, -0.15181813, 0.03307712, 0.08523183, -0.0378058, -0.06312874, 0.01485772, 0.06307944, -0.0632959, 0.18963205, 0.0369126, -0.01485447, 0.04037722, 0.0085057, -0.06912903, 0.02073998, 0.1174351, 0.17599277, -0.06842139, 0.12587608, 0.07698113, -0.0032394, -0.12045792, -0.03132877, 0.05047314, 0.02013453, 0.04080741, 0.00158392, 0.10237899, -0.09069682, 0.09242174, -0.15445323, 0.09190278, 0.07138498, 0.03002497, 0.02495252, 0.01286942, 0.06449978, 0.03031802, 0.11754861, -0.02322272, 0.00455867, -0.02132251, 0.09119446, -0.03210086, -0.06509545, 0.07306443, 0.04330647, 0.078111, -0.04146907, 0.05705476, 0.02492201, -0.03200572, -0.02859788, -0.05893749, 0.00089538, 0.0432551, 0.04001474, 0.04888828, -0.17708392, 0.16478644, 0.1171006, 0.11664846, 0.01410477, -0.12458953, -0.11692081, 0.0413047, -0.09292439, -0.07042327, 0.14119701, -0.05114335, 0.04994696, -0.09520663, 0.04829406, -0.01603065, -0.1933216, 0.19352763, 0.11819496, 0.04567619, -0.08348306, 0.00812816, -0.00908206, 0.14528945, 0.02901065]) x = np.linspace(0, 1, N) va_1 = 0.3 ** 2 va_2 = 0.7 ** 2 y0 = np.exp(-x ** 2 / (2 * va_1)) + 1.3 * np.exp(-(x - 1) ** 2 / (2 * va_2)) y = y0 + ei kernel = Kernel('gauss', fun=fun) hopt = kernel.hisj(x) kreg = KRegression( x, y, p=0, hs=hs, kernel=kernel, xmin=-2 * hopt, xmax=1 + 2 * hopt) if fast: kreg.__call__ = kreg.eval_grid_fast f = kreg(x, output='plot', title='Kernel regression', plotflag=1) plt.figure(0) f.plot(label='p=0') kreg.p = 1 f1 = kreg(x, output='plot', title='Kernel regression', plotflag=1) f1.plot(label='p=1') # print(f1.data) plt.plot(x, y, '.', label='data') plt.plot(x, y0, 'k', label='True model') from statsmodels.nonparametric.kernel_regression import KernelReg kreg2 = KernelReg(y, x, ('c')) y2 = kreg2.fit(x) plt.plot(x, y2[0], 'm', label='statsmodel') plt.legend()
def _plot_error(neval, err_dic, plot_error): if plot_error: plt.figure(0) for name in err_dic: plt.loglog(neval, err_dic[name], label=name) plt.xlabel('number of function evaluations') plt.ylabel('error') plt.legend()
def demo_kalman_sine(): """Kalman Filter demonstration with sine signal. Examples -------- >>> demo_kalman_sine() >>> plt.close() """ sd = 0.5 dt = 0.1 w = 1 T = np.arange(0, 30 + dt / 2, dt) n = len(T) X = 3 * np.sin(w * T) Y = X + sd * np.random.randn(n) ''' Initialize KF to values x = 0 dx/dt = 0 with great uncertainty in derivative ''' M = np.zeros((2, 1)) P = np.diag([0.1, 2]) R = sd**2 H = np.atleast_2d([1, 0]) q = 0.1 F = np.atleast_2d([[0, 1], [0, 0]]) A, Q = lti_disc(F, L=None, Q=np.diag([0, q]), dt=dt) # Track and animate m = M.shape[0] _MM = np.zeros((m, n)) _PP = np.zeros((m, m, n)) '''In this demonstration we estimate a stationary sine signal from noisy measurements by using the classical Kalman filter.' ''' filt = Kalman(R=R, x=M, P=P, A=A, Q=Q, H=H, B=0) # Generate random voltages and watch the filter operate. # n = 50 # truth = np.random.randn(n) * np.sqrt(q) + V0 # z = truth + np.random.randn(n) * np.sqrt(r) # measurement truth = X z = Y x = np.zeros((n, m)) for i, zi in enumerate(z): x[i] = np.ravel(filt(zi, u=0)) _hz = plt.plot(z, 'r.', label='observations') # a-posteriori state estimates: _hx = plt.plot(x[:, 0], 'b-', label='Kalman output') _ht = plt.plot(truth, 'g-', label='true voltage') plt.legend() plt.title('Automobile Voltimeter Example')
def test_kalman_sine(): """Kalman Filter demonstration with sine signal.""" sd = 1. dt = 0.1 w = 1 T = np.arange(0, 30 + dt / 2, dt) n = len(T) X = np.sin(w * T) Y = X + sd * np.random.randn(n) ''' Initialize KF to values x = 0 dx/dt = 0 with great uncertainty in derivative ''' M = np.zeros((2, 1)) P = np.diag([0.1, 2]) R = sd ** 2 H = np.atleast_2d([1, 0]) q = 0.1 F = np.atleast_2d([[0, 1], [0, 0]]) A, Q = lti_disc(F, L=None, Q=np.diag([0, q]), dt=dt) # Track and animate m = M.shape[0] _MM = np.zeros((m, n)) _PP = np.zeros((m, m, n)) '''In this demonstration we estimate a stationary sine signal from noisy measurements by using the classical Kalman filter.' ''' filt = Kalman(R=R, x=M, P=P, A=A, Q=Q, H=H, B=0) # Generate random voltages and watch the filter operate. # n = 50 # truth = np.random.randn(n) * np.sqrt(q) + V0 # z = truth + np.random.randn(n) * np.sqrt(r) # measurement truth = X z = Y x = np.zeros((n, m)) for i, zi in enumerate(z): x[i] = filt(zi, u=0).ravel() _hz = plt.plot(z, 'r.', label='observations') # a-posteriori state estimates: _hx = plt.plot(x[:, 0], 'b-', label='Kalman output') _ht = plt.plot(truth, 'g-', label='true voltage') plt.legend() plt.title('Automobile Voltimeter Example') plt.show()
def qdemo(f, a, b): ''' Compares different quadrature rules. Parameters ---------- f : callable function a,b : scalars lower and upper integration limits Details ------- qdemo(f,a,b) computes and compares various approximations to the integral of f from a to b. Three approximations are used, the composite trapezoid, Simpson's, and Boole's rules, all with equal length subintervals. In a case like qdemo(exp,0,3) one can see the expected convergence rates for each of the three methods. In a case like qdemo(sqrt,0,3), the convergence rate is limited not by the method, but by the singularity of the integrand. Example ------- >>> import numpy as np >>> qdemo(np.exp,0,3) true value = 19.08553692 ftn Trapezoid Simpsons Booles evals approx error approx error approx error 3, 22.5366862979, 3.4511493747, 19.5061466023, 0.4206096791, 19.4008539142, 0.3153169910 5, 19.9718950387, 0.8863581155, 19.1169646189, 0.0314276957, 19.0910191534, 0.0054822302 9, 19.3086731081, 0.2231361849, 19.0875991312, 0.0020622080, 19.0856414320, 0.0001045088 17, 19.1414188470, 0.0558819239, 19.0856674267, 0.0001305035, 19.0855386464, 0.0000017232 33, 19.0995135407, 0.0139766175, 19.0855451052, 0.0000081821, 19.0855369505, 0.0000000273 65, 19.0890314614, 0.0034945382, 19.0855374350, 0.0000005118, 19.0855369236, 0.0000000004 129, 19.0864105817, 0.0008736585, 19.0855369552, 0.0000000320, 19.0855369232, 0.0000000000 257, 19.0857553393, 0.0002184161, 19.0855369252, 0.0000000020, 19.0855369232, 0.0000000000 513, 19.0855915273, 0.0000546041, 19.0855369233, 0.0000000001, 19.0855369232, 0.0000000000 ftn Clenshaw Chebychev Gauss-L evals approx error approx error approx error 3, 19.5061466023, 0.4206096791, 0.0000000000, 1.0000000000, 19.0803304585, 0.0052064647 5, 19.0834145766, 0.0021223465, 0.0000000000, 1.0000000000, 19.0855365951, 0.0000003281 9, 19.0855369150, 0.0000000082, 0.0000000000, 1.0000000000, 19.0855369232, 0.0000000000 17, 19.0855369232, 0.0000000000, 0.0000000000, 1.0000000000, 19.0855369232, 0.0000000000 33, 19.0855369232, 0.0000000000, 0.0000000000, 1.0000000000, 19.0855369232, 0.0000000000 65, 19.0855369232, 0.0000000000, 0.0000000000, 1.0000000000, 19.0855369232, 0.0000000000 129, 19.0855369232, 0.0000000000, 0.0000000000, 1.0000000000, 19.0855369232, 0.0000000000 257, 19.0855369232, 0.0000000000, 0.0000000000, 1.0000000000, 19.0855369232, 0.0000000000 513, 19.0855369232, 0.0000000000, 0.0000000000, 1.0000000000, 19.0855369232, 0.0000000000 ''' # use quad8 with small tolerance to get "true" value #true1 = quad8(f,a,b,1e-10) #[true tol]= gaussq(f,a,b,1e-12) #[true tol] = agakron(f,a,b,1e-13) true_val, _tol = intg.quad(f, a, b) print('true value = %12.8f' % (true_val,)) kmax = 9 neval = zeros(kmax, dtype=int) qt = zeros(kmax) qs = zeros(kmax) qb = zeros(kmax) qc = zeros(kmax) qc2 = zeros(kmax) qg = zeros(kmax) et = ones(kmax) es = ones(kmax) eb = ones(kmax) ec = ones(kmax) ec2 = ones(kmax) ec3 = ones(kmax) eg = ones(kmax) # try various approximations for k in xrange(kmax): n = 2 ** (k + 1) + 1 neval[k] = n h = (b - a) / (n - 1) x = np.linspace(a, b, n) y = f(x) # trapezoid approximation q = np.trapz(y, x) #h*( (y(1)+y(n))/2 + sum(y(2:n-1)) ) qt[k] = q et[k] = abs(q - true_val) # Simpson approximation q = intg.simps(y, x) #(h/3)*( y(1)+y(n) + 4*sum(y(2:2:n-1)) + 2*sum(y(3:2:n-2)) ) qs[k] = q es[k] = abs(q - true_val) # Boole's rule #q = boole(x,y) q = (2 * h / 45) * (7 * (y[0] + y[-1]) + 12 * np.sum(y[2:n - 1:4]) + 32 * np.sum(y[1:n - 1:2]) + 14 * np.sum(y[4:n - 3:4])) qb[k] = q eb[k] = abs(q - true_val) # Clenshaw-Curtis [q, ec3[k]] = clencurt(f, a, b, (n - 1) / 2) qc[k] = q ec[k] = abs(q - true_val) # Chebychev #ck = chebfit(f,n,a,b) #q = chebval(b,chebint(ck,a,b),a,b) #qc2[k] = q; ec2[k] = abs(q - true) # Gauss-Legendre quadrature q = intg.fixed_quad(f, a, b, n=n)[0] #[x, w]=qrule(n,1) #x = (b-a)/2*x + (a+b)/2 % Transform base points X. #w = (b-a)/2*w % Adjust weigths. #q = sum(feval(f,x)*w) qg[k] = q eg[k] = abs(q - true_val) #% display results formats = ['%4.0f, ', ] + ['%10.10f, ', ]*6 formats[-1] = formats[-1].split(',')[0] data = np.vstack((neval, qt, et, qs, es, qb, eb)).T print(' ftn Trapezoid Simpson''s Boole''s') print('evals approx error approx error approx error') for k in xrange(kmax): tmp = data[k].tolist() print(''.join(fi % t for fi, t in zip(formats, tmp))) # display results data = np.vstack((neval, qc, ec, qc2, ec2, qg, eg)).T print(' ftn Clenshaw Chebychev Gauss-L') print('evals approx error approx error approx error') for k in xrange(kmax): tmp = data[k].tolist() print(''.join(fi % t for fi, t in zip(formats, tmp))) plt.loglog(neval, np.vstack((et, es, eb, ec, ec2, eg)).T) plt.xlabel('number of function evaluations') plt.ylabel('error') plt.legend(('Trapezoid', 'Simpsons', 'Booles', 'Clenshaw', 'Chebychev', 'Gauss-L'))
def demo_tide_filter(): """ Examples -------- >>> demo_tide_filter() >>> plt.close() """ # import statsmodels.api as sa import wafo.spectrum.models as sm sd = 10 Sj = sm.Jonswap(Hm0=4. * sd) S = Sj.tospecdata() q = (0.1 * sd)**2 # variance of process noise s the car operates r = (100 * sd)**2 # variance of measurement error b = 0 # no system input u = 0 # no system input from scipy.signal import butter, filtfilt, lfilter_zi # lfilter, freq_tide = 1. / (12 * 60 * 60) freq_wave = 1. / 10 freq_filt = freq_wave / 10 dt = 1. freq = 1. / dt fn = (freq / 2) P = 10 * np.diag([1, 0.01]) R = r H = np.atleast_2d([1, 0]) F = np.atleast_2d([[0, 1], [0, 0]]) A, Q = lti_disc(F, L=None, Q=np.diag([0, q]), dt=dt) t = np.arange(0, 60 * 12, 1. / freq) w = 2 * np.pi * freq # 1 Hz tide = 100 * np.sin(freq_tide * w * t + 2 * np.pi / 4) + 100 y = tide + S.sim(len(t), dt=1. / freq)[:, 1].ravel() # lowess = sa.nonparametric.lowess # y2 = lowess(y, t, frac=0.5)[:,1] filt = Kalman(R=R, x=np.array([[tide[0]], [0]]), P=P, A=A, Q=Q, H=H, B=b) filt2 = Kalman(R=R, x=np.array([[tide[0]], [0]]), P=P, A=A, Q=Q, H=H, B=b) # y = tide + 0.5 * np.sin(freq_wave * w * t) # Butterworth filter b, a = butter(9, (freq_filt / fn), btype='low') # y2 = [lowess(y[max(i-60,0):i + 1], t[max(i-60,0):i + 1], frac=.3)[-1,1] # for i in range(len(y))] # y2 = [lfilter(b, a, y[:i + 1])[i] for i in range(len(y))] # y3 = filtfilt(b, a, y[:16]).tolist() + [filtfilt(b, a, y[:i + 1])[i] # for i in range(16, len(y))] # y0 = medfilt(y, 41) _zi = lfilter_zi(b, a) # y2 = lfilter(b, a, y)#, zi=y[0]*zi) # standard filter y3 = filtfilt(b, a, y) # filter with phase shift correction y4 = [] y5 = [] for _i, j in enumerate(y): tmp = np.ravel(filt(j, u=u)) tmp = np.ravel(filt2(tmp[0], u=u)) # if i==0: # print(filt.x) # print(filt2.x) y4.append(tmp[0]) y5.append(tmp[1]) _y0 = medfilt(y4, 41) # print(filt.P) # plot plt.plot(t, y, 'r.-', linewidth=2, label='raw data') # plt.plot(t, y2, 'b.-', linewidth=2, label='lowess @ %g Hz' % freq_filt) # plt.plot(t, y2, 'b.-', linewidth=2, label='filter @ %g Hz' % freq_filt) plt.plot(t, y3, 'g.-', linewidth=2, label='filtfilt @ %g Hz' % freq_filt) plt.plot(t, y4, 'k.-', linewidth=2, label='kalman') # plt.plot(t, y5, 'k.', linewidth=2, label='kalman2') plt.plot(t, tide, 'y-', linewidth=2, label='True tide') plt.legend(frameon=False, fontsize=14) plt.xlabel("Time [s]") plt.ylabel("Amplitude")
def qdemo(f, a, b, kmax=9, plot_error=False): ''' Compares different quadrature rules. Parameters ---------- f : callable function a,b : scalars lower and upper integration limits Details ------- qdemo(f,a,b) computes and compares various approximations to the integral of f from a to b. Three approximations are used, the composite trapezoid, Simpson's, and Boole's rules, all with equal length subintervals. In a case like qdemo(exp,0,3) one can see the expected convergence rates for each of the three methods. In a case like qdemo(sqrt,0,3), the convergence rate is limited not by the method, but by the singularity of the integrand. Example ------- >>> import numpy as np >>> qdemo(np.exp,0,3) true value = 19.08553692 ftn, Boole, Chebychev evals approx error approx error 3, 19.4008539142, 0.3153169910, 19.5061466023, 0.4206096791 5, 19.0910191534, 0.0054822302, 19.0910191534, 0.0054822302 9, 19.0856414320, 0.0001045088, 19.0855374134, 0.0000004902 17, 19.0855386464, 0.0000017232, 19.0855369232, 0.0000000000 33, 19.0855369505, 0.0000000273, 19.0855369232, 0.0000000000 65, 19.0855369236, 0.0000000004, 19.0855369232, 0.0000000000 129, 19.0855369232, 0.0000000000, 19.0855369232, 0.0000000000 257, 19.0855369232, 0.0000000000, 19.0855369232, 0.0000000000 513, 19.0855369232, 0.0000000000, 19.0855369232, 0.0000000000 ftn, Clenshaw-Curtis, Gauss-Legendre evals approx error approx error 3, 19.5061466023, 0.4206096791, 19.0803304585, 0.0052064647 5, 19.0834145766, 0.0021223465, 19.0855365951, 0.0000003281 9, 19.0855369150, 0.0000000082, 19.0855369232, 0.0000000000 17, 19.0855369232, 0.0000000000, 19.0855369232, 0.0000000000 33, 19.0855369232, 0.0000000000, 19.0855369232, 0.0000000000 65, 19.0855369232, 0.0000000000, 19.0855369232, 0.0000000000 129, 19.0855369232, 0.0000000000, 19.0855369232, 0.0000000000 257, 19.0855369232, 0.0000000000, 19.0855369232, 0.0000000000 513, 19.0855369232, 0.0000000000, 19.0855369232, 0.0000000000 ftn, Simps, Trapz evals approx error approx error 3, 19.5061466023, 0.4206096791, 22.5366862979, 3.4511493747 5, 19.1169646189, 0.0314276957, 19.9718950387, 0.8863581155 9, 19.0875991312, 0.0020622080, 19.3086731081, 0.2231361849 17, 19.0856674267, 0.0001305035, 19.1414188470, 0.0558819239 33, 19.0855451052, 0.0000081821, 19.0995135407, 0.0139766175 65, 19.0855374350, 0.0000005118, 19.0890314614, 0.0034945382 129, 19.0855369552, 0.0000000320, 19.0864105817, 0.0008736585 257, 19.0855369252, 0.0000000020, 19.0857553393, 0.0002184161 513, 19.0855369233, 0.0000000001, 19.0855915273, 0.0000546041 ''' true_val, _tol = intg.quad(f, a, b) print('true value = %12.8f' % (true_val,)) neval = zeros(kmax, dtype=int) vals_dic = {} err_dic = {} # try various approximations methods = [trapz, simps, boole, ] for k in xrange(kmax): n = 2 ** (k + 1) + 1 neval[k] = n x = np.linspace(a, b, n) y = f(x) for method in methods: name = method.__name__.title() q = method(y, x) vals_dic.setdefault(name, []).append(q) err_dic.setdefault(name, []).append(abs(q - true_val)) name = 'Clenshaw-Curtis' q, _ec3 = clencurt(f, a, b, (n - 1) / 2) vals_dic.setdefault(name, []).append(q[0]) err_dic.setdefault(name, []).append(abs(q[0] - true_val)) name = 'Chebychev' ck = np.polynomial.chebyshev.chebfit(x, y, deg=min(n-1, 36)) cki = np.polynomial.chebyshev.chebint(ck) q = np.polynomial.chebyshev.chebval(x[-1], cki) vals_dic.setdefault(name, []).append(q) err_dic.setdefault(name, []).append(abs(q - true_val)) # ck = chebfit(f,n,a,b) # q = chebval(b,chebint(ck,a,b),a,b) # qc2[k] = q; ec2[k] = abs(q - true) name = 'Gauss-Legendre' # quadrature q = intg.fixed_quad(f, a, b, n=n)[0] # [x, w]=qrule(n,1) # x = (b-a)/2*x + (a+b)/2 % Transform base points X. # w = (b-a)/2*w % Adjust weigths. # q = sum(feval(f,x)*w) vals_dic.setdefault(name, []).append(q) err_dic.setdefault(name, []).append(abs(q - true_val)) # display results names = sorted(vals_dic.keys()) num_cols = 2 formats = ['%4.0f, ', ] + ['%10.10f, ', ] * num_cols * 2 formats[-1] = formats[-1].split(',')[0] formats_h = ['%4s, ', ] + ['%20s, ', ] * num_cols formats_h[-1] = formats_h[-1].split(',')[0] headers = ['evals'] + ['%12s %12s' % ('approx', 'error')] * num_cols while len(names) > 0: print(''.join(fi % t for fi, t in zip(formats_h, ['ftn'] + names[:num_cols]))) print(' '.join(headers)) data = [neval] for name in names[:num_cols]: data.append(vals_dic[name]) data.append(err_dic[name]) data = np.vstack(tuple(data)).T for k in xrange(kmax): tmp = data[k].tolist() print(''.join(fi % t for fi, t in zip(formats, tmp))) if plot_error: plt.figure(0) for name in names[:num_cols]: plt.loglog(neval, err_dic[name], label=name) names = names[num_cols:] if plot_error: plt.xlabel('number of function evaluations') plt.ylabel('error') plt.legend() plt.show('hold')
def test_tide_filter(): # import statsmodels.api as sa import wafo.spectrum.models as sm sd = 10 Sj = sm.Jonswap(Hm0=4.*sd) S = Sj.tospecdata() q = (0.1 * sd) ** 2 # variance of process noise s the car operates r = (100 * sd) ** 2 # variance of measurement error b = 0 # no system input u = 0 # no system input from scipy.signal import butter, filtfilt, lfilter_zi # lfilter, freq_tide = 1. / (12 * 60 * 60) freq_wave = 1. / 10 freq_filt = freq_wave / 10 dt = 1. freq = 1. / dt fn = (freq / 2) P = 10 * np.diag([1, 0.01]) R = r H = np.atleast_2d([1, 0]) F = np.atleast_2d([[0, 1], [0, 0]]) A, Q = lti_disc(F, L=None, Q=np.diag([0, q]), dt=dt) t = np.arange(0, 60 * 12, 1. / freq) w = 2 * np.pi * freq # 1 Hz tide = 100 * np.sin(freq_tide * w * t + 2 * np.pi / 4) + 100 y = tide + S.sim(len(t), dt=1. / freq)[:, 1].ravel() # lowess = sa.nonparametric.lowess # y2 = lowess(y, t, frac=0.5)[:,1] filt = Kalman(R=R, x=np.array([[tide[0]], [0]]), P=P, A=A, Q=Q, H=H, B=b) filt2 = Kalman(R=R, x=np.array([[tide[0]], [0]]), P=P, A=A, Q=Q, H=H, B=b) # y = tide + 0.5 * np.sin(freq_wave * w * t) # Butterworth filter b, a = butter(9, (freq_filt / fn), btype='low') # y2 = [lowess(y[max(i-60,0):i + 1], t[max(i-60,0):i + 1], frac=.3)[-1,1] # for i in range(len(y))] # y2 = [lfilter(b, a, y[:i + 1])[i] for i in range(len(y))] # y3 = filtfilt(b, a, y[:16]).tolist() + [filtfilt(b, a, y[:i + 1])[i] # for i in range(16, len(y))] # y0 = medfilt(y, 41) _zi = lfilter_zi(b, a) # y2 = lfilter(b, a, y)#, zi=y[0]*zi) # standard filter y3 = filtfilt(b, a, y) # filter with phase shift correction y4 = [] y5 = [] for _i, j in enumerate(y): tmp = filt(j, u=u).ravel() tmp = filt2(tmp[0], u=u).ravel() # if i==0: # print(filt.x) # print(filt2.x) y4.append(tmp[0]) y5.append(tmp[1]) _y0 = medfilt(y4, 41) print(filt.P) # plot plt.plot(t, y, 'r.-', linewidth=2, label='raw data') # plt.plot(t, y2, 'b.-', linewidth=2, label='lowess @ %g Hz' % freq_filt) # plt.plot(t, y2, 'b.-', linewidth=2, label='filter @ %g Hz' % freq_filt) plt.plot(t, y3, 'g.-', linewidth=2, label='filtfilt @ %g Hz' % freq_filt) plt.plot(t, y4, 'k.-', linewidth=2, label='kalman') # plt.plot(t, y5, 'k.', linewidth=2, label='kalman2') plt.plot(t, tide, 'y-', linewidth=2, label='True tide') plt.legend(frameon=False, fontsize=14) plt.xlabel("Time [s]") plt.ylabel("Amplitude") plt.show('hold')