Warning
You are reading a version of the website built against the unstable main branch. This content is liable to change without notice and may be inappropriate for your use case.
You can find the documentation for the current stable release here.
firedrake.ensemble package¶
Submodules¶
firedrake.ensemble.ensemble module¶
- class firedrake.ensemble.ensemble.Ensemble(comm, M, **kwargs)[source]¶
Bases:
objectCreate a set of space and ensemble subcommunicators.
- Parameters:
comm – The communicator to split.
M – the size of the communicators used for spatial parallelism.
- Keyword Arguments:
ensemble_name – string used as communicator name prefix, for debugging.
- Raises:
ValueError – if
Mdoes not dividecomm.sizeexactly.
- allreduce(f, f_reduced, op=<mpi4py.MPI.Op object>)[source]¶
Allreduce a function f into f_reduced over
ensemble_comm.- Parameters:
f – The a
Functionto allreduce.f_reduced – the result of the reduction.
op – MPI reduction operator. Defaults to MPI.SUM.
- Raises:
ValueError – if function communicators mismatch each other or the ensemble spatial communicator, or if the functions are in different spaces
- bcast(f, root=0)[source]¶
Broadcast a function f over
ensemble_commfrom rank root- Parameters:
f – The
Functionto broadcast.root – rank to broadcast from. Defaults to 0.
- Raises:
ValueError – if function communicator mismatches the ensemble spatial communicator.
- property ensemble_rank¶
The rank of the local ensemble member.
- property ensemble_size¶
The number of ensemble members.
- iallreduce(f, f_reduced, op=<mpi4py.MPI.Op object>)[source]¶
Allreduce (non-blocking) a function f into f_reduced over
ensemble_comm.- Parameters:
f – The a
Functionto allreduce.f_reduced – the result of the reduction.
op – MPI reduction operator. Defaults to MPI.SUM.
- Returns:
list of MPI.Request objects (one for each of f.subfunctions).
- Raises:
ValueError – if function communicators mismatch each other or the ensemble spatial communicator, or if the functions are in different spaces
- ibcast(f, root=0)[source]¶
Broadcast (non-blocking) a function f over
ensemble_commfrom rank root- Parameters:
f – The
Functionto broadcast.root – rank to broadcast from. Defaults to 0.
- Returns:
list of MPI.Request objects (one for each of f.subfunctions).
- Raises:
ValueError – if function communicator mismatches the ensemble spatial communicator.
- irecv(f, source=-1, tag=-1)[source]¶
Receive (non-blocking) a function f over
ensemble_commfrom another ensemble rank.- Parameters:
f – The a
Functionto receive intosource – the rank to receive from. Defaults to MPI.ANY_SOURCE.
tag – the tag of the message. Defaults to MPI.ANY_TAG.
- Returns:
list of MPI.Request objects (one for each of f.subfunctions).
- Raises:
ValueError – if function communicator mismatches the ensemble spatial communicator.
- ireduce(f, f_reduced, op=<mpi4py.MPI.Op object>, root=0)[source]¶
Reduce (non-blocking) a function f into f_reduced over
ensemble_commto rank root- Parameters:
f – The a
Functionto reduce.f_reduced – the result of the reduction on rank root.
op – MPI reduction operator. Defaults to MPI.SUM.
root – rank to reduce to. Defaults to 0.
- Returns:
list of MPI.Request objects (one for each of f.subfunctions).
- Raises:
ValueError – if function communicators mismatch each other or the ensemble spatial communicator, or is the functions are in different spaces
- isend(f, dest, tag=0)[source]¶
Send (non-blocking) a function f over
ensemble_commto another ensemble rank.- Parameters:
f – The a
Functionto senddest – the rank to send to
tag – the tag of the message. Defaults to 0.
- Returns:
list of MPI.Request objects (one for each of f.subfunctions).
- Raises:
ValueError – if function communicator mismatches the ensemble spatial communicator.
- isendrecv(fsend, dest, sendtag=0, frecv=None, source=-1, recvtag=-1)[source]¶
Send a function fsend and receive a function frecv over
ensemble_commto another ensemble rank.- Parameters:
- Returns:
list of MPI.Request objects (one for each of fsend.subfunctions and frecv.subfunctions).
- Raises:
ValueError – if function communicator mismatches the ensemble spatial communicator.
- recv(f, source=-1, tag=-1, statuses=None)[source]¶
Receive (blocking) a function f over
ensemble_commfrom another ensemble rank.- Parameters:
f – The a
Functionto receive intosource – the rank to receive from. Defaults to MPI.ANY_SOURCE.
tag – the tag of the message. Defaults to MPI.ANY_TAG.
statuses – MPI.Status objects (one for each of f.subfunctions or None).
- Raises:
ValueError – if function communicator mismatches the ensemble spatial communicator.
- reduce(f, f_reduced, op=<mpi4py.MPI.Op object>, root=0)[source]¶
Reduce a function f into f_reduced over
ensemble_commto rank root- Parameters:
f – The a
Functionto reduce.f_reduced – the result of the reduction on rank root.
op – MPI reduction operator. Defaults to MPI.SUM.
root – rank to reduce to. Defaults to 0.
- Raises:
ValueError – if function communicators mismatch each other or the ensemble spatial communicator, or is the functions are in different spaces
- send(f, dest, tag=0)[source]¶
Send (blocking) a function f over
ensemble_commto another ensemble rank.- Parameters:
f – The a
Functionto senddest – the rank to send to
tag – the tag of the message. Defaults to 0
- Raises:
ValueError – if function communicator mismatches the ensemble spatial communicator.
- sendrecv(fsend, dest, sendtag=0, frecv=None, source=-1, recvtag=-1, status=None)[source]¶
Send (blocking) a function fsend and receive a function frecv over
ensemble_commto another ensemble rank.- Parameters:
fsend – The a
Functionto send.dest – the rank to send to.
sendtag – the tag of the send message. Defaults to 0.
frecv – The a
Functionto receive into.source – the rank to receive from. Defaults to MPI.ANY_SOURCE.
recvtag – the tag of the received message. Defaults to MPI.ANY_TAG.
status – MPI.Status object or None.
- Raises:
ValueError – if function communicator mismatches the ensemble spatial communicator.
firedrake.ensemble.ensemble_function module¶
- class firedrake.ensemble.ensemble_function.EnsembleCofunction(function_space: EnsembleDualSpace)[source]¶
Bases:
EnsembleFunctionBaseA mixed finite element Cofunction distributed over an ensemble.
- Parameters:
function_space (
EnsembleDualSpace) – The function space of the cofunction.
- class firedrake.ensemble.ensemble_function.EnsembleFunction(function_space: EnsembleFunctionSpaceBase)[source]¶
Bases:
EnsembleFunctionBaseA mixed Function defined on a
Ensemble. The subcomponents are distributed over the ensemble members, and are specified locally in anEnsembleFunctionSpace.- Parameters:
function_space (
EnsembleFunctionSpace.) – The function space of the Function.
Notes
Passing an
EnsembleDualSpacetoEnsembleFunctionwill return an instance ofEnsembleCofunction.This class does not carry UFL symbolic information, unlike a
Function. UFL expressions can only be defined locally on each ensemble member using aFunctionfromEnsembleFunction.subfunctions.
- class firedrake.ensemble.ensemble_function.EnsembleFunctionBase(function_space: EnsembleFunctionSpaceBase)[source]¶
Bases:
EnsembleFunctionMixinA mixed (co)function defined on a
Ensemble. The subcomponents are distributed over the ensemble members, and are specified locally in anEnsembleFunctionSpace.- Parameters:
function_space – The function space of the (co)function.
Notes
Passing an
EnsembleDualSpacetoEnsembleFunctionwill return an instance ofEnsembleCofunction.This class does not carry UFL symbolic information, unlike a
Function. UFL expressions can only be defined locally on each ensemble member using a ~firedrake.function.Function fromEnsembleFunction.subfunctions.See also
ensemble_functionspace.EnsembleFunctionSpace,ensemble_functionspace.EnsembleDualSpace,EnsembleFunction,EnsembleCofunction- assign(other, subsets=None)[source]¶
Set the
EnsembleFunctionto the value of anotherEnsembleFunctionother.- Parameters:
other (
EnsembleFunction) – The value to assign from.subsets (Collection[Optional[
pyop2.types.set.Subset]]) – One subset for each localfiredrake.function.Function. None elements will be ignored. The values of each local function will only be assigned on the nodes on the corresponding subset.
- copy()[source]¶
Return a deep copy of the
EnsembleFunction.
- riesz_representation(**kwargs)[source]¶
Return the Riesz representation of this
EnsembleFunctionwith respect to the given Riesz map.Internally delegates to the
firedrake.function.Function.riesz_representation()of each component.- Parameters:
riesz_map – The Riesz map to use (l2, L2, or H1). This can also be a callable.
kwargs – other arguments to be passed to the firedrake.riesz_map.
- property subfunctions¶
The (co)functions on the local ensemble member.
- vec()[source]¶
Context manager for the global
PETSc.Vecwith read/write access.It is invalid to access the
Vecoutside of a context manager.
- vec_ro()[source]¶
Context manager for the global
PETSc.Vecwith read only access.It is invalid to access the
Vecoutside of a context manager.
- vec_wo()[source]¶
Context manager for the global
PETSc.Vecwith write only access.It is invalid to access the
Vecoutside of a context manager.
- zero(subsets=None)[source]¶
Set values to zero.
- Parameters:
subsets (Collection[Optional[
pyop2.types.set.Subset]]) – One subset for each localfiredrake.function.Function. None elements will be ignored. The values of each local function will only be zeroed on the nodes on the corresponding subset.
firedrake.ensemble.ensemble_functionspace module¶
- class firedrake.ensemble.ensemble_functionspace.EnsembleDualSpace(local_spaces: Collection, ensemble: Ensemble)[source]¶
Bases:
EnsembleFunctionSpaceBaseA mixed dual function space defined on an
ensemble.Ensemble. The subcomponents are distributed over the ensemble members, but are specified locally on each ensemble member.- Parameters:
local_spaces (Collection) – The list of dual function spaces on the local ensemble.comm.
ensemble (.ensemble.Ensemble) – The communicator that the function space is defined over.
Notes
Passing a list of dual local_spaces to
EnsembleFunctionSpacewill return an instance ofEnsembleDualSpace.This class does not carry UFL symbolic information, unlike a
FiredrakeDualSpace. UFL expressions can only be defined locally on each ensemble member using aFiredrakeDualSpacefrom EnsembleDualSpace.local_spaces.Examples
If U=CG1, V=DG0, and W=U*V, we can have the nested mixed dual space U*xV*xV*xW*xU*. This can be distributed over an
ensemble.Ensemblewith two ensemble members by splitting into [U*xV*]x[V*xW*xU*]. The following code creates the correspondingEnsembleDualSpace:ensemble = Ensemble(COMM_WORLD, COMM_WORLD.size//2) mesh = UnitIntervalMesh(8, comm=ensemble.comm) U = FunctionSpace(mesh, "CG", 1) V = FunctionSpace(mesh, "DG", 0) W = U*V if ensemble.ensemble_rank == 0: local_spaces = [U.dual(), V.dual()] else: local_spaces = [V.dual(), W.dual(), U.dual()] efs0 = EnsembleDualSpace(local_spaces, ensemble) efs1 = EnsembleFunctionSpace(local_spaces, ensemble)
- class firedrake.ensemble.ensemble_functionspace.EnsembleFunctionSpace(local_spaces: Collection, ensemble: Ensemble)[source]¶
Bases:
EnsembleFunctionSpaceBaseA mixed primal function space defined on an
Ensemble. The subcomponents are distributed over the ensemble members, but are specified locally on each ensemble member.- Parameters:
local_spaces (Collection) – The list of primal function spaces on the local
Ensemble.comm.ensemble (
Ensemble) – The communicator that the function space is defined over.
Notes
Passing a list of dual local_spaces to
EnsembleFunctionSpacewill return an instance ofEnsembleDualSpace.This class does not carry UFL symbolic information, unlike a
FunctionSpace. UFL expressions can only be defined locally on each ensemble member using aFunctionSpacefrom EnsembleFunctionSpace.local_spaces.Examples
If U=CG1, V=DG0, and W=UxV, we can have the nested mixed space UxVxVxWxU. This can be distributed over an
ensemble.Ensemblewith two ensemble members by splitting into [UxV]x[VxWxU]. The following code creates the correspondingEnsembleFunctionSpace:ensemble = Ensemble(COMM_WORLD, COMM_WORLD.size//2) mesh = UnitIntervalMesh(8, comm=ensemble.comm) U = FunctionSpace(mesh, "CG", 1) V = FunctionSpace(mesh, "DG", 0) W = U*V if ensemble.ensemble_rank == 0: local_spaces = [U, V] else: local_spaces = [V, W, U] efs = EnsembleFunctionSpace(local_spaces, ensemble)
- class firedrake.ensemble.ensemble_functionspace.EnsembleFunctionSpaceBase(local_spaces: Collection, ensemble: Ensemble)[source]¶
Bases:
objectBase class for mixed function spaces defined on an
Ensemble. The subcomponents are distributed over the ensemble members, and are specified locally.- Parameters:
local_spaces (Collection) – The list of function spaces on the local ensemble.comm.
ensemble (~.ensemble.Ensemble) – The communicator that the function space is defined over.
Notes
Passing a list of dual local_spaces to
EnsembleFunctionSpacewill return an instance ofEnsembleDualSpace.This class does not carry UFL symbolic information, unlike a
FunctionSpace. UFL expressions can only be defined locally on each ensemble member using aFunctionSpacefrom EnsembleFunctionSpace.local_spaces.See also
EnsembleFunctionSpace,EnsembleDualSpace,ensemble_function.EnsembleFunction,ensemble_function.EnsembleCofunction- create_vec()[source]¶
Return a PETSc Vec on the
Ensemble.global_commwith the same layout as anEnsembleFunctionorEnsembleCofunctionin this function space.
- dual()[source]¶
The dual to this function space. A
EnsembleDualSpaceif self is aEnsembleFunctionSpace, and vice-versa.
- property global_spaces_offset¶
Index of the first local subspace in the global mixed space.
- property layout_vec¶
- property local_spaces¶
The
FunctionSpaceon the local ensemble.comm.
- property nglobal_dofs¶
The total number of dofs across all subspaces on all ensemble ranks.
- property nglobal_spaces¶
The total number of subspaces across all ensemble ranks.
- property nlocal_comm_dofs¶
The total number of dofs across all subspaces on the local ensemble.comm.
- property nlocal_rank_dofs¶
The total number of dofs across all subspaces on the local MPI rank.
- property nlocal_spaces¶
The number of subspaces on this ensemble rank.
firedrake.ensemble.ensemble_mat module¶
- firedrake.ensemble.ensemble_mat.EnsembleBlockDiagonalMat(block_mats: Iterable, row_space: EnsembleFunctionSpaceBase, col_space: EnsembleFunctionSpaceBase)[source]¶
A Mat for a block diagonal matrix defined over an
Ensemble. Each block acts on a single subspace of anEnsembleFunctionSpace. This is a convenience function to create a PETSc.Mat with aEnsembleBlockDiagonalMatCtxPython context.- Parameters:
block_mats (Iterable[PETSc.Mat]) – The PETSc Mats for each block. On each ensemble rank there must be as many Mats as there are local subspaces of
row_spaceandcol_space, and the Mat sizes must match the sizes of the corresponding subspaces.row_space – The function space that the matrix acts on. Must have the same number of subspaces on each ensemble rank as col_space.
col_space – The function space for the result of the matrix action. Must have the same number of subspaces on each ensemble rank as row_space.
- Returns:
The PETSc.Mat with an
EnsembleBlockDiagonalMatCtxPython context.- Return type:
PETSc.Mat
See also
- class firedrake.ensemble.ensemble_mat.EnsembleBlockDiagonalMatCtx(block_mats: Iterable, row_space: EnsembleFunctionSpaceBase, col_space: EnsembleFunctionSpaceBase)[source]¶
Bases:
EnsembleMatCtxBaseA python Mat context for a block diagonal matrix defined over an
Ensemble. Each block acts on a single subspace of anEnsembleFunctionSpace.- Parameters:
block_mats (Iterable[PETSc.Mat]) – The PETSc Mats for each block. On each ensemble rank there must be as many Mats as there are local subspaces of
row_spaceandcol_space, and the Mat sizes must match the sizes of the corresponding subspaces.row_space – The function space that the matrix acts on. Must have the same number of subspaces on each ensemble rank as col_space.
col_space – The function space for the result of the matrix action. Must have the same number of subspaces on each ensemble rank as row_space.
Notes
This is a python context, not an actual PETSc.Mat. To create the corresponding PETSc.Mat users should call
EnsembleBlockDiagonalMat().See also
- mult_impl(A, x, y)[source]¶
Apply the action of the matrix to x, putting the result in y.
y is not guaranteed to be zero on entry. This is a convenience method allowing the matrix action to be implemented in terms of EnsembleFunction input and outputs by inheriting classes.
- Parameters:
A (PETSc.Mat) – The PETSc matrix that self is the python context of.
x – The vector acted on by the matrix.
y – The result of the matrix action.
See also
- class firedrake.ensemble.ensemble_mat.EnsembleMatCtxBase(row_space: EnsembleFunctionSpaceBase, col_space: EnsembleFunctionSpaceBase)[source]¶
Bases:
objectBase class for python type Mats defined over an
Ensemble.- Parameters:
row_space – The function space that the matrix acts on. Must have the same number of subspaces on each ensemble rank as col_space.
col_space – The function space for the result of the matrix action. Must have the same number of subspaces on each ensemble rank as row_space.
Notes
The main use of this base class is to enable users to implement the matrix action as acting on and resulting in an
EnsembleFunction. This is done by implementing themult_implmethod.See also
- mult(A, x, y)[source]¶
Apply the action of the matrix to x, putting the result in y.
This method will be called by PETSc with x and y as Vecs, and acts as a wrapper around the
mult_implmethod which has x and y as EnsembleFunction for convenience. y is not guaranteed to be zero on entry.- Parameters:
A (PETSc.Mat) – The PETSc matrix that self is the python context of.
x (PETSc.Vec) – The vector acted on by the matrix.
y (PETSc.Vec) – The result of the matrix action.
See also
- mult_impl(A, x: EnsembleFunctionBase, y: EnsembleFunctionBase)[source]¶
Apply the action of the matrix to x, putting the result in y.
y is not guaranteed to be zero on entry. This is a convenience method allowing the matrix action to be implemented in terms of EnsembleFunction input and outputs by inheriting classes.
- Parameters:
A (PETSc.Mat) – The PETSc matrix that self is the python context of.
x – The vector acted on by the matrix.
y – The result of the matrix action.
See also
firedrake.ensemble.ensemble_pc module¶
- class firedrake.ensemble.ensemble_pc.EnsembleBJacobiPC[source]¶
Bases:
EnsemblePCBaseA python PC context for a block Jacobi method defined over an
Ensemble. Each block acts on a single subspace of anEnsembleFunctionSpaceand is (approximately) solved with its own KSP, which defaults to -ksp_type preonly.Available options:
-pc_use_amat- use Amat to apply block of operator in inner Krylov method-sub_%d- set options for the%d’th block, numbered from ensemble rank 0.-sub_- set default options for all blocks.
Notes
Currently this is only implemented for
EnsembleBlockDiagonalMatCtxmatrices.- initialize(pc)[source]¶
Initialize any state in this preconditioner.
This method is only called on the first time that the
setUpmethod is called.
- prefix = 'ebjacobi_'¶
The options prefix of this PC.
- class firedrake.ensemble.ensemble_pc.EnsemblePCBase[source]¶
Bases:
PCBaseBase class for python type PCs defined over an
Ensemble.The pc operators must be python Mats with
EnsembleMatCtxBase.Notes
The main use of this base class is to enable users to implement the preconditioner action as acting on and resulting in an
EnsembleFunction. This is done by implementing theapply_implmethod.See also
- apply(pc, x, y)[source]¶
Apply the preconditioner to x, putting the result in y.
Both x and y are PETSc Vecs, y is not guaranteed to be zero on entry.
- initialize(pc)[source]¶
Initialize any state in this preconditioner.
This method is only called on the first time that the
setUpmethod is called.
- needs_python_pmat = True¶
Set this to False if the P matrix needs to be Python (matfree).
