xref: /petsc/src/ts/tutorials/ex34.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static const char help[] = "An elastic wave equation driven by Dieterich-Ruina friction\n";
2*c4762a1bSJed Brown /*
3*c4762a1bSJed Brown This whole derivation comes from Erickson, Birnir, and Lavallee [2010]. The model comes from the continuum limit in Carlson and Langer [1989],
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown   u_{tt}   = c^2 u_{xx} - \tilde\gamma^2 u - (\gamma^2 / \xi) (\theta + \ln(u_t + 1))
6*c4762a1bSJed Brown   \theta_t = -(u_t + 1) (\theta + (1 + \epsilon) \ln(u_t +1))
7*c4762a1bSJed Brown 
8*c4762a1bSJed Brown which can be reduced to a first order system,
9*c4762a1bSJed Brown 
10*c4762a1bSJed Brown   u_t      = v
11*c4762a1bSJed Brown   v_t      = c^2 u_{xx} - \tilde\gamma^2 u - (\gamma^2 / \xi)(\theta + ln(v + 1)))
12*c4762a1bSJed Brown   \theta_t = -(v + 1) (\theta + (1 + \epsilon) \ln(v+1))
13*c4762a1bSJed Brown */
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown #include <petscdm.h>
16*c4762a1bSJed Brown #include <petscdmda.h>
17*c4762a1bSJed Brown #include <petscts.h>
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown typedef struct {
20*c4762a1bSJed Brown   PetscScalar u,v, th;
21*c4762a1bSJed Brown } Field;
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown typedef struct _User *User;
24*c4762a1bSJed Brown struct _User {
25*c4762a1bSJed Brown   PetscReal epsilon;    /* inverse of seismic ratio, B-A / A */
26*c4762a1bSJed Brown   PetscReal gamma;      /* wave frequency for interblock coupling */
27*c4762a1bSJed Brown   PetscReal gammaTilde; /* wave frequency for coupling to plate */
28*c4762a1bSJed Brown   PetscReal xi;         /* interblock spring constant */
29*c4762a1bSJed Brown   PetscReal c;          /* wavespeed */
30*c4762a1bSJed Brown };
31*c4762a1bSJed Brown 
32*c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts, PetscReal t, Vec U, Vec F, void *ctx)
33*c4762a1bSJed Brown {
34*c4762a1bSJed Brown   User              user = (User) ctx;
35*c4762a1bSJed Brown   DM                dm, cdm;
36*c4762a1bSJed Brown   DMDALocalInfo     info;
37*c4762a1bSJed Brown   Vec               C;
38*c4762a1bSJed Brown   Field             *f;
39*c4762a1bSJed Brown   const Field       *u;
40*c4762a1bSJed Brown   const PetscScalar *x;
41*c4762a1bSJed Brown   PetscInt          i;
42*c4762a1bSJed Brown   PetscErrorCode    ierr;
43*c4762a1bSJed Brown 
44*c4762a1bSJed Brown   PetscFunctionBeginUser;
45*c4762a1bSJed Brown   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
47*c4762a1bSJed Brown   ierr = DMGetCoordinatesLocal(dm, &C);CHKERRQ(ierr);
48*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(dm, &info);CHKERRQ(ierr);
49*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(dm,  U, (void*)&u);CHKERRQ(ierr);
50*c4762a1bSJed Brown   ierr = DMDAVecGetArray(dm,  F, &f);CHKERRQ(ierr);
51*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(cdm, C, (void*)&x);CHKERRQ(ierr);
52*c4762a1bSJed Brown   for (i = info.xs; i < info.xs+info.xm; ++i) {
53*c4762a1bSJed Brown     const PetscScalar hx = i+1 == info.xs+info.xm ? x[i] - x[i-1] : x[i+1] - x[i];
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown     f[i].u  =  hx*(u[i].v);
56*c4762a1bSJed Brown     f[i].v  = -hx*(PetscSqr(user->gammaTilde)*u[i].u + (PetscSqr(user->gamma) / user->xi)*(u[i].th + PetscLogScalar(u[i].v + 1)));
57*c4762a1bSJed Brown     f[i].th = -hx*(u[i].v + 1)*(u[i].th + (1 + user->epsilon)*PetscLogScalar(u[i].v + 1));
58*c4762a1bSJed Brown   }
59*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(dm,  U, (void*)&u);CHKERRQ(ierr);
60*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(dm,  F, &f);CHKERRQ(ierr);
61*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(cdm, C, (void*)&x);CHKERRQ(ierr);
62*c4762a1bSJed Brown   PetscFunctionReturn(0);
63*c4762a1bSJed Brown }
64*c4762a1bSJed Brown 
65*c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx)
66*c4762a1bSJed Brown {
67*c4762a1bSJed Brown   User           user = (User) ctx;
68*c4762a1bSJed Brown   DM             dm, cdm;
69*c4762a1bSJed Brown   DMDALocalInfo  info;
70*c4762a1bSJed Brown   Vec            Uloc, C;
71*c4762a1bSJed Brown   Field         *u, *udot, *f;
72*c4762a1bSJed Brown   PetscScalar   *x;
73*c4762a1bSJed Brown   PetscInt       i;
74*c4762a1bSJed Brown   PetscErrorCode ierr;
75*c4762a1bSJed Brown 
76*c4762a1bSJed Brown   PetscFunctionBeginUser;
77*c4762a1bSJed Brown   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
78*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(dm, &info);CHKERRQ(ierr);
79*c4762a1bSJed Brown   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
80*c4762a1bSJed Brown   ierr = DMGetCoordinatesLocal(dm, &C);CHKERRQ(ierr);
81*c4762a1bSJed Brown   ierr = DMGetLocalVector(dm, &Uloc);CHKERRQ(ierr);
82*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(dm, U, INSERT_VALUES, Uloc);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(dm, U, INSERT_VALUES, Uloc);CHKERRQ(ierr);
84*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(dm,  Uloc, &u);CHKERRQ(ierr);
85*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(dm,  Udot, &udot);CHKERRQ(ierr);
86*c4762a1bSJed Brown   ierr = DMDAVecGetArray(dm,  F,    &f);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(cdm, C,    &x);CHKERRQ(ierr);
88*c4762a1bSJed Brown   for (i = info.xs; i < info.xs+info.xm; ++i) {
89*c4762a1bSJed Brown     if (i == 0) {
90*c4762a1bSJed Brown       const PetscScalar hx = x[i+1] - x[i];
91*c4762a1bSJed Brown       f[i].u  = hx * udot[i].u;
92*c4762a1bSJed Brown       f[i].v  = hx * udot[i].v - PetscSqr(user->c) * (u[i+1].u - u[i].u) / hx;
93*c4762a1bSJed Brown       f[i].th = hx * udot[i].th;
94*c4762a1bSJed Brown     } else if (i == info.mx-1) {
95*c4762a1bSJed Brown       const PetscScalar hx = x[i] - x[i-1];
96*c4762a1bSJed Brown       f[i].u  = hx * udot[i].u;
97*c4762a1bSJed Brown       f[i].v  = hx * udot[i].v - PetscSqr(user->c) * (u[i-1].u - u[i].u) / hx;
98*c4762a1bSJed Brown       f[i].th = hx * udot[i].th;
99*c4762a1bSJed Brown     } else {
100*c4762a1bSJed Brown       const PetscScalar hx = x[i+1] - x[i];
101*c4762a1bSJed Brown       f[i].u  = hx * udot[i].u;
102*c4762a1bSJed Brown       f[i].v  = hx * udot[i].v - PetscSqr(user->c) * (u[i-1].u - 2.*u[i].u + u[i+1].u) / hx;
103*c4762a1bSJed Brown       f[i].th = hx * udot[i].th;
104*c4762a1bSJed Brown     }
105*c4762a1bSJed Brown   }
106*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(dm,  Uloc, &u);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(dm,  Udot, &udot);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(dm,  F,    &f);CHKERRQ(ierr);
109*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(cdm, C,    &x);CHKERRQ(ierr);
110*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(dm, &Uloc);CHKERRQ(ierr);
111*c4762a1bSJed Brown   PetscFunctionReturn(0);
112*c4762a1bSJed Brown }
113*c4762a1bSJed Brown 
114*c4762a1bSJed Brown /* IJacobian - Compute IJacobian = dF/dU + a dF/dUdot */
115*c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat J, Mat Jpre, void *ctx)
116*c4762a1bSJed Brown {
117*c4762a1bSJed Brown   User           user = (User) ctx;
118*c4762a1bSJed Brown   DM             dm, cdm;
119*c4762a1bSJed Brown   DMDALocalInfo  info;
120*c4762a1bSJed Brown   Vec            C;
121*c4762a1bSJed Brown   Field         *u, *udot;
122*c4762a1bSJed Brown   PetscScalar   *x;
123*c4762a1bSJed Brown   PetscInt       i;
124*c4762a1bSJed Brown   PetscErrorCode ierr;
125*c4762a1bSJed Brown 
126*c4762a1bSJed Brown   PetscFunctionBeginUser;
127*c4762a1bSJed Brown   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
128*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(dm, &info);CHKERRQ(ierr);
129*c4762a1bSJed Brown   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
130*c4762a1bSJed Brown   ierr = DMGetCoordinatesLocal(dm, &C);CHKERRQ(ierr);
131*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(dm,  U,    &u);CHKERRQ(ierr);
132*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(dm,  Udot, &udot);CHKERRQ(ierr);
133*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(cdm, C,    &x);CHKERRQ(ierr);
134*c4762a1bSJed Brown   for (i = info.xs; i < info.xs+info.xm; ++i) {
135*c4762a1bSJed Brown     if (i == 0) {
136*c4762a1bSJed Brown       const PetscScalar hx            = x[i+1] - x[i];
137*c4762a1bSJed Brown       const PetscInt    row           = i, col[] = {i,i+1};
138*c4762a1bSJed Brown       const PetscScalar dxx0          = PetscSqr(user->c)/hx,dxxR = -PetscSqr(user->c)/hx;
139*c4762a1bSJed Brown       const PetscScalar vals[3][2][3] = {{{a*hx,     0,0},{0,0,   0}},
140*c4762a1bSJed Brown                                          {{0,a*hx+dxx0,0},{0,dxxR,0}},
141*c4762a1bSJed Brown                                          {{0,0,     a*hx},{0,0,   0}}};
142*c4762a1bSJed Brown 
143*c4762a1bSJed Brown       ierr = MatSetValuesBlocked(Jpre, 1, &row, 2, col, &vals[0][0][0], INSERT_VALUES);CHKERRQ(ierr);
144*c4762a1bSJed Brown     } else if (i == info.mx-1) {
145*c4762a1bSJed Brown       const PetscScalar hx            = x[i+1] - x[i];
146*c4762a1bSJed Brown       const PetscInt    row           = i, col[] = {i-1,i};
147*c4762a1bSJed Brown       const PetscScalar dxxL          = -PetscSqr(user->c)/hx, dxx0 = PetscSqr(user->c)/hx;
148*c4762a1bSJed Brown       const PetscScalar vals[3][2][3] = {{{0,0,   0},{a*hx,     0,0}},
149*c4762a1bSJed Brown                                          {{0,dxxL,0},{0,a*hx+dxx0,0}},
150*c4762a1bSJed Brown                                          {{0,0,   0},{0,0,     a*hx}}};
151*c4762a1bSJed Brown 
152*c4762a1bSJed Brown       ierr = MatSetValuesBlocked(Jpre, 1, &row, 2, col, &vals[0][0][0], INSERT_VALUES);CHKERRQ(ierr);
153*c4762a1bSJed Brown     } else {
154*c4762a1bSJed Brown       const PetscScalar hx            = x[i+1] - x[i];
155*c4762a1bSJed Brown       const PetscInt    row           = i, col[] = {i-1,i,i+1};
156*c4762a1bSJed Brown       const PetscScalar dxxL          = -PetscSqr(user->c)/hx, dxx0 = 2.*PetscSqr(user->c)/hx,dxxR = -PetscSqr(user->c)/hx;
157*c4762a1bSJed Brown       const PetscScalar vals[3][3][3] = {{{0,0,   0},{a*hx,     0,0},{0,0,   0}},
158*c4762a1bSJed Brown                                          {{0,dxxL,0},{0,a*hx+dxx0,0},{0,dxxR,0}},
159*c4762a1bSJed Brown                                          {{0,0,   0},{0,0,     a*hx},{0,0,   0}}};
160*c4762a1bSJed Brown 
161*c4762a1bSJed Brown       ierr = MatSetValuesBlocked(Jpre, 1, &row, 3, col, &vals[0][0][0], INSERT_VALUES);CHKERRQ(ierr);
162*c4762a1bSJed Brown     }
163*c4762a1bSJed Brown   }
164*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(dm,  U,    &u);CHKERRQ(ierr);
165*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(dm,  Udot, &udot);CHKERRQ(ierr);
166*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(cdm, C,    &x);CHKERRQ(ierr);
167*c4762a1bSJed Brown   ierr = MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
168*c4762a1bSJed Brown   ierr = MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
169*c4762a1bSJed Brown   if (J != Jpre) {
170*c4762a1bSJed Brown     ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
171*c4762a1bSJed Brown     ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
172*c4762a1bSJed Brown   }
173*c4762a1bSJed Brown   PetscFunctionReturn(0);
174*c4762a1bSJed Brown }
175*c4762a1bSJed Brown 
176*c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts, Vec U, void *ctx)
177*c4762a1bSJed Brown {
178*c4762a1bSJed Brown   /* User            user = (User) ctx; */
179*c4762a1bSJed Brown   DM              dm, cdm;
180*c4762a1bSJed Brown   DMDALocalInfo   info;
181*c4762a1bSJed Brown   Vec             C;
182*c4762a1bSJed Brown   Field          *u;
183*c4762a1bSJed Brown   PetscScalar    *x;
184*c4762a1bSJed Brown   const PetscReal sigma = 1.0;
185*c4762a1bSJed Brown   PetscInt        i;
186*c4762a1bSJed Brown   PetscErrorCode  ierr;
187*c4762a1bSJed Brown 
188*c4762a1bSJed Brown   PetscFunctionBeginUser;
189*c4762a1bSJed Brown   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
190*c4762a1bSJed Brown   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
191*c4762a1bSJed Brown   ierr = DMGetCoordinatesLocal(dm, &C);CHKERRQ(ierr);
192*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(dm, &info);CHKERRQ(ierr);
193*c4762a1bSJed Brown   ierr = DMDAVecGetArray(dm,  U, &u);CHKERRQ(ierr);
194*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(cdm, C, &x);CHKERRQ(ierr);
195*c4762a1bSJed Brown   for (i = info.xs; i < info.xs+info.xm; ++i) {
196*c4762a1bSJed Brown     u[i].u  = 1.5 * PetscExpScalar(-PetscSqr(x[i] - 10)/PetscSqr(sigma));
197*c4762a1bSJed Brown     u[i].v  = 0.0;
198*c4762a1bSJed Brown     u[i].th = 0.0;
199*c4762a1bSJed Brown   }
200*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(dm,  U, &u);CHKERRQ(ierr);
201*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(cdm, C, &x);CHKERRQ(ierr);
202*c4762a1bSJed Brown   PetscFunctionReturn(0);
203*c4762a1bSJed Brown }
204*c4762a1bSJed Brown 
205*c4762a1bSJed Brown int main(int argc, char **argv)
206*c4762a1bSJed Brown {
207*c4762a1bSJed Brown   DM                dm;
208*c4762a1bSJed Brown   TS                ts;
209*c4762a1bSJed Brown   Vec               X;
210*c4762a1bSJed Brown   Mat               J;
211*c4762a1bSJed Brown   PetscInt          steps, mx;
212*c4762a1bSJed Brown   PetscReal         ftime, hx, dt;
213*c4762a1bSJed Brown   TSConvergedReason reason;
214*c4762a1bSJed Brown   struct _User      user;
215*c4762a1bSJed Brown   PetscErrorCode    ierr;
216*c4762a1bSJed Brown 
217*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
218*c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 11, 3, 1, NULL, &dm);CHKERRQ(ierr);
219*c4762a1bSJed Brown   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
220*c4762a1bSJed Brown   ierr = DMSetUp(dm);CHKERRQ(ierr);
221*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(dm, 0.0, 20.0, 0.0, 0.0, 0.0, 0.0);CHKERRQ(ierr);
222*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm, &X);CHKERRQ(ierr);
223*c4762a1bSJed Brown 
224*c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Dynamic Friction Options", "");
225*c4762a1bSJed Brown   {
226*c4762a1bSJed Brown     user.epsilon    = 0.1;
227*c4762a1bSJed Brown     user.gamma      = 0.5;
228*c4762a1bSJed Brown     user.gammaTilde = 0.5;
229*c4762a1bSJed Brown     user.xi         = 0.5;
230*c4762a1bSJed Brown     user.c          = 0.5;
231*c4762a1bSJed Brown     ierr = PetscOptionsReal("-epsilon", "Inverse of seismic ratio", "", user.epsilon, &user.epsilon, NULL);CHKERRQ(ierr);
232*c4762a1bSJed Brown     ierr = PetscOptionsReal("-gamma", "Wave frequency for interblock coupling", "", user.gamma, &user.gamma, NULL);CHKERRQ(ierr);
233*c4762a1bSJed Brown     ierr = PetscOptionsReal("-gamma_tilde", "Wave frequency for plate coupling", "", user.gammaTilde, &user.gammaTilde, NULL);CHKERRQ(ierr);
234*c4762a1bSJed Brown     ierr = PetscOptionsReal("-xi", "Interblock spring constant", "", user.xi, &user.xi, NULL);CHKERRQ(ierr);
235*c4762a1bSJed Brown     ierr = PetscOptionsReal("-c", "Wavespeed", "", user.c, &user.c, NULL);CHKERRQ(ierr);
236*c4762a1bSJed Brown   }
237*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
238*c4762a1bSJed Brown 
239*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr);
240*c4762a1bSJed Brown   ierr = TSSetDM(ts, dm);CHKERRQ(ierr);
241*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts, NULL, FormRHSFunction, &user);CHKERRQ(ierr);
242*c4762a1bSJed Brown   ierr = TSSetIFunction(ts, NULL, FormIFunction, &user);CHKERRQ(ierr);
243*c4762a1bSJed Brown   ierr = DMSetMatType(dm, MATAIJ);CHKERRQ(ierr);
244*c4762a1bSJed Brown   ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr);
245*c4762a1bSJed Brown   ierr = TSSetIJacobian(ts, J, J, FormIJacobian, &user);CHKERRQ(ierr);
246*c4762a1bSJed Brown 
247*c4762a1bSJed Brown   ftime = 800.0;
248*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr);
249*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
250*c4762a1bSJed Brown   ierr = FormInitialSolution(ts, X, &user);CHKERRQ(ierr);
251*c4762a1bSJed Brown   ierr = TSSetSolution(ts, X);CHKERRQ(ierr);
252*c4762a1bSJed Brown   ierr = VecGetSize(X, &mx);CHKERRQ(ierr);
253*c4762a1bSJed Brown   hx   = 20.0/(PetscReal)(mx-1);
254*c4762a1bSJed Brown   dt   = 0.4 * PetscSqr(hx) / PetscSqr(user.c); /* Diffusive stability limit */
255*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
256*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
257*c4762a1bSJed Brown 
258*c4762a1bSJed Brown   ierr = TSSolve(ts, X);CHKERRQ(ierr);
259*c4762a1bSJed Brown   ierr = TSGetSolveTime(ts, &ftime);CHKERRQ(ierr);
260*c4762a1bSJed Brown   ierr = TSGetStepNumber(ts, &steps);CHKERRQ(ierr);
261*c4762a1bSJed Brown   ierr = TSGetConvergedReason(ts, &reason);CHKERRQ(ierr);
262*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %D steps\n", TSConvergedReasons[reason], (double)ftime, steps);CHKERRQ(ierr);
263*c4762a1bSJed Brown 
264*c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
265*c4762a1bSJed Brown   ierr = VecDestroy(&X);CHKERRQ(ierr);
266*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
267*c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
268*c4762a1bSJed Brown   ierr = PetscFinalize();
269*c4762a1bSJed Brown   return ierr;
270*c4762a1bSJed Brown }
271*c4762a1bSJed Brown 
272*c4762a1bSJed Brown /*TEST
273*c4762a1bSJed Brown 
274*c4762a1bSJed Brown     build:
275*c4762a1bSJed Brown       requires: !single  !complex c99
276*c4762a1bSJed Brown 
277*c4762a1bSJed Brown     test:
278*c4762a1bSJed Brown       TODO: broken, was not nightly tested, SNES solve eventually fails, -snes_test_jacobian indicates Jacobian is wrong, but even -snes_mf_operator fails
279*c4762a1bSJed Brown 
280*c4762a1bSJed Brown TEST*/
281