fem.splines#

Classes#

Inheritance diagram of psydac.fem.splines

SplineSpace(degree[, knots, grid, ...])

a 1D Splines Finite Element space

Details#

class SplineSpace(degree, knots=None, grid=None, multiplicity=None, parent_multiplicity=None, periodic=False, dirichlet=(False, False), basis='B', pads=None)[source]#

Bases: FemSpace

a 1D Splines Finite Element space

Parameters:
degreeint

Polynomial degree.

knotsarray_like

Coordinates of knots (clamped or extended by periodicity).

grid: array_like

Coordinates of the grid. Used to construct the knots sequence, if not given.

multiplicity: int

Multiplicity of the knots in the knot sequence.

parent_multiplicity: int

Multiplicity of the parent knot sequence, if the space is reduced space.

periodicbool

True if domain is periodic, False otherwise. Default: False

dirichlettuple, list

True if using homogeneous dirichlet boundary conditions, False otherwise. Must be specified for each bound Default: (False, False)

basisstr

Set to “B” for B-splines (have partition of unity) Set to “M” for M-splines (have unit integrals)

property histopolation_grid#

Coordinates of the N+1 points x[i] that define the N 1D edges (x[i], x[i+1]) for histopolation, where N is equal to the number of basis functions (i.e. the cardinality of the space).

In the non-periodic case x is simply the array of extended Greville points. In the periodic case we “unroll” the 1D edges to ensure that they correspond to positive, well-defined intervals with x[i] < x[i+1].

init_interpolation(dtype=<class 'float'>)[source]#

Compute the 1D collocation matrix and factorize it, in preparation for the calculation of a spline interpolant given the values at the Greville points.

init_histopolation(dtype=<class 'float'>)[source]#

Compute the 1D histopolation matrix and factorize it, in preparation for the calculation of a spline interpolant given the integrals within the cells defined by the extended Greville points.

property ldim#

Parametric dimension.

property periodic#

True if domain is periodic, False otherwise.

property pads#

Padding for potential parallel assembly.

property mapping#

Assume identity mapping for now.

property vector_space#

Returns the topological associated vector space.

property is_product#

Boolean flag that describes whether the space is a product space. If True, an element of this space can be decomposed into separate fields.

property symbolic_space#

Symbolic space.

eval_field(field, *eta, weights=None)[source]#

Evaluate field at location(s) eta.

Parameters:
fieldFemField

Field object (element of FemSpace) to be evaluated.

etalist of float or numpy.ndarray

Evaluation point(s) in logical domain.

weightsStencilVector, optional

Weights of the basis functions, such that weights.space == field.coeffs.space.

Returns:
valuefloat or numpy.ndarray

Field value(s) at location(s) eta.

eval_field_gradient(field, *eta, weights=None)[source]#

Evaluate field gradient at location(s) eta.

Parameters:
fieldFemField

Field object (element of FemSpace) to be evaluated.

etalist of float or numpy.ndarray

Evaluation point(s) in logical domain.

weightsStencilVector, optional

Weights of the basis functions, such that weights.space == field.coeffs.space.

Returns:
valuefloat or numpy.ndarray

Value(s) of field gradient at location(s) eta.

property is_scalar#

Only scalar field is implemented for now.

property basis#
property interpolation_grid#
property nbasis#

Number of basis functions, i.e. cardinality of spline space.

property degree#

Spline degree.

property ncells#

Number of cells in domain.

property dirichlet#

True if using homogeneous dirichlet boundary conditions, False otherwise.

property knots#

Knot sequence.

property multiplicity#
property parent_multiplicity#
property breaks#

List of breakpoints.

property domain#

Domain boundaries [a,b].

property greville#

Coordinates of all Greville points. Used for interpolation.

property ext_greville#

Greville coordinates of ‘extended’ space with degree p+1. Used for histopolation.

property scaling_array#

If self.basis==’M’, return array used to rescale B-splines to M-splines If self.basis==’B’, return None.

The length of the scaling array is (len(knots)-degree-1).

compute_interpolant(values, field)[source]#

Compute field (i.e. update its spline coefficients) such that it interpolates a certain function \(f(x)\) at the Greville points.

Parameters:
valuesarray_like (nbasis,)

Function values \(f(x_i)\) at the ‘nbasis’ Greville points \(x_i\), to be interpolated.

fieldFemField

Input/output argument: spline that has to interpolate the given values.

compute_histopolant(values, field)[source]#

Compute field (i.e. update its spline coefficients) such that its integrals between the extended Greville points match the given values.

Parameters:
valuesarray_like (nbasis,)

Integral values between the ‘nbasis’ extended Greville cells \([x_i, x_{i+1}]\), to be matched by the spline.

fieldFemField

Input/output argument: spline that has to match the given integral values.

refine(ncells)[source]#

Create a refined 1D spline space with the given number of cells.

Parameters:
ncellsint

Number of cells of refined space. Must be multiple of self.ncells.

Returns:
SplineSpace

Refined 1D spline space which contains the original space.

draw()[source]#