-
Notifications
You must be signed in to change notification settings - Fork 178
BDDC: cellwise subdomains, matfree MatIS, and user-defined primal vertices #4757
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Open
pbrubeck
wants to merge
40
commits into
main
Choose a base branch
from
pbrubeck/bddc-primal-vertices
base: main
Could not load branches
Branch not found: {{ refName }}
Loading
Could not load tags
Nothing to show
Loading
Are you sure you want to change the base?
Some commits from the old base branch may be removed from the timeline,
and old review comments may become outdated.
+214
−49
Open
Changes from all commits
Commits
Show all changes
40 commits
Select commit
Hold shift + click to select a range
9e9cb95
Submesh on subcommunicator
pbrubeck a1cca74
test
pbrubeck 586c545
PETSc version
pbrubeck 9e69d99
BDDC: create_matis
pbrubeck 5e017f0
Merge branch 'pbrubeck/submesh-comm-self' into pbrubeck/bddc-primal-v…
pbrubeck ed2f1ab
reorder kwarg
pbrubeck 60a6a8e
Merge branch 'pbrubeck/submesh-comm-self' into pbrubeck/bddc-primal-v…
pbrubeck 26e505c
Optional comm kwarg
pbrubeck af8b9e8
Merge branch 'pbrubeck/submesh-comm-self' into HEAD
pbrubeck 6fb20e0
FML: fix subtraction and support BaseForm
pbrubeck 4a91dee
Cleanup
pbrubeck 2df0024
Cleanup
pbrubeck d5df9d1
Cleanup
pbrubeck fb52209
fixup
pbrubeck cec8083
fixup
pbrubeck 5e1fdb3
Merge branch 'pbrubeck/submesh-comm-self' into pbrubeck/bddc-primal-v…
pbrubeck fdae167
fixup
pbrubeck c125fd3
fix
pbrubeck 5bd3312
Merge branch 'pbrubeck/submesh-comm-self' into pbrubeck/bddc-primal-v…
pbrubeck 2f2e77d
add assemble test
pbrubeck 9961010
add assemble test
pbrubeck 2cf9427
merge conflict
pbrubeck 825d2ef
add assemble test
pbrubeck 7124a31
Merge branch 'pbrubeck/submesh-comm-self' into pbrubeck/bddc-primal-v…
pbrubeck 6134381
use correct COMM_SELF
pbrubeck d78a670
merge conflict
pbrubeck 554d875
Apply suggestions from code review
pbrubeck 2588e75
merge conflict
pbrubeck 9bdf6ca
Merge branch 'pbrubeck/bddc-primal-vertices' of github.com:firedrakep…
pbrubeck 090ff86
do not reorder
pbrubeck 49a9d9e
Merge branch 'main' into pbrubeck/bddc-primal-vertices
pbrubeck 7a27812
create_matis: support reordered meshes
pbrubeck 4a5564f
cleanup
pbrubeck 91d0786
add a test
pbrubeck da5f10e
Merge branch 'main' into pbrubeck/bddc-primal-vertices
pbrubeck d88d80d
fixup
pbrubeck b261d17
Merge branch 'main' into pbrubeck/bddc-primal-vertices
pbrubeck 5c43f94
lint
pbrubeck 9cfe641
remove is_lagrange
pbrubeck e77fe20
Merge branch 'main' into pbrubeck/bddc-primal-vertices
pbrubeck File filter
Filter by extension
Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
Some comments aren't visible on the classic Files Changed page.
There are no files selected for viewing
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
| Original file line number | Diff line number | Diff line change | ||||||
|---|---|---|---|---|---|---|---|---|
| @@ -1,3 +1,5 @@ | ||||||||
| from itertools import repeat | ||||||||
|
|
||||||||
| from firedrake.preconditioners.base import PCBase | ||||||||
| from firedrake.preconditioners.patch import bcdofs | ||||||||
| from firedrake.preconditioners.facet_split import get_restriction_indices | ||||||||
|
|
@@ -6,11 +8,16 @@ | |||||||
| from firedrake.ufl_expr import TestFunction, TrialFunction | ||||||||
| from firedrake.function import Function | ||||||||
| from firedrake.functionspace import FunctionSpace, VectorFunctionSpace, TensorFunctionSpace | ||||||||
| from firedrake.preconditioners.fdm import tabulate_exterior_derivative | ||||||||
| from firedrake.preconditioners.fdm import broken_function, tabulate_exterior_derivative | ||||||||
| from firedrake.preconditioners.hiptmair import curl_to_grad | ||||||||
| from ufl import H1, H2, inner, dx, JacobianDeterminant | ||||||||
| from firedrake.parloops import par_loop, INC, READ | ||||||||
| from firedrake.utils import cached_property | ||||||||
| from firedrake.bcs import DirichletBC | ||||||||
| from firedrake.mesh import Submesh | ||||||||
| from ufl import Form, H1, H2, JacobianDeterminant, dx, inner, replace | ||||||||
| from finat.ufl import BrokenElement | ||||||||
| from pyop2.mpi import COMM_SELF | ||||||||
| from pyop2.utils import as_tuple | ||||||||
| import gem | ||||||||
| import numpy | ||||||||
|
|
||||||||
| __all__ = ("BDDCPC",) | ||||||||
|
|
@@ -23,6 +30,8 @@ class BDDCPC(PCBase): | |||||||
|
|
||||||||
| Internally, this PC creates a PETSc PCBDDC object that can be controlled by | ||||||||
| the options: | ||||||||
| - ``'bddc_cellwise'`` to set up a MatIS on cellwise subdomains if P.type == python, | ||||||||
| - ``'bddc_matfree'`` to set up a matrix-free MatIS if A.type == python, | ||||||||
| - ``'bddc_pc_bddc_neumann'`` to set sub-KSPs on subdomains excluding corners, | ||||||||
| - ``'bddc_pc_bddc_dirichlet'`` to set sub-KSPs on subdomain interiors, | ||||||||
| - ``'bddc_pc_bddc_coarse'`` to set the coarse solver KSP. | ||||||||
|
|
@@ -34,36 +43,58 @@ class BDDCPC(PCBase): | |||||||
| - ``'get_divergence_mat'`` for problems in H(div) (resp. 2D H(curl)), this is | ||||||||
| provide the arguments (a Mat with the assembled bilinear form testing the divergence | ||||||||
| (curl) against an L2 space) and keyword arguments supplied to ``PETSc.PC.setDivergenceMat``. | ||||||||
| - ``'primal_markers'`` a Function marking degrees of freedom to be included in the | ||||||||
| coarse space. If a DG(0) Function is provided, then all degrees of freedom on the cell | ||||||||
| are marked as primal. Alternatively, primal_markers can be a list with the global | ||||||||
| degrees of freedom to be included in the coarse space. | ||||||||
| """ | ||||||||
|
|
||||||||
| _prefix = "bddc_" | ||||||||
|
|
||||||||
| def initialize(self, pc): | ||||||||
| # Get context from pc | ||||||||
| _, P = pc.getOperators() | ||||||||
| dm = pc.getDM() | ||||||||
| self.prefix = (pc.getOptionsPrefix() or "") + self._prefix | ||||||||
| prefix = (pc.getOptionsPrefix() or "") + self._prefix | ||||||||
|
|
||||||||
| dm = pc.getDM() | ||||||||
| V = get_function_space(dm) | ||||||||
| variant = V.ufl_element().variant() | ||||||||
| sobolev_space = V.ufl_element().sobolev_space | ||||||||
|
|
||||||||
| # Create new PC object as BDDC type | ||||||||
| bddcpc = PETSc.PC().create(comm=pc.comm) | ||||||||
| bddcpc.incrementTabLevel(1, parent=pc) | ||||||||
| bddcpc.setOptionsPrefix(self.prefix) | ||||||||
| bddcpc.setOperators(*pc.getOperators()) | ||||||||
| bddcpc.setOptionsPrefix(prefix) | ||||||||
| bddcpc.setType(PETSc.PC.Type.BDDC) | ||||||||
|
|
||||||||
| opts = PETSc.Options(bddcpc.getOptionsPrefix()) | ||||||||
| matfree = opts.getBool("matfree", False) | ||||||||
|
|
||||||||
| # Set operators | ||||||||
| assemblers = [] | ||||||||
| A, P = pc.getOperators() | ||||||||
| if P.type == "python": | ||||||||
| # Reconstruct P as MatIS | ||||||||
| cellwise = opts.getBool("cellwise", False) | ||||||||
| P, assembleP = create_matis(P, "aij", cellwise=cellwise) | ||||||||
| assemblers.append(assembleP) | ||||||||
|
|
||||||||
| if P.type != "is": | ||||||||
| raise ValueError(f"Expecting P to be either 'matfree' or 'is', not {P.type}.") | ||||||||
|
Member
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Wouldn't this currently error even for
Suggested change
|
||||||||
|
|
||||||||
| if A.type == "python" and matfree: | ||||||||
| # Reconstruct A as MatIS | ||||||||
| A, assembleA = create_matis(A, "matfree", cellwise=P.getISAllowRepeated()) | ||||||||
| assemblers.append(assembleA) | ||||||||
| bddcpc.setOperators(A, P) | ||||||||
| self.assemblers = assemblers | ||||||||
|
|
||||||||
| # Do not use CSR of local matrix to define dofs connectivity unless requested | ||||||||
| # Using the CSR only makes sense for H1/H2 problems | ||||||||
| is_h1h2 = sobolev_space in [H1, H2] | ||||||||
| if "pc_bddc_use_local_mat_graph" not in opts and (not is_h1h2 or not is_lagrange(V.finat_element)): | ||||||||
| is_h1h2 = V.ufl_element().sobolev_space in {H1, H2} | ||||||||
| if "pc_bddc_use_local_mat_graph" not in opts and (not is_h1h2 or not V.finat_element.has_pointwise_dual_basis): | ||||||||
| opts["pc_bddc_use_local_mat_graph"] = False | ||||||||
|
|
||||||||
| # Handle boundary dofs | ||||||||
| # Get context from DM | ||||||||
| ctx = get_appctx(dm) | ||||||||
|
|
||||||||
| # Handle boundary dofs | ||||||||
| bcs = tuple(ctx._problem.dirichlet_bcs()) | ||||||||
| mesh = V.mesh().unique() | ||||||||
| if mesh.extruded and not mesh.extruded_periodic: | ||||||||
|
|
@@ -85,16 +116,17 @@ def initialize(self, pc): | |||||||
| bddcpc.setBDDCNeumannBoundaries(neu_bndr) | ||||||||
|
|
||||||||
| appctx = self.get_appctx(pc) | ||||||||
| degree = max(as_tuple(V.ufl_element().degree())) | ||||||||
|
|
||||||||
| # Set coordinates only if corner selection is requested | ||||||||
| # There's no API to query from PC | ||||||||
| if "pc_bddc_corner_selection" in opts: | ||||||||
| W = VectorFunctionSpace(V.mesh(), "Lagrange", degree, variant=variant) | ||||||||
| coords = Function(W).interpolate(V.mesh().coordinates) | ||||||||
| degree = max(as_tuple(V.ufl_element().degree())) | ||||||||
| variant = V.ufl_element().variant() | ||||||||
| W = VectorFunctionSpace(mesh, "Lagrange", degree, variant=variant) | ||||||||
| coords = Function(W).interpolate(mesh.coordinates) | ||||||||
| bddcpc.setCoordinates(coords.dat.data_ro.repeat(V.block_size, axis=0)) | ||||||||
|
|
||||||||
| tdim = V.mesh().topological_dimension | ||||||||
| tdim = mesh.topological_dimension | ||||||||
| if tdim >= 2 and V.finat_element.formdegree == tdim-1: | ||||||||
| allow_repeated = P.getISAllowRepeated() | ||||||||
| get_divergence = appctx.get("get_divergence_mat", get_divergence_mat) | ||||||||
|
|
@@ -116,14 +148,22 @@ def initialize(self, pc): | |||||||
| grad_kwargs = dict() | ||||||||
| bddcpc.setBDDCDiscreteGradient(*grad_args, **grad_kwargs) | ||||||||
|
|
||||||||
| # Set the user-defined primal (coarse) degrees of freedom | ||||||||
| primal_markers = appctx.get("primal_markers") | ||||||||
| if primal_markers is not None: | ||||||||
| primal_indices = get_primal_indices(V, primal_markers) | ||||||||
| primal_is = PETSc.IS().createGeneral(primal_indices.astype(PETSc.IntType), comm=pc.comm) | ||||||||
| bddcpc.setBDDCPrimalVerticesIS(primal_is) | ||||||||
|
|
||||||||
| bddcpc.setFromOptions() | ||||||||
| self.pc = bddcpc | ||||||||
|
|
||||||||
| def view(self, pc, viewer=None): | ||||||||
| self.pc.view(viewer=viewer) | ||||||||
|
|
||||||||
| def update(self, pc): | ||||||||
| pass | ||||||||
| for c in self.assemblers: | ||||||||
| c() | ||||||||
|
|
||||||||
| def apply(self, pc, x, y): | ||||||||
| self.pc.apply(x, y) | ||||||||
|
|
@@ -132,6 +172,104 @@ def applyTranspose(self, pc, x, y): | |||||||
| self.pc.applyTranspose(x, y) | ||||||||
|
|
||||||||
|
|
||||||||
| class BrokenDirichletBC(DirichletBC): | ||||||||
| def __init__(self, bc): | ||||||||
| self.bc = bc | ||||||||
| V = bc.function_space().broken_space() | ||||||||
| g = bc._original_arg | ||||||||
| super().__init__(V, g, bc.sub_domain) | ||||||||
|
|
||||||||
| @cached_property | ||||||||
| def nodes(self): | ||||||||
| u = Function(self.bc.function_space()) | ||||||||
| self.bc.set(u, 1) | ||||||||
| u = broken_function(u.function_space(), val=u.dat) | ||||||||
| return numpy.flatnonzero(u.dat.data) | ||||||||
|
|
||||||||
|
|
||||||||
| def create_matis(Amat, local_mat_type, cellwise=False): | ||||||||
| from firedrake.assemble import get_assembler | ||||||||
|
|
||||||||
| def local_mesh(mesh): | ||||||||
| key = "local_submesh" | ||||||||
| cache = mesh._shared_data_cache["local_submesh_cache"] | ||||||||
| try: | ||||||||
| return cache[key] | ||||||||
| except KeyError: | ||||||||
| if mesh.comm.size > 1: | ||||||||
| submesh = Submesh(mesh, mesh.topological_dimension, None, ignore_halo=True, reorder=False, comm=COMM_SELF) | ||||||||
| else: | ||||||||
| submesh = None | ||||||||
| return cache.setdefault(key, submesh) | ||||||||
|
|
||||||||
| def local_space(V, cellwise): | ||||||||
| mesh = local_mesh(V.mesh().unique()) | ||||||||
| element = BrokenElement(V.ufl_element()) if cellwise else None | ||||||||
| return V.reconstruct(mesh=mesh, element=element) | ||||||||
|
|
||||||||
| def local_argument(arg, cellwise): | ||||||||
| return arg.reconstruct(function_space=local_space(arg.function_space(), cellwise)) | ||||||||
|
|
||||||||
| def local_integral(it): | ||||||||
| extra_domain_integral_type_map = dict(it.extra_domain_integral_type_map()) | ||||||||
| extra_domain_integral_type_map[it.ufl_domain()] = it.integral_type() | ||||||||
| return it.reconstruct(domain=local_mesh(it.ufl_domain()), | ||||||||
| extra_domain_integral_type_map=extra_domain_integral_type_map) | ||||||||
|
|
||||||||
| def local_bc(bc, cellwise): | ||||||||
| V = bc.function_space() | ||||||||
| Vsub = local_space(V, False) | ||||||||
| sub_domain = list(bc.sub_domain) | ||||||||
| if "on_boundary" in sub_domain: | ||||||||
| sub_domain.remove("on_boundary") | ||||||||
| sub_domain.extend(V.mesh().unique().exterior_facets.unique_markers) | ||||||||
|
|
||||||||
| valid_markers = Vsub.mesh().unique().exterior_facets.unique_markers | ||||||||
| sub_domain = list(set(sub_domain) & set(valid_markers)) | ||||||||
| bc = bc.reconstruct(V=Vsub, g=0, sub_domain=sub_domain) | ||||||||
| if cellwise: | ||||||||
| bc = BrokenDirichletBC(bc) | ||||||||
| return bc | ||||||||
|
|
||||||||
| def local_to_global_map(V, cellwise): | ||||||||
| u = Function(V) | ||||||||
| u.dat.data_wo[:] = numpy.arange(*V.dof_dset.layout_vec.getOwnershipRange()) | ||||||||
|
|
||||||||
| Vsub = local_space(V, False) | ||||||||
| usub = Function(Vsub).assign(u) | ||||||||
| if cellwise: | ||||||||
| usub = broken_function(usub.function_space(), val=usub.dat) | ||||||||
| indices = usub.dat.data_ro.astype(PETSc.IntType) | ||||||||
| return PETSc.LGMap().create(indices, comm=V.comm) | ||||||||
|
|
||||||||
| assert Amat.type == "python" | ||||||||
| ctx = Amat.getPythonContext() | ||||||||
| form = ctx.a | ||||||||
| bcs = ctx.bcs | ||||||||
|
|
||||||||
| local_form = replace(form, {arg: local_argument(arg, cellwise) for arg in form.arguments()}) | ||||||||
| local_form = Form(list(map(local_integral, local_form.integrals()))) | ||||||||
| local_bcs = tuple(map(local_bc, bcs, repeat(cellwise))) | ||||||||
|
|
||||||||
| assembler = get_assembler(local_form, bcs=local_bcs, mat_type=local_mat_type) | ||||||||
| tensor = assembler.assemble() | ||||||||
|
|
||||||||
| rmap = local_to_global_map(form.arguments()[0].function_space(), cellwise) | ||||||||
| cmap = local_to_global_map(form.arguments()[1].function_space(), cellwise) | ||||||||
|
|
||||||||
| Amatis = PETSc.Mat().createIS(Amat.getSizes(), comm=Amat.getComm()) | ||||||||
| Amatis.setISAllowRepeated(cellwise) | ||||||||
| Amatis.setLGMap(rmap, cmap) | ||||||||
| Amatis.setISLocalMat(tensor.petscmat) | ||||||||
| Amatis.setUp() | ||||||||
| Amatis.assemble() | ||||||||
|
|
||||||||
| def update(): | ||||||||
| assembler.assemble(tensor=tensor) | ||||||||
| Amatis.assemble() | ||||||||
| return Amatis, update | ||||||||
|
|
||||||||
|
|
||||||||
| def get_restricted_dofs(V, domain): | ||||||||
| W = FunctionSpace(V.mesh(), V.ufl_element()[domain]) | ||||||||
| indices = get_restriction_indices(V, W) | ||||||||
|
|
@@ -166,7 +304,7 @@ def get_discrete_gradient(V): | |||||||
| nsp = VectorSpaceBasis([basis]) | ||||||||
| nsp.orthonormalize() | ||||||||
| gradient.setNullSpace(nsp.nullspace()) | ||||||||
| if not is_lagrange(Q.finat_element): | ||||||||
| if not Q.finat_element.has_pointwise_dual_basis: | ||||||||
| vdofs = get_restricted_dofs(Q, "vertex") | ||||||||
| gradient.compose('_elements_corners', vdofs) | ||||||||
|
|
||||||||
|
|
@@ -176,23 +314,25 @@ def get_discrete_gradient(V): | |||||||
| return grad_args, grad_kwargs | ||||||||
|
|
||||||||
|
|
||||||||
| def is_lagrange(finat_element): | ||||||||
| """Returns whether finat_element.dual_basis consists only of point evaluation dofs.""" | ||||||||
| try: | ||||||||
| Q, ps = finat_element.dual_basis | ||||||||
| except NotImplementedError: | ||||||||
| return False | ||||||||
| # Inspect the weight matrix | ||||||||
| # Lagrange elements have gem.Delta as the only terminal nodes | ||||||||
| children = [Q] | ||||||||
| while children: | ||||||||
| nodes = [] | ||||||||
| for c in children: | ||||||||
| if isinstance(c, gem.Delta): | ||||||||
| pass | ||||||||
| elif isinstance(c, gem.gem.Terminal): | ||||||||
| return False | ||||||||
| else: | ||||||||
| nodes.extend(c.children) | ||||||||
| children = nodes | ||||||||
| return True | ||||||||
| def get_primal_indices(V, primal_markers): | ||||||||
| if isinstance(primal_markers, Function): | ||||||||
| marker_space = primal_markers.function_space() | ||||||||
| if marker_space == V: | ||||||||
| markers = primal_markers | ||||||||
| elif marker_space.finat_element.space_dimension() == 1: | ||||||||
| shapes = (V.finat_element.space_dimension(), V.block_size) | ||||||||
| domain = "{[i,j]: 0 <= i < %d and 0 <= j < %d}" % shapes | ||||||||
| instructions = """ | ||||||||
| for i, j | ||||||||
| w[i,j] = w[i,j] + t[0] | ||||||||
| end | ||||||||
| """ | ||||||||
| markers = Function(V) | ||||||||
| par_loop((domain, instructions), dx, {"w": (markers, INC), "t": (primal_markers, READ)}) | ||||||||
| else: | ||||||||
| raise ValueError(f"Expecting markers in either {V.ufl_element()} or DG(0).") | ||||||||
| primal_indices = numpy.flatnonzero(markers.dat.data >= 1E-12) | ||||||||
| primal_indices += V.dof_dset.layout_vec.getOwnershipRange()[0] | ||||||||
| else: | ||||||||
| primal_indices = numpy.asarray(primal_markers) | ||||||||
| return primal_indices | ||||||||
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Thanks for adding these to the docstring. I just have a couple of questions.
Functioncan/should be is a bit unclear to me. What values should be in theDG(0)Function? Can I provide Functions in other spaces?