UFL user-facing API#

UFL: The Unified Form Language.

The Unified Form Language is an embedded domain specific language for definition of variational forms intended for finite element discretization. More precisely, it defines a fixed interface for choosing finite element spaces and defining expressions for weak forms in a notation close to the mathematical one.

This Python module contains the language as well as algorithms to work with it.

  • To import the language, type:

    import ufl
    
  • To import the underlying classes an UFL expression tree is built from, type

    import ufl.classes
    
  • Various algorithms for working with UFL expression trees can be accessed by

Classes and algorithms are considered implementation details and should not be used in form definitions.

For more details on the language, see

and

The development version can be found in the repository at

A very brief overview of the language contents follows:

  • Cells:

    -AbstractCell
    -Cell
    -TensorProductCell
    -vertex
    -interval
    -triangle
    -tetrahedron
    -quadrilateral
    -hexahedron
    -prism
    -pyramid
    -pentatope
    -tesseract
    
  • Domains:

    -AbstractDomain
    -Mesh
    -MeshSequence
    -MeshView
    
  • Sobolev spaces:

    -L2
    -H1
    -H2
    -HInf
    -HDiv
    -HCurl
    -HEin
    -HDivDiv
    -HCurlDiv
    
  • Pull backs:

    -identity_pullback
    -contravariant_piola
    -covariant_piola
    -l2_piola
    -double_contravariant_piola
    -double_covariant_piola
    -covariant_contravariant_piola
    
  • Function spaces:

    -FunctionSpace
    -MixedFunctionSpace
    
  • Arguments:

    -Argument
    -TestFunction
    -TrialFunction
    -Arguments
    -TestFunctions
    -TrialFunctions
    
  • Coefficients:

    -Coefficient
    -Constant
    -VectorConstant
    -TensorConstant
    
  • Splitting form arguments in mixed spaces:

    -split
    
  • Literal constants:

    -Identity
    -PermutationSymbol
    
  • Geometric quantities:

    -SpatialCoordinate
    -FacetNormal
    -CellNormal
    -CellVolume
    -CellDiameter
    -Circumradius
    -MinCellEdgeLength
    -MaxCellEdgeLength
    -FacetArea
    -MinFacetEdgeLength
    -MaxFacetEdgeLength
    -Jacobian
    -JacobianDeterminant
    -JacobianInverse
    
  • Indices:

    -Index
    -indices
    -i, j, k, l
    -p, q, r, s
    
  • Scalar to tensor expression conversion:

    -as_tensor
    -as_vector
    -as_matrix
    
  • Unit vectors and matrices:

    -unit_vector
    -unit_vectors
    -unit_matrix
    -unit_matrices
    
  • Tensor algebra operators:

    -outer, inner, dot, cross, perp
    -det, inv, cofac
    -transpose, tr, diag, diag_vector
    -dev, skew, sym
    
  • Elementwise tensor operators:

    -elem_mult
    -elem_div
    -elem_pow
    -elem_op
    
  • Differential operators:

    -variable
    (-diff,)
    -grad, nabla_grad
    -div, nabla_div
    -curl, rot
    -Dx, Dn
    
  • Nonlinear functions:

    -max_value, min_value
    -abs, sign
    -sqrt
    -exp, ln, erf
    -cos, sin, tan
    -acos, asin, atan, atan2
    -cosh, sinh, tanh
    -bessel_J, bessel_Y, bessel_I, bessel_K
    
  • Complex operations:

    - conj, real, imag
    conjugate is an alias for conj
    
  • Discontinuous Galerkin operators:

    -v("+"), v("-")
    -jump
    -avg
    -cell_avg, facet_avg
    
  • Conditional operators:

    - eq, ne, le, ge, lt, gt
    - <, >, <=, >=
    - And, Or, Not
    - conditional
    
  • Integral measures:

    -dx, ds, dS, dP, dl
    -dc, dC, dO, dI, dX
    -ds_b, ds_t, ds_tb, ds_v, dS_h, dS_v
    
  • Form transformations:

    -rhs, lhs
    -system
    -functional
    -replace
    -adjoint
    -action
    (-energy_norm,)
    -sensitivity_rhs
    -derivative
    
class ufl.AbstractCell[source]#

Bases: UFLObject

A base class for all cells.

abstract property cellname: str#

Return the cellname of the cell.

property edge_types: tuple[AbstractCell, ...]#

Get the unique edge types.

property edges: tuple[AbstractCell, ...]#

Get the edges.

property face_types: tuple[AbstractCell, ...]#

Get the unique face types.

property faces: tuple[AbstractCell, ...]#

Get the faces.

property facet_types: tuple[AbstractCell, ...]#

Get the unique facet types.

Facets are entities of dimension tdim-1.

property facets: tuple[AbstractCell, ...]#

Get the facets.

Facets are entities of dimension tdim-1.

abstract property has_simplex_facets: bool#

Return True if all the facets of this cell are simplex cells.

abstract property is_simplex: bool#

Return True if this is a simplex cell.

property num_edges: int#

Get the number of edges.

property num_faces: int#

Get the number of faces.

property num_facets: int#

Get the number of facets.

Facets are entities of dimension tdim-1.

property num_peaks: int#

Get the number of peaks.

Peaks are entities of dimension tdim-3.

property num_ridges: int#

Get the number of ridges.

Ridges are entities of dimension tdim-2.

abstractmethod num_sub_entities(dim: int) int[source]#

Get the number of sub-entities of the given dimension.

property num_vertices: int#

Get the number of vertices.

property peak_types: tuple[AbstractCell, ...]#

Get the unique peak types.

Peaks are entities of dimension tdim-3.

property peaks: tuple[AbstractCell, ...]#

Get the peaks.

Peaks are entities of dimension tdim-3.

abstractmethod reconstruct(**kwargs: Any) AbstractCell[source]#

Reconstruct this cell, overwriting properties by those in kwargs.

property ridge_types: tuple[AbstractCell, ...]#

Get the unique ridge types.

Ridges are entities of dimension tdim-2.

property ridges: tuple[AbstractCell, ...]#

Get the ridges.

Ridges are entities of dimension tdim-2.

abstractmethod sub_entities(dim: int) tuple[AbstractCell, ...][source]#

Get the sub-entities of the given dimension.

abstractmethod sub_entity_types(dim: int) tuple[AbstractCell, ...][source]#

Get the unique sub-entity types of the given dimension.

abstract property topological_dimension: int#

Return the dimension of the topology of this cell.

property vertex_types: tuple[AbstractCell, ...]#

Get the unique vertices types.

property vertices: tuple[AbstractCell, ...]#

Get the vertices.

class ufl.AbstractDomain(topological_dimension, geometric_dimension)[source]#

Bases: object

Symbolic representation of a geometric domain.

Domain has only a geometric and a topological dimension.

Initialise.

can_make_function_space(element: AbstractFiniteElement) bool[source]#

Check whether this mesh can make a function space with element.

property geometric_dimension#

Return the dimension of the space this domain is embedded in.

iterable_like(element: AbstractFiniteElement) Iterable[Mesh] | MeshSequence[source]#

Return iterable object that is iterable like element.

property meshes#

Return the component meshes.

property topological_dimension#

Return the dimension of the topology of this domain.

class ufl.AbstractFiniteElement[source]#

Bases: ABC

Base class for all finite elements.

To make your element library compatible with UFL, you should make a subclass of AbstractFiniteElement and provide implementions of all the abstract methods and properties. All methods and properties that are not marked as abstract are implemented here and should not need to be overwritten in your subclass.

An example of how the methods in your subclass could be implemented can be found in Basix; see FEniCS/basix

abstract property cell: Cell#

Return the cell of the finite element.

abstract property embedded_subdegree: int#

Degree of the maximum degree Lagrange space that is spanned by this element.

This returns the degree of the highest degree Lagrange space such that the polynomial space of the Lagrange space is a subspace of this element’s polynomial space. If this element’s polynomial space does not include the constant function, this function should return -1.

Note that on a simplex cells, the polynomial space of Lagrange space is a complete polynomial space, but on other cells this is not true. For example, on quadrilateral cells, the degree 1 Lagrange space includes the degree 2 polynomial xy.

abstract property embedded_superdegree: int | None#

Degree of the minimum degree Lagrange space that spans this element.

This returns the degree of the lowest degree Lagrange space such that the polynomial space of the Lagrange space is a superspace of this element’s polynomial space. If this element contains basis functions that are not in any Lagrange space, this function should return None.

Note that on a simplex cells, the polynomial space of Lagrange space is a complete polynomial space, but on other cells this is not true. For example, on quadrilateral cells, the degree 1 Lagrange space includes the degree 2 polynomial xy.

is_cellwise_constant() bool[source]#

Check whether this element is spatially constant over each cell.

property num_sub_elements: int#

Return number of sub-elements.

This function does not recurse: ie it does not count the sub-elements of sub-elements.

abstract property pullback: AbstractPullback#

Return the pullback for this element.

abstract property reference_value_shape: tuple[int, ...]#

Return the shape of the value space on the reference cell.

property reference_value_size: int#

Return the integer product of the reference value shape.

abstract property sobolev_space: SobolevSpace#

Return the underlying Sobolev space.

abstract property sub_elements: Sequence#

Return list of sub-elements.

This function does not recurse: ie it does not extract the sub-elements of sub-elements.

class ufl.AbstractPullback[source]#

Bases: ABC

An abstract pull back.

apply(expr: Expr, domain: AbstractDomain | None = None) Expr[source]#

Apply the pull back.

Parameters:
  • expr – A function on a physical cell

  • domain – The domain on which the function is defined

Returns: The function pulled back to the reference cell

abstract property is_identity: bool#

Is this pull back the identity (or the identity applied to mutliple components).

abstractmethod physical_value_shape(element, domain) tuple[int, ...][source]#

Get the physical value shape when this pull back is applied to an element on a domain.

Parameters:
  • element – The element that the pull back is applied to

  • domain – The domain

Returns:

The value shape when the pull back is applied to the given element

class ufl.Action(*args, **kw)[source]#

Bases: BaseForm

UFL base form type: respresents the action of an object on another.

For example:

res = Ax

A would be the first argument, left and x would be the second argument, right.

Action objects will result when the action of an assembled object (e.g. a Matrix) is taken. This delays the evaluation of the action until assembly occurs.

Initialise.

equals(other)[source]#

Check if two Actions are equal.

left()[source]#

Get left.

right()[source]#

Get right.

ufl_domains()[source]#

Return all domains found in the base form.

ufl_function_spaces()[source]#

Get the tuple of function spaces of the underlying form.

ufl_operands: tuple[FormArgument, ...]#
class ufl.Adjoint(*args, **kw)[source]#

Bases: BaseForm

UFL base form type: represents the adjoint of an object.

Adjoint objects will result when the adjoint of an assembled object (e.g. a Matrix) is taken. This delays the evaluation of the adjoint until assembly occurs.

Initialise.

equals(other)[source]#

Check if two Adjoints are equal.

form()[source]#

Return the form.

ufl_domains()[source]#

Return all domains found in the base form.

ufl_function_spaces()[source]#

Get the tuple of function spaces of the underlying form.

ufl_operands: tuple[FormArgument, ...]#
ufl.And(left, right)[source]#

A boolean expression (left and right) for use with conditional.

class ufl.Argument(*args, **kw)[source]#

Bases: FormArgument, BaseArgument

UFL value: Representation of an argument to a form.

Initialise.

ufl_domains()[source]#

Return UFL domains.

ufl.Arguments(function_space, number)[source]#

Create an Argument in a mixed space.

Returns a tuple with the function components corresponding to the subelements.

class ufl.BaseForm[source]#

Bases: object

Description of an object containing arguments.

Initialise.

arguments()[source]#

Return all Argument objects found in form.

coefficients()[source]#

Return all Coefficient objects found in form.

empty()[source]#

Returns whether the BaseForm has no components.

geometric_quantities()[source]#

Return all GeometricQuantity objects found in form.

ufl_domain()[source]#

Return the single geometric integration domain occuring in the base form.

Fails if multiple domains are found.

ufl_operands: tuple[FormArgument, ...]#
class ufl.Cell(cellname: str)[source]#

Bases: AbstractCell

Representation of a named finite element cell with known structure.

Initialise.

Parameters:

cellname – Name of the cell

property cellname: str#

Return the cellname of the cell.

property has_simplex_facets: bool#

Return True if all the facets of this cell are simplex cells.

property is_simplex: bool#

Return True if this is a simplex cell.

num_sub_entities(dim: int) int[source]#

Get the number of sub-entities of the given dimension.

reconstruct(**kwargs: Any) Cell[source]#

Reconstruct this cell, overwriting properties by those in kwargs.

sub_entities(dim: int) tuple[AbstractCell, ...][source]#

Get the sub-entities of the given dimension.

sub_entity_types(dim: int) tuple[AbstractCell, ...][source]#

Get the unique sub-entity types of the given dimension.

property topological_dimension: int#

Return the dimension of the topology of this cell.

class ufl.CellDiameter(domain)[source]#

Bases: GeometricCellQuantity

The diameter of the cell, i.e., maximal distance of two points in the cell.

Initialise.

name = 'diameter'#
class ufl.CellNormal(domain)[source]#

Bases: GeometricCellQuantity

The upwards pointing normal vector of the current manifold cell.

Initialise.

is_cellwise_constant()[source]#

Return whether this expression is spatially constant over each cell.

name = 'cell_normal'#
property ufl_shape#

Return the number of coordinates defined (i.e. the geometric dimension of the domain).

class ufl.CellVolume(domain)[source]#

Bases: GeometricCellQuantity

The volume of the cell.

Initialise.

name = 'volume'#
class ufl.Circumradius(domain)[source]#

Bases: GeometricCellQuantity

The circumradius of the cell.

Initialise.

name = 'circumradius'#
class ufl.Coargument(*args, **kw)[source]#

Bases: BaseForm, BaseArgument

UFL value: Representation of an argument to a form in a dual space.

Initialise.

arguments(outer_form=None)[source]#

Return all Argument objects found in form.

equals(other)[source]#

Check equality.

ufl_domain()[source]#

Return the UFL domain.

ufl_operands: tuple[FormArgument, ...]#
class ufl.Coefficient(*args, **kw)[source]#

Bases: FormArgument, BaseCoefficient

UFL form argument type: Representation of a form coefficient.

Initialise.

ufl_domains()[source]#

Get the UFL domains.

ufl.Coefficients(function_space)[source]#

Create a Coefficient in a mixed space.

Returns a tuple with the function components corresponding to the subelements.

class ufl.Cofunction(*args, **kw)[source]#

Bases: BaseCoefficient, BaseForm

UFL form argument type: Representation of a form coefficient from a dual space.

Initialise.

equals(other)[source]#

Check equality.

ufl_operands: tuple[FormArgument, ...]#
class ufl.Constant(domain, shape=(), count=None)[source]#

Bases: Terminal, Counted

Constant.

Initalise.

is_cellwise_constant()[source]#

Return True if the function is cellwise constant.

ufl_domain()[source]#

Get the UFL domain.

ufl_domains()[source]#

Get the UFL domains.

property ufl_shape#

Get the UFL shape.

class ufl.DSIntegralDomain(coordinate_element_positive_side: AbstractFiniteElement, coordinate_element_negative_side: AbstractFiniteElement, facet_type: AbstractCell)[source]#

Bases: AbstractIntegralDomain

Integral domain for dS.

Initialise.

property integral_type: str#

Integral type.

ufl.Dn(f)[source]#

Take the directional derivative of f in the facet normal direction.

The facet normal is Dn(f) := dot(grad(f), n).

class ufl.DsIntegralDomain(coordinate_element: AbstractFiniteElement, facet_type: AbstractCell)[source]#

Bases: AbstractIntegralDomain

Integral domain for ds.

Initialise.

property integral_type: str#

Integral type.

ufl.Dx(f, *i)[source]#

Take the partial derivative of f with respect to spatial variable number i.

Equivalent to f.dx(*i).

class ufl.DxIntegralDomain(coordinate_element: AbstractFiniteElement)[source]#

Bases: AbstractIntegralDomain

Integral domain for dx.

Initialise.

property integral_type: str#

Integral type.

class ufl.ExternalOperator(*operands, function_space, derivatives=None, argument_slots=())[source]#

Bases: BaseFormOperator

External operator.

Initialise.

Parameters:
  • operands – operands on which acts the ExternalOperator.

  • function_space – the FunctionSpace, or MixedFunctionSpace on which to build this Function.

  • derivatives – tuple specifying the derivative multiindex.

  • argument_slots – tuple composed containing expressions with ufl.Argument or ufl.Coefficient objects.

assemble(*args, **kwargs)[source]#

Assemble the external operator.

grad()[source]#

Returns the symbolic grad of the external operator.

ufl_element()[source]#

Shortcut to get the finite element of the function space of the external operator.

ufl_operands: tuple[ufl.core.terminal.FormArgument, ...]#
class ufl.FacetArea(domain)[source]#

Bases: GeometricFacetQuantity

The area of the facet.

Initialise.

name = 'facetarea'#
class ufl.FacetNormal(domain)[source]#

Bases: GeometricFacetQuantity

The outwards pointing normal vector of the current facet.

Initialise.

is_cellwise_constant()[source]#

Return whether this expression is spatially constant over each cell.

name = 'n'#
property ufl_shape#

Return the number of coordinates defined (i.e. the geometric dimension of the domain).

class ufl.Form(integrals)[source]#

Bases: BaseForm

Description of a weak form consisting of a sum of integrals over subdomains.

Initialise.

base_form_operators()[source]#

Return all BaseFormOperator objects found in form.

coefficient_numbering()[source]#

Return a contiguous numbering of coefficients in a mapping {coefficient:number}.

coefficients()[source]#

Return all Coefficient objects found in form.

constant_numbering()[source]#

Return a contiguous numbering of constants in a mapping {constant:number}.

constants()[source]#

Get constants.

domain_numbering()[source]#

Return a contiguous numbering of domains in a mapping {domain:number}.

empty()[source]#

Returns whether the form has no integrals.

equals(other)[source]#

Evaluate bool(lhs_form == rhs_form).

geometric_dimension()[source]#

Return the geometric dimension shared by all domains and functions in this form.

integrals()[source]#

Return a sequence of all integrals in form.

integrals_by_domain(domain)[source]#

Return a sequence of all integrals with a particular integration domain.

integrals_by_type(integral_type)[source]#

Return a sequence of all integrals with a particular domain type.

signature()[source]#

Signature for use with jit cache (independent of incidental numbering of indices etc).

subdomain_data()[source]#

Returns a mapping on the form {domain:{integral_type: subdomain_data}}.

terminal_numbering()[source]#

Return a contiguous numbering for all counted objects in the form.

The returned object is mapping from terminal to its number (an integer).

The numbering is computed per type so Coefficient, Constant, etc will each be numbered from zero.

ufl_cell()[source]#

Return the single cell this form is defined on.

Fails if multiple cells are found.

ufl_domains()[source]#

Return the geometric integration domains occuring in the form.

NB! This does not include domains of coefficients defined on other meshes.

The return type is a tuple even if only a single domain exists.

ufl_operands: tuple[FormArgument, ...]#
class ufl.FormSum(*args, **kwargs)[source]#

Bases: BaseForm

Form sum.

Description of a weighted sum of variational forms and form-like objects components is the list of Forms to be summed arg_weights is a list of tuples of component index and weight

Initialise.

components()[source]#

Get components.

empty()[source]#

Returns whether the FormSum has no components.

equals(other)[source]#

Evaluate bool(lhs_form == rhs_form).

ufl_domains()[source]#

Return all domains found in the base form.

ufl_operands: tuple[FormArgument, ...]#
weights()[source]#

Get weights.

class ufl.FunctionSpace(domain, element, label='')[source]#

Bases: BaseFunctionSpace, UFLObject

Representation of a Function space.

Initialise.

dual()[source]#

Get the dual of the space.

class ufl.Identity(dim)[source]#

Bases: ConstantValue

Representation of an identity matrix.

Initialise.

evaluate(x, mapping, component, index_values)[source]#

Evaluate.

ufl_shape: tuple[int, ...]#
class ufl.Index(count=None)[source]#

Bases: IndexBase, Counted

UFL value: An index with no value assigned.

Used to represent free indices in Einstein indexing notation.

Initialise.

class ufl.Integral(integrand: Expr, integral_type: str, domain, subdomain_id: int | tuple[int, ...], metadata: dict, subdomain_data: Any, extra_domain_integral_type_map: dict[Mesh, str] | None = None)[source]#

Bases: object

An integral over a single domain.

Initialise.

Parameters:
  • integrand – Integrand of the integral

  • integral_type – Type of integral, see: ufl.measure._integral_types for available types.

  • domain – Integration domain as an ufl.Mesh.

  • subdomain_id – Integer or tuple of integer restricting the integral to a subset of entities, marked in some way through subdomain_data.

  • metadata – Metadata that is usually passed to the form compiler

  • subdomain_data – Data associated with the subdomain, can be anything (depends on the compiler and user-facing API)

  • extra_domain_integral_type_map – Mapping from other ufl.Mesh objects present in integrand to specific integral types on this domain; see: ufl.measure._integral_types for possible types.

extra_domain_integral_type_map()[source]#

Return the extra domain-integral_type map.

integral_type()[source]#

Return the domain type of this integral.

integrand()[source]#

Return the integrand expression, which is an Expr instance.

metadata()[source]#

Return the compiler metadata this integral has been annotated with.

reconstruct(integrand=None, integral_type=None, domain=None, subdomain_id=None, metadata=None, subdomain_data=None, extra_domain_integral_type_map=None)[source]#

Construct a new Integral object with some properties replaced with new values.

Example

<a = Integral instance> b = a.reconstruct(expand_compounds(a.integrand())) c = a.reconstruct(metadata={‘quadrature_degree’:2})

subdomain_data()[source]#

Return the domain data of this integral.

subdomain_id()[source]#

Return the subdomain id of this integral.

ufl_domain()[source]#

Return the integration domain of this integral.

class ufl.Interpolate(expr, v)[source]#

Bases: BaseFormOperator

Symbolic representation of the interpolation operator.

Initialise.

Parameters:
  • expr – a UFL expression to interpolate.

  • v – the FunctionSpace to interpolate into or the Coargument defined on the dual of the FunctionSpace to interpolate into.

ufl_operands: tuple[ufl.core.terminal.FormArgument, ...]#
class ufl.Jacobian(domain)[source]#

Bases: GeometricCellQuantity

The Jacobian of the mapping from reference cell to spatial coordinates.

\[\begin{split}J_{ij} = \\frac{dx_i}{dX_j}\end{split}\]

Initialise.

is_cellwise_constant()[source]#

Return whether this expression is spatially constant over each cell.

name = 'J'#
property ufl_shape#

Return the number of coordinates defined (i.e. the geometric dimension of the domain).

class ufl.JacobianDeterminant(domain)[source]#

Bases: GeometricCellQuantity

The determinant of the Jacobian.

Represents the signed determinant of a square Jacobian or the pseudo-determinant of a non-square Jacobian.

Initialise.

is_cellwise_constant()[source]#

Return whether this expression is spatially constant over each cell.

name = 'detJ'#
class ufl.JacobianInverse(domain)[source]#

Bases: GeometricCellQuantity

The inverse of the Jacobian.

Represents the inverse of a square Jacobian or the pseudo-inverse of a non-square Jacobian.

Initialise.

is_cellwise_constant()[source]#

Return whether this expression is spatially constant over each cell.

name = 'K'#
property ufl_shape#

Return the number of coordinates defined (i.e. the geometric dimension of the domain).

class ufl.Label(count=None)[source]#

Bases: Terminal, Counted

Label.

Initialise.

is_cellwise_constant()[source]#

Return true if the object is constant on each cell.

ufl_domains()[source]#

Return tuple of domains related to this terminal object.

property ufl_free_indices#

Get the UFL free indices.

property ufl_index_dimensions#

Get the UFL index dimensions.

property ufl_shape#

Get the UFL shape.

class ufl.Matrix(row_space, column_space, count=None)[source]#

Bases: BaseForm, Counted

An assemble linear operator between two function spaces.

Initialise.

equals(other)[source]#

Check equality.

ufl_function_spaces()[source]#

Get the tuple of function spaces of this coefficient.

ufl_operands: tuple[FormArgument, ...]#
class ufl.MaxCellEdgeLength(domain)[source]#

Bases: GeometricCellQuantity

The maximum edge length of the cell.

Initialise.

name = 'maxcelledgelength'#
class ufl.MaxFacetEdgeLength(domain)[source]#

Bases: GeometricFacetQuantity

The maximum edge length of the facet.

Initialise.

name = 'maxfacetedgelength'#
class ufl.Measure(integral_type, domain=None, subdomain_id='everywhere', metadata=None, subdomain_data=None, intersect_measures=None)[source]#

Bases: object

Representation of an integration measure.

The Measure object holds information about integration properties to be transferred to a Form on multiplication with a scalar expression.

Initialise.

Parameters:
  • integral_type – one of “cell”, etc, or short form “dx”, etc; the integral_type on the primal domain in a multi-domain problem.

  • domain – an AbstractDomain object (most often a Mesh); the primal domain in a multi-domain problem.

  • subdomain_id – either string “everywhere”, a single subdomain id int, or tuple of ints

  • metadata – dict, with additional compiler-specific parameters affecting how code is generated, including parameters for optimization or debugging of generated code

  • subdomain_data – object representing data to interpret subdomain_id with

  • intersect_measures – a tuple of intersect measures that defines domain-integral_type map for each domain in a multi-domain problem. For instance, if the integral_type is “dS” on the primal domain, while integral_type is “ds” on domainZ, say, one must pass intersect_measures=(Measure(“ds”, domainZ),). Currently, the intersect measures must have “everywhere” subdomain_id; i.e., subdomain_id can only be specified for the primal domain.

integral_type()[source]#

Return the domain type.

Valid domain types are “cell”, “exterior_facet”, “interior_facet”, etc.

intersect_measures()[source]#

Return the intersect measures.

metadata()[source]#

Return the integral metadata.

This data is not interpreted by UFL. It is passed to the form compiler which can ignore it or use it to compile each integral of a form in a different way.

reconstruct(integral_type=None, subdomain_id=None, domain=None, metadata=None, subdomain_data=None, intersect_measures=None)[source]#

Construct a new Measure object with some properties replaced with new values.

Example

<dm = Measure instance> b = dm.reconstruct(subdomain_id=2) c = dm.reconstruct(metadata={ “quadrature_degree”: 3 })

Used by the call operator, so this is equivalent:

b = dm(2) c = dm(0, { “quadrature_degree”: 3 })

subdomain_data()[source]#

Return the integral subdomain_data.

This data is not interpreted by UFL. Its intension is to give a context in which the domain id is interpreted.

subdomain_id()[source]#

Return the domain id of this measure (integer).

ufl_domain()[source]#

Return the domain associated with this measure.

This may be None or a Domain object.

class ufl.Mesh(coordinate_element, ufl_id=None, cargo=None)[source]#

Bases: AbstractDomain, UFLObject

Symbolic representation of a mesh.

Initialise.

can_make_function_space(element: AbstractFiniteElement) bool[source]#

Check whether this mesh can make a function space with element.

is_piecewise_linear_simplex_domain()[source]#

Check if the domain is a piecewise linear simplex.

iterable_like(element: AbstractFiniteElement) Iterable[Mesh][source]#

Return iterable object that is iterable like element.

property meshes#

Return the component meshes.

ufl_cargo()[source]#

Return carried object that will not be used by UFL.

ufl_cell()[source]#

Get the cell.

ufl_coordinate_element()[source]#

Get the coordinate element.

ufl_id()#

Return the ufl_id of this object.

class ufl.MeshSequence(meshes: Sequence[Mesh])[source]#

Bases: AbstractDomain, UFLObject

Symbolic representation of a mixed mesh.

This class represents a collection of meshes that, along with a MixedElement, represent a mixed function space defined on multiple domains. This abstraction allows for defining the mixed function space with the conventional FunctionSpace class and integrating multi-domain problems seamlessly.

Currently, all component meshes must have the same cell type (and thus the same topological dimension).

Currently, one can only perform cell integrations when MeshSequence is used.

cell = triangle
mesh0 = Mesh(FiniteElement("Lagrange", cell, 1, (2,), identity_pullback, H1))
mesh1 = Mesh(FiniteElement("Lagrange", cell, 1, (2,), identity_pullback, H1))
domain = MeshSequence([mesh0, mesh1])
elem0 = FiniteElement("Lagrange", cell, 1, (), identity_pullback, H1)
elem1 = FiniteElement("Lagrange", cell, 2, (), identity_pullback, H1)
elem = MixedElement([elem0, elem1])
V = FunctionSpace(domain, elem)
v = TestFunction(V)
v0, v1 = split(v)

Initialise.

can_make_function_space(element: AbstractFiniteElement) bool[source]#

Check whether this mesh can make a function space with element.

iterable_like(element: AbstractFiniteElement) MeshSequence[source]#

Return iterable object that is iterable like element.

property meshes#

Return the component meshes.

ufl_cell()[source]#

Get the cell.

class ufl.MeshView(mesh, topological_dimension, ufl_id=None)[source]#

Bases: AbstractDomain, UFLObject

Symbolic representation of a mesh.

Initialise.

is_piecewise_linear_simplex_domain()[source]#

Check if the domain is a piecewise linear simplex.

ufl_cell()[source]#

Get the cell.

ufl_id()#

Return the ufl_id of this object.

ufl_mesh()[source]#

Get the mesh.

class ufl.MinCellEdgeLength(domain)[source]#

Bases: GeometricCellQuantity

The minimum edge length of the cell.

Initialise.

name = 'mincelledgelength'#
class ufl.MinFacetEdgeLength(domain)[source]#

Bases: GeometricFacetQuantity

The minimum edge length of the facet.

Initialise.

name = 'minfacetedgelength'#
class ufl.MixedFunctionSpace(*args)[source]#

Bases: AbstractFunctionSpace, UFLObject

Mixed function space.

Initialise.

dual(*args)[source]#

Return the dual to this function space.

If no additional arguments are passed then a MixedFunctionSpace is returned whose components are the duals of the originals.

If additional arguments are passed, these must be integers. In this case, the MixedFunctionSpace which is returned will have dual components in the positions corresponding to the arguments passed, and the original components in the other positions.

num_sub_spaces()[source]#

Return number of subspaces.

ufl_domain()[source]#

Return ufl domain.

ufl_domains()[source]#

Return ufl domains.

ufl_element()[source]#

Return ufl element.

ufl_elements()[source]#

Return ufl elements.

ufl_sub_space(i)[source]#

Return i-th ufl sub space.

ufl_sub_spaces()[source]#

Return ufl sub spaces.

class ufl.MixedPullback(element: _AbstractFiniteElement)[source]#

Bases: AbstractPullback

Pull back for a mixed element.

Initalise.

Parameters:

element – The mixed element

apply(expr, domain=None)[source]#

Apply the pull back.

Parameters:
  • expr – A function on a physical cell

  • domain – The domain on which the function is defined

Returns: The function pulled back to the reference cell

property is_identity: bool#

Is this pull back the identity (or the identity applied to mutliple components).

physical_value_shape(element, domain) tuple[int, ...][source]#

Get the physical value shape when this pull back is applied to an element on a domain.

Parameters:
  • element – The element that the pull back is applied to

  • domain – The domain

Returns:

The value shape when the pull back is applied to the given element

ufl.Not(condition)[source]#

A boolean expression (not condition) for use with conditional.

ufl.Or(left, right)[source]#

A boolean expression (left or right) for use with conditional.

class ufl.PermutationSymbol(dim)[source]#

Bases: ConstantValue

Representation of a permutation symbol.

This is also known as the Levi-Civita symbol, antisymmetric symbol, or alternating symbol.

Initialise.

evaluate(x, mapping, component, index_values)[source]#

Evaluate.

ufl_shape: tuple[int, ...]#
class ufl.RidgeJacobian(domain)[source]#

Bases: GeometricRidgeQuantity

The Jacobian of the mapping from reference ridge to spatial coordinates.

RJ_ij = dx_i/dXe_j

The RidgeJacobian is the product of the Jacobian and CellRidgeJacobian:

RJ = dx/dXe = dx/dX dX/dXe = J * CRJ

Initialise.

is_cellwise_constant()[source]#

Return whether this expression is spatially constant over each cell.

name = 'RJ'#
property ufl_shape#

Get the UFL shape.

class ufl.RidgeJacobianDeterminant(domain)[source]#

Bases: GeometricRidgeQuantity

The pseudo-determinant of the RidgeJacobian.

Initialise.

is_cellwise_constant()[source]#

Return whether this expression is spatially constant over each cell.

name = 'detRJ'#
class ufl.RidgeJacobianInverse(domain)[source]#

Bases: GeometricRidgeQuantity

The pseudo-inverse of the RidgeJacobian.

Initialise.

is_cellwise_constant()[source]#

Return whether this expression is spatially constant over each cell.

name = 'EK'#
property ufl_shape#

Get the UFL shape.

class ufl.SpatialCoordinate(domain)[source]#

Bases: GeometricCellQuantity

The coordinate in a domain.

In the context of expression integration, represents the domain coordinate of each quadrature point.

In the context of expression evaluation in a point, represents the value of that point.

Initialise.

count()[source]#

Count.

evaluate(x, mapping, component, index_values)[source]#

Return the value of the coordinate.

is_cellwise_constant()[source]#

Return whether this expression is spatially constant over each cell.

name = 'x'#
property ufl_shape#

Return the number of coordinates defined (i.e. the geometric dimension of the domain).

class ufl.SymmetricPullback(element: _AbstractFiniteElement, symmetry: dict[tuple[int, ...], int])[source]#

Bases: AbstractPullback

Pull back for an element with symmetry.

Initalise.

Parameters:
  • element – The element

  • symmetry – A dictionary mapping from the component in physical space to the local component

apply(expr, domain=None)[source]#

Apply the pull back.

Parameters:
  • expr – A function on a physical cell

  • domain – The domain on which the function is defined

Returns: The function pulled back to the reference cell

property is_identity: bool#

Is this pull back the identity (or the identity applied to mutliple components).

physical_value_shape(element, domain) tuple[int, ...][source]#

Get the physical value shape when this pull back is applied to an element on a domain.

Parameters:
  • element – The element that the pull back is applied to

  • domain – The domain

Returns:

The value shape when the pull back is applied to the given element

ufl.TensorConstant(domain, count=None)[source]#

Tensor constant.

class ufl.TensorProductCell(*cells: AbstractCell)[source]#

Bases: AbstractCell

Tensor product cell.

Initialise.

Parameters:

cells – Cells to take the tensor product of

property cellname: str#

Return the cellname of the cell.

property has_simplex_facets: bool#

Return True if all the facets of this cell are simplex cells.

property is_simplex: bool#

Return True if this is a simplex cell.

num_sub_entities(dim: int) int[source]#

Get the number of sub-entities of the given dimension.

reconstruct(**kwargs: Any) AbstractCell[source]#

Reconstruct this cell, overwriting properties by those in kwargs.

property sub_cells: tuple[AbstractCell, ...]#

Return list of cell factors.

sub_entities(dim: int) tuple[AbstractCell, ...][source]#

Get the sub-entities of the given dimension.

sub_entity_types(dim: int) tuple[AbstractCell, ...][source]#

Get the unique sub-entity types of the given dimension.

property topological_dimension: int#

Return the dimension of the topology of this cell.

ufl.TestFunction(function_space, part=None)[source]#

UFL value: Create a test function argument to a form.

ufl.TestFunctions(function_space)[source]#

Create a TestFunction in a mixed space.

Returns a tuple with the function components corresponding to the subelements.

ufl.TrialFunction(function_space, part=None)[source]#

UFL value: Create a trial function argument to a form.

ufl.TrialFunctions(function_space)[source]#

Create a TrialFunction in a mixed space.

Returns a tuple with the function components corresponding to the subelements.

ufl.VectorConstant(domain, count=None)[source]#

Vector constant.

class ufl.ZeroBaseForm(arguments)[source]#

Bases: BaseForm

Description of a zero base form.

ZeroBaseForm is idempotent with respect to assembly and is mostly used for sake of simplifying base-form expressions.

Initialise.

form#
ufl_domains()[source]#

Return all domains found in the base form.

ufl_operands: tuple[FormArgument, ...]#
ufl.acos(f)[source]#

Take the inverse cosine of f.

ufl.action(form, coefficient=None, derivatives_expanded=None)[source]#

Get the action.

Given a bilinear form, return a linear form with an additional coefficient, representing the action of the form on the coefficient. This can be used for matrix-free methods. For formbase objects,coefficient can be any object of the correct type, and this function returns an Action object.

When action is being called multiple times on the same form, expanding derivatives become expensive -> derivatives_expanded enables to use caching mechanisms to avoid that.

ufl.adjoint(form, reordered_arguments=None, derivatives_expanded=None)[source]#

Get the adjoint.

Given a combined bilinear form, compute the adjoint form by changing the ordering (count) of the test and trial functions, and taking the complex conjugate of the result.

By default, new ufl.Argument objects will be created with opposite ordering. However, if the adjoint form is to be added to other forms later, their arguments must match. In that case, the user must provide a tuple reordered_arguments=(u2,v2).

If the form is a baseform instance instead of a Form object, we return an Adjoint object instructing the adjoint to be computed at a later point.

When adjoint() is being called multiple times on the same form, expanding derivatives become expensive -> derivatives_expanded enables to use caching mechanisms to avoid that.

ufl.as_cell(cell: AbstractCell | str | tuple[AbstractCell, ...]) AbstractCell[source]#

Convert any valid object to a Cell or return cell if it is already a Cell.

Allows an already valid cell, a known cellname string, or a tuple of cells for a product cell.

ufl.as_matrix(expressions, indices=None)[source]#

As as_tensor(), but limited to rank 2 tensors.

ufl.as_tensor(expressions, indices=None)[source]#

Make a tensor valued expression.

This works in two different ways, by using indices or lists.

1) Returns \(A\) such that \(A\) [indices] = expressions. If indices are provided, expressions must be a scalar valued expression with all the provided indices among its free indices. This operator will then map each of these indices to a tensor axis, thereby making a tensor valued expression from a scalar valued expression with free indices.

2) Returns \(A\) such that \(A[k,...]\) = expressions*[k]. If no indices are provided, *expressions must be a list or tuple of expressions. The expressions can also consist of recursively nested lists to build higher rank tensors.

ufl.as_ufl(expression)[source]#

Converts expression to an Expr if possible.

ufl.as_vector(expressions, index=None)[source]#

As as_tensor(), but limited to rank 1 tensors.

ufl.asin(f)[source]#

Take the inverse sine of f.

ufl.atan(f)[source]#

Take the inverse tangent of f.

ufl.atan2(f1, f2)[source]#

Take the inverse tangent with two the arguments f1 and f2.

ufl.avg(v)[source]#

Take the average of v across a facet.

ufl.bessel_I(nu, f)[source]#

Regular modified cylindrical Bessel function.

ufl.bessel_J(nu, f)[source]#

Cylindrical Bessel function of the first kind.

ufl.bessel_K(nu, f)[source]#

Irregular modified cylindrical Bessel function.

ufl.bessel_Y(nu, f)[source]#

Cylindrical Bessel function of the second kind.

ufl.cell_avg(f)[source]#

Take the average of v over a cell.

ufl.cofac(A)[source]#

Take the cofactor of A.

ufl.conditional(condition, true_value, false_value)[source]#

A conditional expression.

This takes the value of true_value when condition evaluates to true and false_value otherwise.

ufl.conj(f)[source]#

The complex conjugate of f.

ufl.cos(f)[source]#

Take the cosine of f.

ufl.cosh(f)[source]#

Take the hyperbolic cosine of f.

ufl.cross(a, b)[source]#

Take the cross product of a and b.

ufl.curl(f)[source]#

Take the curl of f.

ufl.derivative(form, coefficient, argument=None, coefficient_derivatives=None)[source]#

Compute the Gateaux derivative of form w.r.t. coefficient in direction of argument.

If the argument is omitted, a new Argument is created in the same space as the coefficient, with argument number one higher than the highest one in the form.

The resulting form has one additional Argument in the same finite element space as the coefficient.

A tuple of Coefficient s may be provided in place of a single Coefficient, in which case the new Argument argument is based on a MixedElement created from this tuple.

An indexed Coefficient from a mixed space may be provided, in which case the argument should be in the corresponding subspace of the coefficient space.

If provided, coefficient_derivatives should be a mapping from Coefficient instances to their derivatives w.r.t. coefficient.

ufl.det(A)[source]#

Take the determinant of A.

ufl.dev(A)[source]#

Take the deviatoric part of A.

ufl.diag(A)[source]#

Diagonal ranl-2 tensor.

Take the diagonal part of rank 2 tensor A or make a diagonal rank 2 tensor from a rank 1 tensor.

Always returns a rank 2 tensor. See also diag_vector.

ufl.diag_vector(A)[source]#

Take the diagonal part of rank 2 tensor A and return as a vector.

See also diag.

ufl.diff(f, v)[source]#

Take the derivative of f with respect to the variable v.

If f is a form, diff is applied to each integrand.

ufl.div(f)[source]#

Take the divergence of f.

This operator follows the div convention where

div(v) = v[i].dx(i)
div(T)[:] = T[:,i].dx(i)

for vector expressions v, and arbitrary rank tensor expressions T.

See also: nabla_div()

ufl.dot(a, b)[source]#

Take the dot product of a and b.

This won’t take the complex conjugate of the second argument.

ufl.elem_div(A, B)[source]#

Take the elementwise division of tensors A and B with the same shape.

ufl.elem_mult(A, B)[source]#

Take the elementwise multiplication of tensors A and B with the same shape.

ufl.elem_op(op, *args)[source]#

Apply element-wise operations.

Take the element-wise application of operator op on scalar values from one or more tensor arguments.

ufl.elem_pow(A, B)[source]#

Take the elementwise power of tensors A and B with the same shape.

ufl.energy_norm(form, coefficient=None)[source]#

Get the energy norm.

Given a bilinear form a and a coefficient f, return the functional \(a(f,f)\).

ufl.eq(left, right)[source]#

A boolean expression (left == right) for use with conditional.

ufl.erf(f)[source]#

Take the error function of f.

ufl.exp(f)[source]#

Take the exponential of f.

ufl.exterior_derivative(f)[source]#

Take the exterior derivative of f.

The exterior derivative uses the element Sobolev space to determine whether id, grad, curl or div should be used.

Note that this uses the grad and div operators, as opposed to nabla_grad and nabla_div.

ufl.extract_blocks(form: Form, i: int | None = None, j: int | None = None, replace_argument: bool = True)[source]#

Extract blocks.

Given a linear or bilinear form on a mixed space, extract the block corresponding to the indices ix, iy.

Example

a = inner(grad(u), grad(v))*dx + div(u)*q*dx + div(v)*p*dx extract_blocks(a, 0, 0) -> inner(grad(u), grad(v))*dx extract_blocks(a) -> [inner(grad(u), grad(v))*dx, div(v)*p*dx, div(u)*q*dx, 0]

ufl.facet_avg(f)[source]#

Take the average of v over a facet.

ufl.functional(form)[source]#

Extract the functional part of form.

ufl.ge(left, right)[source]#

A boolean expression (left >= right) for use with conditional.

ufl.grad(f)[source]#

Take the gradient of f.

This operator follows the grad convention where

grad(s)[i] = s.dx(i)
grad(v)[i, j] = v[i].dx(j)
grad(T)[:, i] = T[:].dx(i)

for scalar expressions s, vector expressions v, and arbitrary rank tensor expressions T.

See also: nabla_grad()

ufl.gt(left, right)[source]#

A boolean expression (left > right) for use with conditional.

ufl.imag(f)[source]#

The imaginary part of f.

ufl.indices(n)[source]#

Return a tuple of n new Index objects.

ufl.inner(a, b)[source]#

Take the inner product of a and b.

The complex conjugate of the second argument is taken.

ufl.integral_types()[source]#

Return a tuple of all domain type strings.

ufl.interpolate(expr, v)[source]#

Create symbolic representation of the interpolation operator.

Parameters:
  • expr – a UFL expression to interpolate.

  • v – the FunctionSpace to interpolate into or the Coargument defined on the dual of the FunctionSpace to interpolate into.

ufl.inv(A)[source]#

Take the inverse of A.

ufl.jump(v, n=None)[source]#

Take the jump of v across a facet.

ufl.le(left, right)[source]#

A boolean expression (left <= right) for use with conditional.

ufl.lhs(form)[source]#

Get the left hand side.

Given a combined bilinear and linear form, extract the left hand side (bilinear form part).

Example

a = u*v*dx + f*v*dx a = lhs(a) -> u*v*dx

ufl.ln(f)[source]#

Take the natural logarithm of f.

ufl.lt(left, right)[source]#

A boolean expression (left < right) for use with conditional.

ufl.max_value(x, y)[source]#

Take the maximum of x and y.

ufl.min_value(x, y)[source]#

Take the minimum of x and y.

ufl.nabla_div(f)[source]#

Take the divergence of f.

This operator follows the div convention where

nabla_div(v) = v[i].dx(i)
nabla_div(T)[:] = T[i,:].dx(i)

for vector expressions v, and arbitrary rank tensor expressions T.

See also: div()

ufl.nabla_grad(f)[source]#

Take the gradient of f.

This operator follows the grad convention where

nabla_grad(s)[i] = s.dx(i)
nabla_grad(v)[i,j] = v[j].dx(i)
nabla_grad(T)[i,:] = T[:].dx(i)

for scalar expressions s, vector expressions v, and arbitrary rank tensor expressions T.

See also: grad()

ufl.ne(left, right)[source]#

A boolean expression (left != right) for use with conditional.

ufl.outer(*operands)[source]#

Take the outer product of two or more operands.

The complex conjugate of the first argument is taken.

ufl.perp(v)[source]#

Take the perp of v.

I.e. \((-v_1, +v_0)\).

ufl.product(sequence)[source]#

Return the product of all elements in a sequence.

ufl.rank(f)[source]#

The rank of f.

ufl.real(f)[source]#

The real part of f.

ufl.register_integral_type(integral_type, measure_name)[source]#

Register an integral type.

ufl.replace(e, mapping)[source]#

Replace subexpressions in expression.

Params:

e: An Expr or Form mapping: A dict with from:to replacements to perform

Returns:

The expression with replacements performed

ufl.rhs(form)[source]#

Get the right hand side.

Given a combined bilinear and linear form, extract the right hand side (negated linear form part).

Example

a = u*v*dx + f*v*dx L = rhs(a) -> -f*v*dx

ufl.rot(f)#

Take the curl of f.

ufl.sensitivity_rhs(a, u, L, v)[source]#

Compute the right hand side for a sensitivity calculation system.

The derivation behind this computation is as follows. Assume a, L to be bilinear and linear forms corresponding to the assembled linear system

\[Ax = b.\]

Where x is the vector of the discrete function corresponding to u. Let v be some scalar variable this equation depends on. Then we can write

\[ \begin{align}\begin{aligned}\begin{split}0 = \\frac{d}{dv}(Ax-b) = \\frac{dA}{dv} x + A \\frac{dx}{dv} - \\frac{db}{dv},\end{split}\\\begin{split}A \\frac{dx}{dv} = \\frac{db}{dv} - \\frac{dA}{dv} x,\end{split}\end{aligned}\end{align} \]

and solve this system for \(\\frac{dx}{dv}\), using the same bilinear form a and matrix A from the original system. Assume the forms are written

v = variable(v_expression)
L = IL(v) * dx
a = Ia(v) * dx

where IL and Ia are integrand expressions. Define a Coefficient u representing the solution to the equations. Then we can compute \(\\frac{db}{dv}\) and \(\\frac{dA}{dv}\) from the forms

da = diff(a, v)
dL = diff(L, v)

and the action of da on u by

dau = action(da, u)

In total, we can build the right hand side of the system to compute \(\\frac{du}{dv}\) with the single line

dL = diff(L, v) - action(diff(a, v), u)

or, using this function,

dL = sensitivity_rhs(a, u, L, v)
ufl.shape(f)[source]#

The shape of f.

ufl.sign(x)[source]#

Return the sign of x.

This returns +1 if x is positive, -1 if x is negative, and 0 if x is 0.

ufl.sin(f)[source]#

Take the sine of f.

ufl.sinh(f)[source]#

Take the hyperbolic sine of f.

ufl.skew(A)[source]#

Take the skew symmetric part of A.

ufl.split(v)[source]#

Split a coefficient or argument.

If v is a Coefficient or Argument in a mixed space, returns a tuple with the function components corresponding to the subelements.

ufl.sqrt(f)[source]#

Take the square root of f.

ufl.sym(A)[source]#

Take the symmetric part of A.

ufl.system(form)[source]#

Split a form into the left hand side and right hand side.

See lhs and rhs.

ufl.tan(f)[source]#

Take the tangent of f.

ufl.tanh(f)[source]#

Take the hyperbolic tangent of f.

ufl.tr(A)[source]#

Take the trace of A.

ufl.transpose(A)[source]#

Take the transposed of tensor A.

ufl.unit_matrices(d)[source]#

A tuple of constant unit matrices in all directions with dimension d.

ufl.unit_matrix(i, j, d)[source]#

A constant unit matrix in direction i,*j* with dimension d.

ufl.unit_vector(i, d)[source]#

A constant unit vector in direction i with dimension d.

ufl.unit_vectors(d)[source]#

A tuple of constant unit vectors in all directions with dimension d.

ufl.variable(e)[source]#

Define a variable representing the given expression.

See also diff().

ufl.zero(*shape)[source]#

UFL literal constant: Return a zero tensor with the given shape.