xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex5adj.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "Demonstrates adjoint sensitivity analysis for Reaction-Diffusion Equations.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown   See ex5.c for details on the equation.
5c4762a1bSJed Brown   This code demonestrates the TSAdjoint interface to a system of time-dependent partial differential equations.
6c4762a1bSJed Brown   It computes the sensitivity of a component in the final solution, which locates in the center of the 2D domain, w.r.t. the initial conditions.
7c4762a1bSJed Brown   The user does not need to provide any additional functions. The required functions in the original simulation are reused in the adjoint run.
8c4762a1bSJed Brown 
9c4762a1bSJed Brown   Runtime options:
10c4762a1bSJed Brown     -forwardonly  - run the forward simulation without adjoint
11c4762a1bSJed Brown     -implicitform - provide IFunction and IJacobian to TS, if not set, RHSFunction and RHSJacobian will be used
12c4762a1bSJed Brown     -aijpc        - set the preconditioner matrix to be aij (the Jacobian matrix can be of a different type such as ELL)
13c4762a1bSJed Brown */
1460f0b76eSHong Zhang #include "reaction_diffusion.h"
15c4762a1bSJed Brown #include <petscdm.h>
16c4762a1bSJed Brown #include <petscdmda.h>
17c4762a1bSJed Brown 
18c4762a1bSJed Brown PetscErrorCode InitialConditions(DM, Vec);
19c4762a1bSJed Brown 
20*9371c9d4SSatish Balay PetscErrorCode InitializeLambda(DM da, Vec lambda, PetscReal x, PetscReal y) {
21c4762a1bSJed Brown   PetscInt i, j, Mx, My, xs, ys, xm, ym;
22c4762a1bSJed Brown   Field  **l;
23c4762a1bSJed Brown   PetscFunctionBegin;
24c4762a1bSJed Brown 
259566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
26c4762a1bSJed Brown   /* locate the global i index for x and j index for y */
27c4762a1bSJed Brown   i = (PetscInt)(x * (Mx - 1));
28c4762a1bSJed Brown   j = (PetscInt)(y * (My - 1));
299566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   if (xs <= i && i < xs + xm && ys <= j && j < ys + ym) {
32c4762a1bSJed Brown     /* the i,j vertex is on this process */
339566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(da, lambda, &l));
34c4762a1bSJed Brown     l[j][i].u = 1.0;
35c4762a1bSJed Brown     l[j][i].v = 1.0;
369566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(da, lambda, &l));
37c4762a1bSJed Brown   }
38c4762a1bSJed Brown   PetscFunctionReturn(0);
39c4762a1bSJed Brown }
40c4762a1bSJed Brown 
41*9371c9d4SSatish Balay int main(int argc, char **argv) {
42c4762a1bSJed Brown   TS        ts; /* ODE integrator */
43c4762a1bSJed Brown   Vec       x;  /* solution */
44c4762a1bSJed Brown   DM        da;
45c4762a1bSJed Brown   AppCtx    appctx;
46c4762a1bSJed Brown   Vec       lambda[1];
47c4762a1bSJed Brown   PetscBool forwardonly = PETSC_FALSE, implicitform = PETSC_TRUE;
48c4762a1bSJed Brown 
49327415f7SBarry Smith   PetscFunctionBeginUser;
509566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-forwardonly", &forwardonly, NULL));
529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-implicitform", &implicitform, NULL));
53c4762a1bSJed Brown   appctx.aijpc = PETSC_FALSE;
549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-aijpc", &appctx.aijpc, NULL));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   appctx.D1    = 8.0e-5;
57c4762a1bSJed Brown   appctx.D2    = 4.0e-5;
58c4762a1bSJed Brown   appctx.gamma = .024;
59c4762a1bSJed Brown   appctx.kappa = .06;
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
63c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
649566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, DM_BOUNDARY_PERIODIC, DMDA_STENCIL_STAR, 64, 64, PETSC_DECIDE, PETSC_DECIDE, 2, 1, NULL, NULL, &da));
659566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
669566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
679566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 0, "u"));
689566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 1, "v"));
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
72c4762a1bSJed Brown      vectors that are the same types
73c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
749566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &x));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77c4762a1bSJed Brown      Create timestepping solver context
78c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
799566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
809566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
819566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
829566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
83c4762a1bSJed Brown   if (!implicitform) {
849566063dSJacob Faibussowitsch     PetscCall(TSSetType(ts, TSRK));
859566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &appctx));
869566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ts, NULL, NULL, RHSJacobian, &appctx));
87c4762a1bSJed Brown   } else {
889566063dSJacob Faibussowitsch     PetscCall(TSSetType(ts, TSCN));
899566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, IFunction, &appctx));
90c4762a1bSJed Brown     if (appctx.aijpc) {
91c4762a1bSJed Brown       Mat A, B;
92c4762a1bSJed Brown 
939566063dSJacob Faibussowitsch       PetscCall(DMSetMatType(da, MATSELL));
949566063dSJacob Faibussowitsch       PetscCall(DMCreateMatrix(da, &A));
959566063dSJacob Faibussowitsch       PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &B));
96c4762a1bSJed Brown       /* FIXME do we need to change viewer to display matrix in natural ordering as DMCreateMatrix_DA does? */
979566063dSJacob Faibussowitsch       PetscCall(TSSetIJacobian(ts, A, B, IJacobian, &appctx));
989566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A));
999566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&B));
100c4762a1bSJed Brown     } else {
1019566063dSJacob Faibussowitsch       PetscCall(TSSetIJacobian(ts, NULL, NULL, IJacobian, &appctx));
102c4762a1bSJed Brown     }
103c4762a1bSJed Brown   }
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106c4762a1bSJed Brown      Set initial conditions
107c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1089566063dSJacob Faibussowitsch   PetscCall(InitialConditions(da, x));
1099566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, x));
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /*
112c4762a1bSJed Brown     Have the TS save its trajectory so that TSAdjointSolve() may be used
113c4762a1bSJed Brown   */
1149566063dSJacob Faibussowitsch   if (!forwardonly) PetscCall(TSSetSaveTrajectory(ts));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117c4762a1bSJed Brown      Set solver options
118c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1199566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 200.0));
1209566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.5));
1219566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
1229566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125c4762a1bSJed Brown      Solve ODE system
126c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1279566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
128c4762a1bSJed Brown   if (!forwardonly) {
129c4762a1bSJed Brown     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130c4762a1bSJed Brown        Start the Adjoint model
131c4762a1bSJed Brown        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1329566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &lambda[0]));
133c4762a1bSJed Brown     /*   Reset initial conditions for the adjoint integration */
1349566063dSJacob Faibussowitsch     PetscCall(InitializeLambda(da, lambda[0], 0.5, 0.5));
1359566063dSJacob Faibussowitsch     PetscCall(TSSetCostGradients(ts, 1, lambda, NULL));
1369566063dSJacob Faibussowitsch     PetscCall(TSAdjointSolve(ts));
1379566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lambda[0]));
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
141c4762a1bSJed Brown      are no longer needed.
142c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1439566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1449566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1459566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1469566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
147b122ec5aSJacob Faibussowitsch   return 0;
148c4762a1bSJed Brown }
149c4762a1bSJed Brown 
150c4762a1bSJed Brown /* ------------------------------------------------------------------- */
151*9371c9d4SSatish Balay PetscErrorCode InitialConditions(DM da, Vec U) {
152c4762a1bSJed Brown   PetscInt  i, j, xs, ys, xm, ym, Mx, My;
153c4762a1bSJed Brown   Field   **u;
154c4762a1bSJed Brown   PetscReal hx, hy, x, y;
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   PetscFunctionBegin;
1579566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   hx = 2.5 / (PetscReal)Mx;
160c4762a1bSJed Brown   hy = 2.5 / (PetscReal)My;
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   /*
163c4762a1bSJed Brown      Get pointers to vector data
164c4762a1bSJed Brown   */
1659566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, U, &u));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /*
168c4762a1bSJed Brown      Get local grid boundaries
169c4762a1bSJed Brown   */
1709566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /*
173c4762a1bSJed Brown      Compute function over the locally owned part of the grid
174c4762a1bSJed Brown   */
175c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
176c4762a1bSJed Brown     y = j * hy;
177c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
178c4762a1bSJed Brown       x = i * hx;
179*9371c9d4SSatish Balay       if (PetscApproximateGTE(x, 1.0) && PetscApproximateLTE(x, 1.5) && PetscApproximateGTE(y, 1.0) && PetscApproximateLTE(y, 1.5))
180*9371c9d4SSatish Balay         u[j][i].v = PetscPowReal(PetscSinReal(4.0 * PETSC_PI * x), 2.0) * PetscPowReal(PetscSinReal(4.0 * PETSC_PI * y), 2.0) / 4.0;
181c4762a1bSJed Brown       else u[j][i].v = 0.0;
182c4762a1bSJed Brown 
183c4762a1bSJed Brown       u[j][i].u = 1.0 - 2.0 * u[j][i].v;
184c4762a1bSJed Brown     }
185c4762a1bSJed Brown   }
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   /*
188c4762a1bSJed Brown      Restore vectors
189c4762a1bSJed Brown   */
1909566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, U, &u));
191c4762a1bSJed Brown   PetscFunctionReturn(0);
192c4762a1bSJed Brown }
193c4762a1bSJed Brown 
194c4762a1bSJed Brown /*TEST
195c4762a1bSJed Brown 
196c4762a1bSJed Brown    build:
19760f0b76eSHong Zhang       depends: reaction_diffusion.c
198c4762a1bSJed Brown       requires: !complex !single
199c4762a1bSJed Brown 
200c4762a1bSJed Brown    test:
201c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -da_grid_x 20 -da_grid_y 20
202c4762a1bSJed Brown       output_file: output/ex5adj_1.out
203c4762a1bSJed Brown 
204c4762a1bSJed Brown    test:
205c4762a1bSJed Brown       suffix: 2
206c4762a1bSJed Brown       nsize: 2
207c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_dt 10 -ts_monitor -ts_adjoint_monitor -ksp_monitor_short -da_grid_x 20 -da_grid_y 20 -ts_trajectory_dirname Test-dir -ts_trajectory_file_template test-%06D.cp
208c4762a1bSJed Brown 
209c4762a1bSJed Brown    test:
210c4762a1bSJed Brown       suffix: 3
211c4762a1bSJed Brown       nsize: 2
212c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_dt 10 -ts_adjoint_monitor_draw_sensi
213c4762a1bSJed Brown 
214c4762a1bSJed Brown    test:
215c4762a1bSJed Brown       suffix: 4
216c4762a1bSJed Brown       nsize: 2
217c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_dt 10 -ts_monitor -ts_adjoint_monitor -ksp_monitor_short -da_grid_x 20 -da_grid_y 20 -snes_fd_color
218c4762a1bSJed Brown       output_file: output/ex5adj_2.out
219c4762a1bSJed Brown 
220c4762a1bSJed Brown    test:
221c4762a1bSJed Brown       suffix: 5
222c4762a1bSJed Brown       nsize: 2
223c4762a1bSJed Brown       args: -ts_max_steps 10 -implicitform 0 -ts_type rk -ts_rk_type 4 -ts_monitor -ts_adjoint_monitor -da_grid_x 20 -da_grid_y 20 -snes_fd_color
224c4762a1bSJed Brown       output_file: output/ex5adj_1.out
225c4762a1bSJed Brown 
226c4762a1bSJed Brown    test:
227c4762a1bSJed Brown       suffix: knl
228c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -ts_trajectory_type memory -ts_trajectory_solution_only 0 -malloc_hbw -ts_trajectory_use_dram 1
229c4762a1bSJed Brown       output_file: output/ex5adj_3.out
230c4762a1bSJed Brown       requires: knl
231c4762a1bSJed Brown 
232c4762a1bSJed Brown    test:
233c4762a1bSJed Brown       suffix: sell
234c4762a1bSJed Brown       nsize: 4
235c4762a1bSJed Brown       args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type sell -pc_type none
236c4762a1bSJed Brown       output_file: output/ex5adj_sell_1.out
237c4762a1bSJed Brown 
238c4762a1bSJed Brown    test:
239c4762a1bSJed Brown       suffix: aijsell
240c4762a1bSJed Brown       nsize: 4
241c4762a1bSJed Brown       args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type aijsell -pc_type none
242c4762a1bSJed Brown       output_file: output/ex5adj_sell_1.out
243c4762a1bSJed Brown 
244c4762a1bSJed Brown    test:
245c4762a1bSJed Brown       suffix: sell2
246c4762a1bSJed Brown       nsize: 4
247c4762a1bSJed Brown       args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type sell -pc_type mg -pc_mg_levels 2 -mg_coarse_pc_type sor
248c4762a1bSJed Brown       output_file: output/ex5adj_sell_2.out
249c4762a1bSJed Brown 
250c4762a1bSJed Brown    test:
251c4762a1bSJed Brown       suffix: aijsell2
252c4762a1bSJed Brown       nsize: 4
253c4762a1bSJed Brown       args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type aijsell -pc_type mg -pc_mg_levels 2 -mg_coarse_pc_type sor
254c4762a1bSJed Brown       output_file: output/ex5adj_sell_2.out
255c4762a1bSJed Brown 
256c4762a1bSJed Brown    test:
257c4762a1bSJed Brown       suffix: sell3
258c4762a1bSJed Brown       nsize: 4
259c4762a1bSJed Brown       args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type sell -pc_type mg -pc_mg_levels 2 -mg_coarse_pc_type bjacobi -mg_levels_pc_type bjacobi
260c4762a1bSJed Brown       output_file: output/ex5adj_sell_3.out
261c4762a1bSJed Brown 
262c4762a1bSJed Brown    test:
263c4762a1bSJed Brown       suffix: sell4
264c4762a1bSJed Brown       nsize: 4
265c4762a1bSJed Brown       args: -forwardonly -implicitform -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type sell -pc_type mg -pc_mg_levels 2 -mg_coarse_pc_type bjacobi -mg_levels_pc_type bjacobi
266c4762a1bSJed Brown       output_file: output/ex5adj_sell_4.out
267c4762a1bSJed Brown 
268c4762a1bSJed Brown    test:
269c4762a1bSJed Brown       suffix: sell5
270c4762a1bSJed Brown       nsize: 4
271c4762a1bSJed Brown       args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type sell -aijpc
272c4762a1bSJed Brown       output_file: output/ex5adj_sell_5.out
273c4762a1bSJed Brown 
274c4762a1bSJed Brown    test:
275c4762a1bSJed Brown       suffix: aijsell5
276c4762a1bSJed Brown       nsize: 4
277c4762a1bSJed Brown       args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type aijsell
278c4762a1bSJed Brown       output_file: output/ex5adj_sell_5.out
279c4762a1bSJed Brown 
280c4762a1bSJed Brown    test:
281c4762a1bSJed Brown       suffix: sell6
282c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -ts_trajectory_type memory -ts_trajectory_solution_only 0 -dm_mat_type sell -pc_type jacobi
283c4762a1bSJed Brown       output_file: output/ex5adj_sell_6.out
284c4762a1bSJed Brown 
285fcb023d4SJed Brown    test:
286fcb023d4SJed Brown       suffix: gamg1
28773f7197eSJed Brown       args: -pc_type gamg -pc_gamg_esteig_ksp_type gmres -pc_gamg_esteig_ksp_max_it 10 -pc_mg_levels 2 -ksp_monitor_short -ts_max_steps 2 -ts_monitor -ts_adjoint_monitor -ts_trajectory_type memory -ksp_rtol 1e-2 -pc_gamg_use_sa_esteig 0
28873f7197eSJed Brown       output_file: output/ex5adj_gamg.out
289fcb023d4SJed Brown 
29075f8e9bdSHong Zhang    test:
29175f8e9bdSHong Zhang       suffix: gamg2
29273f7197eSJed Brown       args: -pc_type gamg -pc_gamg_esteig_ksp_type gmres -pc_gamg_esteig_ksp_max_it 10 -pc_mg_levels 2 -ksp_monitor_short -ts_max_steps 2 -ts_monitor -ts_adjoint_monitor -ts_trajectory_type memory -ksp_use_explicittranspose -ksp_rtol 1e-2 -pc_gamg_use_sa_esteig 0
29373f7197eSJed Brown       output_file: output/ex5adj_gamg.out
294c4762a1bSJed Brown TEST*/
295