Basic classes

BSplineBasis

class splipy.BSplineBasis[source]

Represents a one-dimensional B-Spline basis.

BSplineBasis objects support basic arithmetic operators, which are interpreted as acting on the parametric domain.

__getitem__(i)[source]

Returns the knot at a given index.

__init__(order=2, knots=None, periodic=-1)[source]

Construct a B-Spline basis with a given order and knot vector.

Parameters:
  • order (int) – Spline order, i.e. one greater than the polynomial degree.
  • knots ([float]) – Knot vector of non-decreasing components. Defaults to open knot vector on domain [0,1].
  • periodic (int) – Number of continuous derivatives at start and end. –1 is not periodic, 0 is continuous, etc.
Raises:

ValueError – for inapproriate knot vectors

__len__()[source]

Returns the number of knots in this basis.

clone()[source]

Clone the object.

continuity(knot)[source]

Get the continuity of the basis functions at a given point.

Returns:pm–1 at a knot with multiplicity m, or inf between knots.
Return type:int or float
end()[source]

End point of parametric domain. For open knot vectors, this is the last knot.

Returns:Knot number np, where p is the spline order and n is the number of knots
Return type:Float
evaluate(t, d=0, from_right=True, sparse=False)[source]

Evaluate all basis functions in a given set of points. :param t: The parametric coordinate(s) in which to evaluate :type t: float or [float] :param int d: Number of derivatives to compute :param bool from_right: True if evaluation should be done in the limit

from above
Parameters:sparse (bool) – True if computed matrix should be returned as sparse
Returns:A matrix N[i,j] of all basis functions j evaluated in all points i
Return type:numpy.array
greville(index=None)[source]

Fetch greville points, also known as knot averages:

\[\sum_{j=i+1}^{i+p-1} \frac{t_j}{p-1}\]
Returns:One, or all of the Greville points
Return type:[float] (if index is None) or float
insert_knot(new_knot)[source]

Inserts a knot in the knot vector.

The return value is a sparse matrix C (actually, a dense matrix with lots of zeros), such that N_new = N_old x C, where N are row vectors of basis functions.

Parameters:new_knot (float) – The parametric coordinate of the point to insert
Returns:Transformation matrix C
Return type:numpy.array
Raises:ValueError – If the new knot is outside the domain
integrate(t0, t1)[source]

Integrate all basis functions over a given domain

Parameters:
  • t0 (float) – The parametric starting point
  • t1 (float) – The parametric end point
Returns:

The integration of all functions over the input domain

Return type:

list

knot_spans(include_ghost_knots=False)[source]

Return the set of unique knots in the knot vector.

Parameters:include_ghost_knots (bool) – if knots outside start/end are to be included. These knots are used by periodic basis.
Returns:List of unique knots
Return type:[float]
lower_order(amount)[source]

Create a knot vector with lower order.

The continuity at the knots are kept unchanged by decreasing their multiplicities.

Returns:

New knot vector

Return type:

[float]

Raises:
  • TypeError – If amount is not an int
  • ValueError – If amount is negative
make_periodic(continuity)[source]

Create a periodic basis with a given continuity.

matches(bspline, reverse=False)[source]

Checks if this basis equals another basis, when disregarding scaling and translation of the knots vector. I.e. will this basis and bspline yield the same spline object if paired with identical controlpoints

normalize()[source]

Set the parametric domain to be (0,1).

num_functions()[source]

Returns the number of basis functions in the basis.

Warning

This is different from splipy.BSplineBasis.__len__().

raise_order(amount)[source]

Create a knot vector with higher order.

The continuity at the knots are kept unchanged by increasing their multiplicities.

Returns:

New knot vector

Return type:

[float]

Raises:
  • TypeError – If amount is not an int
  • ValueError – If amount is negative
reparam(start=0, end=1)[source]

Set the parametric domain to be (start, end)

Raises:ValueError – If endstart
reverse()[source]

Reverse parametric domain, keeping start/end values unchanged.

roll(new_start)[source]

rotate a periodic knot vector by setting a new starting index.

Parameters:new_start (int) – The index of to the new first knot
snap(t)[source]

Snap evaluation points to knots if they are sufficiently close as given in by state.state.knot_tolerance. This will modify the input vector t

Parameters:t ([float]) – evaluation points
Returns:none
start()[source]

Start point of parametric domain. For open knot vectors, this is the first knot.

Returns:Knot number p, where p is the spline order
Return type:float

SplineObject

class splipy.SplineObject(bases=None, controlpoints=None, rational=False, raw=False)[source]

Master class for spline objects with arbitrary dimensions.

This class should be subclassed instead of used directly.

All SplineObjects support basic arithmetic operators, which are interpreted as translation and scaling. In-place operators (e.g. +=) mutate the object, while infix operators (e.g. +) create new objects.

__init__(bases=None, controlpoints=None, rational=False, raw=False)[source]

Construct a spline object with the given bases and control points.

The default is to create a linear one-element mapping from and to the unit (hyper)cube.

Parameters:
  • bases ([BSplineBasis]) – The basis of each parameter direction
  • controlpoints (array-like) – An n1 × n2 × … × d matrix of control points
  • rational (bool) – Whether the object is rational (in which case the control points are interpreted as pre-multiplied with the weight, which is the last coordinate)
  • raw (bool) – If True, skip any control point reordering. (For internal use.)
bounding_box()[source]

Gets the bounding box of a spline object, computed from the control-point values. Could be inaccurate for rational splines.

Returns the minima and maxima for each direction: [(xmin, xmax), (ymin, ymax), …]

Returns:Bounding box
Return type:[(float)]
center()[source]

Gets the center of the domain

For curves this will return \((\tilde{x}, \tilde{y},...)\), where

\[\tilde{x} = \frac{1}{L} \int_{t_0}^{t_1} x(t) \; dt\]

and \(L=t_1-t_0\) is the length of the parametric domain \([t_0,t_1]\).

For surfaces this will return \((\tilde{x}, \tilde{y},...)\), where

\[\tilde{x} = \frac{1}{A} \int_{v_0}^{v_1} \int_{u_0}^{u_1} x(u,v) \; du \; dv\]

and \(A=(u_1-u_0)(v_1-v_0)\) is the area of the parametric domain \([u_0,u_1]\times[v_0,v_1]\).

Warning

For rational splines, this will integrate in projective coordinates, then project the centerpoint. This is as opposed to integrating the rational functions \(\frac{N_i(t)w_i}{\sum_j N_j(t)w_j}\).

clone()[source]

Clone the object.

corners(order='C')[source]

Return the corner control points.

The order parameter determines which order to use, either 'F' or 'C', for row-major or column-major ordering. E.g. for a volume, in parametric coordinates,

  • 'C' gives (0,0,0), (1,0,0), (0,1,0), (1,1,0), (0,0,1), etc.
  • 'F' gives (0,0,0), (0,0,1), (0,1,0), (0,1,1), (1,0,0), etc.
Parameters:order (str) – The ordering to use
Returns:Corners
Return type:np.array

Warning

For rational splines, this will return the corners in projective coordinates, including weights.

derivative(*params, **kwargs)[source]

Evaluate the derivative of the object at the given parametric values.

If tensor is true, evaluation will take place on a tensor product grid, i.e. it will return an n1 × n2 × … × dim array, where ni is the number of evaluation points in direction i, and dim is the physical dimension of the object.

If tensor is false, there must be an equal number n of evaluation points in all directions, and the return value will be an n × dim array.

If there is only one evaluation point, a vector of length dim is returned instead.

Examples:

# Tangent of curve at single point
curve.derivative(1.0)

# Double derivative of curve at single point:
curve.derivative(1.0, d=2)

# Third derivative of curve at several points:
curve.derivative([0.0, 1.0, 2.0], d=3)

# Tangents of surface:
surface.derivative(0.5, 0.7, d=(1,0))
surface.derivative(0.5, 0.7, d=(0,1))

# Cross-derivative of surface:
surface.derivative(0.5, 0.7, d=(1,1))
Parameters:
  • u,v,.. (float or [float]) – Parametric coordinates in which to evaluate
  • d ((int)) – Order of derivative to compute
  • above ((bool)) – Evaluation in the limit from above
  • tensor (bool) – Whether to evaluate on a tensor product grid
Returns:

Derivatives

Return type:

numpy.array

end(direction=None)[source]

Return the end of the parametric domain.

If direction is given, returns the end of that direction, as a float. If it is not given, returns the end of all directions, as a tuple.

Parameters:direction (int) – Direction in which to get the end.
Raises:ValueError – For invalid direction
evaluate(*params, **kwargs)[source]

Evaluate the object at given parametric values.

If tensor is true, evaluation will take place on a tensor product grid, i.e. it will return an n1 × n2 × … × dim array, where ni is the number of evaluation points in direction i, and dim is the physical dimension of the object.

If tensor is false, there must be an equal number n of evaluation points in all directions, and the return value will be an n × dim array.

If there is only one evaluation point, a vector of length dim is returned instead.

Parameters:
  • u,v,.. (float or [float]) – Parametric coordinates in which to evaluate
  • tensor (bool) – Whether to evaluate on a tensor product grid
Returns:

Geometry coordinates

Return type:

numpy.array

force_rational()[source]

Force a rational representation of the object.

The weights of a non-rational object will be set to 1.

Returns:self
get_derivative_spline(direction=None)[source]

Compute the controlpoints associated with the derivative spline object

If direction is given, only the derivatives in that direction are returned.

If direction is not given, this function returns a tuple of all partial derivatives

# Create a 4x4 element cubic spline surface
surf = Surface()
surf.raise_order(2,2)
surf.refine(3,3)
surf[1:4,1:4,:] += 0.1 # make the surface non-trivial by moving controlpoints

# Create the derivative surface
du = surf.get_derivative_spline(direction='u')

# evaluation is identical
print(du.evaluate(0.3, 0.4))
print(surf.derivative(0.3, 0.4, d=(1,0)))

print(surf.order()) # prints (3,3)
print(du.order())   # prints (2,3)
Parameters:direction (int) – The tangential direction
Returns:Derivative spline
Return type:SplineObject
insert_knot(knot, direction=0)[source]

Insert a new knot into the spline.

Parameters:
  • direction (int) – The direction to insert in
  • knot (float or [float]) – The new knot(s) to insert
Raises:

ValueError – For invalid direction

Returns:

self

knots(direction=None, with_multiplicities=False)[source]

Return knots vector

If direction is given, returns the knots in that direction, as a list. If it is not given, returns the knots of all directions, as a tuple.

Parameters:
  • direction (int) – Direction in which to get the knots.
  • with_multiplicities (bool) – If true, return knots with multiplicities (i.e. repeated).
Raises:

ValueError – For invalid direction

lower_order(*lowers)[source]

Lower the polynomial order of the object. If only one argument is given, the order is lowered equally over all directions.

Parameters:u,v,.. (int) – Number of times to lower the order in a given direction.
Return SplineObject:
 Approximation of the current object on a lower order basis
lower_periodic(periodic, direction=0)[source]

Sets the periodicity of the spline object in the given direction, keeping the geometry unchanged.

Parameters:
  • periodic (int) – new periodicity, i.e. the basis is C^k over the start/end
  • direction (int) – the parametric direction of the basis to modify
Returns:

self

make_periodic(continuity=None, direction=0)[source]

Make the spline object periodic in a given parametric direction.

Parameters:
  • continuity (int) – The continuity along the boundary (default max).
  • direction (int) – The direction to ensure continuity in.
classmethod make_splines_compatible(spline1, spline2)[source]

Ensure that two splines are compatible.

This will manipulate one or both to ensure that they are both rational or nonrational, and that they lie in the same physical space.

Parameters:
classmethod make_splines_identical(spline1, spline2)[source]

Ensure that two splines have identical discretization.

This will first make them compatible (see splipy.SplineObject.make_curves_compatible()), reparametrize them, and possibly raise the order and insert knots as required.

Parameters:
mirror(normal)[source]

Mirror the object around a plane through the origin.

Parameters:normal (array-like) – The plane normal to mirror about.
Raises:RuntimeError – If the physical dimension is not 2 or 3
Returns:self
order(direction=None)[source]

Return polynomial order (degree + 1).

If direction is given, returns the order of that direction, as an int. If it is not given, returns the order of all directions, as a tuple.

Parameters:direction (int) – Direction in which to get the order.
Raises:ValueError – For invalid direction
pardim

The number of parametric dimensions: 1 for curves, 2 for surfaces, 3 for volumes, etc.

periodic(direction=0)[source]

Returns true if the spline object is periodic in the given parametric direction

project(plane)[source]

Projects the geometry onto a plane or axis.

  • project(‘xy’) will project the object onto the xy plane, setting all z components to zero.
  • project(‘y’) will project the object onto the y axis, setting all x and z components to zero.
Parameters:plane (string) – Any combination of ‘x’, ‘y’ and ‘z’
Returns:self
raise_order(*raises)[source]

Raise the polynomial order of the object. If only one argument is given, the order is raised equally over all directions.

Parameters:u,v,.. (int) – Number of times to raise the order in a given direction.
Returns:self
refine(*ns, **kwargs)[source]

Enrich the spline space by inserting knots into each existing knot span.

This method supports three different usage patterns:

# Refine each direction by a factor n
obj.refine(n)

# Refine a single direction by a factor n
obj.refine(n, direction='v')

# Refine all directions by given factors
obj.refine(nu, nv, ...)
Parameters:
  • nu,nv,.. (int) – Number of new knots to insert into each span
  • direction (int) – Direction to refine in
Returns:

self

reparam(*args, **kwargs)[source]

Redefine the parametric domain. This function accepts two calling conventions:

reparametrize(u, v, …) reparametrizes each direction to the domains given by the tuples u, v, etc. It is equivalent to calling reparametrize(u[0], u[1]) on each basis. The default domain for directions not given is (0,1). In particular, if no arguments are given, the new parametric domain will be the unit (hyper)cube.

reparametrize(u, direction=d) reparametrizes just the direction given by d and leaves the others untouched.

Parameters:
  • u, v, .. (tuple) – New parametric domains, default to (0,1)
  • direction (int) – The direction to reparametrize
Returns:

self

reverse(direction=0)[source]

Swap the direction of a parameter by making it go in the reverse direction. The parametric domain remains unchanged.

Parameters:direction (int) – The direction to flip.
Returns:self
rotate(theta, normal=(0, 0, 1))[source]

Rotate the object around an axis.

Parameters:
  • theta (float) – Angle to rotate about, measured in radians
  • normal (array-like) – The normal axis (if 3D) to rotate about
Raises:

RuntimeError – If the physical dimension is not 2 or 3

Returns:

self

scale(*args)[source]

Scale, or magnify the object by a given amount.

In case of one input argument, the scaling is uniform.

Parameters:args (array-like or float) – Scaling factors, possibly different in each direction.
Returns:self
section(*args, **kwargs)[source]

Returns a section from the object. A section can be any sub-object of parametric dimension not exceeding that of the object. E.g. for a volume, sections include vertices, edges, faces, etc.

The arguments are control point indices for each direction. None means that direction should be variable in the returned object.

# Get the face with u=max
vol.section(-1, None, None)

# Keyword arguments are supported for u, v and w
# This is the same as the above
vol.section(u=-1)

# Get the edge with u=min, v=max
vol.section(0, -1, None)

# This is equivalent to vol.clone()
vol.section()

If a specific subclass of SplineObject is found that handles the requested number of variable directions (parametric dimension), then the return value is of that type. If not, it will be a generic SplineObject.

If the section has no variable directions (it is a point), then the return value will be an array, unless the keyword argument unwrap_points is true, in which case it will return a zero-dimensional SplineObject.

Parameters:u,v,.. (int or None) – Control point indices
Returns:Section
Return type:SplineObject or np.array
set_dimension(new_dim)[source]

Sets the physical dimension of the object. If increased, the new components are set to zero.

Parameters:new_dim (int) – New dimension.
Returns:self
set_order(*order)[source]

Set the polynomial order of the object. If only one argument is given, the order is set uniformly over all directions.

Parameters:u,v,.. (int) – The new order in a given direction.
Raises:ValueError – If the order is reduced in any direction.
Returns:self
shape

The dimensions of the control point array.

split(knots, direction=0)[source]

Split an object into two or more separate representations with C0 continuity between them.

Parameters:
  • knots (float or [float]) – The splitting points
  • direction (int) – Parametric direction
Returns:

The new objects

Return type:

[SplineObject]

start(direction=None)[source]

Return the start of the parametric domain.

If direction is given, returns the start of that direction, as a float. If it is not given, returns the start of all directions, as a tuple.

Parameters:direction (int) – Direction in which to get the start.
Raises:ValueError – For invalid direction
swap(dir1=0, dir2=1)[source]

Swaps two parameter directions.

This function silently passes for curves.

Parameters:
  • dir1 (direction) – The first direction (default u)
  • dir2 (direction) – The second direction (default v)
Returns:

self

tangent(*params, **kwargs)[source]

Evaluate the tangents of the object at the given parametric values.

If direction is given, only the derivatives in that direction are evaluated. This is equivalent to calling splipy.SplineObject.derivative() with d=(0,…,0,1,0,…,0), the unit vector corresponding to the given direction.

If direction is not given, this function returns a tuple of all tangential derivatives at the given points.

Parameters:
  • u,v,.. (float or [float]) – Parametric coordinates in which to evaluate
  • direction (int) – The tangential direction
  • above ((bool)) – Evaluation in the limit from above
  • tensor (bool) – Whether to evaluate on a tensor product grid
Returns:

Tangents

Return type:

tuple<numpy.array>

translate(x)[source]

Translate (i.e. move) the object by a given distance.

Parameters:x (array-like) – The vector to translate by.
Returns:self

Curve

class splipy.Curve[source]

Bases: splipy.SplineObject.SplineObject

Represents a curve: an object with a one-dimensional parameter space.

evaluate(u)[source]

Evaluate the curve at the given parametric values.

This function returns an n × dim array, where n is the number of evaluation points, and dim is the physical dimension of the curve.

If there is only one evaluation point, a vector of length dim is returned instead.

Parameters:u (float or [float]) – Parametric coordinates in which to evaluate
Returns:Geometry coordinates
Return type:numpy.array
__getitem__(i)

Get the control point at a given index.

Indexing is in column-major order. Examples of supported indexing modes:

# Flat indexing with an int
obj[4]

# Flat indexing from the end
obj[-1]

# Flat indexing with a slice
obj[2:5]

# Multi-indexing with ints, negative ints and slices
obj[0,-1,:]
Return type:numpy.array
__init__(basis=None, controlpoints=None, rational=False, **kwargs)[source]

Construct a curve with the given basis and control points.

The default is to create a linear one-element mapping from (0,1) to the unit interval.

Parameters:
  • basis (BSplineBasis) – The underlying B-Spline basis
  • controlpoints (array-like) – An n × d matrix of control points
  • rational (bool) – Whether the curve is rational (in which case the control points are interpreted as pre-multiplied with the weight, which is the last coordinate)
__len__()

Return the number of control points (basis functions) for the object.

__setitem__(i, cp)

Set the control points at given indices.

This function supports the same indexing modes as SplineObject.__getitem__()

Parameters:
  • i (int) – Index or indices
  • cp (numpy.array) – New control point(s)
append(curve)[source]

Extend the curve by merging another curve to the end of it.

The curves are glued together in a C0 fashion with enough repeated knots. The function assumes that the end of this curve perfectly matches the start of the input curve.

Parameters:curve (Curve) – Another curve
Raises:RuntimeError – If either curve is periodic
Returns:self
binormal(t, above=True)[source]

Evaluate the normalized binormal of the curve at the given parametric value(s).

This function returns an n × 3 array, where n is the number of evaluation points.

The binormal is computed as the normalized cross product between the velocity and acceleration of the curve.

Parameters:
  • t (float or [float]) – Parametric coordinates in which to evaluate
  • above (bool) – Evaluation in the limit from above
Returns:

Derivative array

Return type:

numpy.array

bounding_box()

Gets the bounding box of a spline object, computed from the control-point values. Could be inaccurate for rational splines.

Returns the minima and maxima for each direction: [(xmin, xmax), (ymin, ymax), …]

Returns:Bounding box
Return type:[(float)]
center()

Gets the center of the domain

For curves this will return \((\tilde{x}, \tilde{y},...)\), where

\[\tilde{x} = \frac{1}{L} \int_{t_0}^{t_1} x(t) \; dt\]

and \(L=t_1-t_0\) is the length of the parametric domain \([t_0,t_1]\).

For surfaces this will return \((\tilde{x}, \tilde{y},...)\), where

\[\tilde{x} = \frac{1}{A} \int_{v_0}^{v_1} \int_{u_0}^{u_1} x(u,v) \; du \; dv\]

and \(A=(u_1-u_0)(v_1-v_0)\) is the area of the parametric domain \([u_0,u_1]\times[v_0,v_1]\).

Warning

For rational splines, this will integrate in projective coordinates, then project the centerpoint. This is as opposed to integrating the rational functions \(\frac{N_i(t)w_i}{\sum_j N_j(t)w_j}\).

clone()

Clone the object.

continuity(knot)[source]

Get the parametric continuity of the curve at a given point. Will return p-1-m, where m is the knot multiplicity and inf between knots

corners(order='C')

Return the corner control points.

The order parameter determines which order to use, either 'F' or 'C', for row-major or column-major ordering. E.g. for a volume, in parametric coordinates,

  • 'C' gives (0,0,0), (1,0,0), (0,1,0), (1,1,0), (0,0,1), etc.
  • 'F' gives (0,0,0), (0,0,1), (0,1,0), (0,1,1), (1,0,0), etc.
Parameters:order (str) – The ordering to use
Returns:Corners
Return type:np.array

Warning

For rational splines, this will return the corners in projective coordinates, including weights.

curvature(t, above=True)[source]

Evaluate the curvaure at specified point(s). The curvature is defined as

\[\frac{|\boldsymbol{v}\times \boldsymbol{a}|}{|\boldsymbol{v}|^3}\]
Parameters:
  • t (float or [float]) – Parametric coordinates in which to evaluate
  • above (bool) – Evaluation in the limit from above
Returns:

Derivative array

Return type:

numpy.array

derivative(t, d=1, above=True, tensor=True)[source]

Evaluate the derivative of the curve at the given parametric values.

This function returns an n × dim array, where n is the number of evaluation points, and dim is the physical dimension of the curve.

If there is only one evaluation point, a vector of length dim is returned instead.

Parameters:
  • t (float or [float]) – Parametric coordinates in which to evaluate
  • d (int) – Number of derivatives to compute
  • above (bool) – Evaluation in the limit from above
  • tensor (bool) – Not used in this method
Returns:

Derivative array

Return type:

numpy.array

end(direction=None)

Return the end of the parametric domain.

If direction is given, returns the end of that direction, as a float. If it is not given, returns the end of all directions, as a tuple.

Parameters:direction (int) – Direction in which to get the end.
Raises:ValueError – For invalid direction
error(target)[source]

Computes the L2 (squared and per knot span) and max error between this curve and a target curve

\[||\boldsymbol{x_h}(t)-\boldsymbol{x}(t)||_{L^2(t_1,t_2)}^2 = \int_{t_1}^{t_2} |\boldsymbol{x_h}(t)-\boldsymbol{x}(t)|^2 dt, \quad \forall \;\text{knots}\;t_1 < t_2\]
\[||\boldsymbol{x_h}(t)-\boldsymbol{x}(t)||_{L^\infty} = \max_t |\boldsymbol{x_h}(t)-\boldsymbol{x}(t)|\]
Parameters:target (function) – callable function which takes as input a vector of evaluation points t and gives as output a matrix x where x[i,j] is component j evaluated at point t[i]
Returns:L2 error per knot-span and the maximum error
Return type:tuple(list(float), float)

Examples:

import numpy as np

def arclength_circle(t):
    return np.array( [np.cos(t), np.sin(t)] ).T

crv = curve_factory.circle(r=1)
(err2, maxerr) = crv.error(arclength_circle)
print '|| e ||_L2  = ', np.sqrt(np.sum(err2))
print '|| e ||_max = ', maxerr
force_rational()

Force a rational representation of the object.

The weights of a non-rational object will be set to 1.

Returns:self
get_derivative_spline(direction=None)

Compute the controlpoints associated with the derivative spline object

If direction is given, only the derivatives in that direction are returned.

If direction is not given, this function returns a tuple of all partial derivatives

# Create a 4x4 element cubic spline surface
surf = Surface()
surf.raise_order(2,2)
surf.refine(3,3)
surf[1:4,1:4,:] += 0.1 # make the surface non-trivial by moving controlpoints

# Create the derivative surface
du = surf.get_derivative_spline(direction='u')

# evaluation is identical
print(du.evaluate(0.3, 0.4))
print(surf.derivative(0.3, 0.4, d=(1,0)))

print(surf.order()) # prints (3,3)
print(du.order())   # prints (2,3)
Parameters:direction (int) – The tangential direction
Returns:Derivative spline
Return type:SplineObject
get_kinks()[source]

Get the parametric coordinates at all points which have C0- continuity

insert_knot(knot, direction=0)

Insert a new knot into the spline.

Parameters:
  • direction (int) – The direction to insert in
  • knot (float or [float]) – The new knot(s) to insert
Raises:

ValueError – For invalid direction

Returns:

self

knots(direction=None, with_multiplicities=False)

Return knots vector

If direction is given, returns the knots in that direction, as a list. If it is not given, returns the knots of all directions, as a tuple.

Parameters:
  • direction (int) – Direction in which to get the knots.
  • with_multiplicities (bool) – If true, return knots with multiplicities (i.e. repeated).
Raises:

ValueError – For invalid direction

length(t0=None, t1=None)[source]

Computes the euclidian length of the curve in geometric space

\[\int_{t_0}^{t_1}\sqrt{x(t)^2 + y(t)^2 + z(t)^2} dt\]
lower_order(*lowers)

Lower the polynomial order of the object. If only one argument is given, the order is lowered equally over all directions.

Parameters:u,v,.. (int) – Number of times to lower the order in a given direction.
Return SplineObject:
 Approximation of the current object on a lower order basis
lower_periodic(periodic, direction=0)

Sets the periodicity of the spline object in the given direction, keeping the geometry unchanged.

Parameters:
  • periodic (int) – new periodicity, i.e. the basis is C^k over the start/end
  • direction (int) – the parametric direction of the basis to modify
Returns:

self

make_periodic(continuity=None, direction=0)

Make the spline object periodic in a given parametric direction.

Parameters:
  • continuity (int) – The continuity along the boundary (default max).
  • direction (int) – The direction to ensure continuity in.
classmethod make_splines_compatible(spline1, spline2)

Ensure that two splines are compatible.

This will manipulate one or both to ensure that they are both rational or nonrational, and that they lie in the same physical space.

Parameters:
classmethod make_splines_identical(spline1, spline2)

Ensure that two splines have identical discretization.

This will first make them compatible (see splipy.SplineObject.make_curves_compatible()), reparametrize them, and possibly raise the order and insert knots as required.

Parameters:
mirror(normal)

Mirror the object around a plane through the origin.

Parameters:normal (array-like) – The plane normal to mirror about.
Raises:RuntimeError – If the physical dimension is not 2 or 3
Returns:self
normal(t, above=True)[source]

Evaluate the normal of the curve at the given parametric value(s).

This function returns an n × 3 array, where n is the number of evaluation points.

The normal is computed as the cross product between the binormal and the tangent of the curve.

Parameters:
  • t (float or [float]) – Parametric coordinates in which to evaluate
  • above (bool) – Evaluation in the limit from above
Returns:

Derivative array

Return type:

numpy.array

order(direction=None)

Return polynomial order (degree + 1).

If direction is given, returns the order of that direction, as an int. If it is not given, returns the order of all directions, as a tuple.

Parameters:direction (int) – Direction in which to get the order.
Raises:ValueError – For invalid direction
pardim

The number of parametric dimensions: 1 for curves, 2 for surfaces, 3 for volumes, etc.

periodic(direction=0)

Returns true if the spline object is periodic in the given parametric direction

project(plane)

Projects the geometry onto a plane or axis.

  • project(‘xy’) will project the object onto the xy plane, setting all z components to zero.
  • project(‘y’) will project the object onto the y axis, setting all x and z components to zero.
Parameters:plane (string) – Any combination of ‘x’, ‘y’ and ‘z’
Returns:self
raise_order(amount)[source]

Raise the polynomial order of the curve.

Parameters:amount (int) – Number of times to raise the order
Returns:self
rebuild(p, n)[source]

Creates an approximation to this curve by resampling it using a uniform knot vector of order p with n control points.

Parameters:
  • p (int) – Polynomial discretization order
  • n (int) – Number of control points
Returns:

A new approximate curve

Return type:

Curve

refine(*ns, **kwargs)

Enrich the spline space by inserting knots into each existing knot span.

This method supports three different usage patterns:

# Refine each direction by a factor n
obj.refine(n)

# Refine a single direction by a factor n
obj.refine(n, direction='v')

# Refine all directions by given factors
obj.refine(nu, nv, ...)
Parameters:
  • nu,nv,.. (int) – Number of new knots to insert into each span
  • direction (int) – Direction to refine in
Returns:

self

reparam(*args, **kwargs)

Redefine the parametric domain. This function accepts two calling conventions:

reparametrize(u, v, …) reparametrizes each direction to the domains given by the tuples u, v, etc. It is equivalent to calling reparametrize(u[0], u[1]) on each basis. The default domain for directions not given is (0,1). In particular, if no arguments are given, the new parametric domain will be the unit (hyper)cube.

reparametrize(u, direction=d) reparametrizes just the direction given by d and leaves the others untouched.

Parameters:
  • u, v, .. (tuple) – New parametric domains, default to (0,1)
  • direction (int) – The direction to reparametrize
Returns:

self

reverse(direction=0)

Swap the direction of a parameter by making it go in the reverse direction. The parametric domain remains unchanged.

Parameters:direction (int) – The direction to flip.
Returns:self
rotate(theta, normal=(0, 0, 1))

Rotate the object around an axis.

Parameters:
  • theta (float) – Angle to rotate about, measured in radians
  • normal (array-like) – The normal axis (if 3D) to rotate about
Raises:

RuntimeError – If the physical dimension is not 2 or 3

Returns:

self

scale(*args)

Scale, or magnify the object by a given amount.

In case of one input argument, the scaling is uniform.

Parameters:args (array-like or float) – Scaling factors, possibly different in each direction.
Returns:self
section(*args, **kwargs)

Returns a section from the object. A section can be any sub-object of parametric dimension not exceeding that of the object. E.g. for a volume, sections include vertices, edges, faces, etc.

The arguments are control point indices for each direction. None means that direction should be variable in the returned object.

# Get the face with u=max
vol.section(-1, None, None)

# Keyword arguments are supported for u, v and w
# This is the same as the above
vol.section(u=-1)

# Get the edge with u=min, v=max
vol.section(0, -1, None)

# This is equivalent to vol.clone()
vol.section()

If a specific subclass of SplineObject is found that handles the requested number of variable directions (parametric dimension), then the return value is of that type. If not, it will be a generic SplineObject.

If the section has no variable directions (it is a point), then the return value will be an array, unless the keyword argument unwrap_points is true, in which case it will return a zero-dimensional SplineObject.

Parameters:u,v,.. (int or None) – Control point indices
Returns:Section
Return type:SplineObject or np.array
set_dimension(new_dim)

Sets the physical dimension of the object. If increased, the new components are set to zero.

Parameters:new_dim (int) – New dimension.
Returns:self
set_order(*order)

Set the polynomial order of the object. If only one argument is given, the order is set uniformly over all directions.

Parameters:u,v,.. (int) – The new order in a given direction.
Raises:ValueError – If the order is reduced in any direction.
Returns:self
shape

The dimensions of the control point array.

split(knots, direction=0)

Split an object into two or more separate representations with C0 continuity between them.

Parameters:
  • knots (float or [float]) – The splitting points
  • direction (int) – Parametric direction
Returns:

The new objects

Return type:

[SplineObject]

start(direction=None)

Return the start of the parametric domain.

If direction is given, returns the start of that direction, as a float. If it is not given, returns the start of all directions, as a tuple.

Parameters:direction (int) – Direction in which to get the start.
Raises:ValueError – For invalid direction
swap(dir1=0, dir2=1)

Swaps two parameter directions.

This function silently passes for curves.

Parameters:
  • dir1 (direction) – The first direction (default u)
  • dir2 (direction) – The second direction (default v)
Returns:

self

tangent(*params, **kwargs)

Evaluate the tangents of the object at the given parametric values.

If direction is given, only the derivatives in that direction are evaluated. This is equivalent to calling splipy.SplineObject.derivative() with d=(0,…,0,1,0,…,0), the unit vector corresponding to the given direction.

If direction is not given, this function returns a tuple of all tangential derivatives at the given points.

Parameters:
  • u,v,.. (float or [float]) – Parametric coordinates in which to evaluate
  • direction (int) – The tangential direction
  • above ((bool)) – Evaluation in the limit from above
  • tensor (bool) – Whether to evaluate on a tensor product grid
Returns:

Tangents

Return type:

tuple<numpy.array>

torsion(t, above=True)[source]

Evaluate the torsion for a 3D curve at specified point(s). The torsion is defined as

\[\frac{(\boldsymbol{v}\times \boldsymbol{a})\cdot (d\boldsymbol{a}/dt)}{|\boldsymbol{v}\times \boldsymbol{a}|^2}\]
Parameters:
  • t (float or [float]) – Parametric coordinates in which to evaluate
  • above (bool) – Evaluation in the limit from above
Returns:

Derivative array

Return type:

numpy.array

translate(x)

Translate (i.e. move) the object by a given distance.

Parameters:x (array-like) – The vector to translate by.
Returns:self

Surface

class splipy.Surface[source]

Bases: splipy.SplineObject.SplineObject

Represents a surface: an object with a two-dimensional parameter space.

evaluate(u, v)

Evaluate the surface at the given parametric values.

This function returns an n1 × n2 × dim array, where ni is the number of evaluation points in each direction, and dim is the dimension of the surface.

If there is only one evaluation point, a vector of length dim is returned instead.

Parameters:
  • u (float or [float]) – Parametric coordinates in the first direction
  • u – Parametric coordinates in the second direction
Returns:

Geometry coordinates

Return type:

numpy.array

evalute_derivative(u, v[, d=(1, 1)])

Evaluate the derivative of the surface at the given parametric values.

This function returns an n1 × n2 × dim array, where ni is the number of evaluation points in direction i, and dim is the dimension of the surface.

If there is only one evaluation point, a vector of length dim is returned instead.

Parameters:
  • u (float or [float]) – Parametric coordinates in the first direction
  • u – Parametric coordinates in the second direction
  • d ((int)) – Order of derivative to compute
Returns:

Derivatives

Return type:

numpy.array

__getitem__(i)

Get the control point at a given index.

Indexing is in column-major order. Examples of supported indexing modes:

# Flat indexing with an int
obj[4]

# Flat indexing from the end
obj[-1]

# Flat indexing with a slice
obj[2:5]

# Multi-indexing with ints, negative ints and slices
obj[0,-1,:]
Return type:numpy.array
__init__(basis1=None, basis2=None, controlpoints=None, rational=False, **kwargs)[source]

Construct a surface with the given basis and control points.

The default is to create a linear one-element mapping from and to the unit square.

Parameters:
  • basis1 (BSplineBasis) – The basis of the first parameter direction
  • basis2 (BSplineBasis) – The basis of the second parameter direction
  • controlpoints (array-like) – An n1 × n2 × d matrix of control points
  • rational (bool) – Whether the surface is rational (in which case the control points are interpreted as pre-multiplied with the weight, which is the last coordinate)
__len__()

Return the number of control points (basis functions) for the object.

__setitem__(i, cp)

Set the control points at given indices.

This function supports the same indexing modes as SplineObject.__getitem__()

Parameters:
  • i (int) – Index or indices
  • cp (numpy.array) – New control point(s)
area()[source]

Computes the area of the surface in geometric space

bounding_box()

Gets the bounding box of a spline object, computed from the control-point values. Could be inaccurate for rational splines.

Returns the minima and maxima for each direction: [(xmin, xmax), (ymin, ymax), …]

Returns:Bounding box
Return type:[(float)]
center()

Gets the center of the domain

For curves this will return \((\tilde{x}, \tilde{y},...)\), where

\[\tilde{x} = \frac{1}{L} \int_{t_0}^{t_1} x(t) \; dt\]

and \(L=t_1-t_0\) is the length of the parametric domain \([t_0,t_1]\).

For surfaces this will return \((\tilde{x}, \tilde{y},...)\), where

\[\tilde{x} = \frac{1}{A} \int_{v_0}^{v_1} \int_{u_0}^{u_1} x(u,v) \; du \; dv\]

and \(A=(u_1-u_0)(v_1-v_0)\) is the area of the parametric domain \([u_0,u_1]\times[v_0,v_1]\).

Warning

For rational splines, this will integrate in projective coordinates, then project the centerpoint. This is as opposed to integrating the rational functions \(\frac{N_i(t)w_i}{\sum_j N_j(t)w_j}\).

clone()

Clone the object.

const_par_curve(knot, direction)[source]

Get a Curve representation of the parametric line of some constant knot value. :param float knot: The constant knot value to sample the surface :param int direction: The parametric direction for the constant value :return: curve on this surface :rtype: Curve

corners(order='C')

Return the corner control points.

The order parameter determines which order to use, either 'F' or 'C', for row-major or column-major ordering. E.g. for a volume, in parametric coordinates,

  • 'C' gives (0,0,0), (1,0,0), (0,1,0), (1,1,0), (0,0,1), etc.
  • 'F' gives (0,0,0), (0,0,1), (0,1,0), (0,1,1), (1,0,0), etc.
Parameters:order (str) – The ordering to use
Returns:Corners
Return type:np.array

Warning

For rational splines, this will return the corners in projective coordinates, including weights.

derivative(*params, **kwargs)

Evaluate the derivative of the object at the given parametric values.

If tensor is true, evaluation will take place on a tensor product grid, i.e. it will return an n1 × n2 × … × dim array, where ni is the number of evaluation points in direction i, and dim is the physical dimension of the object.

If tensor is false, there must be an equal number n of evaluation points in all directions, and the return value will be an n × dim array.

If there is only one evaluation point, a vector of length dim is returned instead.

Examples:

# Tangent of curve at single point
curve.derivative(1.0)

# Double derivative of curve at single point:
curve.derivative(1.0, d=2)

# Third derivative of curve at several points:
curve.derivative([0.0, 1.0, 2.0], d=3)

# Tangents of surface:
surface.derivative(0.5, 0.7, d=(1,0))
surface.derivative(0.5, 0.7, d=(0,1))

# Cross-derivative of surface:
surface.derivative(0.5, 0.7, d=(1,1))
Parameters:
  • u,v,.. (float or [float]) – Parametric coordinates in which to evaluate
  • d ((int)) – Order of derivative to compute
  • above ((bool)) – Evaluation in the limit from above
  • tensor (bool) – Whether to evaluate on a tensor product grid
Returns:

Derivatives

Return type:

numpy.array

edges()[source]

Return the four edge curves in (parametric) order: umin, umax, vmin, vmax

Returns:Edge curves
Return type:(Curve)
end(direction=None)

Return the end of the parametric domain.

If direction is given, returns the end of that direction, as a float. If it is not given, returns the end of all directions, as a tuple.

Parameters:direction (int) – Direction in which to get the end.
Raises:ValueError – For invalid direction
force_rational()

Force a rational representation of the object.

The weights of a non-rational object will be set to 1.

Returns:self
get_derivative_spline(direction=None)

Compute the controlpoints associated with the derivative spline object

If direction is given, only the derivatives in that direction are returned.

If direction is not given, this function returns a tuple of all partial derivatives

# Create a 4x4 element cubic spline surface
surf = Surface()
surf.raise_order(2,2)
surf.refine(3,3)
surf[1:4,1:4,:] += 0.1 # make the surface non-trivial by moving controlpoints

# Create the derivative surface
du = surf.get_derivative_spline(direction='u')

# evaluation is identical
print(du.evaluate(0.3, 0.4))
print(surf.derivative(0.3, 0.4, d=(1,0)))

print(surf.order()) # prints (3,3)
print(du.order())   # prints (2,3)
Parameters:direction (int) – The tangential direction
Returns:Derivative spline
Return type:SplineObject
get_derivative_surface(direction=None)

Compute the controlpoints associated with the derivative spline object

If direction is given, only the derivatives in that direction are returned.

If direction is not given, this function returns a tuple of all partial derivatives

# Create a 4x4 element cubic spline surface
surf = Surface()
surf.raise_order(2,2)
surf.refine(3,3)
surf[1:4,1:4,:] += 0.1 # make the surface non-trivial by moving controlpoints

# Create the derivative surface
du = surf.get_derivative_spline(direction='u')

# evaluation is identical
print(du.evaluate(0.3, 0.4))
print(surf.derivative(0.3, 0.4, d=(1,0)))

print(surf.order()) # prints (3,3)
print(du.order())   # prints (2,3)
Parameters:direction (int) – The tangential direction
Returns:Derivative spline
Return type:SplineObject
insert_knot(knot, direction=0)

Insert a new knot into the spline.

Parameters:
  • direction (int) – The direction to insert in
  • knot (float or [float]) – The new knot(s) to insert
Raises:

ValueError – For invalid direction

Returns:

self

knots(direction=None, with_multiplicities=False)

Return knots vector

If direction is given, returns the knots in that direction, as a list. If it is not given, returns the knots of all directions, as a tuple.

Parameters:
  • direction (int) – Direction in which to get the knots.
  • with_multiplicities (bool) – If true, return knots with multiplicities (i.e. repeated).
Raises:

ValueError – For invalid direction

lower_order(*lowers)

Lower the polynomial order of the object. If only one argument is given, the order is lowered equally over all directions.

Parameters:u,v,.. (int) – Number of times to lower the order in a given direction.
Return SplineObject:
 Approximation of the current object on a lower order basis
lower_periodic(periodic, direction=0)

Sets the periodicity of the spline object in the given direction, keeping the geometry unchanged.

Parameters:
  • periodic (int) – new periodicity, i.e. the basis is C^k over the start/end
  • direction (int) – the parametric direction of the basis to modify
Returns:

self

make_periodic(continuity=None, direction=0)

Make the spline object periodic in a given parametric direction.

Parameters:
  • continuity (int) – The continuity along the boundary (default max).
  • direction (int) – The direction to ensure continuity in.
classmethod make_splines_compatible(spline1, spline2)

Ensure that two splines are compatible.

This will manipulate one or both to ensure that they are both rational or nonrational, and that they lie in the same physical space.

Parameters:
classmethod make_splines_identical(spline1, spline2)

Ensure that two splines have identical discretization.

This will first make them compatible (see splipy.SplineObject.make_curves_compatible()), reparametrize them, and possibly raise the order and insert knots as required.

Parameters:
mirror(normal)

Mirror the object around a plane through the origin.

Parameters:normal (array-like) – The plane normal to mirror about.
Raises:RuntimeError – If the physical dimension is not 2 or 3
Returns:self
normal(u, v, above=(True, True), tensor=True)[source]

Evaluate the normal of the surface at given parametric values.

This is equal to the cross-product between tangents. The return value is normalized.

If tensor is true, evaluation will take place on a tensor product grid, i.e. it will return an n1 × n2 × … × dim array, where ni is the number of evaluation points in direction i, and dim is the physical dimension of the object.

If tensor is false, there must be an equal number n of evaluation points in all directions, and the return value will be an n × dim array.

Parameters:
  • u (float or [float]) – Parametric coordinate(s) in the first direction
  • v (float or [float]) – Parametric coordinate(s) in the second direction
  • above ((bool)) – Evaluation in the limit from above
  • tensor (bool) – Whether to evaluate on a tensor product grid
Returns:

Normal array X[i,j,k] of component xk evaluated at (u[i], v[j])

Return type:

numpy.array

Raises:

RuntimeError – If the physical dimension is not 2 or 3

order(direction=None)

Return polynomial order (degree + 1).

If direction is given, returns the order of that direction, as an int. If it is not given, returns the order of all directions, as a tuple.

Parameters:direction (int) – Direction in which to get the order.
Raises:ValueError – For invalid direction
pardim

The number of parametric dimensions: 1 for curves, 2 for surfaces, 3 for volumes, etc.

periodic(direction=0)

Returns true if the spline object is periodic in the given parametric direction

project(plane)

Projects the geometry onto a plane or axis.

  • project(‘xy’) will project the object onto the xy plane, setting all z components to zero.
  • project(‘y’) will project the object onto the y axis, setting all x and z components to zero.
Parameters:plane (string) – Any combination of ‘x’, ‘y’ and ‘z’
Returns:self
raise_order(*raises)

Raise the polynomial order of the object. If only one argument is given, the order is raised equally over all directions.

Parameters:u,v,.. (int) – Number of times to raise the order in a given direction.
Returns:self
rebuild(p, n)[source]

Creates an approximation to this surface by resampling it using uniform knot vectors of order p with n control points.

Parameters:
  • p ((int)) – Tuple of polynomial discretization order in each direction
  • n ((int)) – Tuple of number of control points in each direction
Returns:

A new approximate surface

Return type:

Surface

refine(*ns, **kwargs)

Enrich the spline space by inserting knots into each existing knot span.

This method supports three different usage patterns:

# Refine each direction by a factor n
obj.refine(n)

# Refine a single direction by a factor n
obj.refine(n, direction='v')

# Refine all directions by given factors
obj.refine(nu, nv, ...)
Parameters:
  • nu,nv,.. (int) – Number of new knots to insert into each span
  • direction (int) – Direction to refine in
Returns:

self

reparam(*args, **kwargs)

Redefine the parametric domain. This function accepts two calling conventions:

reparametrize(u, v, …) reparametrizes each direction to the domains given by the tuples u, v, etc. It is equivalent to calling reparametrize(u[0], u[1]) on each basis. The default domain for directions not given is (0,1). In particular, if no arguments are given, the new parametric domain will be the unit (hyper)cube.

reparametrize(u, direction=d) reparametrizes just the direction given by d and leaves the others untouched.

Parameters:
  • u, v, .. (tuple) – New parametric domains, default to (0,1)
  • direction (int) – The direction to reparametrize
Returns:

self

reverse(direction=0)

Swap the direction of a parameter by making it go in the reverse direction. The parametric domain remains unchanged.

Parameters:direction (int) – The direction to flip.
Returns:self
rotate(theta, normal=(0, 0, 1))

Rotate the object around an axis.

Parameters:
  • theta (float) – Angle to rotate about, measured in radians
  • normal (array-like) – The normal axis (if 3D) to rotate about
Raises:

RuntimeError – If the physical dimension is not 2 or 3

Returns:

self

scale(*args)

Scale, or magnify the object by a given amount.

In case of one input argument, the scaling is uniform.

Parameters:args (array-like or float) – Scaling factors, possibly different in each direction.
Returns:self
section(*args, **kwargs)

Returns a section from the object. A section can be any sub-object of parametric dimension not exceeding that of the object. E.g. for a volume, sections include vertices, edges, faces, etc.

The arguments are control point indices for each direction. None means that direction should be variable in the returned object.

# Get the face with u=max
vol.section(-1, None, None)

# Keyword arguments are supported for u, v and w
# This is the same as the above
vol.section(u=-1)

# Get the edge with u=min, v=max
vol.section(0, -1, None)

# This is equivalent to vol.clone()
vol.section()

If a specific subclass of SplineObject is found that handles the requested number of variable directions (parametric dimension), then the return value is of that type. If not, it will be a generic SplineObject.

If the section has no variable directions (it is a point), then the return value will be an array, unless the keyword argument unwrap_points is true, in which case it will return a zero-dimensional SplineObject.

Parameters:u,v,.. (int or None) – Control point indices
Returns:Section
Return type:SplineObject or np.array
set_dimension(new_dim)

Sets the physical dimension of the object. If increased, the new components are set to zero.

Parameters:new_dim (int) – New dimension.
Returns:self
set_order(*order)

Set the polynomial order of the object. If only one argument is given, the order is set uniformly over all directions.

Parameters:u,v,.. (int) – The new order in a given direction.
Raises:ValueError – If the order is reduced in any direction.
Returns:self
shape

The dimensions of the control point array.

split(knots, direction=0)

Split an object into two or more separate representations with C0 continuity between them.

Parameters:
  • knots (float or [float]) – The splitting points
  • direction (int) – Parametric direction
Returns:

The new objects

Return type:

[SplineObject]

start(direction=None)

Return the start of the parametric domain.

If direction is given, returns the start of that direction, as a float. If it is not given, returns the start of all directions, as a tuple.

Parameters:direction (int) – Direction in which to get the start.
Raises:ValueError – For invalid direction
swap(dir1=0, dir2=1)

Swaps two parameter directions.

This function silently passes for curves.

Parameters:
  • dir1 (direction) – The first direction (default u)
  • dir2 (direction) – The second direction (default v)
Returns:

self

tangent(*params, **kwargs)

Evaluate the tangents of the object at the given parametric values.

If direction is given, only the derivatives in that direction are evaluated. This is equivalent to calling splipy.SplineObject.derivative() with d=(0,…,0,1,0,…,0), the unit vector corresponding to the given direction.

If direction is not given, this function returns a tuple of all tangential derivatives at the given points.

Parameters:
  • u,v,.. (float or [float]) – Parametric coordinates in which to evaluate
  • direction (int) – The tangential direction
  • above ((bool)) – Evaluation in the limit from above
  • tensor (bool) – Whether to evaluate on a tensor product grid
Returns:

Tangents

Return type:

tuple<numpy.array>

translate(x)

Translate (i.e. move) the object by a given distance.

Parameters:x (array-like) – The vector to translate by.
Returns:self

Volume

class splipy.Volume[source]

Bases: splipy.SplineObject.SplineObject

Represents a volume: an object with a three-dimensional parameter space.

evaluate(u, v, w)

Evaluate the volume at the given parametric values.

This function returns an n1 × n2 × n3 × dim array, where ni is the number of evaluation points in each direction, and dim is the dimension of the volume.

If there is only one evaluation point, a vector of length dim is returned instead.

Parameters:
  • u (float or [float]) – Parametric coordinates in the first direction
  • v (float or [float]) – Parametric coordinates in the second direction
  • w (float or [float]) – Parametric coordinates in the third direction
Returns:

Geometry coordinates

Return type:

numpy.array

evalute_derivative(u, v, w[, d=(1, 1, 1)])

Evaluate the derivative of the volume at the given parametric values.

This function returns an n1 × n2 × n3 × dim array, where ni is the number of evaluation points in direction i, and dim is the dimension of the volume.

If there is only one evaluation point, a vector of length dim is returned instead.

Parameters:
  • u (float or [float]) – Parametric coordinates in the first direction
  • v (float or [float]) – Parametric coordinates in the second direction
  • w (float or [float]) – Parametric coordinates in the third direction
  • d ((int)) – Order of derivative to compute
Returns:

Derivatives

Return type:

numpy.array

__getitem__(i)

Get the control point at a given index.

Indexing is in column-major order. Examples of supported indexing modes:

# Flat indexing with an int
obj[4]

# Flat indexing from the end
obj[-1]

# Flat indexing with a slice
obj[2:5]

# Multi-indexing with ints, negative ints and slices
obj[0,-1,:]
Return type:numpy.array
__init__(basis1=None, basis2=None, basis3=None, controlpoints=None, rational=False, **kwargs)[source]

Construct a volume with the given basis and control points.

The default is to create a linear one-element mapping from and to the unit cube.

Parameters:
  • basis1 (BSplineBasis) – The basis of the first parameter direction
  • basis2 (BSplineBasis) – The basis of the second parameter direction
  • basis3 (BSplineBasis) – The basis of the third parameter direction
  • controlpoints (array-like) – An n1 × n2 × n3 × d matrix of control points
  • rational (bool) – Whether the volume is rational (in which case the control points are interpreted as pre-multiplied with the weight, which is the last coordinate)
__len__()

Return the number of control points (basis functions) for the object.

__setitem__(i, cp)

Set the control points at given indices.

This function supports the same indexing modes as SplineObject.__getitem__()

Parameters:
  • i (int) – Index or indices
  • cp (numpy.array) – New control point(s)
bounding_box()

Gets the bounding box of a spline object, computed from the control-point values. Could be inaccurate for rational splines.

Returns the minima and maxima for each direction: [(xmin, xmax), (ymin, ymax), …]

Returns:Bounding box
Return type:[(float)]
center()

Gets the center of the domain

For curves this will return \((\tilde{x}, \tilde{y},...)\), where

\[\tilde{x} = \frac{1}{L} \int_{t_0}^{t_1} x(t) \; dt\]

and \(L=t_1-t_0\) is the length of the parametric domain \([t_0,t_1]\).

For surfaces this will return \((\tilde{x}, \tilde{y},...)\), where

\[\tilde{x} = \frac{1}{A} \int_{v_0}^{v_1} \int_{u_0}^{u_1} x(u,v) \; du \; dv\]

and \(A=(u_1-u_0)(v_1-v_0)\) is the area of the parametric domain \([u_0,u_1]\times[v_0,v_1]\).

Warning

For rational splines, this will integrate in projective coordinates, then project the centerpoint. This is as opposed to integrating the rational functions \(\frac{N_i(t)w_i}{\sum_j N_j(t)w_j}\).

clone()

Clone the object.

corners(order='C')

Return the corner control points.

The order parameter determines which order to use, either 'F' or 'C', for row-major or column-major ordering. E.g. for a volume, in parametric coordinates,

  • 'C' gives (0,0,0), (1,0,0), (0,1,0), (1,1,0), (0,0,1), etc.
  • 'F' gives (0,0,0), (0,0,1), (0,1,0), (0,1,1), (1,0,0), etc.
Parameters:order (str) – The ordering to use
Returns:Corners
Return type:np.array

Warning

For rational splines, this will return the corners in projective coordinates, including weights.

derivative(*params, **kwargs)

Evaluate the derivative of the object at the given parametric values.

If tensor is true, evaluation will take place on a tensor product grid, i.e. it will return an n1 × n2 × … × dim array, where ni is the number of evaluation points in direction i, and dim is the physical dimension of the object.

If tensor is false, there must be an equal number n of evaluation points in all directions, and the return value will be an n × dim array.

If there is only one evaluation point, a vector of length dim is returned instead.

Examples:

# Tangent of curve at single point
curve.derivative(1.0)

# Double derivative of curve at single point:
curve.derivative(1.0, d=2)

# Third derivative of curve at several points:
curve.derivative([0.0, 1.0, 2.0], d=3)

# Tangents of surface:
surface.derivative(0.5, 0.7, d=(1,0))
surface.derivative(0.5, 0.7, d=(0,1))

# Cross-derivative of surface:
surface.derivative(0.5, 0.7, d=(1,1))
Parameters:
  • u,v,.. (float or [float]) – Parametric coordinates in which to evaluate
  • d ((int)) – Order of derivative to compute
  • above ((bool)) – Evaluation in the limit from above
  • tensor (bool) – Whether to evaluate on a tensor product grid
Returns:

Derivatives

Return type:

numpy.array

edges()[source]

Return the twelve edges of this volume in order:

  • umin, vmin
  • umax, vmin
  • umin, vmax
  • umax, vmax
  • umin, wmin
  • umax, wmin
  • umin, wmax
  • umax, wmax
  • vmin, wmin
  • vmax, wmin
  • vmin, wmax
  • vmax, wmax
Returns:Edges
Return type:(Curve)
end(direction=None)

Return the end of the parametric domain.

If direction is given, returns the end of that direction, as a float. If it is not given, returns the end of all directions, as a tuple.

Parameters:direction (int) – Direction in which to get the end.
Raises:ValueError – For invalid direction
faces()[source]

Return the six faces of this volume in order: umin, umax, vmin, vmax, wmin, wmax.

Returns:Boundary faces
Return type:(Surface)
force_rational()

Force a rational representation of the object.

The weights of a non-rational object will be set to 1.

Returns:self
get_derivative_spline(direction=None)

Compute the controlpoints associated with the derivative spline object

If direction is given, only the derivatives in that direction are returned.

If direction is not given, this function returns a tuple of all partial derivatives

# Create a 4x4 element cubic spline surface
surf = Surface()
surf.raise_order(2,2)
surf.refine(3,3)
surf[1:4,1:4,:] += 0.1 # make the surface non-trivial by moving controlpoints

# Create the derivative surface
du = surf.get_derivative_spline(direction='u')

# evaluation is identical
print(du.evaluate(0.3, 0.4))
print(surf.derivative(0.3, 0.4, d=(1,0)))

print(surf.order()) # prints (3,3)
print(du.order())   # prints (2,3)
Parameters:direction (int) – The tangential direction
Returns:Derivative spline
Return type:SplineObject
get_derivative_volume(direction=None)

Compute the controlpoints associated with the derivative spline object

If direction is given, only the derivatives in that direction are returned.

If direction is not given, this function returns a tuple of all partial derivatives

# Create a 4x4 element cubic spline surface
surf = Surface()
surf.raise_order(2,2)
surf.refine(3,3)
surf[1:4,1:4,:] += 0.1 # make the surface non-trivial by moving controlpoints

# Create the derivative surface
du = surf.get_derivative_spline(direction='u')

# evaluation is identical
print(du.evaluate(0.3, 0.4))
print(surf.derivative(0.3, 0.4, d=(1,0)))

print(surf.order()) # prints (3,3)
print(du.order())   # prints (2,3)
Parameters:direction (int) – The tangential direction
Returns:Derivative spline
Return type:SplineObject
insert_knot(knot, direction=0)

Insert a new knot into the spline.

Parameters:
  • direction (int) – The direction to insert in
  • knot (float or [float]) – The new knot(s) to insert
Raises:

ValueError – For invalid direction

Returns:

self

knots(direction=None, with_multiplicities=False)

Return knots vector

If direction is given, returns the knots in that direction, as a list. If it is not given, returns the knots of all directions, as a tuple.

Parameters:
  • direction (int) – Direction in which to get the knots.
  • with_multiplicities (bool) – If true, return knots with multiplicities (i.e. repeated).
Raises:

ValueError – For invalid direction

lower_order(*lowers)

Lower the polynomial order of the object. If only one argument is given, the order is lowered equally over all directions.

Parameters:u,v,.. (int) – Number of times to lower the order in a given direction.
Return SplineObject:
 Approximation of the current object on a lower order basis
lower_periodic(periodic, direction=0)

Sets the periodicity of the spline object in the given direction, keeping the geometry unchanged.

Parameters:
  • periodic (int) – new periodicity, i.e. the basis is C^k over the start/end
  • direction (int) – the parametric direction of the basis to modify
Returns:

self

make_periodic(continuity=None, direction=0)

Make the spline object periodic in a given parametric direction.

Parameters:
  • continuity (int) – The continuity along the boundary (default max).
  • direction (int) – The direction to ensure continuity in.
classmethod make_splines_compatible(spline1, spline2)

Ensure that two splines are compatible.

This will manipulate one or both to ensure that they are both rational or nonrational, and that they lie in the same physical space.

Parameters:
classmethod make_splines_identical(spline1, spline2)

Ensure that two splines have identical discretization.

This will first make them compatible (see splipy.SplineObject.make_curves_compatible()), reparametrize them, and possibly raise the order and insert knots as required.

Parameters:
mirror(normal)

Mirror the object around a plane through the origin.

Parameters:normal (array-like) – The plane normal to mirror about.
Raises:RuntimeError – If the physical dimension is not 2 or 3
Returns:self
order(direction=None)

Return polynomial order (degree + 1).

If direction is given, returns the order of that direction, as an int. If it is not given, returns the order of all directions, as a tuple.

Parameters:direction (int) – Direction in which to get the order.
Raises:ValueError – For invalid direction
pardim

The number of parametric dimensions: 1 for curves, 2 for surfaces, 3 for volumes, etc.

periodic(direction=0)

Returns true if the spline object is periodic in the given parametric direction

project(plane)

Projects the geometry onto a plane or axis.

  • project(‘xy’) will project the object onto the xy plane, setting all z components to zero.
  • project(‘y’) will project the object onto the y axis, setting all x and z components to zero.
Parameters:plane (string) – Any combination of ‘x’, ‘y’ and ‘z’
Returns:self
raise_order(*raises)

Raise the polynomial order of the object. If only one argument is given, the order is raised equally over all directions.

Parameters:u,v,.. (int) – Number of times to raise the order in a given direction.
Returns:self
rebuild(p, n)[source]

Creates an approximation to this volume by resampling it using uniform knot vectors of order p with n control points.

Parameters:
  • p ((int)) – Tuple of polynomial discretization order in each direction
  • n ((int)) – Tuple of number of control points in each direction
Returns:

A new approximate volume

Return type:

Volume

refine(*ns, **kwargs)

Enrich the spline space by inserting knots into each existing knot span.

This method supports three different usage patterns:

# Refine each direction by a factor n
obj.refine(n)

# Refine a single direction by a factor n
obj.refine(n, direction='v')

# Refine all directions by given factors
obj.refine(nu, nv, ...)
Parameters:
  • nu,nv,.. (int) – Number of new knots to insert into each span
  • direction (int) – Direction to refine in
Returns:

self

reparam(*args, **kwargs)

Redefine the parametric domain. This function accepts two calling conventions:

reparametrize(u, v, …) reparametrizes each direction to the domains given by the tuples u, v, etc. It is equivalent to calling reparametrize(u[0], u[1]) on each basis. The default domain for directions not given is (0,1). In particular, if no arguments are given, the new parametric domain will be the unit (hyper)cube.

reparametrize(u, direction=d) reparametrizes just the direction given by d and leaves the others untouched.

Parameters:
  • u, v, .. (tuple) – New parametric domains, default to (0,1)
  • direction (int) – The direction to reparametrize
Returns:

self

reverse(direction=0)

Swap the direction of a parameter by making it go in the reverse direction. The parametric domain remains unchanged.

Parameters:direction (int) – The direction to flip.
Returns:self
rotate(theta, normal=(0, 0, 1))

Rotate the object around an axis.

Parameters:
  • theta (float) – Angle to rotate about, measured in radians
  • normal (array-like) – The normal axis (if 3D) to rotate about
Raises:

RuntimeError – If the physical dimension is not 2 or 3

Returns:

self

scale(*args)

Scale, or magnify the object by a given amount.

In case of one input argument, the scaling is uniform.

Parameters:args (array-like or float) – Scaling factors, possibly different in each direction.
Returns:self
section(*args, **kwargs)

Returns a section from the object. A section can be any sub-object of parametric dimension not exceeding that of the object. E.g. for a volume, sections include vertices, edges, faces, etc.

The arguments are control point indices for each direction. None means that direction should be variable in the returned object.

# Get the face with u=max
vol.section(-1, None, None)

# Keyword arguments are supported for u, v and w
# This is the same as the above
vol.section(u=-1)

# Get the edge with u=min, v=max
vol.section(0, -1, None)

# This is equivalent to vol.clone()
vol.section()

If a specific subclass of SplineObject is found that handles the requested number of variable directions (parametric dimension), then the return value is of that type. If not, it will be a generic SplineObject.

If the section has no variable directions (it is a point), then the return value will be an array, unless the keyword argument unwrap_points is true, in which case it will return a zero-dimensional SplineObject.

Parameters:u,v,.. (int or None) – Control point indices
Returns:Section
Return type:SplineObject or np.array
set_dimension(new_dim)

Sets the physical dimension of the object. If increased, the new components are set to zero.

Parameters:new_dim (int) – New dimension.
Returns:self
set_order(*order)

Set the polynomial order of the object. If only one argument is given, the order is set uniformly over all directions.

Parameters:u,v,.. (int) – The new order in a given direction.
Raises:ValueError – If the order is reduced in any direction.
Returns:self
shape

The dimensions of the control point array.

split(knots, direction=0)

Split an object into two or more separate representations with C0 continuity between them.

Parameters:
  • knots (float or [float]) – The splitting points
  • direction (int) – Parametric direction
Returns:

The new objects

Return type:

[SplineObject]

start(direction=None)

Return the start of the parametric domain.

If direction is given, returns the start of that direction, as a float. If it is not given, returns the start of all directions, as a tuple.

Parameters:direction (int) – Direction in which to get the start.
Raises:ValueError – For invalid direction
swap(dir1=0, dir2=1)

Swaps two parameter directions.

This function silently passes for curves.

Parameters:
  • dir1 (direction) – The first direction (default u)
  • dir2 (direction) – The second direction (default v)
Returns:

self

tangent(*params, **kwargs)

Evaluate the tangents of the object at the given parametric values.

If direction is given, only the derivatives in that direction are evaluated. This is equivalent to calling splipy.SplineObject.derivative() with d=(0,…,0,1,0,…,0), the unit vector corresponding to the given direction.

If direction is not given, this function returns a tuple of all tangential derivatives at the given points.

Parameters:
  • u,v,.. (float or [float]) – Parametric coordinates in which to evaluate
  • direction (int) – The tangential direction
  • above ((bool)) – Evaluation in the limit from above
  • tensor (bool) – Whether to evaluate on a tensor product grid
Returns:

Tangents

Return type:

tuple<numpy.array>

translate(x)

Translate (i.e. move) the object by a given distance.

Parameters:x (array-like) – The vector to translate by.
Returns:self
volume()[source]

Computes the volume of the object in geometric space