core.bsplines_kernels#

Functions#

basis_ders_on_irregular_grid_p(knots, ...)

Evaluate B-Splines and their derivatives on an irregular_grid.

basis_ders_on_quad_grid_p(knots, degree, ...)

Evaluate B-Splines and their derivatives on the quadrature grid.

basis_funs_1st_der_p(knots, degree, x, span, out)

Compute the first derivative of the non-vanishing B-splines at a location.

basis_funs_all_ders_p(knots, degree, x, ...)

Evaluate value and n derivatives at x of all basis functions with support in interval \([x_{span-1}, x_{span}]\).

basis_funs_array_p(knots, degree, x, span, out)

Compute the non-vanishing B-splines at locations in x, given the knot sequence, polynomial degree and knot span.

basis_funs_p(knots, degree, x, span, out)

Compute the non-vanishing B-splines at a unique location.

basis_integrals_p(knots, degree, out)

Return the integral of each B-spline basis function over the real line:

breakpoints_p(knots, degree, out[, tol])

Determine breakpoints' coordinates.

cell_index_p(breaks, i_grid, tol, out)

Computes in which cells a given sorted array of locations belong.

collocation_matrix_p(knots, degree, ...[, ...])

Compute the collocation matrix \(C_ij = B_j(x_i)\), which contains the values of each B-spline basis function \(B_j\) at all locations \(x_i\).

elements_spans_p(knots, degree, out)

Compute the index of the last non-vanishing spline on each grid element (cell).

elevate_knots_p(knots, degree, periodic, out)

Given the knot sequence of a spline space S of degree p, compute the knot sequence of a spline space S_0 of degree p+1 such that u' is in S for all u in S_0.

find_span_p(knots, degree, x)

Determine the knot span index at location x, given the B-Splines' knot sequence and polynomial degree.

find_spans_p(knots, degree, x, out)

Determine the knot span index at a set of locations x, given the B-Splines' knot sequence and polynomial degree.

greville_p(knots, degree, periodic, out[, ...])

Compute coordinates of all Greville points.

histopolation_matrix_p(knots, degree, ...[, ...])

Computes the histopolation matrix.

make_knots_p(breaks, degree, periodic, out)

Create spline knots from breakpoints, with appropriate boundary conditions.

merge_sort(a)

Performs a 'in place' merge sort of the input list

quadrature_grid_p(breaks, quad_rule_x, ...)

Compute the quadrature points and weights for performing integrals over each element (interval) of the 1D domain, given a certain Gaussian quadrature rule.

Details#

find_span_p(knots: float[:], degree: int, x: float)[source]#

Determine the knot span index at location x, given the B-Splines’ knot sequence and polynomial degree. See Algorithm A2.1 in [1].

For a degree p, the knot span index i identifies the indices [i-p:i] of all p+1 non-zero basis functions at a given location x.

Parameters:
knotsarray_like

Knots sequence.

degreeint

Polynomial degree of B-splines.

xfloat

Location of interest.

Returns:
spanint

Knot span index.

References

[1]

L. Piegl and W. Tiller. The NURBS Book, 2nd ed., Springer-Verlag Berlin Heidelberg GmbH, 1997.

find_spans_p(knots: float[:], degree: int, x: float[:], out: int[:])[source]#

Determine the knot span index at a set of locations x, given the B-Splines’ knot sequence and polynomial degree. See Algorithm A2.1 in [1].

For a degree p, the knot span index i identifies the indices [i-p:i] of all p+1 non-zero basis functions at a given location x.

Parameters:
knotsarray_like

Knots sequence.

degreeint

Polynomial degree of B-splines.

xarray_like of floats

Locations of interest.

outarray

The result will be inserted into this array. It should be of the appropriate shape and dtype.

basis_funs_p(knots: float[:], degree: int, x: float, span: int, out: float[:])[source]#

Compute the non-vanishing B-splines at a unique location.

Parameters:
knotsarray_like of floats

Knots sequence.

degreeint

Polynomial degree of B-splines.

xfloat

Evaluation point.

spanint

Knot span index.

outarray

The result will be inserted into this array. It should be of the appropriate shape and dtype.

Notes

The original Algorithm A2.2 in The NURBS Book [1] is here slightly improved by using ‘left’ and ‘right’ temporary arrays that are one element shorter.

References

[1]

L. Piegl and W. Tiller. The NURBS Book, 2nd ed., Springer-Verlag Berlin Heidelberg GmbH, 1997.

basis_funs_array_p(knots: float[:], degree: int, x: float[:], span: int[:], out: float[:, :])[source]#

Compute the non-vanishing B-splines at locations in x, given the knot sequence, polynomial degree and knot span. See Algorithm A2.2 in [1].

Parameters:
knotsarray_like of floats

Knots sequence.

degreeint

Polynomial degree of B-splines.

xarray_like of floats

Evaluation points.

spanarray_like of int

Knot span indexes.

outarray

The result will be inserted into this array. It should be of the appropriate shape and dtype.

References

[1]

L. Piegl and W. Tiller. The NURBS Book, 2nd ed., Springer-Verlag Berlin Heidelberg GmbH, 1997.

basis_funs_1st_der_p(knots: float[:], degree: int, x: float, span: int, out: float[:])[source]#

Compute the first derivative of the non-vanishing B-splines at a location.

Parameters:
knotsarray_like

Knots sequence.

degreeint

Polynomial degree of B-splines.

xfloat

Evaluation point.

spanint

Knot span index.

outarray

The result will be inserted into this array. It should be of the appropriate shape and dtype.

Notes

See function ‘s_bsplines_non_uniform__eval_deriv’ in Selalib’s ([2]) source file ‘src/splines/sll_m_bsplines_non_uniform.F90’.

References

[2]

SELALIB, Semi-Lagrangian Library. http://selalib.gforge.inria.fr

basis_funs_all_ders_p(knots: float[:], degree: int, x: float, span: int, n: int, normalization: bool, out: float[:, :])[source]#

Evaluate value and n derivatives at x of all basis functions with support in interval \([x_{span-1}, x_{span}]\).

If called with normalization=True, this uses M-splines instead of B-splines.

Fills a 2D array with n+1 (from 0-th to n-th) derivatives at x of all (degree+1) non-vanishing basis functions in given span.

\[ders[i,j] = \frac{d^i}{dx^i} B_k(x) \, \text{with} k=(span-degree+j), \forall (i,j), 0 \leq i \leq n \, 0 \leq j \leq \text{degree}+1.\]
Parameters:
knotsarray_like

Knots sequence.

degreeint

Polynomial degree of B-splines.

xfloat

Evaluation point.

spanint

Knot span index.

nint

Max derivative of interest.

normalization: bool

Set to False to get B-Splines and True to get M-Splines

outarray

The result will be inserted into this array. It should be of the appropriate shape and dtype.

Notes

The original Algorithm A2.3 in The NURBS Book [1] is here improved:
  • ‘left’ and ‘right’ arrays are 1 element shorter;

  • inverse of knot differences are saved to avoid unnecessary divisions;

  • innermost loops are replaced with vector operations on slices.

References

[1]

L. Piegl and W. Tiller. The NURBS Book, 2nd ed., Springer-Verlag Berlin Heidelberg GmbH, 1997.

basis_integrals_p(knots: float[:], degree: int, out: float[:])[source]#

Return the integral of each B-spline basis function over the real line:

Math:

K[i] := int_{-infty}^{+infty} B[i](x) dx = (T[i+p+1]-T[i]) / (p+1).

This array can be used to convert B-splines to M-splines, which have unit integral over the real line but no partition-of-unity property.

Parameters:
knotsarray_like

Knots sequence.

degreeint

Polynomial degree of B-splines.

outarray

The result will be inserted into this array. It should be of the appropriate shape and dtype.

Notes

For convenience, this function does not distinguish between periodic and non-periodic spaces, hence the length of the output array is always equal to (len(knots)-degree-1). In the periodic case the last (degree) values in the array are redundant, as they are a copy of the first (degree) values.

collocation_matrix_p(knots: float[:], degree: int, periodic: bool, normalization: bool, xgrid: float[:], out: float[:, :], multiplicity: int = 1)[source]#

Compute the collocation matrix \(C_ij = B_j(x_i)\), which contains the values of each B-spline basis function \(B_j\) at all locations \(x_i\).

If called with normalization=True, this uses M-splines instead of B-splines.

Parameters:
knotsarray_like

Knots sequence.

degreeint

Polynomial degree of spline space.

periodicbool

True if domain is periodic, False otherwise.

normalizationbool

Set to False for B-splines, and True for M-splines.

xgridarray_like

Evaluation points.

outarray

The result will be inserted into this array. It should be of the appropriate shape and dtype.

multiplicityint

Multiplicity of the knots in the knot sequence, we assume that the same multiplicity applies to each interior knot.

histopolation_matrix_p(knots: float[:], degree: int, periodic: bool, normalization: bool, xgrid: float[:], check_boundary: bool, elevated_knots: float[:], out: float[:, :], multiplicity: int = 1)[source]#

Computes the histopolation matrix.

If called with normalization=True, this uses M-splines instead of B-splines.

Parameters:
knotsarray_like

Knots sequence.

degreeint

Polynomial degree of spline space.

periodicbool

True if domain is periodic, False otherwise.

normalizationstr

Set to False for B-splines, and True for M-splines.

xgridarray_like

Grid points.

check_boundarybool, default=True

If true and periodic, will check the boundaries of xgrid.

outarray

The result will be inserted into this array. It should be of the appropriate shape and dtype.

multiplicityint

Multiplicity of the knots in the knot sequence, we assume that the same multiplicity applies to each interior knot.

Notes

The histopolation matrix \(H_{ij} = \int_{x_i}^{x_{i+1}}B_j(x)\,dx\) contains the integrals of each B-spline basis function \(B_j\) between two successive grid points.

merge_sort(a: float[:]) float[:][source]#

Performs a ‘in place’ merge sort of the input list

Parameters:
aarray_like

1D array/list to sort using merge sort

breakpoints_p(knots: float[:], degree: int, out: float[:], tol: float = 1e-15)[source]#

Determine breakpoints’ coordinates.

Parameters:
knotsarray_like

Knots sequence.

degreeint

Polynomial degree of B-splines.

tol: float

If the distance between two knots is less than tol, we assume that they are repeated knots which correspond to the same break point.

outarray

The result will be inserted into this array. It should be of the appropriate shape and dtype.

Returns:
last_indexint

Last meaningful index + 1, e.g. the actual interesting result is out[:last_index].

greville_p(knots: float[:], degree: int, periodic: bool, out: float[:], multiplicity: int = 1)[source]#

Compute coordinates of all Greville points.

Parameters:
knotsarray_like

Knots sequence.

degreeint

Polynomial degree of B-splines.

periodicbool

True if domain is periodic, False otherwise.

outarray

The result will be inserted into this array. It should be of the appropriate shape and dtype.

multiplicityint

Multiplicity of the knots in the knot sequence, we assume that the same multiplicity applies to each interior knot.

elements_spans_p(knots: float[:], degree: int, out: int[:])[source]#

Compute the index of the last non-vanishing spline on each grid element (cell). The length of the returned array is the number of cells.

Parameters:
knotsarray_like

Knots sequence.

degreeint

Polynomial degree of B-splines.

outarray

The result will be inserted into this array. It should be of the appropriate shape and dtype.

Returns:
last_indexint

Last meaningful index + 1, e.g. the actual interesting result is out[:last_index].

Notes

  1. Numbering of basis functions starts from 0, not 1;

  2. This function could be written in two lines:

    breaks = breakpoints( knots, degree ) spans = np.searchsorted( knots, breaks[:-1], side=’right’ ) - 1

make_knots_p(breaks: float[:], degree: int, periodic: bool, out: float[:], multiplicity: int = 1)[source]#

Create spline knots from breakpoints, with appropriate boundary conditions.

If domain is periodic, knot sequence is extended by periodicity to have a total of (n_cells-1)*mult+2p+2 knots (all break points are repeated mult time and we add p+1-mult knots by periodicity at each side).

Otherwise, knot sequence is clamped (i.e. endpoints have multiplicity p+1).

Parameters:
breaksarray_like

Coordinates of breakpoints (= cell edges); given in increasing order and with no duplicates.

degreeint

Spline degree (= polynomial degree within each interval).

periodicbool

True if domain is periodic, False otherwise.

multiplicity: int

Multiplicity of the knots in the knot sequence, we assume that the same multiplicity applies to each interior knot.

outarray

The result will be inserted into this array. It should be of the appropriate shape and dtype.

elevate_knots_p(knots: float[:], degree: int, periodic: bool, out: float[:], multiplicity: int = 1, tol: float = 1e-15)[source]#

Given the knot sequence of a spline space S of degree p, compute the knot sequence of a spline space S_0 of degree p+1 such that u’ is in S for all u in S_0.

Specifically, on bounded domains the first and last knots are repeated in the sequence, and in the periodic case the knot sequence is extended by periodicity.

Parameters:
knotsarray_like

Knots sequence of spline space of degree p.

degreeint

Spline degree (= polynomial degree within each interval).

periodicbool

True if domain is periodic, False otherwise.

multiplicityint

Multiplicity of the knots in the knot sequence, we assume that the same multiplicity applies to each interior knot.

tol: float

If the distance between two knots is less than tol, we assume that they are repeated knots which correspond to the same break point.

outarray

The result will be inserted into this array. It should be of the appropriate shape and dtype.

quadrature_grid_p(breaks: float[:], quad_rule_x: float[:], quad_rule_w: float[:], out1: float[:, :], out2: float[:, :])[source]#

Compute the quadrature points and weights for performing integrals over each element (interval) of the 1D domain, given a certain Gaussian quadrature rule.

An n-point Gaussian quadrature rule for the canonical interval \([-1,+1]\) and trivial weighting function \(\omega(x)=1\) is defined by the n abscissas \(x_i\) and n weights \(w_i\) that satisfy the following identity for polynomial functions \(f(x)\) of degree \(2n-1\) or less:

\[\int_{-1}^{+1} f(x) dx = \sum_{i=0}^{n-1} w_i f(x_i)\]
Parameters:
breaksarray_like of floats

Coordinates of spline breakpoints.

quad_rule_xarray_like of ints

Coordinates of quadrature points on canonical interval [-1,1].

quad_rule_warray_like of ints

Weights assigned to quadrature points on canonical interval [-1,1].

out1array

2D output array where the quadrature points will be inserted. It should be of the appropriate shape and dtype.

out2array

2D output array where the quadrature weights will be inserted. It should be of the appropriate shape and dtype.

Notes

Contents of 2D output arrays ‘out1’ and ‘out2’ are accessed with two indices (ie,iq) where:

  • ie is the global element index;

  • iq is the local index of a quadrature point within the element.

basis_ders_on_quad_grid_p(knots: float[:], degree: int, quad_grid: float[:, :], nders: int, normalization: bool, offset: int, out: float[:, :, :, :])[source]#

Evaluate B-Splines and their derivatives on the quadrature grid.

If called with normalization=True, this uses M-splines instead of B-splines.

Parameters:
knotsarray_like

Knots sequence.

degreeint

Polynomial degree of B-splines.

quad_grid: ndarray

2D array of shape (ne, nq). Coordinates of quadrature points of each element in 1D domain, which can be given by quadrature_grid() or chosen arbitrarily.

ndersint

Maximum derivative of interest.

normalizationbool

Set to False for B-splines, and True for M-splines.

offsetint, default=0

Assumes that the quadrature grid starts from cell number offset.

outarray

The result will be inserted into this array. It should be of the appropriate shape and dtype.

Notes

4D output array ‘out’ contains values of B-Splines and their derivatives at quadrature points in each element of 1D domain. Indices are . ie: global element (0 <= ie < ne ) . il: local basis function (0 <= il <= degree) . id: derivative (0 <= id <= nders ) . iq: local quadrature point (0 <= iq < nq )

cell_index_p(breaks: float[:], i_grid: float[:], tol: float, out: int[:])[source]#

Computes in which cells a given sorted array of locations belong.

Locations close to a interior breakpoint will be assumed to be present twice in the grid, once of for each cell. Boundary breakpoints are snapped to the interior of the domain.

Parameters:
breaksarray_like

Coordinates of breakpoints (= cell edges); given in increasing order and with no duplicates.

i_gridndarray

1D array of locations.

tolfloat

If the distance between a given point in i_grid and a breakpoint is less than tol then it is considered to be the breakpoint.

outarray

The result will be inserted into this array. It should be of the appropriate shape and dtype.

Returns:
statusint

0 if everything worked as intended, 1 if not.

basis_ders_on_irregular_grid_p(knots: float[:], degree: int, i_grid: float[:], cell_index: int[:], nders: int, normalization: bool, out: float[:, :, :])[source]#

Evaluate B-Splines and their derivatives on an irregular_grid.

If called with normalization=True, this uses M-splines instead of B-splines.

Parameters:
knotsarray_like

Knots sequence.

degreeint

Polynomial degree of B-splines.

i_gridndarray

1D array of all of the points on which to evaluate the basis functions. The points do not need to be sorted

cell_indexndarray

1D array of the same shape as i_grid. cell_index[i] is the index of the cell in which i_grid[i] belong.

ndersint

Maximum derivative of interest.

normalizationbool

Set to False for B-splines, and True for M-splines.

outarray

The result will be inserted into this array. It should be of the appropriate shape and dtype.

Notes

3D output array ‘out’ contains values of B-Splines and their derivatives at quadrature points in each element of 1D domain. Indices are . ie: location (0 <= ie < nx ) . il: local basis function (0 <= il <= degree) . id: derivative (0 <= id <= nders )