15808f684SSatish Balayfrom petsc4py import PETSc 25808f684SSatish Balayimport unittest 35808f684SSatish Balayimport numpy as np 45808f684SSatish Balay 55808f684SSatish Balay 65808f684SSatish Balayclass TestDMShell(unittest.TestCase): 75808f684SSatish Balay COMM = PETSc.COMM_WORLD 85808f684SSatish Balay 95808f684SSatish Balay def setUp(self): 105808f684SSatish Balay self.dm = PETSc.DMShell().create(comm=self.COMM) 115808f684SSatish Balay 125808f684SSatish Balay def tearDown(self): 1362e5d2d2SJDBetteridge self.dm.destroy() 145808f684SSatish Balay self.dm = None 1562e5d2d2SJDBetteridge PETSc.garbage_cleanup() 165808f684SSatish Balay 175808f684SSatish Balay def testSetGlobalVector(self): 185808f684SSatish Balay vec = PETSc.Vec().create(comm=self.COMM) 195808f684SSatish Balay vec.setSizes((10, None)) 205808f684SSatish Balay vec.setUp() 215808f684SSatish Balay self.dm.setGlobalVector(vec) 225808f684SSatish Balay gvec = self.dm.createGlobalVector() 235808f684SSatish Balay self.assertEqual(vec.getSizes(), gvec.getSizes()) 245808f684SSatish Balay self.assertEqual(vec.comm, gvec.comm) 255808f684SSatish Balay 265808f684SSatish Balay def testSetCreateGlobalVector(self): 275808f684SSatish Balay def create_vec(dm): 285808f684SSatish Balay v = PETSc.Vec().create(comm=dm.comm) 295808f684SSatish Balay v.setSizes((10, None)) 305808f684SSatish Balay v.setUp() 315808f684SSatish Balay return v 32*6f336411SStefano Zampini 335808f684SSatish Balay self.dm.setCreateGlobalVector(create_vec) 345808f684SSatish Balay gvec = self.dm.createGlobalVector() 355808f684SSatish Balay self.assertEqual(gvec.comm, self.dm.comm) 365808f684SSatish Balay self.assertEqual(gvec.getLocalSize(), 10) 375808f684SSatish Balay 385808f684SSatish Balay def testSetLocalVector(self): 395808f684SSatish Balay vec = PETSc.Vec().create(comm=PETSc.COMM_SELF) 405808f684SSatish Balay vec.setSizes((1 + 10 * self.COMM.rank, None)) 415808f684SSatish Balay vec.setUp() 425808f684SSatish Balay self.dm.setLocalVector(vec) 435808f684SSatish Balay lvec = self.dm.createLocalVector() 445808f684SSatish Balay self.assertEqual(vec.getSizes(), lvec.getSizes()) 455808f684SSatish Balay lsize, gsize = lvec.getSizes() 465808f684SSatish Balay self.assertEqual(lsize, gsize) 475808f684SSatish Balay self.assertEqual(lvec.comm, PETSc.COMM_SELF) 485808f684SSatish Balay 495808f684SSatish Balay def testSetCreateLocalVector(self): 505808f684SSatish Balay def create_vec(dm): 515808f684SSatish Balay v = PETSc.Vec().create(comm=PETSc.COMM_SELF) 525808f684SSatish Balay v.setSizes((1 + 10 * dm.comm.rank, None)) 535808f684SSatish Balay v.setUp() 545808f684SSatish Balay return v 55*6f336411SStefano Zampini 565808f684SSatish Balay self.dm.setCreateLocalVector(create_vec) 575808f684SSatish Balay lvec = self.dm.createLocalVector() 585808f684SSatish Balay lsize, gsize = lvec.getSizes() 595808f684SSatish Balay self.assertEqual(lsize, gsize) 605808f684SSatish Balay self.assertEqual(lsize, 1 + 10 * self.dm.comm.rank) 615808f684SSatish Balay self.assertEqual(lvec.comm, PETSc.COMM_SELF) 625808f684SSatish Balay 635808f684SSatish Balay def testSetMatrix(self): 645808f684SSatish Balay mat = PETSc.Mat().create(comm=self.COMM) 655808f684SSatish Balay mat.setSizes(((10, None), (2, None))) 665808f684SSatish Balay mat.setUp() 675808f684SSatish Balay mat.assemble() 685808f684SSatish Balay self.dm.setMatrix(mat) 695808f684SSatish Balay nmat = self.dm.createMatrix() 705808f684SSatish Balay self.assertEqual(nmat.getSizes(), mat.getSizes()) 715808f684SSatish Balay 725808f684SSatish Balay def testSetCreateMatrix(self): 735808f684SSatish Balay def create_mat(dm): 745808f684SSatish Balay mat = PETSc.Mat().create(comm=self.COMM) 755808f684SSatish Balay mat.setSizes(((10, None), (2, None))) 765808f684SSatish Balay mat.setUp() 775808f684SSatish Balay return mat 78*6f336411SStefano Zampini 795808f684SSatish Balay self.dm.setCreateMatrix(create_mat) 805808f684SSatish Balay nmat = self.dm.createMatrix() 815808f684SSatish Balay self.assertEqual(nmat.getSizes(), create_mat(self.dm).getSizes()) 825808f684SSatish Balay 835808f684SSatish Balay def testGlobalToLocal(self): 845808f684SSatish Balay def begin(dm, ivec, mode, ovec): 855808f684SSatish Balay if mode == PETSc.InsertMode.INSERT_VALUES: 865808f684SSatish Balay ovec[...] = ivec[...] 875808f684SSatish Balay elif mode == PETSc.InsertMode.ADD_VALUES: 885808f684SSatish Balay ovec[...] += ivec[...] 89*6f336411SStefano Zampini 905808f684SSatish Balay def end(dm, ivec, mode, ovec): 915808f684SSatish Balay pass 92*6f336411SStefano Zampini 935808f684SSatish Balay vec = PETSc.Vec().create(comm=self.COMM) 945808f684SSatish Balay vec.setSizes((10, None)) 955808f684SSatish Balay vec.setUp() 965808f684SSatish Balay vec[...] = self.dm.comm.rank + 1 975808f684SSatish Balay ovec = PETSc.Vec().create(comm=PETSc.COMM_SELF) 985808f684SSatish Balay ovec.setSizes((10, None)) 995808f684SSatish Balay ovec.setUp() 1005808f684SSatish Balay self.dm.setGlobalToLocal(begin, end) 1015808f684SSatish Balay self.dm.globalToLocal(vec, ovec, addv=PETSc.InsertMode.INSERT_VALUES) 1025808f684SSatish Balay self.assertTrue(np.allclose(vec.getArray(), ovec.getArray())) 1035808f684SSatish Balay self.dm.globalToLocal(vec, ovec, addv=PETSc.InsertMode.ADD_VALUES) 1045808f684SSatish Balay self.assertTrue(np.allclose(2 * vec.getArray(), ovec.getArray())) 1055808f684SSatish Balay 1065808f684SSatish Balay def testLocalToGlobal(self): 1075808f684SSatish Balay def begin(dm, ivec, mode, ovec): 1085808f684SSatish Balay if mode == PETSc.InsertMode.INSERT_VALUES: 1095808f684SSatish Balay ovec[...] = ivec[...] 1105808f684SSatish Balay elif mode == PETSc.InsertMode.ADD_VALUES: 1115808f684SSatish Balay ovec[...] += ivec[...] 112*6f336411SStefano Zampini 1135808f684SSatish Balay def end(dm, ivec, mode, ovec): 1145808f684SSatish Balay pass 115*6f336411SStefano Zampini 1165808f684SSatish Balay vec = PETSc.Vec().create(comm=PETSc.COMM_SELF) 1175808f684SSatish Balay vec.setSizes((10, None)) 1185808f684SSatish Balay vec.setUp() 1195808f684SSatish Balay vec[...] = self.dm.comm.rank + 1 1205808f684SSatish Balay ovec = PETSc.Vec().create(comm=self.COMM) 1215808f684SSatish Balay ovec.setSizes((10, None)) 1225808f684SSatish Balay ovec.setUp() 1235808f684SSatish Balay self.dm.setLocalToGlobal(begin, end) 1245808f684SSatish Balay self.dm.localToGlobal(vec, ovec, addv=PETSc.InsertMode.INSERT_VALUES) 1255808f684SSatish Balay self.assertTrue(np.allclose(vec.getArray(), ovec.getArray())) 1265808f684SSatish Balay self.dm.localToGlobal(vec, ovec, addv=PETSc.InsertMode.ADD_VALUES) 1275808f684SSatish Balay self.assertTrue(np.allclose(2 * vec.getArray(), ovec.getArray())) 1285808f684SSatish Balay 1295808f684SSatish Balay def testLocalToLocal(self): 1305808f684SSatish Balay def begin(dm, ivec, mode, ovec): 1315808f684SSatish Balay if mode == PETSc.InsertMode.INSERT_VALUES: 1325808f684SSatish Balay ovec[...] = ivec[...] 1335808f684SSatish Balay elif mode == PETSc.InsertMode.ADD_VALUES: 1345808f684SSatish Balay ovec[...] += ivec[...] 135*6f336411SStefano Zampini 1365808f684SSatish Balay def end(dm, ivec, mode, ovec): 1375808f684SSatish Balay pass 138*6f336411SStefano Zampini 1395808f684SSatish Balay vec = PETSc.Vec().create(comm=PETSc.COMM_SELF) 1405808f684SSatish Balay vec.setSizes((10, None)) 1415808f684SSatish Balay vec.setUp() 1425808f684SSatish Balay vec[...] = self.dm.comm.rank + 1 1435808f684SSatish Balay ovec = vec.duplicate() 1445808f684SSatish Balay self.dm.setLocalToLocal(begin, end) 1455808f684SSatish Balay self.dm.localToLocal(vec, ovec, addv=PETSc.InsertMode.INSERT_VALUES) 1465808f684SSatish Balay self.assertTrue(np.allclose(vec.getArray(), ovec.getArray())) 1475808f684SSatish Balay self.dm.localToLocal(vec, ovec, addv=PETSc.InsertMode.ADD_VALUES) 1485808f684SSatish Balay self.assertTrue(np.allclose(2 * vec.getArray(), ovec.getArray())) 1495808f684SSatish Balay 1505808f684SSatish Balay def testGlobalToLocalVecScatter(self): 1515808f684SSatish Balay vec = PETSc.Vec().create() 1525808f684SSatish Balay vec.setSizes((10, None)) 1535808f684SSatish Balay vec.setUp() 1545808f684SSatish Balay sct, ovec = PETSc.Scatter.toAll(vec) 1555808f684SSatish Balay self.dm.setGlobalToLocalVecScatter(sct) 1565808f684SSatish Balay self.dm.globalToLocal(vec, ovec, addv=PETSc.InsertMode.INSERT_VALUES) 1575808f684SSatish Balay 1585808f684SSatish Balay def testLocalToGlobalVecScatter(self): 1595808f684SSatish Balay vec = PETSc.Vec().create() 1605808f684SSatish Balay vec.setSizes((10, None)) 1615808f684SSatish Balay vec.setUp() 1625808f684SSatish Balay sct, ovec = PETSc.Scatter.toAll(vec) 1635808f684SSatish Balay self.dm.setLocalToGlobalVecScatter(sct) 1645808f684SSatish Balay self.dm.localToGlobal(vec, ovec, addv=PETSc.InsertMode.INSERT_VALUES) 1655808f684SSatish Balay 1665808f684SSatish Balay def testLocalToLocalVecScatter(self): 1675808f684SSatish Balay vec = PETSc.Vec().create() 1685808f684SSatish Balay vec.setSizes((10, None)) 1695808f684SSatish Balay vec.setUp() 1705808f684SSatish Balay sct, ovec = PETSc.Scatter.toAll(vec) 1715808f684SSatish Balay self.dm.setLocalToLocalVecScatter(sct) 1725808f684SSatish Balay self.dm.localToLocal(vec, ovec, addv=PETSc.InsertMode.INSERT_VALUES) 1735808f684SSatish Balay 1745808f684SSatish Balay def testCoarsenRefine(self): 1755808f684SSatish Balay cdm = PETSc.DMShell().create(comm=self.COMM) 176*6f336411SStefano Zampini 1775808f684SSatish Balay def coarsen(dm, comm): 1785808f684SSatish Balay return cdm 179*6f336411SStefano Zampini 1805808f684SSatish Balay def refine(dm, comm): 1815808f684SSatish Balay return self.dm 182*6f336411SStefano Zampini 1835808f684SSatish Balay cdm.setRefine(refine) 1845808f684SSatish Balay self.dm.setCoarsen(coarsen) 1855808f684SSatish Balay coarsened = self.dm.coarsen() 1865808f684SSatish Balay self.assertEqual(coarsened, cdm) 1875808f684SSatish Balay refined = coarsened.refine() 1885808f684SSatish Balay self.assertEqual(refined, self.dm) 1895808f684SSatish Balay 1905808f684SSatish Balay def testCreateInterpolation(self): 1915808f684SSatish Balay mat = PETSc.Mat().create() 1925808f684SSatish Balay mat.setSizes(((10, None), (10, None))) 1935808f684SSatish Balay mat.setUp() 1945808f684SSatish Balay vec = PETSc.Vec().create() 1955808f684SSatish Balay vec.setSizes((10, None)) 1965808f684SSatish Balay vec.setUp() 197*6f336411SStefano Zampini 1985808f684SSatish Balay def create_interp(dm, dmf): 1995808f684SSatish Balay return mat, vec 200*6f336411SStefano Zampini 2015808f684SSatish Balay self.dm.setCreateInterpolation(create_interp) 2025808f684SSatish Balay m, v = self.dm.createInterpolation(self.dm) 2035808f684SSatish Balay self.assertEqual(m, mat) 2045808f684SSatish Balay self.assertEqual(v, vec) 2055808f684SSatish Balay 2065808f684SSatish Balay def testCreateInjection(self): 2075808f684SSatish Balay mat = PETSc.Mat().create() 2085808f684SSatish Balay mat.setSizes(((10, None), (10, None))) 2095808f684SSatish Balay mat.setUp() 210*6f336411SStefano Zampini 2115808f684SSatish Balay def create_inject(dm, dmf): 2125808f684SSatish Balay return mat 213*6f336411SStefano Zampini 2145808f684SSatish Balay self.dm.setCreateInjection(create_inject) 2155808f684SSatish Balay m = self.dm.createInjection(self.dm) 2165808f684SSatish Balay self.assertEqual(m, mat) 2175808f684SSatish Balay 2185808f684SSatish Balay 2195808f684SSatish Balayif __name__ == '__main__': 2205808f684SSatish Balay unittest.main() 221