Hierarchical Spline Spaces

This module implements hierarchical spline spaces and supports the hierarchical B-spline (HB-spline) and truncated hierarchical B-spline (THB-spline) bases.

The main user-facing class is HSpace, which describes a hierarchical spline space and supports both HB- and THB-spline representations. Individual functions living in such a spline space are represented by the class HSplineFunc, which follows the interface of BSplineFunc. L2 projection into a hierarchical spline space can be done using the pyiga.approx.project_L2() function.

In order to compute the stiffness matrix and right-hand side vector for the Galerkin discretization of a variational problem in a hierarchical spline space, use pyiga.assemble.assemble(). Internally, this uses the HDiscretization class.

A tensor product B-spline basis function is usually referred to by a multi-index represented as a tuple (i_1, …, i_d), where d is the space dimension and i_k is the index of the univariate B-spline function used in the k-th coordinate direction. Similarly, cells in the underlying tensor product mesh are indexed by multi-indices (j_1, .., j_d), where j_k is the index of the knot span along the k-th axis.

Whenever an ordering of the degrees of freedom in a hierarchical spline space is required, for instance when assembling a stiffness matrix, we use the following canonical order: first, all active basis function on the coarsest level, then all active basis functions on the next finer level, and so on until the finest level. Within each level, the functions are ordered lexicographically with respect to their tensor product multi-indices (i_1, …, i_d).

The canonical order on active cells is defined in the same way.

The implementation of HSpace is loosely based on the approach described in [GV2018] and the corresponding implementation in [GeoPDEs].

[GV2018]Garau, Vázquez: “Algorithms for the implementation of adaptive isogeometric methods using hierarchical B-splines”, 2018.
[GeoPDEs]http://rafavzqz.github.io/geopdes/

Hierarchical spline spaces

class pyiga.hierarchical.HSpace(kvs, truncate=False, disparity=inf, bdspecs=None)

Represents a HB-spline or THB-spline space over an adaptively refined mesh.

Parameters:
  • kvs – a sequence of d KnotVector instances, representing the tensor product B-spline space on the coarsest level
  • truncate (bool) – if True, the space represents a THB-spline space, otherwise an HB-spline space
  • disparity (int) –

    the desired mesh level disparity. This means that an active basis function on level lv may have interactions at most with coarse functions from level lv - disparity, but not from any coarser levels.

    This disparity is respected when calling refine(). Lower disparity leads to more gradual changes in the sizes of neighboring active cells.

    If no restriction on the number of overlapping mesh levels is desired, pass numpy.inf (which is the default).

  • bdspecs – optionally, a list of boundary specifications on which degrees of freedom should be eliminated (usually for treating Dirichlet boundary conditions). See assemble.compute_dirichlet_bc() for the precise format.
active_cells(lv=None, flat=False)

If lv is specified, return the set of active cells on that level. Otherwise, return a list containing, for each level, the set of active cells.

If lv=None and flat=True, return a flat list of (lv, (j_1, …, j_d)) pairs of all active cells in canonical order, where the first entry is the level and the second entry is the multi-index of the cell on that level. The length of the returned list is total_active_cells.

active_functions(lv=None, flat=False)

If lv is specified, return the set of active functions on that level. Otherwise, return a list containing, for each level, the set of active functions.

If lv=None and flat=True, return a flat list of (lv, (i_1, …, i_d)) pairs of all active functions in canonical order, where the first entry is the level and the second entry is the multi-index of the function on that level. The length of the returned list is numdofs.

active_indices()

Return a tuple which contains, per level, the raveled (sequential) indices of active basis functions.

cell_extents(lv, c)

Return the extents (as a tuple of min/max pairs) of the cell c on level lv.

cell_supp_indices(remove_dirichlet=True)

Return a tuple which contains tuples which contain, per level, the raveled (sequential) indices of cell_supp basis functions in the virtual hierarchy per level.

coeffs_to_levelwise_funcs(coeffs, truncate=None)

Compute the levelwise contributions of the hierarchical spline function given by the coefficient vector coeffs as a list containing one BSplineFunc per level.

If truncate=True, the coefficients are interpreted as THB-spline coefficients; if False, as HB-splines; if None (default), HSpace.truncate is used.

compute_supports(functions)

Compute the union of the supports, in terms of active mesh cells, of the given list-of-seqs of functions and return the active cells as a dict-of-sets.

deactivated_cells(lv=None)

If lv is specified, return the set of deactivated cells on that level. Otherwise, return a list containing, for each level, the set of deactivated cells.

deactivated_indices()

Return a tuple which contains, per level, the raveled (sequential) indices of deactivated basis functions.

dirichlet_dofs(lv=None)

Matrix indices which lie on the specified Dirichlet boundaries.

func_supp_indices()

Return a tuple which contains tuples which contain, per level, the raveled (sequential) indices of func_supp basis functions in the virtual hierarchy per level.

function_support(lv, jj)

Return the support (as a tuple of pairs) of the function on level lv with multi-index jj.

global_indices()

Return a tuple which contains tuples which contain, per level, the raveled (sequential) indices of global basis functions in the virtual hierarchy per level.

grid_eval(coeffs, gridaxes, truncate=None)

Evaluate an HB-spline function with the given coefficients over a tensor product grid.

hb_to_thb()

Return a sparse square matrix of size numdofs which transforms HB-spline coefficients into the corresponding THB-spline coefficients.

incidence_matrix()

Compute the incidence matrix which contains one row per active basis function and one column per active cell in the hierarchical mesh. An entry (i,j) is 1 if the function i is nonzero in the cell j, and 0 otherwise.

knotvectors(lv)

Return a tuple of KnotVector instances describing the tensor product space on level lv.

mesh(lv)

Return the underlying TPMesh on the given level.

new_indices()

Return a tuple which contains tuples which contain, per level, the raveled (sequential) indices of newly added basis functions in the virtual hierarchy per level.

non_dirichlet_dofs()

Matrix indices which do not lie on the specified Dirichlet boundaries.

numactive

A tuple containing the number of active basis functions per level.

numdofs

The total number of active basis functions in this hierarchical space.

numlevels

The number of levels in this hierarchical space.

ravel_indices(indices)

Given a list indices which contains, per level, a list or set of function multi-indices on that level, return a list of arrays with the corresponding raveled indices, i.e., the sequential indices corresponding to the multi-indices in lexicographic order.

These raveled indices are useful, e.g., for indexing into vectors or matrices which are defined over the full tensor product space on a level.

See numpy.ravel_multi_index() for details on the raveling operation.

refine(marked, truncate=False)

Refine the given cells; marked is a dictionary which has the levels as indices and the list of marked cells on that level as values.

The refinement procedure preserves the mesh level disparity, following the method described in [Bracco, Giannelli, Vázquez, 2018].

Returns:the actually refined cells in the same format as marked; if disparity is less than infinity, this is a superset of the input cells
refine_region(lv, region_function)

Refine all active cells on level lv whose cell center satisfies region_function.

region_function should be a function of dim scalar arguments (e.g., (x,y)) which returns True if the point is within the refinement region.

represent_fine(lv=None, truncate=None, rows=None, restrict=False)

Compute a matrix which represents HB- or THB-spline basis functions in terms of their coefficients in the finest tensor product spline space.

By default, all active HB-spline functions are represented on the finest tensor product mesh. The returned matrix has size N_fine × N_act, where N_fine is the number of degrees of freedom in the tensor product mesh of the finest level and N_act = numdofs is the total number of active basis functions.

If lv is specified, only active functions up to level lv are represented in terms of their coefficients on level lv.

If truncate is True, the representation of the THB-spline (truncated) basis functions is computed instead of that of the HB-splines. If truncate is None (default), the attribute HSpace.truncate is used.

If rows is given, only those rows are kept in the output. In other words, only the representation with respect to these tensor product basis functions on the fine level is computed. If restrict=False, then the shape of the matrix is not changed, but only the corresponding rows are filled. If restrict=True, the matrix is restricted to the submatrix having only these rows.

split_coeffs(x)

Given a coefficient vector x of length numdofs in canonical order, split it into numlevels vectors which contain the contributions from each individual level.

thb_to_hb()

Return a sparse square matrix of size numdofs which transforms THB-spline coefficients into the corresponding HB-spline coefficients.

total_active_cells

The total number of active cells in the hierarchical mesh.

tp_prolongation(lv, kron=False)

Return the prolongation operator for the underlying tensor product mesh from level lv to lv+1.

If kron is True, the prolongation is returned as a sparse matrix. Otherwise, the prolongation is returned as a tuple of sparse matrices, one per space dimension, whose Kronecker product represents the prolongation operator.

Note

This method, particularly with kron=True, is inherently inefficient since it deals with the full tensor product spaces, not merely the active basis functions.

trunc_indices()

Return a tuple which contains tuples which contain, per level, the raveled (sequential) indices of trunc basis functions in the virtual hierarchy per level.

truncate_one_level(k, num_rows=None, inverse=False)

Compute the matrix which realizes truncation from level k to k+1.

class pyiga.hierarchical.HSplineFunc(hspace, u, truncate=None)

A function living in a hierarchical spline space.

Parameters:
  • hspace (HSpace) – the hierarchical spline space
  • u (array) – the vector of coefficients corresponding to the active basis functions, with length HSpace.numdofs, in canonical order
  • truncate (bool) – if True, the coefficients are interpreted as THB-spline coefficients; if False, as HB-spline coefficients. If it is None (default), the attribute HSpace.truncate of hspace is used.
grid_eval(gridaxes)

Evaluate the function on a tensor product grid.

See BSplineFunc.grid_eval() for details.

grid_hessian(gridaxes)

Evaluate the Hessian matrix on a tensor product grid.

See BSplineFunc.grid_hessian() for details.

grid_jacobian(gridaxes)

Evaluate the Jacobian on a tensor product grid.

See BSplineFunc.grid_jacobian() for details.

support

Return a sequence of pairs (lower,upper), one per source dimension, which describe the extent of the support in the parameter space.

class pyiga.hierarchical.TPMesh(kvs)

A tensor product mesh described by a sequence of knot vectors.

cell_extents(c)

Return the extents (as a tuple of min/max pairs) of the cell c.

cells()

Return a list of all cells in this mesh.

functions()

Return a list of all basis functions defined on this mesh.

neighbors(indices)

Return all functions which have nontrivial support intersection with the given ones.

support(indices)

Return the set of cells where any of the given functions does not vanish.

supported_in(cells)

Return the set of functions whose support intersects the given cells.

Discretization

class pyiga.hierarchical.HDiscretization(hspace, vform, asm_args)

Represents the discretization of a variational problem over a hierarchical spline space.

Parameters:
  • hspace (HSpace) – the HB- or THB-spline space in which to discretize
  • vform (VForm) – the variational form describing the problem
  • asm_args (dict) –

    a dictionary which provides named inputs for the assembler. Most problems will require at least a geometry map; this can be given in the form {'geo': geo}, where geo is a geometry function defined using the pyiga.geometry module. Further inputs declared via the VForm.input() method must be included in this dict.

    The assemblers both for the matrix and any linear functionals will draw their input arguments from this dict.

assemble_functional(vf)

Assemble a linear functional described by the VForm vf over the hierarchical spline space.

Returns:a vector whose length is equal to the HSpace.numdofs attribute of hspace, corresponding to the active basis functions in canonical order
assemble_matrix()

Assemble the stiffness matrix for the hierarchical discretization and return it.

Returns:a sparse matrix whose size corresponds to the HSpace.numdofs attribute of hspace
assemble_rhs(vf=None)

Assemble the right-hand side vector for the hierarchical discretization and return it.

By default (if vf=None), a standard L2 inner product <f, v> is used for computing the right-hand side, and the function f is taken from the key 'f' of the asm_args dict. It is assumed to be given in physical coordinates.

A different functional can be specified by passing a VForm with arity=1 as the vf parameter.

Returns:a vector whose length is equal to the HSpace.numdofs attribute of hspace, corresponding to the active basis functions in canonical order