firedrake.slate package¶
Subpackages¶
- firedrake.slate.slac package
- Submodules
- firedrake.slate.slac.compiler module
- firedrake.slate.slac.kernel_builder module
CellFacetKernelArg
CoefficientInfo
IndexCreator
LayerCountKernelArg
LocalLoopyKernelBuilder
LocalLoopyKernelBuilder.cell_facets_arg_name
LocalLoopyKernelBuilder.cell_orientations_arg_name
LocalLoopyKernelBuilder.cell_sizes_arg_name
LocalLoopyKernelBuilder.collect_coefficients()
LocalLoopyKernelBuilder.collect_constants()
LocalLoopyKernelBuilder.collect_tsfc_kernel_data()
LocalLoopyKernelBuilder.coordinates_arg_name
LocalLoopyKernelBuilder.extent()
LocalLoopyKernelBuilder.facet_integral_predicates()
LocalLoopyKernelBuilder.generate_lhs()
LocalLoopyKernelBuilder.generate_tsfc_calls()
LocalLoopyKernelBuilder.generate_wrapper_kernel_args()
LocalLoopyKernelBuilder.initialise_terminals()
LocalLoopyKernelBuilder.is_integral_type()
LocalLoopyKernelBuilder.layer_arg_name
LocalLoopyKernelBuilder.layer_count_name
LocalLoopyKernelBuilder.layer_integral_predicates()
LocalLoopyKernelBuilder.local_facet_array_arg_name
LocalLoopyKernelBuilder.loopify_tsfc_kernel_data()
LocalLoopyKernelBuilder.shape()
LocalLoopyKernelBuilder.slate_call()
LocalLoopyKernelBuilder.supported_integral_types
LocalLoopyKernelBuilder.supported_subdomain_types
LocalLoopyKernelBuilder.tsfc_cxt_kernels()
SlateWrapperBag
- firedrake.slate.slac.optimise module
- firedrake.slate.slac.tsfc_driver module
- firedrake.slate.slac.utils module
- Module contents
- firedrake.slate.static_condensation package
- Submodules
- firedrake.slate.static_condensation.hybridization module
- firedrake.slate.static_condensation.la_utils module
- firedrake.slate.static_condensation.sc_base module
- firedrake.slate.static_condensation.scpc module
- Module contents
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.
- property arg_function_spaces¶
Returns a tuple of function spaces that the tensor is defined on.
- 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
.- property arg_function_spaces¶
Returns a tuple of function spaces that the tensor is defined on.
- assembled = True¶
- property form¶
- property integrals¶
- operands = ()¶
- prec = 0¶
- subdomain_data()[source]¶
Returns a mapping on the tensor:
{domain:{integral_type: subdomain_data}}
.
- terminal = True¶
- 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.
- property arg_function_spaces¶
Returns a tuple of function spaces that the tensor is defined on.
- property assembled¶
bool(x) -> bool
Returns True when the argument x is true, False otherwise. The builtins True and False are the only two instances of the class bool. The class bool is a subclass of the class int, and cannot be subclassed.
- property form¶
- prec = 0¶
- subdomain_data()[source]¶
Returns a mapping on the tensor:
{domain:{integral_type: subdomain_data}}
.
- property terminal¶
Blocks are only terminal when they sit on Tensors or AssembledVectors
- 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
.- property arg_function_spaces¶
Returns a tuple of function spaces associated to the corresponding block.
- property form¶
- 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.
- property arg_function_spaces¶
Returns a tuple of function spaces that the tensor is defined on.
- 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:
LU with full or partial pivoting (‘FullPivLU’ and ‘PartialPivLU’);
QR using Householder reflectors (‘HouseholderQR’) with the option to use column pivoting (‘ColPivHouseholderQR’) or full pivoting (‘FullPivHouseholderQR’);
standard Cholesky (‘LLT’) and stabilized Cholesky factorizations with pivoting (‘LDLT’);
a rank-revealing complete orthogonal decomposition using Householder transformations (‘CompleteOrthogonalDecomposition’); and
singular-valued decompositions (‘JacobiSVD’ and ‘BDCSVD’). For larger matrices, ‘BDCSVD’ is recommended.
Constructor for the Factorization class.
- property arg_function_spaces¶
Returns a tuple of function spaces that the tensor is defined on.
- prec = 0¶
- 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.
- property arg_function_spaces¶
Returns a tuple of function spaces that the tensor is defined on.
- 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.
- property arg_function_spaces¶
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.
- property arg_function_spaces¶
Returns a tuple of function spaces that the tensor is defined on.
- class firedrake.slate.slate.Reciprocal(A)[source]¶
Bases:
UnaryOp
An abstract Slate class representing the reciprocal of a vector.
Constructor for the Inverse class.
- property arg_function_spaces¶
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.
- property arg_function_spaces¶
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: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.If the form has one ufl.Argument as in the case of a typical linear form, then this will create a rank-1 Vector.
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.
- property arg_function_spaces¶
Returns a tuple of function spaces that the tensor is defined on.
- operands = ()¶
- prec = 0¶
- subdomain_data()[source]¶
Returns a mapping on the tensor:
{domain:{integral_type: subdomain_data}}
.
- terminal = True¶