nurbspy.jax.nurbs_curve module#

class nurbspy.jax.nurbs_curve.NurbsCurve(control_points=None, weights=None, degree=None, knots=None)[source]#

Bases: Module

Create a NURBS (Non-Uniform Rational Basis Spline) curve object

Parameters:
control_pointsndarray with shape (ndim, n+1)

Array containing the coordinates of the control points The first dimension of P spans the coordinates of the control points (any number of dimensions) The second dimension of P spans the u-direction control points (0,1,…,n)

weightsndarray with shape (n+1,)

Array containing the weight of the control points

degreeint

Degree of the basis polynomials

knotsndarray with shape (r+1=n+p+2,)

The knot vector in the u-direction Set the multiplicity of the first and last entries equal to p+1 to obtain a clamped NURBS

Methods

get_arclength([u1, u2, n_points])

Compute the arc length of a parametric curve in the interval [u1,u2] using numerical quadrature

get_binormal(u)

Evaluate the unit binormal vector along the curve for the given u-parameterization.

get_curvature(u)

Evaluate the curvature of the curve for the input u-parametrization

get_derivative(u, order)

Evaluate the derivative of the curve for the input u-parametrization

get_normal(u)

Evaluate the unit normal vector along the curve for the given u-parameterization.

get_tangent(u)

Evaluate the unit tangent vector along the curve for the given u-parameterization.

get_torsion(u)

Evaluate the torsion of the curve for the input u-parametrization

get_value(u)

Evaluate the coordinates of the curve for the input u parametrization

plot([fig, ax, curve, control_points, ...])

Create a plot and return the figure and axes handles

plot_control_points(fig, ax[, linewidth, ...])

Plot the control points of the NURBS curve

plot_curve(fig, ax[, linewidth, linestyle, ...])

Plot the coordinates of the NURBS curve

plot_frenet_serret(fig, ax[, frame_number, ...])

Plot some Frenet-Serret reference frames along the NURBS curve

project_point_to_curve(Q[, max_iters])

Project a point onto the NURBS curve by solving the orthogonality condition.

project_points_to_curve(Q_all)

Vectorized projection of multiple points onto the NURBS curve.

rescale_plot(fig, ax)

Adjust the aspect ratio of the figure

plot_curvature

plot_torsion

Notes

This class includes methods to compute:

  • Curve coordinates for any number of dimensions

  • Analytic curve derivatives of any order and number of dimensions

  • The unitary tangent, normal and binormal vectors (Frenet-Serret reference frame) in 2D and 3D

  • The analytic curvature and torsion in 2D and 3D

  • The arc length of the curve in any number of dimensions.

    The arc length is compute by numerical quadrature using analytic derivative information

The class can be used to represent polynomial and rational Bézier, B-Spline and NURBS curves The type of curve depends on the initialization arguments

  • Polymnomial Bézier: Provide the array of control points

  • Rational Bézier: Provide the arrays of control points and weights

  • B-Spline: Provide the array of control points, degree and knot vector

  • NURBS: Provide the arrays of control points and weights, degree and knot vector

In addition, this class supports operations with real and complex numbers The data type used for the computations is detected from the data type of the arguments Using complex numbers can be useful to compute the derivative of the shape using the complex step method

References

The NURBS Book. See references to equations and algorithms throughout the code L. Piegl and W. Tiller Springer, second edition

Curves and Surfaces for CADGD. See references to equations the source code G. Farin Morgan Kaufmann Publishers, fifth edition

All references correspond to The NURBS book unless it is explicitly stated that they come from Farin’s book

P: Array#
U: Array#
U_mults: Array#
U_values: Array#
W: Array#
curve_type: str#
get_arclength(u1=0.0, u2=1.0, n_points=41)[source]#

Compute the arc length of a parametric curve in the interval [u1,u2] using numerical quadrature

The definition of the arc length is given by equation 10.3 (Farin’s textbook)

Parameters:
u1scalar

Lower limit of integration for the arc length computation

u2scalar

Upper limit of integration for the arc length computation

Returns:
Lscalar

Arc length of NURBS curve in the interval [u1, u2]

get_binormal(u)[source]#

Evaluate the unit binormal vector along the curve for the given u-parameterization.

The binormal is defined as:

b(u) = (C’(u) x C’’(u)) / ||C’(u) x C’’(u)||

Parameters:
uscalar or ndarray (N,)

Parametric coordinates.

Returns:
binormalndarray (ndim, N)

Unit binormal vectors.

get_curvature(u)[source]#

Evaluate the curvature of the curve for the input u-parametrization

The definition of the curvature is given by equation 10.7 (Farin’s textbook)

Parameters:
uscalar or ndarray with shape (N,)

Scalar or array containing the u-parameter used to evaluate the curvature

Returns:
curvaturescalar or ndarray with shape (N, )

Scalar or array containing the curvature of the curve

get_derivative(u, order)[source]#

Evaluate the derivative of the curve for the input u-parametrization

Parameters:
uscalar or ndarray with shape (N,)

Scalar or array containing the u-parameter used to evaluate the curve

orderinteger

Order of the derivative

Returns:
dCndarray with shape (ndim, N)

Array containing the derivative of the desired order The first dimension of dC spans the (x,y,z) coordinates The second dimension of dC spans the u parametrization sample points

get_normal(u)[source]#

Evaluate the unit normal vector along the curve for the given u-parameterization.

For 2D or 3D curves, the normal is defined as:

n(u) = (C’(u) x (C’’(u) x C’(u))) / ||C’(u) x (C’’(u) x C’(u))||

Parameters:
uscalar or ndarray (N,)

Parametric coordinates.

Returns:
normalndarray (ndim, N)

Unit normal vectors.

get_tangent(u)[source]#

Evaluate the unit tangent vector along the curve for the given u-parameterization.

The tangent is defined as:

t(u) = C’(u) / ||C’(u)||

Parameters:
uscalar or ndarray (N,)

Parametric coordinates.

Returns:
tangentndarray (ndim, N)

Unit tangent vectors.

get_torsion(u)[source]#

Evaluate the torsion of the curve for the input u-parametrization

The definition of the torsion is given by equation 10.8 (Farin’s textbook)

Parameters:
uscalar or ndarray with shape (N,)

Scalar or array containing the u-parameter used to evaluate the torsion

Returns:
torsionscalar or ndarray with shape (N, )

Scalar or array containing the torsion of the curve

get_value(u)[source]#

Evaluate the coordinates of the curve for the input u parametrization

Parameters:
uscalar or ndarray with shape (N,)

Parameter used to evaluate the curve

Returns:
Cndarray with shape (ndim, N)

Array containing the coordinates of the curve The first dimension of C spans the (x,y,z) coordinates The second dimension of C spans the u parametrization sample points

ndim: int#
p: int#
plot(fig=None, ax=None, curve=True, control_points=True, frenet_serret=False, axis_off=False, ticks_off=False)[source]#

Create a plot and return the figure and axes handles

plot_control_points(fig, ax, linewidth=1.25, linestyle='-.', color='red', markersize=5, markerstyle='o')[source]#

Plot the control points of the NURBS curve

plot_curvature(fig=None, ax=None, color='black', linestyle='-')[source]#
plot_curve(fig, ax, linewidth=1.5, linestyle='-', color='black', u1=0.0, u2=1.0)[source]#

Plot the coordinates of the NURBS curve

plot_frenet_serret(fig, ax, frame_number=5, frame_scale=0.1)[source]#

Plot some Frenet-Serret reference frames along the NURBS curve

plot_torsion(fig=None, ax=None, color='black', linestyle='-')[source]#
project_point_to_curve(Q, max_iters=32)[source]#

Project a point onto the NURBS curve by solving the orthogonality condition.

The projection point C(u*) minimizes the squared Euclidean distance to Q:

f(u) = ||C(u) - Q||²

The stationary condition is obtained from:

f’(u) = 2 (C(u) - Q) · C’(u) = 0

which ensures that the vector from the curve to the point is orthogonal to the tangent of the curve at the projection point.

The nonlinear equation is solved with a bounded Newton method (optimistix.Newton), ensuring that u_star remains within [0, 1].

The initial guess is determined heuristically from the control polygon (see _initial_guess_from_polygon).

Parameters:
Qarray_like, shape (ndim,)

Coordinates of the point to be projected onto the curve.

max_itersint, optional

Maximum number of Newton iterations used by the solver. Default is 32.

Returns:
u_starfloat

Parameter value in [0, 1] corresponding to the orthogonal projection of Q onto the curve.

project_points_to_curve(Q_all)[source]#

Vectorized projection of multiple points onto the NURBS curve.

Parameters:
Q_allarray_like, shape (ndim, n_points)

Points to project onto the curve.

Returns:
u_allndarray, shape (n_points,)

Parameter values of the orthogonal projections.

rescale_plot(fig, ax)[source]#

Adjust the aspect ratio of the figure

nurbspy.jax.nurbs_curve.binomial_coeff(n, k)[source]#

JAX-compatible binomial coefficient C(n, k) using the gamma function.

nurbspy.jax.nurbs_curve.compute_all_bspline_derivatives(P, p, U, u, up_to_order)[source]#

Compute all analytic derivatives of a polynomial B-spline curve up to a specified order.

This function evaluates the curve and its parametric derivatives using the analytic derivatives of the B-spline basis functions. For each derivative order k, the derivative of the curve is given by

C^(k)(u) = Σ_i P_i * d^k N_{i,p}(u) / du^k

Parameters:
Pndarray (ndim, n+1)

Control point coordinates.

pint

Degree of the B-spline.

Undarray (n+p+2,)

Knot vector.

undarray (Nu,)

Parametric evaluation points.

up_to_orderint

Maximum derivative order to compute (0 ≤ order ≤ p).

Returns:
bspline_derivativesndarray (up_to_order+1, ndim, Nu)

Derivatives of the B-spline curve, where bspline_derivatives[k, :, :] = d^k C(u) / du^k.

nurbspy.jax.nurbs_curve.compute_all_nurbs_derivatives(P, W, p, U, u, up_to_order)[source]#

Compute all analytic derivatives of a NURBS curve up to a specified order.

This function extends the polynomial B-spline derivative computation to rational NURBS curves by applying the quotient rule recursively, following Algorithm A4.2 from The NURBS Book (Piegl & Tiller, 2nd ed.):

Parameters:
Pndarray (ndim, n+1)

Control point coordinates.

Wndarray (n+1,)

Control point weights.

pint

Degree of the NURBS.

Undarray (n+p+2,)

Knot vector.

undarray (Nu,)

Parametric evaluation points.

up_to_orderint

Maximum derivative order to compute (0 ≤ order ≤ p).

Returns:
nurbs_derivativesndarray (up_to_order+1, ndim, Nu)

Derivatives of the NURBS curve

nurbspy.jax.nurbs_curve.compute_bspline_coordinates(P, p, U, u)[source]#

Evaluate the coordinates of a B-spline curve for a given parameter u.

This function computes the coordinates of a polynomial B-spline curve using the standard basis expansion (Equation 3.1 in The NURBS Book). The implementation corresponds to Algorithm A3.1.

Parameters:
Pndarray (ndim, n+1)

Array of control point coordinates. The first dimension spans spatial coordinates (x, y, z, …), and the second spans the control points along the curve (0, 1, …, n).

pint

Degree of the B-spline basis functions.

Undarray (r+1 = n + p + 2,)

Knot vector in the u-direction. The first and last entries typically have multiplicity p + 1 for a clamped spline.

uscalar or ndarray (N,)

Parametric coordinate(s) at which to evaluate the curve.

Returns:
Cndarray (ndim, N)

Coordinates of the evaluated B-spline curve points. The first dimension spans the spatial coordinates, and the second spans the parametric evaluation points u.

Notes

  • The computation is vectorized and uses matrix multiplication for efficiency.

  • For rational curves (NURBS), use compute_nurbs_coordinates instead.

nurbspy.jax.nurbs_curve.compute_nurbs_coordinates(P, W, p, U, u)[source]#

Evaluate the coordinates of a NURBS (Non-Uniform Rational B-Spline) curve for a given parameter u.

This function computes the coordinates of the NURBS curve in homogeneous space using the standard B-spline basis and the control point weights (Equation 4.5 in The NURBS Book), and then maps them back to ordinary space via the rational perspective division (Equation 1.16). The implementation corresponds to Algorithm A4.1 from The NURBS Book.

Parameters:
Pndarray (ndim, n+1)

Array of control point coordinates. The first dimension spans spatial coordinates (x, y, z, …), and the second spans the control points along the curve (0, 1, …, n).

Wndarray (n+1,)

Weights associated with each control point.

pint

Degree of the B-spline basis functions.

Undarray (r+1 = n + p + 2,)

Knot vector in the u-direction. The first and last entries typically have multiplicity p + 1 for a clamped spline.

uscalar or ndarray (N,)

Parametric coordinate(s) at which to evaluate the curve.

Returns:
Cndarray (ndim, N)

Coordinates of the evaluated curve points. The first dimension spans the spatial coordinates, and the second spans the parametric evaluation points u.

Notes

  • The computation is fully vectorized and uses matrix multiplication for efficiency.

  • The function supports arbitrary spatial dimensions.

nurbspy.jax.nurbs_curve.merge_nurbs_curves(nurbs_list)[source]#

Merge multiple NURBS curves into a single continuous curve with C⁰ continuity.

Parameters:
nurbs_listlist of NurbsCurve

List of NURBS curve instances to merge. All must have the same polynomial degree.

Returns:
mergedNurbsCurve

A new NURBS curve representing the concatenation of all input curves.

Notes

  • Each curve is mapped to a subinterval of [0, 1] with equal length.

  • Adjacent curves are connected with multiplicity p + 1 (C⁰ continuity).

  • A small epsilon offset is applied at internal joins to avoid Gmsh issues related to repeated knots with multiplicity exactly equal to degree.