firedrake.slate package

Subpackages

Submodules

firedrake.slate.slate module

Slate is a symbolic language defining a framework for performing linear algebra operations on finite element tensors. It is similar in principle to most linear algebra libraries in notation.

The design of Slate was heavily influenced by UFL, and utilizes much of UFL’s functionality for FEM-specific form manipulation.

Unlike UFL, however, once forms are assembled into Slate Tensor objects, one can utilize the operations defined in Slate to express complicated linear algebra operations (such as the Schur-complement reduction of a block-matrix system).

All Slate expressions are handled by a specialized linear algebra compiler, which interprets expressions and produces C++ kernel functions to be executed within the Firedrake architecture.

class firedrake.slate.slate.Add(A, B)[source]

Bases: BinaryOp

Abstract Slate class representing matrix-matrix, vector-vector

or scalar-scalar addition.

Parameters:
  • A – a TensorBase object.

  • B – another TensorBase object.

Constructor for the Add class.

arg_function_spaces[source]

Returns a tuple of function spaces that the tensor is defined on.

arguments()[source]

Returns a tuple of arguments associated with the tensor.

prec = 1
class firedrake.slate.slate.AssembledVector(function)[source]

Bases: TensorBase

This class is a symbolic representation of an assembled vector of data contained in a Function.

Parameters:

function – A firedrake function.

Initialise a cache for stashing results.

Mirrors Form.

arg_function_spaces[source]

Returns a tuple of function spaces that the tensor is defined on.

arguments()[source]

Returns a tuple of arguments associated with the tensor.

assembled = True
coefficients()[source]

Returns a tuple of coefficients associated with the tensor.

constants()[source]

Returns a tuple of constants associated with the tensor.

form[source]
property integrals
operands = ()
prec = 0
slate_coefficients()[source]

Returns a tuple of coefficients associated with the tensor.

subdomain_data()[source]

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

terminal = True
ufl_domains()[source]

Returns the integration domains of the integrals associated with the tensor.

class firedrake.slate.slate.Block(tensor, indices)[source]

Bases: TensorBase

This class represents a tensor corresponding to particular block of a mixed tensor. Depending on the indices provided, the subblocks can span multiple test/trial spaces.

Parameters:
  • tensor – A (mixed) tensor.

  • indices – Indices of the test and trial function spaces to extract. This should be a 0-, 1-, or 2-tuple (whose length is equal to the rank of the tensor.) The entries should be an iterable of integer indices.

For example, consider the mixed tensor defined by:

n = FacetNormal(m)
U = FunctionSpace(m, "DRT", 1)
V = FunctionSpace(m, "DG", 0)
M = FunctionSpace(m, "DGT", 0)
W = U * V * M
u, p, r = TrialFunctions(W)
w, q, s = TestFunctions(W)
A = Tensor(dot(u, w)*dx + p*div(w)*dx + r*dot(w, n)*dS
           + div(u)*q*dx + p*q*dx + r*s*ds)

This describes a block 3x3 mixed tensor of the form:

\[\begin{split}\begin{bmatrix} A & B & C \\ D & E & F \\ G & H & J \end{bmatrix}\end{split}\]

Providing the 2-tuple ((0, 1), (0, 1)) returns a tensor corresponding to the upper 2x2 block:

\[\begin{split}\begin{bmatrix} A & B \\ D & E \end{bmatrix}\end{split}\]

More generally, argument indices of the form (idr, idc) produces a tensor of block-size len(idr) x len(idc) spanning the specified test/trial spaces.

Constructor for the Block class.

arg_function_spaces[source]

Returns a tuple of function spaces that the tensor is defined on.

arguments()[source]

Returns a tuple of arguments associated with the tensor.

assembled[source]
coefficients()[source]

Returns a tuple of coefficients associated with the tensor.

constants()[source]

Returns a tuple of constants associated with the tensor.

form[source]
prec = 0
slate_coefficients()[source]

Returns a tuple of coefficients associated with the tensor.

subdomain_data()[source]

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

terminal[source]

Blocks are only terminal when they sit on Tensors or AssembledVectors

ufl_domains()[source]

Returns the integration domains of the integrals associated with the tensor.

class firedrake.slate.slate.BlockAssembledVector(function, expr, indices)[source]

Bases: AssembledVector

This class is a symbolic representation of an assembled vector of data contained in a set of Function s defined on pieces of a split mixed function space.

Parameters:

functions – A tuple of firedrake functions.

Initialise a cache for stashing results.

Mirrors Form.

arg_function_spaces[source]

Returns a tuple of function spaces associated to the corresponding block.

arguments()[source]

Returns a tuple of arguments associated with the corresponding block.

coefficients()[source]

Returns a tuple of coefficients associated with the tensor.

form[source]
slate_coefficients()[source]

Returns a BlockFunction in a tuple which carries all information to generate the right coefficients and maps.

subdomain_data()[source]

Returns mappings on the tensor: {domain:{integral_type: subdomain_data}}.

ufl_domains()[source]

Returns the integration domains of the integrals associated with the tensor.

class firedrake.slate.slate.DiagonalTensor(A)[source]

Bases: UnaryOp

An abstract Slate class representing the diagonal of a tensor.

Warning

This class will raise an error if the tensor is not square.

Constructor for the Diagonal class.

arg_function_spaces[source]

Returns a tuple of function spaces that the tensor is defined on.

arguments()[source]

Returns a tuple of arguments associated with the tensor.

diagonal = True
prec = 0
class firedrake.slate.slate.Factorization(tensor, decomposition=None)[source]

Bases: TensorBase

An abstract Slate class for the factorization of matrices. The factorizations available are the following:

  1. LU with full or partial pivoting (‘FullPivLU’ and ‘PartialPivLU’);

  2. QR using Householder reflectors (‘HouseholderQR’) with the option to use column pivoting (‘ColPivHouseholderQR’) or full pivoting (‘FullPivHouseholderQR’);

  3. standard Cholesky (‘LLT’) and stabilized Cholesky factorizations with pivoting (‘LDLT’);

  4. a rank-revealing complete orthogonal decomposition using Householder transformations (‘CompleteOrthogonalDecomposition’); and

  5. singular-valued decompositions (‘JacobiSVD’ and ‘BDCSVD’). For larger matrices, ‘BDCSVD’ is recommended.

Constructor for the Factorization class.

arg_function_spaces[source]

Returns a tuple of function spaces that the tensor is defined on.

arguments()[source]

Returns a tuple of arguments associated with the tensor.

coefficients()[source]

Returns a tuple of coefficients associated with the tensor.

constants()[source]

Returns a tuple of constants associated with the tensor.

prec = 0
slate_coefficients()[source]

Returns a tuple of coefficients associated with the tensor.

subdomain_data()[source]

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

ufl_domains()[source]

Returns the integration domains of the integrals associated with the tensor.

class firedrake.slate.slate.Inverse(A)[source]

Bases: UnaryOp

An abstract Slate class representing the inverse of a tensor.

Warning

This class will raise an error if the tensor is not square.

Constructor for the Inverse class.

arg_function_spaces[source]

Returns a tuple of function spaces that the tensor is defined on.

arguments()[source]

Returns the expected arguments of the resulting tensor of performing a specific unary operation on a tensor.

class firedrake.slate.slate.Mul(A, B)[source]

Bases: BinaryOp

Abstract Slate class representing the interior product or two tensors. By interior product, we mean an operation that results in a tensor of equal or lower rank via performing a contraction on arguments. This includes Matrix-Matrix and Matrix-Vector multiplication.

Parameters:
  • A – a TensorBase object.

  • B – another TensorBase object.

Constructor for the Mul class.

arg_function_spaces[source]

Returns a tuple of function spaces that the tensor is defined on.

arguments()[source]

Returns the arguments of a tensor resulting from multiplying two tensors A and B.

prec = 2
class firedrake.slate.slate.Negative(*operands)[source]

Bases: UnaryOp

Abstract Slate class representing the negation of a tensor object.

Constructor for the TensorOp class.

arg_function_spaces[source]

Returns a tuple of function spaces that the tensor is defined on.

arguments()[source]

Returns the expected arguments of the resulting tensor of performing a specific unary operation on a tensor.

class firedrake.slate.slate.Reciprocal(A)[source]

Bases: UnaryOp

An abstract Slate class representing the reciprocal of a vector.

Constructor for the Inverse class.

arg_function_spaces[source]

Returns a tuple of function spaces that the tensor is defined on.

arguments()[source]

Returns the expected arguments of the resulting tensor of performing a specific unary operation on a tensor.

prec = 0
class firedrake.slate.slate.Solve(A, B, decomposition=None)[source]

Bases: BinaryOp

Abstract Slate class describing a local linear system of equations. This object is a direct solver, utilizing the application of the inverse of matrix in a decomposed form.

Parameters:
  • A – The left-hand side operator.

  • B – The right-hand side.

  • decomposition – A string denoting the type of matrix decomposition to used. The factorizations available are detailed in the Factorization documentation.

Constructor for the Solve class.

arg_function_spaces[source]

Returns a tuple of function spaces that the tensor is defined on.

arguments()[source]

Returns the arguments of a tensor resulting from applying the inverse of A onto B.

prec = 3
class firedrake.slate.slate.Tensor(form, diagonal=False)[source]

Bases: TensorBase

This class is a symbolic representation of a finite element tensor derived from a bilinear or linear form. This class implements all supported ranks of general tensor (rank-0, rank-1 and rank-2 tensor objects). This class is the primary user-facing class that the Slate symbolic algebra supports.

Parameters:

form – a ufl.Form object.

A ufl.Form is currently the only supported input of creating a slate.Tensor object:

  1. If the form is a bilinear form, namely a form with two ufl.Argument objects, then the Slate Tensor will be a rank-2 Matrix.

  2. If the form has one ufl.Argument as in the case of a typical linear form, then this will create a rank-1 Vector.

  3. A zero-form will create a rank-0 Scalar.

These are all under the same type slate.Tensor. The attribute self.rank is used to determine what kind of tensor object is being handled.

Constructor for the Tensor class.

arg_function_spaces[source]

Returns a tuple of function spaces that the tensor is defined on.

arguments()[source]

Returns a tuple of arguments associated with the tensor.

coefficients()[source]

Returns a tuple of coefficients associated with the tensor.

constants()[source]

Returns a tuple of constants associated with the tensor.

operands = ()
prec = 0
slate_coefficients()[source]

Returns a tuple of coefficients associated with the tensor.

subdomain_data()[source]

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

terminal = True
ufl_domains()[source]

Returns the integration domains of the integrals associated with the tensor.

class firedrake.slate.slate.Transpose(*operands)[source]

Bases: UnaryOp

An abstract Slate class representing the transpose of a tensor.

Constructor for the TensorOp class.

arg_function_spaces[source]

Returns a tuple of function spaces that the tensor is defined on.

arguments()[source]

Returns the expected arguments of the resulting tensor of performing a specific unary operation on a tensor.

Module contents