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
127addb90fSBarry Smith -aijpc - set the matrix used to compute the preconditioner 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
InitializeLambda(DM da,Vec lambda,PetscReal x,PetscReal y)20d71ae5a4SJacob Faibussowitsch PetscErrorCode InitializeLambda(DM da, Vec lambda, PetscReal x, PetscReal y)
21d71ae5a4SJacob Faibussowitsch {
22c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym;
23c4762a1bSJed Brown Field **l;
24c4762a1bSJed Brown
254d86920dSPierre Jolivet PetscFunctionBegin;
269566063dSJacob 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));
27c4762a1bSJed Brown /* locate the global i index for x and j index for y */
28c4762a1bSJed Brown i = (PetscInt)(x * (Mx - 1));
29c4762a1bSJed Brown j = (PetscInt)(y * (My - 1));
309566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
31c4762a1bSJed Brown
32c4762a1bSJed Brown if (xs <= i && i < xs + xm && ys <= j && j < ys + ym) {
33c4762a1bSJed Brown /* the i,j vertex is on this process */
349566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, lambda, &l));
35c4762a1bSJed Brown l[j][i].u = 1.0;
36c4762a1bSJed Brown l[j][i].v = 1.0;
379566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, lambda, &l));
38c4762a1bSJed Brown }
393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
40c4762a1bSJed Brown }
41c4762a1bSJed Brown
main(int argc,char ** argv)42d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
43d71ae5a4SJacob Faibussowitsch {
44c4762a1bSJed Brown TS ts; /* ODE integrator */
45c4762a1bSJed Brown Vec x; /* solution */
46c4762a1bSJed Brown DM da;
47c4762a1bSJed Brown AppCtx appctx;
48c4762a1bSJed Brown Vec lambda[1];
49c4762a1bSJed Brown PetscBool forwardonly = PETSC_FALSE, implicitform = PETSC_TRUE;
50c4762a1bSJed Brown
51327415f7SBarry Smith PetscFunctionBeginUser;
52c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-forwardonly", &forwardonly, NULL));
549566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-implicitform", &implicitform, NULL));
55c4762a1bSJed Brown appctx.aijpc = PETSC_FALSE;
569566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-aijpc", &appctx.aijpc, NULL));
57c4762a1bSJed Brown
58c4762a1bSJed Brown appctx.D1 = 8.0e-5;
59c4762a1bSJed Brown appctx.D2 = 4.0e-5;
60c4762a1bSJed Brown appctx.gamma = .024;
61c4762a1bSJed Brown appctx.kappa = .06;
62c4762a1bSJed Brown
63c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors
65c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
669566063dSJacob 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));
679566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da));
689566063dSJacob Faibussowitsch PetscCall(DMSetUp(da));
699566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "u"));
709566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 1, "v"));
71c4762a1bSJed Brown
72c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining
74c4762a1bSJed Brown vectors that are the same types
75c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
769566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x));
77c4762a1bSJed Brown
78c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79c4762a1bSJed Brown Create timestepping solver context
80c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
819566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
829566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da));
839566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
849566063dSJacob Faibussowitsch PetscCall(TSSetEquationType(ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
85c4762a1bSJed Brown if (!implicitform) {
869566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSRK));
879566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &appctx));
889566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, NULL, NULL, RHSJacobian, &appctx));
89c4762a1bSJed Brown } else {
909566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN));
919566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, IFunction, &appctx));
92c4762a1bSJed Brown if (appctx.aijpc) {
93c4762a1bSJed Brown Mat A, B;
94c4762a1bSJed Brown
959566063dSJacob Faibussowitsch PetscCall(DMSetMatType(da, MATSELL));
969566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da, &A));
979566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &B));
98c4762a1bSJed Brown /* FIXME do we need to change viewer to display matrix in natural ordering as DMCreateMatrix_DA does? */
999566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, B, IJacobian, &appctx));
1009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
1019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B));
102c4762a1bSJed Brown } else {
1039566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, NULL, NULL, IJacobian, &appctx));
104c4762a1bSJed Brown }
105c4762a1bSJed Brown }
106c4762a1bSJed Brown
107c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108c4762a1bSJed Brown Set initial conditions
109c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1109566063dSJacob Faibussowitsch PetscCall(InitialConditions(da, x));
1119566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x));
112c4762a1bSJed Brown
113c4762a1bSJed Brown /*
114c4762a1bSJed Brown Have the TS save its trajectory so that TSAdjointSolve() may be used
115c4762a1bSJed Brown */
1169566063dSJacob Faibussowitsch if (!forwardonly) PetscCall(TSSetSaveTrajectory(ts));
117c4762a1bSJed Brown
118c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119c4762a1bSJed Brown Set solver options
120c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1219566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 200.0));
1229566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.5));
1239566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
1249566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
125c4762a1bSJed Brown
126c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127c4762a1bSJed Brown Solve ODE system
128c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1299566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x));
130c4762a1bSJed Brown if (!forwardonly) {
131c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132c4762a1bSJed Brown Start the Adjoint model
133c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &lambda[0]));
135c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */
1369566063dSJacob Faibussowitsch PetscCall(InitializeLambda(da, lambda[0], 0.5, 0.5));
1379566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 1, lambda, NULL));
1389566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts));
1399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0]));
140c4762a1bSJed Brown }
141c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they
143c4762a1bSJed Brown are no longer needed.
144c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1469566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
1479566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da));
1489566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
149b122ec5aSJacob Faibussowitsch return 0;
150c4762a1bSJed Brown }
151c4762a1bSJed Brown
152c4762a1bSJed Brown /* ------------------------------------------------------------------- */
InitialConditions(DM da,Vec U)153d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(DM da, Vec U)
154d71ae5a4SJacob Faibussowitsch {
155c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My;
156c4762a1bSJed Brown Field **u;
157c4762a1bSJed Brown PetscReal hx, hy, x, y;
158c4762a1bSJed Brown
159c4762a1bSJed Brown PetscFunctionBegin;
1609566063dSJacob 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));
161c4762a1bSJed Brown
162c4762a1bSJed Brown hx = 2.5 / (PetscReal)Mx;
163c4762a1bSJed Brown hy = 2.5 / (PetscReal)My;
164c4762a1bSJed Brown
165c4762a1bSJed Brown /*
166c4762a1bSJed Brown Get pointers to vector data
167c4762a1bSJed Brown */
1689566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u));
169c4762a1bSJed Brown
170c4762a1bSJed Brown /*
171c4762a1bSJed Brown Get local grid boundaries
172c4762a1bSJed Brown */
1739566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
174c4762a1bSJed Brown
175c4762a1bSJed Brown /*
176c4762a1bSJed Brown Compute function over the locally owned part of the grid
177c4762a1bSJed Brown */
178c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
179c4762a1bSJed Brown y = j * hy;
180c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
181c4762a1bSJed Brown x = i * hx;
1829371c9d4SSatish Balay if (PetscApproximateGTE(x, 1.0) && PetscApproximateLTE(x, 1.5) && PetscApproximateGTE(y, 1.0) && PetscApproximateLTE(y, 1.5))
1839371c9d4SSatish 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;
184c4762a1bSJed Brown else u[j][i].v = 0.0;
185c4762a1bSJed Brown
186c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0 * u[j][i].v;
187c4762a1bSJed Brown }
188c4762a1bSJed Brown }
189c4762a1bSJed Brown
190c4762a1bSJed Brown /*
191c4762a1bSJed Brown Restore vectors
192c4762a1bSJed Brown */
1939566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u));
1943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
195c4762a1bSJed Brown }
196c4762a1bSJed Brown
197c4762a1bSJed Brown /*TEST
198c4762a1bSJed Brown
199c4762a1bSJed Brown build:
20060f0b76eSHong Zhang depends: reaction_diffusion.c
201c4762a1bSJed Brown requires: !complex !single
202c4762a1bSJed Brown
203c4762a1bSJed Brown test:
204c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -da_grid_x 20 -da_grid_y 20
205c4762a1bSJed Brown output_file: output/ex5adj_1.out
206c4762a1bSJed Brown
207c4762a1bSJed Brown test:
208c4762a1bSJed Brown suffix: 2
209c4762a1bSJed Brown nsize: 2
210*188af4bfSBarry Smith args: -ts_max_steps 10 -ts_time_step 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
211c4762a1bSJed Brown
212c4762a1bSJed Brown test:
213c4762a1bSJed Brown suffix: 3
214c4762a1bSJed Brown nsize: 2
215*188af4bfSBarry Smith args: -ts_max_steps 10 -ts_time_step 10 -ts_adjoint_monitor_draw_sensi
2163886731fSPierre Jolivet output_file: output/empty.out
217c4762a1bSJed Brown
218c4762a1bSJed Brown test:
219c4762a1bSJed Brown suffix: 4
220c4762a1bSJed Brown nsize: 2
221*188af4bfSBarry Smith args: -ts_max_steps 10 -ts_time_step 10 -ts_monitor -ts_adjoint_monitor -ksp_monitor_short -da_grid_x 20 -da_grid_y 20 -snes_fd_color
222c4762a1bSJed Brown output_file: output/ex5adj_2.out
223c4762a1bSJed Brown
224c4762a1bSJed Brown test:
225c4762a1bSJed Brown suffix: 5
226c4762a1bSJed Brown nsize: 2
227c4762a1bSJed 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
228c4762a1bSJed Brown output_file: output/ex5adj_1.out
229c4762a1bSJed Brown
230c4762a1bSJed Brown test:
231c4762a1bSJed Brown suffix: knl
232c4762a1bSJed 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
2333886731fSPierre Jolivet output_file: output/empty.out
234c4762a1bSJed Brown requires: knl
235c4762a1bSJed Brown
236c4762a1bSJed Brown test:
237c4762a1bSJed Brown suffix: sell
238c4762a1bSJed Brown nsize: 4
239c4762a1bSJed Brown args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type sell -pc_type none
240c4762a1bSJed Brown output_file: output/ex5adj_sell_1.out
241c4762a1bSJed Brown
242c4762a1bSJed Brown test:
243c4762a1bSJed Brown suffix: aijsell
244c4762a1bSJed Brown nsize: 4
245c4762a1bSJed Brown args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type aijsell -pc_type none
246c4762a1bSJed Brown output_file: output/ex5adj_sell_1.out
247c4762a1bSJed Brown
248c4762a1bSJed Brown test:
249c4762a1bSJed Brown suffix: sell2
250c4762a1bSJed Brown nsize: 4
251c4762a1bSJed 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
252c4762a1bSJed Brown output_file: output/ex5adj_sell_2.out
253c4762a1bSJed Brown
254c4762a1bSJed Brown test:
255c4762a1bSJed Brown suffix: aijsell2
256c4762a1bSJed Brown nsize: 4
257c4762a1bSJed 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
258c4762a1bSJed Brown output_file: output/ex5adj_sell_2.out
259c4762a1bSJed Brown
260c4762a1bSJed Brown test:
261c4762a1bSJed Brown suffix: sell3
262c4762a1bSJed Brown nsize: 4
263c4762a1bSJed 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
264c4762a1bSJed Brown output_file: output/ex5adj_sell_3.out
265c4762a1bSJed Brown
266c4762a1bSJed Brown test:
267c4762a1bSJed Brown suffix: sell4
268c4762a1bSJed Brown nsize: 4
269c4762a1bSJed 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
270c4762a1bSJed Brown output_file: output/ex5adj_sell_4.out
271c4762a1bSJed Brown
272c4762a1bSJed Brown test:
273c4762a1bSJed Brown suffix: sell5
274c4762a1bSJed Brown nsize: 4
275c4762a1bSJed Brown args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type sell -aijpc
276c4762a1bSJed Brown output_file: output/ex5adj_sell_5.out
277c4762a1bSJed Brown
278c4762a1bSJed Brown test:
279c4762a1bSJed Brown suffix: aijsell5
280c4762a1bSJed Brown nsize: 4
281c4762a1bSJed Brown args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type aijsell
282c4762a1bSJed Brown output_file: output/ex5adj_sell_5.out
283c4762a1bSJed Brown
284c4762a1bSJed Brown test:
285c4762a1bSJed Brown suffix: sell6
286c4762a1bSJed 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
287c4762a1bSJed Brown output_file: output/ex5adj_sell_6.out
288c4762a1bSJed Brown
289fcb023d4SJed Brown test:
2905f9962eeSHong Zhang suffix: sell7
2915f9962eeSHong Zhang args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -ts_trajectory_type memory -ts_trajectory_solution_only 0 -dm_mat_type sellcuda -dm_vec_type cuda -pc_type jacobi
2925f9962eeSHong Zhang output_file: output/ex5adj_sell_6.out
2935f9962eeSHong Zhang requires: cuda
2945f9962eeSHong Zhang
2955f9962eeSHong Zhang test:
296fcb023d4SJed Brown suffix: gamg1
2975fb461b6SMatthew Knepley args: -pc_type gamg -pc_gamg_esteig_ksp_type gmres -pc_gamg_esteig_ksp_max_it 10 -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
29873f7197eSJed Brown output_file: output/ex5adj_gamg.out
299fcb023d4SJed Brown
30075f8e9bdSHong Zhang test:
30175f8e9bdSHong Zhang suffix: gamg2
3025fb461b6SMatthew Knepley args: -pc_type gamg -pc_gamg_esteig_ksp_type gmres -pc_gamg_esteig_ksp_max_it 10 -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
30373f7197eSJed Brown output_file: output/ex5adj_gamg.out
304c4762a1bSJed Brown TEST*/
305