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 20c4762a1bSJed Brown PetscErrorCode InitializeLambda(DM da,Vec lambda,PetscReal x,PetscReal y) 21c4762a1bSJed Brown { 22c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 23c4762a1bSJed Brown Field **l; 24c4762a1bSJed Brown PetscFunctionBegin; 25c4762a1bSJed Brown 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,lambda,&l)); 35c4762a1bSJed Brown l[j][i].u = 1.0; 36c4762a1bSJed Brown l[j][i].v = 1.0; 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,lambda,&l)); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown PetscFunctionReturn(0); 40c4762a1bSJed Brown } 41c4762a1bSJed Brown 42c4762a1bSJed Brown int main(int argc,char **argv) 43c4762a1bSJed Brown { 44c4762a1bSJed Brown TS ts; /* ODE integrator */ 45c4762a1bSJed Brown Vec x; /* solution */ 46c4762a1bSJed Brown PetscErrorCode ierr; 47c4762a1bSJed Brown DM da; 48c4762a1bSJed Brown AppCtx appctx; 49c4762a1bSJed Brown Vec lambda[1]; 50c4762a1bSJed Brown PetscBool forwardonly=PETSC_FALSE,implicitform=PETSC_TRUE; 51c4762a1bSJed Brown 52c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-forwardonly",&forwardonly,NULL)); 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL)); 55c4762a1bSJed Brown appctx.aijpc = PETSC_FALSE; 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,64,64,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da)); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetFieldName(da,0,"u")); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&x)); 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 79c4762a1bSJed Brown Create timestepping solver context 80c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts,da)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetEquationType(ts,TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 85c4762a1bSJed Brown if (!implicitform) { 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSRK)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunction,&appctx)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,NULL,NULL,RHSJacobian,&appctx)); 89c4762a1bSJed Brown } else { 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSCN)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(ts,NULL,IFunction,&appctx)); 92c4762a1bSJed Brown if (appctx.aijpc) { 93c4762a1bSJed Brown Mat A,B; 94c4762a1bSJed Brown 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(da,MATSELL)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(da,&A)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(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? */ 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIJacobian(ts,A,B,IJacobian,&appctx)); 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 102c4762a1bSJed Brown } else { 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIJacobian(ts,NULL,NULL,IJacobian,&appctx)); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown } 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 108c4762a1bSJed Brown Set initial conditions 109c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(InitialConditions(da,x)); 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,x)); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* 114c4762a1bSJed Brown Have the TS save its trajectory so that TSAdjointSolve() may be used 115c4762a1bSJed Brown */ 116*5f80ce2aSJacob Faibussowitsch if (!forwardonly) CHKERRQ(TSSetSaveTrajectory(ts)); 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 119c4762a1bSJed Brown Set solver options 120c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,200.0)); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,0.5)); 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 124*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 127c4762a1bSJed Brown Solve ODE system 128c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,x)); 130c4762a1bSJed Brown if (!forwardonly) { 131c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 132c4762a1bSJed Brown Start the Adjoint model 133c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&lambda[0])); 135c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */ 136*5f80ce2aSJacob Faibussowitsch CHKERRQ(InitializeLambda(da,lambda[0],0.5,0.5)); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetCostGradients(ts,1,lambda,NULL)); 138*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSAdjointSolve(ts)); 139*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 147*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 148c4762a1bSJed Brown ierr = PetscFinalize(); 149c4762a1bSJed Brown return ierr; 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 153c4762a1bSJed Brown PetscErrorCode InitialConditions(DM da,Vec U) 154c4762a1bSJed Brown { 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; 160*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 168*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,U,&u)); 169c4762a1bSJed Brown 170c4762a1bSJed Brown /* 171c4762a1bSJed Brown Get local grid boundaries 172c4762a1bSJed Brown */ 173*5f80ce2aSJacob Faibussowitsch CHKERRQ(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; 18266baab88SBarry Smith if (PetscApproximateGTE(x,1.0) && PetscApproximateLTE(x,1.5) && PetscApproximateGTE(y,1.0) && PetscApproximateLTE(y,1.5)) u[j][i].v = PetscPowReal(PetscSinReal(4.0*PETSC_PI*x),2.0)*PetscPowReal(PetscSinReal(4.0*PETSC_PI*y),2.0)/4.0; 183c4762a1bSJed Brown else u[j][i].v = 0.0; 184c4762a1bSJed Brown 185c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0*u[j][i].v; 186c4762a1bSJed Brown } 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* 190c4762a1bSJed Brown Restore vectors 191c4762a1bSJed Brown */ 192*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,U,&u)); 193c4762a1bSJed Brown PetscFunctionReturn(0); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196c4762a1bSJed Brown /*TEST 197c4762a1bSJed Brown 198c4762a1bSJed Brown build: 19960f0b76eSHong Zhang depends: reaction_diffusion.c 200c4762a1bSJed Brown requires: !complex !single 201c4762a1bSJed Brown 202c4762a1bSJed Brown test: 203c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -da_grid_x 20 -da_grid_y 20 204c4762a1bSJed Brown output_file: output/ex5adj_1.out 205c4762a1bSJed Brown 206c4762a1bSJed Brown test: 207c4762a1bSJed Brown suffix: 2 208c4762a1bSJed Brown nsize: 2 209c4762a1bSJed 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 210c4762a1bSJed Brown 211c4762a1bSJed Brown test: 212c4762a1bSJed Brown suffix: 3 213c4762a1bSJed Brown nsize: 2 214c4762a1bSJed Brown args: -ts_max_steps 10 -ts_dt 10 -ts_adjoint_monitor_draw_sensi 215c4762a1bSJed Brown 216c4762a1bSJed Brown test: 217c4762a1bSJed Brown suffix: 4 218c4762a1bSJed Brown nsize: 2 219c4762a1bSJed 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 220c4762a1bSJed Brown output_file: output/ex5adj_2.out 221c4762a1bSJed Brown 222c4762a1bSJed Brown test: 223c4762a1bSJed Brown suffix: 5 224c4762a1bSJed Brown nsize: 2 225c4762a1bSJed 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 226c4762a1bSJed Brown output_file: output/ex5adj_1.out 227c4762a1bSJed Brown 228c4762a1bSJed Brown test: 229c4762a1bSJed Brown suffix: knl 230c4762a1bSJed 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 231c4762a1bSJed Brown output_file: output/ex5adj_3.out 232c4762a1bSJed Brown requires: knl 233c4762a1bSJed Brown 234c4762a1bSJed Brown test: 235c4762a1bSJed Brown suffix: sell 236c4762a1bSJed Brown nsize: 4 237c4762a1bSJed Brown args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type sell -pc_type none 238c4762a1bSJed Brown output_file: output/ex5adj_sell_1.out 239c4762a1bSJed Brown 240c4762a1bSJed Brown test: 241c4762a1bSJed Brown suffix: aijsell 242c4762a1bSJed Brown nsize: 4 243c4762a1bSJed Brown args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type aijsell -pc_type none 244c4762a1bSJed Brown output_file: output/ex5adj_sell_1.out 245c4762a1bSJed Brown 246c4762a1bSJed Brown test: 247c4762a1bSJed Brown suffix: sell2 248c4762a1bSJed Brown nsize: 4 249c4762a1bSJed 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 250c4762a1bSJed Brown output_file: output/ex5adj_sell_2.out 251c4762a1bSJed Brown 252c4762a1bSJed Brown test: 253c4762a1bSJed Brown suffix: aijsell2 254c4762a1bSJed Brown nsize: 4 255c4762a1bSJed 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 256c4762a1bSJed Brown output_file: output/ex5adj_sell_2.out 257c4762a1bSJed Brown 258c4762a1bSJed Brown test: 259c4762a1bSJed Brown suffix: sell3 260c4762a1bSJed Brown nsize: 4 261c4762a1bSJed 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 262c4762a1bSJed Brown output_file: output/ex5adj_sell_3.out 263c4762a1bSJed Brown 264c4762a1bSJed Brown test: 265c4762a1bSJed Brown suffix: sell4 266c4762a1bSJed Brown nsize: 4 267c4762a1bSJed 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 268c4762a1bSJed Brown output_file: output/ex5adj_sell_4.out 269c4762a1bSJed Brown 270c4762a1bSJed Brown test: 271c4762a1bSJed Brown suffix: sell5 272c4762a1bSJed Brown nsize: 4 273c4762a1bSJed Brown args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type sell -aijpc 274c4762a1bSJed Brown output_file: output/ex5adj_sell_5.out 275c4762a1bSJed Brown 276c4762a1bSJed Brown test: 277c4762a1bSJed Brown suffix: aijsell5 278c4762a1bSJed Brown nsize: 4 279c4762a1bSJed Brown args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type aijsell 280c4762a1bSJed Brown output_file: output/ex5adj_sell_5.out 281c4762a1bSJed Brown 282c4762a1bSJed Brown test: 283c4762a1bSJed Brown suffix: sell6 284c4762a1bSJed 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 285c4762a1bSJed Brown output_file: output/ex5adj_sell_6.out 286c4762a1bSJed Brown 287fcb023d4SJed Brown test: 288fcb023d4SJed Brown suffix: gamg1 28973f7197eSJed 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 29073f7197eSJed Brown output_file: output/ex5adj_gamg.out 291fcb023d4SJed Brown 29275f8e9bdSHong Zhang test: 29375f8e9bdSHong Zhang suffix: gamg2 29473f7197eSJed 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 29573f7197eSJed Brown output_file: output/ex5adj_gamg.out 296c4762a1bSJed Brown TEST*/ 297