2
2
import pytest
3
3
4
4
5
- def run_gtmg_mixed_poisson ():
5
+ def run_gtmg_mixed_poisson (custom_transfer = False ):
6
6
7
7
m = UnitSquareMesh (10 , 10 )
8
8
nlevels = 2
@@ -14,7 +14,7 @@ def get_p1_space():
14
14
return FunctionSpace (mesh , "CG" , 1 )
15
15
16
16
def get_p1_prb_bcs ():
17
- return DirichletBC (get_p1_space (), Constant ( 0.0 ) , "on_boundary" )
17
+ return DirichletBC (get_p1_space (), 0 , "on_boundary" )
18
18
19
19
def p1_callback ():
20
20
P1 = get_p1_space ()
@@ -57,13 +57,28 @@ def p1_callback():
57
57
appctx = {'get_coarse_operator' : p1_callback ,
58
58
'get_coarse_space' : get_p1_space ,
59
59
'coarse_space_bcs' : get_p1_prb_bcs ()}
60
-
61
- solve (a == L , w , solver_parameters = params , appctx = appctx )
60
+ if custom_transfer :
61
+ P1 = get_p1_space ()
62
+ V = FunctionSpace (mesh , "DGT" , degree - 1 )
63
+ I = assemble (Interpolate (TrialFunction (P1 ), V )).petscmat
64
+ R = PETSc .Mat ().createTranspose (I )
65
+ appctx ['interpolation_matrix' ] = I
66
+ appctx ['restriction_matrix' ] = R
67
+
68
+ problem = LinearVariationalProblem (a , L , w )
69
+ solver = LinearVariationalSolver (problem , solver_parameters = params , appctx = appctx )
70
+ solver .solve ()
62
71
_ , uh = w .subfunctions
63
72
64
73
# Analytical solution
65
74
f .interpolate (x [0 ]* (1 - x [0 ])* x [1 ]* (1 - x [1 ]))
66
75
76
+ if custom_transfer :
77
+ hyb = solver .snes .ksp .pc .getPythonContext ()
78
+ gtmg = hyb .trace_ksp .pc .getPythonContext ()
79
+ assert I .handle != R .handle
80
+ assert gtmg .pc .getMGInterpolation (1 ).handle == I .handle
81
+ assert gtmg .pc .getMGRestriction (1 ).handle == R .handle
67
82
return errornorm (f , uh , norm_type = "L2" )
68
83
69
84
@@ -144,8 +159,9 @@ def p1_callback():
144
159
145
160
146
161
@pytest .mark .skipcomplexnoslate
147
- def test_mixed_poisson_gtmg ():
148
- assert run_gtmg_mixed_poisson () < 1e-5
162
+ @pytest .mark .parametrize ("custom_transfer" , [False , True ])
163
+ def test_mixed_poisson_gtmg (custom_transfer ):
164
+ assert run_gtmg_mixed_poisson (custom_transfer ) < 1e-5
149
165
150
166
151
167
@pytest .mark .skipcomplexnoslate
0 commit comments