xref: /petsc/src/binding/petsc4py/test/test_mat_fact.py (revision 62e5d2d2208a68fdbd23dac44110783534664de8)
15808f684SSatish Balayfrom petsc4py import PETSc
25808f684SSatish Balayimport unittest
35808f684SSatish Balay
45808f684SSatish Balayimport numpy as N
55808f684SSatish Balay
65808f684SSatish Balaydef mkmat(n, mtype, opts):
75808f684SSatish Balay    A = PETSc.Mat().create(PETSc.COMM_SELF)
85808f684SSatish Balay    A.setSizes([n,n])
95808f684SSatish Balay    A.setType(mtype)
105808f684SSatish Balay    A.setUp()
115808f684SSatish Balay    for o in opts:
125808f684SSatish Balay        A.setOption(o, True)
135808f684SSatish Balay    return A
145808f684SSatish Balay
155808f684SSatish Balaydef mksys_diag(n, mtype, opts):
165808f684SSatish Balay    A = mkmat(n, mtype, opts)
175808f684SSatish Balay    x, b = A.createVecs()
185808f684SSatish Balay    for i in range(n):
195808f684SSatish Balay        A[i,i] = i+1
205808f684SSatish Balay        x[i]   = 1.0/(i+1)
215808f684SSatish Balay        b[i]   = 1
225808f684SSatish Balay    A.assemble()
235808f684SSatish Balay    x.assemble()
245808f684SSatish Balay    b.assemble()
255808f684SSatish Balay    return A, x, b
265808f684SSatish Balay
275808f684SSatish Balaydef mksys_poi2(n, mtype, opts):
285808f684SSatish Balay    A = mkmat(n, mtype, opts)
295808f684SSatish Balay    x, b = A.createVecs()
305808f684SSatish Balay    for i in range(n):
315808f684SSatish Balay        if i == 0:
325808f684SSatish Balay            cols = [i, i+1]
335808f684SSatish Balay            vals = [2, -1]
345808f684SSatish Balay        elif i == n-1:
355808f684SSatish Balay            cols = [i-1, i]
365808f684SSatish Balay            vals = [-1,  2]
375808f684SSatish Balay        else:
385808f684SSatish Balay            cols = [i-1, i, i+1]
395808f684SSatish Balay            vals = [-1,  2, -1]
405808f684SSatish Balay        A[i,cols] = vals
415808f684SSatish Balay        x[i]   = i+1
425808f684SSatish Balay        b[i]   = 0
435808f684SSatish Balay    A.assemble()
445808f684SSatish Balay    x.assemble()
455808f684SSatish Balay    b.assemble()
465808f684SSatish Balay    A.mult(x,b)
475808f684SSatish Balay    return A, x, b
485808f684SSatish Balay
495808f684SSatish Balayclass BaseTestMatFactor(object):
505808f684SSatish Balay
515808f684SSatish Balay    MKSYS = None
525808f684SSatish Balay    MTYPE = None
535808f684SSatish Balay    MOPTS = ()
545808f684SSatish Balay
555808f684SSatish Balay    def setUp(self):
565808f684SSatish Balay        A, x, b = self.MKSYS(10, self.MTYPE, self.MOPTS)
575808f684SSatish Balay        self.A = A
585808f684SSatish Balay        self.x = x
595808f684SSatish Balay        self.b = b
605808f684SSatish Balay
615808f684SSatish Balay    def tearDown(self):
625808f684SSatish Balay        self.A.setUnfactored()
635808f684SSatish Balay        self.A.destroy(); self.A = None
645808f684SSatish Balay        self.x.destroy(); self.x = None
655808f684SSatish Balay        self.b.destroy(); self.b = None
66*62e5d2d2SJDBetteridge        PETSc.garbage_cleanup()
675808f684SSatish Balay
685808f684SSatish Balayclass BaseTestMatFactorLU(BaseTestMatFactor):
695808f684SSatish Balay
705808f684SSatish Balay    def testFactorLU(self):
715808f684SSatish Balay        r, c = self.A.getOrdering("nd")
725808f684SSatish Balay        self.A.reorderForNonzeroDiagonal(r, c)
735808f684SSatish Balay        self.A.factorLU(r,c,{'zeropivot':1e-5})
745808f684SSatish Balay        x = self.x.duplicate()
755808f684SSatish Balay        self.A.solve(self.b, x)
765808f684SSatish Balay        x.axpy(-1, self.x)
775808f684SSatish Balay        self.assertTrue(x.norm() < 1e-3)
785808f684SSatish Balay
795808f684SSatish Balayclass BaseTestMatFactorILU(BaseTestMatFactor):
805808f684SSatish Balay
815808f684SSatish Balay    def testFactorILU(self):
825808f684SSatish Balay        r, c = self.A.getOrdering("natural")
835808f684SSatish Balay        self.A.factorILU(r,c,{'levels':0})
845808f684SSatish Balay        x = self.x.duplicate()
855808f684SSatish Balay        self.A.solve(self.b, x)
865808f684SSatish Balay        x.axpy(-1, self.x)
875808f684SSatish Balay        self.assertTrue(x.norm() < 1e-3)
885808f684SSatish Balay
895808f684SSatish Balay## class BaseTestMatFactorILUDT(BaseTestMatFactor):
905808f684SSatish Balay##
915808f684SSatish Balay##     def testFactorILUDT(self):
925808f684SSatish Balay##         r, c = self.A.getOrdering("natural")
935808f684SSatish Balay##         self.A = self.A.factorILUDT(r,c)
945808f684SSatish Balay##         x = self.x.duplicate()
955808f684SSatish Balay##         self.A.solve(self.b, x)
965808f684SSatish Balay##         x.axpy(-1, self.x)
975808f684SSatish Balay##         self.assertTrue(x.norm() < 1e-3)
985808f684SSatish Balay##
995808f684SSatish Balayclass BaseTestMatFactorChol(BaseTestMatFactor):
1005808f684SSatish Balay
1015808f684SSatish Balay    def testFactorChol(self):
1025808f684SSatish Balay        r, c = self.A.getOrdering("natural")
1035808f684SSatish Balay        self.A.factorCholesky(r)
1045808f684SSatish Balay        x = self.x.duplicate()
1055808f684SSatish Balay        self.A.solve(self.b, x)
1065808f684SSatish Balay        x.axpy(-1, self.x)
1075808f684SSatish Balay        self.assertTrue(x.norm() < 1e-3)
1085808f684SSatish Balay
1095808f684SSatish Balayclass BaseTestMatFactorICC(BaseTestMatFactor):
1105808f684SSatish Balay
1115808f684SSatish Balay    def testFactorICC(self):
1125808f684SSatish Balay        r, c = self.A.getOrdering("natural")
1135808f684SSatish Balay        self.A.factorICC(r)
1145808f684SSatish Balay        x = self.x.duplicate()
1155808f684SSatish Balay        self.A.solve(self.b, x)
1165808f684SSatish Balay        x.axpy(-1, self.x)
1175808f684SSatish Balay        self.assertTrue(x.norm() < 1e-3)
1185808f684SSatish Balay
1195808f684SSatish Balay
1205808f684SSatish Balay# --------------------------------------------------------------------
1215808f684SSatish Balay
1225808f684SSatish Balayclass TestMatFactorA1(BaseTestMatFactorLU,
1235808f684SSatish Balay                      BaseTestMatFactorChol,
1245808f684SSatish Balay                      unittest.TestCase):
1255808f684SSatish Balay    MKSYS = staticmethod(mksys_diag)
1265808f684SSatish Balay    MTYPE = PETSc.Mat.Type.SEQDENSE
1275808f684SSatish Balay
1285808f684SSatish Balayclass TestMatFactorA2(BaseTestMatFactorLU,
1295808f684SSatish Balay                      BaseTestMatFactorChol,
1305808f684SSatish Balay                      unittest.TestCase):
1315808f684SSatish Balay    MKSYS = staticmethod(mksys_poi2)
1325808f684SSatish Balay    MTYPE = PETSc.Mat.Type.SEQDENSE
1335808f684SSatish Balay
1345808f684SSatish Balay# ---
1355808f684SSatish Balay
1365808f684SSatish Balayclass TestMatFactorB1(BaseTestMatFactorLU,
1375808f684SSatish Balay                      BaseTestMatFactorILU,
1385808f684SSatish Balay                      ## BaseTestMatFactorILUDT,
1395808f684SSatish Balay                      unittest.TestCase):
1405808f684SSatish Balay    MKSYS = staticmethod(mksys_diag)
1415808f684SSatish Balay    MTYPE = PETSc.Mat.Type.SEQAIJ
1425808f684SSatish Balay
1435808f684SSatish Balayclass TestMatFactorB2(BaseTestMatFactorLU,
1445808f684SSatish Balay                      BaseTestMatFactorILU,
1455808f684SSatish Balay                      ## BaseTestMatFactorILUDT,
1465808f684SSatish Balay                      unittest.TestCase):
1475808f684SSatish Balay    MKSYS = staticmethod(mksys_poi2)
1485808f684SSatish Balay    MTYPE = PETSc.Mat.Type.SEQAIJ
1495808f684SSatish Balay
1505808f684SSatish Balay# ---
1515808f684SSatish Balay
1525808f684SSatish Balayclass TestMatFactorC1(BaseTestMatFactorLU,
1535808f684SSatish Balay                      BaseTestMatFactorILU,
1545808f684SSatish Balay                      unittest.TestCase):
1555808f684SSatish Balay    MKSYS = staticmethod(mksys_diag)
1565808f684SSatish Balay    MTYPE = PETSc.Mat.Type.SEQBAIJ
1575808f684SSatish Balay
1585808f684SSatish Balayclass TestMatFactorC2(BaseTestMatFactorLU,
1595808f684SSatish Balay                      BaseTestMatFactorILU,
1605808f684SSatish Balay                      unittest.TestCase):
1615808f684SSatish Balay    MKSYS = staticmethod(mksys_poi2)
1625808f684SSatish Balay    MTYPE = PETSc.Mat.Type.SEQBAIJ
1635808f684SSatish Balay
1645808f684SSatish Balay# ---
1655808f684SSatish Balay
1665808f684SSatish Balayclass TestMatFactorD1(BaseTestMatFactorChol,
1675808f684SSatish Balay                      BaseTestMatFactorICC,
1685808f684SSatish Balay                      unittest.TestCase):
1695808f684SSatish Balay    MKSYS = staticmethod(mksys_diag)
1705808f684SSatish Balay    MTYPE = PETSc.Mat.Type.SEQSBAIJ
1715808f684SSatish Balay    MOPTS = [PETSc.Mat.Option.IGNORE_LOWER_TRIANGULAR]
1725808f684SSatish Balay
1735808f684SSatish Balayclass TestMatFactorD2(BaseTestMatFactorChol,
1745808f684SSatish Balay                      BaseTestMatFactorICC,
1755808f684SSatish Balay                      unittest.TestCase):
1765808f684SSatish Balay    MKSYS = staticmethod(mksys_poi2)
1775808f684SSatish Balay    MTYPE = PETSc.Mat.Type.SEQSBAIJ
1785808f684SSatish Balay    MOPTS = [PETSc.Mat.Option.IGNORE_LOWER_TRIANGULAR]
1795808f684SSatish Balay
1805808f684SSatish Balay# --------------------------------------------------------------------
1815808f684SSatish Balay
1825808f684SSatish Balayif __name__ == '__main__':
1835808f684SSatish Balay    unittest.main()
184