bezier.curve module¶
Helper for Bézier Curves.
See Curve-Curve Intersection for examples using the
Curve
class to find intersections.
-
bezier._intersection_helpers.
linearization_error
(curve)¶ Compute the maximum error of a linear approximation.
We use the line
\[L(s) = v_0 (1 - s) + v_n s\]and compute a bound on the maximum error
\[\max_{s \in \left[0, 1\right]} \|B(s) - L(s)\|_2.\]Rather than computing the actual maximum (a tight bound), we use an upper bound via the remainder from Lagrange interpolation in each component. This leaves us with \(\frac{s(s - 1)}{2!}\) times the second derivative in each component.
The second derivative curve is degree \(d = n - 2\) and is given by
\[B''(s) = n(n - 1) \sum_{j = 0}^{d} \binom{d}{j} s^j (1 - s)^{d - j} \cdot \Delta^2 v_j\]Due to this form (and the convex combination property of Bézier Curves) we know each component of the second derivative will be bounded by the maximum of that component among the \(\Delta^2 v_j\).
For example, the curve
\[\begin{split}B(s) = \left[\begin{array}{c} 0 \\ 0 \end{array}\right] (1 - s)^2 + \left[\begin{array}{c} 1 \\ 3 \end{array}\right] 2s(1 - s) + \left[\begin{array}{c} -2 \\ 9 \end{array}\right] s^2\end{split}\]has \(B''(s) \equiv \left[\begin{array}{c} -8 \\ 6 \end{array}\right]\) which has norm \(10\) everywhere, hence the maximum error is
\[\left.\frac{s(1 - s)}{2!} \cdot 10\right|_{s = \frac{1}{2}} = \frac{5}{4}.\]>>> nodes = np.array([ ... [ 0.0, 0.0], ... [ 1.0, 3.0], ... [-2.0, 9.0], ... ]) >>> curve = bezier.Curve(nodes) >>> linearization_error(curve) 1.25
Parameters: curve (Curve) – A curve to be approximated by a line. Returns: The maximum error between the curve and the linear approximation. Return type: float
-
bezier._intersection_helpers.
newton_refine
(s, curve1, t, curve2)¶ Apply one step of 2D Newton’s method.
We want to use Newton’s method on the function
\[F(s, t) = B_1(s) - B_2(t)\]to refine \(\left(s_{\ast}, t_{\ast}\right)\). Using this, and the Jacobian \(DF\), we “solve”
\[\begin{split}\left[\begin{array}{c} 0 \\ 0 \end{array}\right] \approx F\left(s_{\ast} + \Delta s, t_{\ast} + \Delta t\right) \approx F\left(s_{\ast}, t_{\ast}\right) + \left[\begin{array}{c c} B_1'\left(s_{\ast}\right) & - B_2'\left(t_{\ast}\right) \end{array}\right] \left[\begin{array}{c} \Delta s \\ \Delta t \end{array}\right]\end{split}\]and refine with the component updates \(\Delta s\) and \(\Delta t\).
Note
This implementation assumes
curve1
andcurve2
live in \(\mathbf{R}^2\).For example, the curves
\[\begin{split}\begin{align*} B_1(s) &= \left[\begin{array}{c} 0 \\ 0 \end{array}\right] (1 - s)^2 + \left[\begin{array}{c} 2 \\ 4 \end{array}\right] 2s(1 - s) + \left[\begin{array}{c} 4 \\ 0 \end{array}\right] s^2 \\ B_2(t) &= \left[\begin{array}{c} 2 \\ 0 \end{array}\right] (1 - t) + \left[\begin{array}{c} 0 \\ 3 \end{array}\right] t \end{align*}\end{split}\]intersect at the point \(B_1\left(\frac{1}{4}\right) = B_2\left(\frac{1}{2}\right) = \frac{1}{2} \left[\begin{array}{c} 2 \\ 3 \end{array}\right]\). However, starting from the wrong point we have
\[\begin{split}\begin{align*} F\left(\frac{3}{8}, \frac{1}{4}\right) &= \frac{1}{8} \left[\begin{array}{c} 0 \\ 9 \end{array}\right] \\ DF\left(\frac{3}{8}, \frac{1}{4}\right) &= \left[\begin{array}{c c} 4 & 2 \\ 2 & -3 \end{array}\right] \\ \Longrightarrow \left[\begin{array}{c} \Delta s \\ \Delta t \end{array}\right] &= \frac{9}{64} \left[\begin{array}{c} -1 \\ 2 \end{array}\right]. \end{align*}\end{split}\]>>> nodes1 = np.array([ ... [0.0, 0.0], ... [2.0, 4.0], ... [4.0, 0.0], ... ]) >>> curve1 = bezier.Curve(nodes1) >>> nodes2 = np.array([ ... [2.0, 0.0], ... [0.0, 3.0], ... ]) >>> curve2 = bezier.Curve(nodes2) >>> s, t = 0.375, 0.25 >>> new_s, new_t = newton_refine(s, curve1, t, curve2) >>> 64.0 * (new_s - s) -9.0 >>> 64.0 * (new_t - t) 18.0
For “typical” curves, we converge to a solution quadratically. This means that the number of correct digits doubles every iteration (until machine precision is reached).
>>> curve1 = bezier.Curve(np.array([ ... [0.0, 0.0], ... [0.25, 2.0], ... [0.5, -2.0], ... [0.75, 2.0], ... [1.0, 0.0], ... ])) >>> curve2 = bezier.Curve(np.array([ ... [0.0, 1.0], ... [0.25, 0.5], ... [0.5, 0.5], ... [0.75, 0.5], ... [1.0, 0.0], ... ])) >>> # The expected intersection is the only real root of >>> # 28 s^3 - 30 s^2 + 9 s - 1. >>> rts = np.roots([28, -30, 9, -1]) >>> expected, = rts[rts.imag == 0.0].real >>> s, t = 0.625, 0.625 >>> np.log2(abs(expected - s)) -4.399... >>> s, t = newton_refine(s, curve1, t, curve2) >>> np.log2(abs(expected - s)) -7.901... >>> s, t = newton_refine(s, curve1, t, curve2) >>> np.log2(abs(expected - s)) -16.010... >>> s, t = newton_refine(s, curve1, t, curve2) >>> np.log2(abs(expected - s)) -32.110... >>> s, t = newton_refine(s, curve1, t, curve2) >>> np.log2(abs(expected - s)) -50.0
However, when the intersection occurs at a point of tangency, the convergence becomes linear. This means that the number of correct digits added each iteration is roughly constant.
>>> curve1 = bezier.Curve(np.array([ ... [0.0, 0.0], ... [0.5, 1.0], ... [1.0, 0.0], ... ])) >>> curve2 = bezier.Curve(np.array([ ... [0.0, 0.5], ... [1.0, 0.5], ... ])) >>> expected = 0.5 >>> s, t = 0.375, 0.375 >>> np.log2(abs(expected - s)) -3.0 >>> s, t = newton_refine(s, curve1, t, curve2) >>> np.log2(abs(expected - s)) -4.0 >>> s, t = newton_refine(s, curve1, t, curve2) >>> np.log2(abs(expected - s)) -5.0 >>> s, t = newton_refine(s, curve1, t, curve2) >>> np.log2(abs(expected - s)) -6.0 >>> s, t = newton_refine(s, curve1, t, curve2) >>> np.log2(abs(expected - s)) -7.0 >>> s, t = newton_refine(s, curve1, t, curve2) >>> np.log2(abs(expected - s)) -8.0
Parameters: Returns: The refined parameters from a single Newton step.
Return type:
-
bezier._intersection_helpers.
segment_intersection
(start0, end0, start1, end1)¶ Determine the intersection of two line segments.
Assumes each line is parametric
\[\begin{split}\begin{alignat*}{2} L_0(s) &= S_0 (1 - s) + E_0 s &&= S_0 + s \Delta_0 \\ L_1(t) &= S_1 (1 - t) + E_1 t &&= S_1 + t \Delta_1. \end{alignat*}\end{split}\]To solve \(S_0 + s \Delta_0 = S_1 + t \Delta_1\), we use the cross product:
\[\left(S_0 + s \Delta_0\right) \times \Delta_1 = \left(S_1 + t \Delta_1\right) \times \Delta_1 \Longrightarrow s \left(\Delta_0 \times \Delta_1\right) = \left(S_1 - S_0\right) \times \Delta_1.\]Similarly
\[\Delta_0 \times \left(S_0 + s \Delta_0\right) = \Delta_0 \times \left(S_1 + t \Delta_1\right) \Longrightarrow \left(S_1 - S_0\right) \times \Delta_0 = \Delta_0 \times \left(S_0 - S_1\right) = t \left(\Delta_0 \times \Delta_1\right).\]Note
Since our points are in \(\mathbf{R}^2\), the “traditional” cross-product in \(\mathbf{R}^3\) will always point in the \(z\) direction, so in the above we mean the \(z\) component of the cross product, rather than the entire vector.
For example, the diagonal lines
\[\begin{split}\begin{align*} L_0(s) &= \left[\begin{array}{c} 0 \\ 0 \end{array}\right] (1 - s) + \left[\begin{array}{c} 2 \\ 2 \end{array}\right] s \\ L_1(t) &= \left[\begin{array}{c} -1 \\ 2 \end{array}\right] (1 - t) + \left[\begin{array}{c} 1 \\ 0 \end{array}\right] t \end{align*}\end{split}\]intersect at \(L_0\left(\frac{1}{4}\right) = L_1\left(\frac{3}{4}\right) = \frac{1}{2} \left[\begin{array}{c} 1 \\ 1 \end{array}\right]\).
>>> start0 = np.array([[0.0, 0.0]]) >>> end0 = np.array([[2.0, 2.0]]) >>> start1 = np.array([[-1.0, 2.0]]) >>> end1 = np.array([[1.0, 0.0]]) >>> s, t, _ = segment_intersection(start0, end0, start1, end1) >>> s 0.25 >>> t 0.75
Taking the parallel (but different) lines
\[\begin{split}\begin{align*} L_0(s) &= \left[\begin{array}{c} 1 \\ 0 \end{array}\right] (1 - s) + \left[\begin{array}{c} 0 \\ 1 \end{array}\right] s \\ L_1(t) &= \left[\begin{array}{c} -1 \\ 3 \end{array}\right] (1 - t) + \left[\begin{array}{c} 3 \\ -1 \end{array}\right] t \end{align*}\end{split}\]we should be able to determine that the lines don’t intersect, but this function is not meant for that check:
>>> start0 = np.array([[1.0, 0.0]]) >>> end0 = np.array([[0.0, 1.0]]) >>> start1 = np.array([[-1.0, 3.0]]) >>> end1 = np.array([[3.0, -1.0]]) >>> _, _, success = segment_intersection(start0, end0, start1, end1) >>> success False
Instead, we use
parallel_different()
:>>> parallel_different(start0, end0, start1, end1) True
Parameters: - start0 (numpy.ndarray) – A 1x2 NumPy array that is the start vector \(S_0\) of the parametric line \(L_0(s)\).
- end0 (numpy.ndarray) – A 1x2 NumPy array that is the end vector \(E_0\) of the parametric line \(L_0(s)\).
- start1 (numpy.ndarray) – A 1x2 NumPy array that is the start vector \(S_1\) of the parametric line \(L_1(s)\).
- end1 (numpy.ndarray) – A 1x2 NumPy array that is the end vector \(E_1\) of the parametric line \(L_1(s)\).
Returns: Pair of \(s_{\ast}\) and \(t_{\ast}\) such that the lines intersect: \(L_0\left(s_{\ast}\right) = L_1\left(t_{\ast}\right)\) and then a boolean indicating if an intersection was found.
Return type:
-
bezier._intersection_helpers.
parallel_different
(start0, end0, start1, end1)¶ Checks if two parallel lines ever meet.
Meant as a back-up when
segment_intersection()
fails.Note
This function assumes but never verifies that the lines are parallel.
In the case that the segments are parallel and lie on the exact same line, finding a unique intersection is not possible. However, if they are parallel but on different lines, then there is a guarantee of no intersection.
In
segment_intersection()
, we utilized the normal form of the lines (via the cross product):\[\begin{split}\begin{align*} L_0(s) \times \Delta_0 &\equiv S_0 \times \Delta_0 \\ L_1(t) \times \Delta_1 &\equiv S_1 \times \Delta_1 \end{align*}\end{split}\]So, we can detect if \(S_1\) is on the first line by checking if
\[S_0 \times \Delta_0 \stackrel{?}{=} S_1 \times \Delta_0.\]If it is not on the first line, then we are done, the segments don’t meet:
>>> # Line: y = 1 >>> start0 = np.array([[0.0, 1.0]]) >>> end0 = np.array([[1.0, 1.0]]) >>> # Vertical shift up: y = 2 >>> start1 = np.array([[-1.0, 2.0]]) >>> end1 = np.array([[3.0, 2.0]]) >>> parallel_different(start0, end0, start1, end1) True
If \(S_1\) is on the first line, we want to check that \(S_1\) and \(E_1\) define parameters outside of \(\left[0, 1\right]\). To compute these parameters:
\[L_1(t) = S_0 + s_{\ast} \Delta_0 \Longrightarrow s_{\ast} = \frac{\Delta_0^T \left( L_1(t) - S_0\right)}{\Delta_0^T \Delta_0}.\]For example, the segments \(\left[0, 1\right]\) and \(\left[\frac{3}{2}, 2\right]\) (via \(S_1 = S_0 + \frac{3}{2} \Delta_0\) and \(E_1 = S_0 + 2 \Delta_0\)) don’t meet:
>>> start0 = np.array([[0.0, 1.0]]) >>> delta0 = np.array([[-1.0, 2.0]]) >>> end0 = start0 + 1.0 * delta0 >>> start1 = start0 + 1.5 * delta0 >>> end1 = start0 + 2.0 * delta0 >>> parallel_different(start0, end0, start1, end1) True
but if the segments, overlap, like \(\left[0, 1\right]\) and \(\left[-1, \frac{1}{2}\right]\), the segments meet:
>>> start1 = start0 - 1.0 * delta0 >>> end1 = start0 + 0.5 * delta0 >>> parallel_different(start0, end0, start1, end1) False
Similarly, if the second interval completely contains the first, the segments meet:
>>> start1 = start0 + 3.0 * delta0 >>> end1 = start0 - 2.0 * delta0 >>> parallel_different(start0, end0, start1, end1) False
Note
This method doesn’t currently allow wiggle room around the desired value, i.e. the two values must be bitwise identical. However, the most “correct” version of this function likely should allow for some round off.
Parameters: - start0 (numpy.ndarray) – A 1x2 NumPy array that is the start vector \(S_0\) of the parametric line \(L_0(s)\).
- end0 (numpy.ndarray) – A 1x2 NumPy array that is the end vector \(E_0\) of the parametric line \(L_0(s)\).
- start1 (numpy.ndarray) – A 1x2 NumPy array that is the start vector \(S_1\) of the parametric line \(L_1(s)\).
- end1 (numpy.ndarray) – A 1x2 NumPy array that is the end vector \(E_1\) of the parametric line \(L_1(s)\).
Returns: Indicating if the lines are different.
Return type:
-
class
bezier.curve.
Curve
(nodes, start=0.0, end=1.0, root=None, _copy=True)¶ Bases:
bezier._base.Base
Represents a Bézier curve.
We take the traditional definition: a Bézier curve is a mapping from \(s \in \left[0, 1\right]\) to convex combinations of points \(v_0, v_1, \ldots, v_n\) in some vector space:
\[B(s) = \sum_{j = 0}^n \binom{n}{j} s^j (1 - s)^{n - j} \cdot v_j\]>>> import bezier >>> nodes = np.array([ ... [0.0 , 0.0], ... [0.625, 0.5], ... [1.0 , 0.5], ... ]) >>> curve = bezier.Curve(nodes) >>> curve <Curve (degree=2, dimension=2)>
Parameters: - nodes (numpy.ndarray) – The nodes in the curve. The rows represent each node while the columns are the dimension of the ambient space.
- start (
Optional
[float
]) – The beginning of the sub-interval that this curve represents. - end (
Optional
[float
]) – The end of the sub-interval that this curve represents. - root (
Optional
[Curve
]) – The root curve that contains this current curve. - _copy (bool) – Flag indicating if the nodes should be copied before
being stored. Defaults to
True
since callers may freely mutatenodes
after passing in.
-
length
¶ float: The length of the current curve.
-
start
¶ float: Start of sub-interval this curve represents.
This value is used to track the current curve in the re-parameterization / subdivision process. The curve is still defined on the unit interval, but this value illustrates how this curve relates to a “parent” curve. For example:
>>> nodes = np.array([ ... [0.0, 0.0], ... [1.0, 2.0], ... ]) >>> curve = bezier.Curve(nodes) >>> curve <Curve (degree=1, dimension=2)> >>> left, right = curve.subdivide() >>> left <Curve (degree=1, dimension=2, start=0, end=0.5)> >>> right <Curve (degree=1, dimension=2, start=0.5, end=1)> >>> _, mid_right = left.subdivide() >>> mid_right <Curve (degree=1, dimension=2, start=0.25, end=0.5)> >>> mid_right.nodes array([[ 0.25, 0.5 ], [ 0.5 , 1. ]])
-
root
¶ Curve: The “root” curve that contains the current curve.
This indicates that the current curve is a section of the “root” curve. For example:
>>> _, right = curve.subdivide() >>> right <Curve (degree=2, dimension=2, start=0.5, end=1)> >>> right.root is curve True >>> right.evaluate(0.0) == curve.evaluate(0.5) array([ True, True], dtype=bool) >>> >>> mid_left, _ = right.subdivide() >>> mid_left <Curve (degree=2, dimension=2, start=0.5, end=0.75)> >>> mid_left.root is curve True >>> mid_left.evaluate(1.0) == curve.evaluate(0.75) array([ True, True], dtype=bool)
-
evaluate
(s)¶ Evaluate \(B(s)\) along the curve.
See
evaluate_multi()
for more details.>>> nodes = np.array([ ... [0.0 , 0.0], ... [0.625, 0.5], ... [1.0 , 0.5], ... ]) >>> curve = bezier.Curve(nodes) >>> curve.evaluate(0.75) array([ 0.796875, 0.46875 ])
Parameters: s (float) – Parameter along the curve. Returns: The point on the curve (as a one dimensional NumPy array). Return type: numpy.ndarray
-
evaluate_multi
(s_vals)¶ Evaluate \(B(s)\) for multiple points along the curve.
This is done by first evaluating each member of the Bernstein basis at each value in
s_vals
and then applying those to the control points for the current curve.This is done instead of using de Casteljau’s algorithm. Implementing de Casteljau is problematic because it requires a choice between one of two methods:
- vectorize operations of the form \((1 - s)v + s w\),
which requires a copy of the curve’s control points for
each value in
s_vals
- avoid vectorization and compute each point in serial
Instead, we can use vectorized operations to build up the Bernstein basis values.
>>> nodes = np.array([ ... [0.0, 0.0, 0.0], ... [1.0, 2.0, 3.0], ... ]) >>> curve = bezier.Curve(nodes) >>> curve <Curve (degree=1, dimension=3)> >>> s_vals = np.linspace(0.0, 1.0, 5) >>> curve.evaluate_multi(s_vals) array([[ 0. , 0. , 0. ], [ 0.25, 0.5 , 0.75], [ 0.5 , 1. , 1.5 ], [ 0.75, 1.5 , 2.25], [ 1. , 2. , 3. ]])
Parameters: s_vals (numpy.ndarray) – Parameters along the curve (as a 1D array). Returns: The points on the curve. As a two dimensional NumPy array, with the rows corresponding to each s
value and the columns to the dimension.Return type: numpy.ndarray - vectorize operations of the form \((1 - s)v + s w\),
which requires a copy of the curve’s control points for
each value in
-
plot
(num_pts, ax=None, show=False)¶ Plot the current curve.
Parameters: - num_pts (int) – Number of points to plot.
- ax (
Optional
[matplotlib.artist.Artist
]) – matplotlib axis object to add plot to. - show (
Optional
[bool
]) – Flag indicating if the plot should be shown.
Returns: The axis containing the plot. This may be a newly created axis.
Return type: Raises: NotImplementedError
– If the curve’s dimension is not2
.
-
subdivide
()¶ Split the curve \(\gamma(s)\) into a left and right half.
Takes the interval \(\left[0, 1\right]\) and splits the curve into \(\gamma_1 = \gamma\left(\left[0, \frac{1}{2}\right]\right)\) and \(\gamma_2 = \gamma\left(\left[\frac{1}{2}, 1\right]\right)\). In order to do this, also reparameterizes the curve, hence the resulting left and right halves have new nodes.
>>> nodes = np.array([ ... [0.0 , 0.0], ... [1.25, 3.0], ... [2.0 , 1.0], ... ]) >>> curve = bezier.Curve(nodes) >>> left, right = curve.subdivide() >>> left <Curve (degree=2, dimension=2, start=0, end=0.5)> >>> left.nodes array([[ 0. , 0. ], [ 0.625, 1.5 ], [ 1.125, 1.75 ]]) >>> right <Curve (degree=2, dimension=2, start=0.5, end=1)> >>> right.nodes array([[ 1.125, 1.75 ], [ 1.625, 2. ], [ 2. , 1. ]])
Returns: The left and right sub-curves. Return type: Tuple
[Curve
,Curve
]
-
intersect
(other)¶ Find the points of intersection with another curve.
>>> nodes1 = np.array([ ... [0.0 , 0.0 ], ... [0.375, 0.75 ], ... [0.75 , 0.375], ... ]) >>> curve1 = bezier.Curve(nodes1) >>> nodes2 = np.array([ ... [0.5, 0.0], ... [0.5, 3.0], ... ]) >>> curve2 = bezier.Curve(nodes2) >>> curve1.intersect(curve2) array([[ 0.5, 0.5]])
Parameters: other (Curve) – Other curve to intersect with.
Returns: Array of intersection points (possibly empty).
Return type: Raises: TypeError
– Ifother
is not a curve.NotImplementedError
– If both curves aren’t two-dimensional.
-
__eq__
(other)¶ Check equality against another shape.
Returns: Boolean indicating if the shapes are the same. Return type: bool
-
__ne__
(other)¶ Check inequality against another shape.
Returns: Boolean indicating if the shapes are not the same. Return type: bool
-
copy
()¶ Make a copy of the current shape.
Returns: Instance of the current shape.
-
degree
¶ int: The degree of the current shape.
-
dimension
¶ int: The dimension that the shape lives in.
For example, if the shape lives in \(\mathbf{R}^3\), then the dimension is
3
.
-
elevate
()¶ Return a degree-elevated version of the current curve.
Does this by converting the current nodes \(v_0, \ldots, v_n\) to new nodes \(w_0, \ldots, w_{n + 1}\) where
\[\begin{split}\begin{align*} w_0 &= v_0 \\ w_j &= \frac{j}{n + 1} v_{j - 1} + \frac{n + 1 - j}{n + 1} v_j \\ w_{n + 1} &= v_n \end{align*}\end{split}\]Returns: The degree-elevated curve. Return type: Curve
-
nodes
¶ numpy.ndarray: The nodes that define the current shape.