xref: /petsc/src/ts/tutorials/ex20opt_ic.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "Solves a ODE-constrained optimization problem -- finding the optimal initial conditions for the van der Pol equation.\n";
2c4762a1bSJed Brown 
30e3d61c9SBarry Smith /*
4c4762a1bSJed Brown   Notes:
5c4762a1bSJed Brown   This code demonstrates how to solve an ODE-constrained optimization problem with TAO, TSAdjoint and TS.
6c4762a1bSJed Brown   The nonlinear problem is written in an ODE equivalent form.
7c4762a1bSJed Brown   The objective is to minimize the difference between observation and model prediction by finding optimal values for initial conditions.
8c4762a1bSJed Brown   The gradient is computed with the discrete adjoint of an implicit method or an explicit method, see ex20adj.c for details.
9c4762a1bSJed Brown */
10c4762a1bSJed Brown 
11c4762a1bSJed Brown #include <petsctao.h>
12c4762a1bSJed Brown #include <petscts.h>
13c4762a1bSJed Brown 
14c4762a1bSJed Brown typedef struct _n_User *User;
15c4762a1bSJed Brown struct _n_User {
16c4762a1bSJed Brown   TS        ts;
17c4762a1bSJed Brown   PetscReal mu;
18c4762a1bSJed Brown   PetscReal next_output;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   /* Sensitivity analysis support */
21c4762a1bSJed Brown   PetscInt  steps;
22c4762a1bSJed Brown   PetscReal ftime;
23c4762a1bSJed Brown   Mat       A;                    /* Jacobian matrix for ODE */
24c4762a1bSJed Brown   Mat       Jacp;                 /* JacobianP matrix for ODE*/
25c4762a1bSJed Brown   Mat       H;                    /* Hessian matrix for optimization */
26c4762a1bSJed Brown   Vec       U, Lambda[1], Mup[1]; /* first-order adjoint variables */
27c4762a1bSJed Brown   Vec       Lambda2[2];           /* second-order adjoint variables */
28c4762a1bSJed Brown   Vec       Ihp1[1];              /* working space for Hessian evaluations */
29c4762a1bSJed Brown   Vec       Dir;                  /* direction vector */
30c4762a1bSJed Brown   PetscReal ob[2];                /* observation used by the cost function */
31c4762a1bSJed Brown   PetscBool implicitform;         /* implicit ODE? */
32c4762a1bSJed Brown };
33c4762a1bSJed Brown PetscErrorCode Adjoint2(Vec, PetscScalar[], User);
34c4762a1bSJed Brown 
35c4762a1bSJed Brown /* ----------------------- Explicit form of the ODE  -------------------- */
36c4762a1bSJed Brown 
RHSFunction(TS ts,PetscReal t,Vec U,Vec F,PetscCtx ctx)37*2a8381b2SBarry Smith static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, PetscCtx ctx)
38d71ae5a4SJacob Faibussowitsch {
39c4762a1bSJed Brown   User               user = (User)ctx;
40c4762a1bSJed Brown   PetscScalar       *f;
41c4762a1bSJed Brown   const PetscScalar *u;
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   PetscFunctionBeginUser;
449566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
459566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
46c4762a1bSJed Brown   f[0] = u[1];
47c4762a1bSJed Brown   f[1] = user->mu * ((1. - u[0] * u[0]) * u[1] - u[0]);
489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
499566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
51c4762a1bSJed Brown }
52c4762a1bSJed Brown 
RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,PetscCtx ctx)53*2a8381b2SBarry Smith static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, PetscCtx ctx)
54d71ae5a4SJacob Faibussowitsch {
55c4762a1bSJed Brown   User               user     = (User)ctx;
56c4762a1bSJed Brown   PetscReal          mu       = user->mu;
57c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
58c4762a1bSJed Brown   PetscScalar        J[2][2];
59c4762a1bSJed Brown   const PetscScalar *u;
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   PetscFunctionBeginUser;
629566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
63c4762a1bSJed Brown   J[0][0] = 0;
64c4762a1bSJed Brown   J[1][0] = -mu * (2.0 * u[1] * u[0] + 1.);
65c4762a1bSJed Brown   J[0][1] = 1.0;
66c4762a1bSJed Brown   J[1][1] = mu * (1.0 - u[0] * u[0]);
679566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
70c4762a1bSJed Brown   if (A != B) {
719566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
729566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
73c4762a1bSJed Brown   }
749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76c4762a1bSJed Brown }
77c4762a1bSJed Brown 
RHSHessianProductUU(TS ts,PetscReal t,Vec U,Vec * Vl,Vec Vr,Vec * VHV,PetscCtx ctx)78*2a8381b2SBarry Smith static PetscErrorCode RHSHessianProductUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, PetscCtx ctx)
79d71ae5a4SJacob Faibussowitsch {
80c4762a1bSJed Brown   const PetscScalar *vl, *vr, *u;
81c4762a1bSJed Brown   PetscScalar       *vhv;
82c4762a1bSJed Brown   PetscScalar        dJdU[2][2][2] = {{{0}}};
83c4762a1bSJed Brown   PetscInt           i, j, k;
84c4762a1bSJed Brown   User               user = (User)ctx;
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   PetscFunctionBeginUser;
879566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
889566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
899566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
909566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   dJdU[1][0][0] = -2. * user->mu * u[1];
93c4762a1bSJed Brown   dJdU[1][1][0] = -2. * user->mu * u[0];
94c4762a1bSJed Brown   dJdU[1][0][1] = -2. * user->mu * u[0];
95c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
96c4762a1bSJed Brown     vhv[j] = 0;
97c4762a1bSJed Brown     for (k = 0; k < 2; k++)
989371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdU[i][j][k] * vr[k];
99c4762a1bSJed Brown   }
1009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1019566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
1029566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
1039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
1043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
105c4762a1bSJed Brown }
106c4762a1bSJed Brown 
107c4762a1bSJed Brown /* ----------------------- Implicit form of the ODE  -------------------- */
108c4762a1bSJed Brown 
IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,PetscCtx ctx)109*2a8381b2SBarry Smith static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, PetscCtx ctx)
110d71ae5a4SJacob Faibussowitsch {
111c4762a1bSJed Brown   User               user = (User)ctx;
112c4762a1bSJed Brown   const PetscScalar *u, *udot;
113c4762a1bSJed Brown   PetscScalar       *f;
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   PetscFunctionBeginUser;
1169566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1179566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
1189566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
119c4762a1bSJed Brown   f[0] = udot[0] - u[1];
120c4762a1bSJed Brown   f[1] = udot[1] - user->mu * ((1.0 - u[0] * u[0]) * u[1] - u[0]);
1219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
1239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
125c4762a1bSJed Brown }
126c4762a1bSJed Brown 
IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,PetscCtx ctx)127*2a8381b2SBarry Smith static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, PetscCtx ctx)
128d71ae5a4SJacob Faibussowitsch {
129c4762a1bSJed Brown   User               user     = (User)ctx;
130c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
131c4762a1bSJed Brown   PetscScalar        J[2][2];
132c4762a1bSJed Brown   const PetscScalar *u;
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   PetscFunctionBeginUser;
1359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1369371c9d4SSatish Balay   J[0][0] = a;
1379371c9d4SSatish Balay   J[0][1] = -1.0;
1389371c9d4SSatish Balay   J[1][0] = user->mu * (1.0 + 2.0 * u[0] * u[1]);
1399371c9d4SSatish Balay   J[1][1] = a - user->mu * (1.0 - u[0] * u[0]);
1409566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
1419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
144c4762a1bSJed Brown   if (A != B) {
1459566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1469566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
147c4762a1bSJed Brown   }
1483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
149c4762a1bSJed Brown }
150c4762a1bSJed Brown 
151c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
Monitor(TS ts,PetscInt step,PetscReal t,Vec U,PetscCtx ctx)152*2a8381b2SBarry Smith static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec U, PetscCtx ctx)
153d71ae5a4SJacob Faibussowitsch {
154c4762a1bSJed Brown   const PetscScalar *u;
155c4762a1bSJed Brown   PetscReal          tfinal, dt;
156c4762a1bSJed Brown   User               user = (User)ctx;
157c4762a1bSJed Brown   Vec                interpolatedU;
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   PetscFunctionBeginUser;
1609566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
1619566063dSJacob Faibussowitsch   PetscCall(TSGetMaxTime(ts, &tfinal));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
1649566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(U, &interpolatedU));
1659566063dSJacob Faibussowitsch     PetscCall(TSInterpolate(ts, user->next_output, interpolatedU));
1669566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(interpolatedU, &u));
1679371c9d4SSatish Balay     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%g] %" PetscInt_FMT " TS %g (dt = %g) X %g %g\n", (double)user->next_output, step, (double)t, (double)dt, (double)PetscRealPart(u[0]), (double)PetscRealPart(u[1])));
1689566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(interpolatedU, &u));
1699566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&interpolatedU));
170c4762a1bSJed Brown     user->next_output += 0.1;
171c4762a1bSJed Brown   }
1723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
173c4762a1bSJed Brown }
174c4762a1bSJed Brown 
IHessianProductUU(TS ts,PetscReal t,Vec U,Vec * Vl,Vec Vr,Vec * VHV,PetscCtx ctx)175*2a8381b2SBarry Smith static PetscErrorCode IHessianProductUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, PetscCtx ctx)
176d71ae5a4SJacob Faibussowitsch {
177c4762a1bSJed Brown   const PetscScalar *vl, *vr, *u;
178c4762a1bSJed Brown   PetscScalar       *vhv;
179c4762a1bSJed Brown   PetscScalar        dJdU[2][2][2] = {{{0}}};
180c4762a1bSJed Brown   PetscInt           i, j, k;
181c4762a1bSJed Brown   User               user = (User)ctx;
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   PetscFunctionBeginUser;
1849566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1859566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
1869566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
1879566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
188c4762a1bSJed Brown   dJdU[1][0][0] = 2. * user->mu * u[1];
189c4762a1bSJed Brown   dJdU[1][1][0] = 2. * user->mu * u[0];
190c4762a1bSJed Brown   dJdU[1][0][1] = 2. * user->mu * u[0];
191c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
192c4762a1bSJed Brown     vhv[j] = 0;
193c4762a1bSJed Brown     for (k = 0; k < 2; k++)
1949371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdU[i][j][k] * vr[k];
195c4762a1bSJed Brown   }
1969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
1989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
1999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
2003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
201c4762a1bSJed Brown }
202c4762a1bSJed Brown 
203c4762a1bSJed Brown /* ------------------ User-defined routines for TAO -------------------------- */
204c4762a1bSJed Brown 
FormFunctionGradient(Tao tao,Vec IC,PetscReal * f,Vec G,PetscCtx ctx)205*2a8381b2SBarry Smith static PetscErrorCode FormFunctionGradient(Tao tao, Vec IC, PetscReal *f, Vec G, PetscCtx ctx)
206d71ae5a4SJacob Faibussowitsch {
207c4762a1bSJed Brown   User               user_ptr = (User)ctx;
208c4762a1bSJed Brown   TS                 ts       = user_ptr->ts;
209c4762a1bSJed Brown   const PetscScalar *x_ptr;
210c4762a1bSJed Brown   PetscScalar       *y_ptr;
211c4762a1bSJed Brown 
212c4762a1bSJed Brown   PetscFunctionBeginUser;
2139566063dSJacob Faibussowitsch   PetscCall(VecCopy(IC, user_ptr->U)); /* set up the initial condition */
214c4762a1bSJed Brown 
2159566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
2169566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
2179566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.001)); /* can be overwritten by command line options */
2189566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
2199566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, user_ptr->U));
220c4762a1bSJed Brown 
2219566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user_ptr->U, &x_ptr));
2229566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->Lambda[0], &y_ptr));
223c4762a1bSJed Brown   *f       = (x_ptr[0] - user_ptr->ob[0]) * (x_ptr[0] - user_ptr->ob[0]) + (x_ptr[1] - user_ptr->ob[1]) * (x_ptr[1] - user_ptr->ob[1]);
224c4762a1bSJed Brown   y_ptr[0] = 2. * (x_ptr[0] - user_ptr->ob[0]);
225c4762a1bSJed Brown   y_ptr[1] = 2. * (x_ptr[1] - user_ptr->ob[1]);
2269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->Lambda[0], &y_ptr));
2279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user_ptr->U, &x_ptr));
228c4762a1bSJed Brown 
2299566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 1, user_ptr->Lambda, NULL));
2309566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
2319566063dSJacob Faibussowitsch   PetscCall(VecCopy(user_ptr->Lambda[0], G));
2323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
233c4762a1bSJed Brown }
234c4762a1bSJed Brown 
FormHessian(Tao tao,Vec U,Mat H,Mat Hpre,PetscCtx ctx)235*2a8381b2SBarry Smith static PetscErrorCode FormHessian(Tao tao, Vec U, Mat H, Mat Hpre, PetscCtx ctx)
236d71ae5a4SJacob Faibussowitsch {
237c4762a1bSJed Brown   User           user_ptr = (User)ctx;
238c4762a1bSJed Brown   PetscScalar    harr[2];
239c4762a1bSJed Brown   PetscScalar   *x_ptr;
240c4762a1bSJed Brown   const PetscInt rows[2] = {0, 1};
241c4762a1bSJed Brown   PetscInt       col;
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   PetscFunctionBeginUser;
2449566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, user_ptr->U));
2459566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->Dir, &x_ptr));
246c4762a1bSJed Brown   x_ptr[0] = 1.;
247c4762a1bSJed Brown   x_ptr[1] = 0.;
2489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->Dir, &x_ptr));
2499566063dSJacob Faibussowitsch   PetscCall(Adjoint2(user_ptr->U, harr, user_ptr));
250c4762a1bSJed Brown   col = 0;
2519566063dSJacob Faibussowitsch   PetscCall(MatSetValues(H, 2, rows, 1, &col, harr, INSERT_VALUES));
252c4762a1bSJed Brown 
2539566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, user_ptr->U));
2549566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->Dir, &x_ptr));
255c4762a1bSJed Brown   x_ptr[0] = 0.;
256c4762a1bSJed Brown   x_ptr[1] = 1.;
2579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->Dir, &x_ptr));
2589566063dSJacob Faibussowitsch   PetscCall(Adjoint2(user_ptr->U, harr, user_ptr));
259c4762a1bSJed Brown   col = 1;
2609566063dSJacob Faibussowitsch   PetscCall(MatSetValues(H, 2, rows, 1, &col, harr, INSERT_VALUES));
261c4762a1bSJed Brown 
2629566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
2639566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
264c4762a1bSJed Brown   if (H != Hpre) {
2659566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Hpre, MAT_FINAL_ASSEMBLY));
2669566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Hpre, MAT_FINAL_ASSEMBLY));
267c4762a1bSJed Brown   }
2683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
269c4762a1bSJed Brown }
270c4762a1bSJed Brown 
MatrixFreeHessian(Tao tao,Vec U,Mat H,Mat Hpre,PetscCtx ctx)271*2a8381b2SBarry Smith static PetscErrorCode MatrixFreeHessian(Tao tao, Vec U, Mat H, Mat Hpre, PetscCtx ctx)
272d71ae5a4SJacob Faibussowitsch {
273c4762a1bSJed Brown   User user_ptr = (User)ctx;
274c4762a1bSJed Brown 
275c4762a1bSJed Brown   PetscFunctionBeginUser;
2769566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, user_ptr->U));
2773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
278c4762a1bSJed Brown }
279c4762a1bSJed Brown 
280c4762a1bSJed Brown /* ------------ Routines calculating second-order derivatives -------------- */
281c4762a1bSJed Brown 
2820e3d61c9SBarry Smith /*
283c4762a1bSJed Brown   Compute the Hessian-vector product for the cost function using Second-order adjoint
284c4762a1bSJed Brown */
Adjoint2(Vec U,PetscScalar arr[],User ctx)285d71ae5a4SJacob Faibussowitsch PetscErrorCode Adjoint2(Vec U, PetscScalar arr[], User ctx)
286d71ae5a4SJacob Faibussowitsch {
287c4762a1bSJed Brown   TS           ts = ctx->ts;
288c4762a1bSJed Brown   PetscScalar *x_ptr, *y_ptr;
289c4762a1bSJed Brown   Mat          tlmsen;
290c4762a1bSJed Brown 
291c4762a1bSJed Brown   PetscFunctionBeginUser;
2929566063dSJacob Faibussowitsch   PetscCall(TSAdjointReset(ts));
293c4762a1bSJed Brown 
2949566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
2959566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
2969566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.001));
2979566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
2989566063dSJacob Faibussowitsch   PetscCall(TSSetCostHessianProducts(ts, 1, ctx->Lambda2, NULL, ctx->Dir));
2999566063dSJacob Faibussowitsch   PetscCall(TSAdjointSetForward(ts, NULL));
3009566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
301c4762a1bSJed Brown 
302c4762a1bSJed Brown   /* Set terminal conditions for first- and second-order adjonts */
3039566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &x_ptr));
3049566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda[0], &y_ptr));
305c4762a1bSJed Brown   y_ptr[0] = 2. * (x_ptr[0] - ctx->ob[0]);
306c4762a1bSJed Brown   y_ptr[1] = 2. * (x_ptr[1] - ctx->ob[1]);
3079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda[0], &y_ptr));
3089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &x_ptr));
3099566063dSJacob Faibussowitsch   PetscCall(TSForwardGetSensitivities(ts, NULL, &tlmsen));
3109566063dSJacob Faibussowitsch   PetscCall(MatDenseGetColumn(tlmsen, 0, &x_ptr));
3119566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda2[0], &y_ptr));
312c4762a1bSJed Brown   y_ptr[0] = 2. * x_ptr[0];
313c4762a1bSJed Brown   y_ptr[1] = 2. * x_ptr[1];
3149566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda2[0], &y_ptr));
3159566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreColumn(tlmsen, &x_ptr));
316c4762a1bSJed Brown 
3179566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 1, ctx->Lambda, NULL));
318c4762a1bSJed Brown   if (ctx->implicitform) {
3199566063dSJacob Faibussowitsch     PetscCall(TSSetIHessianProduct(ts, ctx->Ihp1, IHessianProductUU, NULL, NULL, NULL, NULL, NULL, NULL, ctx));
320c4762a1bSJed Brown   } else {
3219566063dSJacob Faibussowitsch     PetscCall(TSSetRHSHessianProduct(ts, ctx->Ihp1, RHSHessianProductUU, NULL, NULL, NULL, NULL, NULL, NULL, ctx));
322c4762a1bSJed Brown   }
3239566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
324c4762a1bSJed Brown 
3259566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda2[0], &x_ptr));
326c4762a1bSJed Brown   arr[0] = x_ptr[0];
327c4762a1bSJed Brown   arr[1] = x_ptr[1];
3289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda2[0], &x_ptr));
329c4762a1bSJed Brown 
3309566063dSJacob Faibussowitsch   PetscCall(TSAdjointReset(ts));
3319566063dSJacob Faibussowitsch   PetscCall(TSAdjointResetForward(ts));
3323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
333c4762a1bSJed Brown }
334c4762a1bSJed Brown 
FiniteDiff(Vec U,PetscScalar arr[],User ctx)335d71ae5a4SJacob Faibussowitsch PetscErrorCode FiniteDiff(Vec U, PetscScalar arr[], User ctx)
336d71ae5a4SJacob Faibussowitsch {
337c4762a1bSJed Brown   Vec               Up, G, Gp;
338c4762a1bSJed Brown   const PetscScalar eps = PetscRealConstant(1e-7);
339c4762a1bSJed Brown   PetscScalar      *u;
340c4762a1bSJed Brown   Tao               tao = NULL;
341c4762a1bSJed Brown   PetscReal         f;
342c4762a1bSJed Brown 
343c4762a1bSJed Brown   PetscFunctionBeginUser;
3449566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &Up));
3459566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &G));
3469566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &Gp));
347c4762a1bSJed Brown 
3489566063dSJacob Faibussowitsch   PetscCall(FormFunctionGradient(tao, U, &f, G, ctx));
349c4762a1bSJed Brown 
3509566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, Up));
3519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Up, &u));
352c4762a1bSJed Brown   u[0] += eps;
3539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Up, &u));
3549566063dSJacob Faibussowitsch   PetscCall(FormFunctionGradient(tao, Up, &f, Gp, ctx));
3559566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Gp, -1, G));
3569566063dSJacob Faibussowitsch   PetscCall(VecScale(Gp, 1. / eps));
3579566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Gp, &u));
358c4762a1bSJed Brown   arr[0] = u[0];
359c4762a1bSJed Brown   arr[1] = u[1];
3609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Gp, &u));
361c4762a1bSJed Brown 
3629566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, Up));
3639566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Up, &u));
364c4762a1bSJed Brown   u[1] += eps;
3659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Up, &u));
3669566063dSJacob Faibussowitsch   PetscCall(FormFunctionGradient(tao, Up, &f, Gp, ctx));
3679566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Gp, -1, G));
3689566063dSJacob Faibussowitsch   PetscCall(VecScale(Gp, 1. / eps));
3699566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Gp, &u));
370c4762a1bSJed Brown   arr[2] = u[0];
371c4762a1bSJed Brown   arr[3] = u[1];
3729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Gp, &u));
373c4762a1bSJed Brown 
3749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&G));
3759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Gp));
3769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Up));
3773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
378c4762a1bSJed Brown }
379c4762a1bSJed Brown 
HessianProductMat(Mat mat,Vec svec,Vec y)380d71ae5a4SJacob Faibussowitsch static PetscErrorCode HessianProductMat(Mat mat, Vec svec, Vec y)
381d71ae5a4SJacob Faibussowitsch {
382c4762a1bSJed Brown   User         user_ptr;
383c4762a1bSJed Brown   PetscScalar *y_ptr;
384c4762a1bSJed Brown 
385c4762a1bSJed Brown   PetscFunctionBeginUser;
3869566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &user_ptr));
3879566063dSJacob Faibussowitsch   PetscCall(VecCopy(svec, user_ptr->Dir));
3889566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_ptr));
3899566063dSJacob Faibussowitsch   PetscCall(Adjoint2(user_ptr->U, y_ptr, user_ptr));
3909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_ptr));
3913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
392c4762a1bSJed Brown }
393c4762a1bSJed Brown 
main(int argc,char ** argv)394d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
395d71ae5a4SJacob Faibussowitsch {
396c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE, mf = PETSC_TRUE;
397c4762a1bSJed Brown   PetscInt       mode = 0;
398c4762a1bSJed Brown   PetscMPIInt    size;
399c4762a1bSJed Brown   struct _n_User user;
400c4762a1bSJed Brown   Vec            x; /* working vector for TAO */
401c4762a1bSJed Brown   PetscScalar   *x_ptr, arr[4];
402c4762a1bSJed Brown   PetscScalar    ic1 = 2.2, ic2 = -0.7; /* initial guess for TAO */
403c4762a1bSJed Brown   Tao            tao;
404c4762a1bSJed Brown   KSP            ksp;
405c4762a1bSJed Brown   PC             pc;
406c4762a1bSJed Brown 
407c4762a1bSJed Brown   /* Initialize program */
408327415f7SBarry Smith   PetscFunctionBeginUser;
4099566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
4109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
4113c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
412c4762a1bSJed Brown 
413c4762a1bSJed Brown   /* Set runtime options */
414c4762a1bSJed Brown   user.next_output  = 0.0;
415c4762a1bSJed Brown   user.mu           = 1.0e3;
416c4762a1bSJed Brown   user.steps        = 0;
417c4762a1bSJed Brown   user.ftime        = 0.5;
418c4762a1bSJed Brown   user.implicitform = PETSC_TRUE;
419c4762a1bSJed Brown 
4209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
4219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL));
4229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mode", &mode, NULL));
4239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-ic1", &ic1, NULL));
4249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-ic2", &ic2, NULL));
4259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-my_tao_mf", &mf, NULL)); /* matrix-free hessian for optimization */
4269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-implicitform", &user.implicitform, NULL));
427c4762a1bSJed Brown 
428c4762a1bSJed Brown   /* Create necessary matrix and vectors, solve same ODE on every process */
4299566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.A));
4309566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
4319566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.A));
4329566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.A));
4339566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.U, NULL));
4349566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.Dir, NULL));
4359566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.Lambda[0], NULL));
4369566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.Lambda2[0], NULL));
4379566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.Ihp1[0], NULL));
438c4762a1bSJed Brown 
439c4762a1bSJed Brown   /* Create timestepping solver context */
4409566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &user.ts));
4419566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(user.ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
442c4762a1bSJed Brown   if (user.implicitform) {
4439566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(user.ts, NULL, IFunction, &user));
4449566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(user.ts, user.A, user.A, IJacobian, &user));
4459566063dSJacob Faibussowitsch     PetscCall(TSSetType(user.ts, TSCN));
446c4762a1bSJed Brown   } else {
4479566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(user.ts, NULL, RHSFunction, &user));
4489566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(user.ts, user.A, user.A, RHSJacobian, &user));
4499566063dSJacob Faibussowitsch     PetscCall(TSSetType(user.ts, TSRK));
450c4762a1bSJed Brown   }
4519566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(user.ts, user.ftime));
4529566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(user.ts, TS_EXACTFINALTIME_MATCHSTEP));
453c4762a1bSJed Brown 
45448a46eb9SPierre Jolivet   if (monitor) PetscCall(TSMonitorSet(user.ts, Monitor, &user, NULL));
455c4762a1bSJed Brown 
456c4762a1bSJed Brown   /* Set ODE initial conditions */
4579566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user.U, &x_ptr));
458c4762a1bSJed Brown   x_ptr[0] = 2.0;
459c4762a1bSJed Brown   x_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * user.mu) - 292.0 / (2187.0 * user.mu * user.mu);
4609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user.U, &x_ptr));
461c4762a1bSJed Brown 
462c4762a1bSJed Brown   /* Set runtime options */
4639566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(user.ts));
464c4762a1bSJed Brown 
465c4762a1bSJed Brown   /* Obtain the observation */
4669566063dSJacob Faibussowitsch   PetscCall(TSSolve(user.ts, user.U));
4679566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user.U, &x_ptr));
468c4762a1bSJed Brown   user.ob[0] = x_ptr[0];
469c4762a1bSJed Brown   user.ob[1] = x_ptr[1];
4709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user.U, &x_ptr));
471c4762a1bSJed Brown 
4729566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user.U, &x));
4739566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &x_ptr));
474c4762a1bSJed Brown   x_ptr[0] = ic1;
475c4762a1bSJed Brown   x_ptr[1] = ic2;
4769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &x_ptr));
477c4762a1bSJed Brown 
478c4762a1bSJed Brown   /* Save trajectory of solution so that TSAdjointSolve() may be used */
4799566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(user.ts));
480c4762a1bSJed Brown 
481c4762a1bSJed Brown   /* Compare finite difference and second-order adjoint. */
482c4762a1bSJed Brown   switch (mode) {
483c4762a1bSJed Brown   case 2:
4849566063dSJacob Faibussowitsch     PetscCall(FiniteDiff(x, arr, &user));
4859566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n Finite difference approximation of the Hessian\n"));
4869566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%g %g\n%g %g\n", (double)arr[0], (double)arr[1], (double)arr[2], (double)arr[3]));
487c4762a1bSJed Brown     break;
488c4762a1bSJed Brown   case 1: /* Compute the Hessian column by column */
4899566063dSJacob Faibussowitsch     PetscCall(VecCopy(x, user.U));
4909566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user.Dir, &x_ptr));
491c4762a1bSJed Brown     x_ptr[0] = 1.;
492c4762a1bSJed Brown     x_ptr[1] = 0.;
4939566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user.Dir, &x_ptr));
4949566063dSJacob Faibussowitsch     PetscCall(Adjoint2(user.U, arr, &user));
4959566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nFirst column of the Hessian\n"));
4969566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%g\n%g\n", (double)arr[0], (double)arr[1]));
4979566063dSJacob Faibussowitsch     PetscCall(VecCopy(x, user.U));
4989566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user.Dir, &x_ptr));
499c4762a1bSJed Brown     x_ptr[0] = 0.;
500c4762a1bSJed Brown     x_ptr[1] = 1.;
5019566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user.Dir, &x_ptr));
5029566063dSJacob Faibussowitsch     PetscCall(Adjoint2(user.U, arr, &user));
5039566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nSecond column of the Hessian\n"));
5049566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%g\n%g\n", (double)arr[0], (double)arr[1]));
505c4762a1bSJed Brown     break;
506c4762a1bSJed Brown   case 0:
507c4762a1bSJed Brown     /* Create optimization context and set up */
5089566063dSJacob Faibussowitsch     PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
5099566063dSJacob Faibussowitsch     PetscCall(TaoSetType(tao, TAOBLMVM));
5109566063dSJacob Faibussowitsch     PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
511c4762a1bSJed Brown 
512c4762a1bSJed Brown     if (mf) {
5139566063dSJacob Faibussowitsch       PetscCall(MatCreateShell(PETSC_COMM_SELF, 2, 2, 2, 2, (void *)&user, &user.H));
51457d50842SBarry Smith       PetscCall(MatShellSetOperation(user.H, MATOP_MULT, (PetscErrorCodeFn *)HessianProductMat));
5159566063dSJacob Faibussowitsch       PetscCall(MatSetOption(user.H, MAT_SYMMETRIC, PETSC_TRUE));
5169566063dSJacob Faibussowitsch       PetscCall(TaoSetHessian(tao, user.H, user.H, MatrixFreeHessian, (void *)&user));
517c4762a1bSJed Brown     } else { /* Create Hessian matrix */
5189566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_WORLD, &user.H));
5199566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(user.H, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
5209566063dSJacob Faibussowitsch       PetscCall(MatSetUp(user.H));
5219566063dSJacob Faibussowitsch       PetscCall(TaoSetHessian(tao, user.H, user.H, FormHessian, (void *)&user));
522c4762a1bSJed Brown     }
523c4762a1bSJed Brown 
524c4762a1bSJed Brown     /* Not use any preconditioner */
5259566063dSJacob Faibussowitsch     PetscCall(TaoGetKSP(tao, &ksp));
526c4762a1bSJed Brown     if (ksp) {
5279566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(ksp, &pc));
5289566063dSJacob Faibussowitsch       PetscCall(PCSetType(pc, PCNONE));
529c4762a1bSJed Brown     }
530c4762a1bSJed Brown 
531c4762a1bSJed Brown     /* Set initial solution guess */
5329566063dSJacob Faibussowitsch     PetscCall(TaoSetSolution(tao, x));
5339566063dSJacob Faibussowitsch     PetscCall(TaoSetFromOptions(tao));
5349566063dSJacob Faibussowitsch     PetscCall(TaoSolve(tao));
5359566063dSJacob Faibussowitsch     PetscCall(TaoDestroy(&tao));
5369566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&user.H));
537c4762a1bSJed Brown     break;
538d71ae5a4SJacob Faibussowitsch   default:
539d71ae5a4SJacob Faibussowitsch     break;
540c4762a1bSJed Brown   }
541c4762a1bSJed Brown 
542c4762a1bSJed Brown   /* Free work space.  All PETSc objects should be destroyed when they are no longer needed. */
5439566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.A));
5449566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.U));
5459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Lambda[0]));
5469566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Lambda2[0]));
5479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Ihp1[0]));
5489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Dir));
5499566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&user.ts));
5509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
5519566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
552b122ec5aSJacob Faibussowitsch   return 0;
553c4762a1bSJed Brown }
554c4762a1bSJed Brown 
555c4762a1bSJed Brown /*TEST
556c4762a1bSJed Brown     build:
557c4762a1bSJed Brown       requires: !complex !single
558c4762a1bSJed Brown 
559c4762a1bSJed Brown     test:
560188af4bfSBarry Smith       args: -ts_type cn -viewer_binary_skip_info -tao_monitor -tao_view -mu 1000 -ts_time_step 0.03125
561c4762a1bSJed Brown       output_file: output/ex20opt_ic_1.out
562c4762a1bSJed Brown 
563c4762a1bSJed Brown     test:
564c4762a1bSJed Brown       suffix: 2
565188af4bfSBarry Smith       args: -ts_type beuler -viewer_binary_skip_info -tao_monitor -tao_view -mu 100 -ts_time_step 0.01 -tao_type bntr -tao_bnk_pc_type none
566c4762a1bSJed Brown       output_file: output/ex20opt_ic_2.out
567c4762a1bSJed Brown 
568c4762a1bSJed Brown     test:
569c4762a1bSJed Brown       suffix: 3
570188af4bfSBarry Smith       args: -ts_type cn -viewer_binary_skip_info -tao_monitor -tao_view -mu 100 -ts_time_step 0.01 -tao_type bntr -tao_bnk_pc_type none
571c4762a1bSJed Brown       output_file: output/ex20opt_ic_3.out
572c4762a1bSJed Brown 
573c4762a1bSJed Brown     test:
574c4762a1bSJed Brown       suffix: 4
575188af4bfSBarry Smith       args: -implicitform 0 -ts_time_step 0.01 -ts_max_steps 2 -ts_rhs_jacobian_test_mult_transpose -mat_shell_test_mult_transpose_view -ts_rhs_jacobian_test_mult -mat_shell_test_mult_view -mode 1 -my_tao_mf
576c4762a1bSJed Brown TEST*/
577