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 */ 14*60f0b76eSHong 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 PetscErrorCode ierr; 24c4762a1bSJed Brown Field **l; 25c4762a1bSJed Brown PetscFunctionBegin; 26c4762a1bSJed Brown 27c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 28c4762a1bSJed Brown /* locate the global i index for x and j index for y */ 29c4762a1bSJed Brown i = (PetscInt)(x*(Mx-1)); 30c4762a1bSJed Brown j = (PetscInt)(y*(My-1)); 31c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 32c4762a1bSJed Brown 33c4762a1bSJed Brown if (xs <= i && i < xs+xm && ys <= j && j < ys+ym) { 34c4762a1bSJed Brown /* the i,j vertex is on this process */ 35c4762a1bSJed Brown ierr = DMDAVecGetArray(da,lambda,&l);CHKERRQ(ierr); 36c4762a1bSJed Brown l[j][i].u = 1.0; 37c4762a1bSJed Brown l[j][i].v = 1.0; 38c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,lambda,&l);CHKERRQ(ierr); 39c4762a1bSJed Brown } 40c4762a1bSJed Brown PetscFunctionReturn(0); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 43c4762a1bSJed Brown int main(int argc,char **argv) 44c4762a1bSJed Brown { 45c4762a1bSJed Brown TS ts; /* ODE integrator */ 46c4762a1bSJed Brown Vec x; /* solution */ 47c4762a1bSJed Brown PetscErrorCode ierr; 48c4762a1bSJed Brown DM da; 49c4762a1bSJed Brown AppCtx appctx; 50c4762a1bSJed Brown Vec lambda[1]; 51c4762a1bSJed Brown PetscBool forwardonly=PETSC_FALSE,implicitform=PETSC_TRUE; 52c4762a1bSJed Brown 53c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 54c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-forwardonly",&forwardonly,NULL);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL);CHKERRQ(ierr); 56c4762a1bSJed Brown appctx.aijpc = PETSC_FALSE; 57c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-aijpc",&appctx.aijpc,NULL);CHKERRQ(ierr); 58c4762a1bSJed Brown 59c4762a1bSJed Brown appctx.D1 = 8.0e-5; 60c4762a1bSJed Brown appctx.D2 = 4.0e-5; 61c4762a1bSJed Brown appctx.gamma = .024; 62c4762a1bSJed Brown appctx.kappa = .06; 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 65c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 66c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 67c4762a1bSJed Brown ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,64,64,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da);CHKERRQ(ierr); 68c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 69c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 70c4762a1bSJed Brown ierr = DMDASetFieldName(da,0,"u");CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = DMDASetFieldName(da,1,"v");CHKERRQ(ierr); 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 74c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 75c4762a1bSJed Brown vectors that are the same types 76c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 77c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 80c4762a1bSJed Brown Create timestepping solver context 81c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 82c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 83c4762a1bSJed Brown ierr = TSSetDM(ts,da);CHKERRQ(ierr); 84c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 85c4762a1bSJed Brown ierr = TSSetEquationType(ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 86c4762a1bSJed Brown if (!implicitform) { 87c4762a1bSJed Brown ierr = TSSetType(ts,TSRK);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&appctx);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,NULL,NULL,RHSJacobian,&appctx);CHKERRQ(ierr); 90c4762a1bSJed Brown } else { 91c4762a1bSJed Brown ierr = TSSetType(ts,TSCN);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,IFunction,&appctx);CHKERRQ(ierr); 93c4762a1bSJed Brown if (appctx.aijpc) { 94c4762a1bSJed Brown Mat A,B; 95c4762a1bSJed Brown 96c4762a1bSJed Brown ierr = DMSetMatType(da,MATSELL);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = DMCreateMatrix(da,&A);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 99c4762a1bSJed Brown /* FIXME do we need to change viewer to display matrix in natural ordering as DMCreateMatrix_DA does? */ 100c4762a1bSJed Brown ierr = TSSetIJacobian(ts,A,B,IJacobian,&appctx);CHKERRQ(ierr); 101c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 102c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 103c4762a1bSJed Brown } else { 104c4762a1bSJed Brown ierr = TSSetIJacobian(ts,NULL,NULL,IJacobian,&appctx);CHKERRQ(ierr); 105c4762a1bSJed Brown } 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 109c4762a1bSJed Brown Set initial conditions 110c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 111c4762a1bSJed Brown ierr = InitialConditions(da,x);CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = TSSetSolution(ts,x);CHKERRQ(ierr); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* 115c4762a1bSJed Brown Have the TS save its trajectory so that TSAdjointSolve() may be used 116c4762a1bSJed Brown */ 117c4762a1bSJed Brown if (!forwardonly) { ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); } 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 120c4762a1bSJed Brown Set solver options 121c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 122c4762a1bSJed Brown ierr = TSSetMaxTime(ts,200.0);CHKERRQ(ierr); 123c4762a1bSJed Brown ierr = TSSetTimeStep(ts,0.5);CHKERRQ(ierr); 124c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 128c4762a1bSJed Brown Solve ODE system 129c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 130c4762a1bSJed Brown ierr = TSSolve(ts,x);CHKERRQ(ierr); 131c4762a1bSJed Brown if (!forwardonly) { 132c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 133c4762a1bSJed Brown Start the Adjoint model 134c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 135c4762a1bSJed Brown ierr = VecDuplicate(x,&lambda[0]);CHKERRQ(ierr); 136c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */ 137c4762a1bSJed Brown ierr = InitializeLambda(da,lambda[0],0.5,0.5);CHKERRQ(ierr); 138c4762a1bSJed Brown ierr = TSSetCostGradients(ts,1,lambda,NULL);CHKERRQ(ierr); 139c4762a1bSJed Brown ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr); 141c4762a1bSJed Brown } 142c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 143c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 144c4762a1bSJed Brown are no longer needed. 145c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 146c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 148c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = PetscFinalize(); 150c4762a1bSJed Brown return ierr; 151c4762a1bSJed Brown } 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 154c4762a1bSJed Brown PetscErrorCode InitialConditions(DM da,Vec U) 155c4762a1bSJed Brown { 156c4762a1bSJed Brown PetscErrorCode ierr; 157c4762a1bSJed Brown PetscInt i,j,xs,ys,xm,ym,Mx,My; 158c4762a1bSJed Brown Field **u; 159c4762a1bSJed Brown PetscReal hx,hy,x,y; 160c4762a1bSJed Brown 161c4762a1bSJed Brown PetscFunctionBegin; 162c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 163c4762a1bSJed Brown 164c4762a1bSJed Brown hx = 2.5/(PetscReal)Mx; 165c4762a1bSJed Brown hy = 2.5/(PetscReal)My; 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* 168c4762a1bSJed Brown Get pointers to vector data 169c4762a1bSJed Brown */ 170c4762a1bSJed Brown ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr); 171c4762a1bSJed Brown 172c4762a1bSJed Brown /* 173c4762a1bSJed Brown Get local grid boundaries 174c4762a1bSJed Brown */ 175c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* 178c4762a1bSJed Brown Compute function over the locally owned part of the grid 179c4762a1bSJed Brown */ 180c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 181c4762a1bSJed Brown y = j*hy; 182c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 183c4762a1bSJed Brown x = i*hx; 184c803a25aSBarry Smith if (PetscGTE(x,1.0) && PetscLTE(x,1.5) && PetscGTE(y,1.0) && PetscLTE(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; 185c4762a1bSJed Brown else u[j][i].v = 0.0; 186c4762a1bSJed Brown 187c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0*u[j][i].v; 188c4762a1bSJed Brown } 189c4762a1bSJed Brown } 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* 192c4762a1bSJed Brown Restore vectors 193c4762a1bSJed Brown */ 194c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr); 195c4762a1bSJed Brown PetscFunctionReturn(0); 196c4762a1bSJed Brown } 197c4762a1bSJed Brown 198c4762a1bSJed Brown /*TEST 199c4762a1bSJed Brown 200c4762a1bSJed Brown build: 201*60f0b76eSHong Zhang depends: reaction_diffusion.c 202c4762a1bSJed Brown requires: !complex !single 203c4762a1bSJed Brown 204c4762a1bSJed Brown test: 205c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -da_grid_x 20 -da_grid_y 20 206c4762a1bSJed Brown output_file: output/ex5adj_1.out 207c4762a1bSJed Brown 208c4762a1bSJed Brown test: 209c4762a1bSJed Brown suffix: 2 210c4762a1bSJed Brown nsize: 2 211c4762a1bSJed 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 212c4762a1bSJed Brown 213c4762a1bSJed Brown test: 214c4762a1bSJed Brown suffix: 3 215c4762a1bSJed Brown nsize: 2 216c4762a1bSJed Brown args: -ts_max_steps 10 -ts_dt 10 -ts_adjoint_monitor_draw_sensi 217c4762a1bSJed Brown 218c4762a1bSJed Brown test: 219c4762a1bSJed Brown suffix: 4 220c4762a1bSJed Brown nsize: 2 221c4762a1bSJed 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 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 233c4762a1bSJed Brown output_file: output/ex5adj_3.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: 290fcb023d4SJed Brown suffix: gamg1 291fcb023d4SJed Brown args: -pc_type gamg -pc_mg_levels 2 -mg_levels_pc_type jacobi -ksp_monitor_short -ts_max_steps 2 -ts_monitor -ts_adjoint_monitor -ts_trajectory_type memory 292fcb023d4SJed Brown output_file: output/ex5adj_gamg_1.out 293fcb023d4SJed Brown 29475f8e9bdSHong Zhang test: 29575f8e9bdSHong Zhang suffix: gamg2 29675f8e9bdSHong Zhang args: -pc_type gamg -pc_mg_levels 2 -mg_levels_pc_type jacobi -ksp_monitor_short -ts_max_steps 2 -ts_monitor -ts_adjoint_monitor -ts_trajectory_type memory -ksp_use_explicittranspose 29775f8e9bdSHong Zhang output_file: output/ex5adj_gamg_2.out 298c4762a1bSJed Brown TEST*/ 299