Skip to content

Commit 1115dc3

Browse files
committed
add a test
1 parent de74328 commit 1115dc3

File tree

2 files changed

+43
-8
lines changed

2 files changed

+43
-8
lines changed

firedrake/preconditioners/facet_split.py

+3-3
Original file line numberDiff line numberDiff line change
@@ -88,11 +88,11 @@ def initialize(self, pc):
8888
self.set_nullspaces(pc)
8989
self.work_vecs = self.mixed_opmat.createVecs()
9090
elif self.subset:
91-
self.P = PETSc.Mat()
91+
self.mixed_opmat = PETSc.Mat()
9292
global_indices = V.dof_dset.lgmap.apply(self.subset.indices)
9393
self._global_iperm = PETSc.IS().createGeneral(global_indices, comm=pc.comm)
94-
self._permute_op = partial(self.P.createSubMatrixVirtual, P, self._global_iperm, self._global_iperm)
95-
self.mixed_opmat = self._permute_op()
94+
self._permute_op = partial(self.mixed_opmat.createSubMatrixVirtual, P, self._global_iperm, self._global_iperm)
95+
self._permute_op()
9696
self.set_nullspaces(pc)
9797
self.work_vecs = self.mixed_opmat.createVecs()
9898
else:

tests/firedrake/regression/test_facet_split.py

+40-5
Original file line numberDiff line numberDiff line change
@@ -2,7 +2,9 @@
22
from firedrake import *
33

44

5-
def run_facet_split(quadrilateral, pc_type, refine=2):
5+
def run_facet_split(quadrilateral, pc_type, refine=2, betas=None):
6+
if betas is None:
7+
betas = [0]
68
if pc_type == "lu":
79
parameters = {
810
"mat_type": "matfree",
@@ -17,6 +19,7 @@ def run_facet_split(quadrilateral, pc_type, refine=2):
1719
"pc_fieldsplit_schur_precondition": "selfp",
1820
"fieldsplit_0_pc_type": "jacobi",
1921
"fieldsplit_1_pc_type": "lu",
22+
"fieldsplit_ksp_type": "preonly",
2023
},
2124
}
2225
elif pc_type == "jacobi":
@@ -37,6 +40,26 @@ def run_facet_split(quadrilateral, pc_type, refine=2):
3740
"fieldsplit_1_ksp_rtol": 1E-12,
3841
},
3942
}
43+
elif pc_type == "fdm":
44+
parameters = {
45+
"mat_type": "matfree",
46+
"ksp_type": "preonly",
47+
"pc_type": "python",
48+
"pc_python_type": "firedrake.FacetSplitPC",
49+
"facet": {
50+
"mat_type": "submatrix",
51+
"pc_type": "fieldsplit",
52+
"pc_fieldsplit_type": "schur",
53+
"pc_fieldsplit_schur_fact_type": "full",
54+
"pc_fieldsplit_schur_precondition": "a11",
55+
"fieldsplit_0_pc_type": "jacobi",
56+
"fieldsplit_1_pc_type": "python",
57+
"fieldsplit_1_pc_python_type": "firedrake.FDMPC",
58+
"fieldsplit_1_fdm_static_condensation": True,
59+
"fieldsplit_1_fdm_pc_type": "lu",
60+
"fieldsplit_ksp_type": "preonly",
61+
},
62+
}
4063

4164
r = refine
4265
variant = "fdm" if quadrilateral else None
@@ -46,15 +69,21 @@ def run_facet_split(quadrilateral, pc_type, refine=2):
4669
u = TrialFunction(V)
4770
v = TestFunction(V)
4871
uh = Function(V)
49-
50-
a = inner(grad(u), grad(v)) * dx
51-
L = inner(Constant(0), v) * dx
5272
x = SpatialCoordinate(mesh)
5373
u_exact = 42 * x[1]
5474
bcs = [DirichletBC(V, Constant(0), 3),
5575
DirichletBC(V, Constant(42), 4)]
5676

57-
solve(a == L, uh, bcs=bcs, solver_parameters=parameters)
77+
beta = Constant(0)
78+
a = inner(grad(u), grad(v)) * dx + inner(u*beta, v) * dx
79+
L = inner(u_exact * beta, v) * dx
80+
problem = LinearVariationalProblem(a, L, uh, bcs=bcs)
81+
solver = LinearVariationalSolver(problem, solver_parameters=parameters)
82+
for val in betas:
83+
beta.assign(val)
84+
uh.assign(0)
85+
solver.solve()
86+
5887
return sqrt(assemble(inner(uh - u_exact, uh - u_exact) * dx))
5988

6089

@@ -68,3 +97,9 @@ def test_facet_split(quadrilateral, pc_type):
6897
@pytest.mark.parametrize("pc_type", ["lu", "jacobi"])
6998
def test_facet_split_parallel(pc_type):
7099
assert run_facet_split(True, pc_type, refine=3) < 1E-10
100+
101+
102+
@pytest.mark.parametrize("quadrilateral", [True])
103+
@pytest.mark.parametrize("pc_type", ["fdm"])
104+
def test_facet_split_update(quadrilateral, pc_type):
105+
assert run_facet_split(quadrilateral, pc_type, refine=4, betas=[1E4, 0]) < 1E-10

0 commit comments

Comments
 (0)