1*5808f684SSatish Balayfrom petsc4py import PETSc 2*5808f684SSatish Balayimport unittest 3*5808f684SSatish Balay 4*5808f684SSatish Balayimport numpy as N 5*5808f684SSatish Balay 6*5808f684SSatish Balaydef mkmat(n, mtype, opts): 7*5808f684SSatish Balay A = PETSc.Mat().create(PETSc.COMM_SELF) 8*5808f684SSatish Balay A.setSizes([n,n]) 9*5808f684SSatish Balay A.setType(mtype) 10*5808f684SSatish Balay A.setUp() 11*5808f684SSatish Balay for o in opts: 12*5808f684SSatish Balay A.setOption(o, True) 13*5808f684SSatish Balay return A 14*5808f684SSatish Balay 15*5808f684SSatish Balaydef mksys_diag(n, mtype, opts): 16*5808f684SSatish Balay A = mkmat(n, mtype, opts) 17*5808f684SSatish Balay x, b = A.createVecs() 18*5808f684SSatish Balay for i in range(n): 19*5808f684SSatish Balay A[i,i] = i+1 20*5808f684SSatish Balay x[i] = 1.0/(i+1) 21*5808f684SSatish Balay b[i] = 1 22*5808f684SSatish Balay A.assemble() 23*5808f684SSatish Balay x.assemble() 24*5808f684SSatish Balay b.assemble() 25*5808f684SSatish Balay return A, x, b 26*5808f684SSatish Balay 27*5808f684SSatish Balaydef mksys_poi2(n, mtype, opts): 28*5808f684SSatish Balay A = mkmat(n, mtype, opts) 29*5808f684SSatish Balay x, b = A.createVecs() 30*5808f684SSatish Balay for i in range(n): 31*5808f684SSatish Balay if i == 0: 32*5808f684SSatish Balay cols = [i, i+1] 33*5808f684SSatish Balay vals = [2, -1] 34*5808f684SSatish Balay elif i == n-1: 35*5808f684SSatish Balay cols = [i-1, i] 36*5808f684SSatish Balay vals = [-1, 2] 37*5808f684SSatish Balay else: 38*5808f684SSatish Balay cols = [i-1, i, i+1] 39*5808f684SSatish Balay vals = [-1, 2, -1] 40*5808f684SSatish Balay A[i,cols] = vals 41*5808f684SSatish Balay x[i] = i+1 42*5808f684SSatish Balay b[i] = 0 43*5808f684SSatish Balay A.assemble() 44*5808f684SSatish Balay x.assemble() 45*5808f684SSatish Balay b.assemble() 46*5808f684SSatish Balay A.mult(x,b) 47*5808f684SSatish Balay return A, x, b 48*5808f684SSatish Balay 49*5808f684SSatish Balayclass BaseTestMatFactor(object): 50*5808f684SSatish Balay 51*5808f684SSatish Balay MKSYS = None 52*5808f684SSatish Balay MTYPE = None 53*5808f684SSatish Balay MOPTS = () 54*5808f684SSatish Balay 55*5808f684SSatish Balay def setUp(self): 56*5808f684SSatish Balay A, x, b = self.MKSYS(10, self.MTYPE, self.MOPTS) 57*5808f684SSatish Balay self.A = A 58*5808f684SSatish Balay self.x = x 59*5808f684SSatish Balay self.b = b 60*5808f684SSatish Balay 61*5808f684SSatish Balay def tearDown(self): 62*5808f684SSatish Balay self.A.setUnfactored() 63*5808f684SSatish Balay self.A.destroy(); self.A = None 64*5808f684SSatish Balay self.x.destroy(); self.x = None 65*5808f684SSatish Balay self.b.destroy(); self.b = None 66*5808f684SSatish Balay 67*5808f684SSatish Balayclass BaseTestMatFactorLU(BaseTestMatFactor): 68*5808f684SSatish Balay 69*5808f684SSatish Balay def testFactorLU(self): 70*5808f684SSatish Balay r, c = self.A.getOrdering("nd") 71*5808f684SSatish Balay self.A.reorderForNonzeroDiagonal(r, c) 72*5808f684SSatish Balay self.A.factorLU(r,c,{'zeropivot':1e-5}) 73*5808f684SSatish Balay x = self.x.duplicate() 74*5808f684SSatish Balay self.A.solve(self.b, x) 75*5808f684SSatish Balay x.axpy(-1, self.x) 76*5808f684SSatish Balay self.assertTrue(x.norm() < 1e-3) 77*5808f684SSatish Balay 78*5808f684SSatish Balayclass BaseTestMatFactorILU(BaseTestMatFactor): 79*5808f684SSatish Balay 80*5808f684SSatish Balay def testFactorILU(self): 81*5808f684SSatish Balay r, c = self.A.getOrdering("natural") 82*5808f684SSatish Balay self.A.factorILU(r,c,{'levels':0}) 83*5808f684SSatish Balay x = self.x.duplicate() 84*5808f684SSatish Balay self.A.solve(self.b, x) 85*5808f684SSatish Balay x.axpy(-1, self.x) 86*5808f684SSatish Balay self.assertTrue(x.norm() < 1e-3) 87*5808f684SSatish Balay 88*5808f684SSatish Balay## class BaseTestMatFactorILUDT(BaseTestMatFactor): 89*5808f684SSatish Balay## 90*5808f684SSatish Balay## def testFactorILUDT(self): 91*5808f684SSatish Balay## r, c = self.A.getOrdering("natural") 92*5808f684SSatish Balay## self.A = self.A.factorILUDT(r,c) 93*5808f684SSatish Balay## x = self.x.duplicate() 94*5808f684SSatish Balay## self.A.solve(self.b, x) 95*5808f684SSatish Balay## x.axpy(-1, self.x) 96*5808f684SSatish Balay## self.assertTrue(x.norm() < 1e-3) 97*5808f684SSatish Balay## 98*5808f684SSatish Balayclass BaseTestMatFactorChol(BaseTestMatFactor): 99*5808f684SSatish Balay 100*5808f684SSatish Balay def testFactorChol(self): 101*5808f684SSatish Balay r, c = self.A.getOrdering("natural") 102*5808f684SSatish Balay self.A.factorCholesky(r) 103*5808f684SSatish Balay x = self.x.duplicate() 104*5808f684SSatish Balay self.A.solve(self.b, x) 105*5808f684SSatish Balay x.axpy(-1, self.x) 106*5808f684SSatish Balay self.assertTrue(x.norm() < 1e-3) 107*5808f684SSatish Balay 108*5808f684SSatish Balayclass BaseTestMatFactorICC(BaseTestMatFactor): 109*5808f684SSatish Balay 110*5808f684SSatish Balay def testFactorICC(self): 111*5808f684SSatish Balay r, c = self.A.getOrdering("natural") 112*5808f684SSatish Balay self.A.factorICC(r) 113*5808f684SSatish Balay x = self.x.duplicate() 114*5808f684SSatish Balay self.A.solve(self.b, x) 115*5808f684SSatish Balay x.axpy(-1, self.x) 116*5808f684SSatish Balay self.assertTrue(x.norm() < 1e-3) 117*5808f684SSatish Balay 118*5808f684SSatish Balay 119*5808f684SSatish Balay# -------------------------------------------------------------------- 120*5808f684SSatish Balay 121*5808f684SSatish Balayclass TestMatFactorA1(BaseTestMatFactorLU, 122*5808f684SSatish Balay BaseTestMatFactorChol, 123*5808f684SSatish Balay unittest.TestCase): 124*5808f684SSatish Balay MKSYS = staticmethod(mksys_diag) 125*5808f684SSatish Balay MTYPE = PETSc.Mat.Type.SEQDENSE 126*5808f684SSatish Balay 127*5808f684SSatish Balayclass TestMatFactorA2(BaseTestMatFactorLU, 128*5808f684SSatish Balay BaseTestMatFactorChol, 129*5808f684SSatish Balay unittest.TestCase): 130*5808f684SSatish Balay MKSYS = staticmethod(mksys_poi2) 131*5808f684SSatish Balay MTYPE = PETSc.Mat.Type.SEQDENSE 132*5808f684SSatish Balay 133*5808f684SSatish Balay# --- 134*5808f684SSatish Balay 135*5808f684SSatish Balayclass TestMatFactorB1(BaseTestMatFactorLU, 136*5808f684SSatish Balay BaseTestMatFactorILU, 137*5808f684SSatish Balay ## BaseTestMatFactorILUDT, 138*5808f684SSatish Balay unittest.TestCase): 139*5808f684SSatish Balay MKSYS = staticmethod(mksys_diag) 140*5808f684SSatish Balay MTYPE = PETSc.Mat.Type.SEQAIJ 141*5808f684SSatish Balay 142*5808f684SSatish Balayclass TestMatFactorB2(BaseTestMatFactorLU, 143*5808f684SSatish Balay BaseTestMatFactorILU, 144*5808f684SSatish Balay ## BaseTestMatFactorILUDT, 145*5808f684SSatish Balay unittest.TestCase): 146*5808f684SSatish Balay MKSYS = staticmethod(mksys_poi2) 147*5808f684SSatish Balay MTYPE = PETSc.Mat.Type.SEQAIJ 148*5808f684SSatish Balay 149*5808f684SSatish Balay# --- 150*5808f684SSatish Balay 151*5808f684SSatish Balayclass TestMatFactorC1(BaseTestMatFactorLU, 152*5808f684SSatish Balay BaseTestMatFactorILU, 153*5808f684SSatish Balay unittest.TestCase): 154*5808f684SSatish Balay MKSYS = staticmethod(mksys_diag) 155*5808f684SSatish Balay MTYPE = PETSc.Mat.Type.SEQBAIJ 156*5808f684SSatish Balay 157*5808f684SSatish Balayclass TestMatFactorC2(BaseTestMatFactorLU, 158*5808f684SSatish Balay BaseTestMatFactorILU, 159*5808f684SSatish Balay unittest.TestCase): 160*5808f684SSatish Balay MKSYS = staticmethod(mksys_poi2) 161*5808f684SSatish Balay MTYPE = PETSc.Mat.Type.SEQBAIJ 162*5808f684SSatish Balay 163*5808f684SSatish Balay# --- 164*5808f684SSatish Balay 165*5808f684SSatish Balayclass TestMatFactorD1(BaseTestMatFactorChol, 166*5808f684SSatish Balay BaseTestMatFactorICC, 167*5808f684SSatish Balay unittest.TestCase): 168*5808f684SSatish Balay MKSYS = staticmethod(mksys_diag) 169*5808f684SSatish Balay MTYPE = PETSc.Mat.Type.SEQSBAIJ 170*5808f684SSatish Balay MOPTS = [PETSc.Mat.Option.IGNORE_LOWER_TRIANGULAR] 171*5808f684SSatish Balay 172*5808f684SSatish Balayclass TestMatFactorD2(BaseTestMatFactorChol, 173*5808f684SSatish Balay BaseTestMatFactorICC, 174*5808f684SSatish Balay unittest.TestCase): 175*5808f684SSatish Balay MKSYS = staticmethod(mksys_poi2) 176*5808f684SSatish Balay MTYPE = PETSc.Mat.Type.SEQSBAIJ 177*5808f684SSatish Balay MOPTS = [PETSc.Mat.Option.IGNORE_LOWER_TRIANGULAR] 178*5808f684SSatish Balay 179*5808f684SSatish Balay# -------------------------------------------------------------------- 180*5808f684SSatish Balay 181*5808f684SSatish Balayif __name__ == '__main__': 182*5808f684SSatish Balay unittest.main() 183