Skip to content

Fix FacetSplit.update #4238

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

Draft
wants to merge 6 commits into
base: master
Choose a base branch
from
Draft
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
34 changes: 18 additions & 16 deletions firedrake/preconditioners/facet_split.py
Original file line number Diff line number Diff line change
Expand Up @@ -75,27 +75,30 @@ def initialize(self, pc):
self.needs_zeroing = len(indices) < V.dof_count
self.subset = PETSc.IS().createGeneral(indices, comm=PETSc.COMM_SELF)

self.reset_operators = False
if mat_type != "submatrix":
from firedrake.assemble import get_assembler
form_assembler = get_assembler(mixed_operator, bcs=mixed_bcs,
form_compiler_parameters=fcp,
mat_type=mat_type,
options_prefix=options_prefix)
self.P = form_assembler.allocate()
self._assemble_mixed_op = form_assembler.assemble
self._assemble_mixed_op(tensor=self.P)
self.mixed_opmat = self.P.petscmat
self.set_nullspaces(pc)
self.work_vecs = self.mixed_opmat.createVecs()
self._permute_op = partial(form_assembler.assemble, tensor=self.P)
elif self.subset:
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Configuration member requested below would probably be set before this if-else block and then used for this block, and in the update method below.

self.mixed_opmat = PETSc.Mat()
global_indices = V.dof_dset.lgmap.apply(self.subset.indices)
self._global_iperm = PETSc.IS().createGeneral(global_indices, comm=pc.comm)
self._permute_op = partial(PETSc.Mat().createSubMatrixVirtual, P, self._global_iperm, self._global_iperm)
self.mixed_opmat = self._permute_op()
self.set_nullspaces(pc)
self.work_vecs = self.mixed_opmat.createVecs()
self._permute_op = partial(self.mixed_opmat.createSubMatrixVirtual, P, self._global_iperm, self._global_iperm)
self.reset_operators = True
else:
self.mixed_opmat = P
self._permute_op = lambda *args: None

aux = self._permute_op()
if aux is not None:
self.set_nullspaces(pc)
self.work_vecs = self.mixed_opmat.createVecs()

# Internally, we just set up a PC object that the user can configure
# however from the PETSc command line. Since PC allows the user to specify
Expand Down Expand Up @@ -142,14 +145,13 @@ def _restrict_nullspace(nsp):
Pmat.setTransposeNullSpace(_restrict_nullspace(P.getTransposeNullSpace()))

def update(self, pc):
if hasattr(self, "_permute_op"):
for mat in self.pc.getOperators():
mat.destroy()
P = self._permute_op()
self.pc.setOperators(A=P, P=P)
self.mixed_opmat = P
elif hasattr(self, "P"):
self._assemble_mixed_op(tensor=self.P)
self._permute_op()
if self.reset_operators:
self.pc.setOperators(A=self.mixed_opmat, P=self.mixed_opmat)

dm = self.pc.getDM()
with dmhooks.add_hooks(dm, self, appctx=self._ctx_ref):
self.pc.setFromOptions()

def prolong(self, x, y):
if x is not y:
Expand Down
45 changes: 40 additions & 5 deletions tests/firedrake/regression/test_facet_split.py
Original file line number Diff line number Diff line change
Expand Up @@ -2,7 +2,9 @@
from firedrake import *


def run_facet_split(quadrilateral, pc_type, refine=2):
def run_facet_split(quadrilateral, pc_type, refine=2, betas=None):
if betas is None:
betas = [0]
if pc_type == "lu":
parameters = {
"mat_type": "matfree",
Expand All @@ -17,6 +19,7 @@ def run_facet_split(quadrilateral, pc_type, refine=2):
"pc_fieldsplit_schur_precondition": "selfp",
"fieldsplit_0_pc_type": "jacobi",
"fieldsplit_1_pc_type": "lu",
"fieldsplit_ksp_type": "preonly",
},
}
elif pc_type == "jacobi":
Expand All @@ -37,6 +40,26 @@ def run_facet_split(quadrilateral, pc_type, refine=2):
"fieldsplit_1_ksp_rtol": 1E-12,
},
}
elif pc_type == "fdm":
parameters = {
"mat_type": "matfree",
"ksp_type": "preonly",
"pc_type": "python",
"pc_python_type": "firedrake.FacetSplitPC",
"facet": {
"mat_type": "submatrix",
"pc_type": "fieldsplit",
"pc_fieldsplit_type": "schur",
"pc_fieldsplit_schur_fact_type": "full",
"pc_fieldsplit_schur_precondition": "a11",
"fieldsplit_0_pc_type": "jacobi",
"fieldsplit_1_pc_type": "python",
"fieldsplit_1_pc_python_type": "firedrake.FDMPC",
"fieldsplit_1_fdm_static_condensation": True,
"fieldsplit_1_fdm_pc_type": "lu",
"fieldsplit_ksp_type": "preonly",
},
}

r = refine
variant = "fdm" if quadrilateral else None
Expand All @@ -46,15 +69,21 @@ def run_facet_split(quadrilateral, pc_type, refine=2):
u = TrialFunction(V)
v = TestFunction(V)
uh = Function(V)

a = inner(grad(u), grad(v)) * dx
L = inner(Constant(0), v) * dx
x = SpatialCoordinate(mesh)
u_exact = 42 * x[1]
bcs = [DirichletBC(V, Constant(0), 3),
DirichletBC(V, Constant(42), 4)]

solve(a == L, uh, bcs=bcs, solver_parameters=parameters)
beta = Constant(0)
a = inner(grad(u), grad(v)) * dx + inner(u*beta, v) * dx
L = inner(u_exact * beta, v) * dx
problem = LinearVariationalProblem(a, L, uh, bcs=bcs)
solver = LinearVariationalSolver(problem, solver_parameters=parameters)
for val in betas:
beta.assign(val)
uh.assign(0)
solver.solve()

return sqrt(assemble(inner(uh - u_exact, uh - u_exact) * dx))


Expand All @@ -68,3 +97,9 @@ def test_facet_split(quadrilateral, pc_type):
@pytest.mark.parametrize("pc_type", ["lu", "jacobi"])
def test_facet_split_parallel(pc_type):
assert run_facet_split(True, pc_type, refine=3) < 1E-10


@pytest.mark.parametrize("quadrilateral", [True])
@pytest.mark.parametrize("pc_type", ["fdm"])
def test_facet_split_update(quadrilateral, pc_type):
assert run_facet_split(quadrilateral, pc_type, betas=[1E4, 0]) < 1E-10