15808f684SSatish Balayimport unittest 25808f684SSatish Balayfrom petsc4py import PETSc 35808f684SSatish Balayfrom sys import getrefcount 45808f684SSatish Balay 55808f684SSatish Balay# -------------------------------------------------------------------- 65808f684SSatish Balay 75808f684SSatish Balayclass MyODE: 85808f684SSatish Balay """ 95808f684SSatish Balay du/dt + u**2 = 0; 105808f684SSatish Balay u0,u1,u2 = 1,2,3 115808f684SSatish Balay """ 125808f684SSatish Balay def __init__(self): 135808f684SSatish Balay self.rhsfunction_calls = 0 145808f684SSatish Balay self.rhsjacobian_calls = 0 155808f684SSatish Balay self.ifunction_calls = 0 165808f684SSatish Balay self.ijacobian_calls = 0 175808f684SSatish Balay self.presolve_calls = 0 185808f684SSatish Balay self.update_calls = 0 195808f684SSatish Balay self.postsolve_calls = 0 205808f684SSatish Balay self.monitor_calls = 0 215808f684SSatish Balay 225808f684SSatish Balay def rhsfunction(self,ts,t,u,F): 235808f684SSatish Balay # print ('MyODE.rhsfunction()') 245808f684SSatish Balay self.rhsfunction_calls += 1 255808f684SSatish Balay f = -(u * u) 265808f684SSatish Balay f.copy(F) 275808f684SSatish Balay 285808f684SSatish Balay def rhsjacobian(self,ts,t,u,J,P): 295808f684SSatish Balay # print ('MyODE.rhsjacobian()') 305808f684SSatish Balay self.rhsjacobian_calls += 1 315808f684SSatish Balay P.zeroEntries() 325808f684SSatish Balay diag = -2 * u 335808f684SSatish Balay P.setDiagonal(diag) 345808f684SSatish Balay P.assemble() 355808f684SSatish Balay if J != P: J.assemble() 365808f684SSatish Balay return True # same_nz 375808f684SSatish Balay 385808f684SSatish Balay def ifunction(self,ts,t,u,du,F): 395808f684SSatish Balay # print ('MyODE.ifunction()') 405808f684SSatish Balay self.ifunction_calls += 1 415808f684SSatish Balay f = du + u * u 425808f684SSatish Balay f.copy(F) 435808f684SSatish Balay 445808f684SSatish Balay def ijacobian(self,ts,t,u,du,a,J,P): 455808f684SSatish Balay # print ('MyODE.ijacobian()') 465808f684SSatish Balay self.ijacobian_calls += 1 475808f684SSatish Balay P.zeroEntries() 485808f684SSatish Balay diag = a + 2 * u 495808f684SSatish Balay P.setDiagonal(diag) 505808f684SSatish Balay P.assemble() 515808f684SSatish Balay if J != P: J.assemble() 525808f684SSatish Balay return True # same_nz 535808f684SSatish Balay 545808f684SSatish Balay def monitor(self, ts, s, t, u): 555808f684SSatish Balay self.monitor_calls += 1 565808f684SSatish Balay dt = ts.time_step 575808f684SSatish Balay ut = ts.vec_sol.norm() 585808f684SSatish Balay #prn = PETSc.Sys.Print 595808f684SSatish Balay #prn('TS: step %2d, T:%f, dT:%f, u:%f' % (s,t,dt,ut)) 605808f684SSatish Balay 615808f684SSatish Balay 625808f684SSatish Balayclass BaseTestTSNonlinear(object): 635808f684SSatish Balay 645808f684SSatish Balay TYPE = None 655808f684SSatish Balay 665808f684SSatish Balay def setUp(self): 675808f684SSatish Balay self.ts = PETSc.TS().create(PETSc.COMM_SELF) 685808f684SSatish Balay eft = PETSc.TS.ExactFinalTime.STEPOVER 695808f684SSatish Balay self.ts.setExactFinalTime(eft) 705808f684SSatish Balay ptype = PETSc.TS.ProblemType.NONLINEAR 715808f684SSatish Balay self.ts.setProblemType(ptype) 725808f684SSatish Balay self.ts.setType(self.TYPE) 735808f684SSatish Balay if PETSc.ScalarType().dtype.char in 'fF': 745808f684SSatish Balay snes = self.ts.getSNES() 755808f684SSatish Balay snes.setTolerances(rtol=1e-6) 765808f684SSatish Balay 775808f684SSatish Balay def tearDown(self): 785808f684SSatish Balay self.ts = None 795808f684SSatish Balay 805808f684SSatish Balay 815808f684SSatish Balayclass BaseTestTSNonlinearRHS(BaseTestTSNonlinear): 825808f684SSatish Balay 835808f684SSatish Balay def testSolveRHS(self): 845808f684SSatish Balay ts = self.ts 855808f684SSatish Balay dct = self.ts.getDict() 865808f684SSatish Balay self.assertTrue(dct is not None) 875808f684SSatish Balay self.assertTrue(type(dct) is dict) 885808f684SSatish Balay 895808f684SSatish Balay ode = MyODE() 905808f684SSatish Balay J = PETSc.Mat().create(ts.comm) 915808f684SSatish Balay J.setSizes(3); 925808f684SSatish Balay J.setFromOptions() 935808f684SSatish Balay J.setUp() 945808f684SSatish Balay u, f = J.createVecs() 955808f684SSatish Balay 965808f684SSatish Balay ts.setAppCtx(ode) 975808f684SSatish Balay ts.setRHSFunction(ode.rhsfunction, f) 985808f684SSatish Balay ts.setRHSJacobian(ode.rhsjacobian, J, J) 995808f684SSatish Balay ts.setMonitor(ode.monitor) 1005808f684SSatish Balay 1015808f684SSatish Balay ts.snes.ksp.pc.setType('none') 1025808f684SSatish Balay 1035808f684SSatish Balay T0, dT, nT = 0.00, 0.1, 10 1045808f684SSatish Balay T = T0 + nT*dT 1055808f684SSatish Balay ts.setTime(T0) 1065808f684SSatish Balay ts.setTimeStep(dT) 1075808f684SSatish Balay ts.setMaxTime(T) 1085808f684SSatish Balay ts.setMaxSteps(nT) 1095808f684SSatish Balay ts.setFromOptions() 1105808f684SSatish Balay u[0], u[1], u[2] = 1, 2, 3 1115808f684SSatish Balay ts.solve(u) 1125808f684SSatish Balay 1135808f684SSatish Balay self.assertTrue(ode.rhsfunction_calls > 0) 1145808f684SSatish Balay self.assertTrue(ode.rhsjacobian_calls > 0) 1155808f684SSatish Balay 1165808f684SSatish Balay dct = self.ts.getDict() 1175808f684SSatish Balay self.assertTrue('__appctx__' in dct) 1185808f684SSatish Balay self.assertTrue('__rhsfunction__' in dct) 1195808f684SSatish Balay self.assertTrue('__rhsjacobian__' in dct) 1205808f684SSatish Balay self.assertTrue('__monitor__' in dct) 1215808f684SSatish Balay 1225808f684SSatish Balay n = ode.monitor_calls 1235808f684SSatish Balay ts.monitor(ts.step_number, ts.time) 1245808f684SSatish Balay self.assertEqual(ode.monitor_calls, n+1) 1255808f684SSatish Balay n = ode.monitor_calls 126*1dbd64e7SPierre Jolivet ts.monitorCancel() 1275808f684SSatish Balay ts.monitor(ts.step_number, ts.time) 1285808f684SSatish Balay self.assertEqual(ode.monitor_calls, n) 1295808f684SSatish Balay 1305808f684SSatish Balay def testFDColorRHS(self): 1315808f684SSatish Balay ts = self.ts 1325808f684SSatish Balay ode = MyODE() 1335808f684SSatish Balay J = PETSc.Mat().create(ts.comm) 1345808f684SSatish Balay J.setSizes(5); J.setType('aij') 1355808f684SSatish Balay J.setPreallocationNNZ(nnz=1) 1365808f684SSatish Balay u, f = J.createVecs() 1375808f684SSatish Balay 1385808f684SSatish Balay ts.setAppCtx(ode) 1395808f684SSatish Balay ts.setRHSFunction(ode.rhsfunction, f) 1405808f684SSatish Balay ts.setRHSJacobian(ode.rhsjacobian, J, J) 1415808f684SSatish Balay ts.setMonitor(ode.monitor) 1425808f684SSatish Balay 1435808f684SSatish Balay T0, dT, nT = 0.00, 0.1, 10 1445808f684SSatish Balay T = T0 + nT*dT 1455808f684SSatish Balay ts.setTime(T0) 1465808f684SSatish Balay ts.setTimeStep(dT) 1475808f684SSatish Balay ts.setMaxTime(T) 1485808f684SSatish Balay ts.setMaxSteps(nT) 1495808f684SSatish Balay ts.setFromOptions() 1505808f684SSatish Balay u[0], u[1], u[2] = 1, 2, 3 1515808f684SSatish Balay 1525808f684SSatish Balay ts.setSolution(u) 1535808f684SSatish Balay ode.rhsjacobian(ts,0,u,J,J) 1545808f684SSatish Balay ts.setUp() 1555808f684SSatish Balay ts.snes.setUseFD(True) 1565808f684SSatish Balay ts.solve(u) 1575808f684SSatish Balay 1585808f684SSatish Balay def testResetAndSolveRHS(self): 1595808f684SSatish Balay self.ts.reset() 1605808f684SSatish Balay self.ts.setStepNumber(0) 1615808f684SSatish Balay self.testSolveRHS() 1625808f684SSatish Balay self.ts.reset() 1635808f684SSatish Balay self.ts.setStepNumber(0) 1645808f684SSatish Balay self.testSolveRHS() 1655808f684SSatish Balay self.ts.reset() 1665808f684SSatish Balay 1675808f684SSatish Balayclass BaseTestTSNonlinearI(BaseTestTSNonlinear): 1685808f684SSatish Balay 1695808f684SSatish Balay def testSolveI(self): 1705808f684SSatish Balay ts = self.ts 1715808f684SSatish Balay dct = self.ts.getDict() 1725808f684SSatish Balay self.assertTrue(dct is not None) 1735808f684SSatish Balay self.assertTrue(type(dct) is dict) 1745808f684SSatish Balay 1755808f684SSatish Balay ode = MyODE() 1765808f684SSatish Balay J = PETSc.Mat().create(ts.comm) 1775808f684SSatish Balay J.setSizes(3); 1785808f684SSatish Balay J.setFromOptions() 1795808f684SSatish Balay J.setUp() 1805808f684SSatish Balay u, f = J.createVecs() 1815808f684SSatish Balay 1825808f684SSatish Balay ts.setAppCtx(ode) 1835808f684SSatish Balay ts.setIFunction(ode.ifunction, f) 1845808f684SSatish Balay ts.setIJacobian(ode.ijacobian, J, J) 1855808f684SSatish Balay ts.setMonitor(ode.monitor) 1865808f684SSatish Balay 1875808f684SSatish Balay ts.snes.ksp.pc.setType('none') 1885808f684SSatish Balay 1895808f684SSatish Balay T0, dT, nT = 0.00, 0.1, 10 1905808f684SSatish Balay T = T0 + nT*dT 1915808f684SSatish Balay ts.setTime(T0) 1925808f684SSatish Balay ts.setTimeStep(dT) 1935808f684SSatish Balay ts.setMaxTime(T) 1945808f684SSatish Balay ts.setMaxSteps(nT) 1955808f684SSatish Balay ts.setFromOptions() 1965808f684SSatish Balay u[0], u[1], u[2] = 1, 2, 3 1975808f684SSatish Balay ts.solve(u) 1985808f684SSatish Balay 1995808f684SSatish Balay self.assertTrue(ode.ifunction_calls > 0) 2005808f684SSatish Balay self.assertTrue(ode.ijacobian_calls > 0) 2015808f684SSatish Balay 2025808f684SSatish Balay dct = self.ts.getDict() 2035808f684SSatish Balay self.assertTrue('__appctx__' in dct) 2045808f684SSatish Balay self.assertTrue('__ifunction__' in dct) 2055808f684SSatish Balay self.assertTrue('__ijacobian__' in dct) 2065808f684SSatish Balay self.assertTrue('__monitor__' in dct) 2075808f684SSatish Balay 2085808f684SSatish Balay n = ode.monitor_calls 2095808f684SSatish Balay ts.monitor(ts.step_number, ts.time) 2105808f684SSatish Balay self.assertEqual(ode.monitor_calls, n+1) 2115808f684SSatish Balay n = ode.monitor_calls 212*1dbd64e7SPierre Jolivet ts.monitorCancel() 2135808f684SSatish Balay ts.monitor(ts.step_number, ts.time) 2145808f684SSatish Balay self.assertEqual(ode.monitor_calls, n) 2155808f684SSatish Balay 2165808f684SSatish Balay def testFDColorI(self): 2175808f684SSatish Balay ts = self.ts 2185808f684SSatish Balay ode = MyODE() 2195808f684SSatish Balay J = PETSc.Mat().create(ts.comm) 2205808f684SSatish Balay J.setSizes(5); J.setType('aij') 2215808f684SSatish Balay J.setPreallocationNNZ(nnz=1) 2225808f684SSatish Balay J.setFromOptions() 2235808f684SSatish Balay u, f = J.createVecs() 2245808f684SSatish Balay 2255808f684SSatish Balay ts.setAppCtx(ode) 2265808f684SSatish Balay ts.setIFunction(ode.ifunction, f) 2275808f684SSatish Balay ts.setIJacobian(ode.ijacobian, J, J) 2285808f684SSatish Balay ts.setMonitor(ode.monitor) 2295808f684SSatish Balay 2305808f684SSatish Balay T0, dT, nT = 0.00, 0.1, 10 2315808f684SSatish Balay T = T0 + nT*dT 2325808f684SSatish Balay ts.setTime(T0) 2335808f684SSatish Balay ts.setTimeStep(dT) 2345808f684SSatish Balay ts.setMaxTime(T) 2355808f684SSatish Balay ts.setMaxSteps(nT) 2365808f684SSatish Balay ts.setFromOptions() 2375808f684SSatish Balay u[0], u[1], u[2] = 1, 2, 3 2385808f684SSatish Balay 2395808f684SSatish Balay ts.setSolution(u) 2405808f684SSatish Balay ode.ijacobian(ts,0,u,0*u,1,J,J) 2415808f684SSatish Balay ts.setUp() 2425808f684SSatish Balay ts.snes.setUseFD(True) 2435808f684SSatish Balay ts.solve(u) 2445808f684SSatish Balay 2455808f684SSatish Balay def testResetAndSolveI(self): 2465808f684SSatish Balay self.ts.reset() 2475808f684SSatish Balay self.ts.setStepNumber(0) 2485808f684SSatish Balay self.testSolveI() 2495808f684SSatish Balay self.ts.reset() 2505808f684SSatish Balay self.ts.setStepNumber(0) 2515808f684SSatish Balay self.testSolveI() 2525808f684SSatish Balay self.ts.reset() 2535808f684SSatish Balay 2545808f684SSatish Balayclass TestTSBeuler(BaseTestTSNonlinearRHS,BaseTestTSNonlinearI, 2555808f684SSatish Balay unittest.TestCase): 2565808f684SSatish Balay TYPE = PETSc.TS.Type.BEULER 2575808f684SSatish Balay 2585808f684SSatish Balayclass TestTSCN(BaseTestTSNonlinearRHS,BaseTestTSNonlinearI, 2595808f684SSatish Balay unittest.TestCase): 2605808f684SSatish Balay TYPE = PETSc.TS.Type.CN 2615808f684SSatish Balay 2625808f684SSatish Balayclass TestTSTheta(BaseTestTSNonlinearRHS, BaseTestTSNonlinearI, 2635808f684SSatish Balay unittest.TestCase): 2645808f684SSatish Balay TYPE = PETSc.TS.Type.THETA 2655808f684SSatish Balay 2665808f684SSatish Balayclass TestTSAlpha(BaseTestTSNonlinearRHS, BaseTestTSNonlinearI, 2675808f684SSatish Balay unittest.TestCase): 2685808f684SSatish Balay TYPE = PETSc.TS.Type.ALPHA 2695808f684SSatish Balay 2705808f684SSatish Balay# -------------------------------------------------------------------- 2715808f684SSatish Balay 2725808f684SSatish Balayif __name__ == '__main__': 2735808f684SSatish Balay unittest.main() 2745808f684SSatish Balay 2755808f684SSatish Balay# -------------------------------------------------------------------- 276