xref: /petsc/src/binding/petsc4py/demo/legacy/ode/rober.py (revision d53c7dcfa6b151504f6ca0aeb7a446ad68fe4d68)
155a74a43SLisandro Dalcin# Stiff 3-variable ODE system from chemical reactions,
255a74a43SLisandro Dalcin# due to Robertson (1966),
355a74a43SLisandro Dalcin# problem ROBER in Hairer&Wanner, ODE 2, 1996
455a74a43SLisandro Dalcin
5*69777137SStefano Zampiniimport sys
6*69777137SStefano Zampiniimport petsc4py
7*69777137SStefano Zampini
855a74a43SLisandro Dalcinpetsc4py.init(sys.argv)
955a74a43SLisandro Dalcin
1055a74a43SLisandro Dalcinfrom petsc4py import PETSc
1155a74a43SLisandro Dalcin
12*69777137SStefano Zampini
13*69777137SStefano Zampiniclass Rober:
1455a74a43SLisandro Dalcin    n = 3
1555a74a43SLisandro Dalcin    comm = PETSc.COMM_SELF
16*69777137SStefano Zampini
1755a74a43SLisandro Dalcin    def evalSolution(self, t, x):
18*69777137SStefano Zampini        if t != 0.0:
19*69777137SStefano Zampini            raise ValueError('Only for t=0')
2055a74a43SLisandro Dalcin        x[:] = [1, 0, 0]
2155a74a43SLisandro Dalcin        x.assemble()
22*69777137SStefano Zampini
2355a74a43SLisandro Dalcin    def evalFunction(self, ts, t, x, xdot, f):
24*69777137SStefano Zampini        f[:] = [
25*69777137SStefano Zampini            xdot[0] + 0.04 * x[0] - 1e4 * x[1] * x[2],
2655a74a43SLisandro Dalcin            xdot[1] - 0.04 * x[0] + 1e4 * x[1] * x[2] + 3e7 * x[1] ** 2,
27*69777137SStefano Zampini            xdot[2] - 3e7 * x[1] ** 2,
28*69777137SStefano Zampini        ]
2955a74a43SLisandro Dalcin        f.assemble()
30*69777137SStefano Zampini
3155a74a43SLisandro Dalcin    def evalJacobian(self, ts, t, x, xdot, a, A, B):
3255a74a43SLisandro Dalcin        J = B
33*69777137SStefano Zampini        J[:, :] = [
34*69777137SStefano Zampini            [a + 0.04, -1e4 * x[2], -1e4 * x[1]],
3555a74a43SLisandro Dalcin            [-0.04, a + 1e4 * x[2] + 3e7 * 2 * x[1], 1e4 * x[1]],
36*69777137SStefano Zampini            [0, -3e7 * 2 * x[1], a],
37*69777137SStefano Zampini        ]
3855a74a43SLisandro Dalcin        J.assemble()
39*69777137SStefano Zampini        if A != B:
40*69777137SStefano Zampini            A.assemble()
41*69777137SStefano Zampini
4255a74a43SLisandro Dalcin
4355a74a43SLisandro DalcinOptDB = PETSc.Options()
4455a74a43SLisandro Dalcinode = Rober()
4555a74a43SLisandro Dalcin
4655a74a43SLisandro DalcinJ = PETSc.Mat().createDense([ode.n, ode.n], comm=ode.comm)
4755a74a43SLisandro DalcinJ.setUp()
4855a74a43SLisandro Dalcinx = PETSc.Vec().createSeq(ode.n, comm=ode.comm)
4955a74a43SLisandro Dalcinf = x.duplicate()
5055a74a43SLisandro Dalcin
5155a74a43SLisandro Dalcints = PETSc.TS().create(comm=ode.comm)
5255a74a43SLisandro Dalcints.setProblemType(ts.ProblemType.NONLINEAR)
5355a74a43SLisandro Dalcints.setType(ts.Type.THETA)
5455a74a43SLisandro Dalcin
5555a74a43SLisandro Dalcints.setIFunction(ode.evalFunction, f)
5655a74a43SLisandro Dalcints.setIJacobian(ode.evalJacobian, J)
5755a74a43SLisandro Dalcin
5855a74a43SLisandro Dalcinhistory = []
59*69777137SStefano Zampini
60*69777137SStefano Zampini
6155a74a43SLisandro Dalcindef monitor(ts, i, t, x):
6255a74a43SLisandro Dalcin    xx = x[:].tolist()
6355a74a43SLisandro Dalcin    history.append((i, t, xx))
64*69777137SStefano Zampini
65*69777137SStefano Zampini
6655a74a43SLisandro Dalcints.setMonitor(monitor)
6755a74a43SLisandro Dalcin
6855a74a43SLisandro Dalcints.setTime(0.0)
69*69777137SStefano Zampinits.setTimeStep(0.001)
7055a74a43SLisandro Dalcints.setMaxTime(1e30)
7155a74a43SLisandro Dalcints.setMaxSteps(100)
7255a74a43SLisandro Dalcints.setExactFinalTime(PETSc.TS.ExactFinalTime.INTERPOLATE)
7355a74a43SLisandro Dalcin
7455a74a43SLisandro Dalcints.setFromOptions()
7555a74a43SLisandro Dalcinode.evalSolution(0.0, x)
7655a74a43SLisandro Dalcints.solve(x)
7755a74a43SLisandro Dalcin
7855a74a43SLisandro Dalcindel ode, J, x, ts
7955a74a43SLisandro Dalcin
8055a74a43SLisandro Dalcinfrom matplotlib import pylab
8155a74a43SLisandro Dalcinimport numpy as np
82*69777137SStefano Zampini
8355a74a43SLisandro Dalcinii = np.asarray([v[0] for v in history])
8455a74a43SLisandro Dalcintt = np.asarray([v[1] for v in history])
8555a74a43SLisandro Dalcinxx = np.asarray([v[2] for v in history])
8655a74a43SLisandro Dalcin
8755a74a43SLisandro Dalcinpylab.suptitle('Rober')
8855a74a43SLisandro Dalcinpylab.subplot(2, 2, 1)
8955a74a43SLisandro Dalcinpylab.subplots_adjust(wspace=0.3)
90*69777137SStefano Zampinipylab.semilogy(
91*69777137SStefano Zampini    ii[:-1],
92*69777137SStefano Zampini    np.diff(tt),
93*69777137SStefano Zampini)
9455a74a43SLisandro Dalcinpylab.xlabel('step number')
9555a74a43SLisandro Dalcinpylab.ylabel('timestep')
9655a74a43SLisandro Dalcin
97*69777137SStefano Zampinifor i in range(3):
9855a74a43SLisandro Dalcin    pylab.subplot(2, 2, i + 2)
99*69777137SStefano Zampini    pylab.semilogx(tt, xx[:, i], 'rgb'[i])
10055a74a43SLisandro Dalcin    pylab.xlabel('t')
10155a74a43SLisandro Dalcin    pylab.ylabel('x%d' % i)
10255a74a43SLisandro Dalcin
10355a74a43SLisandro Dalcinpylab.show()
104