B-splines#

Basic modules#

Basic module that provides the means for evaluating the B-Splines basis functions and their derivatives. In order to simplify automatic Fortran code generation with Pyccel, no object-oriented features are employed.

References

[1] L. Piegl and W. Tiller. The NURBS Book, 2nd ed.,

Springer-Verlag Berlin Heidelberg GmbH, 1997.

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

struphy.bsplines.bsplines.find_span(knots, degree, x)[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:
  • knots (array_like) – Knots sequence.

  • degree (int) – Polynomial degree of B-splines.

  • x (float) – Location of interest.

Returns:

span – Knot span index.

Return type:

int

struphy.bsplines.bsplines.scaling_vector(knots, degree, span)[source]#

Returns the scaling array for M-splines.

Parameters:
  • knots (array_like) – Knots sequence.

  • degree (int) – Polynomial degree of B-splines.

  • span (int) – Knot span index.

Returns:

x – Scaling vector with elements (p + 1)/(t[i + p + 1] - t[i])

Return type:

array_like

struphy.bsplines.bsplines.basis_funs(knots, degree, x, span, normalize=False)[source]#

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

Parameters:
  • knots (array_like) – Knots sequence.

  • degree (int) – Polynomial degree of B-splines.

  • x (float) – Evaluation point.

  • span (int) – Knot span index.

  • normalize (boolean) – Scaling for M-splines.

Returns:

values – Values of p+1 non-vanishing B-Splines at location x.

Return type:

numpy.ndarray

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.

struphy.bsplines.bsplines.basis_funs_1st_der(knots, degree, x, span)[source]#

Compute the first derivative of the non-vanishing B-splines at location x, given the knot sequence, polynomial degree and knot span.

See function ‘s_bsplines_non_uniform__eval_deriv’ in Selalib’s source file ‘src/splines/sll_m_bsplines_non_uniform.F90’.

Parameters:
  • knots (array_like) – Knots sequence.

  • degree (int) – Polynomial degree of B-splines.

  • x (float) – Evaluation point.

  • span (int) – Knot span index.

  • Results

  • -------

  • ders (numpy.ndarray) – Derivatives of p + 1 non-vanishing B-Splines at location x.

struphy.bsplines.bsplines.basis_funs_all_ders(knots, degree, x, span, n)[source]#

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

ders[i,j] = (d/dx)^i B_k(x) with k=(span-degree+j),

for 0 <= i <= n and 0 <= j <= degree+1.

Parameters:
  • knots (array_like) – Knots sequence.

  • degree (int) – Polynomial degree of B-splines.

  • x (float) – Evaluation point.

  • span (int) – Knot span index.

  • n (int) – Max derivative of interest.

  • Results

  • -------

  • ders (numpy.ndarray (n+1,degree+1)) – 2D array of n+1 (from 0-th to n-th) derivatives at x of all (degree+1) non-vanishing basis functions in given span.

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.

struphy.bsplines.bsplines.collocation_matrix(knots, degree, xgrid, periodic, normalize=False)[source]#

Computes 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$.

Parameters:
  • knots (1D array_like) – Knots sequence.

  • degree (int) – Polynomial degree of B-splines.

  • xgrid (1D array_like) – Evaluation points.

  • periodic (boolean) – True if domain is periodic, False otherwise.

  • normalize (boolean) – Scaling for M-splines.

Returns:

mat – Collocation matrix: values of all basis functions on each point in xgrid.

Return type:

2D numpy.ndarray

struphy.bsplines.bsplines.histopolation_matrix(knots, degree, xgrid, periodic)[source]#

Computes the histopolation matrix $H_ij = int_{x_i}^x_{i+1} D_j(x) dx$ of the M_splines.

Parameters:
  • knots (1D array_like) – Knots sequence.

  • degree (int) – Polynomial degree of B-splines.

  • xgrid (1D array_like) – Evaluation points.

  • periodic (boolean) – True if domain is periodic, False otherwise.

Returns:

mat – Histopolation matrix: values of all integrals between two successive point of the given grid of all M-splines.

Return type:

2D numpy.ndarray

struphy.bsplines.bsplines.breakpoints(knots, degree)[source]#

Determine breakpoints’ coordinates.

Parameters:
  • knots (1D array_like) – Knots sequence.

  • degree (int) – Polynomial degree of B-splines.

Returns:

breaks – Abscissas of all breakpoints.

Return type:

numpy.ndarray (1D)

struphy.bsplines.bsplines.greville(knots, degree, periodic)[source]#

Compute coordinates of all Greville points.

Parameters:
  • knots (1D array_like) – Knots sequence.

  • degree (int) – Polynomial degree of B-splines.

  • periodic (bool) – True if domain is periodic, False otherwise.

Returns:

xg – Abscissas of all Greville points.

Return type:

numpy.ndarray (1D)

struphy.bsplines.bsplines.elements_spans(knots, degree)[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:
  • knots (1D array_like) – Knots sequence.

  • degree (int) – Polynomial degree of B-splines.

Returns:

spans – Index of last non-vanishing spline on each grid element.

Return type:

numpy.ndarray (1D)

Examples

>>> import numpy as np
>>> from psydac.core.bsplines import make_knots, elements_spans
>>> p = 3 ; n = 8
>>> grid  = xp.arange( n-p+1 )
>>> knots = make_knots( breaks=grid, degree=p, periodic=False )
>>> spans = elements_spans( knots=knots, degree=p )
>>> spans
array([3, 4, 5, 6, 7])

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 = xp.searchsorted( knots, breaks[:-1], side=’right’ ) - 1

struphy.bsplines.bsplines.make_knots(breaks, degree, periodic)[source]#

Create spline knots from breakpoints, with appropriate boundary conditions. Let p be spline degree. If domain is periodic, knot sequence is extended by periodicity so that first p basis functions are identical to last p. Otherwise, knot sequence is clamped (i.e. endpoints are repeated p times).

Parameters:
  • breaks (array_like) – Coordinates of breakpoints (= cell edges); given in increasing order and with no duplicates.

  • degree (int) – Spline degree (= polynomial degree within each interval).

  • periodic (bool) – True if domain is periodic, False otherwise.

  • Result

  • ------

  • T (numpy.ndarray (1D)) – Coordinates of spline knots.

struphy.bsplines.bsplines.quadrature_grid(breaks, quad_rule_x, quad_rule_w)[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:
  • breaks (1D array_like) – Coordinates of spline breakpoints.

  • quad_rule_x (1D array_like) – Coordinates of quadrature points on canonical interval [-1,1].

  • quad_rule_w (1D array_like) – Weights assigned to quadrature points on canonical interval [-1,1].

Returns:

  • quad_x (2D numpy.ndarray) – Abscissas of quadrature points on each element (interval) of the 1D domain. See notes below.

  • quad_w (2D numpy.ndarray) – Weights assigned to the quadrature points on each element (interval) of the 1D domain. See notes below.

Notes

Contents of 2D output arrays ‘quad_x’ and ‘quad_w’ 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.

struphy.bsplines.bsplines.basis_ders_on_quad_grid(knots, degree, quad_grid, nders, normalize=False)[source]#

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

Parameters:
  • knots (array_like) – Knots sequence.

  • degree (int) – Polynomial degree of B-splines.

  • quad_grid (2D numpy.ndarray (ne, nq)) – Coordinates of quadrature points of each element in 1D domain, given by quadrature_grid() function.

  • nders (int) – Maximum derivative of interest.

Returns:

basis – 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 )

Return type:

4D numpy.ndarray

Basic functions for point-wise B-spline evaluation

struphy.bsplines.bsplines_kernels.scaling(t_d: Final[float[:]], p_d: int, span_d: int, values: float[:])[source]#

Scaling coefficients for M-splines.

Parameters:
  • knots (array_like) – Knots sequence.

  • degree (int) – Polynomial degree of B-splines.

  • span (int) – Knot span index.

  • values (array[float]) – Output: scaling vector with elements (p + 1)/(t[i + p + 1] - t[i])

struphy.bsplines.bsplines_kernels.find_span(t: Final[float[:]], p: int, eta: float) int[source]#

Computes the knot span index i for which the B-splines i-p until i are non-vanishing at point eta.

Parameters:#

tarray

knot sequence

pinteger

degree of the basis splines

etafloat

Evaluation point

Returns:#

span-index

struphy.bsplines.bsplines_kernels.basis_funs(t: Final[float[:]], p: int, eta: float, span: int, left: float[:], right: float[:], values: float[:])[source]#
Parameters:
  • t (array[float]) – Knots sequence.

  • p (int) – Polynomial degree of B-splines.

  • eta (float) – Evaluation point.

  • span (int) – Knot span index.

  • left (array[float]) – Internal array of size p.

  • right (array[float]) – Internal array of size p.

  • values (array[float]) – Output: values of p + 1 non-vanishing B-Splines at location eta.

struphy.bsplines.bsplines_kernels.basis_funs_all(t: float[:], p: int, eta: float, span: int, left: float[:], right: float[:], values: float[:, :], diff: float[:])[source]#
Parameters:
  • t (array[float]) – Knots sequence.

  • p (int) – Polynomial degree of B-splines.

  • eta (float) – Evaluation point.

  • span (int) – Knot span index.

  • left (array[float]) – Internal array of size p.

  • right (array[float]) – Internal array of size p.

  • values (array[float]) – Outout: values of (p + 1, p + 1) non-vanishing B-Splines at location eta.

  • diff (array[float]) – Output: scaling array (p) for M-splines.

struphy.bsplines.bsplines_kernels.basis_funs_all_ders(knots: Final[float[:]], degree: int, eta: float, span: int, left: float[:], right: float[:], n: int, ders: float[:, :])[source]#

Copied from gvec_to_python, where it is tested.

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

ders[i,j] = (d/deta)^i B_k(eta) with k=(span-degree+j),

for 0 <= i <= n and 0 <= j <= degree+1.

Parameters:
  • knots (array_like) – Knots sequence.

  • degree (int) – Polynomial degree of B-splines.

  • eta (float) – Evaluation point.

  • span (int) – Knot span index.

  • n (int) – Max derivative of interest.

  • Results

  • -------

  • ders (numpy.ndarray (n+1,degree+1)) – 2D array of n+1 (from 0-th to n-th) derivatives at eta of all (degree+1) non-vanishing basis functions in given span.

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.

struphy.bsplines.bsplines_kernels.b_splines_slim(tn: float[:], pn: int, eta: float, span: int, values: float[:])[source]#

Computes the values of pn + 1 non-vanishing B-splines at position eta.

Parameters:
  • tn (array) – Knots sequence.

  • pn (int) – Polynomial degree of B-splines.

  • eta (float) – Evaluation point.

  • span (int) – Knot span index.

  • values (array[float]) – Outout: values of pn + 1 non-vanishing B-Splines at location eta.

struphy.bsplines.bsplines_kernels.d_splines_slim(tn: float[:], pn: int, eta: float, span: int, values: float[:])[source]#

Computes the values of pn non-vanishing D-splines at position eta.

Parameters:
  • tn (array) – Knot sequence of B-splines.

  • pn (int) – Polynomial degree of B-splines.

  • eta (float) – Evaluation point.

  • span (int) – Knot span index.

Returns:

values – Values of p non-vanishing D-Splines at location eta.

Return type:

array

struphy.bsplines.bsplines_kernels.b_d_splines_slim(tn: float[:], pn: int, eta: float, span: int, bn: float[:], bd: float[:])[source]#

One function to compute the values of non-vanishing B-splines and D-splines.

Parameters:
  • tn (array) – Knot sequence of B-splines.

  • pn (int) – Polynomial degree of B-splines.

  • span (integer) – Knot span index i -> [i-p,i] basis functions are non-vanishing.

  • eta (float) – Evaluation point.

  • bn (array[float]) – Output: pn + 1 non-vanishing B-splines at eta

  • bd (array[float]) – Output: pn non-vanishing D-splines at eta

struphy.bsplines.bsplines_kernels.b_der_splines_slim(tn: float[:], pn: int, eta: float, span: int, bn: float[:], der: float[:])[source]#
Parameters:
  • tn (array_like) – Knots sequence of B-splines.

  • pn (int) – Polynomial degree of B-splines.

  • eta (float) – Evaluation point.

  • span (int) – Knot span index.

  • bn (array[float]) – Out: pn + 1 values of B-splines at eta.

  • der (array[float]) – Out: pn + 1 values of derivatives of B-splines at eta.

struphy.bsplines.bsplines_kernels.basis_funs_and_der(t: float[:], p: int, eta: float, span: int, left: float[:], right: float[:], values: float[:, :], diff: float[:], der: float[:])[source]#
Parameters:
  • t (array_like) – Knots sequence.

  • p (int) – Polynomial degree of B-splines.

  • eta (float) – Evaluation point.

  • span (int) – Knot span index.

  • left (array_like) – p left values

  • right (array_like) – p right values

  • values_all (array_like)

Returns:

values – Values of (2, p + 1) non-vanishing B-Splines and derivatives at location eta.

Return type:

numpy.ndarray

struphy.bsplines.bsplines_kernels.basis_funs_1st_der(t: Final[float[:]], p: int, eta: float, span: int, left: float[:], right: float[:], values: float[:])[source]#
Parameters:
  • t (array_like) – Knots sequence.

  • p (int) – Polynomial degree of B-splines.

  • eta (float) – Evaluation point.

  • span (int) – Knot span index.

Returns:

values – Derivatives of p + 1 non-vanishing B-Splines at location eta.

Return type:

numpy.ndarray

struphy.bsplines.bsplines_kernels.b_spl_1st_der_slim(t: float[:], p: int, eta: float, span: int, values: float[:])[source]#
Parameters:
  • t (array_like) – Knots sequence.

  • p (int) – Polynomial degree of B-splines.

  • eta (float) – Evaluation point.

  • span (int) – Knot span index.

Returns:

values – Derivatives of p + 1 non-vanishing B-Splines at location eta.

Return type:

array

struphy.bsplines.bsplines_kernels.piecewise(p: int, delta: float, eta: float) float[source]#

evaluate a hat function (B-spline) centered at eta0 (the center of the support) at eta1, i.e.

\[1.0 / delta * S((eta1 - eta0)/ delta)\]

where S is the B-spline of degree p, delta is the cell size of a uniform mesh. Here we use the expression of the B-spline in each cell directly. For the moment, p = 0, 1, or 2.

Parameters:
  • p (int) – degree of the hat function (B-spline)

  • delta (float) – the cell size of a uniform mesh.

  • eta (float) – eta = eta1 - eta0

  • Returns

  • -------- – value of the hat function

struphy.bsplines.bsplines_kernels.piecewise_der(p: int, delta: float, eta: float) float[source]#

evaluate the derivative of a hat function (B-spline) centered at eta0 (the center of the support) at eta1, i.e. .. math:

1.0 / delta^2 * S'((eta1 - eta0)/ delta)

where S is the B-spline of degree p, delta is the cell size of a uniform mesh. Here we use the expression of the derivative of the B-spline in each cell directly. For the moment, p = 0, 1, or 2.

Parameters:
  • p (int) – degree of the hat function (B-spline)

  • delta (float) – the cell size of a uniform mesh.

  • eta (float) – eta = eta1 - eta0

  • Returns

  • -------- – value of the derivative of the hat function

struphy.bsplines.bsplines_kernels.convolution(p: int, grids: Final[float[:]], eta: float) float[source]#

evaluate a hat function (B-spline) at eta, i.e.

\[1.0 / delta * S(eta/ delta)\]

where S is the B-spline of degree p. Here we use the definition by convolution of the B-splines.

Parameters:
  • p (int) – degree of the hat function (B-spline)

  • grids (float array) – p + 2 points used in the definition of B-splines.

  • eta (float) – evluation point

  • Returns

  • -------- – value of the hat function at eta

struphy.bsplines.bsplines_kernels.convolution_der(p: int, grids: Final[float[:]], eta: float) float[source]#

evaluate the derivative of a hat function (B-spline) at eta, i.e.

\[1.0 / delta^2 * S'(eta/ delta)\]

where S is the B-spline of degree p. Here we use the definition by convolution of the B-splines.

Parameters:
  • p (int) – degree of the hat function (B-spline)

  • grids (float array) – p + 2 points used in the definition of B-splines.

  • eta (float) – evluation point

  • Returns

  • -------- – value of the derivative of the hat function at eta

Spline evaluation#

Acccelerated functions for point-wise evaluation of tensor product B-splines.

S(eta1, eta2, eta3) = sum_ijk [ c_ijk * B_i(eta1) * B_j(eta2) * B_k(eta3) ] with c_ijk in R.

Possible combinations for tensor product (BBB): * (NNN) * (DNN) * (NDN) * (NND) * (NDD) * (DND) * (DDN) * (DDD) * (dN/deta N N) * (N dN/deta N) * (N N dN/deta)

struphy.bsplines.evaluation_kernels_3d.evaluation_kernel_3d(p1: int, p2: int, p3: int, basis1: float[:], basis2: float[:], basis3: float[:], ind1: int[:], ind2: int[:], ind3: int[:], coeff: float[:, :, :]) float[source]#

Summing non-zero contributions of a spline function (serial, needs global arrays).

Parameters:
  • p1 (int) – Degrees of the univariate splines in each direction.

  • p2 (int) – Degrees of the univariate splines in each direction.

  • p3 (int) – Degrees of the univariate splines in each direction.

  • basis1 (array[float]) – The p + 1 values of non-zero basis splines at one point (eta1, eta2, eta3) in each direction.

  • basis2 (array[float]) – The p + 1 values of non-zero basis splines at one point (eta1, eta2, eta3) in each direction.

  • basis3 (array[float]) – The p + 1 values of non-zero basis splines at one point (eta1, eta2, eta3) in each direction.

  • ind1 (array[int]) – Global indices of non-vanishing splines in the element of the considered point.

  • ind2 (array[int]) – Global indices of non-vanishing splines in the element of the considered point.

  • ind3 (array[int]) – Global indices of non-vanishing splines in the element of the considered point.

  • coeff (array[float]) – The spline coefficients c_ijk.

Returns:

spline_value – Value of tensor-product spline at point (eta1, eta2, eta3).

Return type:

float

struphy.bsplines.evaluation_kernels_3d.evaluate_3d(kind1: int, kind2: int, kind3: int, t1: float[:], t2: float[:], t3: float[:], p1: int, p2: int, p3: int, ind1: int[:, :], ind2: int[:, :], ind3: int[:, :], coeff: float[:, :, :], eta1: float, eta2: float, eta3: float) float[source]#

Point-wise evaluation of a tensor-product spline (serial, needs global arrays).

Parameters:
  • kind1 (int) – Kind of univariate spline. 1 for B-spline, 2 for M-spline and 3 for derivative of B-spline.

  • kind2 (int) – Kind of univariate spline. 1 for B-spline, 2 for M-spline and 3 for derivative of B-spline.

  • kind3 (int) – Kind of univariate spline. 1 for B-spline, 2 for M-spline and 3 for derivative of B-spline.

  • t1 (array[float]) – Knot vectors of univariate splines.

  • t2 (array[float]) – Knot vectors of univariate splines.

  • t3 (array[float]) – Knot vectors of univariate splines.

  • p1 (int) – Degrees of univariate splines.

  • p2 (int) – Degrees of univariate splines.

  • p3 (int) – Degrees of univariate splines.

  • ind1 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • ind2 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • ind3 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • coeff (array[float]) – The spline coefficients c_ijk.

  • eta1 (float) – Point of evaluation.

  • eta2 (float) – Point of evaluation.

  • eta3 (float) – Point of evaluation.

Returns:

spline_value – Value of tensor-product spline at point (eta1, eta2, eta3).

Return type:

float

struphy.bsplines.evaluation_kernels_3d.evaluate_tensor_product(t1: float[:], t2: float[:], t3: float[:], p1: int, p2: int, p3: int, ind1: int[:, :], ind2: int[:, :], ind3: int[:, :], coeff: float[:, :, :], eta1: float[:], eta2: float[:], eta3: float[:], spline_values: float[:, :, :], kind: int)[source]#

Tensor-product evaluation of a tensor-product spline (serial, needs global arrays).

Parameters:
  • t1 (array[float]) – Knot vectors of univariate splines.

  • t2 (array[float]) – Knot vectors of univariate splines.

  • t3 (array[float]) – Knot vectors of univariate splines.

  • p1 (int) – Degrees of univariate splines.

  • p2 (int) – Degrees of univariate splines.

  • p3 (int) – Degrees of univariate splines.

  • ind1 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • ind2 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • ind3 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • coeff (array[float]) – The spline coefficients c_ijk.

  • eta1 (array[float]) – Points of evaluation in 1d arrays.

  • eta2 (array[float]) – Points of evaluation in 1d arrays.

  • eta3 (array[float]) – Points of evaluation in 1d arrays.

  • spline_values (array[float]) – Splines evaluated at points S_ijk = S(eta1_i, eta2_j, eta3_k).

  • kind (int) –

    Kind of spline to evaluate.
    • 0 : NNN

    • 11 : DNN

    • 12 : NDN

    • 13 : NND

    • 21 : NDD

    • 22 : DND

    • 23 : DDN

    • 3 : DDD

    • 41 : dN/deta N N

    • 42 : N dN/deta N

    • 43 : N N dN/deta

struphy.bsplines.evaluation_kernels_3d.evaluate_matrix(t1: float[:], t2: float[:], t3: float[:], p1: int, p2: int, p3: int, ind1: int[:, :], ind2: int[:, :], ind3: int[:, :], coeff: float[:, :, :], eta1: float[:, :, :], eta2: float[:, :, :], eta3: float[:, :, :], spline_values: float[:, :, :], kind: int)[source]#

General evaluation of a tensor-product spline (serial, needs global arrays).

Parameters:
  • t1 (array[float]) – Knot vectors of univariate splines.

  • t2 (array[float]) – Knot vectors of univariate splines.

  • t3 (array[float]) – Knot vectors of univariate splines.

  • p1 (int) – Degrees of univariate splines.

  • p2 (int) – Degrees of univariate splines.

  • p3 (int) – Degrees of univariate splines.

  • ind1 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • ind2 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • ind3 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • coeff (array[float]) – The spline coefficients c_ijk.

  • eta1 (array[float]) – Points of evaluation in 3d arrays such that shape(eta1) == shape(eta2) == shape(eta3).

  • eta2 (array[float]) – Points of evaluation in 3d arrays such that shape(eta1) == shape(eta2) == shape(eta3).

  • eta3 (array[float]) – Points of evaluation in 3d arrays such that shape(eta1) == shape(eta2) == shape(eta3).

  • spline_values (array[float]) – Splines evaluated at points S_ijk = S(eta1_i, eta2_j, eta3_k).

  • kind (int) –

    Kind of spline to evaluate.
    • 0 : NNN

    • 11 : DNN

    • 12 : NDN

    • 13 : NND

    • 21 : NDD

    • 22 : DND

    • 23 : DDN

    • 3 : DDD

    • 41 : dN/deta N N

    • 42 : N dN/deta N

    • 43 : N N dN/deta

struphy.bsplines.evaluation_kernels_3d.evaluate_sparse(t1: float[:], t2: float[:], t3: float[:], p1: int, p2: int, p3: int, ind1: int[:, :], ind2: int[:, :], ind3: int[:, :], coeff: float[:, :, :], eta1: float[:, :, :], eta2: float[:, :, :], eta3: float[:, :, :], spline_values: float[:, :, :], kind: int)[source]#

Evaluation of a tensor-product spline using sparse meshgrids (serial, needs global arrays).

Parameters:
  • t1 (array[float]) – Knot vectors of univariate splines.

  • t2 (array[float]) – Knot vectors of univariate splines.

  • t3 (array[float]) – Knot vectors of univariate splines.

  • p1 (int) – Degrees of univariate splines.

  • p2 (int) – Degrees of univariate splines.

  • p3 (int) – Degrees of univariate splines.

  • ind1 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • ind2 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • ind3 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • coeff (array[float]) – The spline coefficients c_ijk.

  • eta1 (array[float]) – Points of evaluation in 3d arrays such that shape(eta1) = (:,1,1), shape(eta2) = (1,:,1), shape(eta3) = (1,1,:).

  • eta2 (array[float]) – Points of evaluation in 3d arrays such that shape(eta1) = (:,1,1), shape(eta2) = (1,:,1), shape(eta3) = (1,1,:).

  • eta3 (array[float]) – Points of evaluation in 3d arrays such that shape(eta1) = (:,1,1), shape(eta2) = (1,:,1), shape(eta3) = (1,1,:).

  • spline_values (array[float]) – Splines evaluated at points S_ijk = S(eta1_i, eta2_j, eta3_k).

  • kind (int) –

    Kind of spline to evaluate.
    • 0 : NNN

    • 11 : DNN

    • 12 : NDN

    • 13 : NND

    • 21 : NDD

    • 22 : DND

    • 23 : DDN

    • 3 : DDD

    • 41 : dN/deta N N

    • 42 : N dN/deta N

    • 43 : N N dN/deta

struphy.bsplines.evaluation_kernels_3d.eval_spline_mpi_kernel(p1: int, p2: int, p3: int, basis1: float[:], basis2: float[:], basis3: float[:], span1: int, span2: int, span3: int, _data: float[:, :, :], starts: int[:]) float[source]#

Summing non-zero contributions of a spline function with distributed memory (domain decomposition).

Parameters:
  • p1 (int) – Degrees of the univariate splines in each direction.

  • p2 (int) – Degrees of the univariate splines in each direction.

  • p3 (int) – Degrees of the univariate splines in each direction.

  • basis1 (array[float]) – The p + 1 values of non-zero basis splines at one point (eta1, eta2, eta3) in each direction.

  • basis2 (array[float]) – The p + 1 values of non-zero basis splines at one point (eta1, eta2, eta3) in each direction.

  • basis3 (array[float]) – The p + 1 values of non-zero basis splines at one point (eta1, eta2, eta3) in each direction.

  • span1 (int) – Knot span index in each direction.

  • span2 (int) – Knot span index in each direction.

  • span3 (int) – Knot span index in each direction.

  • _data (array[float]) – The spline coefficients c_ijk of the current process, ie. the _data attribute of a StencilVector.

  • starts (array[int]) – Starting indices of current process.

Returns:

spline_value – Value of tensor-product spline at point (eta1, eta2, eta3).

Return type:

float

struphy.bsplines.evaluation_kernels_3d.eval_spline_derivative_mpi_kernel(p1: int, p2: int, p3: int, basis1: float[:], basis2: float[:], basis3: float[:], span1: int, span2: int, span3: int, _data: float[:, :, :], starts: int[:], direction: int) float[source]#

Kernel for the derivative of a spline in one direction (distributed, domain decomp.).

Parameters:
  • p1 (int) – Degrees of the univariate splines in each direction.

  • p2 (int) – Degrees of the univariate splines in each direction.

  • p3 (int) – Degrees of the univariate splines in each direction.

  • basis1 (array[float]) – The p + 1 values of non-zero basis splines at one point (eta1, eta2, eta3) in each direction.

  • basis2 (array[float]) – The p + 1 values of non-zero basis splines at one point (eta1, eta2, eta3) in each direction.

  • basis3 (array[float]) – The p + 1 values of non-zero basis splines at one point (eta1, eta2, eta3) in each direction.

  • span1 (int) – Knot span index in each direction.

  • span2 (int) – Knot span index in each direction.

  • span3 (int) – Knot span index in each direction.

  • _data (array[float]) – The spline coefficients c_ijk in current process, ie. the _data attribute of a StencilVector.

  • starts (array[int]) – Starting indices of current process.

Returns:

spline_value – Derivative in one direction of tensor-product spline at point (eta1, eta2, eta3).

Return type:

float

struphy.bsplines.evaluation_kernels_3d.eval_spline_mpi(eta1: float, eta2: float, eta3: float, _data: float[:, :, :], kind: int[:], pn: int[:], tn1: float[:], tn2: float[:], tn3: float[:], starts: int[:]) float[source]#

Point-wise evaluation of a tensor-product spline, distributed.

Parameters:
  • eta1 (float) – Evaluation point in [0, 1]^3.

  • eta2 (float) – Evaluation point in [0, 1]^3.

  • eta3 (float) – Evaluation point in [0, 1]^3.

  • _data (array[float]) – The spline coefficients c_ijk.

  • kind (array[int]) – Kind of 1d basis in each direction: 0 = N-spline, 1 = D-spline.

  • pn (array[int]) – Spline degrees of V0 in each direction.

  • tn1 (array[float]) – Knot vectors of V0 in each direction.

  • tn2 (array[float]) – Knot vectors of V0 in each direction.

  • tn3 (array[float]) – Knot vectors of V0 in each direction.

  • starts (array[float]) – Start indices of splines on current process.

Returns:

value – Value of tensor-product spline at point (eta_1, eta_2, eta_3).

Return type:

float

struphy.bsplines.evaluation_kernels_3d.eval_spline_mpi_tensor_product(eta1: float[:], eta2: float[:], eta3: float[:], _data: float[:, :, :], kind: int[:], pn: int[:], tn1: float[:], tn2: float[:], tn3: float[:], starts: int[:], values: float[:, :, :])[source]#

Tensor-product evaluation of a tensor-product spline, distributed.

Parameters:
  • eta1 (array[float]) – Evaluation points as 1d arrays; points not on local process domain must be flagged as -1. The out array is global and contains the spline values as S_ijk = S(eta1[i], eta2[j], eta3[k]).

  • eta2 (array[float]) – Evaluation points as 1d arrays; points not on local process domain must be flagged as -1. The out array is global and contains the spline values as S_ijk = S(eta1[i], eta2[j], eta3[k]).

  • eta3 (array[float]) – Evaluation points as 1d arrays; points not on local process domain must be flagged as -1. The out array is global and contains the spline values as S_ijk = S(eta1[i], eta2[j], eta3[k]).

  • _data (array[float]) – The spline coefficients c_ijk.

  • kind (array[int]) – Kind of 1d basis in each direction: 0 = N-spline, 1 = D-spline.

  • pn (array[int]) – Spline degrees of V0 in each direction.

  • tn1 (array[float]) – Knot vectors of V0 in each direction.

  • tn2 (array[float]) – Knot vectors of V0 in each direction.

  • tn3 (array[float]) – Knot vectors of V0 in each direction.

  • starts (array[float]) – Start indices of splines on current process.

  • values (array[float]) – Return array for spline values S_ijk = S(eta1[i], eta2[j], eta3[k]).

struphy.bsplines.evaluation_kernels_3d.eval_spline_mpi_tensor_product_fast(eta1: float[:], eta2: float[:], eta3: float[:], _data: float[:, :, :], kind: int[:], pn: int[:], tn1: float[:], tn2: float[:], tn3: float[:], starts: int[:], values: float[:, :, :])[source]#

Tensor-product evaluation of a tensor-product spline, distributed, and optimized by computing spans and spline values only once.

Parameters:
  • eta1 (array[float]) – Evaluation points as 1d arrays; points not on local process domain must be flagged as -1. Spline values are obtained as S_ijk = S(eta1[i], eta2[j], eta3[k]).

  • eta2 (array[float]) – Evaluation points as 1d arrays; points not on local process domain must be flagged as -1. Spline values are obtained as S_ijk = S(eta1[i], eta2[j], eta3[k]).

  • eta3 (array[float]) – Evaluation points as 1d arrays; points not on local process domain must be flagged as -1. Spline values are obtained as S_ijk = S(eta1[i], eta2[j], eta3[k]).

  • _data (array[float]) – The spline coefficients c_ijk.

  • kind (array[int]) – Kind of 1d basis in each direction: 0 = N-spline, 1 = D-spline.

  • pn (array[int]) – Spline degrees of V0 in each direction.

  • tn1 (array[float]) – Knot vectors of V0 in each direction.

  • tn2 (array[float]) – Knot vectors of V0 in each direction.

  • tn3 (array[float]) – Knot vectors of V0 in each direction.

  • starts (array[float]) – Start indices of splines on current process.

  • values (array[float]) – Return array for spline values S_ijk = S(eta1[i], eta2[j], eta3[k]).

struphy.bsplines.evaluation_kernels_3d.eval_spline_mpi_tensor_product_fixed(span1s: int[:], span2s: int[:], span3s: int[:], b1s: float[:, :], b2s: float[:, :], b3s: float[:, :], _data: float[:, :, :], kind: int[:], pn: int[:], starts: int[:], values: float[:, :, :])[source]#

Tensor-product evaluation of a tensor-product spline, distributed, and optimized for a fixed grid (spans and spline values have been pre-evaluated).

Parameters:
  • span1s (array[int]) – Knot span indices.

  • span2s (array[int]) – Knot span indices.

  • span3s (array[int]) – Knot span indices.

  • b1s (array[float]) – Values of p+1 non-zero basis functions at evaluation points.

  • b2s (array[float]) – Values of p+1 non-zero basis functions at evaluation points.

  • b3s (array[float]) – Values of p+1 non-zero basis functions at evaluation points.

  • _data (array[float]) – The spline coefficients c_ijk.

  • kind (array[int]) – Kind of 1d basis in each direction: 0 = N-spline, 1 = D-spline.

  • pn (array[int]) – Spline degrees of V0 in each direction.

  • starts (array[float]) – Start indices of splines on current process.

  • values (array[float]) – Return array for spline values S_ijk where ijk are the flattened indices corresponding to spans.

struphy.bsplines.evaluation_kernels_3d.eval_spline_mpi_matrix(eta1: float[:, :, :], eta2: float[:, :, :], eta3: float[:, :, :], _data: float[:, :, :], kind: int[:], pn: int[:], tn1: float[:], tn2: float[:], tn3: float[:], starts: int[:], values: float[:, :, :])[source]#

3d array evaluation of a tensor-product spline, distributed.

Parameters:
  • eta1 (array[float]) – Evaluation points as 3d arrays; points not on local process domain must be flagged as -1. Spline values are obtained as S_ijk = S(eta1[i,j,k], eta2[i,j,k], eta3[i,j,k]).

  • eta2 (array[float]) – Evaluation points as 3d arrays; points not on local process domain must be flagged as -1. Spline values are obtained as S_ijk = S(eta1[i,j,k], eta2[i,j,k], eta3[i,j,k]).

  • eta3 (array[float]) – Evaluation points as 3d arrays; points not on local process domain must be flagged as -1. Spline values are obtained as S_ijk = S(eta1[i,j,k], eta2[i,j,k], eta3[i,j,k]).

  • _data (array[float]) – The spline coefficients c_ijk.

  • kind (array[int]) – Kind of 1d basis in each direction: 0 = N-spline, 1 = D-spline.

  • pn (array[int]) – Spline degrees of V0 in each direction.

  • tn1 (array[float]) – Knot vectors of V0 in each direction.

  • tn2 (array[float]) – Knot vectors of V0 in each direction.

  • tn3 (array[float]) – Knot vectors of V0 in each direction.

  • starts (array[float]) – Start indices of splines on current process.

  • values (array[float]) – Return array for spline values S_ijk = S(eta1[i,j,k], eta2[i,j,k], eta3[i,j,k]).

struphy.bsplines.evaluation_kernels_3d.eval_spline_mpi_sparse_meshgrid(eta1: float[:, :, :], eta2: float[:, :, :], eta3: float[:, :, :], _data: float[:, :, :], kind: int[:], pn: int[:], tn1: float[:], tn2: float[:], tn3: float[:], starts: int[:], values: float[:, :, :])[source]#

Sparse meshgrid evaluation of a tensor-product spline, distributed.

Parameters:
  • eta1 (array[float]) – Evaluation points as 3d arrays obtained from sparse meshgrid; points not on local process domain must be flagged as -1. Spline values are obtained as S_ijk = S(eta1[i,0,0], eta2[0,j,0], eta3[0,0,k]).

  • eta2 (array[float]) – Evaluation points as 3d arrays obtained from sparse meshgrid; points not on local process domain must be flagged as -1. Spline values are obtained as S_ijk = S(eta1[i,0,0], eta2[0,j,0], eta3[0,0,k]).

  • eta3 (array[float]) – Evaluation points as 3d arrays obtained from sparse meshgrid; points not on local process domain must be flagged as -1. Spline values are obtained as S_ijk = S(eta1[i,0,0], eta2[0,j,0], eta3[0,0,k]).

  • _data (array[float]) – The spline coefficients c_ijk.

  • kind (array[int]) – Kind of 1d basis in each direction: 0 = N-spline, 1 = D-spline.

  • pn (array[int]) – Spline degrees of V0 in each direction.

  • tn1 (array[float]) – Knot vectors of V0 in each direction.

  • tn2 (array[float]) – Knot vectors of V0 in each direction.

  • tn3 (array[float]) – Knot vectors of V0 in each direction.

  • starts (array[float]) – Start indices of splines on current process.

  • values (array[float]) – Return array for spline values S_ijk = S(eta1[i,0,0], eta2[0,j,0], eta3[0,0,k]).

struphy.bsplines.evaluation_kernels_3d.eval_spline_mpi_markers(markers: float[:, :], _data: float[:, :, :], kind: int[:], pn: int[:], tn1: float[:], tn2: float[:], tn3: float[:], starts: int[:], values: float[:])[source]#

Flat (marker) evaluation of a tensor-product spline, distributed.

Parameters:
  • markers (array[float]) – Marker coordinates in format [Np, 3]; markers not on local process domain must be flagged as -1. Spline values are obtained as S_p = S(*markers[p, :]).

  • _data (array[float]) – The spline coefficients c_ijk.

  • kind (array[int]) – Kind of 1d basis in each direction: 0 = N-spline, 1 = D-spline.

  • pn (array[int]) – Spline degrees of V0 in each direction.

  • tn1 (array[float]) – Knot vectors of V0 in each direction.

  • tn2 (array[float]) – Knot vectors of V0 in each direction.

  • tn3 (array[float]) – Knot vectors of V0 in each direction.

  • starts (array[float]) – Start indices of splines on current process.

  • values (array[float]) – Return 1D array for spline values S_p = S(*markers[p, :]).

struphy.bsplines.evaluation_kernels_3d.get_spans(eta1: float, eta2: float, eta3: float, args_derham: DerhamArguments)[source]#

Compute the knot span index, the N-spline values (in bn) and the D-spline values (in bd) at (eta1, eta2, eta3).

struphy.bsplines.evaluation_kernels_3d.eval_0form_spline_mpi(span1: int, span2: int, span3: int, args_derham: DerhamArguments, form_coeffs: float[:, :, :]) float[source]#

Single-point evaluation of Derham 0-form spline defined by form_coeffs, given N-spline values (in bn) and knot span indices span.

struphy.bsplines.evaluation_kernels_3d.eval_1form_spline_mpi(span1: int, span2: int, span3: int, args_derham: DerhamArguments, form_coeffs_1: float[:, :, :], form_coeffs_2: float[:, :, :], form_coeffs_3: float[:, :, :], out: float[:])[source]#

Single-point evaluation of Derham 1-form spline defined by form_coeffs, given N-spline values (in bn), D-spline values (in bd) and knot span indices span.

struphy.bsplines.evaluation_kernels_3d.eval_2form_spline_mpi(span1: int, span2: int, span3: int, args_derham: DerhamArguments, form_coeffs_1: float[:, :, :], form_coeffs_2: float[:, :, :], form_coeffs_3: float[:, :, :], out: float[:])[source]#

Single-point evaluation of Derham 2-form spline defined by form_coeffs, given N-spline values (in bn), D-spline values (in bd) and knot span indices span.

struphy.bsplines.evaluation_kernels_3d.eval_3form_spline_mpi(span1: int, span2: int, span3: int, args_derham: DerhamArguments, form_coeffs: float[:, :, :]) float[source]#

Single-point evaluation of Derham 0-form spline defined by form_coeffs, given D-spline values (in bd) and knot span indices span.

struphy.bsplines.evaluation_kernels_3d.eval_vectorfield_spline_mpi(span1: int, span2: int, span3: int, args_derham: DerhamArguments, form_coeffs_1: float[:, :, :], form_coeffs_2: float[:, :, :], form_coeffs_3: float[:, :, :], out: float[:])[source]#

Single-point evaluation of vector-field spline (H^1)^3 defined by form_coeffs, given N-spline values (in bn) and knot span indices span.

Acccelerated functions for point-wise evaluation of tensor product B-splines.

S(eta1, eta2) = sum_ij [ c_ij * B_i(eta1) * B_j(eta2) ] with c_ij in R.

Possible combinations for tensor product (BB): * (NN) * (DN) * (ND) * (DD) * (dN/deta N) * (N dN/deta)

struphy.bsplines.evaluation_kernels_2d.evaluation_kernel_2d(p1: int, p2: int, basis1: Final[float[:]], basis2: Final[float[:]], ind1: Final[int[:]], ind2: Final[int[:]], coeff: Final[float[:, :]]) float[source]#

Summing non-zero contributions of a spline function.

Parameters:
  • p1 (int) – Degrees of the univariate splines in each direction.

  • p2 (int) – Degrees of the univariate splines in each direction.

  • basis1 (array[float]) – The p + 1 values of non-zero basis splines at one point (eta1, eta2) in each direction.

  • basis2 (array[float]) – The p + 1 values of non-zero basis splines at one point (eta1, eta2) in each direction.

  • ind1 (array[int]) – Global indices of non-vanishing splines in the element of the considered point.

  • ind2 (array[int]) – Global indices of non-vanishing splines in the element of the considered point.

  • coeff (array[float]) – The spline coefficients c_ij.

Returns:

spline_value – Value of tensor-product spline at point (eta1, eta2).

Return type:

float

struphy.bsplines.evaluation_kernels_2d.evaluate_2d(kind1: int, kind2: int, t1: Final[float[:]], t2: Final[float[:]], p1: int, p2: int, ind1: Final[int[:, :]], ind2: Final[int[:, :]], coeff: Final[float[:, :]], eta1: float, eta2: float) float[source]#

Point-wise evaluation of a tensor-product spline.

Parameters:
  • kind1 (int) – Kind of univariate spline. 1 for B-spline, 2 for M-spline and 3 for derivative of B-spline.

  • kind2 (int) – Kind of univariate spline. 1 for B-spline, 2 for M-spline and 3 for derivative of B-spline.

  • t1 (array[float]) – Knot vectors of univariate splines.

  • t2 (array[float]) – Knot vectors of univariate splines.

  • p1 (int) – Degrees of univariate splines.

  • p2 (int) – Degrees of univariate splines.

  • ind1 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • ind2 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • coeff (array[float]) – The spline coefficients c_ij.

  • eta1 (float) – Point of evaluation.

  • eta2 (float) – Point of evaluation.

Returns:

spline_value – Value of tensor-product spline at point (eta1, eta2).

Return type:

float

struphy.bsplines.evaluation_kernels_2d.evaluate_tensor_product_2d(t1: Final[float[:]], t2: Final[float[:]], p1: int, p2: int, ind1: int[:, :], ind2: int[:, :], coeff: Final[float[:, :]], eta1: Final[float[:]], eta2: Final[float[:]], spline_values: float[:, :], kind: int)[source]#

Tensor-product evaluation of a tensor-product spline.

Parameters:
  • t1 (array[float]) – Knot vectors of univariate splines.

  • t2 (array[float]) – Knot vectors of univariate splines.

  • p1 (int) – Degrees of univariate splines.

  • p2 (int) – Degrees of univariate splines.

  • ind1 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • ind2 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • coeff (array[float]) – The spline coefficients c_ij.

  • eta1 (array[float]) – Points of evaluation in 1d arrays.

  • eta2 (array[float]) – Points of evaluation in 1d arrays.

  • spline_values (array[float]) – Splines evaluated at points S_ij = S(eta1_i, eta2_j).

  • kind (int) –

    Kind of spline to evaluate.
    • 0 : NN

    • 11 : DN

    • 12 : ND

    • 2 : DD

    • 31 : dN/deta N

    • 32 : N dN/deta

    • 41 : ddN/deta^2 N

    • 42 : n ddN/deta^3

    • 43 : dN/deta dN/deta

struphy.bsplines.evaluation_kernels_2d.evaluate_matrix_2d(t1: Final[float[:]], t2: Final[float[:]], p1: int, p2: int, ind1: int[:, :], ind2: int[:, :], coeff: Final[float[:, :]], eta1: Final[float[:, :]], eta2: Final[float[:, :]], spline_values: float[:, :], kind: int)[source]#

General evaluation of a tensor-product spline.

Parameters:
  • t1 (array[float]) – Knot vectors of univariate splines.

  • t2 (array[float]) – Knot vectors of univariate splines.

  • p1 (int) – Degrees of univariate splines.

  • p2 (int) – Degrees of univariate splines.

  • ind1 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • ind2 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • coeff (array[float]) – The spline coefficients c_ij.

  • eta1 (array[float]) – Points of evaluation in 3d arrays such that shape(eta1) == shape(eta2).

  • eta2 (array[float]) – Points of evaluation in 3d arrays such that shape(eta1) == shape(eta2).

  • spline_values (array[float]) – Splines evaluated at points S_ij = S(eta1_i, eta2_j).

  • kind (int) –

    Kind of spline to evaluate.
    • 0 : NN

    • 11 : DN

    • 12 : ND

    • 2 : DD

    • 31 : dN/deta N

    • 32 : N dN/deta

    • 41 : ddN/deta^2 N

    • 42 : n ddN/deta^3

    • 43 : dN/deta dN/deta

struphy.bsplines.evaluation_kernels_2d.evaluate_sparse_2d(t1: Final[float[:]], t2: Final[float[:]], p1: int, p2: int, ind1: int[:, :], ind2: int[:, :], coeff: Final[float[:, :]], eta1: Final[float[:, :]], eta2: Final[float[:, :]], spline_values: float[:, :], kind: int)[source]#

Evaluation of a tensor-product spline using sparse meshgrids.

Parameters:
  • t1 (array[float]) – Knot vectors of univariate splines.

  • t2 (array[float]) – Knot vectors of univariate splines.

  • p1 (int) – Degrees of univariate splines.

  • p2 (int) – Degrees of univariate splines.

  • ind1 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • ind2 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • coeff (array[float]) – The spline coefficients c_ij.

  • eta1 (array[float]) – Points of evaluation in 3d arrays such that shape(eta1) = (:,1), shape(eta2) = (1,:).

  • eta2 (array[float]) – Points of evaluation in 3d arrays such that shape(eta1) = (:,1), shape(eta2) = (1,:).

  • spline_values (array[float]) – Splines evaluated at points S_ij = S(eta1_i, eta2_j).

  • kind (int) –

    Kind of spline to evaluate.
    • 0 : NN

    • 11 : DN

    • 12 : ND

    • 2 : DD

    • 31 : dN/deta N

    • 32 : N dN/deta

    • 41 : ddN/deta^2 N

    • 42 : n ddN/deta^3

    • 43 : dN/deta dN/deta

struphy.bsplines.evaluation_kernels_2d.eval_spline_mpi_2d(p1: int, p2: int, basis1: Final[float[:]], basis2: Final[float[:]], span1: int, span2: int, coeff: Final[float[:, :]], starts: Final[int[:]], pn: Final[int[:]]) float[source]#

Evaluate a spline function on the current process.

Parameters:
  • p1 (int) – Degrees of the univariate splines in each direction.

  • p2 (int) – Degrees of the univariate splines in each direction.

  • basis1 (array[float]) – The p + 1 values of non-zero basis splines at one point (eta1, eta2) in each direction.

  • basis2 (array[float]) – The p + 1 values of non-zero basis splines at one point (eta1, eta2) in each direction.

  • span1 (int) – Particle’s element index in each direction.

  • span2 (int) – Particle’s element index in each direction.

  • coeff (array[float]) – The spline coefficients c_ij of the current process, ie. the _data attribute of a StencilVector.

  • starts (array[int]) – Starting indices of current process.

  • pn (array[int]) – B-spline degrees in each direction (=paddings).

Returns:

spline_value – Value of tensor-product spline at point (eta1, eta2).

Return type:

float

Acccelerated functions for point-wise evaluation of tensor product B-splines.

S(eta1) = sum_i [ c_i * B_i(eta1) ] with c_i in R.

Possible combinations for tensor product (B): * (N) * (D) * (dN/deta)

struphy.bsplines.evaluation_kernels_1d.evaluation_kernel_1d(p1: int, basis1: Final[float[:]], ind1: Final[int[:]], coeff: Final[float[:]]) float[source]#

Summing non-zero contributions.

Parameters:
  • p1 (int) – Degree of the univariate spline.

  • basis1 (array[float]) – The p+1 values of non-zero basis splines at one point (eta1,) from ‘basis_funs’ of shape.

  • ind1 (array[int]) – Global indices of non-vanishing splines in the element of the considered point.

  • coeff (array[float]) – The spline coefficients c_i.

Returns:

spline_value – Value of spline at point (eta1,).

Return type:

float

struphy.bsplines.evaluation_kernels_1d.evaluate(kind1: int, t1: Final[float[:]], p1: int, ind1: Final[int[:, :]], coeff: Final[float[:]], eta1: float) float[source]#

Point-wise evaluation of a spline.

Parameters:
  • kind (int) –

    Kind of spline to evaluate.
    • 0 : N

    • 1 : D

    • 2 : dN/deta

    • 3 : ddN/deta^2

  • t1 (array[float]) – Knot vector of univariate spline.

  • p1 (int) – Degree of univariate spline.

  • ind1 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • coeff (array[float]) – The spline coefficients c_i.

  • eta1 (float) – Point of evaluation.

Returns:

spline_value – Value of spline at point (eta1,).

Return type:

float

struphy.bsplines.evaluation_kernels_1d.evaluate_vector(t1: Final[float[:]], p1: int, ind1: Final[int[:, :]], coeff: Final[float[:]], eta1: float[:], spline_values: float[:], kind: int)[source]#

Vector evaluation of a uni-variate spline.

Parameters:
  • t1 (array[float]) – Knot vector of univariate spline.

  • p1 (int) – Degree of univariate spline.

  • ind1 (array[int]) – Global indices of non-vanishing splines in each element. Can be accessed via (element, local index).

  • coeff (array[float]) – The spline coefficients c_i.

  • eta1 (array[float]) – Points of evaluation in a 1d array.

  • spline_values (array[float]) – Splines evaluated at points S_ij = S(eta1_i, eta2_j).

  • kind (int) –

    Kind of spline to evaluate.
    • 0 : N

    • 1 : D

    • 2 : dN/deta

    • 3 : ddN/deta^2

Shape functions#

Basic functions for point-wise B-spline evaluation

struphy.bsplines.shapefunc_kernels.piecewise(p: int, delta: float, eta: float) float[source]#

evaluate a hat function (B-spline) centered at eta0 (the center of the support) at eta1, i.e.

\[1.0 / delta * S((eta1 - eta0)/ delta)\]

where S is the B-spline of degree p, delta is the cell size of a uniform mesh. Here we use the expression of the B-spline in each cell directly. For the moment, p = 0, 1, or 2.

Parameters:
  • p (int) – degree of the hat function (B-spline)

  • delta (float) – the cell size of a uniform mesh.

  • eta (float) – eta = eta1 - eta0

  • Returns

  • -------- – value of the hat function

struphy.bsplines.shapefunc_kernels.piecewise_der(p: int, delta: float, eta: float) float[source]#

evaluate the derivative of a hat function (B-spline) centered at eta0 (the center of the support) at eta1, i.e. .. math:

1.0 / delta^2 * S'((eta1 - eta0)/ delta)

where S is the B-spline of degree p, delta is the cell size of a uniform mesh. Here we use the expression of the derivative of the B-spline in each cell directly. For the moment, p = 0, 1, or 2.

Parameters:
  • p (int) – degree of the hat function (B-spline)

  • delta (float) – the cell size of a uniform mesh.

  • eta (float) – eta = eta1 - eta0

  • Returns

  • -------- – value of the derivative of the hat function

struphy.bsplines.shapefunc_kernels.convolution(p: int, grids: Final[float[:]], eta: float) float[source]#

evaluate a hat function (B-spline) at eta, i.e.

\[1.0 / delta * S(eta/ delta)\]

where S is the B-spline of degree p. Here we use the definition by convolution of the B-splines.

Parameters:
  • p (int) – degree of the hat function (B-spline)

  • grids (float array) – p + 2 points used in the definition of B-splines.

  • eta (float) – evluation point

  • Returns

  • -------- – value of the hat function at eta

struphy.bsplines.shapefunc_kernels.convolution_der(p: int, grids: Final[float[:]], eta: float) float[source]#

evaluate the derivative of a hat function (B-spline) at eta, i.e.

\[1.0 / delta^2 * S'(eta/ delta)\]

where S is the B-spline of degree p. Here we use the definition by convolution of the B-splines.

Parameters:
  • p (int) – degree of the hat function (B-spline)

  • grids (float array) – p + 2 points used in the definition of B-splines.

  • eta (float) – evluation point

  • Returns

  • -------- – value of the derivative of the hat function at eta