ufl Package

ufl Package

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 mathematical notation.

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

  • To import the language, type:

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

    from ufl.classes import *
    
  • Various algorithms for working with UFL expression trees can be found in:

    from ufl.algorithms import *
    

The 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:

  • Domains:

    Domain, ProductDomain
    
  • Cells:

    Cell, ProductCell, OuterProductCell,
    interval, triangle, tetrahedron,
    quadrilateral, hexahedron,
    cell1D, cell2D, cell3D,
    
  • Sobolev spaces:

    L2, H1, H2, HDiv, HCurl
    
  • Elements:

    FiniteElement,
    MixedElement, VectorElement, TensorElement
    EnrichedElement, RestrictedElement,
    TensorProductElement, OuterProductElement,
    OuterProductVectorElement, HDiv, HCurl
    
  • Arguments:

    Argument, TestFunction, TrialFunction
    
  • Coefficients:

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

    split
    
  • Literal constants:

    Identity, PermutationSymbol
    
  • Geometric quantities:

    SpatialCoordinate,
    FacetNormal, CellNormal,
    CellVolume, Circumradius,
    FacetArea, MinFacetEdgeLength, MaxFacetEdgeLength,
    
  • 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, div, nabla_grad, nabla_div,
    Dx, Dn, curl, rot
    
  • Nonlinear functions:

    Max, Min,
    abs, sign, sqrt,
    exp, ln, erf,
    cos, sin, tan,
    acos, asin, atan, atan_2,
    cosh, sinh, tanh,
    bessel_J, bessel_Y, bessel_I, bessel_K
    
  • Discontinuous Galerkin operators:

    jump, avg, v(‘+’), v(‘-‘), cell_avg, facet_avg

  • Conditional operators:

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

    dx, ds, ds_b, ds_t, ds_tb, ds_v, dS, dS_h, dS_v, dP, dc, dE
    
  • Form transformations:

    rhs, lhs, system, functional,
    replace, replace_integral_domains,
    adjoint, action, energy_norm,
    sensitivity_rhs, derivative
    

algebra Module

Basic algebra operations.

class ufl.algebra.Abs(a)

Bases: ufl.operatorbase.AlgebraOperator

evaluate(x, mapping, component, index_values)
free_indices()
index_dimensions()
operands()
shape()
class ufl.algebra.Division(a, b)

Bases: ufl.operatorbase.AlgebraOperator

evaluate(x, mapping, component, index_values)
free_indices()
index_dimensions()
operands()
shape()
class ufl.algebra.Power(a, b)

Bases: ufl.operatorbase.AlgebraOperator

evaluate(x, mapping, component, index_values)
free_indices()
index_dimensions()
operands()
shape()
class ufl.algebra.Product(a, b)

Bases: ufl.operatorbase.AlgebraOperator

The product of two or more UFL objects.

evaluate(x, mapping, component, index_values)
free_indices()
index_dimensions()
operands()
shape()
class ufl.algebra.Sum(a, b)

Bases: ufl.operatorbase.AlgebraOperator

evaluate(x, mapping, component, index_values)
free_indices()
index_dimensions()
operands()
shape()

argument Module

This module defines the class Argument and a number of related classes (functions), including TestFunction and TrialFunction.

class ufl.argument.Argument(element, number, part=None)

Bases: ufl.terminal.FormArgument

UFL value: Representation of an argument to a form.

count()
domains()

Return tuple of domains related to this terminal object.

element()
is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

number()
part()
reconstruct(element=None, number=None, part=None)
shape()
signature_data(domain_numbering)

Signature data for form arguments depend on the global numbering of the form arguments and domains.

ufl.argument.Arguments(element, number)

UFL value: Create an Argument in a mixed space, and return a tuple with the function components corresponding to the subelements.

ufl.argument.TestFunction(element, part=None)

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

ufl.argument.TestFunctions(element)

UFL value: Create a TestFunction in a mixed space, and return a tuple with the function components corresponding to the subelements.

ufl.argument.TrialFunction(element, part=None)

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

ufl.argument.TrialFunctions(element)

UFL value: Create a TrialFunction in a mixed space, and return a tuple with the function components corresponding to the subelements.

assertions Module

This module provides assertion functions used by the UFL implementation.

ufl.assertions.expecting_expr(v)
ufl.assertions.expecting_instance(v, c)
ufl.assertions.expecting_python_scalar(v)
ufl.assertions.expecting_terminal(v)
ufl.assertions.expecting_true_ufl_scalar(v)
ufl.assertions.ufl_assert(condition, *message)

Assert that condition is true and otherwise issue an error with given message.

cell Module

Types for representing a cell.

class ufl.cell.Cell(cellname, geometric_dimension=None, topological_dimension=None)

Bases: object

Representation of a finite element cell.

Initialize basic cell description.

cellname()

Return the cellname of the cell.

circumradius

UFL geometry value. Deprecated, please use the constructor types instead.

d

The dimension of the cell.

Only valid if the geometric and topological dimensions are the same.

facet_area

UFL geometry value. Deprecated, please use the constructor types instead.

facet_cellname()

Return the cellname of the facet of this cell.

geometric_dimension()

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

n

UFL geometry value. Deprecated, please use the constructor types instead.

topological_dimension()

Return the dimension of the topology of this cell.

volume

UFL geometry value. Deprecated, please use the constructor types instead.

x

UFL geometry value. Deprecated, please use the constructor types instead.

class ufl.cell.OuterProductCell(A, B)

Bases: ufl.cell.Cell

Representation of a cell formed as the Cartesian product of two existing cells

facet_cellname()

Return the cellname of the facet of this cell.

facet_horiz
facet_vert
class ufl.cell.ProductCell(*cells)

Bases: ufl.cell.Cell

facet_cellname()

Return the cellname of the facet of this cell.

sub_cells()

Return list of cell factors.

ufl.cell.as_cell(cell)

Convert any valid object to a Cell (in particular, cellname string), or return cell if it is already a Cell.

ufl.cell.cell2dim(cell)

Maps from UFL cell or cell name to topological dimension

checks Module

Utility functions for checking properties of expressions.

ufl.checks.is_globally_constant(expr)

Check if an expression is globally constant, which includes spatially independent constant coefficients that are not known before assembly time.

ufl.checks.is_python_scalar(expression)

Return True iff expression is of a Python scalar type.

ufl.checks.is_scalar_constant_expression(expr)

Check if an expression is a globally constant scalar expression.

ufl.checks.is_true_ufl_scalar(expression)

Return True iff expression is scalar-valued, with no free indices.

ufl.checks.is_ufl_scalar(expression)

Return True iff expression is scalar-valued, but possibly containing free indices.

classes Module

This file is useful for external code like tests and form compilers, since it enables the syntax “from ufl.classes import FooBar” for getting implementation details not exposed through the default ufl namespace. It also contains functionality used by algorithms for dealing with groups of classes, and for mapping types to different handler functions.

coefficient Module

This module defines the Coefficient class and a number of related classes, including Constant.

class ufl.coefficient.Coefficient(element, count=None)

Bases: ufl.terminal.FormArgument

UFL form argument type: Representation of a form coefficient.

count()
domains()

Return tuple of domains related to this terminal object.

element()
is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

reconstruct(element=None, count=None)
shape()
signature_data(count, domain_numbering)

Signature data for form arguments depend on the global numbering of the form arguments and domains.

ufl.coefficient.Coefficients(element)

UFL value: Create a Coefficient in a mixed space, and return a tuple with the function components corresponding to the subelements.

class ufl.coefficient.Constant(domain, count=None)

Bases: ufl.coefficient.ConstantBase

UFL value: Represents a globally constant scalar valued coefficient.

class ufl.coefficient.ConstantBase(element, count)

Bases: ufl.coefficient.Coefficient

class ufl.coefficient.TensorConstant(domain, shape=None, symmetry=None, count=None)

Bases: ufl.coefficient.ConstantBase

UFL value: Represents a globally constant tensor valued coefficient.

class ufl.coefficient.VectorConstant(domain, dim=None, count=None)

Bases: ufl.coefficient.ConstantBase

UFL value: Represents a globally constant vector valued coefficient.

common Module

This module contains a collection of common utilities.

class ufl.common.EmptyDictType

Bases: dict

update(*args, **kwargs)
class ufl.common.ExampleCounted(count=None)

Bases: object

An example class for classes of objects identified by a global counter.

The old inheritance pattern is deprecated. Mimic this class instead.

count()
class ufl.common.Stack(*args)

Bases: list

A stack datastructure.

peek()
push(v)
class ufl.common.StackDict(*args, **kwargs)

Bases: dict

A dict that can be changed incrementally with ‘d.push(k,v)’ and have changes rolled back with ‘k,v = d.pop()’.

pop()
push(k, v)
class ufl.common.Timer(name)
end()
class ufl.common.UFLTypeDefaultDict(default)

Bases: dict

class ufl.common.UFLTypeDict

Bases: dict

ufl.common.and_tuples(seqa, seqb)

Return ‘and’ of all pairs in two sequences of same length.

ufl.common.camel2underscore(name)

Convert a CamelCaps string to underscore_syntax.

ufl.common.component_to_index(component, shape)
ufl.common.counted_init(self, count=None, countedclass=None)
ufl.common.dict_sum(items)

Construct a dict, in between dict(items) and sum(items), by accumulating items for each key.

ufl.common.dstr(d, colsize=80)

Pretty-print dictionary of key-value pairs.

ufl.common.estr(elements)

Format list of elements for printing.

ufl.common.fast_post_traversal(expr)

Yields o for each tree node o in expr, child before parent.

ufl.common.fast_post_traversal2(expr, visited=None)

Yields o for each tree node o in expr, child before parent.

ufl.common.fast_pre_traversal(expr)

Yields o for each tree node o in expr, parent before child.

ufl.common.fast_pre_traversal2(expr, visited=None)

Yields o for each tree node o in expr, parent before child.

This version only visits each node once!

ufl.common.get_status_output(cmd, input=None, cwd=None, env=None)
ufl.common.index_to_component(index, shape)
ufl.common.istr(o)

Format object as string, inserting ? for None.

ufl.common.iter_tree(tree)

Iterate over all nodes in a tree represented by lists of lists of leaves.

ufl.common.lstr(l)

Pretty-print list or tuple, invoking str() on items instead of repr() like str() does.

ufl.common.mergedicts(dicts)
ufl.common.mergedicts2(d1, d2)
ufl.common.openpdf(pdffilename)

Open PDF file in external pdf viewer.

ufl.common.or_tuples(seqa, seqb)

Return ‘or’ of all pairs in two sequences of same length.

ufl.common.pdflatex(latexfilename, pdffilename, flags='')

Execute pdflatex to compile a latex file into pdf.

ufl.common.product(sequence)

Return the product of all elements in a sequence.

ufl.common.recursive_chain(lists)
ufl.common.slice_dict(dictionary, keys, default=None)
ufl.common.some_key(a_dict)

Return an arbitrary key from a dictionary.

ufl.common.sorted_by_count(seq)
ufl.common.sorted_items(mapping)
ufl.common.split_dict(d, criteria)

Split a dict d into two dicts based on a criteria on the keys.

ufl.common.sstr(s)

Pretty-print set.

ufl.common.strides(shape)
ufl.common.subdict(superdict, keys)
ufl.common.topological_sorting(nodes, edges)

Return a topologically sorted list of the nodes

Implemented algorithm from Wikipedia :P

<http://en.wikipedia.org/wiki/Topological_sorting>

No error for cyclic edges...

ufl.common.tstr(t, colsize=80)

Pretty-print list of tuples of key-value pairs.

ufl.common.unique_post_traversal(expr, visited=None)

Yields o for each node o in expr, child before parent.

Never visits a node twice.

ufl.common.unique_pre_traversal(expr, visited=None)

Yields o for each tree node o in expr, parent before child.

This version only visits each node once!

ufl.common.unzip(seq)

Inverse operation of zip: unzip(zip(a, b)) == (a, b)

ufl.common.write_file(filename, text)
ufl.common.xor(a, b)

compound_expressions Module

Functions implementing compound expressions as equivalent representations using basic operators.

ufl.compound_expressions.adj_expr(A)
ufl.compound_expressions.adj_expr_2x2(A)
ufl.compound_expressions.adj_expr_3x3(A)
ufl.compound_expressions.adj_expr_4x4(A)
ufl.compound_expressions.cofactor_expr(A)
ufl.compound_expressions.cofactor_expr_2x2(A)
ufl.compound_expressions.cofactor_expr_3x3(A)
ufl.compound_expressions.cofactor_expr_4x4(A)
ufl.compound_expressions.cross_expr(a, b)
ufl.compound_expressions.determinant_expr(A)

Compute the determinant of A.

ufl.compound_expressions.determinant_expr_2x2(B)
ufl.compound_expressions.determinant_expr_3x3(A)
ufl.compound_expressions.deviatoric_expr(A)
ufl.compound_expressions.deviatoric_expr_2x2(A)
ufl.compound_expressions.deviatoric_expr_3x3(A)
ufl.compound_expressions.inverse_expr(A)

Compute the inverse of A.

ufl.compound_expressions.pseudo_determinant_expr(A)

Compute the pseudo-determinant of A: sqrt(det(A.T*A)).

ufl.compound_expressions.pseudo_inverse_expr(A)

Compute the Penrose-Moore pseudo-inverse of A: (A.T*A)^-1 * A.T.

conditional Module

This module defines classes for conditional expressions.

class ufl.conditional.AndCondition(left, right)

Bases: ufl.conditional.BinaryCondition

evaluate(x, mapping, component, index_values)
class ufl.conditional.BinaryCondition(name, left, right)

Bases: ufl.conditional.Condition

operands()
class ufl.conditional.Condition

Bases: ufl.operatorbase.Operator

free_indices()
index_dimensions()
shape()
class ufl.conditional.Conditional(condition, true_value, false_value)

Bases: ufl.operatorbase.Operator

evaluate(x, mapping, component, index_values)
free_indices()
index_dimensions()
operands()
shape()
class ufl.conditional.EQ(left, right)

Bases: ufl.conditional.BinaryCondition

evaluate(x, mapping, component, index_values)
class ufl.conditional.GE(left, right)

Bases: ufl.conditional.BinaryCondition

evaluate(x, mapping, component, index_values)
class ufl.conditional.GT(left, right)

Bases: ufl.conditional.BinaryCondition

evaluate(x, mapping, component, index_values)
class ufl.conditional.LE(left, right)

Bases: ufl.conditional.BinaryCondition

evaluate(x, mapping, component, index_values)
class ufl.conditional.LT(left, right)

Bases: ufl.conditional.BinaryCondition

evaluate(x, mapping, component, index_values)
class ufl.conditional.NE(left, right)

Bases: ufl.conditional.BinaryCondition

evaluate(x, mapping, component, index_values)
class ufl.conditional.NotCondition(condition)

Bases: ufl.conditional.Condition

evaluate(x, mapping, component, index_values)
operands()
class ufl.conditional.OrCondition(left, right)

Bases: ufl.conditional.BinaryCondition

evaluate(x, mapping, component, index_values)

constantvalue Module

This module defines classes representing constant values.

class ufl.constantvalue.ConstantValue

Bases: ufl.terminal.Terminal

domains()

Return tuple of domains related to this terminal object.

is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

class ufl.constantvalue.FloatValue(value, shape=(), free_indices=(), index_dimensions=None)

Bases: ufl.constantvalue.ScalarValue

UFL literal type: Representation of a constant scalar floating point value.

class ufl.constantvalue.Identity(dim)

Bases: ufl.constantvalue.ConstantValue

UFL literal type: Representation of an identity matrix.

evaluate(x, mapping, component, index_values)
shape()
class ufl.constantvalue.IndexAnnotated(shape=(), free_indices=(), index_dimensions=None)

Bases: ufl.constantvalue.ConstantValue

Class to annotate expressions that don’t depend on indices with a set of free indices, used internally to keep index properties intact during automatic differentiation.

class ufl.constantvalue.IntValue(value, shape=(), free_indices=(), index_dimensions=None)

Bases: ufl.constantvalue.ScalarValue

UFL literal type: Representation of a constant scalar integer value.

class ufl.constantvalue.PermutationSymbol(dim)

Bases: ufl.constantvalue.ConstantValue

UFL literal type: Representation of a permutation symbol.

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

evaluate(x, mapping, component, index_values)
shape()
class ufl.constantvalue.ScalarValue(value, shape=(), free_indices=(), index_dimensions=None)

Bases: ufl.constantvalue.IndexAnnotated

A constant scalar value.

evaluate(x, mapping, component, index_values)
free_indices()
index_dimensions()
reconstruct(free_indices=None)

Reconstruct with new free indices.

shape()
value()
class ufl.constantvalue.Zero(shape=(), free_indices=(), index_dimensions=None)

Bases: ufl.constantvalue.IndexAnnotated

UFL literal type: Representation of a zero valued expression.

evaluate(x, mapping, component, index_values)
free_indices()
index_dimensions()
reconstruct(free_indices=None)
shape()
ufl.constantvalue.as_ufl(expression)

Converts expression to an Expr if possible.

ufl.constantvalue.format_float(x)

Format float value based on global UFL precision.

ufl.constantvalue.zero(*shape)

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

differentiation Module

Differential operators.

class ufl.differentiation.CoefficientDerivative(integrand, coefficients, arguments, coefficient_derivatives)

Bases: ufl.differentiation.Derivative

Derivative of the integrand of a form w.r.t. the degrees of freedom in a discrete Coefficient.

free_indices()
index_dimensions()
operands()
shape()
class ufl.differentiation.CompoundDerivative

Bases: ufl.differentiation.Derivative

Base class for all compound derivative types.

class ufl.differentiation.Curl(f)

Bases: ufl.differentiation.CompoundDerivative

free_indices()
index_dimensions()
operands()
shape()
class ufl.differentiation.Derivative

Bases: ufl.operatorbase.Operator

Base class for all derivative types.

class ufl.differentiation.Div(f)

Bases: ufl.differentiation.CompoundDerivative

free_indices()
index_dimensions()
operands()
shape()
class ufl.differentiation.Grad(f)

Bases: ufl.differentiation.CompoundDerivative

evaluate(x, mapping, component, index_values, derivatives=())

Get child from mapping and return the component asked for.

free_indices()
index_dimensions()
operands()
reconstruct(op)

Return a new object of the same type with new operands.

shape()
class ufl.differentiation.NablaDiv(f)

Bases: ufl.differentiation.CompoundDerivative

free_indices()
index_dimensions()
operands()
shape()
class ufl.differentiation.NablaGrad(f)

Bases: ufl.differentiation.CompoundDerivative

free_indices()
index_dimensions()
operands()
reconstruct(op)

Return a new object of the same type with new operands.

shape()
class ufl.differentiation.ReferenceGrad(f)

Bases: ufl.differentiation.CompoundDerivative

evaluate(x, mapping, component, index_values, derivatives=())

Get child from mapping and return the component asked for.

free_indices()
index_dimensions()
operands()
reconstruct(op)

Return a new object of the same type with new operands.

shape()
class ufl.differentiation.VariableDerivative(f, v)

Bases: ufl.differentiation.Derivative

free_indices()
index_dimensions()
operands()
shape()
ufl.differentiation.split_indices(expression, idx)

domain Module

Types for representing a geometric domain.

class ufl.domain.Domain(*args, **kwargs)

Bases: object

Symbolic representation of a geometrical domain.

Used in the definition of geometric terminal expressions, finite element spaces, and integration measures.

Takes a single positional argument which is either the cell of the underlying mesh

D = Domain(triangle)

or the coordinate field which is a vector valued Coefficient.

P2 = VectorElement(“CG”, D, 2) x = Coefficient(P2) E = Domain(x)

With the cell variant of the constructor, an optional label can be passed to distinguish two domains from each other.

Da = Domain(cell, label=”a”) Db = Domain(cell, label=”b”)

an optional data argument can also be passed, for integration with problem solver environments (e.g. dolfin), this is typically the underlying mesh.

Da = Domain(cell, label=”a”, data=mesha) Db = Domain(cell, label=”b”, data=meshb)
cell()

Return the cell this domain is defined in terms of.

coordinates()

Return the coordinate vector field this domain is defined in terms of.

data()

Return attached data object.

geometric_dimension()

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

hash_data()
is_piecewise_linear_simplex_domain()
label()

Return the label identifying this domain. None means no label has been set.

reconstruct(cell=None, coordinates=None, label=None, data=None)

Create a new Domain object with possibly changed label or data.

reconstruction_signature()

Format as string for evaluation as Python object.

For use with cross language frameworks, stored in generated code and evaluated later in Python to reconstruct this object.

This differs from repr in that it does not include domain label and data or coordinates, which must be reconstructed or supplied by other means.

signature_data(domain_numbering)

Signature data of domain depend on the global domain numbering.

topological_dimension()

Return the dimension of the topology of this domain.

class ufl.domain.IntersectionDomain(domain1, domain2, label=None, data=None)

Bases: ufl.domain.Domain

WARNING: This is work in progress, design is in no way completed.

child_domains()
class ufl.domain.OverlapDomain(domain1, domain2, label=None, data=None)

Bases: ufl.domain.Domain

WARNING: This is work in progress, design is in no way completed.

child_domains()
class ufl.domain.ProductDomain(domains, data=None)

Bases: ufl.domain.Domain

WARNING: This is work in progress, design is in no way completed.

child_domains()
ufl.domain.as_domain(domain)

Convert any valid object to a Domain (in particular, cell or cellname string), or return domain if it is already a Domain.

ufl.domain.check_domain_compatibility(domains)
ufl.domain.extract_domains(expr)
ufl.domain.join_domains(domains)

Take a list of Domains and return a list with only unique domain objects.

Checks that domains with the same label are compatible, and allows data to be None or

ufl.domain.join_subdomain_data(subdomain_datas)

equation Module

The Equation class, used to express equations like a == L.

class ufl.equation.Equation(lhs, rhs)

This class is used to represent equations expressed by the “==” operator. Examples include a == L and F == 0 where a, L and F are Form objects.

Create equation lhs == rhs

expr Module

This module defines the Expr class, the superclass for all expression tree node types in UFL.

NB! A note about other operators not implemented here:

More operators (special functions) on Exprs are defined in exproperators.py, as well as the transpose “A.T” and spatial derivative “a.dx(i)”. This is to avoid circular dependencies between Expr and its subclasses.

class ufl.expr.Expr

Bases: object

Base class for all UFL objects.

T

Transposed a rank two tensor expression. For more general transpose operations of higher order tensor expressions, use indexing and Tensor.

cell()

Return the cell this expression is defined on.

domain()

Return the single unique domain this expression is defined on or throw an error.

domains()
dx(*ii)

Return the partial derivative with respect to spatial variable number i.

evaluate(x, mapping, component, index_values)

Evaluate expression at given coordinate with given values for terminals.

free_indices()

Return a tuple with the free indices (unassigned) of the expression.

geometric_dimension()

Return the geometric dimension this expression lives in.

index_dimensions()

Return a dict with the free or repeated indices in the expression as keys and the dimensions of those indices as values.

is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

operands()

Return a sequence with all subtree nodes in expression tree.

rank()

Return the tensor rank of the expression.

reconstruct(*operands)

Return a new object of the same type with new operands.

shape()

Return the tensor shape of the expression.

signature_data()

Return data that uniquely identifies form compiler relevant aspects of this object.

x__del__()
ufl.expr.print_expr_statistics()

exprcontainers Module

This module defines special types for representing mapping of expressions to expressions.

class ufl.exprcontainers.ExprList(*operands)

Bases: ufl.operatorbase.WrapperType

List of Expr objects. For internal use, never to be created by end users.

operands()
class ufl.exprcontainers.ExprMapping(*operands)

Bases: ufl.operatorbase.WrapperType

Mapping of Expr objects. For internal use, never to be created by end users.

domains()
operands()

exprequals Module

ufl.exprequals.expr_equals(self, other)

Checks whether the two expressions are represented the exact same way. This does not check if the expressions are mathematically equal or equivalent! Used by sets and dicts.

ufl.exprequals.print_collisions()

exproperators Module

This module attaches special functions to Expr. This way we avoid circular dependencies between e.g. Sum and its superclass Expr.

ufl.exproperators.analyse_key(ii, rank)

Takes something the user might input as an index tuple inside [], which could include complete slices (:) and ellipsis (...), and returns tuples of actual UFL index objects.

The return value is a tuple (indices, axis_indices), each being a tuple of IndexBase instances.

The return value ‘indices’ corresponds to all input objects of these types: - Index - FixedIndex - int => Wrapped in FixedIndex

The return value ‘axis_indices’ corresponds to all input objects of these types: - Complete slice (:) => Replaced by a single new index - Ellipsis (...) => Replaced by multiple new indices

form Module

The Form class.

class ufl.form.Form(integrals)

Bases: object

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

cell()
compute_form_data(object_names=None)

Compute and return form metadata

domain()

Return the geometric integration domain occuring in the form.

NB! This does not include domains of coefficients defined on other meshes, look at form data for that additional information.

domains()

Return the geometric integration domains occuring in the form.

NB! This does not include domains of coefficients defined on other meshes, look at form data for that additional information.

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

empty()
equals(other)

Evaluate ‘bool(lhs_form == rhs_form)’.

form_data()

Return form metadata (None if form has not been preprocessed)

integrals()

Return a sequence of all integrals in form.

integrals_by_type(integral_type)

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

is_preprocessed()

Return true if this form is the result of a preprocessing of another form.

x_repr_latex_()
x_repr_png_()
ufl.form.as_form(form)

Convert to form if not a form, otherwise return form.

ufl.form.integral_sort_key(integral)
ufl.form.replace_integral_domains(form, common_domain)

Given a form and a domain, assign a common integration domain to all integrals.

Does not modify the input form (Form should always be immutable). This is to support ill formed forms with no domain specified, some times occuring in pydolfin, e.g. assemble(1*dx, mesh=mesh).

formoperators Module

Various high level ways to transform a complete Form into a new Form.

ufl.formoperators.action(form, coefficient=None)

UFL form operator: 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.

ufl.formoperators.adjoint(form, reordered_arguments=None)

UFL form operator: Given a combined bilinear form, compute the adjoint form by changing the ordering (count) of the test and trial functions.

By default, new 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).

ufl.formoperators.derivative(form, coefficient, argument=None, coefficient_derivatives=None)

UFL form operator: 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 Coefficients 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.formoperators.energy_norm(form, coefficient=None)

UFL form operator: 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.

ufl.formoperators.functional(form)

UFL form operator: Extract the functional part of form.

ufl.formoperators.lhs(form)

UFL form operator: 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.formoperators.rhs(form)

UFL form operator: 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.formoperators.sensitivity_rhs(a, u, L, v)

UFL form operator: 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

0 = d/dv[Ax-b] = dA/dv x + A dx/dv - db/dv, A dx/dv = db/dv - dA/dv x,

and solve this system for 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 db/dv and 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 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.formoperators.set_list_item(li, i, v)
ufl.formoperators.system(form)

UFL form operator: Split a form into the left hand side and right hand side, see lhs and rhs.

ufl.formoperators.zero_lists(shape)

geometry Module

Types for representing symbolic expressions for geometric quantities.

class ufl.geometry.CellCoordinate(domain)

Bases: ufl.geometry.GeometricCellQuantity

UFL geometry representation: The coordinate in a reference cell.

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

In the context of expression evaluation in a point in a cell, represents that point in the reference coordinate system of the cell.

is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

name = 'X'
shape()
class ufl.geometry.CellFacetJacobian(domain)

Bases: ufl.geometry.GeometricFacetQuantity

UFL geometry representation: The Jacobian of the mapping from reference facet to reference cell coordinates.

CFJ_ij = dX_i/dXf_j

is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

name = 'CFJ'
shape()
class ufl.geometry.CellFacetJacobianDeterminant(domain)

Bases: ufl.geometry.GeometricFacetQuantity

UFL geometry representation: The pseudo-determinant of the CellFacetJacobian.

is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

name = 'detCFJ'
class ufl.geometry.CellFacetJacobianInverse(domain)

Bases: ufl.geometry.GeometricFacetQuantity

UFL geometry representation: The pseudo-inverse of the CellFacetJacobian.

is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

name = 'CFK'
shape()
class ufl.geometry.CellFacetOrigin(domain)

Bases: ufl.geometry.GeometricFacetQuantity

UFL geometry representation: The reference cell coordinate corresponding to origin of a reference facet.

name = 'X0f'
shape()
class ufl.geometry.CellNormal(domain)

Bases: ufl.geometry.GeometricCellQuantity

UFL geometry representation: The upwards pointing normal vector of the current manifold cell.

name = 'cell_normal'
shape()
class ufl.geometry.CellOrientation(domain)

Bases: ufl.geometry.GeometricCellQuantity

UFL geometry representation: The orientation (+1/-1) of the current cell.

For non-manifold cells (tdim == gdim), this equals the sign of the Jacobian determinant, i.e. +1 if the physical cell is oriented the same way as the reference cell and -1 otherwise.

For manifold cells of tdim==gdim-1 this is input data belonging to the mesh, used to distinguish between the sides of the manifold.

name = 'cell_orientation'
class ufl.geometry.CellOrigin(domain)

Bases: ufl.geometry.GeometricCellQuantity

UFL geometry representation: The spatial coordinate corresponding to origin of a reference cell.

is_cellwise_constant()
name = 'x0'
shape()
class ufl.geometry.CellVolume(domain)

Bases: ufl.geometry.GeometricCellQuantity

UFL geometry representation: The volume of the cell.

name = 'volume'
class ufl.geometry.Circumradius(domain)

Bases: ufl.geometry.GeometricCellQuantity

UFL geometry representation: The circumradius of the cell.

name = 'circumradius'
class ufl.geometry.FacetArea(domain)

Bases: ufl.geometry.GeometricFacetQuantity

UFL geometry representation: The area of the facet.

name = 'facetarea'
class ufl.geometry.FacetCoordinate(domain)

Bases: ufl.geometry.GeometricFacetQuantity

UFL geometry representation: The coordinate in a reference cell of a facet.

In the context of expression integration over a facet, represents the reference facet coordinate of each quadrature point.

In the context of expression evaluation in a point on a facet, represents that point in the reference coordinate system of the facet.

is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

name = 'Xf'
shape()
class ufl.geometry.FacetJacobian(domain)

Bases: ufl.geometry.GeometricFacetQuantity

UFL geometry representation: The Jacobian of the mapping from reference facet to spatial coordinates.

FJ_ij = dx_i/dXf_j

The FacetJacobian is the product of the Jacobian and CellFacetJacobian:

FJ = dx/dXf = dx/dX dX/dXf = J * CFJ
is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

name = 'FJ'
shape()
class ufl.geometry.FacetJacobianDeterminant(domain)

Bases: ufl.geometry.GeometricFacetQuantity

UFL geometry representation: The pseudo-determinant of the FacetJacobian.

is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

name = 'detFJ'
class ufl.geometry.FacetJacobianInverse(domain)

Bases: ufl.geometry.GeometricFacetQuantity

UFL geometry representation: The pseudo-inverse of the FacetJacobian.

is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

name = 'FK'
shape()
class ufl.geometry.FacetNormal(domain)

Bases: ufl.geometry.GeometricFacetQuantity

UFL geometry representation: The outwards pointing normal vector of the current facet.

is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

name = 'n'
shape()
class ufl.geometry.FacetOrigin(domain)

Bases: ufl.geometry.GeometricFacetQuantity

UFL geometry representation: The spatial coordinate corresponding to origin of a reference facet.

name = 'x0f'
shape()
class ufl.geometry.GeometricCellQuantity(domain)

Bases: ufl.geometry.GeometricQuantity

class ufl.geometry.GeometricFacetQuantity(domain)

Bases: ufl.geometry.GeometricQuantity

class ufl.geometry.GeometricQuantity(domain)

Bases: ufl.terminal.Terminal

domains()
is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

shape()

Scalar shaped.

signature_data(domain_numbering)

Signature data of geometric quantities depend on the domain numbering.

class ufl.geometry.Jacobian(domain)

Bases: ufl.geometry.GeometricCellQuantity

UFL geometry representation: The Jacobian of the mapping from reference cell to spatial coordinates.

J_ij = dx_i/dX_j

is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

name = 'J'
shape()
class ufl.geometry.JacobianDeterminant(domain)

Bases: ufl.geometry.GeometricCellQuantity

UFL geometry representation: The determinant of the Jacobian.

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

is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

name = 'detJ'
class ufl.geometry.JacobianInverse(domain)

Bases: ufl.geometry.GeometricCellQuantity

UFL geometry representation: The inverse of the Jacobian.

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

is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

name = 'K'
shape()
class ufl.geometry.MaxFacetEdgeLength(domain)

Bases: ufl.geometry.GeometricFacetQuantity

UFL geometry representation: The maximum edge length of the facet.

name = 'maxfacetedgelength'
class ufl.geometry.MinFacetEdgeLength(domain)

Bases: ufl.geometry.GeometricFacetQuantity

UFL geometry representation: The minimum edge length of the facet.

name = 'minfacetedgelength'
class ufl.geometry.QuadratureWeight(domain)

Bases: ufl.geometry.GeometricQuantity

UFL geometry representation: The current quadrature weight.

Only used inside a quadrature context.

is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

name = 'weight'
class ufl.geometry.SpatialCoordinate(domain)

Bases: ufl.geometry.GeometricCellQuantity

UFL geometry representation: 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.

evaluate(x, mapping, component, index_values)
is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

name = 'x'
shape()

indexed Module

This module defines the Indexed class.

class ufl.indexed.Indexed(expression, indices)

Bases: ufl.operatorbase.WrapperType

evaluate(x, mapping, component, index_values, derivatives=())
free_indices()
index_dimensions()
is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

operands()
shape()

indexing Module

This module defines the single index types and some internal index utilities.

class ufl.indexing.FixedIndex(value)

Bases: ufl.indexing.IndexBase

UFL value: An index with a specific value assigned.

class ufl.indexing.Index(count=None)

Bases: ufl.indexing.IndexBase

UFL value: An index with no value assigned.

Used to represent free indices in Einstein indexing notation.

count()
class ufl.indexing.IndexBase

Bases: object

class ufl.indexing.MultiIndex(ii, idims)

Bases: ufl.terminal.UtilityType

Represents a sequence of indices, either fixed or free.

evaluate(x, mapping, component, index_values)
free_indices()
index_dimensions()
ufl.indexing.as_index(i)
ufl.indexing.as_multi_index(ii, shape=None)
ufl.indexing.indices(n)

UFL value: Return a tuple of n new Index objects.

indexsum Module

This module defines the IndexSum class.

class ufl.indexsum.IndexSum(summand, index)

Bases: ufl.operatorbase.AlgebraOperator

dimension()
evaluate(x, mapping, component, index_values)
free_indices()
index()
index_dimensions()
is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

operands()
shape()

indexutils Module

Some utilities for working with index tuples.

ufl.indexutils.complete_shape(shape, default_dim)

Complete shape tuple by replacing non-integers with a default dimension.

ufl.indexutils.repeated_indices(indices)

Return tuple of indices occuring more than once in input.

ufl.indexutils.shared_indices(ai, bi)

Return a tuple of indices occuring in both index tuples ai and bi.

ufl.indexutils.single_indices(indices)

Return a tuple of all indices occuring exactly once in input.

ufl.indexutils.unique_indices(indices)

Return a tuple of all indices from input, with no single index occuring more than once in output.

integral Module

The Integral class.

class ufl.integral.Integral(integrand, integral_type, domain, subdomain_id, metadata, subdomain_data)

Bases: object

An integral over a single domain.

domain()

Return the integration domain of this integral.

integral_type()

Return the domain type of this integral.

integrand()

Return the integrand expression, which is an Expr instance.

metadata()

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)

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()

Return the domain data of this integral.

subdomain_id()

Return the subdomain id of this integral.

log Module

This module provides functions used by the UFL implementation to output messages. These may be redirected by the user of UFL.

class ufl.log.Logger(name, exception_type=<type 'exceptions.Exception'>)

Create logger instance.

add_indent(increment=1)

Add to indentation level.

add_logfile(filename=None, mode='a', level=10)
begin(*message)

Begin task: write message and increase indentation level.

debug(*message)

Write debug message.

deprecate(*message)

Write deprecation message.

end()

End task: write a newline and decrease indentation level.

error(*message)

Write error message and raise an exception.

get_handler()

Get handler for logging.

get_logfile_handler(filename)
get_logger()

Return message logger.

info(*message)

Write info message.

info_blue(*message)

Write info message in blue.

info_green(*message)

Write info message in green.

info_red(*message)

Write info message in red.

log(level, *message)

Write a log message on given log level

pop_level()

Pop log level from the level stack, reverting to before the last push_level.

push_level(level)

Push a log level on the level stack.

set_handler(handler)

Replace handler for logging. To add additional handlers instead of replacing the existing, use log.get_logger().addHandler(myhandler). See the logging module for more details.

set_indent(level)

Set indentation level.

set_level(level)

Set log level.

set_prefix(prefix)

Set prefix for log messages.

warning(*message)

Write warning message.

warning_blue(*message)

Write warning message in blue.

warning_green(*message)

Write warning message in green.

warning_red(*message)

Write warning message in red.

mathfunctions Module

This module provides basic mathematical functions.

class ufl.mathfunctions.Acos(argument)

Bases: ufl.mathfunctions.MathFunction

class ufl.mathfunctions.Asin(argument)

Bases: ufl.mathfunctions.MathFunction

class ufl.mathfunctions.Atan(argument)

Bases: ufl.mathfunctions.MathFunction

class ufl.mathfunctions.Atan2(arg1, arg2)

Bases: ufl.operatorbase.Operator

evaluate(x, mapping, component, index_values)
free_indices()
index_dimensions()
operands()
shape()
class ufl.mathfunctions.BesselFunction(name, classname, nu, argument)

Bases: ufl.operatorbase.Operator

Base class for all bessel functions

evaluate(x, mapping, component, index_values)
free_indices()
index_dimensions()
operands()
shape()
class ufl.mathfunctions.BesselI(nu, argument)

Bases: ufl.mathfunctions.BesselFunction

class ufl.mathfunctions.BesselJ(nu, argument)

Bases: ufl.mathfunctions.BesselFunction

class ufl.mathfunctions.BesselK(nu, argument)

Bases: ufl.mathfunctions.BesselFunction

class ufl.mathfunctions.BesselY(nu, argument)

Bases: ufl.mathfunctions.BesselFunction

class ufl.mathfunctions.Cos(argument)

Bases: ufl.mathfunctions.MathFunction

class ufl.mathfunctions.Cosh(argument)

Bases: ufl.mathfunctions.MathFunction

class ufl.mathfunctions.Erf(argument)

Bases: ufl.mathfunctions.MathFunction

evaluate(x, mapping, component, index_values)
class ufl.mathfunctions.Exp(argument)

Bases: ufl.mathfunctions.MathFunction

class ufl.mathfunctions.Ln(argument)

Bases: ufl.mathfunctions.MathFunction

evaluate(x, mapping, component, index_values)
class ufl.mathfunctions.MathFunction(name, argument)

Bases: ufl.operatorbase.Operator

Base class for all math functions

evaluate(x, mapping, component, index_values)
free_indices()
index_dimensions()
operands()
shape()
class ufl.mathfunctions.Sin(argument)

Bases: ufl.mathfunctions.MathFunction

class ufl.mathfunctions.Sinh(argument)

Bases: ufl.mathfunctions.MathFunction

class ufl.mathfunctions.Sqrt(argument)

Bases: ufl.mathfunctions.MathFunction

class ufl.mathfunctions.Tan(argument)

Bases: ufl.mathfunctions.MathFunction

class ufl.mathfunctions.Tanh(argument)

Bases: ufl.mathfunctions.MathFunction

measure Module

The Measure class.

class ufl.measure.Measure(integral_type, domain=None, subdomain_id='everywhere', metadata=None, subdomain_data=None)

Bases: object

integral_type:
str, one of “cell”, etc., or short form “dx”, etc.
domain:
a Domain object (includes cell, dims, label, domain data)
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.
CELL = 'cell'
EXTERIOR_FACET = 'exterior_facet'
EXTERIOR_FACET_BOTTOM = 'exterior_facet_bottom'
EXTERIOR_FACET_TOP = 'exterior_facet_top'
EXTERIOR_FACET_VERT = 'exterior_facet_vert'
INTERIOR_FACET = 'interior_facet'
INTERIOR_FACET_HORIZ = 'interior_facet_horiz'
INTERIOR_FACET_VERT = 'interior_facet_vert'
MACRO_CELL = 'macro_cell'
POINT = 'point'
QUADRATURE = 'quadrature'
SURFACE = 'surface'
domain()

Return the domain associated with this measure.

This may be None or a Domain object.

integral_type()

Return the domain type.

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

metadata()

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)

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()

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()

Return the domain id of this measure (integer).

class ufl.measure.MeasureProduct(*measures)

Bases: object

Represents a product of measures.

This is a notational intermediate object to handle the notation

f*(dm1*dm2)

This is work in progress and not functional. It needs support in other parts of ufl and the rest of the code generation chain.

Create MeasureProduct from given list of measures.

sub_measures()

Return submeasures.

class ufl.measure.MeasureSum(*measures)

Bases: object

Represents a sum of measures.

This is a notational intermediate object to translate the notation

f*(ds(1)+ds(3))

into

f*ds(1) + f*ds(3)
ufl.measure.as_integral_type(integral_type)

Map short name to long name and require a valid one.

ufl.measure.integral_types()

Return a tuple of all domain type strings.

ufl.measure.measure_names()

Return a tuple of all measure name strings.

ufl.measure.register_integral_type(integral_type, measure_name)

objects Module

Utility objects for pretty syntax in user code.

operatorbase Module

Base class for all operators, i.e. non-terminal expr types.

class ufl.operatorbase.AlgebraOperator

Bases: ufl.operatorbase.Operator

class ufl.operatorbase.Operator

Bases: ufl.expr.Expr

is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

reconstruct(*operands)

Return a new object of the same type with new operands.

signature_data()
class ufl.operatorbase.WrapperType

Bases: ufl.operatorbase.Operator

free_indices()
index_dimensions()
shape()
ufl.operatorbase.compute_hash(expr)
ufl.operatorbase.compute_hash_with_stats(expr)

operators Module

This module extends the form language with free function operators, which are either already available as member functions on UFL objects or defined as compound operators involving basic operations on the UFL objects.

ufl.operators.And(left, right)

UFL operator: A boolean expresion (left and right) for use with conditional.

ufl.operators.Dn(f)

UFL operator: Take the directional derivative of f in the facet normal direction, Dn(f) := dot(grad(f), n).

ufl.operators.Dt(f)

UFL operator: <Not implemented yet!> The partial derivative of f with respect to time.

ufl.operators.Dx(f, *i)

UFL operator: Take the partial derivative of f with respect to spatial variable number i. Equivalent to f.dx(*i).

ufl.operators.Max(x, y)

UFL operator: Take the maximum of x and y.

ufl.operators.Min(x, y)

UFL operator: Take the minimum of x and y.

ufl.operators.Not(condition)

UFL operator: A boolean expresion (not condition) for use with conditional.

ufl.operators.Or(left, right)

UFL operator: A boolean expresion (left or right) for use with conditional.

ufl.operators.acos(f)

UFL operator: Take the inverse cosinus of f.

ufl.operators.asin(f)

UFL operator: Take the inverse sinus of f.

ufl.operators.atan(f)

UFL operator: Take the inverse tangent of f.

ufl.operators.atan_2(f1, f2)

UFL operator: Take the inverse tangent of f.

ufl.operators.avg(v)

UFL operator: Take the average of v across a facet.

ufl.operators.bessel_I(nu, f)

UFL operator: regular modified cylindrical Bessel function.

ufl.operators.bessel_J(nu, f)

UFL operator: cylindrical Bessel function of the first kind.

ufl.operators.bessel_K(nu, f)

UFL operator: irregular modified cylindrical Bessel function.

ufl.operators.bessel_Y(nu, f)

UFL operator: cylindrical Bessel function of the second kind.

ufl.operators.cell_avg(f)

UFL operator: Take the average of v over a cell.

ufl.operators.cofac(A)

UFL operator: Take the cofactor of A.

ufl.operators.conditional(condition, true_value, false_value)

UFL operator: A conditional expression, taking the value of true_value when condition evaluates to true and false_value otherwise.

ufl.operators.contraction(a, a_axes, b, b_axes)

UFL operator: Take the contraction of a and b over given axes.

ufl.operators.cos(f)

UFL operator: Take the cosinus of f.

ufl.operators.cosh(f)

UFL operator: Take the cosinus hyperbolicus of f.

ufl.operators.cross(a, b)

UFL operator: Take the cross product of a and b.

ufl.operators.curl(f)

UFL operator: Take the curl of f.

ufl.operators.det(A)

UFL operator: Take the determinant of A.

ufl.operators.dev(A)

UFL operator: Take the deviatoric part of A.

ufl.operators.diag(A)

UFL operator: 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.operators.diag_vector(A)

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

See also diag.

ufl.operators.diff(f, v)

UFL operator: Take the derivative of f with respect to the variable v.

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

ufl.operators.div(f)

UFL operator: 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.operators.dot(a, b)

UFL operator: Take the dot product of a and b.

ufl.operators.elem_div(A, B)

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

ufl.operators.elem_mult(A, B)

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

ufl.operators.elem_op(op, *args)

UFL operator: Take the elementwise application of operator op on scalar values from one or more tensor arguments.

ufl.operators.elem_op_items(op_ind, indices, *args)
ufl.operators.elem_pow(A, B)

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

ufl.operators.eq(left, right)

UFL operator: A boolean expresion (left == right) for use with conditional.

ufl.operators.erf(f)

UFL operator: Take the error function of f.

ufl.operators.exp(f)

UFL operator: Take the exponential of f.

ufl.operators.exterior_derivative(f)

UFL operator: Take the exterior derivative of f.

The exterior derivative uses the element family 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.operators.facet_avg(f)

UFL operator: Take the average of v over a facet.

ufl.operators.ge(left, right)

UFL operator: A boolean expresion (left >= right) for use with conditional.

ufl.operators.grad(f)

UFL operator: 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.operators.gt(left, right)

UFL operator: A boolean expresion (left > right) for use with conditional.

ufl.operators.inner(a, b)

UFL operator: Take the inner product of a and b.

ufl.operators.inv(A)

UFL operator: Take the inverse of A.

ufl.operators.jump(v, n=None)

UFL operator: Take the jump of v across a facet.

ufl.operators.le(left, right)

UFL operator: A boolean expresion (left <= right) for use with conditional.

ufl.operators.ln(f)

UFL operator: Take the natural logarithm of f.

ufl.operators.lt(left, right)

UFL operator: A boolean expresion (left < right) for use with conditional.

ufl.operators.nabla_div(f)

UFL operator: 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.operators.nabla_grad(f)

UFL operator: 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.operators.ne(left, right)

UFL operator: A boolean expresion (left != right) for use with conditional.

ufl.operators.outer(*operands)

UFL operator: Take the outer product of two or more operands.

ufl.operators.perp(v)

UFL operator: Take the perp of v, i.e. (-v1, +v0).

ufl.operators.rank(f)

UFL operator: The rank of f.

ufl.operators.rot(f)

UFL operator: Take the curl of f.

ufl.operators.shape(f)

UFL operator: The shape of f.

ufl.operators.sign(x)

UFL operator: Take the sign (+1 or -1) of x.

ufl.operators.sin(f)

UFL operator: Take the sinus of f.

ufl.operators.sinh(f)

UFL operator: Take the sinus hyperbolicus of f.

ufl.operators.skew(A)

UFL operator: Take the skew symmetric part of A.

ufl.operators.sqrt(f)

UFL operator: Take the square root of f.

ufl.operators.sym(A)

UFL operator: Take the symmetric part of A.

ufl.operators.tan(f)

UFL operator: Take the tangent of f.

ufl.operators.tanh(f)

UFL operator: Take the tangent hyperbolicus of f.

ufl.operators.tr(A)

UFL operator: Take the trace of A.

ufl.operators.transpose(A)

UFL operator: Take the transposed of tensor A.

ufl.operators.variable(e)

UFL operator: Define a variable representing the given expression, see also diff().

permutation Module

This module provides utility functions for computing permutations and generating index lists.

ufl.permutation.build_component_numbering(shape, symmetry)

Build a numbering of components within the given value shape, taking into consideration a symmetry mapping which leaves the mapping noncontiguous. Returns a dict { component -> numbering } and an ordered list of components [ numbering -> component ]. The dict contains all components while the list only contains the ones not mapped by the symmetry mapping.

ufl.permutation.compute_indices(shape)

Compute all index combinations for given shape

ufl.permutation.compute_indices2(shape)

Compute all index combinations for given shape

ufl.permutation.compute_order_tuples(k, n)

Compute all tuples of n integers such that the sum is k

ufl.permutation.compute_permutation_pairs(j, k)

Compute all permutations of j + k elements from (0, j + k) in rising order within (0, j) and (j, j + k) respectively.

ufl.permutation.compute_permutations(k, n, skip=None)

Compute all permutations of k elements from (0, n) in rising order. Any elements that are contained in the list skip are not included.

ufl.permutation.compute_sign(permutation)

Compute sign by sorting.

precedence Module

Precedence handling.

ufl.precedence.assign_precedences(precedence_list)

Given a precedence list, assign ints to class._precedence.

ufl.precedence.build_precedence_list()
ufl.precedence.build_precedence_mapping(precedence_list)

Given a precedence list, build a dict with class->int mappings. Utility function used by some external code.

ufl.precedence.parstr(child, parent, pre='(', post=')', format=<type 'str'>)

protocols Module

ufl.protocols.id_or_none(obj)

Returns None if the object is None, obj.ufl_id() if available, or id(obj) if not.

This allows external libraries to implement an alternative to id(obj) in the ufl_id() function, such that ufl can identify objects as the same without knowing about their types.

ufl.protocols.metadata_equal(a, b)
ufl.protocols.metadata_hashdata(md)

pullbackof Module

restriction Module

Restriction operations.

class ufl.restriction.CellAvg(f)

Bases: ufl.operatorbase.Operator

evaluate(x, mapping, component, index_values)

Performs an approximate symbolic evaluation, since we dont have a cell.

free_indices()
index_dimensions()
operands()
shape()
class ufl.restriction.FacetAvg(f)

Bases: ufl.operatorbase.Operator

evaluate(x, mapping, component, index_values)

Performs an approximate symbolic evaluation, since we dont have a cell.

free_indices()
index_dimensions()
operands()
shape()
class ufl.restriction.NegativeRestricted(f)

Bases: ufl.restriction.Restricted

class ufl.restriction.PositiveRestricted(f)

Bases: ufl.restriction.Restricted

class ufl.restriction.Restricted(f, side)

Bases: ufl.operatorbase.Operator

evaluate(x, mapping, component, index_values)
free_indices()
index_dimensions()
operands()
shape()

sobolevspace Module

This module defines a symbolic heirarchy of Sobolev spaces to enable symbolic reasoning about the spaces in which finite elements lie.

class ufl.sobolevspace.SobolevSpace(name, parents=None)

Bases: object

Symbolic representation of a Sobolev space. This implements a subset of the methods of a Python set so that finite elements and other Sobolev spaces can be tested for inclusion.

Instantiate a SobolevSpace object.

Parameters:
  • name – The name of this space,
  • parents – A set of Sobolev spaces of which this

space is a subspace.

sorting Module

This module contains a sorting rule for expr objects that is more robust w.r.t. argument numbering than using repr.

class ufl.sorting.ExprKey(x)

Bases: object

x
ufl.sorting.cmp_expr(a, b)

Sorting rule for Expr objects. NB! Do not use to compare for equality!

ufl.sorting.expr_key(expr)
ufl.sorting.sorted_expr(seq)
ufl.sorting.sorted_expr_sum(seq)

split_functions Module

Algorithm for splitting a Coefficient or Argument into subfunctions.

ufl.split_functions.split(v)

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

tensoralgebra Module

Compound tensor algebra operations.

class ufl.tensoralgebra.Cofactor(A)

Bases: ufl.tensoralgebra.CompoundTensorOperator

free_indices()
index_dimensions()
operands()
shape()
class ufl.tensoralgebra.CompoundTensorOperator

Bases: ufl.operatorbase.AlgebraOperator

class ufl.tensoralgebra.Cross(a, b)

Bases: ufl.tensoralgebra.CompoundTensorOperator

free_indices()
index_dimensions()
operands()
shape()
class ufl.tensoralgebra.Determinant(A)

Bases: ufl.tensoralgebra.CompoundTensorOperator

free_indices()
index_dimensions()
operands()
shape()
class ufl.tensoralgebra.Deviatoric(A)

Bases: ufl.tensoralgebra.CompoundTensorOperator

free_indices()
index_dimensions()
operands()
shape()
class ufl.tensoralgebra.Dot(a, b)

Bases: ufl.tensoralgebra.CompoundTensorOperator

free_indices()
index_dimensions()
operands()
shape()
class ufl.tensoralgebra.Inner(a, b)

Bases: ufl.tensoralgebra.CompoundTensorOperator

free_indices()
index_dimensions()
operands()
shape()
class ufl.tensoralgebra.Inverse(A)

Bases: ufl.tensoralgebra.CompoundTensorOperator

free_indices()
index_dimensions()
operands()
shape()
class ufl.tensoralgebra.Outer(a, b)

Bases: ufl.tensoralgebra.CompoundTensorOperator

free_indices()
index_dimensions()
operands()
shape()
class ufl.tensoralgebra.Skew(A)

Bases: ufl.tensoralgebra.CompoundTensorOperator

free_indices()
index_dimensions()
operands()
shape()
class ufl.tensoralgebra.Sym(A)

Bases: ufl.tensoralgebra.CompoundTensorOperator

free_indices()
index_dimensions()
operands()
shape()
class ufl.tensoralgebra.Trace(A)

Bases: ufl.tensoralgebra.CompoundTensorOperator

free_indices()
index_dimensions()
operands()
shape()
class ufl.tensoralgebra.Transposed(A)

Bases: ufl.tensoralgebra.CompoundTensorOperator

free_indices()
index_dimensions()
operands()
shape()
ufl.tensoralgebra.merge_indices(a, b)

tensors Module

Classes used to group scalar expressions into expressions with rank > 0.

class ufl.tensors.ComponentTensor(expression, indices)

Bases: ufl.operatorbase.WrapperType

UFL operator type: Maps the free indices of a scalar valued expression to tensor axes.

evaluate(x, mapping, component, index_values)
free_indices()
index_dimensions()
indices()
is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

operands()
reconstruct(expressions, indices)
shape()
class ufl.tensors.ListTensor(*expressions)

Bases: ufl.operatorbase.WrapperType

UFL operator type: Wraps a list of expressions into a tensor valued expression of one higher rank.

evaluate(x, mapping, component, index_values, derivatives=())
free_indices()
index_dimensions()
is_cellwise_constant()

Return whether this expression is spatially constant over each cell.

operands()
shape()
ufl.tensors.as_matrix(expressions, indices=None)

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

ufl.tensors.as_scalar(expression)

Given a scalar or tensor valued expression A, returns either of the tuples:

(a,b) = (A, ())
(a,b) = (A[indices], indices)

such that a is always a scalar valued expression.

ufl.tensors.as_tensor(expressions, indices=None)

UFL operator: 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.tensors.as_vector(expressions, index=None)

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

ufl.tensors.dyad(d, *iota)

TODO: Develop this concept, can e.g. write A[i,j]*dyad(j,i) for the transpose.

ufl.tensors.from_numpy_to_lists(expressions)
ufl.tensors.numpy2nestedlists(arr)
ufl.tensors.relabel(A, indexmap)

UFL operator: Relabel free indices of A with new indices, using the given mapping.

ufl.tensors.unit_indexed_tensor(shape, component)
ufl.tensors.unit_list(i, n)
ufl.tensors.unit_list2(i, j, n)
ufl.tensors.unit_matrices(d)

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

ufl.tensors.unit_matrix(i, j, d)

UFL value: A constant unit matrix in direction i,j with dimension d.

ufl.tensors.unit_vector(i, d)

UFL value: A constant unit vector in direction i with dimension d.

ufl.tensors.unit_vectors(d)

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

ufl.tensors.unwrap_list_tensor(lt)

terminal Module

This module defines the Terminal class, the superclass for all types that are terminal nodes in the expression trees.

class ufl.terminal.FormArgument

Bases: ufl.terminal.Terminal

class ufl.terminal.Terminal

Bases: ufl.expr.Expr

A terminal node in the UFL expression tree.

domains()

Return tuple of domains related to this terminal object.

evaluate(x, mapping, component, index_values, derivatives=())

Get self from mapping and return the component asked for.

free_indices()

A Terminal object never has free indices.

index_dimensions()

A Terminal object never has free indices.

operands()

A Terminal object never has operands.

reconstruct(*operands)

Return self.

signature_data()

Default signature data for of terminals just return the repr string.

class ufl.terminal.UtilityType

Bases: ufl.terminal.Terminal

domains()

Return tuple of domains related to this terminal object.

free_indices()
index_dimensions()
is_cellwise_constant()
shape()

testobjects Module

Some premade objects useful for quick testing.

variable Module

Defines the Variable and Label classes, used to label expressions as variables for differentiation.

class ufl.variable.Label(count=None)

Bases: ufl.terminal.UtilityType

count()
class ufl.variable.Variable(expression, label=None)

Bases: ufl.operatorbase.WrapperType

A Variable is a representative for another expression.

It will be used by the end-user mainly for defining a quantity to differentiate w.r.t. using diff. Example:

e = <...>
e = variable(e)
f = exp(e**2)
df = diff(f, e)
domains()
evaluate(x, mapping, component, index_values)
expression()
free_indices()
index_dimensions()
is_cellwise_constant()
label()
operands()
shape()