sfepy.discrete.iga.iga module

Isogeometric analysis utilities.

Notes

The functions compute_bezier_extraction_1d() and eval_nurbs_basis_tp() implement the algorithms described in [1].

[1] Michael J. Borden, Michael A. Scott, John A. Evans, Thomas J. R. Hughes:

Isogeometric finite element data structures based on Bezier extraction of NURBS, Institute for Computational Engineering and Sciences, The University of Texas at Austin, Austin, Texas, March 2010.

sfepy.discrete.iga.iga.combine_bezier_extraction(cs)[source]

For a nD B-spline parametric domain, combine the 1D element extraction operators in each parametric dimension into a single operator for each nD element.

Parameters:
cslist of lists of 2D arrays

The element extraction operators in each parametric dimension.

Returns:
ccslist of 2D arrays

The combined element extraction operators.

sfepy.discrete.iga.iga.compute_bezier_control(control_points, weights, ccs, conn, bconn)[source]

Compute the control points and weights of the Bezier mesh.

Parameters:
control_pointsarray

The NURBS control points.

weightsarray

The NURBS weights.

ccslist of 2D arrays

The combined element extraction operators.

connarray

The connectivity of the global NURBS basis.

bconnarray

The connectivity of the Bezier basis.

Returns:
bezier_control_pointsarray

The control points of the Bezier mesh.

bezier_weightsarray

The weights of the Bezier mesh.

sfepy.discrete.iga.iga.compute_bezier_extraction(knots, degrees)[source]

Compute local (element) Bezier extraction operators for a nD B-spline parametric domain.

Parameters:
knotssequence of array or array

The knot vectors.

degreessequence of ints or int

Polynomial degrees in each parametric dimension.

Returns:
cslist of lists of 2D arrays

The element extraction operators in each parametric dimension.

sfepy.discrete.iga.iga.compute_bezier_extraction_1d(knots, degree)[source]

Compute local (element) Bezier extraction operators for a 1D B-spline parametric domain.

Parameters:
knotsarray

The knot vector.

degreeint

The curve degree.

Returns:
csarray of 2D arrays (3D array)

The element extraction operators.

sfepy.discrete.iga.iga.create_boundary_qp(coors, dim)[source]

Create boundary quadrature points from the surface quadrature points.

Uses the Bezier element tensor product structure.

Parameters:
coorsarray, shape (n_qp, d)

The coordinates of the surface quadrature points.

dimint

The topological dimension.

Returns:
bcoorsarray, shape (n_qp, d + 1)

The coordinates of the boundary quadrature points.

sfepy.discrete.iga.iga.create_connectivity(n_els, knots, degrees)[source]

Create connectivity arrays of nD Bezier elements.

Parameters:
n_elssequence of ints

The number of elements in each parametric dimension.

knotssequence of array or array

The knot vectors.

degreessequence of ints or int

The basis degrees in each parametric dimension.

Returns:
connarray

The connectivity of the global NURBS basis.

bconnarray

The connectivity of the Bezier basis.

sfepy.discrete.iga.iga.create_connectivity_1d(n_el, knots, degree)[source]

Create connectivity arrays of 1D Bezier elements.

Parameters:
n_elint

The number of elements.

knotsarray

The knot vector.

degreeint

The basis degree.

Returns:
connarray

The connectivity of the global NURBS basis.

bconnarray

The connectivity of the Bezier basis.

sfepy.discrete.iga.iga.eval_bernstein_basis(x, degree)[source]

Evaluate the Bernstein polynomial basis of the given degree, and its derivatives, in a point x in [0, 1].

Parameters:
xfloat

The point in [0, 1].

degreeint

The basis degree.

Returns:
funsarray

The degree + 1 values of the Bernstein polynomial basis.

dersarray

The degree + 1 values of the Bernstein polynomial basis derivatives.

sfepy.discrete.iga.iga.eval_mapping_data_in_qp(qps, control_points, weights, degrees, cs, conn, cells=None)[source]

Evaluate data required for the isogeometric domain reference mapping in the given quadrature points. The quadrature points are the same for all Bezier elements and should correspond to the Bernstein basis degree.

Parameters:
qpsarray

The quadrature points coordinates with components in [0, 1] reference element domain.

control_pointsarray

The NURBS control points.

weightsarray

The NURBS weights.

degreessequence of ints or int

The basis degrees in each parametric dimension.

cslist of lists of 2D arrays

The element extraction operators in each parametric dimension.

connarray

The connectivity of the global NURBS basis.

cellsarray, optional

If given, use only the given Bezier elements.

Returns:
bfsarray

The NURBS shape functions in the physical quadrature points of all elements.

bfgsarray

The NURBS shape functions derivatives w.r.t. the physical coordinates in the physical quadrature points of all elements.

detsarray

The Jacobians of the mapping to the unit reference element in the physical quadrature points of all elements.

sfepy.discrete.iga.iga.eval_nurbs_basis_tp(qp, ie, control_points, weights, degrees, cs, conn)[source]

Evaluate the tensor-product NURBS shape functions in a quadrature point for a given Bezier element.

Parameters:
qparray

The quadrature point coordinates with components in [0, 1] reference element domain.

ieint

The Bezier element index.

control_pointsarray

The NURBS control points.

weightsarray

The NURBS weights.

degreessequence of ints or int

The basis degrees in each parametric dimension.

cslist of lists of 2D arrays

The element extraction operators in each parametric dimension.

connarray

The connectivity of the global NURBS basis.

Returns:
Rarray

The NURBS shape functions.

dR_dxarray

The NURBS shape functions derivatives w.r.t. the physical coordinates.

detarray

The Jacobian of the mapping to the unit reference element.

sfepy.discrete.iga.iga.eval_variable_in_qp(variable, qps, control_points, weights, degrees, cs, conn, cells=None)[source]

Evaluate a field variable in the given quadrature points. The quadrature points are the same for all Bezier elements and should correspond to the Bernstein basis degree. The field variable is defined by its DOFs - the coefficients of the NURBS basis.

Parameters:
variablearray

The DOF values of the variable with n_c components, shape (:, n_c).

qpsarray

The quadrature points coordinates with components in [0, 1] reference element domain.

control_pointsarray

The NURBS control points.

weightsarray

The NURBS weights.

degreessequence of ints or int

The basis degrees in each parametric dimension.

cslist of lists of 2D arrays

The element extraction operators in each parametric dimension.

connarray

The connectivity of the global NURBS basis.

cellsarray, optional

If given, use only the given Bezier elements.

Returns:
coorsarray

The physical coordinates of the quadrature points of all elements.

valsarray

The field variable values in the physical quadrature points.

detsarray

The Jacobians of the mapping to the unit reference element in the physical quadrature points.

sfepy.discrete.iga.iga.get_bezier_element_entities(degrees)[source]

Get faces and edges of a Bezier mesh element in terms of indices into the element’s connectivity (reference Bezier element entities).

Parameters:
degreessequence of ints or int

Polynomial degrees in each parametric dimension.

Returns:
faceslist of arrays

The indices for each face or None if not 3D.

edgeslist of arrays

The indices for each edge or None if not at least 2D.

verticeslist of arrays

The indices for each vertex.

Notes

The ordering of faces and edges has to be the same as in sfepy.discrete.fem.geometry_element.geometry_data.

sfepy.discrete.iga.iga.get_bezier_topology(bconn, degrees)[source]

Get a topology connectivity corresponding to the Bezier mesh connectivity.

In the referenced Bezier control points the Bezier mesh is interpolatory.

Parameters:
bconnarray

The connectivity of the Bezier basis.

degreessequence of ints or int

The basis degrees in each parametric dimension.

Returns:
tconnarray

The topology connectivity (corner nodes, or vertices, of Bezier elements) with vertex ordering suitable for a FE mesh.

sfepy.discrete.iga.iga.get_facet_axes(dim)[source]

For each reference Bezier element facet return the facet axes followed by the remaining (perpendicular) axis, as well as the remaining axis coordinate of the facet.

Parameters:
dimint

The topological dimension.

Returns:
axesarray

The axes of the reference element facets.

coorsarray

The remaining coordinate of the reference element facets.

sfepy.discrete.iga.iga.get_patch_box_regions(n_els, degrees)[source]

Get box regions of Bezier topological mesh in terms of element corner vertices of Bezier mesh.

Parameters:
n_elssequence of ints

The number of elements in each parametric dimension.

degreessequence of ints or int

Polynomial degrees in each parametric dimension.

Returns:
regionsdict

The Bezier mesh vertices of box regions.

sfepy.discrete.iga.iga.get_raveled_index(indices, shape)[source]

Get a global raveled index corresponding to nD indices into an array of the given shape.

sfepy.discrete.iga.iga.get_surface_degrees(degrees)[source]

Get degrees of the NURBS patch surfaces.

Parameters:
degreessequence of ints or int

Polynomial degrees in each parametric dimension.

Returns:
sdegreeslist of arrays

The degrees of the patch surfaces, in the order of the reference Bezier element facets.

sfepy.discrete.iga.iga.get_unraveled_indices(index, shape)[source]

Get nD indices into an array of the given shape corresponding to a global raveled index.

sfepy.discrete.iga.iga.tensor_product(a, b)[source]

Compute tensor product of two 2D arrays with possibly different shapes. The result has the form:

c = [[a00 b, a01 b, ...],
     [a10 b, a11 b, ...],
      ...
      ...               ]