83. plugins.nurbs — Using NURBS in pyFormex.

The nurbs module defines functions and classes to manipulate NURBS curves and surfaces in pyFormex. The most important classes and functions for the end user are:

83.1. Classes defined in module plugins.nurbs

class plugins.nurbs.NurbsCurve(control, degree=None, wts=None, knots=None, blended=True, closed=False, wrap=True, norm_urange=True)[source]

A NURBS curve.

This class can represent a NURBS curve as well as its special case of a non-rational (or simple) B-spline curve. The difference sits only in the weights attributed to the control points: if they are constant, the interpolation functions are simple polynoms and the curve is a B-spline. If variable weights are attributed, the curve becomes a more general Non-Uniform Rational B-Spline or NURBS.

A B-spline is defined by nctrl control points, a degree (>=1) and a sequence of nknots = nctrl+degree+1 parametric values called knots. The knots divide the curve in parametric regions where the mathematical representation of the curve is constant. They form a non-descending sequence , but values may be repeated to model appropriate discontinuities. Given a set of control points and a knot vector, the degree of the curve is obviously degree = nknots-nctrl-1. In literature often the order of a NURBS is used. This is just order = degree+1 = nknots-nctrl. The order is also the minimum number of control points.

To model a NURBS curve, we add a weight wi to each control point and use 4-dimensional coordinates (xi*wi, yi*wi, zi*wi, wi) for each point i. Internally, the NurbsCurve class uses 4-dim points, with default weights set to 1. Before returning, the coordinates are divided At the endInput and output are always 3-dimensional points.

Parameters:
  • control (Coords-like (nctrl,3)) – The vertices of the control polygon.

  • degree (int) – The degree of the Nurbs curve. If not specified, it is derived from the length of the knot vector (knots).

  • wts (float array (nctrl)) – Weights to be attributed to the control points. Default is to attribute a weight 1.0 to all points. Using different weights allows for more versatile modeling (like perfect circles and arcs.)

  • knots (KnotVector | list of floats) – The nknots knot values to be used. If a list, the values should be in ascending order. Identical values have to be repeated to their multiplicity. The values are only defined upon a multiplicative constant and will be normalized to set the last value to 1. If degree is specified, default values are constructed automatically by calling genKnotVector(). If no knots are given and no degree is specified, the degree is set to the nctrl-1 if the curve is blended. If not blended, the degree is not set larger than 3.

  • closed (bool, optional) – Determines whether the curve is closed. Default is False. The functionality of closed NurbsCurves is currently limited.

  • blended (bool, optional) – Determines that the curve is blended. Default is True. Set blended==False to define a nonblended curve. A nonblended curve is a chain of independent curves, Bezier curves if the weights are all ones. See also decompose(). The number of control points should be a multiple of the degree, plus one. This parameter is only used if no knots are specified.

property coords4

The 4-dim coordinates

property coords

The 3-dim coordinates

property knots

Return the full list of knot values

property nctrl

The number of control points

property nknots

The number of knots

property order

The order of the Nurbs curve = nknots - nctrl = degree + 1

urange()[source]

Return the parameter range on which the curve is defined.

Returns a (2,) float array with the minimum and maximum parameter value for which the curve is defined.

isClamped()[source]

Return True if the NurbsCurve uses a clamped knot vector.

A clamped knot vector has a multiplicity p+1 for the first and last knot. All our generated knot vectors are clamped.

isUniform()[source]

Return True if the NurbsCurve has a uniform knot vector.

A uniform knot vector has a constant spacing between the knot values.

isRational()[source]

Return True if the NurbsCurve is rational.

The curve is rational if the weights are not constant. The curve is polynomial if the weights are constant.

Returns True for a rational curve, False for a polynomial curve.

isBlended()[source]

Return True if the NurbsCurve is blended.

An clamped NurbsCurve is unblended (or decomposed) if it consists of a chain of independent Bezier curves. Such a curve has multiplicity p for all internal knots and p+1 for the end knots of an open curve. Any other NurbsCurve is blended.

Returns True for a blended curve, False for an unblended one.

Note: for testing whether an unclamped curve is blended or not, first clamp it.

bbox()[source]

Return the bounding box of the NURBS curve.

copy()[source]

Return a (deep) copy of self.

Changing the copy will not change the original.

pointsAt(u)[source]

Return the points on the Nurbs curve at given parametric values.

Parameters:

u (float array_like (nu,)) – The parametric values at which a point is to be placed. Note that valid points are only obtained for parameter values in the range().

Returns:

Coords (nu, 3) – The coordinates of nu points on the curve, at the specified parametric values.

Examples

>>> N = NurbsCurve(control=Coords('0121'), degree=3)
>>> N.pointsAt([0.0, 0.5, 1.0])
Coords([[0. , 0. , 0. ],
        [1. , 0.5, 0. ],
        [2. , 1. , 0. ]])
>>> N.pointsAt(0.5)
Coords([[1. , 0.5, 0. ]])
derivs(u, d=1)[source]

Returns the points and derivatives up to d at parameter values u

Parameters:
  • u (float array_like | int) –

    If a float array (nu,), these are the parameter values at which to compute points and derivatives.

    If an int, specifies the number of parameter values (nu) at which to evaluate the points and derivatives of the curve. The points are equally spaced in parameter space.

  • d (int) – The highest derivative to compute.

Returns:

float array (d+1,nu,3) – The coordinates of the points and the derivates up to the order d at those points.

Examples

>>> N = NurbsCurve(control=Coords('0121'), degree=3)
>>> N.derivs([0.0, 0.5, 1.0], 2)
Coords([[[ 0. ,  0. ,  0. ],
         [ 1. ,  0.5,  0. ],
         [ 2. ,  1. ,  0. ]],

        [[ 3. ,  0. ,  0. ],
         [ 1.5,  1.5,  0. ],
         [ 3. ,  0. ,  0. ]],

        [[-6. ,  6. ,  0. ],
         [ 0. ,  0. ,  0. ],
         [ 6. , -6. ,  0. ]]])
frenet(u)[source]

Compute Frenet vectors, curvature and torsion at parameter values u

Parameters:

u (float array_like | int) –

If a float array (nu,), these are the parameter values at which to compute points and derivatives.

If an int, specifies the number of parameter values (nu) at which to evaluate the points and derivatives of the curve. The points are equally spaced in parameter space.

Returns:

  • T (float array(nu,3)) – Normalized tangent vectors (nu,3) at nu points.

  • N (float array(nu,3)) – Normalized normal vectors (nu,3) at nu points.

  • B (float array(nu,3)) – Normalized binormal vectors (nu,3) at nu points.

  • k (float array(nu,3)) – Curvature of the curve (nu) at nu points.

  • t (float array(nu,3)) – Torsion of the curve (nu) at nu points.

curvature(u, torsion=False)[source]

Compute curvature and torsion at parameter values u

This is like frenet() but only returns the curvature and optionally the torsion.

knotPoints(multiple=False)[source]

Returns the points on the curve at the knot values.

If multiple is True, points are returned with their multiplicity. The default is to return the points just once.

Examples

>>> N = NurbsCurve(Coords('012141'), degree=3)
>>> print(N.knotu)
KnotVector(10): 0*4, 0.333333*1, 0.666667*1, 1*4
>>> N.knotPoints()
Coords([[0.    , 0.    , 0.    ],
        [1.1667, 0.75  , 0.    ],
        [1.8333, 0.75  , 0.    ],
        [3.    , 0.    , 0.    ]])
insertKnots(u)[source]

Insert knots in the Nurbs curve.

Parameters:

u (float array_like (nu,)) – A list of parameter values to be inserted into the curve’s knot vector. The control points are adapted to keep the curve unchanged.

Returns:

NurbsCurve – A new curve equivalent with the original but with the specified knot values inserted in the knot vector, and the control points adapted.

requireKnots(val, mul)[source]

Insert knots until the required multiplicity is reached.

Inserts knot values only if they are currently not there or their multiplicity is lower than the required one.

Parameters:
  • val (float array_like (nval,)) – The list of knot values required in the knot vector.

  • mul (int array_like (nval,)) – The list of multiplicities required for the knot values.

Returns:

NurbsCurve – A curve equivalent with the original but where the knot vector is guaranteed to contain the values in val with at least the corresponding multiplicity in mul. If all requirements were already fulfilled at the beginning, returns self.

removeKnot(u, m, tol=1e-05)[source]

Remove a knot from the knot vector of the Nurbs curve.

Parameters:
  • u (float) – The knot value to remove.

  • m (int) – How many times to remove the knot u. If negative, remove maximally.

  • tol (float) – Acceptable error (distance between old and new curve).

Returns:

NurbsCurve – A Nurbs curve equivalent (to the specified tolerance) with the original but with a knot vector where the value u has been removed m times, if possible, or else as many times as possible. The control points are adapted accordingly. Returns self if no value was removed.

Notes

removeKnot currently only works on open curves

removeAllKnots(tol=1e-05)[source]

Remove all removable knots.

Parameters:

tol (float) – Acceptable error (distance between old and new curve).

Returns:

NurbsCurve – An equivalent (if tol is small) NurbsCurve with all extraneous knots removed.

Notes

removeAllKnots() and blend() are aliases.

blend(tol=1e-05)

Remove all removable knots.

Parameters:

tol (float) – Acceptable error (distance between old and new curve).

Returns:

NurbsCurve – An equivalent (if tol is small) NurbsCurve with all extraneous knots removed.

Notes

removeAllKnots() and blend() are aliases.

decompose()[source]

Decompose a curve in subsequent Bezier curves.

Returns:

NurbsCurve – An equivalent unblended NurbsCurve.

See also

toCurve

convert the NurbsCurve to a BezierSpline or PolyLine

Notes

decompose() and unblend() are aliases.

unblend()

Decompose a curve in subsequent Bezier curves.

Returns:

NurbsCurve – An equivalent unblended NurbsCurve.

See also

toCurve

convert the NurbsCurve to a BezierSpline or PolyLine

Notes

decompose() and unblend() are aliases.

subCurve(u1, u2)[source]

Extract the subcurve between parameter values u1 and u2

Parameters:
  • u1 (float) – Parametric value of the start of the subcurve to extract. The value should be in urange().

  • u2 (float) – Parametric value of the end of the subcurve to extract. The value should be in urange() and > u1.

Returns:

NurbsCurve – A NurbsCurve containing only the part between u1 and u2 of the original.

clamp()[source]

Clamp the knot vector of the curve.

A clamped knot vector starts and ends with multiplicities p-1. See also isClamped().

Returns:

NurbsCurve – An equivalent curve with clamped knot vector, or self if the curve was already clamped.

Notes

The use of unclamped knot vectors is deprecated. This method is provided only as a convenient method to import curves from legacy systems using unclamped knot vectors.

unclamp()[source]

Unclamp the knot vector of the curve.

Warning

The use of unclamped knot vectors is deprecated.

Returns:

NurbsCurve – An equivalent curve with unclamped knot vector, or self if the curve was already unclamped.

toCurve(force_Bezier=False)[source]

Convert a (nonrational) NurbsCurve to a BezierSpline or PolyLine.

This decomposes the curve in a chain of Bezier curves and converts the chain to a BezierSpline or PolyLine.

This only works for nonrational NurbsCurves, as the BezierSpline and PolyLine classes do not allow homogeneous coordinates required for rational curves.

Returns:

BezierSpline | PolyLine – A BezierSpline (or PolyLine if degree is 1) that is equivalent with the NurbsCurve.

See also

unblend

decompose both rational and nonrational NurbsCurves

toBezier()[source]

Convert a (nonrational) NurbsCurve to a BezierSpline.

This is equivalent with toCurve(force_Bezier=True) and returns a BezierSpline in all cases.

elevateDegree(t=1)[source]

Elevate the degree of the Nurbs curve.

Parameters:

t (int) – How much to elevate the degree of the curve

Returns:

NurbsCurve – A NurbsCurve equivalent with the original but of a higher degree.

reduceDegree(t=1, tol=inf)[source]

Reduce the degree of the Nurbs curve.

Parameters:

t (int) – How much to reduce the degree (max. = degree-1)

Returns:

NurbsCurve – A NurbsCurve approximating the original but of a lower degree.

projectPoint(P, *, u0=None, nseed=20, eps1=1e-05, eps2=1e-05, maxit=50)[source]

Project a given point on the Nurbs curve.

This can also be used to determine the parameter value of a point lying on the curve.

Parameters:
  • P (coords_like (3,)) – A set of npts points in space.

  • u0 (float, optional) – Start value for the parameter of the projected point. Providing a value close to the correct foot point will speed up the operation. If not provided, a number of values is tried and the one giving the closest point to P is chosen. See nseed.

  • nseed (int) – Number of intervals to divide the parameter range into for guessing a suitable initial value u0. Only used if u0 is not provided.

Returns:

  • u (float) – Parameter value of the base point X of the projection of P on the NurbsCurve.

  • P (Coords (3,)) – The coordinates of the base point of the projection of P on the NurbsCurve.

  • d (float) – The distance of the point P from the NurbsCurve

Notes

Based on section 6.1 of The NurbsBook.

projectPoints(P, **kargs)[source]

Project multiple points

distancePoints(P, **kargs)[source]

Compute the distance of points P to the NurbsCurve

approx(ndiv=None, nseg=None, **kargs)[source]

Return a PolyLine approximation of the Nurbs curve

If no nseg is given, the curve is approximated by a PolyLine through equidistant ndiv+1 point in parameter space. These points may be far from equidistant in Cartesian space.

If nseg is given, a second approximation is computed with nseg straight segments of nearly equal length. The lengths are computed based on the first approximation with ndiv segments.

reverse()[source]

Return the reversed Nurbs curve.

The reversed curve is geometrically identical, but start and en point are interchanged and parameter values increase in the opposite direction.

actor(**kargs)[source]

Graphical representation

class plugins.nurbs.NurbsSurface(control, degree=(None, None), wts=None, knots=(None, None), closed=(False, False), blended=(True, True))[source]

A NURBS surface

The Nurbs surface is defined as a tensor product of NURBS curves in two parametrical directions u and v. The control points form a grid of (nctrlu,nctrlv) points. The other data are like those for a NURBS curve, but need to be specified as a tuple for the (u,v) directions.

The knot values are only defined upon a multiplicative constant, equal to the largest value. Sensible default values are constructed automatically by a call to the genKnotVector() function.

If no knots are given and no degree is specified, the degree is set to the number of control points - 1 if the curve is blended. If not blended, the degree is not set larger than 3.

Warning

This is a class under development!

urange()[source]

Return the u-parameter range on which the curve is defined.

Returns a (2,) float array with the minimum and maximum parameter value u for which the curve is defined.

vrange()[source]

Return the v-parameter range on which the curve is defined.

Returns a (2,) float array with the minimum and maximum parameter value v for which the curve is defined.

bbox()[source]

Return the bounding box of the NURBS surface.

pointsAt(u)[source]

Return the points on the Nurbs surface at given parametric values.

Parameters:

  • u: (nu,2) shaped float array: nu parametric values (u,v) at which a point is to be placed.

Returns (nu,3) shaped Coords with nu points at the specified parametric values.

derivs(u, m)[source]

Return points and derivatives at given parametric values.

Parameters:

  • u: (nu,2) shaped float array: nu parametric values (u,v) at which the points and derivatives are evaluated.

  • m: tuple of two int values (mu,mv). The points and derivatives up to order mu in u direction and mv in v direction are returned.

Returns:

(nu+1,nv+1,nu,3) shaped Coords with nu points at the specified parametric values. The slice (0,0,:,:) contains the points.

approx(ndiv=None, **kargs)[source]

Return a Quad4 Mesh approximation of the Nurbs surface

Parameters:

  • ndiv: number of divisions of the parametric space.

actor(**kargs)[source]

Graphical representation

class plugins.nurbs.Coords4(data=None, w=None, normalize=True, dtyp=<class 'numpy.float32'>, copy=False)[source]

A collection of points represented by their homogeneous coordinates.

While most of the pyFormex implementation is based on the 3D Cartesian coordinates class Coords, some applications may benefit from using 4-dimensional coordinates. The Coords4 class provides some basic functions, including conversion to and from cartesian coordinates. Through the conversion, all pyFormex functions defined on Coords, such as transformations, are possible.

Coords4 is implemented as a float type numpy.ndarray whose last axis has a length equal to 4. Each set of 4 values (x,y,z,w) along the last axis represents a single point in 3D space. The cartesian coordinates of the point are obtained by dividing the first three values by the fourth: (x/w, y/w, z/w). A zero w-value represents a point at infinity. Converting such points to Coords will result in Inf or NaN values in the resulting object.

The float datatype is only checked at creation time. It is the responsibility of the user to keep this consistent throughout the lifetime of the object.

Just like Coords, the class Coords4 is derived from numpy.ndarray.

Parameters:
  • data (array_like, optional) – An array of floats, with the length of its last axis not larger than 4. If equal to four, each tuple along the last axis represents a single point in homogeneous coordinates. If smaller than four, the last axis is expanded to four by adding zero values in the second and third position and unity values in the last position (the w-coordinate). If no data are given, a single point (0.,0.,0.,1.) is created.

  • w (array_like, optional) – If specified, the w values are used to denormalize the homogeneous data such that the last component becomes w.

  • dtyp (data-type, optional) – The datatype to be used. It not specified, the datatype of data is used, or the default Float (which is equivalent to numpy.float32).

  • copy (bool, optional) – If True, the data are copied. By default, the original data are used if possible, e.g. if a correctly shaped and typed numpy.ndarray is provided.

normalize()[source]

Normalize the homogeneous coordinates.

Two sets of homogeneous coordinates that differ only by a multiplicative constant refer to the same points in cartesian space. Normalization of the coordinates is a way to make the representation of a single point unique. Normalization is done so that the last component (w) is equal to 1.

The normalization of the coordinates is done in place.

Warning

Normalizing points at infinity will result in Inf or NaN values.

Examples

>>> X4 = Coords4([2, 3, 4, 2])
>>> X4.normalize()
>>> X4
Coords4([1. , 1.5, 2. , 1. ])
deNormalize(w)[source]

Denormalize the homogeneous coordinates.

This multiplies the homogeneous coordinates with the values w. w normally is a constant or an array with shape self.shape[:-1] + (1,). It then multiplies all 4 coordinates of a point with the same value, thus resulting in a denormalization while keeping the position of the point unchanged.

The denormalization of the coordinates is done in place. If the Coords4 object was normalized, it will have precisely w as its 4-th coordinate value after the call.

Examples

>>> X4 = Coords4([1, 2, 3, 1])
>>> X4.deNormalize(2)
>>> X4
Coords4([2., 4., 6., 2.])
toCoords()[source]

Convert homogeneous coordinates to cartesian coordinates.

Returns:

  • A Coords object with the cartesian coordinates

  • of the points. Points at infinity (w=0) will result in

  • Inf or NaN value. If there are no points at infinity, the

  • resulting Coords point set is equivalent to the

  • Coords4 input.

Examples

>>> Coords4([2, 3, 4, 2]).toCoords()
Coords([1. , 1.5, 2. ])
npoints()[source]

Return the total number of points in the Coords4.

ncoords()

Return the total number of points in the Coords4.

x()[source]

Return the x-plane

y()[source]

Return the y-plane

z()[source]

Return the z-plane

w()[source]

Return the w-plane

bbox()[source]

Return the bounding box of a set of points.

Returns the bounding box of the cartesian coordinates of the object.

actor(**kargs)[source]

Graphical representation

class plugins.nurbs.Geometry4[source]

This is a preliminary class intended to provide some transforms in 4D

class plugins.nurbs.KnotVector(data=None, val=None, mul=None)[source]

A knot vector.

A knot vector is sequence of float values sorted in ascending order. Values can occur multiple times. In the typical use case (NURBS) most values do indeed occur multiple times, and the multiplicity of the values is an important quantity. Therefore, the knot vector is stored in two arrays of the same length: val and mul.

val

The unique float values, in a strictly ascending sequence.

Type:

float array (nval)

mul

The multiplicity of each of the values in val.

Type:

int array (nval)

See also

genKnotVector

generate a knot vector for a NURBS curve

Examples

>>> K = KnotVector([0.,0.,0.,0.5,0.5,1.,1.,1.])
>>> print(K.val)
[0.  0.5 1. ]
>>> print(K.mul)
[3 2 3]
>>> print(K)
KnotVector(8): 0*3, 0.5*2, 1*3
>>> print([(v,m) for v,m in zip(K.val,K.mul)])
[(0.0, 3), (0.5, 2), (1.0, 3)]
>>> print(K.values())
[0.  0.  0.  0.5 0.5 1.  1.  1. ]
>>> print(K.nknots)
8
>>> K.index(0.5)
1
>>> K.span(1.0)
5
>>> K.mult(0.5)
2
>>> K.mult(0.7)
0
>>> K[4],K[-1]
(0.5, 1.0)
>>> K[4:6]
array([0.5, 1. ])
property nknots

Return the total number of knots

values()[source]

Return the full list of knot values

index(u)[source]

Find the index of knot value u.

If the value does not exist, a ValueError is raised.

vindex(j)[source]

Find the index of knot value with index j.

j should be in the range(0,self.nknots), else a ValueError is raised. This is equivalent to (but more efficient than):

self.index(self[j])
mult(u)[source]

Return the multiplicity of knot value u.

Returns:

int – The multiplicity of the knot value u, or 0 if u is not in the KnotVector.

span(u)[source]

Find the (first) index of knot value u in the full knot values vector.

If the value does not exist, a ValueError is raised.

copy()[source]

Return a copy of the KnotVector.

Changing the copy will not change the original.

reverse()[source]

Return the reverse knot vector.

Examples

>>> print(KnotVector([0,0,0,1,3,6,6,8,8,8]).reverse().values())
[0. 0. 0. 2. 2. 5. 7. 8. 8. 8.]

83.2. Functions defined in module plugins.nurbs

plugins.nurbs.genKnotVector(nctrl, degree, blended=True, closed=False)[source]

Compute a sensible knot vector for a Nurbs curve.

A knot vector is a sequence of non-decreasing parametric values. These values define the knots, i.e. the points where the analytical expression of the Nurbs curve may change. The knot values are only meaningful upon a multiplicative constant, and they are usually normalized to the range [0.0..1.0].

A Nurbs curve with nctrl points and of given degree needs a knot vector with nknots = nctrl+degree+1 values. A degree curve needs at least nctrl = degree+1 control points, and thus at least nknots = 2*(degree+1) knot values.

By default, the returned knot vector will produce a blended, clamped, open curve. A blended curve uses overlapping basis functions spanning multiple intervals of the knot vector. This gives the highest smoothness. An unblended (or decomposed) curve is like a chain of separate curves, each formed on a sequence of degree+1 control points.

A clamped (open) curve starts and ends exactly at the first and last control points. Unclamped curves are deprecated and therefore genKnotVector only produces knot vectors for clamped curves. To achieve this, the knot vector should have multiplicity degree+1 for its end values.

Thus, for an open blended curve, the policy is to set the knot values at the ends to 0.0, resp. 1.0, both with multiplicity degree+1, and to spread the remaining nctrl-degree-1 values equally over the interval.

For an open unblended curve, all internal knots get multiplicity degree. This results in a curve that is only one time continuously derivable at the knots, thus the curve is smooth between the knots, but can have crusts at the knots. There is an extra requirement in this case: nctrl should be a multiple of degree plus 1.

For a closed curve, we currently only produce blended curves. The knots are equally spread over the interval, all having a multiplicity 1 for maximum continuity of the curve.

Parameters:
  • nctrl (int) – The number of control points. This should be higher than the degree of the curve.

  • degree (int) – The intended degree of the curve. Should be smaller than nctrl.

Returns:

KnotVector – A KnotVector to create a NurbsCurve with the specified degree and number of control points.

Examples

>>> nctrl = 5
>>> for i in range(1, nctrl): print(genKnotVector(nctrl, i))
KnotVector(7): 0*2, 0.25*1, 0.5*1, 0.75*1, 1*2
KnotVector(8): 0*3, 0.333333*1, 0.666667*1, 1*3
KnotVector(9): 0*4, 0.5*1, 1*4
KnotVector(10): 0*5, 1*5
>>> print(genKnotVector(7,3))
KnotVector(11): 0*4, 0.25*1, 0.5*1, 0.75*1, 1*4
>>> print(genKnotVector(7,3,blended=False))
KnotVector(11): 0*4, 1*3, 2*4
>>> print(genKnotVector(3,2,closed=True))
KnotVector(6): 0*1, 0.2*1, 0.4*1, 0.6*1, 0.8*1, 1*1
plugins.nurbs.akimTangents(Q, corner=0.5, normalize=True)[source]

Estimate tangents to a curve passing through given points.

Parameters:
  • Q (coords_like (npts, 3)) – The points through which the curve should pass.

  • corner (float) – A value in the range 0.0..1.0. A value 1.0 will set the tangent at a corner to that after the corner. A value sets it to that before. A value 0.5 gives an intermediate value.

  • normalize (bool) – If True (default), normalized vectors are returned. If False, the vectors are not normalized and are good approximations for the derivative vectors of a curve interpolating/approximating the points.

Returns:

array (npts, 3) – The normalized tangent vectors at the points Q to a curve that interpolates those points.

Notes

Based on The NURBS Book (9.29) and (9.31).

plugins.nurbs.cubicInterpolate(Q, T=None, return_param=False)[source]

Create a C1 cubic interpolation curve

Parameters:
  • Q (coords_like (npts, 3)) – The points through which the curve should pass.

  • T (coords_like (npts, 3), optional) – The normalized tangent vectors at the points Q to a curve that interpolates those points. If not provided, they are computed with akimTangents().

  • return_param (bool) – If True, also returns also parameter values where the given points occur.

Returns:

  • NurbsCurve – A cubic NurbsCurve interpolating the points Q and having tangents T at those points. The curve is guaranteed to have C1 continuity (including in the speed of the curve) and has a fairly uniform parametrization.

  • u (float array (npts, )) – Only returned if return_param=True: the parametric values where the input points are found on the NurbsCurve. Thus, N.pointsAt(u) produces Q.

Notes

Based on The NURBS Book 9.3.4.

plugins.nurbs.quadraticInterpolate(Q, T=None, corners=False, reduce=True, return_param=False)[source]

Create a G1/C1 quadratic interpolation curve

Parameters:
  • Q (coords_like (npts, 3)) – The points through which the curve should pass.

  • T (coords_like (npts, 3), optional) – The normalized tangent vectors at the points Q to a curve that interpolates those points. If not provided, they are computed with akimTangents(), with a corner value of 1.0 if corners is True, and 0.5 if not.

  • corners (bool, optional) – If True, obvious corners will be retained in the result. The default (False) will round corners.

  • reduce (bool, optional) – If True (default), the curve is reduced by removing the internal oncurve control points. If set False, the oncurve points are retained. Note that there may be more oncurve points than the originally specified set.

  • return_param (bool) – If True, also returns also parameter values where the given points occur.

Returns:

NurbsCurve – A quadratic NurbsCurve interpolating the given points and having tangents equal to the provided or computed ones. If not reduced, and corners is False, the curve has G1 (geometric) continuity and a fairly uniform parametrization. If reduced, it has C1 continuity (thus including the curve speed) and a smaller footprint, at the price of a less uniform parametrization. If corners is True, it has C0 continuity.

Notes

Based on The NURBS Book 9.3.2

plugins.nurbs.globalInterpolationParameters(Q, exp=1.0)[source]

Compute parameters for a global interpolation curve

The global interpolation algorithm computes the control points that produce a NurbsCurve with given points occurring at predefined parameter values. The curve shape depends on the choosen values. This function provides a way to set values that work well under mosr conditions/

Parameters:
  • Q (coords_like (npts, 3)) – An ordered set of points through which the curve should pass. Two consecutive points should not coincide.

  • exp (float) – The exponent to be used in the interpolation algorithm. See Notes.

Returns:

u (float array (npts, )) – Only returned if return_param=True: the parametric values where the input points are found on the NurbsCurve. Thus, N.pointsAt(u) produces Q.

Notes

The algorithm to set these values uses a variable exponent. Different values produce (slighly) different curves. The smaller the value, the more the two spans get curved. Values above 1 will almost straighten the end spans, but intermediate spans become more curved. Typical values are:

0.0: equally spaced (not recommended) 0.5: centripetal (recommended when data set take sharp turns) 0.7: our prefered default 1.0: chord length (widely used)

plugins.nurbs.cubicSpline(Q, D0, D1, *, u=None, exp=0.7, return_param=False)[source]

Cubic spline interpolation.

Computes a traditional C2 cubic spline through the points Q with given end tangents.

Parameters:
  • Q (coords_like (npts, 3)) – An ordered set of points through which the curve should pass. Two consecutive points should not coincide.

  • D0 (coords_like (3,)) – The derivative of the curve at the start point Q[0]. The length of the vector is significant.

  • D1 (coords_like (3,)) – The derivative of the curve at the end point Q[-1]. The length of the vector is significant.

  • u (float array (npts,), optional) – The parameter values where the points Q should be obtained on the curve. If not provided, a default set of parameter values is computed from globalInterpolationParameters() (Q, exp=exp)``.

  • exp (float, optional) – The exponent to be used in computing the parameter values if no u was provided. The default value 0.7 will work well in most cases. See globalInterpolationParameters().

  • return_param (bool) – If True, also returns also parameter values where the given points occur.

Returns:

  • N (NurbsCurve) – A NurbsCurve of the third degree that passes through the given point set Q and has tangents D0 and D1 at its ends. The number of control points of the curve is npts + 2. The parametric values of the input points Q can be got from N.knots[3:-3].

  • u (float array (npts, )) – Only returned if return_param=True: the parametric values where the input points are found on the NurbsCurve. Thus, N.pointsAt(u) produces Q.

See also

globalInterpolationCurve

global interpolation curve of any degree

plugins.nurbs.globalInterpolationCurve(Q, degree=3, *, D=None, D0=None, D1=None, u=None, exp=0.7, return_param=False)[source]

Create a global interpolation NurbsCurve.

An interpolation curve is a curve passing through all the given points.

Parameters:
  • Q (coords_like (npts, 3)) – An ordered set of points through which the curve should pass. Two consecutive points should not coincide.

  • degree (int) – The degree of the resulting curve. Usually 2 or 3 is used. For degree 1, the result is a PolyLine through the points Q. For degrees higher than 4, it is better to create a degree 3 curve and then increase the degree.

  • D (coords_like (npts, 3)) – The derivatives of the curve at the points Q. The longer the derivative vectors, the further the curve will follow the tangent. If not provided, the interpolation is done without imposing tangent directions at the intermediate points.

  • D0 (coords_like (3,)) – The derivative of the curve at the start point. Can only be used if D is not provided, and derivatives at one or both ends need to be imposed. The length of the vector is significant.

  • D1 (coords_like (3,)) – The derivative of the curve at the end point. Can only be used if D is not provided, and derivatives at one or both ends need to be imposed. The length of the vector is significant.

  • u (float array_like (npts,), optional) – The parameter values where the points Q should be obtained on the curve. If not provided, a default set of parameter values is computed from globalInterpolationParameters(Q, exp=exp).

  • exp (float, optional) – The exponent to be used in computing the parameter values if no u was provided. The default value 0.7 works well in most cases. See globalInterpolationParameters().

  • return_param (bool) – If True, also returns also parameter values where the given points occur.

Returns:

  • N (NurbsCurve) – A NurbsCurve of the specified degree that passes through the given point set. The number of control points is equal to the number of input points plus one for every end tangent set.

  • u (float array (npts, )) – Only returned if return_param=True: the parametric values where the input points are found on the NurbsCurve. Thus, N.pointsAt(u) produces Q.

See also

cubicSpline

returns the classical C2 spline interpolate

plugins.nurbs.lsqCurve(Q, p, nctrl, Wq=None, D=None, I=None, Wd=None, return_param=False)[source]

Weighed and constrained least squares curve fit

Creates a NurbsCurve of a given degree and with given number of control points, approximating a sequence of points. Points can be attributed weights to force the curve closer to specific points. Tangents can be specified at any point, and can be weighed to have better approximation at specific points. Specifying negative weights for specific points or tangents turns them into constraints to be met accurately.

Parameters:
  • Q (coords_like (npts, 3)) – The points through which to fit a NurbsCurve

  • p (int) – The degree of the target curve. A typical value is 3.

  • nctrl (int) – The number of control points in the target curve. It should be higher than the degree p, but smaller smaller than the number of points npts. The higher the number, the closer the curve can come to all points, at the cost of higher oscillations though. A value of 2 with p=1 fits a straight line throught the data.

  • Wq (array_like float (npts,), optional) – The weigths of the npts points, where the value is non-negative. Where negative, a point constraint is added to force the curve through the point. If not specified, all weights are set equal to 1.

  • D (coords_like (nder, 3), optional) – The derivatives at nder points. Note that the length of the vectors is significant. See Notes.

  • I (array_like int (nder,), optional) – The index of the points in Q for which the derivatives D are specified. Required if D is provided.

  • Wd (array_like float (nder,), optional) – The weight of the nder derivatives, where the value is non-negative. Where negative, a constraint is added to force the curve to have that derivative. If not specified, all weights are set equal to 1.

  • return_param (bool) – If True, also returns also parameter values where the given points occur.

Returns:

  • N (NurbsCurve) – A NurbsCurve of the specified degree and with the specified number of control points, passing as close as possible to the given (weighed) points and obeying any specified constraints.

  • u (float array (npts, )) – Only returned if return_param=True: the parametric values where the input points are found on the NurbsCurve. Thus, N.pointsAt(u) produces points in the neighborhood of Q.

Notes

If derivatives or constraints are provided, the number of control points is further restricted. Let nc be the number of constraints points and derivatives and nu be the number of unconstrained points and derivatives, then nctrl should satisfy: nc-2 < nctrl < nu-nc+2.

The length of the derivatives has a significant influence on the resulting curve. The recommended length is the polygonal length of the point set. This can be achieved with:

import pyformex.arraytools as at
from pyformex.core import PolyLine
D = at.normalize(D) * PolyLine(Q).length()

This is algorithm A9.6 p.417 of the NurbsBook.

This function requires scipy.

plugins.nurbs.optLsqCurve(Q, p, dtol, nmin=None, nmax=None)[source]

Compute an optimal lsqCurve

Finds the lsqCurve with minimum number of control points that passes not further than dtol from any of the given points.

Parameters:
  • Q (see lsqCurve) –

  • p (see lsqCurve) –

  • dtol (float) – Maximum distance from the curve allowed for any of the points Q. Note that specifying a very small value may lead to high oscillations in the curve.

  • nmin (int) – Minimum number of control points in the curve. If not specified, it is set to p+1.

  • nmax (int) – Minimum number of control points in the curve. If not specified, it is set to npts-1.

Returns:

  • N (NurbsCurve) – The NurbsCurve obtained with lsqCurve with the minimum number of control points guaranteeing that no point is further than dtol from the curve.

  • d (float) – The maximum distance of any point Q from the curve N.

plugins.nurbs.NurbsCircle(C=[0.0, 0.0, 0.0], r=1.0, X=[1.0, 0.0, 0.0], Y=[0.0, 1.0, 0.0], ths=0.0, the=360.0)[source]

Create a NurbsCurve representing a perfect circle or arc.

Parameters:
  • C (float (3,)) – The center of the circle

  • r (float) – The radius of the circle

  • X (float (3,)) – A unit vector in the plane of the circle

  • Y (float (3,)) – A unit vector in the plane of the circle and perpendicular to X

  • ths (float) – Start angle (in degrees) of the arc, measured from the X axis, counterclockwise in X-Y plane

  • the (float) – End angle (in degrees) of the arc, measured from the X axis, counterclockwise in X-Y plane

Returns:

NurbsCurve – A NurbsCurve that is a perfect circle or arc in the plane of the XY-axes.

Notes

NurbsBook algorithm A7.1

plugins.nurbs.toCoords4(x)[source]

Convert cartesian coordinates to homogeneous

x: Coords

Array with cartesian coordinates.

Returns a Coords4 object corresponding to the input cartesian coordinates.

plugins.nurbs.pointsOnBezierCurve(P, u)[source]

Compute points on a Bezier curve

Parameters:

P is an array with n+1 points defining a Bezier curve of degree n. u is a vector with nu parameter values between 0 and 1.

Returns:

An array with the nu points of the Bezier curve corresponding with the specified parametric values. ERROR: currently u is a single paramtric value!

See also: examples BezierCurve, Casteljau

plugins.nurbs.frenet(d1, d2, d3=None)[source]

Compute Frenet vectors, curvature and torsion.

This function computes Frenet vectors, curvature and torsion from the provided first, second, and optional third derivatives of curve. The derivatives can be obtained from NurbsCurve.deriv(). Curvature is computed as abs| d1 x d2 | / |d1|**3

Parameters:
  • d1 (float array_like (npts,3)) – First derivative at npts points of a nurbs curve

  • d2 (float array_like (npts,3)) – Second derivative at npts points of a nurbs curve

  • d3 (float array_like (npts,3) , optional) – Third derivative at npts points of a nurbs curve

Returns:

  • T (float array(npts,3)) – Normalized tangent vector to the curve at npts points.

  • N (float array(npts,3)) – Normalized normal vector to the curve at npts points.

  • B (float array(npts,3)) – Normalized binormal vector to the curve at npts points.

  • k (float array(npts,3)) – Curvature of the curve at npts points.

  • t (float array(npts,3), optional) – Torsion of the curve at npts points. This value is only returned if d3 was provided.

See also

NurbsCurve.frenet

the corresponding NurbsCurve method

NurbsCurve.deriv

computation of the derivatives of a NurbsCurve