Mercurial > repos > public > sbplib
comparison +parametrization/dataSpline.m @ 1098:36d092a00040 feature/dataspline
Bugfix Curve.arcLengthParametrization() by using linear interpolation instead of higher-order spline for parameter as a funciton of arclength. Generalize dataSpline() to curves in higher dimensions.
author | Martin Almquist <malmquist@stanford.edu> |
---|---|
date | Tue, 09 Apr 2019 09:50:44 -0700 |
parents | fdff002ebb32 |
children | 60c875c18de3 |
comparison
equal
deleted
inserted
replaced
1097:eec03e78c6f2 | 1098:36d092a00040 |
---|---|
1 % Accepts data points (t_i, f_i) and returns a Curve object, | 1 % dataSpline calculates a Curve through the points f_i using cubic spline interpolation. |
2 % using spline interpolation. | |
3 % The spline curve is parametrized with the arc length parametrization | 2 % The spline curve is parametrized with the arc length parametrization |
4 % to facilitate better grids. | 3 % to facilitate better grids. |
5 % | 4 % |
6 % t_data - vector | 5 % f - m x D matrix of m points in D dimensions |
7 % f_data - vector | 6 function C = dataSpline(f) |
8 % C - curve object | 7 m = size(f, 1); |
9 function C = dataSpline(t_data, f_data) | |
10 | 8 |
11 assert(length(t_data)==length(f_data),'Vectors must be same length'); | 9 t = linspace(0,1,m); |
12 m_data = length(t_data); | |
13 | 10 |
14 pp_g = spapi(4, t_data, f_data); | 11 pp_g = spapi(4, t, f'); |
15 pp_gp = fnder(pp_g); | 12 pp_gp = fnder(pp_g); |
16 | 13 |
17 % Reparametrize with parameter s from 0 to 1 to use Curve class | 14 g = @(t) fnval(pp_g, t); |
18 tmin = min(t_data); | 15 gp = @(t) fnval(pp_gp, t); |
19 tmax = max(t_data); | |
20 t = @(s) tmin + s*(tmax-tmin); | |
21 dt_ds = @(s) 0*s + (tmax-tmin); | |
22 g = @(s) [t(s); fnval(pp_g, t(s))]; | |
23 gp = @(s) [dt_ds(s); fnval(pp_gp, t(s)).*dt_ds(s)]; | |
24 | 16 |
25 % Create Curve object and reparametrize with arclength parametrization | |
26 C = parametrization.Curve(g, gp); | 17 C = parametrization.Curve(g, gp); |
27 C = C.arcLengthParametrization(m_data); | 18 C = C.arcLengthParametrization(); |
28 end | 19 end |